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