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