xref: /petsc/src/ts/utils/dmplexlandau/tutorials/ex2.c (revision 7a46b59513fac54acb317385fa9864cc0988b1fb)
1 static char help[] = "Runaway electron model with Landau collision operator\n\n";
2 
3 #include <petscdmplex.h>
4 #include <petsclandau.h>
5 #include <petscts.h>
6 #include <petscds.h>
7 #include <petscdmcomposite.h>
8 #include "petsc/private/petscimpl.h"
9 
10 #if defined(PETSC_HAVE_CUDA_NVTX)
11 #include <nvToolsExt.h>
12 #endif
13 
14 /* data for runaway electron model */
15 typedef struct REctx_struct {
16   PetscErrorCode (*test)(TS, Vec, PetscInt, PetscReal, PetscBool,  LandauCtx *, struct REctx_struct *);
17   PetscErrorCode (*impuritySrcRate)(PetscReal, PetscReal *, LandauCtx*);
18   PetscErrorCode (*E)(Vec, Vec, PetscInt, PetscReal, LandauCtx*, PetscReal *);
19   PetscReal     T_cold;        /* temperature of newly ionized electrons and impurity ions */
20   PetscReal     ion_potential; /* ionization potential of impurity */
21   PetscReal     Ne_ion;        /* effective number of electrons shed in ioization of impurity */
22   PetscReal     Ez_initial;
23   PetscReal     L;             /* inductance */
24   Vec           X_0;
25   PetscInt      imp_idx;       /* index for impurity ionizing sink */
26   PetscReal     pulse_start;
27   PetscReal     pulse_width;
28   PetscReal     pulse_rate;
29   PetscReal     current_rate;
30   PetscInt      plotIdx;
31   PetscInt      plotStep;
32   PetscInt      idx; /* cache */
33   PetscReal     j; /* cache */
34   PetscReal     plotDt;
35   PetscBool     plotting;
36   PetscBool     use_spitzer_eta;
37   PetscInt      print_period;
38   PetscInt      grid_view_idx;
39 } REctx;
40 
41 static const PetscReal kev_joul = 6.241506479963235e+15; /* 1/1000e */
42 
43 #define RE_CUT 3.
44 /* < v, u_re * v * q > */
45 static void f0_j_re(PetscInt dim, PetscInt Nf, PetscInt NfAux,
46                     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
47                     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
48                     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 }
65 
66 /* sum < v, u*v*q > */
67 static void f0_jz_sum(PetscInt dim, PetscInt Nf, PetscInt NfAux,
68                       const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
69                       const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
70                       PetscReal t, const PetscReal x[],  PetscInt numConstants, const PetscScalar q[], PetscScalar *f0)
71 {
72   PetscInt ii;
73   f0[0] = 0;
74   if (dim==2) {
75     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 */
76   } else {
77     for (ii=0;ii<Nf;ii++) f0[0] += u[ii]                * x[2] * q[ii]; /* n * v_|| * q  * v_0 */
78   }
79 }
80 
81 /* < v, n_e > */
82 static void f0_n(PetscInt dim, PetscInt Nf, PetscInt NfAux,
83                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
84                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
85                   PetscReal t, const PetscReal x[],  PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0)
86 {
87   PetscInt ii = (PetscInt)PetscRealPart(constants[0]);
88   if (dim==2) f0[0] = 2.*PETSC_PI*x[0]*u[ii];
89   else {
90     f0[0] =                        u[ii];
91   }
92 }
93 
94 /* < v, n_e v_|| > */
95 static void f0_vz(PetscInt dim, PetscInt Nf, PetscInt NfAux,
96                    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
97                    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
98                    PetscReal t, const PetscReal x[],  PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0)
99 {
100   PetscInt ii = (PetscInt)PetscRealPart(constants[0]);
101   if (dim==2) f0[0] = u[ii] * 2.*PETSC_PI*x[0] * x[1]; /* n r v_|| */
102   else {
103     f0[0] =           u[ii] *                x[2]; /* n v_|| */
104   }
105 }
106 
107 /* < v, n_e (v-shift) > */
108 static void f0_ve_shift(PetscInt dim, PetscInt Nf, PetscInt NfAux,
109                          const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
110                          const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
111                          PetscReal t, const PetscReal x[],  PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0)
112 {
113   PetscReal vz = numConstants>0 ? PetscRealPart(constants[0]) : 0;
114   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 */
115   else {
116     *f0 =           u[0] *                PetscSqrtReal(x[0]*x[0] + x[1]*x[1] + (x[2]-vz)*(x[2]-vz)); /* n v */
117   }
118 }
119 
120  /* CalculateE - Calculate the electric field  */
121  /*  T        -- Electron temperature  */
122  /*  n        -- Electron density  */
123  /*  lnLambda --   */
124  /*  eps0     --  */
125  /*  E        -- output E, input \hat E */
126 static PetscReal CalculateE(PetscReal Tev, PetscReal n, PetscReal lnLambda, PetscReal eps0, PetscReal *E)
127 {
128   PetscReal c,e,m;
129 
130   PetscFunctionBegin;
131   c = 299792458.0;
132   e = 1.602176e-19;
133   m = 9.10938e-31;
134   if (1) {
135     double Ec, Ehat = *E, betath = PetscSqrtReal(2*Tev*e/(m*c*c)), j0 = Ehat * 7/(PetscSqrtReal(2)*2) * PetscPowReal(betath,3) * n * e * c;
136     Ec = n*lnLambda*PetscPowReal(e,3) / (4*PETSC_PI*PetscPowReal(eps0,2)*m*c*c);
137     *E = Ec;
138     PetscPrintf(PETSC_COMM_WORLD, "CalculateE j0=%g Ec = %g\n",j0,Ec);
139   } else {
140     PetscReal Ed, vth;
141     vth = PetscSqrtReal(8*Tev*e/(m*PETSC_PI));
142     Ed =  n*lnLambda*PetscPowReal(e,3) / (4*PETSC_PI*PetscPowReal(eps0,2)*m*vth*vth);
143     *E = Ed;
144   }
145   PetscFunctionReturn(0);
146 }
147 
148 static PetscReal Spitzer(PetscReal m_e, PetscReal e, PetscReal Z, PetscReal epsilon0,  PetscReal lnLam, PetscReal kTe_joules)
149 {
150   PetscReal Fz = (1+1.198*Z+0.222*Z*Z)/(1+2.966*Z+0.753*Z*Z), eta;
151   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);
152   return eta;
153 }
154 
155 /*  */
156 static PetscErrorCode testNone(TS ts, Vec X, PetscInt stepi, PetscReal time, PetscBool islast, LandauCtx *ctx, REctx *rectx)
157 {
158   PetscFunctionBeginUser;
159   PetscFunctionReturn(0);
160 }
161 
162 /*  */
163 static PetscErrorCode testSpitzer(TS ts, Vec X, PetscInt stepi, PetscReal time, PetscBool islast, LandauCtx *ctx, REctx *rectx)
164 {
165   PetscInt          ii,nDMs;
166   PetscDS           prob;
167   static PetscReal  old_ratio = 1e10;
168   TSConvergedReason reason;
169   PetscReal         J,J_re,spit_eta,Te_kev=0,E,ratio,Z,n_e,v,v2;
170   PetscScalar       user[2] = {0.,ctx->charges[0]}, q[LANDAU_MAX_SPECIES],tt[LANDAU_MAX_SPECIES],vz;
171   PetscReal         dt;
172   DM                pack, plexe = ctx->plex[0], plexi = (ctx->num_grids==1) ? NULL : ctx->plex[1];
173   Vec               *XsubArray;
174 
175   PetscFunctionBeginUser;
176   PetscCheck(ctx->num_species==2,PETSC_COMM_SELF, PETSC_ERR_PLIB, "ctx->num_species %" PetscInt_FMT " != 2",ctx->num_species);
177   PetscCall(VecGetDM(X, &pack));
178   PetscCheck(pack,PETSC_COMM_SELF, PETSC_ERR_PLIB, "no DM");
179   PetscCall(DMCompositeGetNumberDM(pack,&nDMs));
180   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);
181   PetscCall(PetscMalloc(sizeof(*XsubArray)*nDMs, &XsubArray));
182   PetscCall(DMCompositeGetAccessArray(pack, X, nDMs, NULL, XsubArray)); // read only
183   PetscCall(TSGetTimeStep(ts,&dt));
184   /* get current for each grid */
185   for (ii=0;ii<ctx->num_species;ii++) q[ii] = ctx->charges[ii];
186   PetscCall(DMGetDS(plexe, &prob));
187   PetscCall(PetscDSSetConstants(prob, 2, &q[0]));
188   PetscCall(PetscDSSetObjective(prob, 0, &f0_jz_sum));
189   PetscCall(DMPlexComputeIntegralFEM(plexe,XsubArray[ LAND_PACK_IDX(ctx->batch_view_idx,0) ],tt,NULL));
190   J = -ctx->n_0*ctx->v_0*PetscRealPart(tt[0]);
191   if (plexi) { // add first (only) ion
192     PetscCall(DMGetDS(plexi, &prob));
193     PetscCall(PetscDSSetConstants(prob, 1, &q[1]));
194     PetscCall(PetscDSSetObjective(prob, 0, &f0_jz_sum));
195     PetscCall(DMPlexComputeIntegralFEM(plexi,XsubArray[LAND_PACK_IDX(ctx->batch_view_idx,1)],tt,NULL));
196     J += -ctx->n_0*ctx->v_0*PetscRealPart(tt[0]);
197   }
198   /* get N_e */
199   PetscCall(DMGetDS(plexe, &prob));
200   PetscCall(PetscDSSetConstants(prob, 1, user));
201   PetscCall(PetscDSSetObjective(prob, 0, &f0_n));
202   PetscCall(DMPlexComputeIntegralFEM(plexe,XsubArray[LAND_PACK_IDX(ctx->batch_view_idx,0)],tt,NULL));
203   n_e = PetscRealPart(tt[0])*ctx->n_0;
204   /* Z */
205   Z = -ctx->charges[1]/ctx->charges[0];
206   /* remove drift */
207   if (0) {
208     user[0] = 0; // electrons
209     PetscCall(DMGetDS(plexe, &prob));
210     PetscCall(PetscDSSetConstants(prob, 1, user));
211     PetscCall(PetscDSSetObjective(prob, 0, &f0_vz));
212     PetscCall(DMPlexComputeIntegralFEM(plexe,XsubArray[LAND_PACK_IDX(ctx->batch_view_idx,0)],tt,NULL));
213     vz = ctx->n_0*PetscRealPart(tt[0])/n_e; /* non-dimensional */
214   } else vz = 0;
215   /* thermal velocity */
216   PetscCall(DMGetDS(plexe, &prob));
217   PetscCall(PetscDSSetConstants(prob, 1, &vz));
218   PetscCall(PetscDSSetObjective(prob, 0, &f0_ve_shift));
219   PetscCall(DMPlexComputeIntegralFEM(plexe,XsubArray[LAND_PACK_IDX(ctx->batch_view_idx,0)],tt,NULL));
220   v = ctx->n_0*ctx->v_0*PetscRealPart(tt[0])/n_e;   /* remove number density to get velocity */
221   v2 = PetscSqr(v);                                    /* use real space: m^2 / s^2 */
222   Te_kev = (v2*ctx->masses[0]*PETSC_PI/8)*kev_joul;    /* temperature in kev */
223   spit_eta = Spitzer(ctx->masses[0],-ctx->charges[0],Z,ctx->epsilon0,ctx->lnLam,Te_kev/kev_joul); /* kev --> J (kT) */
224   if (0) {
225     PetscCall(DMGetDS(plexe, &prob));
226     PetscCall(PetscDSSetConstants(prob, 1, q));
227     PetscCall(PetscDSSetObjective(prob, 0, &f0_j_re));
228     PetscCall(DMPlexComputeIntegralFEM(plexe,XsubArray[LAND_PACK_IDX(ctx->batch_view_idx,0)],tt,NULL));
229   } else tt[0] = 0;
230   J_re = -ctx->n_0*ctx->v_0*PetscRealPart(tt[0]);
231   PetscCall(DMCompositeRestoreAccessArray(pack, X, nDMs, NULL, XsubArray)); // read only
232   PetscCall(PetscFree(XsubArray));
233 
234   if (rectx->use_spitzer_eta) {
235     E = ctx->Ez = spit_eta*(rectx->j-J_re);
236   } else {
237     E = ctx->Ez; /* keep real E */
238     rectx->j = J; /* cache */
239   }
240 
241   ratio = E/J/spit_eta;
242   if (stepi>10 && !rectx->use_spitzer_eta && (
243         (old_ratio-ratio < 1.e-6))) {
244     rectx->pulse_start = time + 0.98*dt;
245     rectx->use_spitzer_eta = PETSC_TRUE;
246   }
247   PetscCall(TSGetConvergedReason(ts,&reason));
248   PetscCall(TSGetConvergedReason(ts,&reason));
249   if ((rectx->plotting) || stepi == 0 || reason || rectx->pulse_start == time + 0.98*dt) {
250     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,(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));
251     PetscCheck(rectx->pulse_start != (time + 0.98*dt),PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG,"Spitzer complete ratio=%g",(double)ratio);
252   }
253   old_ratio = ratio;
254   PetscFunctionReturn(0);
255 }
256 
257 static const double ppp = 2;
258 static void f0_0_diff_lp(PetscInt dim, PetscInt Nf, PetscInt NfAux,
259                           const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
260                           const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
261                           PetscReal t, const PetscReal x[],  PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0)
262 {
263   LandauCtx       *ctx = (LandauCtx *)constants;
264   REctx           *rectx = (REctx*)ctx->data;
265   PetscInt        ii = rectx->idx, i;
266   const PetscReal kT_m = ctx->k*ctx->thermal_temps[ii]/ctx->masses[ii]; /* kT/m */
267   const PetscReal n = ctx->n[ii];
268   PetscReal       diff, f_maxwell, v2 = 0, theta = 2*kT_m/(ctx->v_0*ctx->v_0); /* theta = 2kT/mc^2 */
269   for (i = 0; i < dim; ++i) v2 += x[i]*x[i];
270   f_maxwell = n*PetscPowReal(PETSC_PI*theta,-1.5)*(PetscExpReal(-v2/theta));
271   diff = 2.*PETSC_PI*x[0]*(PetscRealPart(u[ii]) - f_maxwell);
272   f0[0] = PetscPowReal(diff,ppp);
273 }
274 static void f0_0_maxwellian_lp(PetscInt dim, PetscInt Nf, PetscInt NfAux,
275                           const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
276                           const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
277                           PetscReal t, const PetscReal x[],  PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0)
278 {
279   LandauCtx       *ctx = (LandauCtx *)constants;
280   REctx           *rectx = (REctx*)ctx->data;
281   PetscInt        ii = rectx->idx, i;
282   const PetscReal kT_m = ctx->k*ctx->thermal_temps[ii]/ctx->masses[ii]; /* kT/m */
283   const PetscReal n = ctx->n[ii];
284   PetscReal       f_maxwell, v2 = 0, theta = 2*kT_m/(ctx->v_0*ctx->v_0); /* theta = 2kT/mc^2 */
285   for (i = 0; i < dim; ++i) v2 += x[i]*x[i];
286   f_maxwell = 2.*PETSC_PI*x[0] * n*PetscPowReal(PETSC_PI*theta,-1.5)*(PetscExpReal(-v2/theta));
287   f0[0] = PetscPowReal(f_maxwell,ppp);
288 }
289 
290 /*  */
291 static PetscErrorCode testStable(TS ts, Vec X, PetscInt stepi, PetscReal time, PetscBool islast, LandauCtx *ctx, REctx *rectx)
292 {
293   PetscDS           prob;
294   Vec               X2;
295   PetscReal         ediff,idiff=0,lpm0,lpm1=1;
296   PetscScalar       tt[LANDAU_MAX_SPECIES];
297   DM                dm, plex = ctx->plex[0];
298 
299   PetscFunctionBeginUser;
300   PetscCall(VecGetDM(X, &dm));
301   PetscCall(DMGetDS(plex, &prob));
302   PetscCall(VecDuplicate(X,&X2));
303   PetscCall(VecCopy(X,X2));
304   if (!rectx->X_0) {
305     PetscCall(VecDuplicate(X,&rectx->X_0));
306     PetscCall(VecCopy(X,rectx->X_0));
307   }
308   PetscCall(VecAXPY(X,-1.0,rectx->X_0));
309   PetscCall(PetscDSSetConstants(prob, sizeof(LandauCtx)/sizeof(PetscScalar), (PetscScalar*)ctx));
310   rectx->idx = 0;
311   PetscCall(PetscDSSetObjective(prob, 0, &f0_0_diff_lp));
312   PetscCall(DMPlexComputeIntegralFEM(plex,X2,tt,NULL));
313   ediff = PetscPowReal(PetscRealPart(tt[0]),1./ppp);
314   PetscCall(PetscDSSetObjective(prob, 0, &f0_0_maxwellian_lp));
315   PetscCall(DMPlexComputeIntegralFEM(plex,X2,tt,NULL));
316   lpm0 = PetscPowReal(PetscRealPart(tt[0]),1./ppp);
317   if (ctx->num_species>1) {
318     rectx->idx = 1;
319     PetscCall(PetscDSSetObjective(prob, 0, &f0_0_diff_lp));
320     PetscCall(DMPlexComputeIntegralFEM(plex,X2,tt,NULL));
321     idiff = PetscPowReal(PetscRealPart(tt[0]),1./ppp);
322     PetscCall(PetscDSSetObjective(prob, 0, &f0_0_maxwellian_lp));
323     PetscCall(DMPlexComputeIntegralFEM(plex,X2,tt,NULL));
324     lpm1 = PetscPowReal(PetscRealPart(tt[0]),1./ppp);
325   }
326   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)));
327   /* view */
328   PetscCall(VecCopy(X2,X));
329   PetscCall(VecDestroy(&X2));
330   if (islast) {
331     PetscCall(VecDestroy(&rectx->X_0));
332     rectx->X_0 = NULL;
333   }
334   PetscFunctionReturn(0);
335 }
336 
337 static PetscErrorCode EInduction(Vec X, Vec X_t, PetscInt step, PetscReal time, LandauCtx *ctx, PetscReal *a_E)
338 {
339   REctx             *rectx = (REctx*)ctx->data;
340   PetscInt          ii;
341   DM                dm,plex;
342   PetscScalar       tt[LANDAU_MAX_SPECIES], qv0[LANDAU_MAX_SPECIES];
343   PetscReal         dJ_dt;
344   PetscDS           prob;
345 
346   PetscFunctionBeginUser;
347   for (ii=0;ii<ctx->num_species;ii++) qv0[ii] = ctx->charges[ii]*ctx->v_0;
348   PetscCall(VecGetDM(X, &dm));
349   PetscCall(DMGetDS(dm, &prob));
350   PetscCall(DMConvert(dm, DMPLEX, &plex));
351   /* get d current / dt */
352   PetscCall(PetscDSSetConstants(prob, ctx->num_species, qv0));
353   PetscCall(PetscDSSetObjective(prob, 0, &f0_jz_sum));
354   PetscCheck(X_t,PETSC_COMM_SELF, PETSC_ERR_PLIB, "X_t");
355   PetscCall(DMPlexComputeIntegralFEM(plex,X_t,tt,NULL));
356   dJ_dt = -ctx->n_0*PetscRealPart(tt[0])/ctx->t_0;
357   /* E induction */
358   *a_E = -rectx->L*dJ_dt + rectx->Ez_initial;
359   PetscCall(DMDestroy(&plex));
360   PetscFunctionReturn(0);
361 }
362 
363 static PetscErrorCode EConstant(Vec X,  Vec X_t, PetscInt step, PetscReal time, LandauCtx *ctx, PetscReal *a_E)
364 {
365   PetscFunctionBeginUser;
366   *a_E = ctx->Ez;
367   PetscFunctionReturn(0);
368 }
369 
370 static PetscErrorCode ENone(Vec X,  Vec X_t, PetscInt step, PetscReal time, LandauCtx *ctx, PetscReal *a_E)
371 {
372   PetscFunctionBeginUser;
373   *a_E = 0;
374   PetscFunctionReturn(0);
375 }
376 
377 /* ------------------------------------------------------------------- */
378 /*
379    FormSource - Evaluates source terms F(t).
380 
381    Input Parameters:
382 .  ts - the TS context
383 .  time -
384 .  X_dummmy - input vector
385 .  dummy - optional user-defined context, as set by SNESSetFunction()
386 
387    Output Parameter:
388 .  F - function vector
389  */
390 static PetscErrorCode FormSource(TS ts, PetscReal ftime, Vec X_dummmy, Vec F, void *dummy)
391 {
392   PetscReal      new_imp_rate;
393   LandauCtx      *ctx;
394   DM             pack;
395   REctx          *rectx;
396 
397   PetscFunctionBeginUser;
398   PetscCall(TSGetDM(ts,&pack));
399   PetscCall(DMGetApplicationContext(pack, &ctx));
400   rectx = (REctx*)ctx->data;
401   /* check for impurities */
402   PetscCall(rectx->impuritySrcRate(ftime,&new_imp_rate,ctx));
403   if (new_imp_rate != 0) {
404     if (new_imp_rate != rectx->current_rate) {
405       PetscInt       ii;
406       PetscReal      dne_dt,dni_dt,tilda_ns[LANDAU_MAX_SPECIES],temps[LANDAU_MAX_SPECIES];
407       Vec            globFarray[LANDAU_MAX_GRIDS*LANDAU_MAX_BATCH_SZ];
408       rectx->current_rate = new_imp_rate;
409       for (ii=1;ii<LANDAU_MAX_SPECIES;ii++) tilda_ns[ii] = 0;
410       for (ii=1;ii<LANDAU_MAX_SPECIES;ii++)    temps[ii] = 1;
411       dni_dt = new_imp_rate               /* *ctx->t_0 */; /* fully ionized immediately, no normalize, stay in non-dim */
412       dne_dt = new_imp_rate*rectx->Ne_ion /* *ctx->t_0 */;
413       tilda_ns[0] = dne_dt;        tilda_ns[rectx->imp_idx] = dni_dt;
414       temps[0]    = rectx->T_cold;    temps[rectx->imp_idx] = rectx->T_cold;
415       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));
416       PetscCall(DMCompositeGetAccessArray(pack, F, ctx->num_grids*ctx->batch_sz, NULL, globFarray));
417       for (PetscInt grid=0 ; grid<ctx->num_grids ; grid++) {
418         /* add it */
419         PetscCall(DMPlexLandauAddMaxwellians(ctx->plex[grid],globFarray[ LAND_PACK_IDX(0,grid) ],ftime,temps,tilda_ns,grid,0,1,ctx));
420       }
421       // Does DMCompositeRestoreAccessArray copy the data back? (no)
422       PetscCall(DMCompositeRestoreAccessArray(pack, F, ctx->num_grids*ctx->batch_sz, NULL, globFarray));
423     }
424   } else {
425     PetscCall(VecZeroEntries(F));
426     rectx->current_rate = 0;
427   }
428   PetscFunctionReturn(0);
429 }
430 
431 PetscErrorCode Monitor(TS ts, PetscInt stepi, PetscReal time, Vec X, void *actx)
432 {
433   LandauCtx         *ctx = (LandauCtx*) actx;   /* user-defined application context */
434   REctx             *rectx = (REctx*)ctx->data;
435   DM                pack=NULL;
436   Vec               globXArray[LANDAU_MAX_GRIDS*LANDAU_MAX_BATCH_SZ];
437   TSConvergedReason reason;
438 
439   PetscFunctionBeginUser;
440   PetscCall(TSGetConvergedReason(ts,&reason));
441   if (rectx->grid_view_idx != -1 || (reason && ctx->verbose > 3)) {
442     PetscCall(VecGetDM(X, &pack));
443     PetscCall(DMCompositeGetAccessArray(pack, X, ctx->num_grids*ctx->batch_sz, NULL, globXArray));
444   }
445   if (stepi > rectx->plotStep && rectx->plotting) {
446     rectx->plotting = PETSC_FALSE; /* was doing diagnostics, now done */
447     rectx->plotIdx++;
448   }
449   /* view */
450   if (time/rectx->plotDt >= (PetscReal)rectx->plotIdx || reason) {
451     if ((reason || stepi==0 || rectx->plotIdx%rectx->print_period==0) && ctx->verbose > 1) {
452       /* print norms */
453       PetscCall(DMPlexLandauPrintNorms(X, stepi));
454     }
455     if (!rectx->plotting) { /* first step of possible backtracks */
456       rectx->plotting = PETSC_TRUE;
457       /* diagnostics + change E field with Sptizer (not just a monitor) */
458       PetscCall(rectx->test(ts,X,stepi,time,reason ? PETSC_TRUE : PETSC_FALSE, ctx, rectx));
459     } else {
460       PetscPrintf(PETSC_COMM_WORLD, "\t\t ERROR SKIP test spit ------\n");
461       rectx->plotting = PETSC_TRUE;
462     }
463     if (rectx->grid_view_idx != -1) {
464       PetscCall(PetscObjectSetName((PetscObject) globXArray[ LAND_PACK_IDX(ctx->batch_view_idx,rectx->grid_view_idx) ], rectx->grid_view_idx==0 ? "ue" : "ui"));
465       /* view, overwrite step when back tracked */
466       PetscCall(DMSetOutputSequenceNumber(pack, rectx->plotIdx, time*ctx->t_0));
467       PetscCall(VecViewFromOptions(globXArray[ LAND_PACK_IDX(ctx->batch_view_idx, rectx->grid_view_idx) ],NULL,"-ex2_vec_view"));
468     }
469     rectx->plotStep = stepi;
470   } else {
471     if (rectx->plotting) PetscPrintf(PETSC_COMM_WORLD," ERROR rectx->plotting=%s step %" PetscInt_FMT "\n",PetscBools[rectx->plotting],stepi);
472     /* diagnostics + change E field with Sptizer (not just a monitor) - can we lag this? */
473     PetscCall(rectx->test(ts,X,stepi,time,reason ? PETSC_TRUE : PETSC_FALSE, ctx, rectx));
474   }
475   /* parallel check that only works of all batches are identical */
476   if (reason && ctx->verbose > 3) {
477     PetscReal    val,rval;
478     PetscMPIInt  rank;
479     PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
480     for (PetscInt grid=0;grid<ctx->num_grids;grid++) {
481       PetscInt nerrors=0;
482       for (PetscInt i=0; i<ctx->batch_sz;i++) {
483         PetscCall(VecNorm(globXArray[ LAND_PACK_IDX(i,grid) ],NORM_2,&val));
484         if (i==0) rval = val;
485         else if ((val=PetscAbs(val-rval)/rval) > 1000*PETSC_MACHINE_EPSILON) {
486           PetscCall(PetscPrintf(PETSC_COMM_SELF, " [%d] Warning %" PetscInt_FMT ".%" PetscInt_FMT ") diff = %2.15e\n",rank,grid,i,(double)val));
487           nerrors++;
488         }
489       }
490       if (nerrors) {
491         PetscCall(PetscPrintf(PETSC_COMM_SELF, " ***** [%d] ERROR max %" PetscInt_FMT " errors\n",rank,nerrors));
492       } else {
493         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "[%d] %" PetscInt_FMT ") batch consistency check OK\n",rank,grid));
494       }
495     }
496   }
497   rectx->idx = 0;
498   if (rectx->grid_view_idx != -1 || (reason && ctx->verbose > 3)) {
499     PetscCall(DMCompositeRestoreAccessArray(pack, X, ctx->num_grids*ctx->batch_sz, NULL, globXArray));
500   }
501   PetscFunctionReturn(0);
502 }
503 
504 PetscErrorCode PreStep(TS ts)
505 {
506   LandauCtx      *ctx;
507   REctx          *rectx;
508   DM             dm;
509   PetscInt       stepi;
510   PetscReal      time;
511   Vec            X;
512 
513   PetscFunctionBeginUser;
514   /* not used */
515   PetscCall(TSGetDM(ts,&dm));
516   PetscCall(TSGetTime(ts,&time));
517   PetscCall(TSGetSolution(ts,&X));
518   PetscCall(DMGetApplicationContext(dm, &ctx));
519   rectx = (REctx*)ctx->data;
520   PetscCall(TSGetStepNumber(ts, &stepi));
521   /* update E */
522   PetscCall(rectx->E(X, NULL, stepi, time, ctx, &ctx->Ez));
523   PetscFunctionReturn(0);
524 }
525 
526 /* 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) */
527 static PetscErrorCode stepSrc(PetscReal time, PetscReal *rho, LandauCtx *ctx)
528 {
529   REctx         *rectx = (REctx*)ctx->data;
530 
531   PetscFunctionBeginUser;
532   if (time >= rectx->pulse_start) *rho = rectx->pulse_rate;
533   else *rho = 0.;
534   PetscFunctionReturn(0);
535 }
536 static PetscErrorCode zeroSrc(PetscReal time, PetscReal *rho, LandauCtx *ctx)
537 {
538   PetscFunctionBeginUser;
539   *rho = 0.;
540   PetscFunctionReturn(0);
541 }
542 static PetscErrorCode pulseSrc(PetscReal time, PetscReal *rho, LandauCtx *ctx)
543 {
544   REctx *rectx = (REctx*)ctx->data;
545 
546   PetscFunctionBeginUser;
547   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'");
548   if (time < rectx->pulse_start || time > rectx->pulse_start + 3*rectx->pulse_width) *rho = 0;
549   else {
550     double x = PetscSinReal((time-rectx->pulse_start)/(3*rectx->pulse_width)*2*PETSC_PI - PETSC_PI/2) + 1; /* 0:2, integrates to 1.0 */
551     *rho = rectx->pulse_rate * x / (3*rectx->pulse_width);
552     if (!rectx->use_spitzer_eta) rectx->use_spitzer_eta = PETSC_TRUE; /* use it next time */
553   }
554   PetscFunctionReturn(0);
555 }
556 
557 #undef __FUNCT__
558 #define __FUNCT__ "ProcessREOptions"
559 static PetscErrorCode ProcessREOptions(REctx *rectx, const LandauCtx *ctx, DM dm, const char prefix[])
560 {
561   PetscFunctionList plist = NULL, testlist = NULL, elist = NULL;
562   char              pname[256],testname[256],ename[256];
563   DM                dm_dummy;
564   PetscBool         Connor_E = PETSC_FALSE;
565 
566   PetscFunctionBeginUser;
567   PetscCall(DMCreate(PETSC_COMM_WORLD,&dm_dummy));
568   rectx->Ne_ion = 1;                 /* number of electrons given up by impurity ion */
569   rectx->T_cold = .005;              /* kev */
570   rectx->ion_potential = 15;         /* ev */
571   rectx->L = 2;
572   rectx->X_0 = NULL;
573   rectx->imp_idx = ctx->num_species - 1; /* default ionized impurity as last one */
574   rectx->pulse_start = PETSC_MAX_REAL;
575   rectx->pulse_width = 1;
576   rectx->plotStep = PETSC_MAX_INT;
577   rectx->pulse_rate = 1.e-1;
578   rectx->current_rate = 0;
579   rectx->plotIdx = 0;
580   rectx->j = 0;
581   rectx->plotDt = 1.0;
582   rectx->plotting = PETSC_FALSE;
583   rectx->use_spitzer_eta = PETSC_FALSE;
584   rectx->idx = 0;
585   rectx->print_period = 10;
586   rectx->grid_view_idx = -1; // do not get if not needed
587   /* Register the available impurity sources */
588   PetscCall(PetscFunctionListAdd(&plist,"step",&stepSrc));
589   PetscCall(PetscFunctionListAdd(&plist,"none",&zeroSrc));
590   PetscCall(PetscFunctionListAdd(&plist,"pulse",&pulseSrc));
591   PetscCall(PetscStrcpy(pname,"none"));
592   PetscCall(PetscFunctionListAdd(&testlist,"none",&testNone));
593   PetscCall(PetscFunctionListAdd(&testlist,"spitzer",&testSpitzer));
594   PetscCall(PetscFunctionListAdd(&testlist,"stable",&testStable));
595   PetscCall(PetscStrcpy(testname,"none"));
596   PetscCall(PetscFunctionListAdd(&elist,"none",&ENone));
597   PetscCall(PetscFunctionListAdd(&elist,"induction",&EInduction));
598   PetscCall(PetscFunctionListAdd(&elist,"constant",&EConstant));
599   PetscCall(PetscStrcpy(ename,"constant"));
600 
601   PetscOptionsBegin(PETSC_COMM_SELF, prefix, "Options for Runaway/seed electron model", "none");
602   PetscCall(PetscOptionsReal("-ex2_plot_dt", "Plotting interval", "ex2.c", rectx->plotDt, &rectx->plotDt, NULL));
603   if (rectx->plotDt < 0) rectx->plotDt = 1e30;
604   if (rectx->plotDt == 0) rectx->plotDt = 1e-30;
605   PetscCall(PetscOptionsInt("-ex2_print_period", "Plotting interval", "ex2.c", rectx->print_period, &rectx->print_period, NULL));
606   PetscCall(PetscOptionsInt("-ex2_grid_view_idx", "grid_view_idx", "ex2.c", rectx->grid_view_idx, &rectx->grid_view_idx, NULL));
607   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);
608   PetscCall(PetscOptionsFList("-ex2_impurity_source_type","Name of impurity source to run","",plist,pname,pname,sizeof(pname),NULL));
609   PetscCall(PetscOptionsFList("-ex2_test_type","Name of test to run","",testlist,testname,testname,sizeof(testname),NULL));
610   PetscCall(PetscOptionsInt("-ex2_impurity_index", "index of sink for impurities", "none", rectx->imp_idx, &rectx->imp_idx, NULL));
611   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);
612   PetscCall(PetscOptionsFList("-ex2_e_field_type","Electric field type","",elist,ename,ename,sizeof(ename),NULL));
613   rectx->Ne_ion = -ctx->charges[rectx->imp_idx]/ctx->charges[0];
614   PetscCall(PetscOptionsReal("-ex2_t_cold","Temperature of cold electron and ions after ionization in keV","none",rectx->T_cold,&rectx->T_cold, NULL));
615   PetscCall(PetscOptionsReal("-ex2_pulse_start_time","Time at which pulse happens for 'pulse' source","none",rectx->pulse_start,&rectx->pulse_start, NULL));
616   PetscCall(PetscOptionsReal("-ex2_pulse_width_time","Width of pulse 'pulse' source","none",rectx->pulse_width,&rectx->pulse_width, NULL));
617   PetscCall(PetscOptionsReal("-ex2_pulse_rate","Number density of pulse for 'pulse' source","none",rectx->pulse_rate,&rectx->pulse_rate, NULL));
618   rectx->T_cold *= 1.16e7; /* convert to Kelvin */
619   PetscCall(PetscOptionsReal("-ex2_ion_potential","Potential to ionize impurity (should be array) in ev","none",rectx->ion_potential,&rectx->ion_potential, NULL));
620   PetscCall(PetscOptionsReal("-ex2_inductance","Inductance E feild","none",rectx->L,&rectx->L, NULL));
621   PetscCall(PetscOptionsBool("-ex2_connor_e_field_units","Scale Ex but Connor-Hastie E_c","none",Connor_E,&Connor_E, NULL));
622   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));
623   PetscOptionsEnd();
624   /* get impurity source rate function */
625   PetscCall(PetscFunctionListFind(plist,pname,&rectx->impuritySrcRate));
626   PetscCheck(rectx->impuritySrcRate,PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG,"No impurity source function found '%s'",pname);
627   PetscCall(PetscFunctionListFind(testlist,testname,&rectx->test));
628   PetscCheck(rectx->test,PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG,"No test found '%s'",testname);
629   PetscCall(PetscFunctionListFind(elist,ename,&rectx->E));
630   PetscCheck(rectx->E,PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG,"No E field function found '%s'",ename);
631   PetscCall(PetscFunctionListDestroy(&plist));
632   PetscCall(PetscFunctionListDestroy(&testlist));
633   PetscCall(PetscFunctionListDestroy(&elist));
634 
635   /* convert E from Connor-Hastie E_c units to real if doing Spitzer E */
636   if (Connor_E) {
637     PetscReal E = ctx->Ez, Tev = ctx->thermal_temps[0]*8.621738e-5, n = ctx->n_0*ctx->n[0];
638     CalculateE(Tev, n, ctx->lnLam, ctx->epsilon0, &E);
639     ((LandauCtx*)ctx)->Ez *= E;
640   }
641   PetscCall(DMDestroy(&dm_dummy));
642   PetscFunctionReturn(0);
643 }
644 
645 int main(int argc, char **argv)
646 {
647   DM             pack;
648   Vec            X;
649   PetscInt       dim = 2, nDMs;
650   TS             ts;
651   Mat            J;
652   PetscDS        prob;
653   LandauCtx      *ctx;
654   REctx          *rectx;
655 #if defined PETSC_USE_LOG
656   PetscLogStage  stage;
657 #endif
658   PetscMPIInt    rank;
659 #if defined(PETSC_HAVE_THREADSAFETY)
660   double         starttime, endtime;
661 #endif
662   PetscCall(PetscInitialize(&argc, &argv, NULL,help));
663   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
664   if (rank) { /* turn off output stuff for duplicate runs */
665     PetscCall(PetscOptionsClearValue(NULL,"-ex2_dm_view"));
666     PetscCall(PetscOptionsClearValue(NULL,"-ex2_vec_view"));
667     PetscCall(PetscOptionsClearValue(NULL,"-ex2_vec_view_init"));
668     PetscCall(PetscOptionsClearValue(NULL,"-ex2_dm_view_init"));
669     PetscCall(PetscOptionsClearValue(NULL,"-info")); /* this does not work */
670   }
671   PetscCall(PetscOptionsGetInt(NULL,NULL, "-dim", &dim, NULL));
672   /* Create a mesh */
673   PetscCall(DMPlexLandauCreateVelocitySpace(PETSC_COMM_WORLD, dim, "", &X, &J, &pack));
674   PetscCall(DMCompositeGetNumberDM(pack,&nDMs));
675   PetscCall(PetscObjectSetName((PetscObject)J, "Jacobian"));
676   PetscCall(PetscObjectSetName((PetscObject)X, "f"));
677   PetscCall(DMGetApplicationContext(pack, &ctx));
678   PetscCall(DMSetUp(pack));
679   /* context */
680   PetscCall(PetscNew(&rectx));
681   ctx->data = rectx;
682   PetscCall(ProcessREOptions(rectx,ctx,pack,""));
683   PetscCall(DMGetDS(pack, &prob));
684   if (rectx->grid_view_idx != -1) {
685     Vec *XsubArray=NULL;
686     PetscCall(PetscMalloc(sizeof(*XsubArray)*nDMs, &XsubArray));
687     PetscCall(DMCompositeGetAccessArray(pack, X, nDMs, NULL, XsubArray)); // read only
688     PetscCall(PetscObjectSetName((PetscObject) XsubArray[ LAND_PACK_IDX(ctx->batch_view_idx, rectx->grid_view_idx) ], rectx->grid_view_idx==0 ? "ue" : "ui"));
689     PetscCall(DMViewFromOptions(ctx->plex[rectx->grid_view_idx],NULL,"-ex2_dm_view"));
690     PetscCall(DMViewFromOptions(ctx->plex[rectx->grid_view_idx], NULL,"-ex2_dm_view_init"));
691     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)
692     PetscCall(DMCompositeRestoreAccessArray(pack, X, nDMs, NULL, XsubArray)); // read only
693     PetscCall(PetscFree(XsubArray));
694   }
695   PetscCall(DMSetOutputSequenceNumber(pack, 0, 0.0));
696   /* Create timestepping solver context */
697   PetscCall(TSCreate(PETSC_COMM_SELF,&ts));
698   PetscCall(TSSetDM(ts,pack));
699   PetscCall(TSSetIFunction(ts,NULL,DMPlexLandauIFunction,NULL));
700   PetscCall(TSSetIJacobian(ts,J,J,DMPlexLandauIJacobian,NULL));
701   PetscCall(TSSetRHSFunction(ts,NULL,FormSource,NULL));
702   PetscCall(TSSetFromOptions(ts));
703   PetscCall(TSSetSolution(ts,X));
704   PetscCall(TSSetApplicationContext(ts, ctx));
705   PetscCall(TSMonitorSet(ts,Monitor,ctx,NULL));
706   PetscCall(TSSetPreStep(ts,PreStep));
707   rectx->Ez_initial = ctx->Ez;       /* cache for induction caclulation - applied E field */
708   if (1) { /* warm up an test just DMPlexLandauIJacobian */
709     Vec           vec;
710     PetscInt      nsteps;
711     PetscReal     dt;
712     PetscCall(PetscLogStageRegister("Warmup", &stage));
713     PetscCall(PetscLogStagePush(stage));
714     PetscCall(VecDuplicate(X,&vec));
715     PetscCall(VecCopy(X,vec));
716     PetscCall(TSGetMaxSteps(ts,&nsteps));
717     PetscCall(TSGetTimeStep(ts,&dt));
718     PetscCall(TSSetMaxSteps(ts,1));
719     PetscCall(TSSolve(ts,X));
720     PetscCall(TSSetMaxSteps(ts,nsteps));
721     PetscCall(TSSetStepNumber(ts,0));
722     PetscCall(TSSetTime(ts,0));
723     PetscCall(TSSetTimeStep(ts,dt));
724     rectx->plotIdx = 0;
725     rectx->plotting = PETSC_FALSE;
726     PetscCall(PetscLogStagePop());
727     PetscCall(VecCopy(vec,X));
728     PetscCall(VecDestroy(&vec));
729     PetscCall(PetscObjectStateIncrease((PetscObject)ctx->J));
730   }
731   /* go */
732   PetscCall(PetscLogStageRegister("Solve", &stage));
733   ctx->stage = 0; // lets not use this stage
734 #if defined(PETSC_HAVE_THREADSAFETY)
735   ctx->stage = 1; // not set with thread safety
736 #endif
737   //PetscCall(TSSetSolution(ts,X));
738   PetscCall(PetscLogStagePush(stage));
739 #if defined(PETSC_HAVE_THREADSAFETY)
740   starttime = MPI_Wtime();
741 #endif
742 #if defined(PETSC_HAVE_CUDA_NVTX)
743   nvtxRangePushA("ex2-TSSolve-warm");
744 #endif
745   PetscCall(TSSolve(ts,X));
746 #if defined(PETSC_HAVE_CUDA_NVTX)
747   nvtxRangePop();
748 #endif
749   PetscCall(PetscLogStagePop());
750 #if defined(PETSC_HAVE_THREADSAFETY)
751   endtime = MPI_Wtime();
752   ctx->times[LANDAU_EX2_TSSOLVE] += (endtime - starttime);
753 #endif
754   /* clean up */
755   PetscCall(DMPlexLandauDestroyVelocitySpace(&pack));
756   PetscCall(TSDestroy(&ts));
757   PetscCall(VecDestroy(&X));
758   PetscCall(PetscFree(rectx));
759   PetscCall(PetscFinalize());
760   return 0;
761 }
762 
763 /*TEST
764 
765   testset:
766     requires: p4est !complex double defined(PETSC_USE_DMLANDAU_2D)
767     output_file: output/ex2_0.out
768     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-10 -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 -1 -ts_rtol 1e-3 -ts_dt 1.e-1 -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 -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
769     test:
770       suffix: cpu
771       args: -dm_landau_device_type cpu -ksp_type bicg -pc_type jacobi
772     test:
773       suffix: kokkos
774       requires: kokkos_kernels
775       args: -dm_landau_device_type kokkos -dm_mat_type aijkokkos -dm_vec_type kokkos -ksp_type bicg -pc_type jacobi
776     test:
777       suffix: cuda
778       requires: cuda
779       args: -dm_landau_device_type cuda -dm_mat_type aijcusparse -dm_vec_type cuda -mat_cusparse_use_cpu_solve -ksp_type bicg -pc_type jacobi
780     test:
781       suffix: kokkos_batch
782       requires: kokkos_kernels
783       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
784     test:
785       suffix: kokkos_batch_coo
786       requires: kokkos_kernels
787       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 -dm_landau_coo_assembly
788 
789 TEST*/
790