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