Actual source code: ex2.c

  1: static char help[] = "Runaway electron model with Landau collision operator\n\n";

  3: #include <petscdmplex.h>
  4: #include <petsclandau.h>
  5: #include <petscts.h>
  6: #include <petscds.h>
  7: #include <petscdmcomposite.h>
  8: #include <petsc/private/petscimpl.h>

 10: #if defined(PETSC_HAVE_CUDA_NVTX)
 11:   #if PETSC_PKG_CUDA_VERSION_GE(10, 0, 0)
 12:     #include <nvtx3/nvToolsExt.h>
 13:   #else
 14:     #include <nvToolsExt.h>
 15:   #endif
 16: #endif

 18: /* data for runaway electron model */
 19: typedef struct REctx_struct {
 20:   PetscErrorCode (*test)(TS, Vec, PetscInt, PetscReal, PetscBool, LandauCtx *, struct REctx_struct *);
 21:   PetscErrorCode (*impuritySrcRate)(PetscReal, PetscReal *, LandauCtx *);
 22:   PetscErrorCode (*E)(Vec, Vec, PetscInt, PetscReal, LandauCtx *, PetscReal *);
 23:   PetscReal T_cold;        /* temperature of newly ionized electrons and impurity ions */
 24:   PetscReal ion_potential; /* ionization potential of impurity */
 25:   PetscReal Ne_ion;        /* effective number of electrons shed in ioization of impurity */
 26:   PetscReal Ez_initial;
 27:   PetscReal L; /* inductance */
 28:   Vec       X_0;
 29:   PetscInt  imp_idx; /* index for impurity ionizing sink */
 30:   PetscReal pulse_start;
 31:   PetscReal pulse_width;
 32:   PetscReal pulse_rate;
 33:   PetscReal current_rate;
 34:   PetscInt  plotIdx;
 35:   PetscInt  plotStep;
 36:   PetscInt  idx; /* cache */
 37:   PetscReal j;   /* cache */
 38:   PetscReal plotDt;
 39:   PetscBool plotting;
 40:   PetscBool use_spitzer_eta;
 41:   PetscInt  print_period;
 42:   PetscInt  grid_view_idx;
 43: } REctx;

 45: static const PetscReal kev_joul = 6.241506479963235e+15; /* 1/1000e */

 47: #define RE_CUT 3.
 48: /* < v, u_re * v * q > */
 49: static void f0_j_re(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0)
 50: {
 51:   PetscReal n_e = PetscRealPart(u[0]);
 52:   if (dim == 2) {
 53:     if (x[1] > RE_CUT || x[1] < -RE_CUT) {                    /* simply a cutoff for REs. v_|| > 3 v(T_e) */
 54:       *f0 = n_e * 2. * PETSC_PI * x[0] * x[1] * constants[0]; /* n * r * v_|| * q */
 55:     } else {
 56:       *f0 = 0;
 57:     }
 58:   } else {
 59:     if (x[2] > RE_CUT || x[2] < -RE_CUT) { /* simply a cutoff for REs. v_|| > 3 v(T_e) */
 60:       *f0 = n_e * x[2] * constants[0];
 61:     } else {
 62:       *f0 = 0;
 63:     }
 64:   }
 65: }

 67: /* sum < v, u*v*q > */
 68: static void f0_jz_sum(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar q[], PetscScalar *f0)
 69: {
 70:   PetscInt ii;
 71:   f0[0] = 0;
 72:   if (dim == 2) {
 73:     for (ii = 0; ii < Nf; ii++) f0[0] += u[ii] * 2. * PETSC_PI * x[0] * x[1] * q[ii]; /* n * r * v_|| * q * v_0 */
 74:   } else {
 75:     for (ii = 0; ii < Nf; ii++) f0[0] += u[ii] * x[2] * q[ii]; /* n * v_|| * q  * v_0 */
 76:   }
 77: }

 79: /* < v, n_e > */
 80: static void f0_n(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0)
 81: {
 82:   PetscInt ii = (PetscInt)PetscRealPart(constants[0]);
 83:   if (dim == 2) f0[0] = 2. * PETSC_PI * x[0] * u[ii];
 84:   else f0[0] = u[ii];
 85: }

 87: /* < v, n_e v_|| > */
 88: static void f0_vz(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0)
 89: {
 90:   PetscInt ii = (PetscInt)PetscRealPart(constants[0]);
 91:   if (dim == 2) f0[0] = u[ii] * 2. * PETSC_PI * x[0] * x[1]; /* n r v_|| */
 92:   else f0[0] = u[ii] * x[2];                                 /* n v_|| */
 93: }

 95: /* < v, n_e (v-shift) > */
 96: static void f0_ve_shift(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0)
 97: {
 98:   PetscReal vz = numConstants > 0 ? PetscRealPart(constants[0]) : 0;
 99:   if (dim == 2) *f0 = u[0] * 2. * PETSC_PI * x[0] * PetscSqrtReal(x[0] * x[0] + (x[1] - vz) * (x[1] - vz)); /* n r v */
100:   else *f0 = u[0] * PetscSqrtReal(x[0] * x[0] + x[1] * x[1] + (x[2] - vz) * (x[2] - vz));                   /* n v */
101: }

103: /* CalculateE - Calculate the electric field  */
104: /*  T        -- Electron temperature  */
105: /*  n        -- Electron density  */
106: /*  lnLambda --   */
107: /*  eps0     --  */
108: /*  E        -- output E, input \hat E */
109: static PetscReal CalculateE(PetscReal Tev, PetscReal n, PetscReal lnLambda, PetscReal eps0, PetscReal *E)
110: {
111:   PetscReal c, e, m;

113:   PetscFunctionBegin;
114:   c = 299792458.0;
115:   e = 1.602176e-19;
116:   m = 9.10938e-31;
117:   if (1) {
118:     double Ec, Ehat = *E, betath = PetscSqrtReal(2 * Tev * e / (m * c * c)), j0 = Ehat * 7 / (PetscSqrtReal(2) * 2) * PetscPowReal(betath, 3) * n * e * c;
119:     Ec = n * lnLambda * PetscPowReal(e, 3) / (4 * PETSC_PI * PetscPowReal(eps0, 2) * m * c * c);
120:     *E = Ec;
121:     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "CalculateE j0=%g Ec = %g\n", j0, Ec));
122:   } else {
123:     PetscReal Ed, vth;
124:     vth = PetscSqrtReal(8 * Tev * e / (m * PETSC_PI));
125:     Ed  = n * lnLambda * PetscPowReal(e, 3) / (4 * PETSC_PI * PetscPowReal(eps0, 2) * m * vth * vth);
126:     *E  = Ed;
127:   }
128:   PetscFunctionReturn(PETSC_SUCCESS);
129: }

131: static PetscReal Spitzer(PetscReal m_e, PetscReal e, PetscReal Z, PetscReal epsilon0, PetscReal lnLam, PetscReal kTe_joules)
132: {
133:   PetscReal Fz = (1 + 1.198 * Z + 0.222 * Z * Z) / (1 + 2.966 * Z + 0.753 * Z * Z), eta;
134:   eta          = Fz * 4. / 3. * PetscSqrtReal(2. * PETSC_PI) * Z * PetscSqrtReal(m_e) * PetscSqr(e) * lnLam * PetscPowReal(4 * PETSC_PI * epsilon0, -2.) * PetscPowReal(kTe_joules, -1.5);
135:   return eta;
136: }

138: static PetscErrorCode testNone(TS ts, Vec X, PetscInt stepi, PetscReal time, PetscBool islast, LandauCtx *ctx, REctx *rectx)
139: {
140:   PetscFunctionBeginUser;
141:   PetscFunctionReturn(PETSC_SUCCESS);
142: }

144: static PetscErrorCode testSpitzer(TS ts, Vec X, PetscInt stepi, PetscReal time, PetscBool islast, LandauCtx *ctx, REctx *rectx)
145: {
146:   PetscInt          nDMs;
147:   PetscDS           prob;
148:   static PetscReal  old_ratio = 1e10;
149:   TSConvergedReason reason;
150:   PetscReal         J, J_re, spit_eta, Te_kev = 0, E, ratio, Z, n_e, v, v2;
151:   PetscScalar       user[2] = {0., ctx->charges[0]}, q[LANDAU_MAX_SPECIES], tt[LANDAU_MAX_SPECIES], vz;
152:   PetscReal         dt;
153:   DM                pack, plexe = ctx->plex[0], plexi = (ctx->num_grids == 1) ? NULL : ctx->plex[1];
154:   Vec              *XsubArray;

156:   PetscFunctionBeginUser;
157:   PetscCheck(ctx->num_species == 2, PETSC_COMM_SELF, PETSC_ERR_PLIB, "ctx->num_species %" PetscInt_FMT " != 2", ctx->num_species);
158:   PetscCall(VecGetDM(X, &pack));
159:   PetscCheck(pack, PETSC_COMM_SELF, PETSC_ERR_PLIB, "no DM");
160:   PetscCall(DMCompositeGetNumberDM(pack, &nDMs));
161:   PetscCheck(nDMs == ctx->num_grids * ctx->batch_sz, PETSC_COMM_SELF, PETSC_ERR_PLIB, "nDMs != ctx->num_grids*ctx->batch_sz %" PetscInt_FMT " != %" PetscInt_FMT, nDMs, ctx->num_grids * ctx->batch_sz);
162:   PetscCall(PetscMalloc(sizeof(*XsubArray) * nDMs, &XsubArray));
163:   PetscCall(DMCompositeGetAccessArray(pack, X, nDMs, NULL, XsubArray)); // read only
164:   PetscCall(TSGetTimeStep(ts, &dt));
165:   /* get current for each grid */
166:   for (PetscInt ii = 0; ii < ctx->num_species; ii++) q[ii] = ctx->charges[ii];
167:   PetscCall(DMGetDS(plexe, &prob));
168:   PetscCall(PetscDSSetConstants(prob, 2, &q[0]));
169:   PetscCall(PetscDSSetObjective(prob, 0, &f0_jz_sum));
170:   PetscCall(DMPlexComputeIntegralFEM(plexe, XsubArray[LAND_PACK_IDX(ctx->batch_view_idx, 0)], tt, NULL));
171:   J = -ctx->n_0 * ctx->v_0 * PetscRealPart(tt[0]);
172:   if (plexi) { // add first (only) ion
173:     PetscCall(DMGetDS(plexi, &prob));
174:     PetscCall(PetscDSSetConstants(prob, 1, &q[1]));
175:     PetscCall(PetscDSSetObjective(prob, 0, &f0_jz_sum));
176:     PetscCall(DMPlexComputeIntegralFEM(plexi, XsubArray[LAND_PACK_IDX(ctx->batch_view_idx, 1)], tt, NULL));
177:     J += -ctx->n_0 * ctx->v_0 * PetscRealPart(tt[0]);
178:   }
179:   /* get N_e */
180:   PetscCall(DMGetDS(plexe, &prob));
181:   PetscCall(PetscDSSetConstants(prob, 1, user));
182:   PetscCall(PetscDSSetObjective(prob, 0, &f0_n));
183:   PetscCall(DMPlexComputeIntegralFEM(plexe, XsubArray[LAND_PACK_IDX(ctx->batch_view_idx, 0)], tt, NULL));
184:   n_e = PetscRealPart(tt[0]) * ctx->n_0;
185:   /* Z */
186:   Z = -ctx->charges[1] / ctx->charges[0];
187:   /* remove drift */
188:   if (0) {
189:     user[0] = 0; // electrons
190:     PetscCall(DMGetDS(plexe, &prob));
191:     PetscCall(PetscDSSetConstants(prob, 1, user));
192:     PetscCall(PetscDSSetObjective(prob, 0, &f0_vz));
193:     PetscCall(DMPlexComputeIntegralFEM(plexe, XsubArray[LAND_PACK_IDX(ctx->batch_view_idx, 0)], tt, NULL));
194:     vz = ctx->n_0 * PetscRealPart(tt[0]) / n_e; /* non-dimensional */
195:   } else vz = 0;
196:   /* thermal velocity */
197:   PetscCall(DMGetDS(plexe, &prob));
198:   PetscCall(PetscDSSetConstants(prob, 1, &vz));
199:   PetscCall(PetscDSSetObjective(prob, 0, &f0_ve_shift));
200:   PetscCall(DMPlexComputeIntegralFEM(plexe, XsubArray[LAND_PACK_IDX(ctx->batch_view_idx, 0)], tt, NULL));
201:   v        = ctx->n_0 * ctx->v_0 * PetscRealPart(tt[0]) / n_e;                                                   /* remove number density to get velocity */
202:   v2       = PetscSqr(v);                                                                                        /* use real space: m^2 / s^2 */
203:   Te_kev   = (v2 * ctx->masses[0] * PETSC_PI / 8) * kev_joul;                                                    /* temperature in kev */
204:   spit_eta = Spitzer(ctx->masses[0], -ctx->charges[0], Z, ctx->epsilon0, ctx->lambdas[0][1], Te_kev / kev_joul); /* kev --> J (kT) */
205:   if (0) {
206:     PetscCall(DMGetDS(plexe, &prob));
207:     PetscCall(PetscDSSetConstants(prob, 1, q));
208:     PetscCall(PetscDSSetObjective(prob, 0, &f0_j_re));
209:     PetscCall(DMPlexComputeIntegralFEM(plexe, XsubArray[LAND_PACK_IDX(ctx->batch_view_idx, 0)], tt, NULL));
210:   } else tt[0] = 0;
211:   J_re = -ctx->n_0 * ctx->v_0 * PetscRealPart(tt[0]);
212:   PetscCall(DMCompositeRestoreAccessArray(pack, X, nDMs, NULL, XsubArray)); // read only
213:   PetscCall(PetscFree(XsubArray));

215:   if (rectx->use_spitzer_eta) {
216:     E = ctx->Ez = spit_eta * (rectx->j - J_re);
217:   } else {
218:     E        = ctx->Ez; /* keep real E */
219:     rectx->j = J;       /* cache */
220:   }

222:   ratio = E / J / spit_eta;
223:   if (stepi > 10 && !rectx->use_spitzer_eta && (old_ratio - ratio < 1.e-6)) {
224:     rectx->pulse_start     = time + 0.98 * dt;
225:     rectx->use_spitzer_eta = PETSC_TRUE;
226:   }
227:   PetscCall(TSGetConvergedReason(ts, &reason));
228:   PetscCall(TSGetConvergedReason(ts, &reason));
229:   if (rectx->plotting || stepi == 0 || reason || rectx->pulse_start == time + 0.98 * dt) {
230:     PetscCall(PetscPrintf(ctx->comm, "testSpitzer: %4" PetscInt_FMT ") time=%11.4e n_e= %10.3e E= %10.3e J= %10.3e J_re= %10.3e %.3g%% Te_kev= %10.3e Z_eff=%g E/J to eta ratio= %g (diff=%g) %s %s spit_eta=%g\n", stepi, (double)time,
231:                           (double)(n_e / ctx->n_0), (double)ctx->Ez, (double)J, (double)J_re, (double)(100 * J_re / J), (double)Te_kev, (double)Z, (double)ratio, (double)(old_ratio - ratio), rectx->use_spitzer_eta ? "using Spitzer eta*J E" : "constant E", rectx->pulse_start != time + 0.98 * dt ? "normal" : "transition", (double)spit_eta));
232:     PetscCheck(rectx->pulse_start != (time + 0.98 * dt), PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Spitzer complete ratio=%g", (double)ratio);
233:   }
234:   old_ratio = ratio;
235:   PetscFunctionReturn(PETSC_SUCCESS);
236: }

238: static const double ppp = 2;
239: static void f0_0_diff_lp(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0)
240: {
241:   LandauCtx      *ctx   = (LandauCtx *)constants;
242:   REctx          *rectx = (REctx *)ctx->data;
243:   PetscInt        ii    = rectx->idx, i;
244:   const PetscReal kT_m  = ctx->k * ctx->thermal_temps[ii] / ctx->masses[ii]; /* kT/m */
245:   const PetscReal n     = ctx->n[ii];
246:   PetscReal       diff, f_maxwell, v2 = 0, theta = 2 * kT_m / (ctx->v_0 * ctx->v_0); /* theta = 2kT/mc^2 */
247:   for (i = 0; i < dim; ++i) v2 += x[i] * x[i];
248:   f_maxwell = n * PetscPowReal(PETSC_PI * theta, -1.5) * (PetscExpReal(-v2 / theta));
249:   diff      = 2. * PETSC_PI * x[0] * (PetscRealPart(u[ii]) - f_maxwell);
250:   f0[0]     = PetscPowReal(diff, ppp);
251: }

253: static void f0_0_maxwellian_lp(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0)
254: {
255:   LandauCtx      *ctx   = (LandauCtx *)constants;
256:   REctx          *rectx = (REctx *)ctx->data;
257:   PetscInt        ii    = rectx->idx, i;
258:   const PetscReal kT_m  = ctx->k * ctx->thermal_temps[ii] / ctx->masses[ii]; /* kT/m */
259:   const PetscReal n     = ctx->n[ii];
260:   PetscReal       f_maxwell, v2 = 0, theta = 2 * kT_m / (ctx->v_0 * ctx->v_0); /* theta = 2kT/mc^2 */
261:   for (i = 0; i < dim; ++i) v2 += x[i] * x[i];
262:   f_maxwell = 2. * PETSC_PI * x[0] * n * PetscPowReal(PETSC_PI * theta, -1.5) * (PetscExpReal(-v2 / theta));
263:   f0[0]     = PetscPowReal(f_maxwell, ppp);
264: }

266: static PetscErrorCode testStable(TS ts, Vec X, PetscInt stepi, PetscReal time, PetscBool islast, LandauCtx *ctx, REctx *rectx)
267: {
268:   PetscDS     prob;
269:   Vec         X2;
270:   PetscReal   ediff, idiff = 0, lpm0, lpm1 = 1;
271:   PetscScalar tt[LANDAU_MAX_SPECIES];
272:   DM          dm, plex = ctx->plex[0];

274:   PetscFunctionBeginUser;
275:   PetscCall(VecGetDM(X, &dm));
276:   PetscCall(DMGetDS(plex, &prob));
277:   PetscCall(VecDuplicate(X, &X2));
278:   PetscCall(VecCopy(X, X2));
279:   if (!rectx->X_0) {
280:     PetscCall(VecDuplicate(X, &rectx->X_0));
281:     PetscCall(VecCopy(X, rectx->X_0));
282:   }
283:   PetscCall(VecAXPY(X, -1.0, rectx->X_0));
284:   PetscCall(PetscDSSetConstants(prob, sizeof(LandauCtx) / sizeof(PetscScalar), (PetscScalar *)ctx));
285:   rectx->idx = 0;
286:   PetscCall(PetscDSSetObjective(prob, 0, &f0_0_diff_lp));
287:   PetscCall(DMPlexComputeIntegralFEM(plex, X2, tt, NULL));
288:   ediff = PetscPowReal(PetscRealPart(tt[0]), 1. / ppp);
289:   PetscCall(PetscDSSetObjective(prob, 0, &f0_0_maxwellian_lp));
290:   PetscCall(DMPlexComputeIntegralFEM(plex, X2, tt, NULL));
291:   lpm0 = PetscPowReal(PetscRealPart(tt[0]), 1. / ppp);
292:   if (ctx->num_species > 1) {
293:     rectx->idx = 1;
294:     PetscCall(PetscDSSetObjective(prob, 0, &f0_0_diff_lp));
295:     PetscCall(DMPlexComputeIntegralFEM(plex, X2, tt, NULL));
296:     idiff = PetscPowReal(PetscRealPart(tt[0]), 1. / ppp);
297:     PetscCall(PetscDSSetObjective(prob, 0, &f0_0_maxwellian_lp));
298:     PetscCall(DMPlexComputeIntegralFEM(plex, X2, tt, NULL));
299:     lpm1 = PetscPowReal(PetscRealPart(tt[0]), 1. / ppp);
300:   }
301:   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s %" PetscInt_FMT ") time=%10.3e n-%d norm electrons/max=%20.13e ions/max=%20.13e\n", "----", stepi, (double)time, (int)ppp, (double)(ediff / lpm0), (double)(idiff / lpm1)));
302:   /* view */
303:   PetscCall(VecCopy(X2, X));
304:   PetscCall(VecDestroy(&X2));
305:   if (islast) {
306:     PetscCall(VecDestroy(&rectx->X_0));
307:     rectx->X_0 = NULL;
308:   }
309:   PetscFunctionReturn(PETSC_SUCCESS);
310: }

312: static PetscErrorCode EInduction(Vec X, Vec X_t, PetscInt step, PetscReal time, LandauCtx *ctx, PetscReal *a_E)
313: {
314:   REctx      *rectx = (REctx *)ctx->data;
315:   DM          dm, plex;
316:   PetscScalar tt[LANDAU_MAX_SPECIES], qv0[LANDAU_MAX_SPECIES];
317:   PetscReal   dJ_dt;
318:   PetscDS     prob;

320:   PetscFunctionBeginUser;
321:   for (PetscInt ii = 0; ii < ctx->num_species; ii++) qv0[ii] = ctx->charges[ii] * ctx->v_0;
322:   PetscCall(VecGetDM(X, &dm));
323:   PetscCall(DMGetDS(dm, &prob));
324:   PetscCall(DMConvert(dm, DMPLEX, &plex));
325:   /* get d current / dt */
326:   PetscCall(PetscDSSetConstants(prob, ctx->num_species, qv0));
327:   PetscCall(PetscDSSetObjective(prob, 0, &f0_jz_sum));
328:   PetscCheck(X_t, PETSC_COMM_SELF, PETSC_ERR_PLIB, "X_t");
329:   PetscCall(DMPlexComputeIntegralFEM(plex, X_t, tt, NULL));
330:   dJ_dt = -ctx->n_0 * PetscRealPart(tt[0]) / ctx->t_0;
331:   /* E induction */
332:   *a_E = -rectx->L * dJ_dt + rectx->Ez_initial;
333:   PetscCall(DMDestroy(&plex));
334:   PetscFunctionReturn(PETSC_SUCCESS);
335: }

337: static PetscErrorCode EConstant(Vec X, Vec X_t, PetscInt step, PetscReal time, LandauCtx *ctx, PetscReal *a_E)
338: {
339:   PetscFunctionBeginUser;
340:   *a_E = ctx->Ez;
341:   PetscFunctionReturn(PETSC_SUCCESS);
342: }

344: static PetscErrorCode ENone(Vec X, Vec X_t, PetscInt step, PetscReal time, LandauCtx *ctx, PetscReal *a_E)
345: {
346:   PetscFunctionBeginUser;
347:   *a_E = 0;
348:   PetscFunctionReturn(PETSC_SUCCESS);
349: }

351: /* ------------------------------------------------------------------- */
352: /*
353:    FormSource - Evaluates source terms F(t).

355:    Input Parameters:
356: .  ts - the TS context
357: .  time -
358: .  X_dummmy - input vector
359: .  dummy - optional user-defined context, as set by SNESSetFunction()

361:    Output Parameter:
362: .  F - function vector
363:  */
364: static PetscErrorCode FormSource(TS ts, PetscReal ftime, Vec X_dummmy, Vec F, void *dummy)
365: {
366:   PetscReal  new_imp_rate;
367:   LandauCtx *ctx;
368:   DM         pack;
369:   REctx     *rectx;

371:   PetscFunctionBeginUser;
372:   PetscCall(TSGetDM(ts, &pack));
373:   PetscCall(DMGetApplicationContext(pack, &ctx));
374:   rectx = (REctx *)ctx->data;
375:   /* check for impurities */
376:   PetscCall(rectx->impuritySrcRate(ftime, &new_imp_rate, ctx));
377:   if (new_imp_rate != 0) {
378:     if (new_imp_rate != rectx->current_rate) {
379:       PetscReal dne_dt, dni_dt, tilda_ns[LANDAU_MAX_SPECIES], temps[LANDAU_MAX_SPECIES];
380:       Vec       globFarray[LANDAU_MAX_GRIDS * LANDAU_MAX_BATCH_SZ];
381:       rectx->current_rate = new_imp_rate;
382:       for (PetscInt ii = 1; ii < LANDAU_MAX_SPECIES; ii++) tilda_ns[ii] = 0;
383:       for (PetscInt ii = 1; ii < LANDAU_MAX_SPECIES; ii++) temps[ii] = 1;
384:       dni_dt                   = new_imp_rate /* *ctx->t_0 */; /* fully ionized immediately, no normalize, stay in non-dim */
385:       dne_dt                   = new_imp_rate * rectx->Ne_ion /* *ctx->t_0 */;
386:       tilda_ns[0]              = dne_dt;
387:       tilda_ns[rectx->imp_idx] = dni_dt;
388:       temps[0]                 = rectx->T_cold;
389:       temps[rectx->imp_idx]    = rectx->T_cold;
390:       PetscCall(PetscInfo(ctx->plex[0], "\tHave new_imp_rate= %10.3e time= %10.3e de/dt= %10.3e di/dt= %10.3e ***\n", (double)new_imp_rate, (double)ftime, (double)dne_dt, (double)dni_dt));
391:       PetscCall(DMCompositeGetAccessArray(pack, F, ctx->num_grids * ctx->batch_sz, NULL, globFarray));
392:       for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
393:         /* add it */
394:         PetscCall(DMPlexLandauAddMaxwellians(ctx->plex[grid], globFarray[LAND_PACK_IDX(0, grid)], ftime, temps, tilda_ns, grid, 0, 1, ctx));
395:       }
396:       // Does DMCompositeRestoreAccessArray copy the data back? (no)
397:       PetscCall(DMCompositeRestoreAccessArray(pack, F, ctx->num_grids * ctx->batch_sz, NULL, globFarray));
398:     }
399:   } else {
400:     PetscCall(VecZeroEntries(F));
401:     rectx->current_rate = 0;
402:   }
403:   PetscFunctionReturn(PETSC_SUCCESS);
404: }

406: static PetscErrorCode Monitor(TS ts, PetscInt stepi, PetscReal time, Vec X, void *actx)
407: {
408:   LandauCtx        *ctx   = (LandauCtx *)actx; /* user-defined application context */
409:   REctx            *rectx = (REctx *)ctx->data;
410:   DM                pack  = NULL;
411:   Vec               globXArray[LANDAU_MAX_GRIDS * LANDAU_MAX_BATCH_SZ];
412:   TSConvergedReason reason;

414:   PetscFunctionBeginUser;
415:   PetscCall(TSGetConvergedReason(ts, &reason));
416:   if (rectx->grid_view_idx != -1 || (reason && ctx->verbose > 3)) {
417:     PetscCall(VecGetDM(X, &pack));
418:     PetscCall(DMCompositeGetAccessArray(pack, X, ctx->num_grids * ctx->batch_sz, NULL, globXArray));
419:   }
420:   if (stepi > rectx->plotStep && rectx->plotting) {
421:     rectx->plotting = PETSC_FALSE; /* was doing diagnostics, now done */
422:     rectx->plotIdx++;
423:   }
424:   /* view */
425:   if (time / rectx->plotDt >= (PetscReal)rectx->plotIdx || reason) {
426:     if ((reason || stepi == 0 || rectx->plotIdx % rectx->print_period == 0) && ctx->verbose > 1) {
427:       /* print norms */
428:       PetscCall(DMPlexLandauPrintNorms(X, stepi));
429:     }
430:     if (!rectx->plotting) { /* first step of possible backtracks */
431:       rectx->plotting = PETSC_TRUE;
432:       /* diagnostics + change E field with Sptizer (not just a monitor) */
433:       PetscCall(rectx->test(ts, X, stepi, time, reason ? PETSC_TRUE : PETSC_FALSE, ctx, rectx));
434:     } else {
435:       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\t\t ERROR SKIP test spit ------\n"));
436:       rectx->plotting = PETSC_TRUE;
437:     }
438:     if (rectx->grid_view_idx != -1) {
439:       PetscCall(PetscObjectSetName((PetscObject)globXArray[LAND_PACK_IDX(ctx->batch_view_idx, rectx->grid_view_idx)], rectx->grid_view_idx == 0 ? "ue" : "ui"));
440:       /* view, overwrite step when back tracked */
441:       PetscCall(DMSetOutputSequenceNumber(ctx->plex[rectx->grid_view_idx], rectx->plotIdx, time * ctx->t_0));
442:       PetscCall(VecViewFromOptions(globXArray[LAND_PACK_IDX(ctx->batch_view_idx, rectx->grid_view_idx)], NULL, "-ex2_vec_view"));
443:     }
444:     rectx->plotStep = stepi;
445:   } else {
446:     if (rectx->plotting) PetscCall(PetscPrintf(PETSC_COMM_WORLD, " ERROR rectx->plotting=%s step %" PetscInt_FMT "\n", PetscBools[rectx->plotting], stepi));
447:     /* diagnostics + change E field with Sptizer (not just a monitor) - can we lag this? */
448:     PetscCall(rectx->test(ts, X, stepi, time, reason ? PETSC_TRUE : PETSC_FALSE, ctx, rectx));
449:   }
450:   /* parallel check that only works of all batches are identical */
451:   if (reason && ctx->verbose > 3 && ctx->batch_sz > 1) {
452:     PetscReal   val, rval;
453:     PetscMPIInt rank;
454:     PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
455:     for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
456:       PetscInt nerrors = 0;
457:       for (PetscInt i = 0; i < ctx->batch_sz; i++) {
458:         PetscCall(VecNorm(globXArray[LAND_PACK_IDX(i, grid)], NORM_2, &val));
459:         if (i == 0) rval = val;
460:         else if ((val = PetscAbs(val - rval) / rval) > 1000 * PETSC_MACHINE_EPSILON) {
461:           PetscCall(PetscPrintf(PETSC_COMM_SELF, " [%d] Warning %" PetscInt_FMT ".%" PetscInt_FMT ") diff = %2.15e\n", rank, grid, i, (double)val));
462:           nerrors++;
463:         }
464:       }
465:       if (nerrors) {
466:         PetscCall(PetscPrintf(PETSC_COMM_SELF, " ***** [%d] ERROR max %" PetscInt_FMT " errors\n", rank, nerrors));
467:       } else {
468:         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "[%d] %" PetscInt_FMT ") batch consistency check OK\n", rank, grid));
469:       }
470:     }
471:   }
472:   rectx->idx = 0;
473:   if (rectx->grid_view_idx != -1 || (reason && ctx->verbose > 3)) PetscCall(DMCompositeRestoreAccessArray(pack, X, ctx->num_grids * ctx->batch_sz, NULL, globXArray));
474:   PetscFunctionReturn(PETSC_SUCCESS);
475: }

477: static PetscErrorCode PreStep(TS ts)
478: {
479:   LandauCtx *ctx;
480:   REctx     *rectx;
481:   DM         dm;
482:   PetscInt   stepi;
483:   PetscReal  time;
484:   Vec        X;

486:   PetscFunctionBeginUser;
487:   /* not used */
488:   PetscCall(TSGetDM(ts, &dm));
489:   PetscCall(TSGetTime(ts, &time));
490:   PetscCall(TSGetSolution(ts, &X));
491:   PetscCall(DMGetApplicationContext(dm, &ctx));
492:   rectx = (REctx *)ctx->data;
493:   PetscCall(TSGetStepNumber(ts, &stepi));
494:   /* update E */
495:   PetscCall(rectx->E(X, NULL, stepi, time, ctx, &ctx->Ez));
496:   PetscFunctionReturn(PETSC_SUCCESS);
497: }

499: /* model for source of non-ionized impurities, profile provided by model, in du/dt form in normalized units (tricky because n_0 is normalized with electrons) */
500: static PetscErrorCode stepSrc(PetscReal time, PetscReal *rho, LandauCtx *ctx)
501: {
502:   REctx *rectx = (REctx *)ctx->data;

504:   PetscFunctionBeginUser;
505:   if (time >= rectx->pulse_start) *rho = rectx->pulse_rate;
506:   else *rho = 0.;
507:   PetscFunctionReturn(PETSC_SUCCESS);
508: }
509: static PetscErrorCode zeroSrc(PetscReal time, PetscReal *rho, LandauCtx *ctx)
510: {
511:   PetscFunctionBeginUser;
512:   *rho = 0.;
513:   PetscFunctionReturn(PETSC_SUCCESS);
514: }

516: static PetscErrorCode pulseSrc(PetscReal time, PetscReal *rho, LandauCtx *ctx)
517: {
518:   REctx *rectx = (REctx *)ctx->data;

520:   PetscFunctionBeginUser;
521:   PetscCheck(rectx->pulse_start != PETSC_MAX_REAL, PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "'-ex2_pulse_start_time X' must be used with '-ex2_impurity_source_type pulse'");
522:   if (time < rectx->pulse_start || time > rectx->pulse_start + 3 * rectx->pulse_width) *rho = 0;
523:   else {
524:     double x = PetscSinReal((time - rectx->pulse_start) / (3 * rectx->pulse_width) * 2 * PETSC_PI - PETSC_PI / 2) + 1; /* 0:2, integrates to 1.0 */
525:     *rho     = rectx->pulse_rate * x / (3 * rectx->pulse_width);
526:     if (!rectx->use_spitzer_eta) rectx->use_spitzer_eta = PETSC_TRUE; /* use it next time */
527:   }
528:   PetscFunctionReturn(PETSC_SUCCESS);
529: }

531: static PetscErrorCode ProcessREOptions(REctx *rectx, const LandauCtx *ctx, DM dm, const char prefix[])
532: {
533:   PetscFunctionList plist = NULL, testlist = NULL, elist = NULL;
534:   char              pname[256], testname[256], ename[256];
535:   DM                dm_dummy;
536:   PetscBool         Connor_E = PETSC_FALSE;

538:   PetscFunctionBeginUser;
539:   PetscCall(DMCreate(PETSC_COMM_WORLD, &dm_dummy));
540:   rectx->Ne_ion          = 1;    /* number of electrons given up by impurity ion */
541:   rectx->T_cold          = .005; /* kev */
542:   rectx->ion_potential   = 15;   /* ev */
543:   rectx->L               = 2;
544:   rectx->X_0             = NULL;
545:   rectx->imp_idx         = ctx->num_species - 1; /* default ionized impurity as last one */
546:   rectx->pulse_start     = PETSC_MAX_REAL;
547:   rectx->pulse_width     = 1;
548:   rectx->plotStep        = PETSC_INT_MAX;
549:   rectx->pulse_rate      = 1.e-1;
550:   rectx->current_rate    = 0;
551:   rectx->plotIdx         = 0;
552:   rectx->j               = 0;
553:   rectx->plotDt          = 1.0;
554:   rectx->plotting        = PETSC_FALSE;
555:   rectx->use_spitzer_eta = PETSC_FALSE;
556:   rectx->idx             = 0;
557:   rectx->print_period    = 10;
558:   rectx->grid_view_idx   = -1; // do not get if not needed
559:   /* Register the available impurity sources */
560:   PetscCall(PetscFunctionListAdd(&plist, "step", &stepSrc));
561:   PetscCall(PetscFunctionListAdd(&plist, "none", &zeroSrc));
562:   PetscCall(PetscFunctionListAdd(&plist, "pulse", &pulseSrc));
563:   PetscCall(PetscStrncpy(pname, "none", sizeof(pname)));
564:   PetscCall(PetscFunctionListAdd(&testlist, "none", &testNone));
565:   PetscCall(PetscFunctionListAdd(&testlist, "spitzer", &testSpitzer));
566:   PetscCall(PetscFunctionListAdd(&testlist, "stable", &testStable));
567:   PetscCall(PetscStrncpy(testname, "none", sizeof(testname)));
568:   PetscCall(PetscFunctionListAdd(&elist, "none", &ENone));
569:   PetscCall(PetscFunctionListAdd(&elist, "induction", &EInduction));
570:   PetscCall(PetscFunctionListAdd(&elist, "constant", &EConstant));
571:   PetscCall(PetscStrncpy(ename, "constant", sizeof(ename)));

573:   PetscOptionsBegin(PETSC_COMM_SELF, prefix, "Options for Runaway/seed electron model", "none");
574:   PetscCall(PetscOptionsReal("-ex2_plot_dt", "Plotting interval", "ex2.c", rectx->plotDt, &rectx->plotDt, NULL));
575:   if (rectx->plotDt < 0) rectx->plotDt = 1e30;
576:   if (rectx->plotDt == 0) rectx->plotDt = 1e-30;
577:   PetscCall(PetscOptionsInt("-ex2_print_period", "Plotting interval", "ex2.c", rectx->print_period, &rectx->print_period, NULL));
578:   PetscCall(PetscOptionsInt("-ex2_grid_view_idx", "grid_view_idx", "ex2.c", rectx->grid_view_idx, &rectx->grid_view_idx, NULL));
579:   PetscCheck(rectx->grid_view_idx < ctx->num_grids || rectx->grid_view_idx == -1, PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "rectx->grid_view_idx (%" PetscInt_FMT ") >= ctx->num_grids (%" PetscInt_FMT ")", rectx->imp_idx, ctx->num_grids);
580:   PetscCall(PetscOptionsFList("-ex2_impurity_source_type", "Name of impurity source to run", "", plist, pname, pname, sizeof(pname), NULL));
581:   PetscCall(PetscOptionsFList("-ex2_test_type", "Name of test to run", "", testlist, testname, testname, sizeof(testname), NULL));
582:   PetscCall(PetscOptionsInt("-ex2_impurity_index", "index of sink for impurities", "none", rectx->imp_idx, &rectx->imp_idx, NULL));
583:   PetscCheck((rectx->imp_idx < ctx->num_species && rectx->imp_idx >= 1) || ctx->num_species <= 1, PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "index of sink for impurities ions is out of range (%" PetscInt_FMT "), must be > 0 && < NS", rectx->imp_idx);
584:   PetscCall(PetscOptionsFList("-ex2_e_field_type", "Electric field type", "", elist, ename, ename, sizeof(ename), NULL));
585:   rectx->Ne_ion = -ctx->charges[rectx->imp_idx] / ctx->charges[0];
586:   PetscCall(PetscOptionsReal("-ex2_t_cold", "Temperature of cold electron and ions after ionization in keV", "none", rectx->T_cold, &rectx->T_cold, NULL));
587:   PetscCall(PetscOptionsReal("-ex2_pulse_start_time", "Time at which pulse happens for 'pulse' source", "none", rectx->pulse_start, &rectx->pulse_start, NULL));
588:   PetscCall(PetscOptionsReal("-ex2_pulse_width_time", "Width of pulse 'pulse' source", "none", rectx->pulse_width, &rectx->pulse_width, NULL));
589:   PetscCall(PetscOptionsReal("-ex2_pulse_rate", "Number density of pulse for 'pulse' source", "none", rectx->pulse_rate, &rectx->pulse_rate, NULL));
590:   rectx->T_cold *= 1.16e7; /* convert to Kelvin */
591:   PetscCall(PetscOptionsReal("-ex2_ion_potential", "Potential to ionize impurity (should be array) in ev", "none", rectx->ion_potential, &rectx->ion_potential, NULL));
592:   PetscCall(PetscOptionsReal("-ex2_inductance", "Inductance E field", "none", rectx->L, &rectx->L, NULL));
593:   PetscCall(PetscOptionsBool("-ex2_connor_e_field_units", "Scale Ex but Connor-Hastie E_c", "none", Connor_E, &Connor_E, NULL));
594:   PetscCall(PetscInfo(dm_dummy, "Num electrons from ions=%g, T_cold=%10.3e, ion potential=%10.3e, E_z=%10.3e v_0=%10.3e\n", (double)rectx->Ne_ion, (double)rectx->T_cold, (double)rectx->ion_potential, (double)ctx->Ez, (double)ctx->v_0));
595:   PetscOptionsEnd();
596:   /* get impurity source rate function */
597:   PetscCall(PetscFunctionListFind(plist, pname, &rectx->impuritySrcRate));
598:   PetscCheck(rectx->impuritySrcRate, PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "No impurity source function found '%s'", pname);
599:   PetscCall(PetscFunctionListFind(testlist, testname, &rectx->test));
600:   PetscCheck(rectx->test, PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "No test found '%s'", testname);
601:   PetscCall(PetscFunctionListFind(elist, ename, &rectx->E));
602:   PetscCheck(rectx->E, PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "No E field function found '%s'", ename);
603:   PetscCall(PetscFunctionListDestroy(&plist));
604:   PetscCall(PetscFunctionListDestroy(&testlist));
605:   PetscCall(PetscFunctionListDestroy(&elist));

607:   /* convert E from Connor-Hastie E_c units to real if doing Spitzer E */
608:   if (Connor_E) {
609:     PetscReal E = ctx->Ez, Tev = ctx->thermal_temps[0] * 8.621738e-5, n = ctx->n_0 * ctx->n[0];
610:     CalculateE(Tev, n, ctx->lambdas[0][1], ctx->epsilon0, &E);
611:     ((LandauCtx *)ctx)->Ez *= E;
612:   }
613:   PetscCall(DMDestroy(&dm_dummy));
614:   PetscFunctionReturn(PETSC_SUCCESS);
615: }

617: int main(int argc, char **argv)
618: {
619:   DM            pack;
620:   Vec           X;
621:   PetscInt      dim = 2, nDMs;
622:   TS            ts;
623:   Mat           J;
624:   PetscDS       prob;
625:   LandauCtx    *ctx;
626:   REctx        *rectx;
627:   PetscMPIInt   rank;
628:   PetscLogStage stage;

630:   PetscFunctionBeginUser;
631:   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
632:   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
633:   if (rank) { /* turn off output stuff for duplicate runs */
634:     PetscCall(PetscOptionsClearValue(NULL, "-ex2_dm_view"));
635:     PetscCall(PetscOptionsClearValue(NULL, "-ex2_vec_view"));
636:     PetscCall(PetscOptionsClearValue(NULL, "-ex2_vec_view_init"));
637:     PetscCall(PetscOptionsClearValue(NULL, "-ex2_dm_view_init"));
638:     PetscCall(PetscOptionsClearValue(NULL, "-info")); /* this does not work */
639:   }
640:   PetscCall(PetscOptionsGetInt(NULL, NULL, "-dim", &dim, NULL));
641:   /* Create a mesh */
642:   PetscCall(DMPlexLandauCreateVelocitySpace(PETSC_COMM_WORLD, dim, "", &X, &J, &pack));
643:   PetscCall(DMCompositeGetNumberDM(pack, &nDMs));
644:   PetscCall(PetscObjectSetName((PetscObject)J, "Jacobian"));
645:   PetscCall(PetscObjectSetName((PetscObject)X, "f"));
646:   PetscCall(DMGetApplicationContext(pack, &ctx));
647:   PetscCall(DMSetUp(pack));
648:   /* context */
649:   PetscCall(PetscNew(&rectx));
650:   ctx->data = rectx;
651:   PetscCall(ProcessREOptions(rectx, ctx, pack, ""));
652:   PetscCall(DMGetDS(pack, &prob));
653:   if (rectx->grid_view_idx != -1) {
654:     Vec *XsubArray = NULL;
655:     PetscCall(PetscMalloc(sizeof(*XsubArray) * nDMs, &XsubArray));
656:     PetscCall(DMCompositeGetAccessArray(pack, X, nDMs, NULL, XsubArray)); // read only
657:     PetscCall(PetscObjectSetName((PetscObject)XsubArray[LAND_PACK_IDX(ctx->batch_view_idx, rectx->grid_view_idx)], rectx->grid_view_idx == 0 ? "ue" : "ui"));
658:     PetscCall(DMSetOutputSequenceNumber(ctx->plex[rectx->grid_view_idx], 0, 0.0));
659:     PetscCall(DMViewFromOptions(ctx->plex[rectx->grid_view_idx], NULL, "-ex2_dm_view"));
660:     PetscCall(DMViewFromOptions(ctx->plex[rectx->grid_view_idx], NULL, "-ex2_dm_view_init"));
661:     PetscCall(VecViewFromOptions(XsubArray[LAND_PACK_IDX(ctx->batch_view_idx, rectx->grid_view_idx)], NULL, "-ex2_vec_view"));      // initial condition (monitor plots after step)
662:     PetscCall(VecViewFromOptions(XsubArray[LAND_PACK_IDX(ctx->batch_view_idx, rectx->grid_view_idx)], NULL, "-ex2_vec_view_init")); // initial condition (monitor plots after step)
663:     PetscCall(DMCompositeRestoreAccessArray(pack, X, nDMs, NULL, XsubArray));                                                       // read only
664:     PetscCall(PetscFree(XsubArray));
665:   }
666:   /* Create timestepping solver context */
667:   PetscCall(TSCreate(PETSC_COMM_SELF, &ts));
668:   PetscCall(TSSetDM(ts, pack));
669:   PetscCall(TSSetIFunction(ts, NULL, DMPlexLandauIFunction, NULL));
670:   PetscCall(TSSetIJacobian(ts, J, J, DMPlexLandauIJacobian, NULL));
671:   PetscCall(TSSetRHSFunction(ts, NULL, FormSource, NULL));
672:   PetscCall(TSSetFromOptions(ts));
673:   PetscCall(TSSetSolution(ts, X));
674:   PetscCall(TSSetApplicationContext(ts, ctx));
675:   PetscCall(TSMonitorSet(ts, Monitor, ctx, NULL));
676:   PetscCall(TSSetPreStep(ts, PreStep));
677:   rectx->Ez_initial = ctx->Ez; /* cache for induction calculation - applied E field */
678:   if (1) {                     /* warm up an test just DMPlexLandauIJacobian */
679:     Vec       vec;
680:     PetscInt  nsteps;
681:     PetscReal dt;
682:     PetscCall(PetscLogStageRegister("Warmup", &stage));
683:     PetscCall(PetscLogStagePush(stage));
684:     PetscCall(VecDuplicate(X, &vec));
685:     PetscCall(VecCopy(X, vec));
686:     PetscCall(TSGetMaxSteps(ts, &nsteps));
687:     PetscCall(TSGetTimeStep(ts, &dt));
688:     PetscCall(TSSetMaxSteps(ts, 1));
689:     PetscCall(TSSolve(ts, X));
690:     PetscCall(TSSetMaxSteps(ts, nsteps));
691:     PetscCall(TSSetStepNumber(ts, 0));
692:     PetscCall(TSSetTime(ts, 0));
693:     PetscCall(TSSetTimeStep(ts, dt));
694:     rectx->plotIdx  = 0;
695:     rectx->plotting = PETSC_FALSE;
696:     PetscCall(PetscLogStagePop());
697:     PetscCall(VecCopy(vec, X));
698:     PetscCall(VecDestroy(&vec));
699:     PetscCall(PetscObjectStateIncrease((PetscObject)ctx->J));
700:   }
701:   /* go */
702:   PetscCall(PetscLogStageRegister("Solve", &stage));
703:   ctx->stage = 0; // lets not use this stage
704:   PetscCall(PetscLogStagePush(stage));
705: #if defined(PETSC_HAVE_CUDA_NVTX)
706:   nvtxRangePushA("ex2-TSSolve-warm");
707: #endif
708:   PetscCall(TSSolve(ts, X));
709: #if defined(PETSC_HAVE_CUDA_NVTX)
710:   nvtxRangePop();
711: #endif
712:   PetscCall(PetscLogStagePop());
713:   /* clean up */
714:   PetscCall(DMPlexLandauDestroyVelocitySpace(&pack));
715:   PetscCall(TSDestroy(&ts));
716:   PetscCall(VecDestroy(&X));
717:   PetscCall(PetscFree(rectx));
718:   PetscCall(PetscFinalize());
719:   return 0;
720: }

722: /*TEST

724:   testset:
725:     requires: p4est !complex double defined(PETSC_USE_DMLANDAU_2D)
726:     output_file: output/ex2_0.out
727:     args: -dm_landau_num_species_grid 1,1 -dm_landau_Ez 0 -petscspace_degree 3 -petscspace_poly_tensor 1 -dm_landau_type p4est -dm_landau_ion_masses 2 \
728:     -dm_landau_ion_charges 1 -dm_landau_thermal_temps 5,5 -dm_landau_n 2,2 -dm_landau_n_0 5e19 -ts_monitor -snes_rtol 1.e-9 -snes_stol 1.e-14 -snes_monitor \
729:     -snes_converged_reason -snes_max_it 10 -ts_type arkimex -ts_arkimex_type 1bee -ts_max_snes_failures unlimited -ts_rtol 1e-3 -ts_time_step 1.e-2 -ts_max_time 1 -ts_adapt_clip .5,1.25 -ts_max_steps 2 -ts_adapt_scale_solve_failed 0.75 -ts_adapt_time_step_increase_delay 5 -dm_landau_amr_levels_max 2,2 -dm_landau_amr_re_levels 2 -dm_landau_re_radius 0 -ex2_impurity_source_type pulse -ex2_pulse_start_time 1e-1 -ex2_pulse_width_time 10 -ex2_pulse_rate 1e-2 -ex2_t_cold .05 -ex2_plot_dt 1e-1 -dm_refine 0 -dm_landau_gpu_assembly true -dm_landau_batch_size 2 -dm_landau_verbose 2 -dm_landau_domain_radius 5.,5.
730:     test:
731:       suffix: cpu
732:       args: -dm_landau_device_type cpu -ksp_type bicg -pc_type jacobi
733:     test:
734:       suffix: kokkos
735:       requires: kokkos_kernels !defined(PETSC_HAVE_CUDA_CLANG)
736:       args: -dm_landau_device_type kokkos -dm_mat_type aijkokkos -dm_vec_type kokkos -ksp_type bicg -pc_type jacobi
737:     test:
738:       suffix: kokkos_batch
739:       requires: kokkos_kernels
740:       args: -dm_landau_device_type kokkos -dm_mat_type aijkokkos -dm_vec_type kokkos -ksp_type preonly -pc_type bjkokkos -pc_bjkokkos_ksp_type bicg -pc_bjkokkos_pc_type jacobi
741:     test:
742:       suffix: kokkos_batch_tfqmr
743:       requires: kokkos_kernels !defined(PETSC_HAVE_CUDA_CLANG)
744:       args: -dm_landau_device_type kokkos -dm_mat_type aijkokkos -dm_vec_type kokkos -ksp_type preonly -pc_type bjkokkos -pc_bjkokkos_ksp_type tfqmr -pc_bjkokkos_pc_type jacobi

746:   testset:
747:     requires: p4est !complex double defined(PETSC_USE_DMLANDAU_2D)
748:     output_file: output/ex2_deg4.out
749:     args: -dm_landau_num_species_grid 1,1 -dm_landau_Ez 0 -petscspace_degree 4 -petscspace_poly_tensor 1 -dm_landau_type p4est -dm_landau_ion_masses 2 \
750:     -dm_landau_ion_charges 1 -dm_landau_thermal_temps 5,5 -dm_landau_n 2,2 -dm_landau_n_0 5e19 -ts_monitor -snes_rtol 1.e-9 -snes_stol 1.e-14 -snes_monitor \
751:     -snes_converged_reason -snes_max_it 10 -ts_type arkimex -ts_arkimex_type 1bee -ts_max_snes_failures unlimited -ts_rtol 1e-3 -ts_time_step 1.e-2 -ts_max_time 1 -ts_adapt_clip .5,1.25 -ts_max_steps 2 -ts_adapt_scale_solve_failed 0.75 -ts_adapt_time_step_increase_delay 5 -dm_landau_amr_levels_max 1,2 -dm_landau_amr_re_levels 2 -dm_landau_re_radius 0 -ex2_impurity_source_type pulse -ex2_pulse_start_time 1e-1 -ex2_pulse_width_time 10 -ex2_pulse_rate 1e-2 -ex2_t_cold .05 -ex2_plot_dt 1e-1 -dm_refine 0 -dm_landau_gpu_assembly true -dm_landau_batch_size 2 -dm_landau_verbose 2 -dm_landau_domain_radius 5.,5.
752:     test:
753:       suffix: cpu_deg4
754:       args: -dm_landau_device_type cpu -ksp_type bicg -pc_type jacobi
755:     test:
756:       suffix: kokkos_deg4
757:       requires: kokkos_kernels !defined(PETSC_HAVE_CUDA_CLANG)
758:       args: -dm_landau_device_type kokkos -dm_mat_type aijkokkos -dm_vec_type kokkos -ksp_type bicg -pc_type jacobi

760:   test:
761:     requires: !complex double defined(PETSC_USE_DMLANDAU_2D) !cuda
762:     suffix: single
763:     nsize: 1
764:     args: -dm_refine 2 -dm_landau_num_species_grid 1 -dm_landau_thermal_temps 1 -dm_landau_electron_shift 1.25 -petscspace_degree 3 -snes_converged_reason -ts_type beuler -ts_time_step .1\
765:      -ex2_plot_dt .1 -ts_max_steps 1 -ex2_grid_view_idx 0 -ex2_dm_view -snes_rtol 1.e-13 -snes_stol 1.e-13 -dm_landau_verbose 2 -ex2_print_period 1 -ksp_type preonly -pc_type lu -dm_landau_device_type cpu -dm_landau_use_relativistic_corrections

767:   testset:
768:     requires: !complex double defined(PETSC_USE_DMLANDAU_2D)
769:     nsize: 1
770:     output_file: output/ex2_simplex.out
771:     args: -dim 2 -dm_landau_num_species_grid 1,1 -petscspace_degree 2 -dm_landau_simplex -dm_landau_ion_masses 2 -dm_landau_ion_charges 1 -dm_landau_thermal_temps 2,1 -dm_landau_n 1,1 -snes_rtol 1e-15\
772:      -snes_stol 1e-15 -snes_monitor -ts_type beuler -snes_converged_reason -ts_exact_final_time stepover -ts_time_step .1 -ts_max_steps 1 -ts_max_snes_failures unlimited -ksp_type preonly\
773:       -pc_type lu -dm_landau_verbose 2 -ex2_grid_view_idx 0 -ex2_dm_view -dm_refine 1 -ksp_type bicg -pc_type jacobi
774:     test:
775:       suffix: simplex
776:       args: -dm_landau_device_type cpu
777:     test:
778:       suffix: simplexkokkos
779:       requires: kokkos_kernels !defined(PETSC_HAVE_CUDA_CLANG) !sycl
780:       args: -dm_landau_device_type kokkos -dm_mat_type aijkokkos -dm_vec_type kokkos

782:   test:
783:     requires: double !defined(PETSC_USE_DMLANDAU_2D)
784:     suffix: sphere_3d
785:     nsize: 1
786:     args: -dim 3 -dm_landau_thermal_temps 2 -ts_type beuler -ts_time_step .1 -ts_max_steps 1 -dm_landau_verbose 2 -ksp_type preonly -pc_type lu -dm_landau_device_type cpu -snes_rtol 1.e-14 -snes_stol 1.e-14 -snes_converged_reason \
787:      -dm_landau_sphere -ex2_grid_view_idx 0 -ex2_dm_view -dm_landau_domain_radius 6 -dm_landau_sphere_inner_radius_90degree_scale .35 -petscspace_degree 3 -dm_refine 0 # -ex2_dm_view hdf5:my.hdf5:hdf5_viz -ex2_vec_view hdf5:my.hdf5:hdf5_viz:append

789: TEST*/