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