xref: /petsc/src/ts/utils/dmplexlandau/tutorials/ex2.c (revision cb359cd464fc40f56c50f7eb5c225ee3d03df42f)
1e0eea495SMark static char help[] = "Runaway electron model with Landau collision operator\n\n";
2e0eea495SMark 
3e0eea495SMark #include <petscdmplex.h>
4e0eea495SMark #include <petsclandau.h>
5e0eea495SMark #include <petscts.h>
6e0eea495SMark #include <petscds.h>
78a6f2e61SMark Adams #include <petscdmcomposite.h>
846233b44SBarry Smith #include <petsc/private/petscimpl.h>
9e0eea495SMark 
10419c2fb8Smarkadams4 #if defined(PETSC_HAVE_CUDA_NVTX)
11*cb359cd4SJunchao Zhang   #if PETSC_PKG_CUDA_VERSION_GE(10, 0, 0)
12*cb359cd4SJunchao Zhang     #include <nvtx3/nvToolsExt.h>
13*cb359cd4SJunchao Zhang   #else
14419c2fb8Smarkadams4     #include <nvToolsExt.h>
15419c2fb8Smarkadams4   #endif
16*cb359cd4SJunchao Zhang #endif
17419c2fb8Smarkadams4 
18e0eea495SMark /* data for runaway electron model */
19e0eea495SMark typedef struct REctx_struct {
208a6f2e61SMark Adams   PetscErrorCode (*test)(TS, Vec, PetscInt, PetscReal, PetscBool, LandauCtx *, struct REctx_struct *);
21e0eea495SMark   PetscErrorCode (*impuritySrcRate)(PetscReal, PetscReal *, LandauCtx *);
22e0eea495SMark   PetscErrorCode (*E)(Vec, Vec, PetscInt, PetscReal, LandauCtx *, PetscReal *);
23e0eea495SMark   PetscReal T_cold;        /* temperature of newly ionized electrons and impurity ions */
24e0eea495SMark   PetscReal ion_potential; /* ionization potential of impurity */
25e0eea495SMark   PetscReal Ne_ion;        /* effective number of electrons shed in ioization of impurity */
26e0eea495SMark   PetscReal Ez_initial;
27e0eea495SMark   PetscReal L; /* inductance */
28e0eea495SMark   Vec       X_0;
29e0eea495SMark   PetscInt  imp_idx; /* index for impurity ionizing sink */
30e0eea495SMark   PetscReal pulse_start;
31e0eea495SMark   PetscReal pulse_width;
32e0eea495SMark   PetscReal pulse_rate;
33e0eea495SMark   PetscReal current_rate;
34e0eea495SMark   PetscInt  plotIdx;
35e0eea495SMark   PetscInt  plotStep;
36e0eea495SMark   PetscInt  idx; /* cache */
37e0eea495SMark   PetscReal j;   /* cache */
38e0eea495SMark   PetscReal plotDt;
39e0eea495SMark   PetscBool plotting;
40e0eea495SMark   PetscBool use_spitzer_eta;
418a6f2e61SMark Adams   PetscInt  print_period;
428fdabdddSMark Adams   PetscInt  grid_view_idx;
43e0eea495SMark } REctx;
44e0eea495SMark 
45e0eea495SMark static const PetscReal kev_joul = 6.241506479963235e+15; /* 1/1000e */
46e0eea495SMark 
4790d163c9SMark Adams #define RE_CUT 3.
48e0eea495SMark /* < v, u_re * v * q > */
49d71ae5a4SJacob Faibussowitsch 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)
50d71ae5a4SJacob Faibussowitsch {
51e0eea495SMark   PetscReal n_e = PetscRealPart(u[0]);
52e0eea495SMark   if (dim == 2) {
53e0eea495SMark     if (x[1] > RE_CUT || x[1] < -RE_CUT) {                    /* simply a cutoff for REs. v_|| > 3 v(T_e) */
54e0eea495SMark       *f0 = n_e * 2. * PETSC_PI * x[0] * x[1] * constants[0]; /* n * r * v_|| * q */
55e0eea495SMark     } else {
56e0eea495SMark       *f0 = 0;
57e0eea495SMark     }
58e0eea495SMark   } else {
59e0eea495SMark     if (x[2] > RE_CUT || x[2] < -RE_CUT) { /* simply a cutoff for REs. v_|| > 3 v(T_e) */
60e0eea495SMark       *f0 = n_e * x[2] * constants[0];
61e0eea495SMark     } else {
62e0eea495SMark       *f0 = 0;
63e0eea495SMark     }
64e0eea495SMark   }
65e0eea495SMark }
66e0eea495SMark 
67e0eea495SMark /* sum < v, u*v*q > */
68d71ae5a4SJacob Faibussowitsch 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)
69d71ae5a4SJacob Faibussowitsch {
70e0eea495SMark   PetscInt ii;
71e0eea495SMark   f0[0] = 0;
72e0eea495SMark   if (dim == 2) {
738a6f2e61SMark Adams     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 */
74e0eea495SMark   } else {
758a6f2e61SMark Adams     for (ii = 0; ii < Nf; ii++) f0[0] += u[ii] * x[2] * q[ii]; /* n * v_|| * q  * v_0 */
76e0eea495SMark   }
77e0eea495SMark }
78e0eea495SMark 
79e0eea495SMark /* < v, n_e > */
80d71ae5a4SJacob Faibussowitsch 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)
81d71ae5a4SJacob Faibussowitsch {
82e0eea495SMark   PetscInt ii = (PetscInt)PetscRealPart(constants[0]);
83e0eea495SMark   if (dim == 2) f0[0] = 2. * PETSC_PI * x[0] * u[ii];
84ad540459SPierre Jolivet   else f0[0] = u[ii];
85e0eea495SMark }
86e0eea495SMark 
87e0eea495SMark /* < v, n_e v_|| > */
88d71ae5a4SJacob Faibussowitsch 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)
89d71ae5a4SJacob Faibussowitsch {
90e0eea495SMark   PetscInt ii = (PetscInt)PetscRealPart(constants[0]);
91e0eea495SMark   if (dim == 2) f0[0] = u[ii] * 2. * PETSC_PI * x[0] * x[1]; /* n r v_|| */
92ad540459SPierre Jolivet   else f0[0] = u[ii] * x[2];                                 /* n v_|| */
93e0eea495SMark }
94e0eea495SMark 
95e0eea495SMark /* < v, n_e (v-shift) > */
96d71ae5a4SJacob Faibussowitsch 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)
97d71ae5a4SJacob Faibussowitsch {
988a6f2e61SMark Adams   PetscReal vz = numConstants > 0 ? PetscRealPart(constants[0]) : 0;
99e0eea495SMark   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 */
100d71ae5a4SJacob Faibussowitsch   else {
101d71ae5a4SJacob Faibussowitsch     *f0 = u[0] * PetscSqrtReal(x[0] * x[0] + x[1] * x[1] + (x[2] - vz) * (x[2] - vz)); /* n v */
102d71ae5a4SJacob Faibussowitsch   }
103e0eea495SMark }
104e0eea495SMark 
105e0eea495SMark /* CalculateE - Calculate the electric field  */
106e0eea495SMark /*  T        -- Electron temperature  */
107e0eea495SMark /*  n        -- Electron density  */
108e0eea495SMark /*  lnLambda --   */
109e0eea495SMark /*  eps0     --  */
110e0eea495SMark /*  E        -- output E, input \hat E */
111d71ae5a4SJacob Faibussowitsch static PetscReal CalculateE(PetscReal Tev, PetscReal n, PetscReal lnLambda, PetscReal eps0, PetscReal *E)
112d71ae5a4SJacob Faibussowitsch {
113e0eea495SMark   PetscReal c, e, m;
114e0eea495SMark 
115e0eea495SMark   PetscFunctionBegin;
116cefb98e8SMark Adams   c = 299792458.0;
117e0eea495SMark   e = 1.602176e-19;
118e0eea495SMark   m = 9.10938e-31;
119e0eea495SMark   if (1) {
120e0eea495SMark     double Ec, Ehat = *E, betath = PetscSqrtReal(2 * Tev * e / (m * c * c)), j0 = Ehat * 7 / (PetscSqrtReal(2) * 2) * PetscPowReal(betath, 3) * n * e * c;
121e0eea495SMark     Ec = n * lnLambda * PetscPowReal(e, 3) / (4 * PETSC_PI * PetscPowReal(eps0, 2) * m * c * c);
122e0eea495SMark     *E = Ec;
1233ba16761SJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "CalculateE j0=%g Ec = %g\n", j0, Ec));
124e0eea495SMark   } else {
125e0eea495SMark     PetscReal Ed, vth;
126e0eea495SMark     vth = PetscSqrtReal(8 * Tev * e / (m * PETSC_PI));
127e0eea495SMark     Ed  = n * lnLambda * PetscPowReal(e, 3) / (4 * PETSC_PI * PetscPowReal(eps0, 2) * m * vth * vth);
128e0eea495SMark     *E  = Ed;
129e0eea495SMark   }
1303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
131e0eea495SMark }
132e0eea495SMark 
133d71ae5a4SJacob Faibussowitsch static PetscReal Spitzer(PetscReal m_e, PetscReal e, PetscReal Z, PetscReal epsilon0, PetscReal lnLam, PetscReal kTe_joules)
134d71ae5a4SJacob Faibussowitsch {
135e0eea495SMark   PetscReal Fz = (1 + 1.198 * Z + 0.222 * Z * Z) / (1 + 2.966 * Z + 0.753 * Z * Z), eta;
136e0eea495SMark   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);
137e0eea495SMark   return eta;
138e0eea495SMark }
139e0eea495SMark 
140d71ae5a4SJacob Faibussowitsch static PetscErrorCode testNone(TS ts, Vec X, PetscInt stepi, PetscReal time, PetscBool islast, LandauCtx *ctx, REctx *rectx)
141d71ae5a4SJacob Faibussowitsch {
142e0eea495SMark   PetscFunctionBeginUser;
1433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
144e0eea495SMark }
145e0eea495SMark 
146d71ae5a4SJacob Faibussowitsch static PetscErrorCode testSpitzer(TS ts, Vec X, PetscInt stepi, PetscReal time, PetscBool islast, LandauCtx *ctx, REctx *rectx)
147d71ae5a4SJacob Faibussowitsch {
148f37ccb5fSMark Adams   PetscInt          ii, nDMs;
149e0eea495SMark   PetscDS           prob;
150a587d139SMark   static PetscReal  old_ratio = 1e10;
151a587d139SMark   TSConvergedReason reason;
152a587d139SMark   PetscReal         J, J_re, spit_eta, Te_kev = 0, E, ratio, Z, n_e, v, v2;
1538a6f2e61SMark Adams   PetscScalar       user[2] = {0., ctx->charges[0]}, q[LANDAU_MAX_SPECIES], tt[LANDAU_MAX_SPECIES], vz;
154930e68a5SMark Adams   PetscReal         dt;
1558a6f2e61SMark Adams   DM                pack, plexe = ctx->plex[0], plexi = (ctx->num_grids == 1) ? NULL : ctx->plex[1];
156f37ccb5fSMark Adams   Vec              *XsubArray;
157e0eea495SMark 
158e0eea495SMark   PetscFunctionBeginUser;
15963a3b9bcSJacob Faibussowitsch   PetscCheck(ctx->num_species == 2, PETSC_COMM_SELF, PETSC_ERR_PLIB, "ctx->num_species %" PetscInt_FMT " != 2", ctx->num_species);
1609566063dSJacob Faibussowitsch   PetscCall(VecGetDM(X, &pack));
1613c633725SBarry Smith   PetscCheck(pack, PETSC_COMM_SELF, PETSC_ERR_PLIB, "no DM");
1629566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetNumberDM(pack, &nDMs));
16363a3b9bcSJacob Faibussowitsch   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);
1649566063dSJacob Faibussowitsch   PetscCall(PetscMalloc(sizeof(*XsubArray) * nDMs, &XsubArray));
1659566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetAccessArray(pack, X, nDMs, NULL, XsubArray)); // read only
1669566063dSJacob Faibussowitsch   PetscCall(TSGetTimeStep(ts, &dt));
1678a6f2e61SMark Adams   /* get current for each grid */
1688a6f2e61SMark Adams   for (ii = 0; ii < ctx->num_species; ii++) q[ii] = ctx->charges[ii];
1699566063dSJacob Faibussowitsch   PetscCall(DMGetDS(plexe, &prob));
1709566063dSJacob Faibussowitsch   PetscCall(PetscDSSetConstants(prob, 2, &q[0]));
1719566063dSJacob Faibussowitsch   PetscCall(PetscDSSetObjective(prob, 0, &f0_jz_sum));
1729566063dSJacob Faibussowitsch   PetscCall(DMPlexComputeIntegralFEM(plexe, XsubArray[LAND_PACK_IDX(ctx->batch_view_idx, 0)], tt, NULL));
173e0eea495SMark   J = -ctx->n_0 * ctx->v_0 * PetscRealPart(tt[0]);
1748fdabdddSMark Adams   if (plexi) { // add first (only) ion
1759566063dSJacob Faibussowitsch     PetscCall(DMGetDS(plexi, &prob));
1769566063dSJacob Faibussowitsch     PetscCall(PetscDSSetConstants(prob, 1, &q[1]));
1779566063dSJacob Faibussowitsch     PetscCall(PetscDSSetObjective(prob, 0, &f0_jz_sum));
1789566063dSJacob Faibussowitsch     PetscCall(DMPlexComputeIntegralFEM(plexi, XsubArray[LAND_PACK_IDX(ctx->batch_view_idx, 1)], tt, NULL));
1798a6f2e61SMark Adams     J += -ctx->n_0 * ctx->v_0 * PetscRealPart(tt[0]);
1808a6f2e61SMark Adams   }
181a587d139SMark   /* get N_e */
1829566063dSJacob Faibussowitsch   PetscCall(DMGetDS(plexe, &prob));
1839566063dSJacob Faibussowitsch   PetscCall(PetscDSSetConstants(prob, 1, user));
1849566063dSJacob Faibussowitsch   PetscCall(PetscDSSetObjective(prob, 0, &f0_n));
1859566063dSJacob Faibussowitsch   PetscCall(DMPlexComputeIntegralFEM(plexe, XsubArray[LAND_PACK_IDX(ctx->batch_view_idx, 0)], tt, NULL));
186a587d139SMark   n_e = PetscRealPart(tt[0]) * ctx->n_0;
187a587d139SMark   /* Z */
188a587d139SMark   Z = -ctx->charges[1] / ctx->charges[0];
189a587d139SMark   /* remove drift */
1908a6f2e61SMark Adams   if (0) {
1918a6f2e61SMark Adams     user[0] = 0; // electrons
1929566063dSJacob Faibussowitsch     PetscCall(DMGetDS(plexe, &prob));
1939566063dSJacob Faibussowitsch     PetscCall(PetscDSSetConstants(prob, 1, user));
1949566063dSJacob Faibussowitsch     PetscCall(PetscDSSetObjective(prob, 0, &f0_vz));
1959566063dSJacob Faibussowitsch     PetscCall(DMPlexComputeIntegralFEM(plexe, XsubArray[LAND_PACK_IDX(ctx->batch_view_idx, 0)], tt, NULL));
196a587d139SMark     vz = ctx->n_0 * PetscRealPart(tt[0]) / n_e; /* non-dimensional */
197a587d139SMark   } else vz = 0;
198a587d139SMark   /* thermal velocity */
1999566063dSJacob Faibussowitsch   PetscCall(DMGetDS(plexe, &prob));
2009566063dSJacob Faibussowitsch   PetscCall(PetscDSSetConstants(prob, 1, &vz));
2019566063dSJacob Faibussowitsch   PetscCall(PetscDSSetObjective(prob, 0, &f0_ve_shift));
2029566063dSJacob Faibussowitsch   PetscCall(DMPlexComputeIntegralFEM(plexe, XsubArray[LAND_PACK_IDX(ctx->batch_view_idx, 0)], tt, NULL));
203a587d139SMark   v        = ctx->n_0 * ctx->v_0 * PetscRealPart(tt[0]) / n_e;                                                   /* remove number density to get velocity */
204a587d139SMark   v2       = PetscSqr(v);                                                                                        /* use real space: m^2 / s^2 */
205a587d139SMark   Te_kev   = (v2 * ctx->masses[0] * PETSC_PI / 8) * kev_joul;                                                    /* temperature in kev */
206ae6bf4a8Smarkadams4   spit_eta = Spitzer(ctx->masses[0], -ctx->charges[0], Z, ctx->epsilon0, ctx->lambdas[0][1], Te_kev / kev_joul); /* kev --> J (kT) */
2078a6f2e61SMark Adams   if (0) {
2089566063dSJacob Faibussowitsch     PetscCall(DMGetDS(plexe, &prob));
2099566063dSJacob Faibussowitsch     PetscCall(PetscDSSetConstants(prob, 1, q));
2109566063dSJacob Faibussowitsch     PetscCall(PetscDSSetObjective(prob, 0, &f0_j_re));
2119566063dSJacob Faibussowitsch     PetscCall(DMPlexComputeIntegralFEM(plexe, XsubArray[LAND_PACK_IDX(ctx->batch_view_idx, 0)], tt, NULL));
212a587d139SMark   } else tt[0] = 0;
213e0eea495SMark   J_re = -ctx->n_0 * ctx->v_0 * PetscRealPart(tt[0]);
2149566063dSJacob Faibussowitsch   PetscCall(DMCompositeRestoreAccessArray(pack, X, nDMs, NULL, XsubArray)); // read only
2159566063dSJacob Faibussowitsch   PetscCall(PetscFree(XsubArray));
216a587d139SMark 
21790d163c9SMark Adams   if (rectx->use_spitzer_eta) {
21890d163c9SMark Adams     E = ctx->Ez = spit_eta * (rectx->j - J_re);
21990d163c9SMark Adams   } else {
22090d163c9SMark Adams     E        = ctx->Ez; /* keep real E */
22190d163c9SMark Adams     rectx->j = J;       /* cache */
22290d163c9SMark Adams   }
22390d163c9SMark Adams 
224e0eea495SMark   ratio = E / J / spit_eta;
225f4f49eeaSPierre Jolivet   if (stepi > 10 && !rectx->use_spitzer_eta && (old_ratio - ratio < 1.e-6)) {
22690d163c9SMark Adams     rectx->pulse_start     = time + 0.98 * dt;
227a587d139SMark     rectx->use_spitzer_eta = PETSC_TRUE;
228a587d139SMark   }
2299566063dSJacob Faibussowitsch   PetscCall(TSGetConvergedReason(ts, &reason));
2309566063dSJacob Faibussowitsch   PetscCall(TSGetConvergedReason(ts, &reason));
231f4f49eeaSPierre Jolivet   if (rectx->plotting || stepi == 0 || reason || rectx->pulse_start == time + 0.98 * dt) {
2329371c9d4SSatish Balay     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,
2339371c9d4SSatish Balay                           (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));
23463a3b9bcSJacob Faibussowitsch     PetscCheck(rectx->pulse_start != (time + 0.98 * dt), PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Spitzer complete ratio=%g", (double)ratio);
235930e68a5SMark Adams   }
236e0eea495SMark   old_ratio = ratio;
2373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
238e0eea495SMark }
239e0eea495SMark 
240e0eea495SMark static const double ppp = 2;
241d71ae5a4SJacob Faibussowitsch 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)
242d71ae5a4SJacob Faibussowitsch {
243e0eea495SMark   LandauCtx      *ctx   = (LandauCtx *)constants;
244e0eea495SMark   REctx          *rectx = (REctx *)ctx->data;
245e0eea495SMark   PetscInt        ii    = rectx->idx, i;
246e0eea495SMark   const PetscReal kT_m  = ctx->k * ctx->thermal_temps[ii] / ctx->masses[ii]; /* kT/m */
247e0eea495SMark   const PetscReal n     = ctx->n[ii];
248e0eea495SMark   PetscReal       diff, f_maxwell, v2 = 0, theta = 2 * kT_m / (ctx->v_0 * ctx->v_0); /* theta = 2kT/mc^2 */
249e0eea495SMark   for (i = 0; i < dim; ++i) v2 += x[i] * x[i];
250e0eea495SMark   f_maxwell = n * PetscPowReal(PETSC_PI * theta, -1.5) * (PetscExpReal(-v2 / theta));
251e0eea495SMark   diff      = 2. * PETSC_PI * x[0] * (PetscRealPart(u[ii]) - f_maxwell);
252e0eea495SMark   f0[0]     = PetscPowReal(diff, ppp);
253e0eea495SMark }
254d71ae5a4SJacob Faibussowitsch 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)
255d71ae5a4SJacob Faibussowitsch {
256e0eea495SMark   LandauCtx      *ctx   = (LandauCtx *)constants;
257e0eea495SMark   REctx          *rectx = (REctx *)ctx->data;
258e0eea495SMark   PetscInt        ii    = rectx->idx, i;
259e0eea495SMark   const PetscReal kT_m  = ctx->k * ctx->thermal_temps[ii] / ctx->masses[ii]; /* kT/m */
260e0eea495SMark   const PetscReal n     = ctx->n[ii];
261e0eea495SMark   PetscReal       f_maxwell, v2 = 0, theta = 2 * kT_m / (ctx->v_0 * ctx->v_0); /* theta = 2kT/mc^2 */
262e0eea495SMark   for (i = 0; i < dim; ++i) v2 += x[i] * x[i];
263e0eea495SMark   f_maxwell = 2. * PETSC_PI * x[0] * n * PetscPowReal(PETSC_PI * theta, -1.5) * (PetscExpReal(-v2 / theta));
264e0eea495SMark   f0[0]     = PetscPowReal(f_maxwell, ppp);
265e0eea495SMark }
266e0eea495SMark 
267d71ae5a4SJacob Faibussowitsch static PetscErrorCode testStable(TS ts, Vec X, PetscInt stepi, PetscReal time, PetscBool islast, LandauCtx *ctx, REctx *rectx)
268d71ae5a4SJacob Faibussowitsch {
269e0eea495SMark   PetscDS     prob;
270e0eea495SMark   Vec         X2;
271e0eea495SMark   PetscReal   ediff, idiff = 0, lpm0, lpm1 = 1;
272e0eea495SMark   PetscScalar tt[LANDAU_MAX_SPECIES];
2738a6f2e61SMark Adams   DM          dm, plex = ctx->plex[0];
274e0eea495SMark 
275e0eea495SMark   PetscFunctionBeginUser;
2769566063dSJacob Faibussowitsch   PetscCall(VecGetDM(X, &dm));
2779566063dSJacob Faibussowitsch   PetscCall(DMGetDS(plex, &prob));
2789566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(X, &X2));
2799566063dSJacob Faibussowitsch   PetscCall(VecCopy(X, X2));
280e0eea495SMark   if (!rectx->X_0) {
2819566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(X, &rectx->X_0));
2829566063dSJacob Faibussowitsch     PetscCall(VecCopy(X, rectx->X_0));
283e0eea495SMark   }
2849566063dSJacob Faibussowitsch   PetscCall(VecAXPY(X, -1.0, rectx->X_0));
2859566063dSJacob Faibussowitsch   PetscCall(PetscDSSetConstants(prob, sizeof(LandauCtx) / sizeof(PetscScalar), (PetscScalar *)ctx));
286e0eea495SMark   rectx->idx = 0;
2879566063dSJacob Faibussowitsch   PetscCall(PetscDSSetObjective(prob, 0, &f0_0_diff_lp));
2889566063dSJacob Faibussowitsch   PetscCall(DMPlexComputeIntegralFEM(plex, X2, tt, NULL));
289e0eea495SMark   ediff = PetscPowReal(PetscRealPart(tt[0]), 1. / ppp);
2909566063dSJacob Faibussowitsch   PetscCall(PetscDSSetObjective(prob, 0, &f0_0_maxwellian_lp));
2919566063dSJacob Faibussowitsch   PetscCall(DMPlexComputeIntegralFEM(plex, X2, tt, NULL));
292e0eea495SMark   lpm0 = PetscPowReal(PetscRealPart(tt[0]), 1. / ppp);
293e0eea495SMark   if (ctx->num_species > 1) {
294e0eea495SMark     rectx->idx = 1;
2959566063dSJacob Faibussowitsch     PetscCall(PetscDSSetObjective(prob, 0, &f0_0_diff_lp));
2969566063dSJacob Faibussowitsch     PetscCall(DMPlexComputeIntegralFEM(plex, X2, tt, NULL));
297e0eea495SMark     idiff = PetscPowReal(PetscRealPart(tt[0]), 1. / ppp);
2989566063dSJacob Faibussowitsch     PetscCall(PetscDSSetObjective(prob, 0, &f0_0_maxwellian_lp));
2999566063dSJacob Faibussowitsch     PetscCall(DMPlexComputeIntegralFEM(plex, X2, tt, NULL));
300e0eea495SMark     lpm1 = PetscPowReal(PetscRealPart(tt[0]), 1. / ppp);
301e0eea495SMark   }
30263a3b9bcSJacob Faibussowitsch   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)));
303e0eea495SMark   /* view */
3049566063dSJacob Faibussowitsch   PetscCall(VecCopy(X2, X));
3059566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&X2));
306e0eea495SMark   if (islast) {
3079566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&rectx->X_0));
308e0eea495SMark     rectx->X_0 = NULL;
309e0eea495SMark   }
3103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
311e0eea495SMark }
312e0eea495SMark 
313d71ae5a4SJacob Faibussowitsch static PetscErrorCode EInduction(Vec X, Vec X_t, PetscInt step, PetscReal time, LandauCtx *ctx, PetscReal *a_E)
314d71ae5a4SJacob Faibussowitsch {
315e0eea495SMark   REctx      *rectx = (REctx *)ctx->data;
316e0eea495SMark   PetscInt    ii;
317e0eea495SMark   DM          dm, plex;
3188a6f2e61SMark Adams   PetscScalar tt[LANDAU_MAX_SPECIES], qv0[LANDAU_MAX_SPECIES];
319e0eea495SMark   PetscReal   dJ_dt;
320e0eea495SMark   PetscDS     prob;
321e0eea495SMark 
322e0eea495SMark   PetscFunctionBeginUser;
3238a6f2e61SMark Adams   for (ii = 0; ii < ctx->num_species; ii++) qv0[ii] = ctx->charges[ii] * ctx->v_0;
3249566063dSJacob Faibussowitsch   PetscCall(VecGetDM(X, &dm));
3259566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &prob));
3269566063dSJacob Faibussowitsch   PetscCall(DMConvert(dm, DMPLEX, &plex));
327e0eea495SMark   /* get d current / dt */
3289566063dSJacob Faibussowitsch   PetscCall(PetscDSSetConstants(prob, ctx->num_species, qv0));
3299566063dSJacob Faibussowitsch   PetscCall(PetscDSSetObjective(prob, 0, &f0_jz_sum));
3303c633725SBarry Smith   PetscCheck(X_t, PETSC_COMM_SELF, PETSC_ERR_PLIB, "X_t");
3319566063dSJacob Faibussowitsch   PetscCall(DMPlexComputeIntegralFEM(plex, X_t, tt, NULL));
3328a6f2e61SMark Adams   dJ_dt = -ctx->n_0 * PetscRealPart(tt[0]) / ctx->t_0;
333e0eea495SMark   /* E induction */
334e0eea495SMark   *a_E = -rectx->L * dJ_dt + rectx->Ez_initial;
3359566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&plex));
3363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
337e0eea495SMark }
338e0eea495SMark 
339d71ae5a4SJacob Faibussowitsch static PetscErrorCode EConstant(Vec X, Vec X_t, PetscInt step, PetscReal time, LandauCtx *ctx, PetscReal *a_E)
340d71ae5a4SJacob Faibussowitsch {
341e0eea495SMark   PetscFunctionBeginUser;
342e0eea495SMark   *a_E = ctx->Ez;
3433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
344e0eea495SMark }
345e0eea495SMark 
346d71ae5a4SJacob Faibussowitsch static PetscErrorCode ENone(Vec X, Vec X_t, PetscInt step, PetscReal time, LandauCtx *ctx, PetscReal *a_E)
347d71ae5a4SJacob Faibussowitsch {
348e0eea495SMark   PetscFunctionBeginUser;
349e0eea495SMark   *a_E = 0;
3503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
351e0eea495SMark }
352e0eea495SMark 
353e0eea495SMark /* ------------------------------------------------------------------- */
354e0eea495SMark /*
355e0eea495SMark    FormSource - Evaluates source terms F(t).
356e0eea495SMark 
357e0eea495SMark    Input Parameters:
358e0eea495SMark .  ts - the TS context
359e0eea495SMark .  time -
360e0eea495SMark .  X_dummmy - input vector
361e0eea495SMark .  dummy - optional user-defined context, as set by SNESSetFunction()
362e0eea495SMark 
363e0eea495SMark    Output Parameter:
364e0eea495SMark .  F - function vector
365e0eea495SMark  */
366d71ae5a4SJacob Faibussowitsch static PetscErrorCode FormSource(TS ts, PetscReal ftime, Vec X_dummmy, Vec F, void *dummy)
367d71ae5a4SJacob Faibussowitsch {
368e0eea495SMark   PetscReal  new_imp_rate;
369e0eea495SMark   LandauCtx *ctx;
3708a6f2e61SMark Adams   DM         pack;
371e0eea495SMark   REctx     *rectx;
372e0eea495SMark 
373e0eea495SMark   PetscFunctionBeginUser;
3749566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &pack));
3759566063dSJacob Faibussowitsch   PetscCall(DMGetApplicationContext(pack, &ctx));
376e0eea495SMark   rectx = (REctx *)ctx->data;
377e0eea495SMark   /* check for impurities */
3789566063dSJacob Faibussowitsch   PetscCall(rectx->impuritySrcRate(ftime, &new_imp_rate, ctx));
379e0eea495SMark   if (new_imp_rate != 0) {
380e0eea495SMark     if (new_imp_rate != rectx->current_rate) {
381e0eea495SMark       PetscInt  ii;
382e0eea495SMark       PetscReal dne_dt, dni_dt, tilda_ns[LANDAU_MAX_SPECIES], temps[LANDAU_MAX_SPECIES];
3838fdabdddSMark Adams       Vec       globFarray[LANDAU_MAX_GRIDS * LANDAU_MAX_BATCH_SZ];
384e0eea495SMark       rectx->current_rate = new_imp_rate;
385e0eea495SMark       for (ii = 1; ii < LANDAU_MAX_SPECIES; ii++) tilda_ns[ii] = 0;
386e0eea495SMark       for (ii = 1; ii < LANDAU_MAX_SPECIES; ii++) temps[ii] = 1;
3878a6f2e61SMark Adams       dni_dt                   = new_imp_rate /* *ctx->t_0 */; /* fully ionized immediately, no normalize, stay in non-dim */
3888a6f2e61SMark Adams       dne_dt                   = new_imp_rate * rectx->Ne_ion /* *ctx->t_0 */;
3899371c9d4SSatish Balay       tilda_ns[0]              = dne_dt;
3909371c9d4SSatish Balay       tilda_ns[rectx->imp_idx] = dni_dt;
3919371c9d4SSatish Balay       temps[0]                 = rectx->T_cold;
3929371c9d4SSatish Balay       temps[rectx->imp_idx]    = rectx->T_cold;
39363a3b9bcSJacob Faibussowitsch       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));
3949566063dSJacob Faibussowitsch       PetscCall(DMCompositeGetAccessArray(pack, F, ctx->num_grids * ctx->batch_sz, NULL, globFarray));
3958a6f2e61SMark Adams       for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
396e0eea495SMark         /* add it */
397c3e4dd79SMark Adams         PetscCall(DMPlexLandauAddMaxwellians(ctx->plex[grid], globFarray[LAND_PACK_IDX(0, grid)], ftime, temps, tilda_ns, grid, 0, 1, ctx));
398e0eea495SMark       }
3998fdabdddSMark Adams       // Does DMCompositeRestoreAccessArray copy the data back? (no)
4009566063dSJacob Faibussowitsch       PetscCall(DMCompositeRestoreAccessArray(pack, F, ctx->num_grids * ctx->batch_sz, NULL, globFarray));
401e0eea495SMark     }
402e0eea495SMark   } else {
4039566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(F));
404e0eea495SMark     rectx->current_rate = 0;
405e0eea495SMark   }
4063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
407e0eea495SMark }
408419c2fb8Smarkadams4 
409d71ae5a4SJacob Faibussowitsch PetscErrorCode Monitor(TS ts, PetscInt stepi, PetscReal time, Vec X, void *actx)
410d71ae5a4SJacob Faibussowitsch {
411e0eea495SMark   LandauCtx        *ctx   = (LandauCtx *)actx; /* user-defined application context */
412e0eea495SMark   REctx            *rectx = (REctx *)ctx->data;
413419c2fb8Smarkadams4   DM                pack  = NULL;
4148fdabdddSMark Adams   Vec               globXArray[LANDAU_MAX_GRIDS * LANDAU_MAX_BATCH_SZ];
415e0eea495SMark   TSConvergedReason reason;
416419c2fb8Smarkadams4 
417e0eea495SMark   PetscFunctionBeginUser;
418419c2fb8Smarkadams4   PetscCall(TSGetConvergedReason(ts, &reason));
419419c2fb8Smarkadams4   if (rectx->grid_view_idx != -1 || (reason && ctx->verbose > 3)) {
4209566063dSJacob Faibussowitsch     PetscCall(VecGetDM(X, &pack));
4219566063dSJacob Faibussowitsch     PetscCall(DMCompositeGetAccessArray(pack, X, ctx->num_grids * ctx->batch_sz, NULL, globXArray));
422419c2fb8Smarkadams4   }
423e0eea495SMark   if (stepi > rectx->plotStep && rectx->plotting) {
424e0eea495SMark     rectx->plotting = PETSC_FALSE; /* was doing diagnostics, now done */
425e0eea495SMark     rectx->plotIdx++;
426e0eea495SMark   }
427e0eea495SMark   /* view */
428e0eea495SMark   if (time / rectx->plotDt >= (PetscReal)rectx->plotIdx || reason) {
429c3e4dd79SMark Adams     if ((reason || stepi == 0 || rectx->plotIdx % rectx->print_period == 0) && ctx->verbose > 1) {
430e0eea495SMark       /* print norms */
4319566063dSJacob Faibussowitsch       PetscCall(DMPlexLandauPrintNorms(X, stepi));
432a587d139SMark     }
433930e68a5SMark Adams     if (!rectx->plotting) { /* first step of possible backtracks */
434930e68a5SMark Adams       rectx->plotting = PETSC_TRUE;
435930e68a5SMark Adams       /* diagnostics + change E field with Sptizer (not just a monitor) */
4369566063dSJacob Faibussowitsch       PetscCall(rectx->test(ts, X, stepi, time, reason ? PETSC_TRUE : PETSC_FALSE, ctx, rectx));
437930e68a5SMark Adams     } else {
4383ba16761SJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\t\t ERROR SKIP test spit ------\n"));
439930e68a5SMark Adams       rectx->plotting = PETSC_TRUE;
440930e68a5SMark Adams     }
441419c2fb8Smarkadams4     if (rectx->grid_view_idx != -1) {
4429566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject)globXArray[LAND_PACK_IDX(ctx->batch_view_idx, rectx->grid_view_idx)], rectx->grid_view_idx == 0 ? "ue" : "ui"));
443930e68a5SMark Adams       /* view, overwrite step when back tracked */
444e5d13be4Smarkadams4       PetscCall(DMSetOutputSequenceNumber(ctx->plex[rectx->grid_view_idx], rectx->plotIdx, time * ctx->t_0));
445419c2fb8Smarkadams4       PetscCall(VecViewFromOptions(globXArray[LAND_PACK_IDX(ctx->batch_view_idx, rectx->grid_view_idx)], NULL, "-ex2_vec_view"));
446419c2fb8Smarkadams4     }
447e0eea495SMark     rectx->plotStep = stepi;
448930e68a5SMark Adams   } else {
4493ba16761SJacob Faibussowitsch     if (rectx->plotting) PetscCall(PetscPrintf(PETSC_COMM_WORLD, " ERROR rectx->plotting=%s step %" PetscInt_FMT "\n", PetscBools[rectx->plotting], stepi));
450930e68a5SMark Adams     /* diagnostics + change E field with Sptizer (not just a monitor) - can we lag this? */
4519566063dSJacob Faibussowitsch     PetscCall(rectx->test(ts, X, stepi, time, reason ? PETSC_TRUE : PETSC_FALSE, ctx, rectx));
452e0eea495SMark   }
453f37ccb5fSMark Adams   /* parallel check that only works of all batches are identical */
454ae6bf4a8Smarkadams4   if (reason && ctx->verbose > 3 && ctx->batch_sz > 1) {
455e0eea495SMark     PetscReal   val, rval;
456e0eea495SMark     PetscMPIInt rank;
4579566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
458f37ccb5fSMark Adams     for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
459f37ccb5fSMark Adams       PetscInt nerrors = 0;
460f37ccb5fSMark Adams       for (PetscInt i = 0; i < ctx->batch_sz; i++) {
4619566063dSJacob Faibussowitsch         PetscCall(VecNorm(globXArray[LAND_PACK_IDX(i, grid)], NORM_2, &val));
462f37ccb5fSMark Adams         if (i == 0) rval = val;
463f37ccb5fSMark Adams         else if ((val = PetscAbs(val - rval) / rval) > 1000 * PETSC_MACHINE_EPSILON) {
46463a3b9bcSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, " [%d] Warning %" PetscInt_FMT ".%" PetscInt_FMT ") diff = %2.15e\n", rank, grid, i, (double)val));
465f37ccb5fSMark Adams           nerrors++;
466f37ccb5fSMark Adams         }
467f37ccb5fSMark Adams       }
468f37ccb5fSMark Adams       if (nerrors) {
46963a3b9bcSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, " ***** [%d] ERROR max %" PetscInt_FMT " errors\n", rank, nerrors));
470e0eea495SMark       } else {
47163a3b9bcSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "[%d] %" PetscInt_FMT ") batch consistency check OK\n", rank, grid));
472f37ccb5fSMark Adams       }
473e0eea495SMark     }
474e0eea495SMark   }
475e0eea495SMark   rectx->idx = 0;
47648a46eb9SPierre Jolivet   if (rectx->grid_view_idx != -1 || (reason && ctx->verbose > 3)) PetscCall(DMCompositeRestoreAccessArray(pack, X, ctx->num_grids * ctx->batch_sz, NULL, globXArray));
4773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
478e0eea495SMark }
479e0eea495SMark 
480d71ae5a4SJacob Faibussowitsch PetscErrorCode PreStep(TS ts)
481d71ae5a4SJacob Faibussowitsch {
482e0eea495SMark   LandauCtx *ctx;
483e0eea495SMark   REctx     *rectx;
484e0eea495SMark   DM         dm;
485e0eea495SMark   PetscInt   stepi;
486e0eea495SMark   PetscReal  time;
487e0eea495SMark   Vec        X;
488e0eea495SMark 
489e0eea495SMark   PetscFunctionBeginUser;
490930e68a5SMark Adams   /* not used */
4919566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
4929566063dSJacob Faibussowitsch   PetscCall(TSGetTime(ts, &time));
4939566063dSJacob Faibussowitsch   PetscCall(TSGetSolution(ts, &X));
4949566063dSJacob Faibussowitsch   PetscCall(DMGetApplicationContext(dm, &ctx));
495e0eea495SMark   rectx = (REctx *)ctx->data;
4969566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts, &stepi));
497e0eea495SMark   /* update E */
4989566063dSJacob Faibussowitsch   PetscCall(rectx->E(X, NULL, stepi, time, ctx, &ctx->Ez));
4993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
500e0eea495SMark }
501e0eea495SMark 
502e0eea495SMark /* 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) */
503d71ae5a4SJacob Faibussowitsch static PetscErrorCode stepSrc(PetscReal time, PetscReal *rho, LandauCtx *ctx)
504d71ae5a4SJacob Faibussowitsch {
505e0eea495SMark   REctx *rectx = (REctx *)ctx->data;
506e0eea495SMark 
507e0eea495SMark   PetscFunctionBeginUser;
508e0eea495SMark   if (time >= rectx->pulse_start) *rho = rectx->pulse_rate;
509e0eea495SMark   else *rho = 0.;
5103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
511e0eea495SMark }
512d71ae5a4SJacob Faibussowitsch static PetscErrorCode zeroSrc(PetscReal time, PetscReal *rho, LandauCtx *ctx)
513d71ae5a4SJacob Faibussowitsch {
514e0eea495SMark   PetscFunctionBeginUser;
515e0eea495SMark   *rho = 0.;
5163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
517e0eea495SMark }
518d71ae5a4SJacob Faibussowitsch static PetscErrorCode pulseSrc(PetscReal time, PetscReal *rho, LandauCtx *ctx)
519d71ae5a4SJacob Faibussowitsch {
520e0eea495SMark   REctx *rectx = (REctx *)ctx->data;
521e0eea495SMark 
522e0eea495SMark   PetscFunctionBeginUser;
52308401ef6SPierre Jolivet   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'");
524e0eea495SMark   if (time < rectx->pulse_start || time > rectx->pulse_start + 3 * rectx->pulse_width) *rho = 0;
525a587d139SMark   else {
526e0eea495SMark     double x = PetscSinReal((time - rectx->pulse_start) / (3 * rectx->pulse_width) * 2 * PETSC_PI - PETSC_PI / 2) + 1; /* 0:2, integrates to 1.0 */
527e0eea495SMark     *rho     = rectx->pulse_rate * x / (3 * rectx->pulse_width);
528a587d139SMark     if (!rectx->use_spitzer_eta) rectx->use_spitzer_eta = PETSC_TRUE; /* use it next time */
529e0eea495SMark   }
5303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
531e0eea495SMark }
532e0eea495SMark 
533e0eea495SMark #undef __FUNCT__
534e0eea495SMark #define __FUNCT__ "ProcessREOptions"
535d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessREOptions(REctx *rectx, const LandauCtx *ctx, DM dm, const char prefix[])
536d71ae5a4SJacob Faibussowitsch {
537e0eea495SMark   PetscFunctionList plist = NULL, testlist = NULL, elist = NULL;
538e0eea495SMark   char              pname[256], testname[256], ename[256];
539930e68a5SMark Adams   DM                dm_dummy;
540a587d139SMark   PetscBool         Connor_E = PETSC_FALSE;
541e0eea495SMark 
542e0eea495SMark   PetscFunctionBeginUser;
5439566063dSJacob Faibussowitsch   PetscCall(DMCreate(PETSC_COMM_WORLD, &dm_dummy));
544e0eea495SMark   rectx->Ne_ion          = 1;    /* number of electrons given up by impurity ion */
545e0eea495SMark   rectx->T_cold          = .005; /* kev */
546e0eea495SMark   rectx->ion_potential   = 15;   /* ev */
547e0eea495SMark   rectx->L               = 2;
548e0eea495SMark   rectx->X_0             = NULL;
549e0eea495SMark   rectx->imp_idx         = ctx->num_species - 1; /* default ionized impurity as last one */
550a587d139SMark   rectx->pulse_start     = PETSC_MAX_REAL;
551e0eea495SMark   rectx->pulse_width     = 1;
5521690c2aeSBarry Smith   rectx->plotStep        = PETSC_INT_MAX;
553e0eea495SMark   rectx->pulse_rate      = 1.e-1;
554e0eea495SMark   rectx->current_rate    = 0;
555e0eea495SMark   rectx->plotIdx         = 0;
556e0eea495SMark   rectx->j               = 0;
557e0eea495SMark   rectx->plotDt          = 1.0;
558930e68a5SMark Adams   rectx->plotting        = PETSC_FALSE;
559e0eea495SMark   rectx->use_spitzer_eta = PETSC_FALSE;
560e0eea495SMark   rectx->idx             = 0;
5618a6f2e61SMark Adams   rectx->print_period    = 10;
562419c2fb8Smarkadams4   rectx->grid_view_idx   = -1; // do not get if not needed
563e0eea495SMark   /* Register the available impurity sources */
5649566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&plist, "step", &stepSrc));
5659566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&plist, "none", &zeroSrc));
5669566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&plist, "pulse", &pulseSrc));
567c6a7a370SJeremy L Thompson   PetscCall(PetscStrncpy(pname, "none", sizeof(pname)));
5689566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&testlist, "none", &testNone));
5699566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&testlist, "spitzer", &testSpitzer));
5709566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&testlist, "stable", &testStable));
571c6a7a370SJeremy L Thompson   PetscCall(PetscStrncpy(testname, "none", sizeof(testname)));
5729566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&elist, "none", &ENone));
5739566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&elist, "induction", &EInduction));
5749566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&elist, "constant", &EConstant));
575c6a7a370SJeremy L Thompson   PetscCall(PetscStrncpy(ename, "constant", sizeof(ename)));
576e0eea495SMark 
577d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_SELF, prefix, "Options for Runaway/seed electron model", "none");
5789566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-ex2_plot_dt", "Plotting interval", "ex2.c", rectx->plotDt, &rectx->plotDt, NULL));
579e0eea495SMark   if (rectx->plotDt < 0) rectx->plotDt = 1e30;
580e0eea495SMark   if (rectx->plotDt == 0) rectx->plotDt = 1e-30;
5819566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-ex2_print_period", "Plotting interval", "ex2.c", rectx->print_period, &rectx->print_period, NULL));
5829566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-ex2_grid_view_idx", "grid_view_idx", "ex2.c", rectx->grid_view_idx, &rectx->grid_view_idx, NULL));
583419c2fb8Smarkadams4   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);
5849566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFList("-ex2_impurity_source_type", "Name of impurity source to run", "", plist, pname, pname, sizeof(pname), NULL));
5859566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFList("-ex2_test_type", "Name of test to run", "", testlist, testname, testname, sizeof(testname), NULL));
5869566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-ex2_impurity_index", "index of sink for impurities", "none", rectx->imp_idx, &rectx->imp_idx, NULL));
587cad9d221SBarry Smith   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);
5889566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFList("-ex2_e_field_type", "Electric field type", "", elist, ename, ename, sizeof(ename), NULL));
589e0eea495SMark   rectx->Ne_ion = -ctx->charges[rectx->imp_idx] / ctx->charges[0];
5909566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-ex2_t_cold", "Temperature of cold electron and ions after ionization in keV", "none", rectx->T_cold, &rectx->T_cold, NULL));
5919566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-ex2_pulse_start_time", "Time at which pulse happens for 'pulse' source", "none", rectx->pulse_start, &rectx->pulse_start, NULL));
5929566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-ex2_pulse_width_time", "Width of pulse 'pulse' source", "none", rectx->pulse_width, &rectx->pulse_width, NULL));
5939566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-ex2_pulse_rate", "Number density of pulse for 'pulse' source", "none", rectx->pulse_rate, &rectx->pulse_rate, NULL));
594e0eea495SMark   rectx->T_cold *= 1.16e7; /* convert to Kelvin */
5959566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-ex2_ion_potential", "Potential to ionize impurity (should be array) in ev", "none", rectx->ion_potential, &rectx->ion_potential, NULL));
596da81f932SPierre Jolivet   PetscCall(PetscOptionsReal("-ex2_inductance", "Inductance E field", "none", rectx->L, &rectx->L, NULL));
5979566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-ex2_connor_e_field_units", "Scale Ex but Connor-Hastie E_c", "none", Connor_E, &Connor_E, NULL));
59863a3b9bcSJacob Faibussowitsch   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));
599d0609cedSBarry Smith   PetscOptionsEnd();
600e0eea495SMark   /* get impurity source rate function */
6019566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListFind(plist, pname, &rectx->impuritySrcRate));
6023c633725SBarry Smith   PetscCheck(rectx->impuritySrcRate, PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "No impurity source function found '%s'", pname);
6039566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListFind(testlist, testname, &rectx->test));
6043c633725SBarry Smith   PetscCheck(rectx->test, PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "No test found '%s'", testname);
6059566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListFind(elist, ename, &rectx->E));
6063c633725SBarry Smith   PetscCheck(rectx->E, PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "No E field function found '%s'", ename);
6079566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListDestroy(&plist));
6089566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListDestroy(&testlist));
6099566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListDestroy(&elist));
610930e68a5SMark Adams 
611a587d139SMark   /* convert E from Connor-Hastie E_c units to real if doing Spitzer E */
612a587d139SMark   if (Connor_E) {
613e0eea495SMark     PetscReal E = ctx->Ez, Tev = ctx->thermal_temps[0] * 8.621738e-5, n = ctx->n_0 * ctx->n[0];
614ae6bf4a8Smarkadams4     CalculateE(Tev, n, ctx->lambdas[0][1], ctx->epsilon0, &E);
615e0eea495SMark     ((LandauCtx *)ctx)->Ez *= E;
616e0eea495SMark   }
6179566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm_dummy));
6183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
619e0eea495SMark }
620e0eea495SMark 
621d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
622d71ae5a4SJacob Faibussowitsch {
6238a6f2e61SMark Adams   DM            pack;
624419c2fb8Smarkadams4   Vec           X;
625f37ccb5fSMark Adams   PetscInt      dim = 2, nDMs;
626e0eea495SMark   TS            ts;
627e0eea495SMark   Mat           J;
628e0eea495SMark   PetscDS       prob;
629e0eea495SMark   LandauCtx    *ctx;
630e0eea495SMark   REctx        *rectx;
6316c2b77d5SStefano Zampini   PetscMPIInt   rank;
632a587d139SMark   PetscLogStage stage;
6336c2b77d5SStefano Zampini 
634327415f7SBarry Smith   PetscFunctionBeginUser;
6359566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
6369566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
637930e68a5SMark Adams   if (rank) { /* turn off output stuff for duplicate runs */
638419c2fb8Smarkadams4     PetscCall(PetscOptionsClearValue(NULL, "-ex2_dm_view"));
639419c2fb8Smarkadams4     PetscCall(PetscOptionsClearValue(NULL, "-ex2_vec_view"));
640419c2fb8Smarkadams4     PetscCall(PetscOptionsClearValue(NULL, "-ex2_vec_view_init"));
641419c2fb8Smarkadams4     PetscCall(PetscOptionsClearValue(NULL, "-ex2_dm_view_init"));
6429566063dSJacob Faibussowitsch     PetscCall(PetscOptionsClearValue(NULL, "-info")); /* this does not work */
643930e68a5SMark Adams   }
6449566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-dim", &dim, NULL));
645e0eea495SMark   /* Create a mesh */
6469566063dSJacob Faibussowitsch   PetscCall(DMPlexLandauCreateVelocitySpace(PETSC_COMM_WORLD, dim, "", &X, &J, &pack));
6479566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetNumberDM(pack, &nDMs));
6489566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)J, "Jacobian"));
6499566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)X, "f"));
6509566063dSJacob Faibussowitsch   PetscCall(DMGetApplicationContext(pack, &ctx));
6519566063dSJacob Faibussowitsch   PetscCall(DMSetUp(pack));
652e0eea495SMark   /* context */
6539566063dSJacob Faibussowitsch   PetscCall(PetscNew(&rectx));
654e0eea495SMark   ctx->data = rectx;
6559566063dSJacob Faibussowitsch   PetscCall(ProcessREOptions(rectx, ctx, pack, ""));
6569566063dSJacob Faibussowitsch   PetscCall(DMGetDS(pack, &prob));
657419c2fb8Smarkadams4   if (rectx->grid_view_idx != -1) {
658419c2fb8Smarkadams4     Vec *XsubArray = NULL;
659419c2fb8Smarkadams4     PetscCall(PetscMalloc(sizeof(*XsubArray) * nDMs, &XsubArray));
6609566063dSJacob Faibussowitsch     PetscCall(DMCompositeGetAccessArray(pack, X, nDMs, NULL, XsubArray)); // read only
6619566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)XsubArray[LAND_PACK_IDX(ctx->batch_view_idx, rectx->grid_view_idx)], rectx->grid_view_idx == 0 ? "ue" : "ui"));
662e5d13be4Smarkadams4     PetscCall(DMSetOutputSequenceNumber(ctx->plex[rectx->grid_view_idx], 0, 0.0));
663419c2fb8Smarkadams4     PetscCall(DMViewFromOptions(ctx->plex[rectx->grid_view_idx], NULL, "-ex2_dm_view"));
664419c2fb8Smarkadams4     PetscCall(DMViewFromOptions(ctx->plex[rectx->grid_view_idx], NULL, "-ex2_dm_view_init"));
665e5d13be4Smarkadams4     PetscCall(VecViewFromOptions(XsubArray[LAND_PACK_IDX(ctx->batch_view_idx, rectx->grid_view_idx)], NULL, "-ex2_vec_view"));      // initial condition (monitor plots after step)
666419c2fb8Smarkadams4     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)
6679566063dSJacob Faibussowitsch     PetscCall(DMCompositeRestoreAccessArray(pack, X, nDMs, NULL, XsubArray));                                                       // read only
6689566063dSJacob Faibussowitsch     PetscCall(PetscFree(XsubArray));
669419c2fb8Smarkadams4   }
670e0eea495SMark   /* Create timestepping solver context */
6719566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_SELF, &ts));
6729566063dSJacob Faibussowitsch   PetscCall(TSSetDM(ts, pack));
6739566063dSJacob Faibussowitsch   PetscCall(TSSetIFunction(ts, NULL, DMPlexLandauIFunction, NULL));
6749566063dSJacob Faibussowitsch   PetscCall(TSSetIJacobian(ts, J, J, DMPlexLandauIJacobian, NULL));
6759566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(ts, NULL, FormSource, NULL));
6769566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
6779566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts, X));
6789566063dSJacob Faibussowitsch   PetscCall(TSSetApplicationContext(ts, ctx));
6799566063dSJacob Faibussowitsch   PetscCall(TSMonitorSet(ts, Monitor, ctx, NULL));
6809566063dSJacob Faibussowitsch   PetscCall(TSSetPreStep(ts, PreStep));
681da81f932SPierre Jolivet   rectx->Ez_initial = ctx->Ez; /* cache for induction calculation - applied E field */
6828594ddcfSMark Adams   if (1) {                     /* warm up an test just DMPlexLandauIJacobian */
683e0eea495SMark     Vec       vec;
684a587d139SMark     PetscInt  nsteps;
685a587d139SMark     PetscReal dt;
6869566063dSJacob Faibussowitsch     PetscCall(PetscLogStageRegister("Warmup", &stage));
6879566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePush(stage));
6889566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(X, &vec));
6899566063dSJacob Faibussowitsch     PetscCall(VecCopy(X, vec));
6909566063dSJacob Faibussowitsch     PetscCall(TSGetMaxSteps(ts, &nsteps));
6919566063dSJacob Faibussowitsch     PetscCall(TSGetTimeStep(ts, &dt));
6929566063dSJacob Faibussowitsch     PetscCall(TSSetMaxSteps(ts, 1));
6939566063dSJacob Faibussowitsch     PetscCall(TSSolve(ts, X));
6949566063dSJacob Faibussowitsch     PetscCall(TSSetMaxSteps(ts, nsteps));
6959566063dSJacob Faibussowitsch     PetscCall(TSSetStepNumber(ts, 0));
6969566063dSJacob Faibussowitsch     PetscCall(TSSetTime(ts, 0));
6979566063dSJacob Faibussowitsch     PetscCall(TSSetTimeStep(ts, dt));
698a587d139SMark     rectx->plotIdx  = 0;
699930e68a5SMark Adams     rectx->plotting = PETSC_FALSE;
7009566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePop());
7019566063dSJacob Faibussowitsch     PetscCall(VecCopy(vec, X));
7029566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&vec));
703c3e4dd79SMark Adams     PetscCall(PetscObjectStateIncrease((PetscObject)ctx->J));
704e0eea495SMark   }
705e0eea495SMark   /* go */
7069566063dSJacob Faibussowitsch   PetscCall(PetscLogStageRegister("Solve", &stage));
707f37ccb5fSMark Adams   ctx->stage = 0; // lets not use this stage
7089566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePush(stage));
709419c2fb8Smarkadams4 #if defined(PETSC_HAVE_CUDA_NVTX)
710419c2fb8Smarkadams4   nvtxRangePushA("ex2-TSSolve-warm");
711419c2fb8Smarkadams4 #endif
7129566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, X));
713419c2fb8Smarkadams4 #if defined(PETSC_HAVE_CUDA_NVTX)
714419c2fb8Smarkadams4   nvtxRangePop();
715419c2fb8Smarkadams4 #endif
7169566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePop());
717e0eea495SMark   /* clean up */
7189566063dSJacob Faibussowitsch   PetscCall(DMPlexLandauDestroyVelocitySpace(&pack));
7199566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
7209566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&X));
7219566063dSJacob Faibussowitsch   PetscCall(PetscFree(rectx));
7229566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
723b122ec5aSJacob Faibussowitsch   return 0;
724e0eea495SMark }
725e0eea495SMark 
726e0eea495SMark /*TEST
727a587d139SMark 
7285dac466eSMark Adams   testset:
729984ed092SMark Adams     requires: p4est !complex double defined(PETSC_USE_DMLANDAU_2D)
7305dac466eSMark Adams     output_file: output/ex2_0.out
73109cb0f53SBarry Smith     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 -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 -snes_converged_reason -snes_max_it 10 -ts_type arkimex -ts_arkimex_type 1bee -ts_max_snes_failures unlimited -ts_rtol 1e-3 -ts_dt 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.
7325dac466eSMark Adams     test:
7335dac466eSMark Adams       suffix: cpu
734c3e4dd79SMark Adams       args: -dm_landau_device_type cpu -ksp_type bicg -pc_type jacobi
735a587d139SMark     test:
736a587d139SMark       suffix: kokkos
737bd4f7539SJunchao Zhang       requires: kokkos_kernels !defined(PETSC_HAVE_CUDA_CLANG)
738c3e4dd79SMark Adams       args: -dm_landau_device_type kokkos -dm_mat_type aijkokkos -dm_vec_type kokkos -ksp_type bicg -pc_type jacobi
73990d163c9SMark Adams     test:
740cb25d741SMark Adams       suffix: kokkos_batch
741cb25d741SMark Adams       requires: kokkos_kernels
742c3e4dd79SMark Adams       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
743bfc784b7SMark Adams     test:
744a4313204SMark Adams       suffix: kokkos_batch_tfqmr
745bd4f7539SJunchao Zhang       requires: kokkos_kernels !defined(PETSC_HAVE_CUDA_CLANG)
746a4313204SMark Adams       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
74790d163c9SMark Adams 
7486b664d00Smarkadams4   test:
7499bdd0a17SMark Adams     requires: !complex double defined(PETSC_USE_DMLANDAU_2D) !cuda
7506b664d00Smarkadams4     suffix: single
7516b664d00Smarkadams4     nsize: 1
7526b664d00Smarkadams4     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_dt .1 -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
7536b664d00Smarkadams4 
754cd27c6deSmarkadams4   testset:
755cd27c6deSmarkadams4     requires: !complex double defined(PETSC_USE_DMLANDAU_2D)
756cd27c6deSmarkadams4     nsize: 1
757cd27c6deSmarkadams4     output_file: output/ex2_simplex.out
75809cb0f53SBarry Smith     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 -snes_stol 1e-15 -snes_monitor -ts_type beuler -snes_converged_reason -ts_exact_final_time stepover -ts_dt .1 -ts_max_steps 1 -ts_max_snes_failures unlimited -ksp_type preonly -pc_type lu -dm_landau_verbose 2 -ex2_grid_view_idx 0 -ex2_dm_view -dm_refine 1 -ksp_type bicg -pc_type jacobi
759cd27c6deSmarkadams4     test:
760cd27c6deSmarkadams4       suffix: simplex
761cd27c6deSmarkadams4       args: -dm_landau_device_type cpu
762cd27c6deSmarkadams4     test:
763cd27c6deSmarkadams4       suffix: simplexkokkos
764bbdffb7fSJunchao Zhang       requires: kokkos_kernels !defined(PETSC_HAVE_CUDA_CLANG) !sycl
765cd27c6deSmarkadams4       args: -dm_landau_device_type kokkos -dm_mat_type aijkokkos -dm_vec_type kokkos
766cd27c6deSmarkadams4 
767d043ef4cSMark Adams   test:
768d043ef4cSMark Adams     requires: double !defined(PETSC_USE_DMLANDAU_2D)
769d043ef4cSMark Adams     suffix: sphere_3d
770d043ef4cSMark Adams     nsize: 1
771d043ef4cSMark Adams     args: -dim 3 -dm_landau_thermal_temps 2 -petscspace_degree 2 -ts_type beuler -ts_dt .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 -dm_landau_sphere -dm_landau_sphere_inner_radius_90degree_scale .55 -dm_landau_sphere_inner_radius_45degree_scale .5
772d043ef4cSMark Adams 
773e0eea495SMark TEST*/
774