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>

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

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

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

 46: #define RE_CUT 3.
 47: /* < v, u_re * v * q > */
 48: 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)
 49: {
 50:   PetscReal n_e = PetscRealPart(u[0]);
 51:   if (dim == 2) {
 52:     if (x[1] > RE_CUT || x[1] < -RE_CUT) {                    /* simply a cutoff for REs. v_|| > 3 v(T_e) */
 53:       *f0 = n_e * 2. * PETSC_PI * x[0] * x[1] * constants[0]; /* n * r * v_|| * q */
 54:     } else {
 55:       *f0 = 0;
 56:     }
 57:   } else {
 58:     if (x[2] > RE_CUT || x[2] < -RE_CUT) { /* simply a cutoff for REs. v_|| > 3 v(T_e) */
 59:       *f0 = n_e * x[2] * constants[0];
 60:     } else {
 61:       *f0 = 0;
 62:     }
 63:   }
 64: }

 66: /* sum < v, u*v*q > */
 67: 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)
 68: {
 69:   PetscInt ii;
 70:   f0[0] = 0;
 71:   if (dim == 2) {
 72:     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 */
 73:   } else {
 74:     for (ii = 0; ii < Nf; ii++) f0[0] += u[ii] * x[2] * q[ii]; /* n * v_|| * q  * v_0 */
 75:   }
 76: }

 78: /* < v, n_e > */
 79: 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)
 80: {
 81:   PetscInt ii = (PetscInt)PetscRealPart(constants[0]);
 82:   if (dim == 2) f0[0] = 2. * PETSC_PI * x[0] * u[ii];
 83:   else f0[0] = u[ii];
 84: }

 86: /* < v, n_e v_|| > */
 87: 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)
 88: {
 89:   PetscInt ii = (PetscInt)PetscRealPart(constants[0]);
 90:   if (dim == 2) f0[0] = u[ii] * 2. * PETSC_PI * x[0] * x[1]; /* n r v_|| */
 91:   else f0[0] = u[ii] * x[2];                                 /* n v_|| */
 92: }

 94: /* < v, n_e (v-shift) > */
 95: 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)
 96: {
 97:   PetscReal vz = numConstants > 0 ? PetscRealPart(constants[0]) : 0;
 98:   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 */
 99:   else *f0 = u[0] * PetscSqrtReal(x[0] * x[0] + x[1] * x[1] + (x[2] - vz) * (x[2] - vz));                   /* n v */
100: }

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

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

130: static PetscReal Spitzer(PetscReal m_e, PetscReal e, PetscReal Z, PetscReal epsilon0, PetscReal lnLam, PetscReal kTe_joules)
131: {
132:   PetscReal Fz = (1 + 1.198 * Z + 0.222 * Z * Z) / (1 + 2.966 * Z + 0.753 * Z * Z), eta;
133:   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);
134:   return eta;
135: }

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

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

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

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

221:   ratio = E / J / spit_eta;
222:   if (stepi > 10 && !rectx->use_spitzer_eta && (old_ratio - ratio < 1.e-6)) {
223:     rectx->pulse_start     = time + 0.98 * dt;
224:     rectx->use_spitzer_eta = PETSC_TRUE;
225:   }
226:   PetscCall(TSGetConvergedReason(ts, &reason));
227:   PetscCall(TSGetConvergedReason(ts, &reason));
228:   if (rectx->plotting || stepi == 0 || reason || rectx->pulse_start == time + 0.98 * dt) {
229:     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,
230:                           (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));
231:     PetscCheck(rectx->pulse_start != (time + 0.98 * dt), PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Spitzer complete ratio=%g", (double)ratio);
232:   }
233:   old_ratio = ratio;
234:   PetscFunctionReturn(PETSC_SUCCESS);
235: }

237: static const double ppp = 2;
238: 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)
239: {
240:   LandauCtx      *ctx   = (LandauCtx *)constants;
241:   REctx          *rectx = (REctx *)ctx->data;
242:   PetscInt        ii    = rectx->idx, i;
243:   const PetscReal kT_m  = ctx->k * ctx->thermal_temps[ii] / ctx->masses[ii]; /* kT/m */
244:   const PetscReal n     = ctx->n[ii];
245:   PetscReal       diff, f_maxwell, v2 = 0, theta = 2 * kT_m / (ctx->v_0 * ctx->v_0); /* theta = 2kT/mc^2 */
246:   for (i = 0; i < dim; ++i) v2 += x[i] * x[i];
247:   f_maxwell = n * PetscPowReal(PETSC_PI * theta, -1.5) * (PetscExpReal(-v2 / theta));
248:   diff      = 2. * PETSC_PI * x[0] * (PetscRealPart(u[ii]) - f_maxwell);
249:   f0[0]     = PetscPowReal(diff, ppp);
250: }

252: 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)
253: {
254:   LandauCtx      *ctx   = (LandauCtx *)constants;
255:   REctx          *rectx = (REctx *)ctx->data;
256:   PetscInt        ii    = rectx->idx, i;
257:   const PetscReal kT_m  = ctx->k * ctx->thermal_temps[ii] / ctx->masses[ii]; /* kT/m */
258:   const PetscReal n     = ctx->n[ii];
259:   PetscReal       f_maxwell, v2 = 0, theta = 2 * kT_m / (ctx->v_0 * ctx->v_0); /* theta = 2kT/mc^2 */
260:   for (i = 0; i < dim; ++i) v2 += x[i] * x[i];
261:   f_maxwell = 2. * PETSC_PI * x[0] * n * PetscPowReal(PETSC_PI * theta, -1.5) * (PetscExpReal(-v2 / theta));
262:   f0[0]     = PetscPowReal(f_maxwell, ppp);
263: }

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

273:   PetscFunctionBeginUser;
274:   PetscCall(VecGetDM(X, &dm));
275:   PetscCall(DMGetDS(plex, &prob));
276:   PetscCall(VecDuplicate(X, &X2));
277:   PetscCall(VecCopy(X, X2));
278:   if (!rectx->X_0) {
279:     PetscCall(VecDuplicate(X, &rectx->X_0));
280:     PetscCall(VecCopy(X, rectx->X_0));
281:   }
282:   PetscCall(VecAXPY(X, -1.0, rectx->X_0));
283:   PetscCall(PetscDSSetConstants(prob, sizeof(LandauCtx) / sizeof(PetscScalar), (PetscScalar *)ctx));
284:   rectx->idx = 0;
285:   PetscCall(PetscDSSetObjective(prob, 0, &f0_0_diff_lp));
286:   PetscCall(DMPlexComputeIntegralFEM(plex, X2, tt, NULL));
287:   ediff = PetscPowReal(PetscRealPart(tt[0]), 1. / ppp);
288:   PetscCall(PetscDSSetObjective(prob, 0, &f0_0_maxwellian_lp));
289:   PetscCall(DMPlexComputeIntegralFEM(plex, X2, tt, NULL));
290:   lpm0 = PetscPowReal(PetscRealPart(tt[0]), 1. / ppp);
291:   if (ctx->num_species > 1) {
292:     rectx->idx = 1;
293:     PetscCall(PetscDSSetObjective(prob, 0, &f0_0_diff_lp));
294:     PetscCall(DMPlexComputeIntegralFEM(plex, X2, tt, NULL));
295:     idiff = PetscPowReal(PetscRealPart(tt[0]), 1. / ppp);
296:     PetscCall(PetscDSSetObjective(prob, 0, &f0_0_maxwellian_lp));
297:     PetscCall(DMPlexComputeIntegralFEM(plex, X2, tt, NULL));
298:     lpm1 = PetscPowReal(PetscRealPart(tt[0]), 1. / ppp);
299:   }
300:   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)));
301:   /* view */
302:   PetscCall(VecCopy(X2, X));
303:   PetscCall(VecDestroy(&X2));
304:   if (islast) {
305:     PetscCall(VecDestroy(&rectx->X_0));
306:     rectx->X_0 = NULL;
307:   }
308:   PetscFunctionReturn(PETSC_SUCCESS);
309: }

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

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

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

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

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

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

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

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

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

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

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

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

498: /* 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) */
499: static PetscErrorCode stepSrc(PetscReal time, PetscReal *rho, LandauCtx *ctx)
500: {
501:   REctx *rectx = (REctx *)ctx->data;

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

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

519:   PetscFunctionBeginUser;
520:   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'");
521:   if (time < rectx->pulse_start || time > rectx->pulse_start + 3 * rectx->pulse_width) *rho = 0;
522:   else {
523:     double x = PetscSinReal((time - rectx->pulse_start) / (3 * rectx->pulse_width) * 2 * PETSC_PI - PETSC_PI / 2) + 1; /* 0:2, integrates to 1.0 */
524:     *rho     = rectx->pulse_rate * x / (3 * rectx->pulse_width);
525:     if (!rectx->use_spitzer_eta) rectx->use_spitzer_eta = PETSC_TRUE; /* use it next time */
526:   }
527:   PetscFunctionReturn(PETSC_SUCCESS);
528: }

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

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

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

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

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

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

721: /*TEST

723:   testset:
724:     requires: p4est !complex double defined(PETSC_USE_DMLANDAU_2D)
725:     output_file: output/ex2_0.out
726:     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 \
727:     -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 \
728:     -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.
729:     test:
730:       suffix: cpu
731:       args: -dm_landau_device_type cpu -ksp_type bicg -pc_type jacobi
732:     test:
733:       suffix: kokkos
734:       requires: kokkos_kernels !defined(PETSC_HAVE_CUDA_CLANG)
735:       args: -dm_landau_device_type kokkos -dm_mat_type aijkokkos -dm_vec_type kokkos -ksp_type bicg -pc_type jacobi
736:     test:
737:       suffix: kokkos_batch
738:       requires: kokkos_kernels
739:       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
740:     test:
741:       suffix: kokkos_batch_tfqmr
742:       requires: kokkos_kernels !defined(PETSC_HAVE_CUDA_CLANG)
743:       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

745:   testset:
746:     requires: p4est !complex double defined(PETSC_USE_DMLANDAU_2D)
747:     output_file: output/ex2_deg4.out
748:     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 \
749:     -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 \
750:     -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.
751:     test:
752:       suffix: cpu_deg4
753:       args: -dm_landau_device_type cpu -ksp_type bicg -pc_type jacobi
754:     test:
755:       suffix: kokkos_deg4
756:       requires: kokkos_kernels !defined(PETSC_HAVE_CUDA_CLANG)
757:       args: -dm_landau_device_type kokkos -dm_mat_type aijkokkos -dm_vec_type kokkos -ksp_type bicg -pc_type jacobi

759:   test:
760:     requires: !complex double defined(PETSC_USE_DMLANDAU_2D) !cuda
761:     suffix: single
762:     nsize: 1
763:     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\
764:      -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

766:   testset:
767:     requires: !complex double defined(PETSC_USE_DMLANDAU_2D)
768:     nsize: 1
769:     output_file: output/ex2_simplex.out
770:     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\
771:      -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\
772:       -pc_type lu -dm_landau_verbose 2 -ex2_grid_view_idx 0 -ex2_dm_view -dm_refine 1 -ksp_type bicg -pc_type jacobi
773:     test:
774:       suffix: simplex
775:       args: -dm_landau_device_type cpu
776:     test:
777:       suffix: simplexkokkos
778:       requires: kokkos_kernels !defined(PETSC_HAVE_CUDA_CLANG) !sycl
779:       args: -dm_landau_device_type kokkos -dm_mat_type aijkokkos -dm_vec_type kokkos

781:   test:
782:     requires: double !defined(PETSC_USE_DMLANDAU_2D)
783:     suffix: sphere_3d
784:     nsize: 1
785:     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 \
786:      -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

788: TEST*/