1e0eea495SMark static char help[] = "Runaway electron model with Landau collision operator\n\n"; 2e0eea495SMark 3e0eea495SMark #include <petscdmplex.h> 4e0eea495SMark #include <petsclandau.h> 5e0eea495SMark #include <petscts.h> 6e0eea495SMark #include <petscds.h> 78a6f2e61SMark Adams #include <petscdmcomposite.h> 846233b44SBarry Smith #include <petsc/private/petscimpl.h> 9e0eea495SMark 10419c2fb8Smarkadams4 #if defined(PETSC_HAVE_CUDA_NVTX) 11419c2fb8Smarkadams4 #include <nvToolsExt.h> 12419c2fb8Smarkadams4 #endif 13419c2fb8Smarkadams4 14e0eea495SMark /* data for runaway electron model */ 15e0eea495SMark typedef struct REctx_struct { 168a6f2e61SMark Adams PetscErrorCode (*test)(TS, Vec, PetscInt, PetscReal, PetscBool, LandauCtx *, struct REctx_struct *); 17e0eea495SMark PetscErrorCode (*impuritySrcRate)(PetscReal, PetscReal *, LandauCtx *); 18e0eea495SMark PetscErrorCode (*E)(Vec, Vec, PetscInt, PetscReal, LandauCtx *, PetscReal *); 19e0eea495SMark PetscReal T_cold; /* temperature of newly ionized electrons and impurity ions */ 20e0eea495SMark PetscReal ion_potential; /* ionization potential of impurity */ 21e0eea495SMark PetscReal Ne_ion; /* effective number of electrons shed in ioization of impurity */ 22e0eea495SMark PetscReal Ez_initial; 23e0eea495SMark PetscReal L; /* inductance */ 24e0eea495SMark Vec X_0; 25e0eea495SMark PetscInt imp_idx; /* index for impurity ionizing sink */ 26e0eea495SMark PetscReal pulse_start; 27e0eea495SMark PetscReal pulse_width; 28e0eea495SMark PetscReal pulse_rate; 29e0eea495SMark PetscReal current_rate; 30e0eea495SMark PetscInt plotIdx; 31e0eea495SMark PetscInt plotStep; 32e0eea495SMark PetscInt idx; /* cache */ 33e0eea495SMark PetscReal j; /* cache */ 34e0eea495SMark PetscReal plotDt; 35e0eea495SMark PetscBool plotting; 36e0eea495SMark PetscBool use_spitzer_eta; 378a6f2e61SMark Adams PetscInt print_period; 388fdabdddSMark Adams PetscInt grid_view_idx; 39e0eea495SMark } REctx; 40e0eea495SMark 41e0eea495SMark static const PetscReal kev_joul = 6.241506479963235e+15; /* 1/1000e */ 42e0eea495SMark 4390d163c9SMark Adams #define RE_CUT 3. 44e0eea495SMark /* < v, u_re * v * q > */ 45d71ae5a4SJacob Faibussowitsch static void f0_j_re(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) 46d71ae5a4SJacob Faibussowitsch { 47e0eea495SMark PetscReal n_e = PetscRealPart(u[0]); 48e0eea495SMark if (dim == 2) { 49e0eea495SMark if (x[1] > RE_CUT || x[1] < -RE_CUT) { /* simply a cutoff for REs. v_|| > 3 v(T_e) */ 50e0eea495SMark *f0 = n_e * 2. * PETSC_PI * x[0] * x[1] * constants[0]; /* n * r * v_|| * q */ 51e0eea495SMark } else { 52e0eea495SMark *f0 = 0; 53e0eea495SMark } 54e0eea495SMark } else { 55e0eea495SMark if (x[2] > RE_CUT || x[2] < -RE_CUT) { /* simply a cutoff for REs. v_|| > 3 v(T_e) */ 56e0eea495SMark *f0 = n_e * x[2] * constants[0]; 57e0eea495SMark } else { 58e0eea495SMark *f0 = 0; 59e0eea495SMark } 60e0eea495SMark } 61e0eea495SMark } 62e0eea495SMark 63e0eea495SMark /* sum < v, u*v*q > */ 64d71ae5a4SJacob Faibussowitsch static void f0_jz_sum(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 q[], PetscScalar *f0) 65d71ae5a4SJacob Faibussowitsch { 66e0eea495SMark PetscInt ii; 67e0eea495SMark f0[0] = 0; 68e0eea495SMark if (dim == 2) { 698a6f2e61SMark Adams for (ii = 0; ii < Nf; ii++) f0[0] += u[ii] * 2. * PETSC_PI * x[0] * x[1] * q[ii]; /* n * r * v_|| * q * v_0 */ 70e0eea495SMark } else { 718a6f2e61SMark Adams for (ii = 0; ii < Nf; ii++) f0[0] += u[ii] * x[2] * q[ii]; /* n * v_|| * q * v_0 */ 72e0eea495SMark } 73e0eea495SMark } 74e0eea495SMark 75e0eea495SMark /* < v, n_e > */ 76d71ae5a4SJacob Faibussowitsch 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) 77d71ae5a4SJacob Faibussowitsch { 78e0eea495SMark PetscInt ii = (PetscInt)PetscRealPart(constants[0]); 79e0eea495SMark if (dim == 2) f0[0] = 2. * PETSC_PI * x[0] * u[ii]; 80ad540459SPierre Jolivet else f0[0] = u[ii]; 81e0eea495SMark } 82e0eea495SMark 83e0eea495SMark /* < v, n_e v_|| > */ 84d71ae5a4SJacob Faibussowitsch 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) 85d71ae5a4SJacob Faibussowitsch { 86e0eea495SMark PetscInt ii = (PetscInt)PetscRealPart(constants[0]); 87e0eea495SMark if (dim == 2) f0[0] = u[ii] * 2. * PETSC_PI * x[0] * x[1]; /* n r v_|| */ 88ad540459SPierre Jolivet else f0[0] = u[ii] * x[2]; /* n v_|| */ 89e0eea495SMark } 90e0eea495SMark 91e0eea495SMark /* < v, n_e (v-shift) > */ 92d71ae5a4SJacob Faibussowitsch static void f0_ve_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) 93d71ae5a4SJacob Faibussowitsch { 948a6f2e61SMark Adams PetscReal vz = numConstants > 0 ? PetscRealPart(constants[0]) : 0; 95e0eea495SMark 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 */ 96d71ae5a4SJacob Faibussowitsch else { 97d71ae5a4SJacob Faibussowitsch *f0 = u[0] * PetscSqrtReal(x[0] * x[0] + x[1] * x[1] + (x[2] - vz) * (x[2] - vz)); /* n v */ 98d71ae5a4SJacob Faibussowitsch } 99e0eea495SMark } 100e0eea495SMark 101e0eea495SMark /* CalculateE - Calculate the electric field */ 102e0eea495SMark /* T -- Electron temperature */ 103e0eea495SMark /* n -- Electron density */ 104e0eea495SMark /* lnLambda -- */ 105e0eea495SMark /* eps0 -- */ 106e0eea495SMark /* E -- output E, input \hat E */ 107d71ae5a4SJacob Faibussowitsch static PetscReal CalculateE(PetscReal Tev, PetscReal n, PetscReal lnLambda, PetscReal eps0, PetscReal *E) 108d71ae5a4SJacob Faibussowitsch { 109e0eea495SMark PetscReal c, e, m; 110e0eea495SMark 111e0eea495SMark PetscFunctionBegin; 112cefb98e8SMark Adams c = 299792458.0; 113e0eea495SMark e = 1.602176e-19; 114e0eea495SMark m = 9.10938e-31; 115e0eea495SMark if (1) { 116e0eea495SMark double Ec, Ehat = *E, betath = PetscSqrtReal(2 * Tev * e / (m * c * c)), j0 = Ehat * 7 / (PetscSqrtReal(2) * 2) * PetscPowReal(betath, 3) * n * e * c; 117e0eea495SMark Ec = n * lnLambda * PetscPowReal(e, 3) / (4 * PETSC_PI * PetscPowReal(eps0, 2) * m * c * c); 118e0eea495SMark *E = Ec; 1193ba16761SJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "CalculateE j0=%g Ec = %g\n", j0, Ec)); 120e0eea495SMark } else { 121e0eea495SMark PetscReal Ed, vth; 122e0eea495SMark vth = PetscSqrtReal(8 * Tev * e / (m * PETSC_PI)); 123e0eea495SMark Ed = n * lnLambda * PetscPowReal(e, 3) / (4 * PETSC_PI * PetscPowReal(eps0, 2) * m * vth * vth); 124e0eea495SMark *E = Ed; 125e0eea495SMark } 1263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 127e0eea495SMark } 128e0eea495SMark 129d71ae5a4SJacob Faibussowitsch static PetscReal Spitzer(PetscReal m_e, PetscReal e, PetscReal Z, PetscReal epsilon0, PetscReal lnLam, PetscReal kTe_joules) 130d71ae5a4SJacob Faibussowitsch { 131e0eea495SMark PetscReal Fz = (1 + 1.198 * Z + 0.222 * Z * Z) / (1 + 2.966 * Z + 0.753 * Z * Z), eta; 132e0eea495SMark 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); 133e0eea495SMark return eta; 134e0eea495SMark } 135e0eea495SMark 136e0eea495SMark /* */ 137d71ae5a4SJacob Faibussowitsch static PetscErrorCode testNone(TS ts, Vec X, PetscInt stepi, PetscReal time, PetscBool islast, LandauCtx *ctx, REctx *rectx) 138d71ae5a4SJacob Faibussowitsch { 139e0eea495SMark PetscFunctionBeginUser; 1403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 141e0eea495SMark } 142e0eea495SMark 143e0eea495SMark /* */ 144d71ae5a4SJacob Faibussowitsch static PetscErrorCode testSpitzer(TS ts, Vec X, PetscInt stepi, PetscReal time, PetscBool islast, LandauCtx *ctx, REctx *rectx) 145d71ae5a4SJacob Faibussowitsch { 146f37ccb5fSMark Adams PetscInt ii, nDMs; 147e0eea495SMark PetscDS prob; 148a587d139SMark static PetscReal old_ratio = 1e10; 149a587d139SMark TSConvergedReason reason; 150a587d139SMark PetscReal J, J_re, spit_eta, Te_kev = 0, E, ratio, Z, n_e, v, v2; 1518a6f2e61SMark Adams PetscScalar user[2] = {0., ctx->charges[0]}, q[LANDAU_MAX_SPECIES], tt[LANDAU_MAX_SPECIES], vz; 152930e68a5SMark Adams PetscReal dt; 1538a6f2e61SMark Adams DM pack, plexe = ctx->plex[0], plexi = (ctx->num_grids == 1) ? NULL : ctx->plex[1]; 154f37ccb5fSMark Adams Vec *XsubArray; 155e0eea495SMark 156e0eea495SMark PetscFunctionBeginUser; 15763a3b9bcSJacob Faibussowitsch PetscCheck(ctx->num_species == 2, PETSC_COMM_SELF, PETSC_ERR_PLIB, "ctx->num_species %" PetscInt_FMT " != 2", ctx->num_species); 1589566063dSJacob Faibussowitsch PetscCall(VecGetDM(X, &pack)); 1593c633725SBarry Smith PetscCheck(pack, PETSC_COMM_SELF, PETSC_ERR_PLIB, "no DM"); 1609566063dSJacob Faibussowitsch PetscCall(DMCompositeGetNumberDM(pack, &nDMs)); 16163a3b9bcSJacob Faibussowitsch PetscCheck(nDMs == ctx->num_grids * ctx->batch_sz, PETSC_COMM_SELF, PETSC_ERR_PLIB, "nDMs != ctx->num_grids*ctx->batch_sz %" PetscInt_FMT " != %" PetscInt_FMT, nDMs, ctx->num_grids * ctx->batch_sz); 1629566063dSJacob Faibussowitsch PetscCall(PetscMalloc(sizeof(*XsubArray) * nDMs, &XsubArray)); 1639566063dSJacob Faibussowitsch PetscCall(DMCompositeGetAccessArray(pack, X, nDMs, NULL, XsubArray)); // read only 1649566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts, &dt)); 1658a6f2e61SMark Adams /* get current for each grid */ 1668a6f2e61SMark Adams for (ii = 0; ii < ctx->num_species; ii++) q[ii] = ctx->charges[ii]; 1679566063dSJacob Faibussowitsch PetscCall(DMGetDS(plexe, &prob)); 1689566063dSJacob Faibussowitsch PetscCall(PetscDSSetConstants(prob, 2, &q[0])); 1699566063dSJacob Faibussowitsch PetscCall(PetscDSSetObjective(prob, 0, &f0_jz_sum)); 1709566063dSJacob Faibussowitsch PetscCall(DMPlexComputeIntegralFEM(plexe, XsubArray[LAND_PACK_IDX(ctx->batch_view_idx, 0)], tt, NULL)); 171e0eea495SMark J = -ctx->n_0 * ctx->v_0 * PetscRealPart(tt[0]); 1728fdabdddSMark Adams if (plexi) { // add first (only) ion 1739566063dSJacob Faibussowitsch PetscCall(DMGetDS(plexi, &prob)); 1749566063dSJacob Faibussowitsch PetscCall(PetscDSSetConstants(prob, 1, &q[1])); 1759566063dSJacob Faibussowitsch PetscCall(PetscDSSetObjective(prob, 0, &f0_jz_sum)); 1769566063dSJacob Faibussowitsch PetscCall(DMPlexComputeIntegralFEM(plexi, XsubArray[LAND_PACK_IDX(ctx->batch_view_idx, 1)], tt, NULL)); 1778a6f2e61SMark Adams J += -ctx->n_0 * ctx->v_0 * PetscRealPart(tt[0]); 1788a6f2e61SMark Adams } 179a587d139SMark /* get N_e */ 1809566063dSJacob Faibussowitsch PetscCall(DMGetDS(plexe, &prob)); 1819566063dSJacob Faibussowitsch PetscCall(PetscDSSetConstants(prob, 1, user)); 1829566063dSJacob Faibussowitsch PetscCall(PetscDSSetObjective(prob, 0, &f0_n)); 1839566063dSJacob Faibussowitsch PetscCall(DMPlexComputeIntegralFEM(plexe, XsubArray[LAND_PACK_IDX(ctx->batch_view_idx, 0)], tt, NULL)); 184a587d139SMark n_e = PetscRealPart(tt[0]) * ctx->n_0; 185a587d139SMark /* Z */ 186a587d139SMark Z = -ctx->charges[1] / ctx->charges[0]; 187a587d139SMark /* remove drift */ 1888a6f2e61SMark Adams if (0) { 1898a6f2e61SMark Adams user[0] = 0; // electrons 1909566063dSJacob Faibussowitsch PetscCall(DMGetDS(plexe, &prob)); 1919566063dSJacob Faibussowitsch PetscCall(PetscDSSetConstants(prob, 1, user)); 1929566063dSJacob Faibussowitsch PetscCall(PetscDSSetObjective(prob, 0, &f0_vz)); 1939566063dSJacob Faibussowitsch PetscCall(DMPlexComputeIntegralFEM(plexe, XsubArray[LAND_PACK_IDX(ctx->batch_view_idx, 0)], tt, NULL)); 194a587d139SMark vz = ctx->n_0 * PetscRealPart(tt[0]) / n_e; /* non-dimensional */ 195a587d139SMark } else vz = 0; 196a587d139SMark /* thermal velocity */ 1979566063dSJacob Faibussowitsch PetscCall(DMGetDS(plexe, &prob)); 1989566063dSJacob Faibussowitsch PetscCall(PetscDSSetConstants(prob, 1, &vz)); 1999566063dSJacob Faibussowitsch PetscCall(PetscDSSetObjective(prob, 0, &f0_ve_shift)); 2009566063dSJacob Faibussowitsch PetscCall(DMPlexComputeIntegralFEM(plexe, XsubArray[LAND_PACK_IDX(ctx->batch_view_idx, 0)], tt, NULL)); 201a587d139SMark v = ctx->n_0 * ctx->v_0 * PetscRealPart(tt[0]) / n_e; /* remove number density to get velocity */ 202a587d139SMark v2 = PetscSqr(v); /* use real space: m^2 / s^2 */ 203a587d139SMark Te_kev = (v2 * ctx->masses[0] * PETSC_PI / 8) * kev_joul; /* temperature in kev */ 204ae6bf4a8Smarkadams4 spit_eta = Spitzer(ctx->masses[0], -ctx->charges[0], Z, ctx->epsilon0, ctx->lambdas[0][1], Te_kev / kev_joul); /* kev --> J (kT) */ 2058a6f2e61SMark Adams if (0) { 2069566063dSJacob Faibussowitsch PetscCall(DMGetDS(plexe, &prob)); 2079566063dSJacob Faibussowitsch PetscCall(PetscDSSetConstants(prob, 1, q)); 2089566063dSJacob Faibussowitsch PetscCall(PetscDSSetObjective(prob, 0, &f0_j_re)); 2099566063dSJacob Faibussowitsch PetscCall(DMPlexComputeIntegralFEM(plexe, XsubArray[LAND_PACK_IDX(ctx->batch_view_idx, 0)], tt, NULL)); 210a587d139SMark } else tt[0] = 0; 211e0eea495SMark J_re = -ctx->n_0 * ctx->v_0 * PetscRealPart(tt[0]); 2129566063dSJacob Faibussowitsch PetscCall(DMCompositeRestoreAccessArray(pack, X, nDMs, NULL, XsubArray)); // read only 2139566063dSJacob Faibussowitsch PetscCall(PetscFree(XsubArray)); 214a587d139SMark 21590d163c9SMark Adams if (rectx->use_spitzer_eta) { 21690d163c9SMark Adams E = ctx->Ez = spit_eta * (rectx->j - J_re); 21790d163c9SMark Adams } else { 21890d163c9SMark Adams E = ctx->Ez; /* keep real E */ 21990d163c9SMark Adams rectx->j = J; /* cache */ 22090d163c9SMark Adams } 22190d163c9SMark Adams 222e0eea495SMark ratio = E / J / spit_eta; 223*f4f49eeaSPierre Jolivet if (stepi > 10 && !rectx->use_spitzer_eta && (old_ratio - ratio < 1.e-6)) { 22490d163c9SMark Adams rectx->pulse_start = time + 0.98 * dt; 225a587d139SMark rectx->use_spitzer_eta = PETSC_TRUE; 226a587d139SMark } 2279566063dSJacob Faibussowitsch PetscCall(TSGetConvergedReason(ts, &reason)); 2289566063dSJacob Faibussowitsch PetscCall(TSGetConvergedReason(ts, &reason)); 229*f4f49eeaSPierre Jolivet if (rectx->plotting || stepi == 0 || reason || rectx->pulse_start == time + 0.98 * dt) { 2309371c9d4SSatish Balay PetscCall(PetscPrintf(ctx->comm, "testSpitzer: %4" PetscInt_FMT ") time=%11.4e n_e= %10.3e E= %10.3e J= %10.3e J_re= %10.3e %.3g%% Te_kev= %10.3e Z_eff=%g E/J to eta ratio= %g (diff=%g) %s %s spit_eta=%g\n", stepi, (double)time, 2319371c9d4SSatish Balay (double)(n_e / ctx->n_0), (double)ctx->Ez, (double)J, (double)J_re, (double)(100 * J_re / J), (double)Te_kev, (double)Z, (double)ratio, (double)(old_ratio - ratio), rectx->use_spitzer_eta ? "using Spitzer eta*J E" : "constant E", rectx->pulse_start != time + 0.98 * dt ? "normal" : "transition", (double)spit_eta)); 23263a3b9bcSJacob Faibussowitsch PetscCheck(rectx->pulse_start != (time + 0.98 * dt), PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Spitzer complete ratio=%g", (double)ratio); 233930e68a5SMark Adams } 234e0eea495SMark old_ratio = ratio; 2353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 236e0eea495SMark } 237e0eea495SMark 238e0eea495SMark static const double ppp = 2; 239d71ae5a4SJacob Faibussowitsch static void f0_0_diff_lp(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) 240d71ae5a4SJacob Faibussowitsch { 241e0eea495SMark LandauCtx *ctx = (LandauCtx *)constants; 242e0eea495SMark REctx *rectx = (REctx *)ctx->data; 243e0eea495SMark PetscInt ii = rectx->idx, i; 244e0eea495SMark const PetscReal kT_m = ctx->k * ctx->thermal_temps[ii] / ctx->masses[ii]; /* kT/m */ 245e0eea495SMark const PetscReal n = ctx->n[ii]; 246e0eea495SMark PetscReal diff, f_maxwell, v2 = 0, theta = 2 * kT_m / (ctx->v_0 * ctx->v_0); /* theta = 2kT/mc^2 */ 247e0eea495SMark for (i = 0; i < dim; ++i) v2 += x[i] * x[i]; 248e0eea495SMark f_maxwell = n * PetscPowReal(PETSC_PI * theta, -1.5) * (PetscExpReal(-v2 / theta)); 249e0eea495SMark diff = 2. * PETSC_PI * x[0] * (PetscRealPart(u[ii]) - f_maxwell); 250e0eea495SMark f0[0] = PetscPowReal(diff, ppp); 251e0eea495SMark } 252d71ae5a4SJacob Faibussowitsch static void f0_0_maxwellian_lp(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) 253d71ae5a4SJacob Faibussowitsch { 254e0eea495SMark LandauCtx *ctx = (LandauCtx *)constants; 255e0eea495SMark REctx *rectx = (REctx *)ctx->data; 256e0eea495SMark PetscInt ii = rectx->idx, i; 257e0eea495SMark const PetscReal kT_m = ctx->k * ctx->thermal_temps[ii] / ctx->masses[ii]; /* kT/m */ 258e0eea495SMark const PetscReal n = ctx->n[ii]; 259e0eea495SMark PetscReal f_maxwell, v2 = 0, theta = 2 * kT_m / (ctx->v_0 * ctx->v_0); /* theta = 2kT/mc^2 */ 260e0eea495SMark for (i = 0; i < dim; ++i) v2 += x[i] * x[i]; 261e0eea495SMark f_maxwell = 2. * PETSC_PI * x[0] * n * PetscPowReal(PETSC_PI * theta, -1.5) * (PetscExpReal(-v2 / theta)); 262e0eea495SMark f0[0] = PetscPowReal(f_maxwell, ppp); 263e0eea495SMark } 264e0eea495SMark 265e0eea495SMark /* */ 266d71ae5a4SJacob Faibussowitsch static PetscErrorCode testStable(TS ts, Vec X, PetscInt stepi, PetscReal time, PetscBool islast, LandauCtx *ctx, REctx *rectx) 267d71ae5a4SJacob Faibussowitsch { 268e0eea495SMark PetscDS prob; 269e0eea495SMark Vec X2; 270e0eea495SMark PetscReal ediff, idiff = 0, lpm0, lpm1 = 1; 271e0eea495SMark PetscScalar tt[LANDAU_MAX_SPECIES]; 2728a6f2e61SMark Adams DM dm, plex = ctx->plex[0]; 273e0eea495SMark 274e0eea495SMark PetscFunctionBeginUser; 2759566063dSJacob Faibussowitsch PetscCall(VecGetDM(X, &dm)); 2769566063dSJacob Faibussowitsch PetscCall(DMGetDS(plex, &prob)); 2779566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X, &X2)); 2789566063dSJacob Faibussowitsch PetscCall(VecCopy(X, X2)); 279e0eea495SMark if (!rectx->X_0) { 2809566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X, &rectx->X_0)); 2819566063dSJacob Faibussowitsch PetscCall(VecCopy(X, rectx->X_0)); 282e0eea495SMark } 2839566063dSJacob Faibussowitsch PetscCall(VecAXPY(X, -1.0, rectx->X_0)); 2849566063dSJacob Faibussowitsch PetscCall(PetscDSSetConstants(prob, sizeof(LandauCtx) / sizeof(PetscScalar), (PetscScalar *)ctx)); 285e0eea495SMark rectx->idx = 0; 2869566063dSJacob Faibussowitsch PetscCall(PetscDSSetObjective(prob, 0, &f0_0_diff_lp)); 2879566063dSJacob Faibussowitsch PetscCall(DMPlexComputeIntegralFEM(plex, X2, tt, NULL)); 288e0eea495SMark ediff = PetscPowReal(PetscRealPart(tt[0]), 1. / ppp); 2899566063dSJacob Faibussowitsch PetscCall(PetscDSSetObjective(prob, 0, &f0_0_maxwellian_lp)); 2909566063dSJacob Faibussowitsch PetscCall(DMPlexComputeIntegralFEM(plex, X2, tt, NULL)); 291e0eea495SMark lpm0 = PetscPowReal(PetscRealPart(tt[0]), 1. / ppp); 292e0eea495SMark if (ctx->num_species > 1) { 293e0eea495SMark rectx->idx = 1; 2949566063dSJacob Faibussowitsch PetscCall(PetscDSSetObjective(prob, 0, &f0_0_diff_lp)); 2959566063dSJacob Faibussowitsch PetscCall(DMPlexComputeIntegralFEM(plex, X2, tt, NULL)); 296e0eea495SMark idiff = PetscPowReal(PetscRealPart(tt[0]), 1. / ppp); 2979566063dSJacob Faibussowitsch PetscCall(PetscDSSetObjective(prob, 0, &f0_0_maxwellian_lp)); 2989566063dSJacob Faibussowitsch PetscCall(DMPlexComputeIntegralFEM(plex, X2, tt, NULL)); 299e0eea495SMark lpm1 = PetscPowReal(PetscRealPart(tt[0]), 1. / ppp); 300e0eea495SMark } 30163a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s %" PetscInt_FMT ") time=%10.3e n-%d norm electrons/max=%20.13e ions/max=%20.13e\n", "----", stepi, (double)time, (int)ppp, (double)(ediff / lpm0), (double)(idiff / lpm1))); 302e0eea495SMark /* view */ 3039566063dSJacob Faibussowitsch PetscCall(VecCopy(X2, X)); 3049566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X2)); 305e0eea495SMark if (islast) { 3069566063dSJacob Faibussowitsch PetscCall(VecDestroy(&rectx->X_0)); 307e0eea495SMark rectx->X_0 = NULL; 308e0eea495SMark } 3093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 310e0eea495SMark } 311e0eea495SMark 312d71ae5a4SJacob Faibussowitsch static PetscErrorCode EInduction(Vec X, Vec X_t, PetscInt step, PetscReal time, LandauCtx *ctx, PetscReal *a_E) 313d71ae5a4SJacob Faibussowitsch { 314e0eea495SMark REctx *rectx = (REctx *)ctx->data; 315e0eea495SMark PetscInt ii; 316e0eea495SMark DM dm, plex; 3178a6f2e61SMark Adams PetscScalar tt[LANDAU_MAX_SPECIES], qv0[LANDAU_MAX_SPECIES]; 318e0eea495SMark PetscReal dJ_dt; 319e0eea495SMark PetscDS prob; 320e0eea495SMark 321e0eea495SMark PetscFunctionBeginUser; 3228a6f2e61SMark Adams for (ii = 0; ii < ctx->num_species; ii++) qv0[ii] = ctx->charges[ii] * ctx->v_0; 3239566063dSJacob Faibussowitsch PetscCall(VecGetDM(X, &dm)); 3249566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &prob)); 3259566063dSJacob Faibussowitsch PetscCall(DMConvert(dm, DMPLEX, &plex)); 326e0eea495SMark /* get d current / dt */ 3279566063dSJacob Faibussowitsch PetscCall(PetscDSSetConstants(prob, ctx->num_species, qv0)); 3289566063dSJacob Faibussowitsch PetscCall(PetscDSSetObjective(prob, 0, &f0_jz_sum)); 3293c633725SBarry Smith PetscCheck(X_t, PETSC_COMM_SELF, PETSC_ERR_PLIB, "X_t"); 3309566063dSJacob Faibussowitsch PetscCall(DMPlexComputeIntegralFEM(plex, X_t, tt, NULL)); 3318a6f2e61SMark Adams dJ_dt = -ctx->n_0 * PetscRealPart(tt[0]) / ctx->t_0; 332e0eea495SMark /* E induction */ 333e0eea495SMark *a_E = -rectx->L * dJ_dt + rectx->Ez_initial; 3349566063dSJacob Faibussowitsch PetscCall(DMDestroy(&plex)); 3353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 336e0eea495SMark } 337e0eea495SMark 338d71ae5a4SJacob Faibussowitsch static PetscErrorCode EConstant(Vec X, Vec X_t, PetscInt step, PetscReal time, LandauCtx *ctx, PetscReal *a_E) 339d71ae5a4SJacob Faibussowitsch { 340e0eea495SMark PetscFunctionBeginUser; 341e0eea495SMark *a_E = ctx->Ez; 3423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 343e0eea495SMark } 344e0eea495SMark 345d71ae5a4SJacob Faibussowitsch static PetscErrorCode ENone(Vec X, Vec X_t, PetscInt step, PetscReal time, LandauCtx *ctx, PetscReal *a_E) 346d71ae5a4SJacob Faibussowitsch { 347e0eea495SMark PetscFunctionBeginUser; 348e0eea495SMark *a_E = 0; 3493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 350e0eea495SMark } 351e0eea495SMark 352e0eea495SMark /* ------------------------------------------------------------------- */ 353e0eea495SMark /* 354e0eea495SMark FormSource - Evaluates source terms F(t). 355e0eea495SMark 356e0eea495SMark Input Parameters: 357e0eea495SMark . ts - the TS context 358e0eea495SMark . time - 359e0eea495SMark . X_dummmy - input vector 360e0eea495SMark . dummy - optional user-defined context, as set by SNESSetFunction() 361e0eea495SMark 362e0eea495SMark Output Parameter: 363e0eea495SMark . F - function vector 364e0eea495SMark */ 365d71ae5a4SJacob Faibussowitsch static PetscErrorCode FormSource(TS ts, PetscReal ftime, Vec X_dummmy, Vec F, void *dummy) 366d71ae5a4SJacob Faibussowitsch { 367e0eea495SMark PetscReal new_imp_rate; 368e0eea495SMark LandauCtx *ctx; 3698a6f2e61SMark Adams DM pack; 370e0eea495SMark REctx *rectx; 371e0eea495SMark 372e0eea495SMark PetscFunctionBeginUser; 3739566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &pack)); 3749566063dSJacob Faibussowitsch PetscCall(DMGetApplicationContext(pack, &ctx)); 375e0eea495SMark rectx = (REctx *)ctx->data; 376e0eea495SMark /* check for impurities */ 3779566063dSJacob Faibussowitsch PetscCall(rectx->impuritySrcRate(ftime, &new_imp_rate, ctx)); 378e0eea495SMark if (new_imp_rate != 0) { 379e0eea495SMark if (new_imp_rate != rectx->current_rate) { 380e0eea495SMark PetscInt ii; 381e0eea495SMark PetscReal dne_dt, dni_dt, tilda_ns[LANDAU_MAX_SPECIES], temps[LANDAU_MAX_SPECIES]; 3828fdabdddSMark Adams Vec globFarray[LANDAU_MAX_GRIDS * LANDAU_MAX_BATCH_SZ]; 383e0eea495SMark rectx->current_rate = new_imp_rate; 384e0eea495SMark for (ii = 1; ii < LANDAU_MAX_SPECIES; ii++) tilda_ns[ii] = 0; 385e0eea495SMark for (ii = 1; ii < LANDAU_MAX_SPECIES; ii++) temps[ii] = 1; 3868a6f2e61SMark Adams dni_dt = new_imp_rate /* *ctx->t_0 */; /* fully ionized immediately, no normalize, stay in non-dim */ 3878a6f2e61SMark Adams dne_dt = new_imp_rate * rectx->Ne_ion /* *ctx->t_0 */; 3889371c9d4SSatish Balay tilda_ns[0] = dne_dt; 3899371c9d4SSatish Balay tilda_ns[rectx->imp_idx] = dni_dt; 3909371c9d4SSatish Balay temps[0] = rectx->T_cold; 3919371c9d4SSatish Balay temps[rectx->imp_idx] = rectx->T_cold; 39263a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ctx->plex[0], "\tHave new_imp_rate= %10.3e time= %10.3e de/dt= %10.3e di/dt= %10.3e ***\n", (double)new_imp_rate, (double)ftime, (double)dne_dt, (double)dni_dt)); 3939566063dSJacob Faibussowitsch PetscCall(DMCompositeGetAccessArray(pack, F, ctx->num_grids * ctx->batch_sz, NULL, globFarray)); 3948a6f2e61SMark Adams for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { 395e0eea495SMark /* add it */ 396c3e4dd79SMark Adams PetscCall(DMPlexLandauAddMaxwellians(ctx->plex[grid], globFarray[LAND_PACK_IDX(0, grid)], ftime, temps, tilda_ns, grid, 0, 1, ctx)); 397e0eea495SMark } 3988fdabdddSMark Adams // Does DMCompositeRestoreAccessArray copy the data back? (no) 3999566063dSJacob Faibussowitsch PetscCall(DMCompositeRestoreAccessArray(pack, F, ctx->num_grids * ctx->batch_sz, NULL, globFarray)); 400e0eea495SMark } 401e0eea495SMark } else { 4029566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F)); 403e0eea495SMark rectx->current_rate = 0; 404e0eea495SMark } 4053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 406e0eea495SMark } 407419c2fb8Smarkadams4 408d71ae5a4SJacob Faibussowitsch PetscErrorCode Monitor(TS ts, PetscInt stepi, PetscReal time, Vec X, void *actx) 409d71ae5a4SJacob Faibussowitsch { 410e0eea495SMark LandauCtx *ctx = (LandauCtx *)actx; /* user-defined application context */ 411e0eea495SMark REctx *rectx = (REctx *)ctx->data; 412419c2fb8Smarkadams4 DM pack = NULL; 4138fdabdddSMark Adams Vec globXArray[LANDAU_MAX_GRIDS * LANDAU_MAX_BATCH_SZ]; 414e0eea495SMark TSConvergedReason reason; 415419c2fb8Smarkadams4 416e0eea495SMark PetscFunctionBeginUser; 417419c2fb8Smarkadams4 PetscCall(TSGetConvergedReason(ts, &reason)); 418419c2fb8Smarkadams4 if (rectx->grid_view_idx != -1 || (reason && ctx->verbose > 3)) { 4199566063dSJacob Faibussowitsch PetscCall(VecGetDM(X, &pack)); 4209566063dSJacob Faibussowitsch PetscCall(DMCompositeGetAccessArray(pack, X, ctx->num_grids * ctx->batch_sz, NULL, globXArray)); 421419c2fb8Smarkadams4 } 422e0eea495SMark if (stepi > rectx->plotStep && rectx->plotting) { 423e0eea495SMark rectx->plotting = PETSC_FALSE; /* was doing diagnostics, now done */ 424e0eea495SMark rectx->plotIdx++; 425e0eea495SMark } 426e0eea495SMark /* view */ 427e0eea495SMark if (time / rectx->plotDt >= (PetscReal)rectx->plotIdx || reason) { 428c3e4dd79SMark Adams if ((reason || stepi == 0 || rectx->plotIdx % rectx->print_period == 0) && ctx->verbose > 1) { 429e0eea495SMark /* print norms */ 4309566063dSJacob Faibussowitsch PetscCall(DMPlexLandauPrintNorms(X, stepi)); 431a587d139SMark } 432930e68a5SMark Adams if (!rectx->plotting) { /* first step of possible backtracks */ 433930e68a5SMark Adams rectx->plotting = PETSC_TRUE; 434930e68a5SMark Adams /* diagnostics + change E field with Sptizer (not just a monitor) */ 4359566063dSJacob Faibussowitsch PetscCall(rectx->test(ts, X, stepi, time, reason ? PETSC_TRUE : PETSC_FALSE, ctx, rectx)); 436930e68a5SMark Adams } else { 4373ba16761SJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\t\t ERROR SKIP test spit ------\n")); 438930e68a5SMark Adams rectx->plotting = PETSC_TRUE; 439930e68a5SMark Adams } 440419c2fb8Smarkadams4 if (rectx->grid_view_idx != -1) { 4419566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)globXArray[LAND_PACK_IDX(ctx->batch_view_idx, rectx->grid_view_idx)], rectx->grid_view_idx == 0 ? "ue" : "ui")); 442930e68a5SMark Adams /* view, overwrite step when back tracked */ 443e5d13be4Smarkadams4 PetscCall(DMSetOutputSequenceNumber(ctx->plex[rectx->grid_view_idx], rectx->plotIdx, time * ctx->t_0)); 444419c2fb8Smarkadams4 PetscCall(VecViewFromOptions(globXArray[LAND_PACK_IDX(ctx->batch_view_idx, rectx->grid_view_idx)], NULL, "-ex2_vec_view")); 445419c2fb8Smarkadams4 } 446e0eea495SMark rectx->plotStep = stepi; 447930e68a5SMark Adams } else { 4483ba16761SJacob Faibussowitsch if (rectx->plotting) PetscCall(PetscPrintf(PETSC_COMM_WORLD, " ERROR rectx->plotting=%s step %" PetscInt_FMT "\n", PetscBools[rectx->plotting], stepi)); 449930e68a5SMark Adams /* diagnostics + change E field with Sptizer (not just a monitor) - can we lag this? */ 4509566063dSJacob Faibussowitsch PetscCall(rectx->test(ts, X, stepi, time, reason ? PETSC_TRUE : PETSC_FALSE, ctx, rectx)); 451e0eea495SMark } 452f37ccb5fSMark Adams /* parallel check that only works of all batches are identical */ 453ae6bf4a8Smarkadams4 if (reason && ctx->verbose > 3 && ctx->batch_sz > 1) { 454e0eea495SMark PetscReal val, rval; 455e0eea495SMark PetscMPIInt rank; 4569566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 457f37ccb5fSMark Adams for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { 458f37ccb5fSMark Adams PetscInt nerrors = 0; 459f37ccb5fSMark Adams for (PetscInt i = 0; i < ctx->batch_sz; i++) { 4609566063dSJacob Faibussowitsch PetscCall(VecNorm(globXArray[LAND_PACK_IDX(i, grid)], NORM_2, &val)); 461f37ccb5fSMark Adams if (i == 0) rval = val; 462f37ccb5fSMark Adams else if ((val = PetscAbs(val - rval) / rval) > 1000 * PETSC_MACHINE_EPSILON) { 46363a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " [%d] Warning %" PetscInt_FMT ".%" PetscInt_FMT ") diff = %2.15e\n", rank, grid, i, (double)val)); 464f37ccb5fSMark Adams nerrors++; 465f37ccb5fSMark Adams } 466f37ccb5fSMark Adams } 467f37ccb5fSMark Adams if (nerrors) { 46863a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " ***** [%d] ERROR max %" PetscInt_FMT " errors\n", rank, nerrors)); 469e0eea495SMark } else { 47063a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "[%d] %" PetscInt_FMT ") batch consistency check OK\n", rank, grid)); 471f37ccb5fSMark Adams } 472e0eea495SMark } 473e0eea495SMark } 474e0eea495SMark rectx->idx = 0; 47548a46eb9SPierre Jolivet if (rectx->grid_view_idx != -1 || (reason && ctx->verbose > 3)) PetscCall(DMCompositeRestoreAccessArray(pack, X, ctx->num_grids * ctx->batch_sz, NULL, globXArray)); 4763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 477e0eea495SMark } 478e0eea495SMark 479d71ae5a4SJacob Faibussowitsch PetscErrorCode PreStep(TS ts) 480d71ae5a4SJacob Faibussowitsch { 481e0eea495SMark LandauCtx *ctx; 482e0eea495SMark REctx *rectx; 483e0eea495SMark DM dm; 484e0eea495SMark PetscInt stepi; 485e0eea495SMark PetscReal time; 486e0eea495SMark Vec X; 487e0eea495SMark 488e0eea495SMark PetscFunctionBeginUser; 489930e68a5SMark Adams /* not used */ 4909566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 4919566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts, &time)); 4929566063dSJacob Faibussowitsch PetscCall(TSGetSolution(ts, &X)); 4939566063dSJacob Faibussowitsch PetscCall(DMGetApplicationContext(dm, &ctx)); 494e0eea495SMark rectx = (REctx *)ctx->data; 4959566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &stepi)); 496e0eea495SMark /* update E */ 4979566063dSJacob Faibussowitsch PetscCall(rectx->E(X, NULL, stepi, time, ctx, &ctx->Ez)); 4983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 499e0eea495SMark } 500e0eea495SMark 501e0eea495SMark /* 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) */ 502d71ae5a4SJacob Faibussowitsch static PetscErrorCode stepSrc(PetscReal time, PetscReal *rho, LandauCtx *ctx) 503d71ae5a4SJacob Faibussowitsch { 504e0eea495SMark REctx *rectx = (REctx *)ctx->data; 505e0eea495SMark 506e0eea495SMark PetscFunctionBeginUser; 507e0eea495SMark if (time >= rectx->pulse_start) *rho = rectx->pulse_rate; 508e0eea495SMark else *rho = 0.; 5093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 510e0eea495SMark } 511d71ae5a4SJacob Faibussowitsch static PetscErrorCode zeroSrc(PetscReal time, PetscReal *rho, LandauCtx *ctx) 512d71ae5a4SJacob Faibussowitsch { 513e0eea495SMark PetscFunctionBeginUser; 514e0eea495SMark *rho = 0.; 5153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 516e0eea495SMark } 517d71ae5a4SJacob Faibussowitsch static PetscErrorCode pulseSrc(PetscReal time, PetscReal *rho, LandauCtx *ctx) 518d71ae5a4SJacob Faibussowitsch { 519e0eea495SMark REctx *rectx = (REctx *)ctx->data; 520e0eea495SMark 521e0eea495SMark PetscFunctionBeginUser; 52208401ef6SPierre Jolivet PetscCheck(rectx->pulse_start != PETSC_MAX_REAL, PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "'-ex2_pulse_start_time X' must be used with '-ex2_impurity_source_type pulse'"); 523e0eea495SMark if (time < rectx->pulse_start || time > rectx->pulse_start + 3 * rectx->pulse_width) *rho = 0; 524a587d139SMark else { 525e0eea495SMark double x = PetscSinReal((time - rectx->pulse_start) / (3 * rectx->pulse_width) * 2 * PETSC_PI - PETSC_PI / 2) + 1; /* 0:2, integrates to 1.0 */ 526e0eea495SMark *rho = rectx->pulse_rate * x / (3 * rectx->pulse_width); 527a587d139SMark if (!rectx->use_spitzer_eta) rectx->use_spitzer_eta = PETSC_TRUE; /* use it next time */ 528e0eea495SMark } 5293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 530e0eea495SMark } 531e0eea495SMark 532e0eea495SMark #undef __FUNCT__ 533e0eea495SMark #define __FUNCT__ "ProcessREOptions" 534d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessREOptions(REctx *rectx, const LandauCtx *ctx, DM dm, const char prefix[]) 535d71ae5a4SJacob Faibussowitsch { 536e0eea495SMark PetscFunctionList plist = NULL, testlist = NULL, elist = NULL; 537e0eea495SMark char pname[256], testname[256], ename[256]; 538930e68a5SMark Adams DM dm_dummy; 539a587d139SMark PetscBool Connor_E = PETSC_FALSE; 540e0eea495SMark 541e0eea495SMark PetscFunctionBeginUser; 5429566063dSJacob Faibussowitsch PetscCall(DMCreate(PETSC_COMM_WORLD, &dm_dummy)); 543e0eea495SMark rectx->Ne_ion = 1; /* number of electrons given up by impurity ion */ 544e0eea495SMark rectx->T_cold = .005; /* kev */ 545e0eea495SMark rectx->ion_potential = 15; /* ev */ 546e0eea495SMark rectx->L = 2; 547e0eea495SMark rectx->X_0 = NULL; 548e0eea495SMark rectx->imp_idx = ctx->num_species - 1; /* default ionized impurity as last one */ 549a587d139SMark rectx->pulse_start = PETSC_MAX_REAL; 550e0eea495SMark rectx->pulse_width = 1; 551e0eea495SMark rectx->plotStep = PETSC_MAX_INT; 552e0eea495SMark rectx->pulse_rate = 1.e-1; 553e0eea495SMark rectx->current_rate = 0; 554e0eea495SMark rectx->plotIdx = 0; 555e0eea495SMark rectx->j = 0; 556e0eea495SMark rectx->plotDt = 1.0; 557930e68a5SMark Adams rectx->plotting = PETSC_FALSE; 558e0eea495SMark rectx->use_spitzer_eta = PETSC_FALSE; 559e0eea495SMark rectx->idx = 0; 5608a6f2e61SMark Adams rectx->print_period = 10; 561419c2fb8Smarkadams4 rectx->grid_view_idx = -1; // do not get if not needed 562e0eea495SMark /* Register the available impurity sources */ 5639566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&plist, "step", &stepSrc)); 5649566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&plist, "none", &zeroSrc)); 5659566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&plist, "pulse", &pulseSrc)); 566c6a7a370SJeremy L Thompson PetscCall(PetscStrncpy(pname, "none", sizeof(pname))); 5679566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&testlist, "none", &testNone)); 5689566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&testlist, "spitzer", &testSpitzer)); 5699566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&testlist, "stable", &testStable)); 570c6a7a370SJeremy L Thompson PetscCall(PetscStrncpy(testname, "none", sizeof(testname))); 5719566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&elist, "none", &ENone)); 5729566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&elist, "induction", &EInduction)); 5739566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&elist, "constant", &EConstant)); 574c6a7a370SJeremy L Thompson PetscCall(PetscStrncpy(ename, "constant", sizeof(ename))); 575e0eea495SMark 576d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_SELF, prefix, "Options for Runaway/seed electron model", "none"); 5779566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ex2_plot_dt", "Plotting interval", "ex2.c", rectx->plotDt, &rectx->plotDt, NULL)); 578e0eea495SMark if (rectx->plotDt < 0) rectx->plotDt = 1e30; 579e0eea495SMark if (rectx->plotDt == 0) rectx->plotDt = 1e-30; 5809566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ex2_print_period", "Plotting interval", "ex2.c", rectx->print_period, &rectx->print_period, NULL)); 5819566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ex2_grid_view_idx", "grid_view_idx", "ex2.c", rectx->grid_view_idx, &rectx->grid_view_idx, NULL)); 582419c2fb8Smarkadams4 PetscCheck(rectx->grid_view_idx < ctx->num_grids || rectx->grid_view_idx == -1, PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "rectx->grid_view_idx (%" PetscInt_FMT ") >= ctx->num_grids (%" PetscInt_FMT ")", rectx->imp_idx, ctx->num_grids); 5839566063dSJacob Faibussowitsch PetscCall(PetscOptionsFList("-ex2_impurity_source_type", "Name of impurity source to run", "", plist, pname, pname, sizeof(pname), NULL)); 5849566063dSJacob Faibussowitsch PetscCall(PetscOptionsFList("-ex2_test_type", "Name of test to run", "", testlist, testname, testname, sizeof(testname), NULL)); 5859566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ex2_impurity_index", "index of sink for impurities", "none", rectx->imp_idx, &rectx->imp_idx, NULL)); 586cad9d221SBarry Smith PetscCheck((rectx->imp_idx < ctx->num_species && rectx->imp_idx >= 1) || ctx->num_species <= 1, PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "index of sink for impurities ions is out of range (%" PetscInt_FMT "), must be > 0 && < NS", rectx->imp_idx); 5879566063dSJacob Faibussowitsch PetscCall(PetscOptionsFList("-ex2_e_field_type", "Electric field type", "", elist, ename, ename, sizeof(ename), NULL)); 588e0eea495SMark rectx->Ne_ion = -ctx->charges[rectx->imp_idx] / ctx->charges[0]; 5899566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ex2_t_cold", "Temperature of cold electron and ions after ionization in keV", "none", rectx->T_cold, &rectx->T_cold, NULL)); 5909566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ex2_pulse_start_time", "Time at which pulse happens for 'pulse' source", "none", rectx->pulse_start, &rectx->pulse_start, NULL)); 5919566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ex2_pulse_width_time", "Width of pulse 'pulse' source", "none", rectx->pulse_width, &rectx->pulse_width, NULL)); 5929566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ex2_pulse_rate", "Number density of pulse for 'pulse' source", "none", rectx->pulse_rate, &rectx->pulse_rate, NULL)); 593e0eea495SMark rectx->T_cold *= 1.16e7; /* convert to Kelvin */ 5949566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ex2_ion_potential", "Potential to ionize impurity (should be array) in ev", "none", rectx->ion_potential, &rectx->ion_potential, NULL)); 595da81f932SPierre Jolivet PetscCall(PetscOptionsReal("-ex2_inductance", "Inductance E field", "none", rectx->L, &rectx->L, NULL)); 5969566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ex2_connor_e_field_units", "Scale Ex but Connor-Hastie E_c", "none", Connor_E, &Connor_E, NULL)); 59763a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(dm_dummy, "Num electrons from ions=%g, T_cold=%10.3e, ion potential=%10.3e, E_z=%10.3e v_0=%10.3e\n", (double)rectx->Ne_ion, (double)rectx->T_cold, (double)rectx->ion_potential, (double)ctx->Ez, (double)ctx->v_0)); 598d0609cedSBarry Smith PetscOptionsEnd(); 599e0eea495SMark /* get impurity source rate function */ 6009566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(plist, pname, &rectx->impuritySrcRate)); 6013c633725SBarry Smith PetscCheck(rectx->impuritySrcRate, PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "No impurity source function found '%s'", pname); 6029566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(testlist, testname, &rectx->test)); 6033c633725SBarry Smith PetscCheck(rectx->test, PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "No test found '%s'", testname); 6049566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(elist, ename, &rectx->E)); 6053c633725SBarry Smith PetscCheck(rectx->E, PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "No E field function found '%s'", ename); 6069566063dSJacob Faibussowitsch PetscCall(PetscFunctionListDestroy(&plist)); 6079566063dSJacob Faibussowitsch PetscCall(PetscFunctionListDestroy(&testlist)); 6089566063dSJacob Faibussowitsch PetscCall(PetscFunctionListDestroy(&elist)); 609930e68a5SMark Adams 610a587d139SMark /* convert E from Connor-Hastie E_c units to real if doing Spitzer E */ 611a587d139SMark if (Connor_E) { 612e0eea495SMark PetscReal E = ctx->Ez, Tev = ctx->thermal_temps[0] * 8.621738e-5, n = ctx->n_0 * ctx->n[0]; 613ae6bf4a8Smarkadams4 CalculateE(Tev, n, ctx->lambdas[0][1], ctx->epsilon0, &E); 614e0eea495SMark ((LandauCtx *)ctx)->Ez *= E; 615e0eea495SMark } 6169566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm_dummy)); 6173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 618e0eea495SMark } 619e0eea495SMark 620d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 621d71ae5a4SJacob Faibussowitsch { 6228a6f2e61SMark Adams DM pack; 623419c2fb8Smarkadams4 Vec X; 624f37ccb5fSMark Adams PetscInt dim = 2, nDMs; 625e0eea495SMark TS ts; 626e0eea495SMark Mat J; 627e0eea495SMark PetscDS prob; 628e0eea495SMark LandauCtx *ctx; 629e0eea495SMark REctx *rectx; 6306c2b77d5SStefano Zampini PetscMPIInt rank; 631a587d139SMark PetscLogStage stage; 6326c2b77d5SStefano Zampini 633327415f7SBarry Smith PetscFunctionBeginUser; 6349566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 6359566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 636930e68a5SMark Adams if (rank) { /* turn off output stuff for duplicate runs */ 637419c2fb8Smarkadams4 PetscCall(PetscOptionsClearValue(NULL, "-ex2_dm_view")); 638419c2fb8Smarkadams4 PetscCall(PetscOptionsClearValue(NULL, "-ex2_vec_view")); 639419c2fb8Smarkadams4 PetscCall(PetscOptionsClearValue(NULL, "-ex2_vec_view_init")); 640419c2fb8Smarkadams4 PetscCall(PetscOptionsClearValue(NULL, "-ex2_dm_view_init")); 6419566063dSJacob Faibussowitsch PetscCall(PetscOptionsClearValue(NULL, "-info")); /* this does not work */ 642930e68a5SMark Adams } 6439566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-dim", &dim, NULL)); 644e0eea495SMark /* Create a mesh */ 6459566063dSJacob Faibussowitsch PetscCall(DMPlexLandauCreateVelocitySpace(PETSC_COMM_WORLD, dim, "", &X, &J, &pack)); 6469566063dSJacob Faibussowitsch PetscCall(DMCompositeGetNumberDM(pack, &nDMs)); 6479566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)J, "Jacobian")); 6489566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)X, "f")); 6499566063dSJacob Faibussowitsch PetscCall(DMGetApplicationContext(pack, &ctx)); 6509566063dSJacob Faibussowitsch PetscCall(DMSetUp(pack)); 651e0eea495SMark /* context */ 6529566063dSJacob Faibussowitsch PetscCall(PetscNew(&rectx)); 653e0eea495SMark ctx->data = rectx; 6549566063dSJacob Faibussowitsch PetscCall(ProcessREOptions(rectx, ctx, pack, "")); 6559566063dSJacob Faibussowitsch PetscCall(DMGetDS(pack, &prob)); 656419c2fb8Smarkadams4 if (rectx->grid_view_idx != -1) { 657419c2fb8Smarkadams4 Vec *XsubArray = NULL; 658419c2fb8Smarkadams4 PetscCall(PetscMalloc(sizeof(*XsubArray) * nDMs, &XsubArray)); 6599566063dSJacob Faibussowitsch PetscCall(DMCompositeGetAccessArray(pack, X, nDMs, NULL, XsubArray)); // read only 6609566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)XsubArray[LAND_PACK_IDX(ctx->batch_view_idx, rectx->grid_view_idx)], rectx->grid_view_idx == 0 ? "ue" : "ui")); 661e5d13be4Smarkadams4 PetscCall(DMSetOutputSequenceNumber(ctx->plex[rectx->grid_view_idx], 0, 0.0)); 662419c2fb8Smarkadams4 PetscCall(DMViewFromOptions(ctx->plex[rectx->grid_view_idx], NULL, "-ex2_dm_view")); 663419c2fb8Smarkadams4 PetscCall(DMViewFromOptions(ctx->plex[rectx->grid_view_idx], NULL, "-ex2_dm_view_init")); 664e5d13be4Smarkadams4 PetscCall(VecViewFromOptions(XsubArray[LAND_PACK_IDX(ctx->batch_view_idx, rectx->grid_view_idx)], NULL, "-ex2_vec_view")); // initial condition (monitor plots after step) 665419c2fb8Smarkadams4 PetscCall(VecViewFromOptions(XsubArray[LAND_PACK_IDX(ctx->batch_view_idx, rectx->grid_view_idx)], NULL, "-ex2_vec_view_init")); // initial condition (monitor plots after step) 6669566063dSJacob Faibussowitsch PetscCall(DMCompositeRestoreAccessArray(pack, X, nDMs, NULL, XsubArray)); // read only 6679566063dSJacob Faibussowitsch PetscCall(PetscFree(XsubArray)); 668419c2fb8Smarkadams4 } 669e0eea495SMark /* Create timestepping solver context */ 6709566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_SELF, &ts)); 6719566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts, pack)); 6729566063dSJacob Faibussowitsch PetscCall(TSSetIFunction(ts, NULL, DMPlexLandauIFunction, NULL)); 6739566063dSJacob Faibussowitsch PetscCall(TSSetIJacobian(ts, J, J, DMPlexLandauIJacobian, NULL)); 6749566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts, NULL, FormSource, NULL)); 6759566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 6769566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts, X)); 6779566063dSJacob Faibussowitsch PetscCall(TSSetApplicationContext(ts, ctx)); 6789566063dSJacob Faibussowitsch PetscCall(TSMonitorSet(ts, Monitor, ctx, NULL)); 6799566063dSJacob Faibussowitsch PetscCall(TSSetPreStep(ts, PreStep)); 680da81f932SPierre Jolivet rectx->Ez_initial = ctx->Ez; /* cache for induction calculation - applied E field */ 6818594ddcfSMark Adams if (1) { /* warm up an test just DMPlexLandauIJacobian */ 682e0eea495SMark Vec vec; 683a587d139SMark PetscInt nsteps; 684a587d139SMark PetscReal dt; 6859566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("Warmup", &stage)); 6869566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stage)); 6879566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X, &vec)); 6889566063dSJacob Faibussowitsch PetscCall(VecCopy(X, vec)); 6899566063dSJacob Faibussowitsch PetscCall(TSGetMaxSteps(ts, &nsteps)); 6909566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts, &dt)); 6919566063dSJacob Faibussowitsch PetscCall(TSSetMaxSteps(ts, 1)); 6929566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, X)); 6939566063dSJacob Faibussowitsch PetscCall(TSSetMaxSteps(ts, nsteps)); 6949566063dSJacob Faibussowitsch PetscCall(TSSetStepNumber(ts, 0)); 6959566063dSJacob Faibussowitsch PetscCall(TSSetTime(ts, 0)); 6969566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, dt)); 697a587d139SMark rectx->plotIdx = 0; 698930e68a5SMark Adams rectx->plotting = PETSC_FALSE; 6999566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 7009566063dSJacob Faibussowitsch PetscCall(VecCopy(vec, X)); 7019566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vec)); 702c3e4dd79SMark Adams PetscCall(PetscObjectStateIncrease((PetscObject)ctx->J)); 703e0eea495SMark } 704e0eea495SMark /* go */ 7059566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("Solve", &stage)); 706f37ccb5fSMark Adams ctx->stage = 0; // lets not use this stage 7079566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stage)); 708419c2fb8Smarkadams4 #if defined(PETSC_HAVE_CUDA_NVTX) 709419c2fb8Smarkadams4 nvtxRangePushA("ex2-TSSolve-warm"); 710419c2fb8Smarkadams4 #endif 7119566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, X)); 712419c2fb8Smarkadams4 #if defined(PETSC_HAVE_CUDA_NVTX) 713419c2fb8Smarkadams4 nvtxRangePop(); 714419c2fb8Smarkadams4 #endif 7159566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 716e0eea495SMark /* clean up */ 7179566063dSJacob Faibussowitsch PetscCall(DMPlexLandauDestroyVelocitySpace(&pack)); 7189566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 7199566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X)); 7209566063dSJacob Faibussowitsch PetscCall(PetscFree(rectx)); 7219566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 722b122ec5aSJacob Faibussowitsch return 0; 723e0eea495SMark } 724e0eea495SMark 725e0eea495SMark /*TEST 726a587d139SMark 7275dac466eSMark Adams testset: 728984ed092SMark Adams requires: p4est !complex double defined(PETSC_USE_DMLANDAU_2D) 7295dac466eSMark Adams output_file: output/ex2_0.out 7306b664d00Smarkadams4 args: -dm_landau_num_species_grid 1,1 -dm_landau_Ez 0 -petscspace_degree 3 -petscspace_poly_tensor 1 -dm_landau_type p4est -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-9 -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-2 -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 -dm_landau_amr_levels_max 2,2 -dm_landau_amr_re_levels 2 -dm_landau_re_radius 0 -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 0 -dm_landau_gpu_assembly true -dm_landau_batch_size 2 -dm_landau_verbose 2 -dm_landau_domain_radius 5.,5. 7315dac466eSMark Adams test: 7325dac466eSMark Adams suffix: cpu 733c3e4dd79SMark Adams args: -dm_landau_device_type cpu -ksp_type bicg -pc_type jacobi 734a587d139SMark test: 735a587d139SMark suffix: kokkos 736c36a59beSJunchao Zhang requires: kokkos_kernels !defined(PETSC_HAVE_CUDA_CLANG) !defined(PETSC_HAVE_SYCL) 737c3e4dd79SMark Adams args: -dm_landau_device_type kokkos -dm_mat_type aijkokkos -dm_vec_type kokkos -ksp_type bicg -pc_type jacobi 73890d163c9SMark Adams test: 739cb25d741SMark Adams suffix: kokkos_batch 740cb25d741SMark Adams requires: kokkos_kernels 741c3e4dd79SMark Adams args: -dm_landau_device_type kokkos -dm_mat_type aijkokkos -dm_vec_type kokkos -ksp_type preonly -pc_type bjkokkos -pc_bjkokkos_ksp_type bicg -pc_bjkokkos_pc_type jacobi 742bfc784b7SMark Adams test: 743a4313204SMark Adams suffix: kokkos_batch_tfqmr 744c36a59beSJunchao Zhang requires: kokkos_kernels !defined(PETSC_HAVE_CUDA_CLANG) !defined(PETSC_HAVE_SYCL) 745a4313204SMark Adams args: -dm_landau_device_type kokkos -dm_mat_type aijkokkos -dm_vec_type kokkos -ksp_type preonly -pc_type bjkokkos -pc_bjkokkos_ksp_type tfqmr -pc_bjkokkos_pc_type jacobi 74690d163c9SMark Adams 7476b664d00Smarkadams4 test: 7486b664d00Smarkadams4 requires: !complex double defined(PETSC_USE_DMLANDAU_2D) !kokkos_kernels !cuda 7496b664d00Smarkadams4 suffix: single 7506b664d00Smarkadams4 nsize: 1 7516b664d00Smarkadams4 args: -dm_refine 2 -dm_landau_num_species_grid 1 -dm_landau_thermal_temps 1 -dm_landau_electron_shift 1.25 -petscspace_degree 3 -snes_converged_reason -ts_type beuler -ts_dt .1 -ex2_plot_dt .1 -ts_max_steps 1 -ex2_grid_view_idx 0 -ex2_dm_view -snes_rtol 1.e-13 -snes_stol 1.e-13 -dm_landau_verbose 2 -ex2_print_period 1 -ksp_type preonly -pc_type lu -dm_landau_device_type cpu -dm_landau_use_relativistic_corrections 7526b664d00Smarkadams4 753cd27c6deSmarkadams4 testset: 754cd27c6deSmarkadams4 requires: !complex double defined(PETSC_USE_DMLANDAU_2D) 755cd27c6deSmarkadams4 nsize: 1 756cd27c6deSmarkadams4 output_file: output/ex2_simplex.out 757cd27c6deSmarkadams4 args: -dim 2 -dm_landau_num_species_grid 1,1 -petscspace_degree 2 -dm_landau_simplex -dm_landau_ion_masses 2 -dm_landau_ion_charges 1 -dm_landau_thermal_temps 2,1 -dm_landau_n 1,1 -snes_rtol 1e-15 -snes_stol 1e-15 -snes_monitor -ts_type beuler -snes_converged_reason -ts_exact_final_time stepover -ts_dt .1 -ts_max_steps 1 -ts_max_snes_failures -1 -ksp_type preonly -pc_type lu -dm_landau_verbose 2 -ex2_grid_view_idx 0 -ex2_dm_view -dm_refine 1 -ksp_type bicg -pc_type jacobi 758cd27c6deSmarkadams4 test: 759cd27c6deSmarkadams4 suffix: simplex 760cd27c6deSmarkadams4 args: -dm_landau_device_type cpu 761cd27c6deSmarkadams4 test: 762cd27c6deSmarkadams4 suffix: simplexkokkos 763c36a59beSJunchao Zhang requires: kokkos_kernels !defined(PETSC_HAVE_CUDA_CLANG) !defined(PETSC_HAVE_SYCL) 764cd27c6deSmarkadams4 args: -dm_landau_device_type kokkos -dm_mat_type aijkokkos -dm_vec_type kokkos 765cd27c6deSmarkadams4 766e0eea495SMark TEST*/ 767