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