xref: /petsc/src/ts/utils/dmplexlandau/tutorials/ex1.c (revision 6b664d0030d045511a188253fb89991de72e26b7)
1*6b664d00Smarkadams4 static char help[] = "Landau collision operator with amnisotropic thermalization verification test as per Hager et al.\n 'A fully non-linear multi-species Fokker-Planck-Landau collision operator for simulation of fusion plasma', and "
2*6b664d00Smarkadams4                      "published as 'A performance portable, fully implicit Landau collision operator with batched linear solvers' https://arxiv.org/abs/2209.03228\n\n";
3e0eea495SMark 
4e0eea495SMark #include <petscts.h>
5e0eea495SMark #include <petsclandau.h>
6*6b664d00Smarkadams4 #include <petscdmcomposite.h>
7*6b664d00Smarkadams4 #include <petscds.h>
824ded41bSmarkadams4 
924ded41bSmarkadams4 /*
10f53e7263SMark Adams  call back method for DMPlexLandauAccess:
1124ded41bSmarkadams4 
1224ded41bSmarkadams4 Input Parameters:
1324ded41bSmarkadams4  .   dm - a DM for this field
1424ded41bSmarkadams4  -   local_field - the local index in the grid for this field
1524ded41bSmarkadams4  .   grid - the grid index
1624ded41bSmarkadams4  +   b_id - the batch index
1724ded41bSmarkadams4  -   vctx - a user context
1824ded41bSmarkadams4 
1924ded41bSmarkadams4  Input/Output Parameters:
2024ded41bSmarkadams4  +   x - Vector to data to
2124ded41bSmarkadams4 
2224ded41bSmarkadams4  */
23f53e7263SMark Adams PetscErrorCode landau_field_print_access_callback(DM dm, Vec x, PetscInt local_field, PetscInt grid, PetscInt b_id, void *vctx)
24d71ae5a4SJacob Faibussowitsch {
25f53e7263SMark Adams   LandauCtx  *ctx;
26*6b664d00Smarkadams4   PetscScalar val;
27*6b664d00Smarkadams4   PetscInt    species;
2824ded41bSmarkadams4 
2924ded41bSmarkadams4   PetscFunctionBegin;
30f53e7263SMark Adams   PetscCall(DMGetApplicationContext(dm, &ctx));
31*6b664d00Smarkadams4   species = ctx->species_offset[grid] + local_field;
32*6b664d00Smarkadams4   val     = (PetscScalar)(LAND_PACK_IDX(b_id, grid) + (species + 1) * 10);
33*6b664d00Smarkadams4   PetscCall(VecSet(x, val));
34*6b664d00Smarkadams4   PetscCall(PetscInfo(dm, "DMPlexLandauAccess user 'add' method to grid %" PetscInt_FMT ", batch %" PetscInt_FMT " and local field %" PetscInt_FMT " with %" PetscInt_FMT " grids\n", grid, b_id, local_field, ctx->num_grids));
35f53e7263SMark Adams 
36*6b664d00Smarkadams4   PetscFunctionReturn(PETSC_SUCCESS);
3724ded41bSmarkadams4 }
38*6b664d00Smarkadams4 
39*6b664d00Smarkadams4 static const PetscReal alphai   = 1 / 1.3;
40*6b664d00Smarkadams4 static const PetscReal kev_joul = 6.241506479963235e+15; /* 1/1000e */
41*6b664d00Smarkadams4 
42*6b664d00Smarkadams4 // constants: [index of (anisotropic) direction of source, z x[1] shift
43*6b664d00Smarkadams4 /* < v, n_s v_|| > */
44*6b664d00Smarkadams4 static void f0_vz(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0)
45*6b664d00Smarkadams4 {
46*6b664d00Smarkadams4   if (dim == 2) f0[0] = u[0] * (2. * PETSC_PI * x[0]) * x[1]; /* n r v_|| */
47*6b664d00Smarkadams4   else f0[0] = u[0] * x[2];
4824ded41bSmarkadams4 }
49*6b664d00Smarkadams4 /* < v, n (v-shift)^2 > */
50*6b664d00Smarkadams4 static void f0_v2_par_shift(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0)
51*6b664d00Smarkadams4 {
52*6b664d00Smarkadams4   PetscReal vz = PetscRealPart(constants[0]);
53*6b664d00Smarkadams4   if (dim == 2) *f0 = u[0] * (2. * PETSC_PI * x[0]) * (x[1] - vz) * (x[1] - vz); /* n r v^2_par|perp */
54*6b664d00Smarkadams4   else *f0 = u[0] * (x[2] - vz) * (x[2] - vz);
55*6b664d00Smarkadams4 }
56*6b664d00Smarkadams4 static void f0_v2_perp(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0)
57*6b664d00Smarkadams4 {
58*6b664d00Smarkadams4   if (dim == 2) *f0 = u[0] * (2. * PETSC_PI * x[0]) * x[0] * x[0]; /* n r v^2_perp */
59*6b664d00Smarkadams4   else *f0 = u[0] * (x[0] * x[0] + x[1] * x[1]);
60*6b664d00Smarkadams4 }
61*6b664d00Smarkadams4 /* < v, n_e > */
62*6b664d00Smarkadams4 static void f0_n(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0)
63*6b664d00Smarkadams4 {
64*6b664d00Smarkadams4   if (dim == 2) f0[0] = 2. * PETSC_PI * x[0] * u[0];
65*6b664d00Smarkadams4   else f0[0] = u[0];
66*6b664d00Smarkadams4 }
67*6b664d00Smarkadams4 static void f0_v2_shift(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0)
68*6b664d00Smarkadams4 {
69*6b664d00Smarkadams4   PetscReal vz = PetscRealPart(constants[0]);
70*6b664d00Smarkadams4   if (dim == 2) f0[0] = u[0] * (2. * PETSC_PI * x[0]) * (x[0] * x[0] + (x[1] - vz) * (x[1] - vz));
71*6b664d00Smarkadams4   else f0[0] = u[0] * (x[0] * x[0] + x[1] * x[1] + (x[2] - vz) * (x[2] - vz));
72*6b664d00Smarkadams4 }
73*6b664d00Smarkadams4 /* Define a Maxwellian function for testing out the operator. */
74*6b664d00Smarkadams4 typedef struct {
75*6b664d00Smarkadams4   PetscReal v_0;
76*6b664d00Smarkadams4   PetscReal kT_m;
77*6b664d00Smarkadams4   PetscReal n;
78*6b664d00Smarkadams4   PetscReal shift;
79*6b664d00Smarkadams4   PetscInt  species;
80*6b664d00Smarkadams4 } MaxwellianCtx;
81*6b664d00Smarkadams4 
82*6b664d00Smarkadams4 static PetscErrorCode maxwellian(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf_dummy, PetscScalar *u, void *actx)
83*6b664d00Smarkadams4 {
84*6b664d00Smarkadams4   MaxwellianCtx *mctx  = (MaxwellianCtx *)actx;
85*6b664d00Smarkadams4   PetscReal      theta = 2 * mctx->kT_m / (mctx->v_0 * mctx->v_0); /* theta = 2kT/mc^2 */
86*6b664d00Smarkadams4   PetscFunctionBegin;
87*6b664d00Smarkadams4   /* evaluate the shifted Maxwellian */
88*6b664d00Smarkadams4   if (dim == 2) u[0] += alphai * mctx->n * PetscPowReal(PETSC_PI * theta, -1.5) * PetscExpReal(-(alphai * x[0] * x[0] + (x[1] - mctx->shift) * (x[1] - mctx->shift)) / theta);
89*6b664d00Smarkadams4   else u[0] += alphai * mctx->n * PetscPowReal(PETSC_PI * theta, -1.5) * PetscExpReal(-(alphai * (x[0] * x[0] + x[1] * x[1]) + (x[2] - mctx->shift) * (x[2] - mctx->shift)) / theta);
90*6b664d00Smarkadams4   PetscFunctionReturn(PETSC_SUCCESS);
91*6b664d00Smarkadams4 }
92*6b664d00Smarkadams4 
93*6b664d00Smarkadams4 static PetscErrorCode SetMaxwellians(DM dm, Vec X, PetscReal time, PetscReal temps[], PetscReal ns[], PetscInt grid, PetscReal shifts[], LandauCtx *ctx)
94*6b664d00Smarkadams4 {
95*6b664d00Smarkadams4   PetscErrorCode (*initu[LANDAU_MAX_SPECIES])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *);
96*6b664d00Smarkadams4   MaxwellianCtx *mctxs[LANDAU_MAX_SPECIES], data[LANDAU_MAX_SPECIES];
97*6b664d00Smarkadams4   PetscFunctionBegin;
98*6b664d00Smarkadams4   if (!ctx) PetscCall(DMGetApplicationContext(dm, &ctx));
99*6b664d00Smarkadams4   for (PetscInt ii = ctx->species_offset[grid], i0 = 0; ii < ctx->species_offset[grid + 1]; ii++, i0++) {
100*6b664d00Smarkadams4     mctxs[i0]        = &data[i0];
101*6b664d00Smarkadams4     data[i0].v_0     = ctx->v_0;                             // v_0 same for all grids
102*6b664d00Smarkadams4     data[i0].kT_m    = ctx->k * temps[ii] / ctx->masses[ii]; /* kT/m = v_th ^ 2*/
103*6b664d00Smarkadams4     data[i0].n       = ns[ii];
104*6b664d00Smarkadams4     initu[i0]        = maxwellian;
105*6b664d00Smarkadams4     data[i0].shift   = 0;
106*6b664d00Smarkadams4     data[i0].species = ii;
107*6b664d00Smarkadams4   }
108*6b664d00Smarkadams4   if (1) {
109*6b664d00Smarkadams4     data[0].shift = -((PetscReal)PetscSign(ctx->charges[ctx->species_offset[grid]])) * ctx->electronShift * ctx->m_0 / ctx->masses[ctx->species_offset[grid]];
110*6b664d00Smarkadams4   } else {
111*6b664d00Smarkadams4     shifts[0]     = 0.5 * PetscSqrtReal(ctx->masses[0] / ctx->masses[1]);
112*6b664d00Smarkadams4     shifts[1]     = 50 * (ctx->masses[0] / ctx->masses[1]);
113*6b664d00Smarkadams4     data[0].shift = ctx->electronShift * shifts[grid] * PetscSqrtReal(data[0].kT_m) / ctx->v_0; // shifts to not matter!!!!
114*6b664d00Smarkadams4   }
115*6b664d00Smarkadams4   PetscCall(DMProjectFunction(dm, time, initu, (void **)mctxs, INSERT_ALL_VALUES, X));
116*6b664d00Smarkadams4   PetscFunctionReturn(PETSC_SUCCESS);
117*6b664d00Smarkadams4 }
118*6b664d00Smarkadams4 
119*6b664d00Smarkadams4 typedef enum {
120*6b664d00Smarkadams4   E_PAR_IDX,
121*6b664d00Smarkadams4   E_PERP_IDX,
122*6b664d00Smarkadams4   I_PAR_IDX,
123*6b664d00Smarkadams4   I_PERP_IDX,
124*6b664d00Smarkadams4   NUM_TEMPS
125*6b664d00Smarkadams4 } TemperatureIDX;
126*6b664d00Smarkadams4 
127*6b664d00Smarkadams4 /* --------------------  Evaluate Function F(x) --------------------- */
128*6b664d00Smarkadams4 static PetscReal n_cm3[2] = {0, 0};
129*6b664d00Smarkadams4 PetscErrorCode   FormFunction(TS ts, PetscReal tdummy, Vec X, Vec F, void *ptr)
130*6b664d00Smarkadams4 {
131*6b664d00Smarkadams4   LandauCtx         *ctx = (LandauCtx *)ptr; /* user-defined application context */
132*6b664d00Smarkadams4   PetscScalar       *f;
133*6b664d00Smarkadams4   const PetscScalar *x;
134*6b664d00Smarkadams4   const PetscReal    k_B = 1.6e-12, e_cgs = 4.8e-10, m_cgs[2] = {9.1094e-28, 9.1094e-28 * ctx->masses[1] / ctx->masses[0]}; // erg/eV, e, m as per NRL;
135*6b664d00Smarkadams4   PetscReal          AA, sqrtA, v_abT, vTe, t1, TeDiff, Te, Ti, Tdiff;
136*6b664d00Smarkadams4 
137*6b664d00Smarkadams4   PetscFunctionBeginUser;
138*6b664d00Smarkadams4   PetscCall(VecGetArrayRead(X, &x));
139*6b664d00Smarkadams4   Te = PetscRealPart(2 * x[E_PERP_IDX] + x[E_PAR_IDX]) / 3, Ti = PetscRealPart(2 * x[I_PERP_IDX] + x[I_PAR_IDX]) / 3;
140*6b664d00Smarkadams4   v_abT = 1.8e-19 * PetscSqrtReal(m_cgs[0] * m_cgs[1]) * n_cm3[0] * ctx->lnLam * PetscPowReal(m_cgs[0] * Ti + m_cgs[1] * Te, -1.5);
141*6b664d00Smarkadams4   PetscCall(VecGetArray(F, &f));
142*6b664d00Smarkadams4   for (PetscInt ii = 0; ii < 2; ii++) {
143*6b664d00Smarkadams4     PetscReal tPerp = PetscRealPart(x[2 * ii + E_PERP_IDX]), tPar = PetscRealPart(x[2 * ii + E_PAR_IDX]);
144*6b664d00Smarkadams4     TeDiff = tPerp - tPar;
145*6b664d00Smarkadams4     AA     = tPerp / tPar - 1;
146*6b664d00Smarkadams4     if (AA < 1e-6) t1 = 0;
147*6b664d00Smarkadams4     else {
148*6b664d00Smarkadams4       sqrtA = PetscSqrtReal(AA);
149*6b664d00Smarkadams4       t1    = (-3 + (AA + 3) * PetscAtanReal(sqrtA) / sqrtA) / PetscSqr(AA);
150*6b664d00Smarkadams4       //PetscReal vTeB = 8.2e-7 * n_cm3[0] * ctx->lnLam * PetscPowReal(Te, -1.5);
151*6b664d00Smarkadams4       vTe = PetscRealPart(2 * PetscSqrtReal(PETSC_PI / m_cgs[ii]) * PetscSqr(PetscSqr(e_cgs)) * n_cm3[0] * ctx->lnLam * PetscPowReal(k_B * x[E_PAR_IDX], -1.5)) * t1;
152*6b664d00Smarkadams4       t1  = vTe * TeDiff * PetscSqrtReal(PETSC_PI); // scaling form NRL that makes it work ???
153*6b664d00Smarkadams4     }
154*6b664d00Smarkadams4     f[2 * ii + E_PAR_IDX]  = 2 * t1; // par
155*6b664d00Smarkadams4     f[2 * ii + E_PERP_IDX] = -t1;    // perp
156*6b664d00Smarkadams4     Tdiff                  = (ii == 0) ? (Ti - Te) : (Te - Ti);
157*6b664d00Smarkadams4     f[2 * ii + E_PAR_IDX] += v_abT * Tdiff;
158*6b664d00Smarkadams4     f[2 * ii + E_PERP_IDX] += v_abT * Tdiff;
159*6b664d00Smarkadams4   }
160*6b664d00Smarkadams4   PetscCall(VecRestoreArrayRead(X, &x));
161*6b664d00Smarkadams4   PetscCall(VecRestoreArray(F, &f));
162*6b664d00Smarkadams4   PetscFunctionReturn(PETSC_SUCCESS);
163*6b664d00Smarkadams4 }
164*6b664d00Smarkadams4 
165*6b664d00Smarkadams4 /* --------------------  Form initial approximation ----------------- */
166*6b664d00Smarkadams4 static PetscReal T0[4] = {300, 390, 200, 260};
167*6b664d00Smarkadams4 PetscErrorCode   createVec_NRL(LandauCtx *ctx, Vec *vec)
168*6b664d00Smarkadams4 {
169*6b664d00Smarkadams4   PetscScalar *x;
170*6b664d00Smarkadams4   Vec          Temps;
171*6b664d00Smarkadams4 
172*6b664d00Smarkadams4   PetscFunctionBeginUser;
173*6b664d00Smarkadams4   PetscCall(VecCreateSeq(PETSC_COMM_SELF, NUM_TEMPS, &Temps));
174*6b664d00Smarkadams4   PetscCall(VecGetArray(Temps, &x));
175*6b664d00Smarkadams4   for (PetscInt i = 0; i < NUM_TEMPS; i++) x[i] = T0[i];
176*6b664d00Smarkadams4   PetscCall(VecRestoreArray(Temps, &x));
177*6b664d00Smarkadams4   *vec = Temps;
178*6b664d00Smarkadams4   PetscFunctionReturn(PETSC_SUCCESS);
179*6b664d00Smarkadams4 }
180*6b664d00Smarkadams4 
181*6b664d00Smarkadams4 PetscErrorCode createTS_NRL(LandauCtx *ctx, Vec Temps)
182*6b664d00Smarkadams4 {
183*6b664d00Smarkadams4   TSAdapt adapt;
184*6b664d00Smarkadams4   TS      ts;
185*6b664d00Smarkadams4 
186*6b664d00Smarkadams4   PetscFunctionBeginUser;
187*6b664d00Smarkadams4   PetscCall(TSCreate(PETSC_COMM_SELF, &ts));
188*6b664d00Smarkadams4   ctx->data = (void *)ts; // 'data' is for applications (eg, monitors)
189*6b664d00Smarkadams4   PetscCall(TSSetApplicationContext(ts, ctx));
190*6b664d00Smarkadams4   PetscCall(TSSetType(ts, TSRK));
191*6b664d00Smarkadams4   PetscCall(TSSetRHSFunction(ts, NULL, FormFunction, ctx));
192*6b664d00Smarkadams4   PetscCall(TSSetSolution(ts, Temps));
193*6b664d00Smarkadams4   PetscCall(TSRKSetType(ts, TSRK2A));
194*6b664d00Smarkadams4   PetscCall(TSSetOptionsPrefix(ts, "nrl_"));
195*6b664d00Smarkadams4   PetscCall(TSSetFromOptions(ts));
196*6b664d00Smarkadams4   PetscCall(TSGetAdapt(ts, &adapt));
197*6b664d00Smarkadams4   PetscCall(TSAdaptSetType(adapt, TSADAPTNONE));
198*6b664d00Smarkadams4   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
199*6b664d00Smarkadams4   PetscCall(TSSetStepNumber(ts, 0));
200*6b664d00Smarkadams4   PetscCall(TSSetMaxSteps(ts, 1));
201*6b664d00Smarkadams4   PetscCall(TSSetTime(ts, 0));
202*6b664d00Smarkadams4 
203*6b664d00Smarkadams4   PetscFunctionReturn(PETSC_SUCCESS);
204*6b664d00Smarkadams4 }
205*6b664d00Smarkadams4 
206*6b664d00Smarkadams4 PetscErrorCode Monitor_nrl(TS ts, PetscInt stepi, PetscReal time, Vec X, void *actx)
207*6b664d00Smarkadams4 {
208*6b664d00Smarkadams4   const PetscScalar *x;
209*6b664d00Smarkadams4   LandauCtx         *ctx = (LandauCtx *)actx; /* user-defined application context */
210*6b664d00Smarkadams4   PetscInt           period, logT;
211*6b664d00Smarkadams4   PetscReal          dt;
212*6b664d00Smarkadams4 
213*6b664d00Smarkadams4   PetscFunctionBeginUser;
214*6b664d00Smarkadams4   PetscCall(TSGetTimeStep(ts, &dt));
215*6b664d00Smarkadams4   logT = (PetscInt)PetscLog2Real(time / dt);
216*6b664d00Smarkadams4   if (logT < 0) logT = 0;
217*6b664d00Smarkadams4   period = PetscPowInt(2, logT) / 2;
218*6b664d00Smarkadams4   if (period == 0) period = 1;
219*6b664d00Smarkadams4   if (stepi % period == 0) {
220*6b664d00Smarkadams4     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "nrl-step %d time= %g ", (int)stepi, (double)(time / ctx->t_0)));
221*6b664d00Smarkadams4     PetscCall(VecGetArrayRead(X, &x));
222*6b664d00Smarkadams4     for (PetscInt i = 0; i < NUM_TEMPS; i++) { PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%g ", (double)PetscRealPart(x[i]))); }
223*6b664d00Smarkadams4     PetscCall(VecRestoreArrayRead(X, &x));
224*6b664d00Smarkadams4     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n"));
225*6b664d00Smarkadams4   }
226*6b664d00Smarkadams4   PetscFunctionReturn(PETSC_SUCCESS);
227*6b664d00Smarkadams4 }
228*6b664d00Smarkadams4 
229*6b664d00Smarkadams4 PetscErrorCode Monitor(TS ts, PetscInt stepi, PetscReal time, Vec X, void *actx)
230*6b664d00Smarkadams4 {
231*6b664d00Smarkadams4   LandauCtx *ctx      = (LandauCtx *)actx; /* user-defined application context */
232*6b664d00Smarkadams4   TS         ts_nrl   = (TS)ctx->data;
233*6b664d00Smarkadams4   PetscInt   printing = 0, logT;
234*6b664d00Smarkadams4 
235*6b664d00Smarkadams4   PetscFunctionBeginUser;
236*6b664d00Smarkadams4   if (ctx->verbose > 0) { // hacks to generate sparse data (eg, use '-dm_landau_verbose 1' and '-dm_landau_verbose -1' to get all steps printed)
237*6b664d00Smarkadams4     PetscReal dt;
238*6b664d00Smarkadams4     PetscCall(TSGetTimeStep(ts, &dt));
239*6b664d00Smarkadams4     logT = (PetscInt)PetscLog2Real(time / dt);
240*6b664d00Smarkadams4     if (logT < 0) logT = 0;
241*6b664d00Smarkadams4     ctx->verbose = PetscPowInt(2, logT) / 2;
242*6b664d00Smarkadams4     if (ctx->verbose == 0) ctx->verbose = 1;
243*6b664d00Smarkadams4   }
244*6b664d00Smarkadams4   if (ctx->verbose) {
245*6b664d00Smarkadams4     TSConvergedReason reason;
246*6b664d00Smarkadams4     PetscCall(TSGetConvergedReason(ts, &reason));
247*6b664d00Smarkadams4     if (stepi % ctx->verbose == 0 || reason || stepi == 1 || ctx->verbose < 0) {
248*6b664d00Smarkadams4       PetscInt nDMs, id;
249*6b664d00Smarkadams4       DM       pack;
250*6b664d00Smarkadams4       Vec     *XsubArray = NULL;
251*6b664d00Smarkadams4       printing           = 1;
252*6b664d00Smarkadams4       PetscCall(TSGetDM(ts, &pack));
253*6b664d00Smarkadams4       PetscCall(DMCompositeGetNumberDM(pack, &nDMs));
254*6b664d00Smarkadams4       PetscCall(DMGetOutputSequenceNumber(ctx->plex[0], &id, NULL));
255*6b664d00Smarkadams4       PetscCall(DMSetOutputSequenceNumber(ctx->plex[0], id + 1, time));
256*6b664d00Smarkadams4       PetscCall(DMSetOutputSequenceNumber(ctx->plex[1], id + 1, time));
257*6b664d00Smarkadams4       PetscCall(PetscInfo(pack, "ex1 plot step %" PetscInt_FMT ", time = %g\n", id, (double)time));
258*6b664d00Smarkadams4       PetscCall(PetscMalloc(sizeof(*XsubArray) * nDMs, &XsubArray));
259*6b664d00Smarkadams4       PetscCall(DMCompositeGetAccessArray(pack, X, nDMs, NULL, XsubArray)); // read only
260*6b664d00Smarkadams4       PetscCall(VecViewFromOptions(XsubArray[LAND_PACK_IDX(ctx->batch_view_idx, 0)], NULL, "-ex1_vec_view_e"));
261*6b664d00Smarkadams4       PetscCall(VecViewFromOptions(XsubArray[LAND_PACK_IDX(ctx->batch_view_idx, 1)], NULL, "-ex1_vec_view_i"));
262*6b664d00Smarkadams4       // temps
263*6b664d00Smarkadams4       for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
264*6b664d00Smarkadams4         PetscDS     prob;
265*6b664d00Smarkadams4         DM          dm      = ctx->plex[grid];
266*6b664d00Smarkadams4         PetscScalar user[2] = {0, 0}, tt[1];
267*6b664d00Smarkadams4         PetscReal   vz_0 = 0, n, energy, e_perp, e_par, m_s = ctx->masses[ctx->species_offset[grid]];
268*6b664d00Smarkadams4         Vec         Xloc = XsubArray[LAND_PACK_IDX(ctx->batch_view_idx, grid)];
269*6b664d00Smarkadams4         PetscCall(DMGetDS(dm, &prob));
270*6b664d00Smarkadams4         /* get n */
271*6b664d00Smarkadams4         PetscCall(PetscDSSetObjective(prob, 0, &f0_n));
272*6b664d00Smarkadams4         PetscCall(DMPlexComputeIntegralFEM(dm, Xloc, tt, NULL));
273*6b664d00Smarkadams4         n = PetscRealPart(tt[0]);
274*6b664d00Smarkadams4         /* get vz */
275*6b664d00Smarkadams4         PetscCall(PetscDSSetObjective(prob, 0, &f0_vz));
276*6b664d00Smarkadams4         PetscCall(DMPlexComputeIntegralFEM(dm, Xloc, tt, NULL));
277*6b664d00Smarkadams4         user[0] = vz_0 = PetscRealPart(tt[0]) / n;
278*6b664d00Smarkadams4         /* energy temp */
279*6b664d00Smarkadams4         PetscCall(PetscDSSetConstants(prob, 2, user));
280*6b664d00Smarkadams4         PetscCall(PetscDSSetObjective(prob, 0, &f0_v2_shift));
281*6b664d00Smarkadams4         PetscCall(DMPlexComputeIntegralFEM(dm, Xloc, tt, ctx));
282*6b664d00Smarkadams4         energy = PetscRealPart(tt[0]) * ctx->v_0 * ctx->v_0 * m_s / n / 3; // scale?
283*6b664d00Smarkadams4         energy *= kev_joul * 1000;                                         // T eV
284*6b664d00Smarkadams4         /* energy temp - par */
285*6b664d00Smarkadams4         PetscCall(PetscDSSetConstants(prob, 2, user));
286*6b664d00Smarkadams4         PetscCall(PetscDSSetObjective(prob, 0, &f0_v2_par_shift));
287*6b664d00Smarkadams4         PetscCall(DMPlexComputeIntegralFEM(dm, Xloc, tt, ctx));
288*6b664d00Smarkadams4         e_par = PetscRealPart(tt[0]) * ctx->v_0 * ctx->v_0 * m_s / n;
289*6b664d00Smarkadams4         e_par *= kev_joul * 1000; // eV
290*6b664d00Smarkadams4         /* energy temp - perp */
291*6b664d00Smarkadams4         PetscCall(PetscDSSetConstants(prob, 2, user));
292*6b664d00Smarkadams4         PetscCall(PetscDSSetObjective(prob, 0, &f0_v2_perp));
293*6b664d00Smarkadams4         PetscCall(DMPlexComputeIntegralFEM(dm, Xloc, tt, ctx));
294*6b664d00Smarkadams4         e_perp = PetscRealPart(tt[0]) * ctx->v_0 * ctx->v_0 * m_s / n / 2;
295*6b664d00Smarkadams4         e_perp *= kev_joul * 1000; // eV
296*6b664d00Smarkadams4         if (grid == 0) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "step %4d) time= %e temperature (eV): ", (int)stepi, (double)time));
297*6b664d00Smarkadams4         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s T= %9.4e T_par= %9.4e T_perp= %9.4e ", (grid == 0) ? "electron:" : ";ion:", (double)energy, (double)e_par, (double)e_perp));
298*6b664d00Smarkadams4         if (n_cm3[grid] == 0) n_cm3[grid] = ctx->n_0 * n * 1e-6; // does not change m^3 --> cm^3
299*6b664d00Smarkadams4       }
300*6b664d00Smarkadams4       // cleanup
301*6b664d00Smarkadams4       PetscCall(DMCompositeRestoreAccessArray(pack, X, nDMs, NULL, XsubArray));
302*6b664d00Smarkadams4       PetscCall(PetscFree(XsubArray));
303*6b664d00Smarkadams4     }
304*6b664d00Smarkadams4   }
305*6b664d00Smarkadams4   /* evolve NRL data, end line */
306*6b664d00Smarkadams4   if (n_cm3[NUM_TEMPS / 2 - 1] < 0 && ts_nrl) {
307*6b664d00Smarkadams4     PetscCall(TSDestroy(&ts_nrl));
308*6b664d00Smarkadams4     ctx->data = NULL;
309*6b664d00Smarkadams4     if (printing) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nSTOP printing NRL Ts\n"));
310*6b664d00Smarkadams4   } else if (ts_nrl) {
311*6b664d00Smarkadams4     const PetscScalar *x;
312*6b664d00Smarkadams4     PetscReal          dt_real, dt;
313*6b664d00Smarkadams4     Vec                U;
314*6b664d00Smarkadams4     PetscCall(TSGetTimeStep(ts, &dt)); // dt for NEXT time step
315*6b664d00Smarkadams4     dt_real = dt * ctx->t_0;
316*6b664d00Smarkadams4     PetscCall(TSGetSolution(ts_nrl, &U));
317*6b664d00Smarkadams4     if (printing) {
318*6b664d00Smarkadams4       PetscCall(VecGetArrayRead(U, &x));
319*6b664d00Smarkadams4       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "NRL_i_par= %9.4e NRL_i_perp= %9.4e ", (double)PetscRealPart(x[I_PAR_IDX]), (double)PetscRealPart(x[I_PERP_IDX])));
320*6b664d00Smarkadams4       if (n_cm3[0] > 0) {
321*6b664d00Smarkadams4         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "NRL_e_par= %9.4e NRL_e_perp= %9.4e\n", (double)PetscRealPart(x[E_PAR_IDX]), (double)PetscRealPart(x[E_PERP_IDX])));
322*6b664d00Smarkadams4       } else PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n"));
323*6b664d00Smarkadams4       PetscCall(VecRestoreArrayRead(U, &x));
324*6b664d00Smarkadams4     }
325*6b664d00Smarkadams4     // we have the next time step, so need to advance now
326*6b664d00Smarkadams4     PetscCall(TSSetTimeStep(ts_nrl, dt_real));
327*6b664d00Smarkadams4     PetscCall(TSSetMaxSteps(ts_nrl, stepi + 1)); // next step
328*6b664d00Smarkadams4     PetscCall(TSSolve(ts_nrl, NULL));
329*6b664d00Smarkadams4   } else if (printing) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n"));
330*6b664d00Smarkadams4   if (printing) { PetscCall(DMPlexLandauPrintNorms(X, stepi /*id + 1*/)); }
331*6b664d00Smarkadams4 
3323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
33324ded41bSmarkadams4 }
33424ded41bSmarkadams4 
335d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
336d71ae5a4SJacob Faibussowitsch {
33724ded41bSmarkadams4   DM          pack;
338*6b664d00Smarkadams4   Vec         X;
339*6b664d00Smarkadams4   PetscInt    dim = 2, nDMs;
340*6b664d00Smarkadams4   TS          ts, ts_nrl = NULL;
341e0eea495SMark   Mat         J;
342*6b664d00Smarkadams4   Vec        *XsubArray = NULL;
343*6b664d00Smarkadams4   LandauCtx  *ctx;
344*6b664d00Smarkadams4   PetscMPIInt rank;
345*6b664d00Smarkadams4   PetscBool   use_nrl   = PETSC_FALSE;
346*6b664d00Smarkadams4   PetscBool   print_nrl = PETSC_TRUE;
347*6b664d00Smarkadams4   PetscReal   dt0;
348327415f7SBarry Smith   PetscFunctionBeginUser;
3499566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
350*6b664d00Smarkadams4   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
351*6b664d00Smarkadams4   if (rank) { /* turn off output stuff for duplicate runs */
352*6b664d00Smarkadams4     PetscCall(PetscOptionsClearValue(NULL, "-ex1_dm_view_e"));
353*6b664d00Smarkadams4     PetscCall(PetscOptionsClearValue(NULL, "-ex1_dm_view_i"));
354*6b664d00Smarkadams4     PetscCall(PetscOptionsClearValue(NULL, "-ex1_vec_view_e"));
355*6b664d00Smarkadams4     PetscCall(PetscOptionsClearValue(NULL, "-ex1_vec_view_i"));
356*6b664d00Smarkadams4     PetscCall(PetscOptionsClearValue(NULL, "-info"));
357*6b664d00Smarkadams4     PetscCall(PetscOptionsClearValue(NULL, "-snes_converged_reason"));
358*6b664d00Smarkadams4     PetscCall(PetscOptionsClearValue(NULL, "-pc_bjkokkos_ksp_converged_reason"));
359*6b664d00Smarkadams4     PetscCall(PetscOptionsClearValue(NULL, "-ksp_converged_reason"));
360*6b664d00Smarkadams4     PetscCall(PetscOptionsClearValue(NULL, "-ts_adapt_monitor"));
361*6b664d00Smarkadams4     PetscCall(PetscOptionsClearValue(NULL, "-ts_monitor"));
362*6b664d00Smarkadams4     PetscCall(PetscOptionsClearValue(NULL, "-snes_monitor"));
363*6b664d00Smarkadams4   }
3649566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-dim", &dim, NULL));
365*6b664d00Smarkadams4   PetscCall(PetscOptionsGetBool(NULL, NULL, "-use_nrl", &use_nrl, NULL));
366*6b664d00Smarkadams4   PetscCall(PetscOptionsGetBool(NULL, NULL, "-print_nrl", &print_nrl, NULL));
367e0eea495SMark   /* Create a mesh */
36824ded41bSmarkadams4   PetscCall(DMPlexLandauCreateVelocitySpace(PETSC_COMM_SELF, dim, "", &X, &J, &pack));
36924ded41bSmarkadams4   PetscCall(DMSetUp(pack));
370*6b664d00Smarkadams4   PetscCall(DMGetApplicationContext(pack, &ctx));
371*6b664d00Smarkadams4   PetscCheck(ctx->num_grids == 2, PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Must have two grids: use '-dm_landau_num_species_grid 1,1'");
372*6b664d00Smarkadams4   PetscCheck(ctx->num_species == 2, PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Must have two species: use '-dm_landau_num_species_grid 1,1'");
373*6b664d00Smarkadams4   PetscCall(DMCompositeGetNumberDM(pack, &nDMs));
374*6b664d00Smarkadams4   /* output plot names */
375*6b664d00Smarkadams4   PetscCall(PetscMalloc(sizeof(*XsubArray) * nDMs, &XsubArray));
376*6b664d00Smarkadams4   PetscCall(DMCompositeGetAccessArray(pack, X, nDMs, NULL, XsubArray)); // read only
377*6b664d00Smarkadams4   PetscCall(PetscObjectSetName((PetscObject)XsubArray[LAND_PACK_IDX(ctx->batch_view_idx, 0)], 0 == 0 ? "ue" : "ui"));
378*6b664d00Smarkadams4   PetscCall(PetscObjectSetName((PetscObject)XsubArray[LAND_PACK_IDX(ctx->batch_view_idx, 1)], 1 == 0 ? "ue" : "ui"));
379*6b664d00Smarkadams4   /* add bimaxwellian anisotropic test */
380*6b664d00Smarkadams4   for (PetscInt b_id = 0; b_id < ctx->batch_sz; b_id++) {
381*6b664d00Smarkadams4     for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
382*6b664d00Smarkadams4       PetscReal shifts[2];
383*6b664d00Smarkadams4       PetscCall(SetMaxwellians(ctx->plex[grid], XsubArray[LAND_PACK_IDX(b_id, grid)], 0.0, ctx->thermal_temps, ctx->n, grid, shifts, ctx));
384*6b664d00Smarkadams4     }
385*6b664d00Smarkadams4   }
386*6b664d00Smarkadams4   PetscCall(DMCompositeRestoreAccessArray(pack, X, nDMs, NULL, XsubArray));
387*6b664d00Smarkadams4   PetscCall(PetscFree(XsubArray));
388*6b664d00Smarkadams4   /* plot */
389*6b664d00Smarkadams4   PetscCall(DMSetOutputSequenceNumber(ctx->plex[0], -1, 0.0));
390*6b664d00Smarkadams4   PetscCall(DMSetOutputSequenceNumber(ctx->plex[1], -1, 0.0));
391*6b664d00Smarkadams4   PetscCall(DMViewFromOptions(ctx->plex[0], NULL, "-ex1_dm_view_e"));
392*6b664d00Smarkadams4   PetscCall(DMViewFromOptions(ctx->plex[1], NULL, "-ex1_dm_view_i"));
393e0eea495SMark   /* Create timestepping solver context */
3949566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_SELF, &ts));
39524ded41bSmarkadams4   PetscCall(TSSetDM(ts, pack));
3969566063dSJacob Faibussowitsch   PetscCall(TSSetIFunction(ts, NULL, DMPlexLandauIFunction, NULL));
3979566063dSJacob Faibussowitsch   PetscCall(TSSetIJacobian(ts, J, J, DMPlexLandauIJacobian, NULL));
3989566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
3999566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
4009566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts, X));
401*6b664d00Smarkadams4   PetscCall(TSMonitorSet(ts, Monitor, ctx, NULL));
402*6b664d00Smarkadams4   /* Create NRL timestepping */
403*6b664d00Smarkadams4   if (use_nrl || print_nrl) {
404*6b664d00Smarkadams4     Vec NRL_vec;
405*6b664d00Smarkadams4     PetscCall(createVec_NRL(ctx, &NRL_vec));
406*6b664d00Smarkadams4     PetscCall(createTS_NRL(ctx, NRL_vec));
407*6b664d00Smarkadams4     PetscCall(VecDestroy(&NRL_vec));
408*6b664d00Smarkadams4   } else ctx->data = NULL;
409*6b664d00Smarkadams4   /* solve */
410*6b664d00Smarkadams4   PetscCall(TSGetTimeStep(ts, &dt0));
411*6b664d00Smarkadams4   PetscCall(TSSetTime(ts, dt0 / 2));
4129566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, X));
413*6b664d00Smarkadams4   /* test add field method & output */
414*6b664d00Smarkadams4   PetscCall(DMPlexLandauAccess(pack, X, landau_field_print_access_callback, NULL));
415*6b664d00Smarkadams4   //PetscCall(Monitor(ts, -1, 1.0, X, ctx));
416e0eea495SMark   /* clean up */
417*6b664d00Smarkadams4   ts_nrl = (TS)ctx->data;
418*6b664d00Smarkadams4   if (print_nrl) {
419*6b664d00Smarkadams4     PetscReal    finalTime, dt_real, tstart = dt0 * ctx->t_0 / 2; // hack
420*6b664d00Smarkadams4     Vec          U;
421*6b664d00Smarkadams4     PetscScalar *x;
422*6b664d00Smarkadams4     PetscInt     nsteps;
423*6b664d00Smarkadams4     dt_real = dt0 * ctx->t_0;
424*6b664d00Smarkadams4     PetscCall(TSSetTimeStep(ts_nrl, dt_real));
425*6b664d00Smarkadams4     PetscCall(TSGetTime(ts, &finalTime));
426*6b664d00Smarkadams4     finalTime *= ctx->t_0;
427*6b664d00Smarkadams4     PetscCall(TSSetMaxTime(ts_nrl, finalTime));
428*6b664d00Smarkadams4     nsteps = (PetscInt)(finalTime / dt_real) + 1;
429*6b664d00Smarkadams4     PetscCall(TSSetMaxSteps(ts_nrl, nsteps));
430*6b664d00Smarkadams4     PetscCall(TSSetStepNumber(ts_nrl, 0));
431*6b664d00Smarkadams4     PetscCall(TSSetTime(ts_nrl, tstart));
432*6b664d00Smarkadams4     PetscCall(TSGetSolution(ts_nrl, &U));
433*6b664d00Smarkadams4     PetscCall(VecGetArray(U, &x));
434*6b664d00Smarkadams4     for (PetscInt i = 0; i < NUM_TEMPS; i++) x[i] = T0[i];
435*6b664d00Smarkadams4     PetscCall(VecRestoreArray(U, &x));
436*6b664d00Smarkadams4     PetscCall(TSMonitorSet(ts_nrl, Monitor_nrl, ctx, NULL));
437*6b664d00Smarkadams4     PetscCall(TSSolve(ts_nrl, NULL));
438*6b664d00Smarkadams4   }
4399566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
440*6b664d00Smarkadams4   PetscCall(TSDestroy(&ts_nrl));
4419566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&X));
442*6b664d00Smarkadams4   PetscCall(DMPlexLandauDestroyVelocitySpace(&pack));
4439566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
444b122ec5aSJacob Faibussowitsch   return 0;
445e0eea495SMark }
446e0eea495SMark 
447e0eea495SMark /*TEST
4485dac466eSMark Adams   testset:
449984ed092SMark Adams     requires: p4est !complex double defined(PETSC_USE_DMLANDAU_2D)
4505dac466eSMark Adams     output_file: output/ex1_0.out
451*6b664d00Smarkadams4     filter: grep -v "DM"
452*6b664d00Smarkadams4     args: -dm_landau_amr_levels_max 0,2 -dm_landau_amr_post_refine 0 -dm_landau_amr_re_levels 2 -dm_landau_domain_radius 6,6 -dm_landau_electron_shift 1.5 -dm_landau_ion_charges 1 -dm_landau_ion_masses 2 -dm_landau_n 1,1 -dm_landau_n_0 1e20 -dm_landau_num_cells 2,4 -dm_landau_num_species_grid 1,1 -dm_landau_re_radius 2 -use_nrl true -print_nrl false -dm_landau_thermal_temps .3,.2 -dm_landau_type p4est -dm_landau_verbose -1 -dm_preallocate_only false -ex1_dm_view_e -ksp_type preonly -pc_type lu -petscspace_degree 3 -snes_converged_reason -snes_rtol 1.e-14 -snes_stol 1.e-14 -ts_adapt_clip .5,1.5 -ts_adapt_dt_max 5 -ts_adapt_monitor -ts_adapt_scale_solve_failed 0.5 -ts_arkimex_type 1bee -ts_dt .01 -ts_max_snes_failures -1 -ts_max_steps 1 -ts_max_time 8 -ts_monitor -ts_rtol 1e-2 -ts_type arkimex
453e0eea495SMark     test:
4545dac466eSMark Adams       suffix: cpu
455*6b664d00Smarkadams4       args: -dm_landau_device_type cpu -dm_landau_use_relativistic_corrections
4565dac466eSMark Adams     test:
4575dac466eSMark Adams       suffix: kokkos
458e0bd6ad4SStefano Zampini       requires: kokkos_kernels !defined(PETSC_HAVE_CUDA_CLANG)
4595dac466eSMark Adams       args: -dm_landau_device_type kokkos -dm_mat_type aijkokkos -dm_vec_type kokkos
4605dac466eSMark Adams     test:
4615dac466eSMark Adams       suffix: cuda
462e0bd6ad4SStefano Zampini       requires: cuda !defined(PETSC_HAVE_CUDA_CLANG)
4635dac466eSMark Adams       args: -dm_landau_device_type cuda -dm_mat_type aijcusparse -dm_vec_type cuda -mat_cusparse_use_cpu_solve
464e0eea495SMark 
465*6b664d00Smarkadams4   testset:
466*6b664d00Smarkadams4     requires: !complex defined(PETSC_USE_DMLANDAU_2D) !cuda p4est
467*6b664d00Smarkadams4     args: -dm_landau_type p4est -dm_landau_num_cells 4,4 -dm_landau_amr_levels_max 3,3 -dm_landau_num_species_grid 1,1 -dm_landau_n 1,1 -dm_landau_thermal_temps 1,1 -dm_landau_ion_charges 1 -dm_landau_ion_masses 2 -petscspace_degree 2 -ts_type beuler -ts_dt .1 -ts_max_steps 0 -dm_landau_verbose 2 -ksp_type preonly -pc_type lu -dm_landau_device_type cpu -use_nrl false -print_nrl -snes_rtol 1.e-14 -snes_stol 1.e-14
468*6b664d00Smarkadams4     nsize: 1
469*6b664d00Smarkadams4     test:
470*6b664d00Smarkadams4       suffix: sphere
471*6b664d00Smarkadams4       args: -dm_landau_sphere -ts_max_steps 0 -dm_landau_amr_post_refine 1
472*6b664d00Smarkadams4     test:
473*6b664d00Smarkadams4       suffix: re
474*6b664d00Smarkadams4       args: -dm_landau_amr_levels_max 0,2 -dm_landau_z_radius_pre 2.5 -dm_landau_z_radius_post 3.75 -dm_landau_amr_z_refine_pre 1 -dm_landau_amr_z_refine_post 1 -dm_landau_electron_shift 1.25 -ts_max_steps 1 -snes_converged_reason -info :vec
475*6b664d00Smarkadams4 
476e0eea495SMark TEST*/
477