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> 78a6f2e61SMark Adams #include <petscdmcomposite.h> 8c3e4dd79SMark Adams #include "petsc/private/petscimpl.h" 9e0eea495SMark 10*419c2fb8Smarkadams4 #if defined(PETSC_HAVE_CUDA_NVTX) 11*419c2fb8Smarkadams4 #include <nvToolsExt.h> 12*419c2fb8Smarkadams4 #endif 13*419c2fb8Smarkadams4 14e0eea495SMark /* data for runaway electron model */ 15e0eea495SMark typedef struct REctx_struct { 168a6f2e61SMark Adams PetscErrorCode (*test)(TS, Vec, PetscInt, PetscReal, PetscBool, LandauCtx *, struct REctx_struct *); 17e0eea495SMark PetscErrorCode (*impuritySrcRate)(PetscReal, PetscReal *, LandauCtx*); 18e0eea495SMark PetscErrorCode (*E)(Vec, Vec, PetscInt, PetscReal, LandauCtx*, PetscReal *); 19e0eea495SMark PetscReal T_cold; /* temperature of newly ionized electrons and impurity ions */ 20e0eea495SMark PetscReal ion_potential; /* ionization potential of impurity */ 21e0eea495SMark PetscReal Ne_ion; /* effective number of electrons shed in ioization of impurity */ 22e0eea495SMark PetscReal Ez_initial; 23e0eea495SMark PetscReal L; /* inductance */ 24e0eea495SMark Vec X_0; 25e0eea495SMark PetscInt imp_idx; /* index for impurity ionizing sink */ 26e0eea495SMark PetscReal pulse_start; 27e0eea495SMark PetscReal pulse_width; 28e0eea495SMark PetscReal pulse_rate; 29e0eea495SMark PetscReal current_rate; 30e0eea495SMark PetscInt plotIdx; 31e0eea495SMark PetscInt plotStep; 32e0eea495SMark PetscInt idx; /* cache */ 33e0eea495SMark PetscReal j; /* cache */ 34e0eea495SMark PetscReal plotDt; 35e0eea495SMark PetscBool plotting; 36e0eea495SMark PetscBool use_spitzer_eta; 378a6f2e61SMark Adams PetscInt print_period; 388fdabdddSMark Adams PetscInt grid_view_idx; 39e0eea495SMark } REctx; 40e0eea495SMark 41e0eea495SMark static const PetscReal kev_joul = 6.241506479963235e+15; /* 1/1000e */ 42e0eea495SMark 4390d163c9SMark Adams #define RE_CUT 3. 44e0eea495SMark /* < v, u_re * v * q > */ 45e0eea495SMark static void f0_j_re(PetscInt dim, PetscInt Nf, PetscInt NfAux, 46e0eea495SMark const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 47e0eea495SMark const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 48e0eea495SMark PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0) 49e0eea495SMark { 50e0eea495SMark PetscReal n_e = PetscRealPart(u[0]); 51e0eea495SMark if (dim==2) { 52e0eea495SMark if (x[1] > RE_CUT || x[1] < -RE_CUT) { /* simply a cutoff for REs. v_|| > 3 v(T_e) */ 53e0eea495SMark *f0 = n_e * 2.*PETSC_PI*x[0] * x[1] * constants[0]; /* n * r * v_|| * q */ 54e0eea495SMark } else { 55e0eea495SMark *f0 = 0; 56e0eea495SMark } 57e0eea495SMark } else { 58e0eea495SMark if (x[2] > RE_CUT || x[2] < -RE_CUT) { /* simply a cutoff for REs. v_|| > 3 v(T_e) */ 59e0eea495SMark *f0 = n_e * x[2] * constants[0]; 60e0eea495SMark } else { 61e0eea495SMark *f0 = 0; 62e0eea495SMark } 63e0eea495SMark } 64e0eea495SMark } 65e0eea495SMark 66e0eea495SMark /* sum < v, u*v*q > */ 67e0eea495SMark static void f0_jz_sum(PetscInt dim, PetscInt Nf, PetscInt NfAux, 68e0eea495SMark const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 69e0eea495SMark const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 708a6f2e61SMark Adams PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar q[], PetscScalar *f0) 71e0eea495SMark { 72e0eea495SMark PetscInt ii; 73e0eea495SMark f0[0] = 0; 74e0eea495SMark if (dim==2) { 758a6f2e61SMark 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 */ 76e0eea495SMark } else { 778a6f2e61SMark Adams for (ii=0;ii<Nf;ii++) f0[0] += u[ii] * x[2] * q[ii]; /* n * v_|| * q * v_0 */ 78e0eea495SMark } 79e0eea495SMark } 80e0eea495SMark 81e0eea495SMark /* < v, n_e > */ 82e0eea495SMark static void f0_n(PetscInt dim, PetscInt Nf, PetscInt NfAux, 83e0eea495SMark const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 84e0eea495SMark const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 85e0eea495SMark PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0) 86e0eea495SMark { 87e0eea495SMark PetscInt ii = (PetscInt)PetscRealPart(constants[0]); 88e0eea495SMark if (dim==2) f0[0] = 2.*PETSC_PI*x[0]*u[ii]; 89e0eea495SMark else { 90e0eea495SMark f0[0] = u[ii]; 91e0eea495SMark } 92e0eea495SMark } 93e0eea495SMark 94e0eea495SMark /* < v, n_e v_|| > */ 95e0eea495SMark static void f0_vz(PetscInt dim, PetscInt Nf, PetscInt NfAux, 96e0eea495SMark const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 97e0eea495SMark const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 98e0eea495SMark PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0) 99e0eea495SMark { 100e0eea495SMark PetscInt ii = (PetscInt)PetscRealPart(constants[0]); 101e0eea495SMark if (dim==2) f0[0] = u[ii] * 2.*PETSC_PI*x[0] * x[1]; /* n r v_|| */ 102e0eea495SMark else { 103e0eea495SMark f0[0] = u[ii] * x[2]; /* n v_|| */ 104e0eea495SMark } 105e0eea495SMark } 106e0eea495SMark 107e0eea495SMark /* < v, n_e (v-shift) > */ 108e0eea495SMark static void f0_ve_shift(PetscInt dim, PetscInt Nf, PetscInt NfAux, 109e0eea495SMark const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 110e0eea495SMark const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 111e0eea495SMark PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0) 112e0eea495SMark { 1138a6f2e61SMark Adams PetscReal vz = numConstants>0 ? PetscRealPart(constants[0]) : 0; 114e0eea495SMark 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 */ 115e0eea495SMark else { 116e0eea495SMark *f0 = u[0] * PetscSqrtReal(x[0]*x[0] + x[1]*x[1] + (x[2]-vz)*(x[2]-vz)); /* n v */ 117e0eea495SMark } 118e0eea495SMark } 119e0eea495SMark 120e0eea495SMark /* CalculateE - Calculate the electric field */ 121e0eea495SMark /* T -- Electron temperature */ 122e0eea495SMark /* n -- Electron density */ 123e0eea495SMark /* lnLambda -- */ 124e0eea495SMark /* eps0 -- */ 125e0eea495SMark /* E -- output E, input \hat E */ 126e0eea495SMark static PetscReal CalculateE(PetscReal Tev, PetscReal n, PetscReal lnLambda, PetscReal eps0, PetscReal *E) 127e0eea495SMark { 128e0eea495SMark PetscReal c,e,m; 129e0eea495SMark 130e0eea495SMark PetscFunctionBegin; 131cefb98e8SMark Adams c = 299792458.0; 132e0eea495SMark e = 1.602176e-19; 133e0eea495SMark m = 9.10938e-31; 134e0eea495SMark if (1) { 135e0eea495SMark double Ec, Ehat = *E, betath = PetscSqrtReal(2*Tev*e/(m*c*c)), j0 = Ehat * 7/(PetscSqrtReal(2)*2) * PetscPowReal(betath,3) * n * e * c; 136e0eea495SMark Ec = n*lnLambda*PetscPowReal(e,3) / (4*PETSC_PI*PetscPowReal(eps0,2)*m*c*c); 137e0eea495SMark *E = Ec; 138e0eea495SMark PetscPrintf(PETSC_COMM_WORLD, "CalculateE j0=%g Ec = %g\n",j0,Ec); 139e0eea495SMark } else { 140e0eea495SMark PetscReal Ed, vth; 141e0eea495SMark vth = PetscSqrtReal(8*Tev*e/(m*PETSC_PI)); 142e0eea495SMark Ed = n*lnLambda*PetscPowReal(e,3) / (4*PETSC_PI*PetscPowReal(eps0,2)*m*vth*vth); 143e0eea495SMark *E = Ed; 144e0eea495SMark } 145e0eea495SMark PetscFunctionReturn(0); 146e0eea495SMark } 147e0eea495SMark 148e0eea495SMark static PetscReal Spitzer(PetscReal m_e, PetscReal e, PetscReal Z, PetscReal epsilon0, PetscReal lnLam, PetscReal kTe_joules) 149e0eea495SMark { 150e0eea495SMark PetscReal Fz = (1+1.198*Z+0.222*Z*Z)/(1+2.966*Z+0.753*Z*Z), eta; 151e0eea495SMark 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); 152e0eea495SMark return eta; 153e0eea495SMark } 154e0eea495SMark 155e0eea495SMark /* */ 1568a6f2e61SMark Adams static PetscErrorCode testNone(TS ts, Vec X, PetscInt stepi, PetscReal time, PetscBool islast, LandauCtx *ctx, REctx *rectx) 157e0eea495SMark { 158e0eea495SMark PetscFunctionBeginUser; 159e0eea495SMark PetscFunctionReturn(0); 160e0eea495SMark } 161e0eea495SMark 162e0eea495SMark /* */ 1638a6f2e61SMark Adams static PetscErrorCode testSpitzer(TS ts, Vec X, PetscInt stepi, PetscReal time, PetscBool islast, LandauCtx *ctx, REctx *rectx) 164e0eea495SMark { 165f37ccb5fSMark Adams PetscInt ii,nDMs; 166e0eea495SMark PetscDS prob; 167a587d139SMark static PetscReal old_ratio = 1e10; 168a587d139SMark TSConvergedReason reason; 169a587d139SMark PetscReal J,J_re,spit_eta,Te_kev=0,E,ratio,Z,n_e,v,v2; 1708a6f2e61SMark Adams PetscScalar user[2] = {0.,ctx->charges[0]}, q[LANDAU_MAX_SPECIES],tt[LANDAU_MAX_SPECIES],vz; 171930e68a5SMark Adams PetscReal dt; 1728a6f2e61SMark Adams DM pack, plexe = ctx->plex[0], plexi = (ctx->num_grids==1) ? NULL : ctx->plex[1]; 173f37ccb5fSMark Adams Vec *XsubArray; 174e0eea495SMark 175e0eea495SMark PetscFunctionBeginUser; 17663a3b9bcSJacob Faibussowitsch PetscCheck(ctx->num_species==2,PETSC_COMM_SELF, PETSC_ERR_PLIB, "ctx->num_species %" PetscInt_FMT " != 2",ctx->num_species); 1779566063dSJacob Faibussowitsch PetscCall(VecGetDM(X, &pack)); 1783c633725SBarry Smith PetscCheck(pack,PETSC_COMM_SELF, PETSC_ERR_PLIB, "no DM"); 1799566063dSJacob Faibussowitsch PetscCall(DMCompositeGetNumberDM(pack,&nDMs)); 18063a3b9bcSJacob Faibussowitsch PetscCheck(nDMs == ctx->num_grids*ctx->batch_sz,PETSC_COMM_SELF, PETSC_ERR_PLIB, "nDMs != ctx->num_grids*ctx->batch_sz %" PetscInt_FMT " != %" PetscInt_FMT,nDMs,ctx->num_grids*ctx->batch_sz); 1819566063dSJacob Faibussowitsch PetscCall(PetscMalloc(sizeof(*XsubArray)*nDMs, &XsubArray)); 1829566063dSJacob Faibussowitsch PetscCall(DMCompositeGetAccessArray(pack, X, nDMs, NULL, XsubArray)); // read only 1839566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts,&dt)); 1848a6f2e61SMark Adams /* get current for each grid */ 1858a6f2e61SMark Adams for (ii=0;ii<ctx->num_species;ii++) q[ii] = ctx->charges[ii]; 1869566063dSJacob Faibussowitsch PetscCall(DMGetDS(plexe, &prob)); 1879566063dSJacob Faibussowitsch PetscCall(PetscDSSetConstants(prob, 2, &q[0])); 1889566063dSJacob Faibussowitsch PetscCall(PetscDSSetObjective(prob, 0, &f0_jz_sum)); 1899566063dSJacob Faibussowitsch PetscCall(DMPlexComputeIntegralFEM(plexe,XsubArray[ LAND_PACK_IDX(ctx->batch_view_idx,0) ],tt,NULL)); 190e0eea495SMark J = -ctx->n_0*ctx->v_0*PetscRealPart(tt[0]); 1918fdabdddSMark Adams if (plexi) { // add first (only) ion 1929566063dSJacob Faibussowitsch PetscCall(DMGetDS(plexi, &prob)); 1939566063dSJacob Faibussowitsch PetscCall(PetscDSSetConstants(prob, 1, &q[1])); 1949566063dSJacob Faibussowitsch PetscCall(PetscDSSetObjective(prob, 0, &f0_jz_sum)); 1959566063dSJacob Faibussowitsch PetscCall(DMPlexComputeIntegralFEM(plexi,XsubArray[LAND_PACK_IDX(ctx->batch_view_idx,1)],tt,NULL)); 1968a6f2e61SMark Adams J += -ctx->n_0*ctx->v_0*PetscRealPart(tt[0]); 1978a6f2e61SMark Adams } 198a587d139SMark /* get N_e */ 1999566063dSJacob Faibussowitsch PetscCall(DMGetDS(plexe, &prob)); 2009566063dSJacob Faibussowitsch PetscCall(PetscDSSetConstants(prob, 1, user)); 2019566063dSJacob Faibussowitsch PetscCall(PetscDSSetObjective(prob, 0, &f0_n)); 2029566063dSJacob Faibussowitsch PetscCall(DMPlexComputeIntegralFEM(plexe,XsubArray[LAND_PACK_IDX(ctx->batch_view_idx,0)],tt,NULL)); 203a587d139SMark n_e = PetscRealPart(tt[0])*ctx->n_0; 204a587d139SMark /* Z */ 205a587d139SMark Z = -ctx->charges[1]/ctx->charges[0]; 206a587d139SMark /* remove drift */ 2078a6f2e61SMark Adams if (0) { 2088a6f2e61SMark Adams user[0] = 0; // electrons 2099566063dSJacob Faibussowitsch PetscCall(DMGetDS(plexe, &prob)); 2109566063dSJacob Faibussowitsch PetscCall(PetscDSSetConstants(prob, 1, user)); 2119566063dSJacob Faibussowitsch PetscCall(PetscDSSetObjective(prob, 0, &f0_vz)); 2129566063dSJacob Faibussowitsch PetscCall(DMPlexComputeIntegralFEM(plexe,XsubArray[LAND_PACK_IDX(ctx->batch_view_idx,0)],tt,NULL)); 213a587d139SMark vz = ctx->n_0*PetscRealPart(tt[0])/n_e; /* non-dimensional */ 214a587d139SMark } else vz = 0; 215a587d139SMark /* thermal velocity */ 2169566063dSJacob Faibussowitsch PetscCall(DMGetDS(plexe, &prob)); 2179566063dSJacob Faibussowitsch PetscCall(PetscDSSetConstants(prob, 1, &vz)); 2189566063dSJacob Faibussowitsch PetscCall(PetscDSSetObjective(prob, 0, &f0_ve_shift)); 2199566063dSJacob Faibussowitsch PetscCall(DMPlexComputeIntegralFEM(plexe,XsubArray[LAND_PACK_IDX(ctx->batch_view_idx,0)],tt,NULL)); 220a587d139SMark v = ctx->n_0*ctx->v_0*PetscRealPart(tt[0])/n_e; /* remove number density to get velocity */ 221a587d139SMark v2 = PetscSqr(v); /* use real space: m^2 / s^2 */ 222a587d139SMark Te_kev = (v2*ctx->masses[0]*PETSC_PI/8)*kev_joul; /* temperature in kev */ 223a587d139SMark spit_eta = Spitzer(ctx->masses[0],-ctx->charges[0],Z,ctx->epsilon0,ctx->lnLam,Te_kev/kev_joul); /* kev --> J (kT) */ 2248a6f2e61SMark Adams if (0) { 2259566063dSJacob Faibussowitsch PetscCall(DMGetDS(plexe, &prob)); 2269566063dSJacob Faibussowitsch PetscCall(PetscDSSetConstants(prob, 1, q)); 2279566063dSJacob Faibussowitsch PetscCall(PetscDSSetObjective(prob, 0, &f0_j_re)); 2289566063dSJacob Faibussowitsch PetscCall(DMPlexComputeIntegralFEM(plexe,XsubArray[LAND_PACK_IDX(ctx->batch_view_idx,0)],tt,NULL)); 229a587d139SMark } else tt[0] = 0; 230e0eea495SMark J_re = -ctx->n_0*ctx->v_0*PetscRealPart(tt[0]); 2319566063dSJacob Faibussowitsch PetscCall(DMCompositeRestoreAccessArray(pack, X, nDMs, NULL, XsubArray)); // read only 2329566063dSJacob Faibussowitsch PetscCall(PetscFree(XsubArray)); 233a587d139SMark 23490d163c9SMark Adams if (rectx->use_spitzer_eta) { 23590d163c9SMark Adams E = ctx->Ez = spit_eta*(rectx->j-J_re); 23690d163c9SMark Adams } else { 23790d163c9SMark Adams E = ctx->Ez; /* keep real E */ 23890d163c9SMark Adams rectx->j = J; /* cache */ 23990d163c9SMark Adams } 24090d163c9SMark Adams 241e0eea495SMark ratio = E/J/spit_eta; 2428a6f2e61SMark Adams if (stepi>10 && !rectx->use_spitzer_eta && ( 2438a6f2e61SMark Adams (old_ratio-ratio < 1.e-6))) { 24490d163c9SMark Adams rectx->pulse_start = time + 0.98*dt; 245a587d139SMark rectx->use_spitzer_eta = PETSC_TRUE; 246a587d139SMark } 2479566063dSJacob Faibussowitsch PetscCall(TSGetConvergedReason(ts,&reason)); 2489566063dSJacob Faibussowitsch PetscCall(TSGetConvergedReason(ts,&reason)); 24990d163c9SMark Adams if ((rectx->plotting) || stepi == 0 || reason || rectx->pulse_start == time + 0.98*dt) { 25063a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(ctx->comm, "testSpitzer: %4" PetscInt_FMT ") time=%11.4e n_e= %10.3e E= %10.3e J= %10.3e J_re= %10.3e %.3g%% Te_kev= %10.3e Z_eff=%g E/J to eta ratio= %g (diff=%g) %s %s spit_eta=%g\n",stepi,(double)time,(double)(n_e/ctx->n_0),(double)ctx->Ez,(double)J,(double)J_re,(double)(100*J_re/J),(double)Te_kev,(double)Z,(double)ratio,(double)(old_ratio-ratio), rectx->use_spitzer_eta ? "using Spitzer eta*J E" : "constant E",rectx->pulse_start != time + 0.98*dt ? "normal" : "transition",(double)spit_eta)); 25163a3b9bcSJacob Faibussowitsch PetscCheck(rectx->pulse_start != (time + 0.98*dt),PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG,"Spitzer complete ratio=%g",(double)ratio); 252930e68a5SMark Adams } 253e0eea495SMark old_ratio = ratio; 254e0eea495SMark PetscFunctionReturn(0); 255e0eea495SMark } 256e0eea495SMark 257e0eea495SMark static const double ppp = 2; 258e0eea495SMark static void f0_0_diff_lp(PetscInt dim, PetscInt Nf, PetscInt NfAux, 259e0eea495SMark const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 260e0eea495SMark const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 261e0eea495SMark PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0) 262e0eea495SMark { 263e0eea495SMark LandauCtx *ctx = (LandauCtx *)constants; 264e0eea495SMark REctx *rectx = (REctx*)ctx->data; 265e0eea495SMark PetscInt ii = rectx->idx, i; 266e0eea495SMark const PetscReal kT_m = ctx->k*ctx->thermal_temps[ii]/ctx->masses[ii]; /* kT/m */ 267e0eea495SMark const PetscReal n = ctx->n[ii]; 268e0eea495SMark PetscReal diff, f_maxwell, v2 = 0, theta = 2*kT_m/(ctx->v_0*ctx->v_0); /* theta = 2kT/mc^2 */ 269e0eea495SMark for (i = 0; i < dim; ++i) v2 += x[i]*x[i]; 270e0eea495SMark f_maxwell = n*PetscPowReal(PETSC_PI*theta,-1.5)*(PetscExpReal(-v2/theta)); 271e0eea495SMark diff = 2.*PETSC_PI*x[0]*(PetscRealPart(u[ii]) - f_maxwell); 272e0eea495SMark f0[0] = PetscPowReal(diff,ppp); 273e0eea495SMark } 274e0eea495SMark static void f0_0_maxwellian_lp(PetscInt dim, PetscInt Nf, PetscInt NfAux, 275e0eea495SMark const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 276e0eea495SMark const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 277e0eea495SMark PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0) 278e0eea495SMark { 279e0eea495SMark LandauCtx *ctx = (LandauCtx *)constants; 280e0eea495SMark REctx *rectx = (REctx*)ctx->data; 281e0eea495SMark PetscInt ii = rectx->idx, i; 282e0eea495SMark const PetscReal kT_m = ctx->k*ctx->thermal_temps[ii]/ctx->masses[ii]; /* kT/m */ 283e0eea495SMark const PetscReal n = ctx->n[ii]; 284e0eea495SMark PetscReal f_maxwell, v2 = 0, theta = 2*kT_m/(ctx->v_0*ctx->v_0); /* theta = 2kT/mc^2 */ 285e0eea495SMark for (i = 0; i < dim; ++i) v2 += x[i]*x[i]; 286e0eea495SMark f_maxwell = 2.*PETSC_PI*x[0] * n*PetscPowReal(PETSC_PI*theta,-1.5)*(PetscExpReal(-v2/theta)); 287e0eea495SMark f0[0] = PetscPowReal(f_maxwell,ppp); 288e0eea495SMark } 289e0eea495SMark 290e0eea495SMark /* */ 2918a6f2e61SMark Adams static PetscErrorCode testStable(TS ts, Vec X, PetscInt stepi, PetscReal time, PetscBool islast, LandauCtx *ctx, REctx *rectx) 292e0eea495SMark { 293e0eea495SMark PetscDS prob; 294e0eea495SMark Vec X2; 295e0eea495SMark PetscReal ediff,idiff=0,lpm0,lpm1=1; 296e0eea495SMark PetscScalar tt[LANDAU_MAX_SPECIES]; 2978a6f2e61SMark Adams DM dm, plex = ctx->plex[0]; 298e0eea495SMark 299e0eea495SMark PetscFunctionBeginUser; 3009566063dSJacob Faibussowitsch PetscCall(VecGetDM(X, &dm)); 3019566063dSJacob Faibussowitsch PetscCall(DMGetDS(plex, &prob)); 3029566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X,&X2)); 3039566063dSJacob Faibussowitsch PetscCall(VecCopy(X,X2)); 304e0eea495SMark if (!rectx->X_0) { 3059566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X,&rectx->X_0)); 3069566063dSJacob Faibussowitsch PetscCall(VecCopy(X,rectx->X_0)); 307e0eea495SMark } 3089566063dSJacob Faibussowitsch PetscCall(VecAXPY(X,-1.0,rectx->X_0)); 3099566063dSJacob Faibussowitsch PetscCall(PetscDSSetConstants(prob, sizeof(LandauCtx)/sizeof(PetscScalar), (PetscScalar*)ctx)); 310e0eea495SMark rectx->idx = 0; 3119566063dSJacob Faibussowitsch PetscCall(PetscDSSetObjective(prob, 0, &f0_0_diff_lp)); 3129566063dSJacob Faibussowitsch PetscCall(DMPlexComputeIntegralFEM(plex,X2,tt,NULL)); 313e0eea495SMark ediff = PetscPowReal(PetscRealPart(tt[0]),1./ppp); 3149566063dSJacob Faibussowitsch PetscCall(PetscDSSetObjective(prob, 0, &f0_0_maxwellian_lp)); 3159566063dSJacob Faibussowitsch PetscCall(DMPlexComputeIntegralFEM(plex,X2,tt,NULL)); 316e0eea495SMark lpm0 = PetscPowReal(PetscRealPart(tt[0]),1./ppp); 317e0eea495SMark if (ctx->num_species>1) { 318e0eea495SMark rectx->idx = 1; 3199566063dSJacob Faibussowitsch PetscCall(PetscDSSetObjective(prob, 0, &f0_0_diff_lp)); 3209566063dSJacob Faibussowitsch PetscCall(DMPlexComputeIntegralFEM(plex,X2,tt,NULL)); 321e0eea495SMark idiff = PetscPowReal(PetscRealPart(tt[0]),1./ppp); 3229566063dSJacob Faibussowitsch PetscCall(PetscDSSetObjective(prob, 0, &f0_0_maxwellian_lp)); 3239566063dSJacob Faibussowitsch PetscCall(DMPlexComputeIntegralFEM(plex,X2,tt,NULL)); 324e0eea495SMark lpm1 = PetscPowReal(PetscRealPart(tt[0]),1./ppp); 325e0eea495SMark } 32663a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s %" PetscInt_FMT ") time=%10.3e n-%d norm electrons/max=%20.13e ions/max=%20.13e\n", "----",stepi,(double)time,(int)ppp,(double)(ediff/lpm0),(double)(idiff/lpm1))); 327e0eea495SMark /* view */ 3289566063dSJacob Faibussowitsch PetscCall(VecCopy(X2,X)); 3299566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X2)); 330e0eea495SMark if (islast) { 3319566063dSJacob Faibussowitsch PetscCall(VecDestroy(&rectx->X_0)); 332e0eea495SMark rectx->X_0 = NULL; 333e0eea495SMark } 334e0eea495SMark PetscFunctionReturn(0); 335e0eea495SMark } 336e0eea495SMark 337e0eea495SMark static PetscErrorCode EInduction(Vec X, Vec X_t, PetscInt step, PetscReal time, LandauCtx *ctx, PetscReal *a_E) 338e0eea495SMark { 339e0eea495SMark REctx *rectx = (REctx*)ctx->data; 340e0eea495SMark PetscInt ii; 341e0eea495SMark DM dm,plex; 3428a6f2e61SMark Adams PetscScalar tt[LANDAU_MAX_SPECIES], qv0[LANDAU_MAX_SPECIES]; 343e0eea495SMark PetscReal dJ_dt; 344e0eea495SMark PetscDS prob; 345e0eea495SMark 346e0eea495SMark PetscFunctionBeginUser; 3478a6f2e61SMark Adams for (ii=0;ii<ctx->num_species;ii++) qv0[ii] = ctx->charges[ii]*ctx->v_0; 3489566063dSJacob Faibussowitsch PetscCall(VecGetDM(X, &dm)); 3499566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &prob)); 3509566063dSJacob Faibussowitsch PetscCall(DMConvert(dm, DMPLEX, &plex)); 351e0eea495SMark /* get d current / dt */ 3529566063dSJacob Faibussowitsch PetscCall(PetscDSSetConstants(prob, ctx->num_species, qv0)); 3539566063dSJacob Faibussowitsch PetscCall(PetscDSSetObjective(prob, 0, &f0_jz_sum)); 3543c633725SBarry Smith PetscCheck(X_t,PETSC_COMM_SELF, PETSC_ERR_PLIB, "X_t"); 3559566063dSJacob Faibussowitsch PetscCall(DMPlexComputeIntegralFEM(plex,X_t,tt,NULL)); 3568a6f2e61SMark Adams dJ_dt = -ctx->n_0*PetscRealPart(tt[0])/ctx->t_0; 357e0eea495SMark /* E induction */ 358e0eea495SMark *a_E = -rectx->L*dJ_dt + rectx->Ez_initial; 3599566063dSJacob Faibussowitsch PetscCall(DMDestroy(&plex)); 360e0eea495SMark PetscFunctionReturn(0); 361e0eea495SMark } 362e0eea495SMark 363e0eea495SMark static PetscErrorCode EConstant(Vec X, Vec X_t, PetscInt step, PetscReal time, LandauCtx *ctx, PetscReal *a_E) 364e0eea495SMark { 365e0eea495SMark PetscFunctionBeginUser; 366e0eea495SMark *a_E = ctx->Ez; 367e0eea495SMark PetscFunctionReturn(0); 368e0eea495SMark } 369e0eea495SMark 370e0eea495SMark static PetscErrorCode ENone(Vec X, Vec X_t, PetscInt step, PetscReal time, LandauCtx *ctx, PetscReal *a_E) 371e0eea495SMark { 372e0eea495SMark PetscFunctionBeginUser; 373e0eea495SMark *a_E = 0; 374e0eea495SMark PetscFunctionReturn(0); 375e0eea495SMark } 376e0eea495SMark 377e0eea495SMark /* ------------------------------------------------------------------- */ 378e0eea495SMark /* 379e0eea495SMark FormSource - Evaluates source terms F(t). 380e0eea495SMark 381e0eea495SMark Input Parameters: 382e0eea495SMark . ts - the TS context 383e0eea495SMark . time - 384e0eea495SMark . X_dummmy - input vector 385e0eea495SMark . dummy - optional user-defined context, as set by SNESSetFunction() 386e0eea495SMark 387e0eea495SMark Output Parameter: 388e0eea495SMark . F - function vector 389e0eea495SMark */ 390e0eea495SMark static PetscErrorCode FormSource(TS ts, PetscReal ftime, Vec X_dummmy, Vec F, void *dummy) 391e0eea495SMark { 392e0eea495SMark PetscReal new_imp_rate; 393e0eea495SMark LandauCtx *ctx; 3948a6f2e61SMark Adams DM pack; 395e0eea495SMark REctx *rectx; 396e0eea495SMark 397e0eea495SMark PetscFunctionBeginUser; 3989566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts,&pack)); 3999566063dSJacob Faibussowitsch PetscCall(DMGetApplicationContext(pack, &ctx)); 400e0eea495SMark rectx = (REctx*)ctx->data; 401e0eea495SMark /* check for impurities */ 4029566063dSJacob Faibussowitsch PetscCall(rectx->impuritySrcRate(ftime,&new_imp_rate,ctx)); 403e0eea495SMark if (new_imp_rate != 0) { 404e0eea495SMark if (new_imp_rate != rectx->current_rate) { 405e0eea495SMark PetscInt ii; 406e0eea495SMark PetscReal dne_dt,dni_dt,tilda_ns[LANDAU_MAX_SPECIES],temps[LANDAU_MAX_SPECIES]; 4078fdabdddSMark Adams Vec globFarray[LANDAU_MAX_GRIDS*LANDAU_MAX_BATCH_SZ]; 408e0eea495SMark rectx->current_rate = new_imp_rate; 409e0eea495SMark for (ii=1;ii<LANDAU_MAX_SPECIES;ii++) tilda_ns[ii] = 0; 410e0eea495SMark for (ii=1;ii<LANDAU_MAX_SPECIES;ii++) temps[ii] = 1; 4118a6f2e61SMark Adams dni_dt = new_imp_rate /* *ctx->t_0 */; /* fully ionized immediately, no normalize, stay in non-dim */ 4128a6f2e61SMark Adams dne_dt = new_imp_rate*rectx->Ne_ion /* *ctx->t_0 */; 413e0eea495SMark tilda_ns[0] = dne_dt; tilda_ns[rectx->imp_idx] = dni_dt; 414e0eea495SMark temps[0] = rectx->T_cold; temps[rectx->imp_idx] = rectx->T_cold; 41563a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(ctx->plex[0], "\tHave new_imp_rate= %10.3e time= %10.3e de/dt= %10.3e di/dt= %10.3e ***\n",(double)new_imp_rate,(double)ftime,(double)dne_dt,(double)dni_dt)); 4169566063dSJacob Faibussowitsch PetscCall(DMCompositeGetAccessArray(pack, F, ctx->num_grids*ctx->batch_sz, NULL, globFarray)); 4178a6f2e61SMark Adams for (PetscInt grid=0 ; grid<ctx->num_grids ; grid++) { 418e0eea495SMark /* add it */ 419c3e4dd79SMark Adams PetscCall(DMPlexLandauAddMaxwellians(ctx->plex[grid],globFarray[ LAND_PACK_IDX(0,grid) ],ftime,temps,tilda_ns,grid,0,1,ctx)); 420e0eea495SMark } 4218fdabdddSMark Adams // Does DMCompositeRestoreAccessArray copy the data back? (no) 4229566063dSJacob Faibussowitsch PetscCall(DMCompositeRestoreAccessArray(pack, F, ctx->num_grids*ctx->batch_sz, NULL, globFarray)); 423e0eea495SMark } 424e0eea495SMark } else { 4259566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F)); 426e0eea495SMark rectx->current_rate = 0; 427e0eea495SMark } 428e0eea495SMark PetscFunctionReturn(0); 429e0eea495SMark } 430*419c2fb8Smarkadams4 431e0eea495SMark PetscErrorCode Monitor(TS ts, PetscInt stepi, PetscReal time, Vec X, void *actx) 432e0eea495SMark { 433e0eea495SMark LandauCtx *ctx = (LandauCtx*) actx; /* user-defined application context */ 434e0eea495SMark REctx *rectx = (REctx*)ctx->data; 435*419c2fb8Smarkadams4 DM pack=NULL; 4368fdabdddSMark Adams Vec globXArray[LANDAU_MAX_GRIDS*LANDAU_MAX_BATCH_SZ]; 437e0eea495SMark TSConvergedReason reason; 438*419c2fb8Smarkadams4 439e0eea495SMark PetscFunctionBeginUser; 440*419c2fb8Smarkadams4 PetscCall(TSGetConvergedReason(ts,&reason)); 441*419c2fb8Smarkadams4 if (rectx->grid_view_idx != -1 || (reason && ctx->verbose > 3)) { 4429566063dSJacob Faibussowitsch PetscCall(VecGetDM(X, &pack)); 4439566063dSJacob Faibussowitsch PetscCall(DMCompositeGetAccessArray(pack, X, ctx->num_grids*ctx->batch_sz, NULL, globXArray)); 444*419c2fb8Smarkadams4 } 445e0eea495SMark if (stepi > rectx->plotStep && rectx->plotting) { 446e0eea495SMark rectx->plotting = PETSC_FALSE; /* was doing diagnostics, now done */ 447e0eea495SMark rectx->plotIdx++; 448e0eea495SMark } 449e0eea495SMark /* view */ 450e0eea495SMark if (time/rectx->plotDt >= (PetscReal)rectx->plotIdx || reason) { 451c3e4dd79SMark Adams if ((reason || stepi==0 || rectx->plotIdx%rectx->print_period==0) && ctx->verbose > 1) { 452e0eea495SMark /* print norms */ 4539566063dSJacob Faibussowitsch PetscCall(DMPlexLandauPrintNorms(X, stepi)); 454a587d139SMark } 455930e68a5SMark Adams if (!rectx->plotting) { /* first step of possible backtracks */ 456930e68a5SMark Adams rectx->plotting = PETSC_TRUE; 457930e68a5SMark Adams /* diagnostics + change E field with Sptizer (not just a monitor) */ 4589566063dSJacob Faibussowitsch PetscCall(rectx->test(ts,X,stepi,time,reason ? PETSC_TRUE : PETSC_FALSE, ctx, rectx)); 459930e68a5SMark Adams } else { 460930e68a5SMark Adams PetscPrintf(PETSC_COMM_WORLD, "\t\t ERROR SKIP test spit ------\n"); 461930e68a5SMark Adams rectx->plotting = PETSC_TRUE; 462930e68a5SMark Adams } 463*419c2fb8Smarkadams4 if (rectx->grid_view_idx != -1) { 4649566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) globXArray[ LAND_PACK_IDX(ctx->batch_view_idx,rectx->grid_view_idx) ], rectx->grid_view_idx==0 ? "ue" : "ui")); 465930e68a5SMark Adams /* view, overwrite step when back tracked */ 4669566063dSJacob Faibussowitsch PetscCall(DMSetOutputSequenceNumber(pack, rectx->plotIdx, time*ctx->t_0)); 467*419c2fb8Smarkadams4 PetscCall(VecViewFromOptions(globXArray[ LAND_PACK_IDX(ctx->batch_view_idx, rectx->grid_view_idx) ],NULL,"-ex2_vec_view")); 468*419c2fb8Smarkadams4 } 469e0eea495SMark rectx->plotStep = stepi; 470930e68a5SMark Adams } else { 47163a3b9bcSJacob Faibussowitsch if (rectx->plotting) PetscPrintf(PETSC_COMM_WORLD," ERROR rectx->plotting=%s step %" PetscInt_FMT "\n",PetscBools[rectx->plotting],stepi); 472930e68a5SMark Adams /* diagnostics + change E field with Sptizer (not just a monitor) - can we lag this? */ 4739566063dSJacob Faibussowitsch PetscCall(rectx->test(ts,X,stepi,time,reason ? PETSC_TRUE : PETSC_FALSE, ctx, rectx)); 474e0eea495SMark } 475f37ccb5fSMark Adams /* parallel check that only works of all batches are identical */ 476f37ccb5fSMark Adams if (reason && ctx->verbose > 3) { 477e0eea495SMark PetscReal val,rval; 478e0eea495SMark PetscMPIInt rank; 4799566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 480f37ccb5fSMark Adams for (PetscInt grid=0;grid<ctx->num_grids;grid++) { 481f37ccb5fSMark Adams PetscInt nerrors=0; 482f37ccb5fSMark Adams for (PetscInt i=0; i<ctx->batch_sz;i++) { 4839566063dSJacob Faibussowitsch PetscCall(VecNorm(globXArray[ LAND_PACK_IDX(i,grid) ],NORM_2,&val)); 484f37ccb5fSMark Adams if (i==0) rval = val; 485f37ccb5fSMark Adams else if ((val=PetscAbs(val-rval)/rval) > 1000*PETSC_MACHINE_EPSILON) { 48663a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " [%d] Warning %" PetscInt_FMT ".%" PetscInt_FMT ") diff = %2.15e\n",rank,grid,i,(double)val)); 487f37ccb5fSMark Adams nerrors++; 488f37ccb5fSMark Adams } 489f37ccb5fSMark Adams } 490f37ccb5fSMark Adams if (nerrors) { 49163a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " ***** [%d] ERROR max %" PetscInt_FMT " errors\n",rank,nerrors)); 492e0eea495SMark } else { 49363a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "[%d] %" PetscInt_FMT ") batch consistency check OK\n",rank,grid)); 494f37ccb5fSMark Adams } 495e0eea495SMark } 496e0eea495SMark } 497e0eea495SMark rectx->idx = 0; 498*419c2fb8Smarkadams4 if (rectx->grid_view_idx != -1 || (reason && ctx->verbose > 3)) { 4999566063dSJacob Faibussowitsch PetscCall(DMCompositeRestoreAccessArray(pack, X, ctx->num_grids*ctx->batch_sz, NULL, globXArray)); 500*419c2fb8Smarkadams4 } 501e0eea495SMark PetscFunctionReturn(0); 502e0eea495SMark } 503e0eea495SMark 504e0eea495SMark PetscErrorCode PreStep(TS ts) 505e0eea495SMark { 506e0eea495SMark LandauCtx *ctx; 507e0eea495SMark REctx *rectx; 508e0eea495SMark DM dm; 509e0eea495SMark PetscInt stepi; 510e0eea495SMark PetscReal time; 511e0eea495SMark Vec X; 512e0eea495SMark 513e0eea495SMark PetscFunctionBeginUser; 514930e68a5SMark Adams /* not used */ 5159566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts,&dm)); 5169566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts,&time)); 5179566063dSJacob Faibussowitsch PetscCall(TSGetSolution(ts,&X)); 5189566063dSJacob Faibussowitsch PetscCall(DMGetApplicationContext(dm, &ctx)); 519e0eea495SMark rectx = (REctx*)ctx->data; 5209566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &stepi)); 521e0eea495SMark /* update E */ 5229566063dSJacob Faibussowitsch PetscCall(rectx->E(X, NULL, stepi, time, ctx, &ctx->Ez)); 523e0eea495SMark PetscFunctionReturn(0); 524e0eea495SMark } 525e0eea495SMark 526e0eea495SMark /* 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) */ 527e0eea495SMark static PetscErrorCode stepSrc(PetscReal time, PetscReal *rho, LandauCtx *ctx) 528e0eea495SMark { 529e0eea495SMark REctx *rectx = (REctx*)ctx->data; 530e0eea495SMark 531e0eea495SMark PetscFunctionBeginUser; 532e0eea495SMark if (time >= rectx->pulse_start) *rho = rectx->pulse_rate; 533e0eea495SMark else *rho = 0.; 534e0eea495SMark PetscFunctionReturn(0); 535e0eea495SMark } 536e0eea495SMark static PetscErrorCode zeroSrc(PetscReal time, PetscReal *rho, LandauCtx *ctx) 537e0eea495SMark { 538e0eea495SMark PetscFunctionBeginUser; 539e0eea495SMark *rho = 0.; 540e0eea495SMark PetscFunctionReturn(0); 541e0eea495SMark } 542e0eea495SMark static PetscErrorCode pulseSrc(PetscReal time, PetscReal *rho, LandauCtx *ctx) 543e0eea495SMark { 544e0eea495SMark REctx *rectx = (REctx*)ctx->data; 545e0eea495SMark 546e0eea495SMark PetscFunctionBeginUser; 54708401ef6SPierre Jolivet 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'"); 548e0eea495SMark if (time < rectx->pulse_start || time > rectx->pulse_start + 3*rectx->pulse_width) *rho = 0; 549a587d139SMark else { 550e0eea495SMark double x = PetscSinReal((time-rectx->pulse_start)/(3*rectx->pulse_width)*2*PETSC_PI - PETSC_PI/2) + 1; /* 0:2, integrates to 1.0 */ 551e0eea495SMark *rho = rectx->pulse_rate * x / (3*rectx->pulse_width); 552a587d139SMark if (!rectx->use_spitzer_eta) rectx->use_spitzer_eta = PETSC_TRUE; /* use it next time */ 553e0eea495SMark } 554e0eea495SMark PetscFunctionReturn(0); 555e0eea495SMark } 556e0eea495SMark 557e0eea495SMark #undef __FUNCT__ 558e0eea495SMark #define __FUNCT__ "ProcessREOptions" 559e0eea495SMark static PetscErrorCode ProcessREOptions(REctx *rectx, const LandauCtx *ctx, DM dm, const char prefix[]) 560e0eea495SMark { 561e0eea495SMark PetscFunctionList plist = NULL, testlist = NULL, elist = NULL; 562e0eea495SMark char pname[256],testname[256],ename[256]; 563930e68a5SMark Adams DM dm_dummy; 564a587d139SMark PetscBool Connor_E = PETSC_FALSE; 565e0eea495SMark 566e0eea495SMark PetscFunctionBeginUser; 5679566063dSJacob Faibussowitsch PetscCall(DMCreate(PETSC_COMM_WORLD,&dm_dummy)); 568e0eea495SMark rectx->Ne_ion = 1; /* number of electrons given up by impurity ion */ 569e0eea495SMark rectx->T_cold = .005; /* kev */ 570e0eea495SMark rectx->ion_potential = 15; /* ev */ 571e0eea495SMark rectx->L = 2; 572e0eea495SMark rectx->X_0 = NULL; 573e0eea495SMark rectx->imp_idx = ctx->num_species - 1; /* default ionized impurity as last one */ 574a587d139SMark rectx->pulse_start = PETSC_MAX_REAL; 575e0eea495SMark rectx->pulse_width = 1; 576e0eea495SMark rectx->plotStep = PETSC_MAX_INT; 577e0eea495SMark rectx->pulse_rate = 1.e-1; 578e0eea495SMark rectx->current_rate = 0; 579e0eea495SMark rectx->plotIdx = 0; 580e0eea495SMark rectx->j = 0; 581e0eea495SMark rectx->plotDt = 1.0; 582930e68a5SMark Adams rectx->plotting = PETSC_FALSE; 583e0eea495SMark rectx->use_spitzer_eta = PETSC_FALSE; 584e0eea495SMark rectx->idx = 0; 5858a6f2e61SMark Adams rectx->print_period = 10; 586*419c2fb8Smarkadams4 rectx->grid_view_idx = -1; // do not get if not needed 587e0eea495SMark /* Register the available impurity sources */ 5889566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&plist,"step",&stepSrc)); 5899566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&plist,"none",&zeroSrc)); 5909566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&plist,"pulse",&pulseSrc)); 5919566063dSJacob Faibussowitsch PetscCall(PetscStrcpy(pname,"none")); 5929566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&testlist,"none",&testNone)); 5939566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&testlist,"spitzer",&testSpitzer)); 5949566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&testlist,"stable",&testStable)); 5959566063dSJacob Faibussowitsch PetscCall(PetscStrcpy(testname,"none")); 5969566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&elist,"none",&ENone)); 5979566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&elist,"induction",&EInduction)); 5989566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&elist,"constant",&EConstant)); 5999566063dSJacob Faibussowitsch PetscCall(PetscStrcpy(ename,"constant")); 600e0eea495SMark 601d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_SELF, prefix, "Options for Runaway/seed electron model", "none"); 6029566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ex2_plot_dt", "Plotting interval", "ex2.c", rectx->plotDt, &rectx->plotDt, NULL)); 603e0eea495SMark if (rectx->plotDt < 0) rectx->plotDt = 1e30; 604e0eea495SMark if (rectx->plotDt == 0) rectx->plotDt = 1e-30; 6059566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ex2_print_period", "Plotting interval", "ex2.c", rectx->print_period, &rectx->print_period, NULL)); 6069566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ex2_grid_view_idx", "grid_view_idx", "ex2.c", rectx->grid_view_idx, &rectx->grid_view_idx, NULL)); 607*419c2fb8Smarkadams4 PetscCheck(rectx->grid_view_idx < ctx->num_grids || rectx->grid_view_idx == -1,PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG,"rectx->grid_view_idx (%" PetscInt_FMT ") >= ctx->num_grids (%" PetscInt_FMT ")",rectx->imp_idx,ctx->num_grids); 6089566063dSJacob Faibussowitsch PetscCall(PetscOptionsFList("-ex2_impurity_source_type","Name of impurity source to run","",plist,pname,pname,sizeof(pname),NULL)); 6099566063dSJacob Faibussowitsch PetscCall(PetscOptionsFList("-ex2_test_type","Name of test to run","",testlist,testname,testname,sizeof(testname),NULL)); 6109566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ex2_impurity_index", "index of sink for impurities", "none", rectx->imp_idx, &rectx->imp_idx, NULL)); 611cad9d221SBarry Smith PetscCheck((rectx->imp_idx < ctx->num_species && rectx->imp_idx >= 1) || ctx->num_species <= 1,PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG,"index of sink for impurities ions is out of range (%" PetscInt_FMT "), must be > 0 && < NS",rectx->imp_idx); 6129566063dSJacob Faibussowitsch PetscCall(PetscOptionsFList("-ex2_e_field_type","Electric field type","",elist,ename,ename,sizeof(ename),NULL)); 613e0eea495SMark rectx->Ne_ion = -ctx->charges[rectx->imp_idx]/ctx->charges[0]; 6149566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ex2_t_cold","Temperature of cold electron and ions after ionization in keV","none",rectx->T_cold,&rectx->T_cold, NULL)); 6159566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ex2_pulse_start_time","Time at which pulse happens for 'pulse' source","none",rectx->pulse_start,&rectx->pulse_start, NULL)); 6169566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ex2_pulse_width_time","Width of pulse 'pulse' source","none",rectx->pulse_width,&rectx->pulse_width, NULL)); 6179566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ex2_pulse_rate","Number density of pulse for 'pulse' source","none",rectx->pulse_rate,&rectx->pulse_rate, NULL)); 618e0eea495SMark rectx->T_cold *= 1.16e7; /* convert to Kelvin */ 6199566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ex2_ion_potential","Potential to ionize impurity (should be array) in ev","none",rectx->ion_potential,&rectx->ion_potential, NULL)); 6209566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-ex2_inductance","Inductance E feild","none",rectx->L,&rectx->L, NULL)); 6219566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ex2_connor_e_field_units","Scale Ex but Connor-Hastie E_c","none",Connor_E,&Connor_E, NULL)); 62263a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(dm_dummy, "Num electrons from ions=%g, T_cold=%10.3e, ion potential=%10.3e, E_z=%10.3e v_0=%10.3e\n",(double)rectx->Ne_ion,(double)rectx->T_cold,(double)rectx->ion_potential,(double)ctx->Ez,(double)ctx->v_0)); 623d0609cedSBarry Smith PetscOptionsEnd(); 624e0eea495SMark /* get impurity source rate function */ 6259566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(plist,pname,&rectx->impuritySrcRate)); 6263c633725SBarry Smith PetscCheck(rectx->impuritySrcRate,PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG,"No impurity source function found '%s'",pname); 6279566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(testlist,testname,&rectx->test)); 6283c633725SBarry Smith PetscCheck(rectx->test,PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG,"No test found '%s'",testname); 6299566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(elist,ename,&rectx->E)); 6303c633725SBarry Smith PetscCheck(rectx->E,PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG,"No E field function found '%s'",ename); 6319566063dSJacob Faibussowitsch PetscCall(PetscFunctionListDestroy(&plist)); 6329566063dSJacob Faibussowitsch PetscCall(PetscFunctionListDestroy(&testlist)); 6339566063dSJacob Faibussowitsch PetscCall(PetscFunctionListDestroy(&elist)); 634930e68a5SMark Adams 635a587d139SMark /* convert E from Connor-Hastie E_c units to real if doing Spitzer E */ 636a587d139SMark if (Connor_E) { 637e0eea495SMark PetscReal E = ctx->Ez, Tev = ctx->thermal_temps[0]*8.621738e-5, n = ctx->n_0*ctx->n[0]; 638e0eea495SMark CalculateE(Tev, n, ctx->lnLam, ctx->epsilon0, &E); 639e0eea495SMark ((LandauCtx*)ctx)->Ez *= E; 640e0eea495SMark } 6419566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm_dummy)); 642e0eea495SMark PetscFunctionReturn(0); 643e0eea495SMark } 644e0eea495SMark 645e0eea495SMark int main(int argc, char **argv) 646e0eea495SMark { 6478a6f2e61SMark Adams DM pack; 648*419c2fb8Smarkadams4 Vec X; 649f37ccb5fSMark Adams PetscInt dim = 2, nDMs; 650e0eea495SMark TS ts; 651e0eea495SMark Mat J; 652e0eea495SMark PetscDS prob; 653e0eea495SMark LandauCtx *ctx; 654e0eea495SMark REctx *rectx; 655a587d139SMark #if defined PETSC_USE_LOG 656a587d139SMark PetscLogStage stage; 657a587d139SMark #endif 658930e68a5SMark Adams PetscMPIInt rank; 6598fdabdddSMark Adams #if defined(PETSC_HAVE_THREADSAFETY) 6608fdabdddSMark Adams double starttime, endtime; 6618fdabdddSMark Adams #endif 6629566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL,help)); 6639566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 664930e68a5SMark Adams if (rank) { /* turn off output stuff for duplicate runs */ 665*419c2fb8Smarkadams4 PetscCall(PetscOptionsClearValue(NULL,"-ex2_dm_view")); 666*419c2fb8Smarkadams4 PetscCall(PetscOptionsClearValue(NULL,"-ex2_vec_view")); 667*419c2fb8Smarkadams4 PetscCall(PetscOptionsClearValue(NULL,"-ex2_vec_view_init")); 668*419c2fb8Smarkadams4 PetscCall(PetscOptionsClearValue(NULL,"-ex2_dm_view_init")); 6699566063dSJacob Faibussowitsch PetscCall(PetscOptionsClearValue(NULL,"-info")); /* this does not work */ 670930e68a5SMark Adams } 6719566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL, "-dim", &dim, NULL)); 672e0eea495SMark /* Create a mesh */ 6739566063dSJacob Faibussowitsch PetscCall(DMPlexLandauCreateVelocitySpace(PETSC_COMM_WORLD, dim, "", &X, &J, &pack)); 6749566063dSJacob Faibussowitsch PetscCall(DMCompositeGetNumberDM(pack,&nDMs)); 6759566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)J, "Jacobian")); 6769566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)X, "f")); 6779566063dSJacob Faibussowitsch PetscCall(DMGetApplicationContext(pack, &ctx)); 6789566063dSJacob Faibussowitsch PetscCall(DMSetUp(pack)); 679e0eea495SMark /* context */ 6809566063dSJacob Faibussowitsch PetscCall(PetscNew(&rectx)); 681e0eea495SMark ctx->data = rectx; 6829566063dSJacob Faibussowitsch PetscCall(ProcessREOptions(rectx,ctx,pack,"")); 6839566063dSJacob Faibussowitsch PetscCall(DMGetDS(pack, &prob)); 684*419c2fb8Smarkadams4 if (rectx->grid_view_idx != -1) { 685*419c2fb8Smarkadams4 Vec *XsubArray=NULL; 686*419c2fb8Smarkadams4 PetscCall(PetscMalloc(sizeof(*XsubArray)*nDMs, &XsubArray)); 6879566063dSJacob Faibussowitsch PetscCall(DMCompositeGetAccessArray(pack, X, nDMs, NULL, XsubArray)); // read only 6889566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) XsubArray[ LAND_PACK_IDX(ctx->batch_view_idx, rectx->grid_view_idx) ], rectx->grid_view_idx==0 ? "ue" : "ui")); 689*419c2fb8Smarkadams4 PetscCall(DMViewFromOptions(ctx->plex[rectx->grid_view_idx],NULL,"-ex2_dm_view")); 690*419c2fb8Smarkadams4 PetscCall(DMViewFromOptions(ctx->plex[rectx->grid_view_idx], NULL,"-ex2_dm_view_init")); 691*419c2fb8Smarkadams4 PetscCall(VecViewFromOptions(XsubArray[ LAND_PACK_IDX(ctx->batch_view_idx,rectx->grid_view_idx) ], NULL,"-ex2_vec_view_init")); // initial condition (monitor plots after step) 6929566063dSJacob Faibussowitsch PetscCall(DMCompositeRestoreAccessArray(pack, X, nDMs, NULL, XsubArray)); // read only 6939566063dSJacob Faibussowitsch PetscCall(PetscFree(XsubArray)); 694*419c2fb8Smarkadams4 } 6959566063dSJacob Faibussowitsch PetscCall(DMSetOutputSequenceNumber(pack, 0, 0.0)); 696e0eea495SMark /* Create timestepping solver context */ 6979566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_SELF,&ts)); 6989566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts,pack)); 6999566063dSJacob Faibussowitsch PetscCall(TSSetIFunction(ts,NULL,DMPlexLandauIFunction,NULL)); 7009566063dSJacob Faibussowitsch PetscCall(TSSetIJacobian(ts,J,J,DMPlexLandauIJacobian,NULL)); 7019566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts,NULL,FormSource,NULL)); 7029566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 7039566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts,X)); 7049566063dSJacob Faibussowitsch PetscCall(TSSetApplicationContext(ts, ctx)); 7059566063dSJacob Faibussowitsch PetscCall(TSMonitorSet(ts,Monitor,ctx,NULL)); 7069566063dSJacob Faibussowitsch PetscCall(TSSetPreStep(ts,PreStep)); 707e0eea495SMark rectx->Ez_initial = ctx->Ez; /* cache for induction caclulation - applied E field */ 7088594ddcfSMark Adams if (1) { /* warm up an test just DMPlexLandauIJacobian */ 709e0eea495SMark Vec vec; 710a587d139SMark PetscInt nsteps; 711a587d139SMark PetscReal dt; 7129566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("Warmup", &stage)); 7139566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stage)); 7149566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X,&vec)); 7159566063dSJacob Faibussowitsch PetscCall(VecCopy(X,vec)); 7169566063dSJacob Faibussowitsch PetscCall(TSGetMaxSteps(ts,&nsteps)); 7179566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts,&dt)); 7189566063dSJacob Faibussowitsch PetscCall(TSSetMaxSteps(ts,1)); 7199566063dSJacob Faibussowitsch PetscCall(TSSolve(ts,X)); 7209566063dSJacob Faibussowitsch PetscCall(TSSetMaxSteps(ts,nsteps)); 7219566063dSJacob Faibussowitsch PetscCall(TSSetStepNumber(ts,0)); 7229566063dSJacob Faibussowitsch PetscCall(TSSetTime(ts,0)); 7239566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts,dt)); 724a587d139SMark rectx->plotIdx = 0; 725930e68a5SMark Adams rectx->plotting = PETSC_FALSE; 7269566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 7279566063dSJacob Faibussowitsch PetscCall(VecCopy(vec,X)); 7289566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vec)); 729c3e4dd79SMark Adams PetscCall(PetscObjectStateIncrease((PetscObject)ctx->J)); 730e0eea495SMark } 731e0eea495SMark /* go */ 7329566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("Solve", &stage)); 733f37ccb5fSMark Adams ctx->stage = 0; // lets not use this stage 7348fdabdddSMark Adams #if defined(PETSC_HAVE_THREADSAFETY) 735e607c864SMark Adams ctx->stage = 1; // not set with thread safety 7368fdabdddSMark Adams #endif 737c3e4dd79SMark Adams //PetscCall(TSSetSolution(ts,X)); 7389566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stage)); 7398fdabdddSMark Adams #if defined(PETSC_HAVE_THREADSAFETY) 7408fdabdddSMark Adams starttime = MPI_Wtime(); 7418fdabdddSMark Adams #endif 742*419c2fb8Smarkadams4 #if defined(PETSC_HAVE_CUDA_NVTX) 743*419c2fb8Smarkadams4 nvtxRangePushA("ex2-TSSolve-warm"); 744*419c2fb8Smarkadams4 #endif 7459566063dSJacob Faibussowitsch PetscCall(TSSolve(ts,X)); 746*419c2fb8Smarkadams4 #if defined(PETSC_HAVE_CUDA_NVTX) 747*419c2fb8Smarkadams4 nvtxRangePop(); 748*419c2fb8Smarkadams4 #endif 7499566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 7508fdabdddSMark Adams #if defined(PETSC_HAVE_THREADSAFETY) 7518fdabdddSMark Adams endtime = MPI_Wtime(); 7528fdabdddSMark Adams ctx->times[LANDAU_EX2_TSSOLVE] += (endtime - starttime); 7538fdabdddSMark Adams #endif 754e0eea495SMark /* clean up */ 7559566063dSJacob Faibussowitsch PetscCall(DMPlexLandauDestroyVelocitySpace(&pack)); 7569566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 7579566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X)); 7589566063dSJacob Faibussowitsch PetscCall(PetscFree(rectx)); 7599566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 760b122ec5aSJacob Faibussowitsch return 0; 761e0eea495SMark } 762e0eea495SMark 763e0eea495SMark /*TEST 764a587d139SMark 7655dac466eSMark Adams testset: 766984ed092SMark Adams requires: p4est !complex double defined(PETSC_USE_DMLANDAU_2D) 7675dac466eSMark Adams output_file: output/ex2_0.out 768c3e4dd79SMark Adams args: -dm_landau_num_species_grid 1,1 -dm_landau_Ez 0 -petscspace_degree 3 -petscspace_poly_tensor 1 -dm_landau_type p4est -dm_landau_ion_masses 2 -dm_landau_ion_charges 1 -dm_landau_thermal_temps 5,5 -dm_landau_n 2,2 -dm_landau_n_0 5e19 -ts_monitor -snes_rtol 1.e-10 -snes_stol 1.e-14 -snes_monitor -snes_converged_reason -snes_max_it 10 -ts_type arkimex -ts_arkimex_type 1bee -ts_max_snes_failures -1 -ts_rtol 1e-3 -ts_dt 1.e-1 -ts_max_time 1 -ts_adapt_clip .5,1.25 -ts_max_steps 2 -ts_adapt_scale_solve_failed 0.75 -ts_adapt_time_step_increase_delay 5 -dm_landau_amr_levels_max 2,2 -ex2_impurity_source_type pulse -ex2_pulse_start_time 1e-1 -ex2_pulse_width_time 10 -ex2_pulse_rate 1e-2 -ex2_t_cold .05 -ex2_plot_dt 1e-1 -dm_refine 0 -dm_landau_gpu_assembly true -dm_landau_batch_size 2 -dm_landau_verbose 2 7695dac466eSMark Adams test: 7705dac466eSMark Adams suffix: cpu 771c3e4dd79SMark Adams args: -dm_landau_device_type cpu -ksp_type bicg -pc_type jacobi 772a587d139SMark test: 773a587d139SMark suffix: kokkos 7745dac466eSMark Adams requires: kokkos_kernels 775c3e4dd79SMark Adams args: -dm_landau_device_type kokkos -dm_mat_type aijkokkos -dm_vec_type kokkos -ksp_type bicg -pc_type jacobi 77690d163c9SMark Adams test: 77790d163c9SMark Adams suffix: cuda 7785dac466eSMark Adams requires: cuda 779c3e4dd79SMark Adams args: -dm_landau_device_type cuda -dm_mat_type aijcusparse -dm_vec_type cuda -mat_cusparse_use_cpu_solve -ksp_type bicg -pc_type jacobi 780cb25d741SMark Adams test: 781cb25d741SMark Adams suffix: kokkos_batch 782cb25d741SMark Adams requires: kokkos_kernels 783c3e4dd79SMark Adams args: -dm_landau_device_type kokkos -dm_mat_type aijkokkos -dm_vec_type kokkos -ksp_type preonly -pc_type bjkokkos -pc_bjkokkos_ksp_type bicg -pc_bjkokkos_pc_type jacobi 784bfc784b7SMark Adams test: 785bfc784b7SMark Adams suffix: kokkos_batch_coo 786bfc784b7SMark Adams requires: kokkos_kernels 787c3e4dd79SMark Adams args: -dm_landau_device_type kokkos -dm_mat_type aijkokkos -dm_vec_type kokkos -ksp_type preonly -pc_type bjkokkos -pc_bjkokkos_ksp_type bicg -pc_bjkokkos_pc_type jacobi -dm_landau_coo_assembly 78890d163c9SMark Adams 789e0eea495SMark TEST*/ 790