1*3ecd1e40SMark Adams static char help[] = "Landau collision operator with anisotropic 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 "
26b664d00Smarkadams4 "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>
66b664d00Smarkadams4 #include <petscdmcomposite.h>
76b664d00Smarkadams4 #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
192fe279fdSBarry Smith Input/Output Parameter:
202fe279fdSBarry Smith . x - Vector to data to
2124ded41bSmarkadams4
2224ded41bSmarkadams4 */
landau_field_print_access_callback(DM dm,Vec x,PetscInt local_field,PetscInt grid,PetscInt b_id,void * vctx)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;
266b664d00Smarkadams4 PetscScalar val;
276b664d00Smarkadams4 PetscInt species;
2824ded41bSmarkadams4
2924ded41bSmarkadams4 PetscFunctionBegin;
30f53e7263SMark Adams PetscCall(DMGetApplicationContext(dm, &ctx));
316b664d00Smarkadams4 species = ctx->species_offset[grid] + local_field;
326b664d00Smarkadams4 val = (PetscScalar)(LAND_PACK_IDX(b_id, grid) + (species + 1) * 10);
336b664d00Smarkadams4 PetscCall(VecSet(x, val));
346b664d00Smarkadams4 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));
356b664d00Smarkadams4 PetscFunctionReturn(PETSC_SUCCESS);
3624ded41bSmarkadams4 }
376b664d00Smarkadams4
386b664d00Smarkadams4 static const PetscReal alphai = 1 / 1.3;
396b664d00Smarkadams4 static const PetscReal kev_joul = 6.241506479963235e+15; /* 1/1000e */
406b664d00Smarkadams4
416b664d00Smarkadams4 // constants: [index of (anisotropic) direction of source, z x[1] shift
426b664d00Smarkadams4 /* < v, n_s v_|| > */
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)436b664d00Smarkadams4 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)
446b664d00Smarkadams4 {
456b664d00Smarkadams4 if (dim == 2) f0[0] = u[0] * (2. * PETSC_PI * x[0]) * x[1]; /* n r v_|| */
466b664d00Smarkadams4 else f0[0] = u[0] * x[2];
4724ded41bSmarkadams4 }
486b664d00Smarkadams4 /* < v, n (v-shift)^2 > */
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)496b664d00Smarkadams4 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)
506b664d00Smarkadams4 {
516b664d00Smarkadams4 PetscReal vz = PetscRealPart(constants[0]);
526b664d00Smarkadams4 if (dim == 2) *f0 = u[0] * (2. * PETSC_PI * x[0]) * (x[1] - vz) * (x[1] - vz); /* n r v^2_par|perp */
536b664d00Smarkadams4 else *f0 = u[0] * (x[2] - vz) * (x[2] - vz);
546b664d00Smarkadams4 }
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)556b664d00Smarkadams4 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)
566b664d00Smarkadams4 {
576b664d00Smarkadams4 if (dim == 2) *f0 = u[0] * (2. * PETSC_PI * x[0]) * x[0] * x[0]; /* n r v^2_perp */
586b664d00Smarkadams4 else *f0 = u[0] * (x[0] * x[0] + x[1] * x[1]);
596b664d00Smarkadams4 }
606b664d00Smarkadams4 /* < v, n_e > */
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)616b664d00Smarkadams4 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)
626b664d00Smarkadams4 {
636b664d00Smarkadams4 if (dim == 2) f0[0] = 2. * PETSC_PI * x[0] * u[0];
646b664d00Smarkadams4 else f0[0] = u[0];
656b664d00Smarkadams4 }
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)666b664d00Smarkadams4 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)
676b664d00Smarkadams4 {
686b664d00Smarkadams4 PetscReal vz = PetscRealPart(constants[0]);
696b664d00Smarkadams4 if (dim == 2) f0[0] = u[0] * (2. * PETSC_PI * x[0]) * (x[0] * x[0] + (x[1] - vz) * (x[1] - vz));
706b664d00Smarkadams4 else f0[0] = u[0] * (x[0] * x[0] + x[1] * x[1] + (x[2] - vz) * (x[2] - vz));
716b664d00Smarkadams4 }
726b664d00Smarkadams4 /* Define a Maxwellian function for testing out the operator. */
736b664d00Smarkadams4 typedef struct {
746b664d00Smarkadams4 PetscReal v_0;
756b664d00Smarkadams4 PetscReal kT_m;
766b664d00Smarkadams4 PetscReal n;
776b664d00Smarkadams4 PetscReal shift;
786b664d00Smarkadams4 PetscInt species;
796b664d00Smarkadams4 } MaxwellianCtx;
806b664d00Smarkadams4
maxwellian(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nf_dummy,PetscScalar * u,void * actx)816b664d00Smarkadams4 static PetscErrorCode maxwellian(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf_dummy, PetscScalar *u, void *actx)
826b664d00Smarkadams4 {
836b664d00Smarkadams4 MaxwellianCtx *mctx = (MaxwellianCtx *)actx;
846b664d00Smarkadams4 PetscReal theta = 2 * mctx->kT_m / (mctx->v_0 * mctx->v_0); /* theta = 2kT/mc^2 */
854d86920dSPierre Jolivet
866b664d00Smarkadams4 PetscFunctionBegin;
876b664d00Smarkadams4 /* evaluate the shifted Maxwellian */
886b664d00Smarkadams4 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);
896b664d00Smarkadams4 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);
906b664d00Smarkadams4 PetscFunctionReturn(PETSC_SUCCESS);
916b664d00Smarkadams4 }
926b664d00Smarkadams4
SetMaxwellians(DM dm,Vec X,PetscReal time,PetscReal temps[],PetscReal ns[],PetscInt grid,PetscReal shifts[],LandauCtx * ctx)936b664d00Smarkadams4 static PetscErrorCode SetMaxwellians(DM dm, Vec X, PetscReal time, PetscReal temps[], PetscReal ns[], PetscInt grid, PetscReal shifts[], LandauCtx *ctx)
946b664d00Smarkadams4 {
956b664d00Smarkadams4 PetscErrorCode (*initu[LANDAU_MAX_SPECIES])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *);
966b664d00Smarkadams4 MaxwellianCtx *mctxs[LANDAU_MAX_SPECIES], data[LANDAU_MAX_SPECIES];
974d86920dSPierre Jolivet
986b664d00Smarkadams4 PetscFunctionBegin;
996b664d00Smarkadams4 if (!ctx) PetscCall(DMGetApplicationContext(dm, &ctx));
1006b664d00Smarkadams4 for (PetscInt ii = ctx->species_offset[grid], i0 = 0; ii < ctx->species_offset[grid + 1]; ii++, i0++) {
1016b664d00Smarkadams4 mctxs[i0] = &data[i0];
1026b664d00Smarkadams4 data[i0].v_0 = ctx->v_0; // v_0 same for all grids
1036b664d00Smarkadams4 data[i0].kT_m = ctx->k * temps[ii] / ctx->masses[ii]; /* kT/m = v_th ^ 2*/
1046b664d00Smarkadams4 data[i0].n = ns[ii];
1056b664d00Smarkadams4 initu[i0] = maxwellian;
1066b664d00Smarkadams4 data[i0].shift = 0;
1076b664d00Smarkadams4 data[i0].species = ii;
1086b664d00Smarkadams4 }
1096b664d00Smarkadams4 if (1) {
11057508eceSPierre Jolivet data[0].shift = (PetscReal)-PetscSign(ctx->charges[ctx->species_offset[grid]]) * ctx->electronShift * ctx->m_0 / ctx->masses[ctx->species_offset[grid]];
1116b664d00Smarkadams4 } else {
1126b664d00Smarkadams4 shifts[0] = 0.5 * PetscSqrtReal(ctx->masses[0] / ctx->masses[1]);
1136b664d00Smarkadams4 shifts[1] = 50 * (ctx->masses[0] / ctx->masses[1]);
1146b664d00Smarkadams4 data[0].shift = ctx->electronShift * shifts[grid] * PetscSqrtReal(data[0].kT_m) / ctx->v_0; // shifts to not matter!!!!
1156b664d00Smarkadams4 }
1166b664d00Smarkadams4 PetscCall(DMProjectFunction(dm, time, initu, (void **)mctxs, INSERT_ALL_VALUES, X));
1176b664d00Smarkadams4 PetscFunctionReturn(PETSC_SUCCESS);
1186b664d00Smarkadams4 }
1196b664d00Smarkadams4
1206b664d00Smarkadams4 typedef enum {
1216b664d00Smarkadams4 E_PAR_IDX,
1226b664d00Smarkadams4 E_PERP_IDX,
1236b664d00Smarkadams4 I_PAR_IDX,
1246b664d00Smarkadams4 I_PERP_IDX,
1256b664d00Smarkadams4 NUM_TEMPS
1266b664d00Smarkadams4 } TemperatureIDX;
1276b664d00Smarkadams4
128cd27c6deSmarkadams4 /* -------------------- Evaluate NRL Function F(x) (analytical solutions exist for this) --------------------- */
1296b664d00Smarkadams4 static PetscReal n_cm3[2] = {0, 0};
FormFunction(TS ts,PetscReal tdummy,Vec X,Vec F,void * ptr)1306b664d00Smarkadams4 PetscErrorCode FormFunction(TS ts, PetscReal tdummy, Vec X, Vec F, void *ptr)
1316b664d00Smarkadams4 {
1326b664d00Smarkadams4 LandauCtx *ctx = (LandauCtx *)ptr; /* user-defined application context */
1336b664d00Smarkadams4 PetscScalar *f;
1346b664d00Smarkadams4 const PetscScalar *x;
1359fa27a79SStefano Zampini const PetscReal k_B = 1.6e-12, e_cgs = 4.8e-10, proton_mass = 9.1094e-28, m_cgs[2] = {proton_mass, proton_mass * ctx->masses[1] / ctx->masses[0]}; // erg/eV, e, m as per NRL;
136ae6bf4a8Smarkadams4 PetscReal AA, v_bar_ab, vTe, t1, TeDiff, Te, Ti, Tdiff;
1376b664d00Smarkadams4
1386b664d00Smarkadams4 PetscFunctionBeginUser;
1396b664d00Smarkadams4 PetscCall(VecGetArrayRead(X, &x));
1406b664d00Smarkadams4 Te = PetscRealPart(2 * x[E_PERP_IDX] + x[E_PAR_IDX]) / 3, Ti = PetscRealPart(2 * x[I_PERP_IDX] + x[I_PAR_IDX]) / 3;
141ae6bf4a8Smarkadams4 // thermalization from NRL Plasma formulary, assume Z = 1, mu = 2, n_i = n_e
142ae6bf4a8Smarkadams4 v_bar_ab = 1.8e-19 * PetscSqrtReal(m_cgs[0] * m_cgs[1]) * n_cm3[0] * ctx->lambdas[0][1] * PetscPowReal(m_cgs[0] * Ti + m_cgs[1] * Te, -1.5);
1436b664d00Smarkadams4 PetscCall(VecGetArray(F, &f));
1446b664d00Smarkadams4 for (PetscInt ii = 0; ii < 2; ii++) {
145ae6bf4a8Smarkadams4 PetscReal tPerp = PetscRealPart(x[2 * ii + E_PERP_IDX]), tPar = PetscRealPart(x[2 * ii + E_PAR_IDX]), ff;
1466b664d00Smarkadams4 TeDiff = tPerp - tPar;
1476b664d00Smarkadams4 AA = tPerp / tPar - 1;
148ae6bf4a8Smarkadams4 if (AA < 0) ff = PetscAtanhReal(PetscSqrtReal(-AA)) / PetscSqrtReal(-AA);
149ae6bf4a8Smarkadams4 else ff = PetscAtanReal(PetscSqrtReal(AA)) / PetscSqrtReal(AA);
150ae6bf4a8Smarkadams4 t1 = (-3 + (AA + 3) * ff) / PetscSqr(AA);
151ae6bf4a8Smarkadams4 //PetscReal vTeB = 8.2e-7 * n_cm3[0] * ctx->lambdas[0][1] * PetscPowReal(Te, -1.5);
1529fa27a79SStefano Zampini vTe = 2 * PetscSqrtReal(PETSC_PI / m_cgs[ii]) * PetscSqr(PetscSqr(e_cgs)) * n_cm3[0] * ctx->lambdas[0][1] * PetscPowReal(PetscRealPart(k_B * x[E_PAR_IDX]), -1.5) * t1;
153ae6bf4a8Smarkadams4 t1 = vTe * TeDiff; // * 2; // scaling from NRL that makes it fit pretty good
154ae6bf4a8Smarkadams4
1556b664d00Smarkadams4 f[2 * ii + E_PAR_IDX] = 2 * t1; // par
1566b664d00Smarkadams4 f[2 * ii + E_PERP_IDX] = -t1; // perp
1576b664d00Smarkadams4 Tdiff = (ii == 0) ? (Ti - Te) : (Te - Ti);
158ae6bf4a8Smarkadams4 f[2 * ii + E_PAR_IDX] += v_bar_ab * Tdiff;
159ae6bf4a8Smarkadams4 f[2 * ii + E_PERP_IDX] += v_bar_ab * Tdiff;
1606b664d00Smarkadams4 }
1616b664d00Smarkadams4 PetscCall(VecRestoreArrayRead(X, &x));
1626b664d00Smarkadams4 PetscCall(VecRestoreArray(F, &f));
1636b664d00Smarkadams4 PetscFunctionReturn(PETSC_SUCCESS);
1646b664d00Smarkadams4 }
1656b664d00Smarkadams4
1666b664d00Smarkadams4 /* -------------------- Form initial approximation ----------------- */
1676b664d00Smarkadams4 static PetscReal T0[4] = {300, 390, 200, 260};
createVec_NRL(LandauCtx * ctx,Vec * vec)1686b664d00Smarkadams4 PetscErrorCode createVec_NRL(LandauCtx *ctx, Vec *vec)
1696b664d00Smarkadams4 {
1706b664d00Smarkadams4 PetscScalar *x;
1716b664d00Smarkadams4 Vec Temps;
1726b664d00Smarkadams4
1736b664d00Smarkadams4 PetscFunctionBeginUser;
1746b664d00Smarkadams4 PetscCall(VecCreateSeq(PETSC_COMM_SELF, NUM_TEMPS, &Temps));
1756b664d00Smarkadams4 PetscCall(VecGetArray(Temps, &x));
1766b664d00Smarkadams4 for (PetscInt i = 0; i < NUM_TEMPS; i++) x[i] = T0[i];
1776b664d00Smarkadams4 PetscCall(VecRestoreArray(Temps, &x));
1786b664d00Smarkadams4 *vec = Temps;
1796b664d00Smarkadams4 PetscFunctionReturn(PETSC_SUCCESS);
1806b664d00Smarkadams4 }
1816b664d00Smarkadams4
createTS_NRL(LandauCtx * ctx,Vec Temps)1826b664d00Smarkadams4 PetscErrorCode createTS_NRL(LandauCtx *ctx, Vec Temps)
1836b664d00Smarkadams4 {
1846b664d00Smarkadams4 TSAdapt adapt;
1856b664d00Smarkadams4 TS ts;
1866b664d00Smarkadams4
1876b664d00Smarkadams4 PetscFunctionBeginUser;
1886b664d00Smarkadams4 PetscCall(TSCreate(PETSC_COMM_SELF, &ts));
1896b664d00Smarkadams4 ctx->data = (void *)ts; // 'data' is for applications (eg, monitors)
1906b664d00Smarkadams4 PetscCall(TSSetApplicationContext(ts, ctx));
1916b664d00Smarkadams4 PetscCall(TSSetType(ts, TSRK));
1926b664d00Smarkadams4 PetscCall(TSSetRHSFunction(ts, NULL, FormFunction, ctx));
1936b664d00Smarkadams4 PetscCall(TSSetSolution(ts, Temps));
1946b664d00Smarkadams4 PetscCall(TSRKSetType(ts, TSRK2A));
1956b664d00Smarkadams4 PetscCall(TSSetOptionsPrefix(ts, "nrl_"));
1966b664d00Smarkadams4 PetscCall(TSSetFromOptions(ts));
1976b664d00Smarkadams4 PetscCall(TSGetAdapt(ts, &adapt));
1986b664d00Smarkadams4 PetscCall(TSAdaptSetType(adapt, TSADAPTNONE));
1996b664d00Smarkadams4 PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
2006b664d00Smarkadams4 PetscCall(TSSetStepNumber(ts, 0));
2016b664d00Smarkadams4 PetscCall(TSSetMaxSteps(ts, 1));
2026b664d00Smarkadams4 PetscCall(TSSetTime(ts, 0));
2036b664d00Smarkadams4 PetscFunctionReturn(PETSC_SUCCESS);
2046b664d00Smarkadams4 }
2056b664d00Smarkadams4
Monitor_nrl(TS ts,PetscInt stepi,PetscReal time,Vec X,void * actx)2066b664d00Smarkadams4 PetscErrorCode Monitor_nrl(TS ts, PetscInt stepi, PetscReal time, Vec X, void *actx)
2076b664d00Smarkadams4 {
2086b664d00Smarkadams4 const PetscScalar *x;
2096b664d00Smarkadams4 LandauCtx *ctx = (LandauCtx *)actx; /* user-defined application context */
2106b664d00Smarkadams4
2116b664d00Smarkadams4 PetscFunctionBeginUser;
212ae6bf4a8Smarkadams4 if (stepi % 100 == 0) {
2136b664d00Smarkadams4 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "nrl-step %d time= %g ", (int)stepi, (double)(time / ctx->t_0)));
2146b664d00Smarkadams4 PetscCall(VecGetArrayRead(X, &x));
2153a7d0413SPierre Jolivet for (PetscInt i = 0; i < NUM_TEMPS; i++) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%g ", (double)PetscRealPart(x[i])));
2166b664d00Smarkadams4 PetscCall(VecRestoreArrayRead(X, &x));
2176b664d00Smarkadams4 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n"));
2186b664d00Smarkadams4 }
2196b664d00Smarkadams4 PetscFunctionReturn(PETSC_SUCCESS);
2206b664d00Smarkadams4 }
2216b664d00Smarkadams4
Monitor(TS ts,PetscInt stepi,PetscReal time,Vec X,void * actx)2226b664d00Smarkadams4 PetscErrorCode Monitor(TS ts, PetscInt stepi, PetscReal time, Vec X, void *actx)
2236b664d00Smarkadams4 {
2246b664d00Smarkadams4 LandauCtx *ctx = (LandauCtx *)actx; /* user-defined application context */
2256b664d00Smarkadams4 TS ts_nrl = (TS)ctx->data;
2266b664d00Smarkadams4 PetscInt printing = 0, logT;
2276b664d00Smarkadams4
2286b664d00Smarkadams4 PetscFunctionBeginUser;
2296b664d00Smarkadams4 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)
2306b664d00Smarkadams4 PetscReal dt;
2316b664d00Smarkadams4 PetscCall(TSGetTimeStep(ts, &dt));
2326b664d00Smarkadams4 logT = (PetscInt)PetscLog2Real(time / dt);
2336b664d00Smarkadams4 if (logT < 0) logT = 0;
2346b664d00Smarkadams4 ctx->verbose = PetscPowInt(2, logT) / 2;
2356b664d00Smarkadams4 if (ctx->verbose == 0) ctx->verbose = 1;
2366b664d00Smarkadams4 }
2376b664d00Smarkadams4 if (ctx->verbose) {
2386b664d00Smarkadams4 TSConvergedReason reason;
2396b664d00Smarkadams4 PetscCall(TSGetConvergedReason(ts, &reason));
2406b664d00Smarkadams4 if (stepi % ctx->verbose == 0 || reason || stepi == 1 || ctx->verbose < 0) {
2416b664d00Smarkadams4 PetscInt nDMs, id;
2426b664d00Smarkadams4 DM pack;
2436b664d00Smarkadams4 Vec *XsubArray = NULL;
2446b664d00Smarkadams4 printing = 1;
2456b664d00Smarkadams4 PetscCall(TSGetDM(ts, &pack));
2466b664d00Smarkadams4 PetscCall(DMCompositeGetNumberDM(pack, &nDMs));
2476b664d00Smarkadams4 PetscCall(DMGetOutputSequenceNumber(ctx->plex[0], &id, NULL));
2486b664d00Smarkadams4 PetscCall(DMSetOutputSequenceNumber(ctx->plex[0], id + 1, time));
2496b664d00Smarkadams4 PetscCall(DMSetOutputSequenceNumber(ctx->plex[1], id + 1, time));
2506b664d00Smarkadams4 PetscCall(PetscInfo(pack, "ex1 plot step %" PetscInt_FMT ", time = %g\n", id, (double)time));
2516b664d00Smarkadams4 PetscCall(PetscMalloc(sizeof(*XsubArray) * nDMs, &XsubArray));
2526b664d00Smarkadams4 PetscCall(DMCompositeGetAccessArray(pack, X, nDMs, NULL, XsubArray)); // read only
2536b664d00Smarkadams4 PetscCall(VecViewFromOptions(XsubArray[LAND_PACK_IDX(ctx->batch_view_idx, 0)], NULL, "-ex1_vec_view_e"));
2546b664d00Smarkadams4 PetscCall(VecViewFromOptions(XsubArray[LAND_PACK_IDX(ctx->batch_view_idx, 1)], NULL, "-ex1_vec_view_i"));
2556b664d00Smarkadams4 // temps
2566b664d00Smarkadams4 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
2576b664d00Smarkadams4 PetscDS prob;
2586b664d00Smarkadams4 DM dm = ctx->plex[grid];
2596b664d00Smarkadams4 PetscScalar user[2] = {0, 0}, tt[1];
2606b664d00Smarkadams4 PetscReal vz_0 = 0, n, energy, e_perp, e_par, m_s = ctx->masses[ctx->species_offset[grid]];
2616b664d00Smarkadams4 Vec Xloc = XsubArray[LAND_PACK_IDX(ctx->batch_view_idx, grid)];
2626b664d00Smarkadams4 PetscCall(DMGetDS(dm, &prob));
2636b664d00Smarkadams4 /* get n */
2646b664d00Smarkadams4 PetscCall(PetscDSSetObjective(prob, 0, &f0_n));
2656b664d00Smarkadams4 PetscCall(DMPlexComputeIntegralFEM(dm, Xloc, tt, NULL));
2666b664d00Smarkadams4 n = PetscRealPart(tt[0]);
2676b664d00Smarkadams4 /* get vz */
2686b664d00Smarkadams4 PetscCall(PetscDSSetObjective(prob, 0, &f0_vz));
2696b664d00Smarkadams4 PetscCall(DMPlexComputeIntegralFEM(dm, Xloc, tt, NULL));
2706b664d00Smarkadams4 user[0] = vz_0 = PetscRealPart(tt[0]) / n;
2716b664d00Smarkadams4 /* energy temp */
2726b664d00Smarkadams4 PetscCall(PetscDSSetConstants(prob, 2, user));
2736b664d00Smarkadams4 PetscCall(PetscDSSetObjective(prob, 0, &f0_v2_shift));
2746b664d00Smarkadams4 PetscCall(DMPlexComputeIntegralFEM(dm, Xloc, tt, ctx));
2756b664d00Smarkadams4 energy = PetscRealPart(tt[0]) * ctx->v_0 * ctx->v_0 * m_s / n / 3; // scale?
2766b664d00Smarkadams4 energy *= kev_joul * 1000; // T eV
2776b664d00Smarkadams4 /* energy temp - par */
2786b664d00Smarkadams4 PetscCall(PetscDSSetConstants(prob, 2, user));
2796b664d00Smarkadams4 PetscCall(PetscDSSetObjective(prob, 0, &f0_v2_par_shift));
2806b664d00Smarkadams4 PetscCall(DMPlexComputeIntegralFEM(dm, Xloc, tt, ctx));
2816b664d00Smarkadams4 e_par = PetscRealPart(tt[0]) * ctx->v_0 * ctx->v_0 * m_s / n;
2826b664d00Smarkadams4 e_par *= kev_joul * 1000; // eV
2836b664d00Smarkadams4 /* energy temp - perp */
2846b664d00Smarkadams4 PetscCall(PetscDSSetConstants(prob, 2, user));
2856b664d00Smarkadams4 PetscCall(PetscDSSetObjective(prob, 0, &f0_v2_perp));
2866b664d00Smarkadams4 PetscCall(DMPlexComputeIntegralFEM(dm, Xloc, tt, ctx));
2876b664d00Smarkadams4 e_perp = PetscRealPart(tt[0]) * ctx->v_0 * ctx->v_0 * m_s / n / 2;
2886b664d00Smarkadams4 e_perp *= kev_joul * 1000; // eV
2896b664d00Smarkadams4 if (grid == 0) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "step %4d) time= %e temperature (eV): ", (int)stepi, (double)time));
2906b664d00Smarkadams4 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));
2916b664d00Smarkadams4 if (n_cm3[grid] == 0) n_cm3[grid] = ctx->n_0 * n * 1e-6; // does not change m^3 --> cm^3
2926b664d00Smarkadams4 }
2936b664d00Smarkadams4 // cleanup
2946b664d00Smarkadams4 PetscCall(DMCompositeRestoreAccessArray(pack, X, nDMs, NULL, XsubArray));
2956b664d00Smarkadams4 PetscCall(PetscFree(XsubArray));
2966b664d00Smarkadams4 }
2976b664d00Smarkadams4 }
2986b664d00Smarkadams4 /* evolve NRL data, end line */
2996b664d00Smarkadams4 if (n_cm3[NUM_TEMPS / 2 - 1] < 0 && ts_nrl) {
3006b664d00Smarkadams4 PetscCall(TSDestroy(&ts_nrl));
3016b664d00Smarkadams4 ctx->data = NULL;
3026b664d00Smarkadams4 if (printing) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nSTOP printing NRL Ts\n"));
3036b664d00Smarkadams4 } else if (ts_nrl) {
3046b664d00Smarkadams4 const PetscScalar *x;
3056b664d00Smarkadams4 PetscReal dt_real, dt;
3066b664d00Smarkadams4 Vec U;
3076b664d00Smarkadams4 PetscCall(TSGetTimeStep(ts, &dt)); // dt for NEXT time step
3086b664d00Smarkadams4 dt_real = dt * ctx->t_0;
3096b664d00Smarkadams4 PetscCall(TSGetSolution(ts_nrl, &U));
3106b664d00Smarkadams4 if (printing) {
3116b664d00Smarkadams4 PetscCall(VecGetArrayRead(U, &x));
312d043ef4cSMark Adams PetscCall(PetscPrintf(PETSC_COMM_WORLD, "NRL_e_par= %9.4e NRL_e_perp= %9.4e ", (double)PetscRealPart(x[E_PAR_IDX]), (double)PetscRealPart(x[E_PERP_IDX])));
313d043ef4cSMark Adams PetscCall(PetscPrintf(PETSC_COMM_WORLD, "NRL_i_par= %9.4e NRL_i_perp= %9.4e\n", (double)PetscRealPart(x[I_PAR_IDX]), (double)PetscRealPart(x[I_PERP_IDX])));
314d043ef4cSMark Adams /* if (n_cm3[0] > 0) { */
315d043ef4cSMark Adams /* } else PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n")); */
3166b664d00Smarkadams4 PetscCall(VecRestoreArrayRead(U, &x));
3176b664d00Smarkadams4 }
3186b664d00Smarkadams4 // we have the next time step, so need to advance now
3196b664d00Smarkadams4 PetscCall(TSSetTimeStep(ts_nrl, dt_real));
3206b664d00Smarkadams4 PetscCall(TSSetMaxSteps(ts_nrl, stepi + 1)); // next step
3216b664d00Smarkadams4 PetscCall(TSSolve(ts_nrl, NULL));
3226b664d00Smarkadams4 } else if (printing) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n"));
3233a7d0413SPierre Jolivet if (printing) PetscCall(DMPlexLandauPrintNorms(X, stepi /*id + 1*/));
3243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
32524ded41bSmarkadams4 }
32624ded41bSmarkadams4
main(int argc,char ** argv)327d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
328d71ae5a4SJacob Faibussowitsch {
32924ded41bSmarkadams4 DM pack;
3306b664d00Smarkadams4 Vec X;
3316b664d00Smarkadams4 PetscInt dim = 2, nDMs;
3326b664d00Smarkadams4 TS ts, ts_nrl = NULL;
333e0eea495SMark Mat J;
3346b664d00Smarkadams4 Vec *XsubArray = NULL;
3356b664d00Smarkadams4 LandauCtx *ctx;
3366b664d00Smarkadams4 PetscMPIInt rank;
337ae6bf4a8Smarkadams4 PetscBool use_nrl = PETSC_TRUE;
338ae6bf4a8Smarkadams4 PetscBool print_nrl = PETSC_FALSE;
3396b664d00Smarkadams4 PetscReal dt0;
3404d86920dSPierre Jolivet
341327415f7SBarry Smith PetscFunctionBeginUser;
3429566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help));
3436b664d00Smarkadams4 PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
3446b664d00Smarkadams4 if (rank) { /* turn off output stuff for duplicate runs */
3456b664d00Smarkadams4 PetscCall(PetscOptionsClearValue(NULL, "-ex1_dm_view_e"));
3466b664d00Smarkadams4 PetscCall(PetscOptionsClearValue(NULL, "-ex1_dm_view_i"));
3476b664d00Smarkadams4 PetscCall(PetscOptionsClearValue(NULL, "-ex1_vec_view_e"));
3486b664d00Smarkadams4 PetscCall(PetscOptionsClearValue(NULL, "-ex1_vec_view_i"));
3496b664d00Smarkadams4 PetscCall(PetscOptionsClearValue(NULL, "-info"));
3506b664d00Smarkadams4 PetscCall(PetscOptionsClearValue(NULL, "-snes_converged_reason"));
3516b664d00Smarkadams4 PetscCall(PetscOptionsClearValue(NULL, "-pc_bjkokkos_ksp_converged_reason"));
3526b664d00Smarkadams4 PetscCall(PetscOptionsClearValue(NULL, "-ksp_converged_reason"));
3536b664d00Smarkadams4 PetscCall(PetscOptionsClearValue(NULL, "-ts_adapt_monitor"));
3546b664d00Smarkadams4 PetscCall(PetscOptionsClearValue(NULL, "-ts_monitor"));
3556b664d00Smarkadams4 PetscCall(PetscOptionsClearValue(NULL, "-snes_monitor"));
3566b664d00Smarkadams4 }
3579566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-dim", &dim, NULL));
3586b664d00Smarkadams4 PetscCall(PetscOptionsGetBool(NULL, NULL, "-use_nrl", &use_nrl, NULL));
3596b664d00Smarkadams4 PetscCall(PetscOptionsGetBool(NULL, NULL, "-print_nrl", &print_nrl, NULL));
360e0eea495SMark /* Create a mesh */
36124ded41bSmarkadams4 PetscCall(DMPlexLandauCreateVelocitySpace(PETSC_COMM_SELF, dim, "", &X, &J, &pack));
36224ded41bSmarkadams4 PetscCall(DMSetUp(pack));
3636b664d00Smarkadams4 PetscCall(DMGetApplicationContext(pack, &ctx));
3646b664d00Smarkadams4 PetscCheck(ctx->num_grids == 2, PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Must have two grids: use '-dm_landau_num_species_grid 1,1'");
3656b664d00Smarkadams4 PetscCheck(ctx->num_species == 2, PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Must have two species: use '-dm_landau_num_species_grid 1,1'");
3666b664d00Smarkadams4 PetscCall(DMCompositeGetNumberDM(pack, &nDMs));
3676b664d00Smarkadams4 /* output plot names */
3686b664d00Smarkadams4 PetscCall(PetscMalloc(sizeof(*XsubArray) * nDMs, &XsubArray));
3696b664d00Smarkadams4 PetscCall(DMCompositeGetAccessArray(pack, X, nDMs, NULL, XsubArray)); // read only
3706b664d00Smarkadams4 PetscCall(PetscObjectSetName((PetscObject)XsubArray[LAND_PACK_IDX(ctx->batch_view_idx, 0)], 0 == 0 ? "ue" : "ui"));
3716b664d00Smarkadams4 PetscCall(PetscObjectSetName((PetscObject)XsubArray[LAND_PACK_IDX(ctx->batch_view_idx, 1)], 1 == 0 ? "ue" : "ui"));
3726b664d00Smarkadams4 /* add bimaxwellian anisotropic test */
3736b664d00Smarkadams4 for (PetscInt b_id = 0; b_id < ctx->batch_sz; b_id++) {
3746b664d00Smarkadams4 for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
3756b664d00Smarkadams4 PetscReal shifts[2];
3766b664d00Smarkadams4 PetscCall(SetMaxwellians(ctx->plex[grid], XsubArray[LAND_PACK_IDX(b_id, grid)], 0.0, ctx->thermal_temps, ctx->n, grid, shifts, ctx));
3776b664d00Smarkadams4 }
3786b664d00Smarkadams4 }
3796b664d00Smarkadams4 PetscCall(DMCompositeRestoreAccessArray(pack, X, nDMs, NULL, XsubArray));
3806b664d00Smarkadams4 PetscCall(PetscFree(XsubArray));
3816b664d00Smarkadams4 /* plot */
3826b664d00Smarkadams4 PetscCall(DMSetOutputSequenceNumber(ctx->plex[0], -1, 0.0));
3836b664d00Smarkadams4 PetscCall(DMSetOutputSequenceNumber(ctx->plex[1], -1, 0.0));
3846b664d00Smarkadams4 PetscCall(DMViewFromOptions(ctx->plex[0], NULL, "-ex1_dm_view_e"));
3856b664d00Smarkadams4 PetscCall(DMViewFromOptions(ctx->plex[1], NULL, "-ex1_dm_view_i"));
386e0eea495SMark /* Create timestepping solver context */
3879566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_SELF, &ts));
38824ded41bSmarkadams4 PetscCall(TSSetDM(ts, pack));
3899566063dSJacob Faibussowitsch PetscCall(TSSetIFunction(ts, NULL, DMPlexLandauIFunction, NULL));
3909566063dSJacob Faibussowitsch PetscCall(TSSetIJacobian(ts, J, J, DMPlexLandauIJacobian, NULL));
3919566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
3929566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts));
3939566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts, X));
3946b664d00Smarkadams4 PetscCall(TSMonitorSet(ts, Monitor, ctx, NULL));
3956b664d00Smarkadams4 /* Create NRL timestepping */
3966b664d00Smarkadams4 if (use_nrl || print_nrl) {
3976b664d00Smarkadams4 Vec NRL_vec;
3986b664d00Smarkadams4 PetscCall(createVec_NRL(ctx, &NRL_vec));
3996b664d00Smarkadams4 PetscCall(createTS_NRL(ctx, NRL_vec));
4006b664d00Smarkadams4 PetscCall(VecDestroy(&NRL_vec));
4016b664d00Smarkadams4 } else ctx->data = NULL;
4026b664d00Smarkadams4 /* solve */
4036b664d00Smarkadams4 PetscCall(TSGetTimeStep(ts, &dt0));
4046b664d00Smarkadams4 PetscCall(TSSetTime(ts, dt0 / 2));
4059566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, X));
4066b664d00Smarkadams4 /* test add field method & output */
4076b664d00Smarkadams4 PetscCall(DMPlexLandauAccess(pack, X, landau_field_print_access_callback, NULL));
408cd27c6deSmarkadams4 // run NRL in separate TS
4096b664d00Smarkadams4 ts_nrl = (TS)ctx->data;
4106b664d00Smarkadams4 if (print_nrl) {
4116b664d00Smarkadams4 PetscReal finalTime, dt_real, tstart = dt0 * ctx->t_0 / 2; // hack
4126b664d00Smarkadams4 Vec U;
4136b664d00Smarkadams4 PetscScalar *x;
4146b664d00Smarkadams4 PetscInt nsteps;
4156b664d00Smarkadams4 dt_real = dt0 * ctx->t_0;
4166b664d00Smarkadams4 PetscCall(TSSetTimeStep(ts_nrl, dt_real));
4176b664d00Smarkadams4 PetscCall(TSGetTime(ts, &finalTime));
4186b664d00Smarkadams4 finalTime *= ctx->t_0;
4196b664d00Smarkadams4 PetscCall(TSSetMaxTime(ts_nrl, finalTime));
4206b664d00Smarkadams4 nsteps = (PetscInt)(finalTime / dt_real) + 1;
4216b664d00Smarkadams4 PetscCall(TSSetMaxSteps(ts_nrl, nsteps));
4226b664d00Smarkadams4 PetscCall(TSSetStepNumber(ts_nrl, 0));
4236b664d00Smarkadams4 PetscCall(TSSetTime(ts_nrl, tstart));
4246b664d00Smarkadams4 PetscCall(TSGetSolution(ts_nrl, &U));
4256b664d00Smarkadams4 PetscCall(VecGetArray(U, &x));
4266b664d00Smarkadams4 for (PetscInt i = 0; i < NUM_TEMPS; i++) x[i] = T0[i];
4276b664d00Smarkadams4 PetscCall(VecRestoreArray(U, &x));
4286b664d00Smarkadams4 PetscCall(TSMonitorSet(ts_nrl, Monitor_nrl, ctx, NULL));
4296b664d00Smarkadams4 PetscCall(TSSolve(ts_nrl, NULL));
4306b664d00Smarkadams4 }
431cd27c6deSmarkadams4 /* clean up */
4329566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts));
4336b664d00Smarkadams4 PetscCall(TSDestroy(&ts_nrl));
4349566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X));
4356b664d00Smarkadams4 PetscCall(DMPlexLandauDestroyVelocitySpace(&pack));
4369566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
437b122ec5aSJacob Faibussowitsch return 0;
438e0eea495SMark }
439e0eea495SMark
440e0eea495SMark /*TEST
4415dac466eSMark Adams testset:
442984ed092SMark Adams requires: p4est !complex double defined(PETSC_USE_DMLANDAU_2D)
4435dac466eSMark Adams output_file: output/ex1_0.out
4446b664d00Smarkadams4 filter: grep -v "DM"
445188af4bfSBarry Smith 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_time_step .01 -ts_max_snes_failures unlimited -ts_max_steps 1 -ts_max_time 8 -ts_monitor -ts_rtol 1e-2 -ts_type arkimex
446e0eea495SMark test:
4475dac466eSMark Adams suffix: cpu
4486b664d00Smarkadams4 args: -dm_landau_device_type cpu -dm_landau_use_relativistic_corrections
4495dac466eSMark Adams test:
4505dac466eSMark Adams suffix: kokkos
451e0bd6ad4SStefano Zampini requires: kokkos_kernels !defined(PETSC_HAVE_CUDA_CLANG)
4525dac466eSMark Adams args: -dm_landau_device_type kokkos -dm_mat_type aijkokkos -dm_vec_type kokkos
453e0eea495SMark
4546b664d00Smarkadams4 testset:
455ae6bf4a8Smarkadams4 requires: !complex defined(PETSC_USE_DMLANDAU_2D) p4est
456188af4bfSBarry Smith args: -dm_landau_type p4est -dm_landau_num_species_grid 1,1 -dm_landau_n 1,1 -dm_landau_thermal_temps 2,1 -dm_landau_ion_charges 1 -dm_landau_ion_masses 2 -petscspace_degree 2 -ts_type beuler -ts_time_step .1 -ts_max_steps 1 -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 -snes_converged_reason -dm_landau_device_type cpu
4576b664d00Smarkadams4 nsize: 1
4586b664d00Smarkadams4 test:
4596b664d00Smarkadams4 suffix: sphere
460cd27c6deSmarkadams4 args: -dm_landau_sphere -dm_landau_amr_levels_max 1,1 -dm_landau_sphere_inner_radius_90degree_scale .55 -dm_landau_sphere_inner_radius_45degree_scale .5
4616b664d00Smarkadams4 test:
4626b664d00Smarkadams4 suffix: re
463cd27c6deSmarkadams4 args: -dm_landau_num_cells 4,4 -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 -info :vec
4646b664d00Smarkadams4
465e0eea495SMark TEST*/
466