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