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