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