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