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> 8e0eea495SMark 9e0eea495SMark /* data for runaway electron model */ 10e0eea495SMark typedef struct REctx_struct { 118a6f2e61SMark 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; 328a6f2e61SMark Adams PetscInt print_period; 338fdabdddSMark Adams PetscInt grid_view_idx; 34e0eea495SMark } REctx; 35e0eea495SMark 36e0eea495SMark static const PetscReal kev_joul = 6.241506479963235e+15; /* 1/1000e */ 37e0eea495SMark 3890d163c9SMark Adams #define RE_CUT 3. 39e0eea495SMark /* < v, u_re * v * q > */ 40e0eea495SMark static void f0_j_re(PetscInt dim, PetscInt Nf, PetscInt NfAux, 41e0eea495SMark const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 42e0eea495SMark const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 43e0eea495SMark PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0) 44e0eea495SMark { 45e0eea495SMark PetscReal n_e = PetscRealPart(u[0]); 46e0eea495SMark if (dim==2) { 47e0eea495SMark if (x[1] > RE_CUT || x[1] < -RE_CUT) { /* simply a cutoff for REs. v_|| > 3 v(T_e) */ 48e0eea495SMark *f0 = n_e * 2.*PETSC_PI*x[0] * x[1] * constants[0]; /* n * r * v_|| * q */ 49e0eea495SMark } else { 50e0eea495SMark *f0 = 0; 51e0eea495SMark } 52e0eea495SMark } else { 53e0eea495SMark if (x[2] > RE_CUT || x[2] < -RE_CUT) { /* simply a cutoff for REs. v_|| > 3 v(T_e) */ 54e0eea495SMark *f0 = n_e * x[2] * constants[0]; 55e0eea495SMark } else { 56e0eea495SMark *f0 = 0; 57e0eea495SMark } 58e0eea495SMark } 59e0eea495SMark } 60e0eea495SMark 61e0eea495SMark /* sum < v, u*v*q > */ 62e0eea495SMark static void f0_jz_sum(PetscInt dim, PetscInt Nf, PetscInt NfAux, 63e0eea495SMark const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 64e0eea495SMark const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 658a6f2e61SMark Adams PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar q[], PetscScalar *f0) 66e0eea495SMark { 67e0eea495SMark PetscInt ii; 68e0eea495SMark f0[0] = 0; 69e0eea495SMark if (dim==2) { 708a6f2e61SMark 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 */ 71e0eea495SMark } else { 728a6f2e61SMark Adams for (ii=0;ii<Nf;ii++) f0[0] += u[ii] * x[2] * q[ii]; /* n * v_|| * q * v_0 */ 73e0eea495SMark } 74e0eea495SMark } 75e0eea495SMark 76e0eea495SMark /* < v, n_e > */ 77e0eea495SMark static void f0_n(PetscInt dim, PetscInt Nf, PetscInt NfAux, 78e0eea495SMark const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 79e0eea495SMark const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 80e0eea495SMark PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0) 81e0eea495SMark { 82e0eea495SMark PetscInt ii = (PetscInt)PetscRealPart(constants[0]); 83e0eea495SMark if (dim==2) f0[0] = 2.*PETSC_PI*x[0]*u[ii]; 84e0eea495SMark else { 85e0eea495SMark f0[0] = u[ii]; 86e0eea495SMark } 87e0eea495SMark } 88e0eea495SMark 89e0eea495SMark /* < v, n_e v_|| > */ 90e0eea495SMark static void f0_vz(PetscInt dim, PetscInt Nf, PetscInt NfAux, 91e0eea495SMark const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 92e0eea495SMark const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 93e0eea495SMark PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0) 94e0eea495SMark { 95e0eea495SMark PetscInt ii = (PetscInt)PetscRealPart(constants[0]); 96e0eea495SMark if (dim==2) f0[0] = u[ii] * 2.*PETSC_PI*x[0] * x[1]; /* n r v_|| */ 97e0eea495SMark else { 98e0eea495SMark f0[0] = u[ii] * x[2]; /* n v_|| */ 99e0eea495SMark } 100e0eea495SMark } 101e0eea495SMark 102e0eea495SMark /* < v, n_e (v-shift) > */ 103e0eea495SMark static void f0_ve_shift(PetscInt dim, PetscInt Nf, PetscInt NfAux, 104e0eea495SMark const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 105e0eea495SMark const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 106e0eea495SMark PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0) 107e0eea495SMark { 1088a6f2e61SMark Adams PetscReal vz = numConstants>0 ? PetscRealPart(constants[0]) : 0; 109e0eea495SMark 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 */ 110e0eea495SMark else { 111e0eea495SMark *f0 = u[0] * PetscSqrtReal(x[0]*x[0] + x[1]*x[1] + (x[2]-vz)*(x[2]-vz)); /* n v */ 112e0eea495SMark } 113e0eea495SMark } 114e0eea495SMark 115e0eea495SMark /* CalculateE - Calculate the electric field */ 116e0eea495SMark /* T -- Electron temperature */ 117e0eea495SMark /* n -- Electron density */ 118e0eea495SMark /* lnLambda -- */ 119e0eea495SMark /* eps0 -- */ 120e0eea495SMark /* E -- output E, input \hat E */ 121e0eea495SMark static PetscReal CalculateE(PetscReal Tev, PetscReal n, PetscReal lnLambda, PetscReal eps0, PetscReal *E) 122e0eea495SMark { 123e0eea495SMark PetscReal c,e,m; 124e0eea495SMark 125e0eea495SMark PetscFunctionBegin; 126cefb98e8SMark Adams c = 299792458.0; 127e0eea495SMark e = 1.602176e-19; 128e0eea495SMark m = 9.10938e-31; 129e0eea495SMark if (1) { 130e0eea495SMark double Ec, Ehat = *E, betath = PetscSqrtReal(2*Tev*e/(m*c*c)), j0 = Ehat * 7/(PetscSqrtReal(2)*2) * PetscPowReal(betath,3) * n * e * c; 131e0eea495SMark Ec = n*lnLambda*PetscPowReal(e,3) / (4*PETSC_PI*PetscPowReal(eps0,2)*m*c*c); 132e0eea495SMark *E = Ec; 133e0eea495SMark PetscPrintf(PETSC_COMM_WORLD, "CalculateE j0=%g Ec = %g\n",j0,Ec); 134e0eea495SMark } else { 135e0eea495SMark PetscReal Ed, vth; 136e0eea495SMark vth = PetscSqrtReal(8*Tev*e/(m*PETSC_PI)); 137e0eea495SMark Ed = n*lnLambda*PetscPowReal(e,3) / (4*PETSC_PI*PetscPowReal(eps0,2)*m*vth*vth); 138e0eea495SMark *E = Ed; 139e0eea495SMark } 140e0eea495SMark PetscFunctionReturn(0); 141e0eea495SMark } 142e0eea495SMark 143e0eea495SMark static PetscReal Spitzer(PetscReal m_e, PetscReal e, PetscReal Z, PetscReal epsilon0, PetscReal lnLam, PetscReal kTe_joules) 144e0eea495SMark { 145e0eea495SMark PetscReal Fz = (1+1.198*Z+0.222*Z*Z)/(1+2.966*Z+0.753*Z*Z), eta; 146e0eea495SMark 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); 147e0eea495SMark return eta; 148e0eea495SMark } 149e0eea495SMark 150e0eea495SMark /* */ 1518a6f2e61SMark Adams static PetscErrorCode testNone(TS ts, Vec X, PetscInt stepi, PetscReal time, PetscBool islast, LandauCtx *ctx, REctx *rectx) 152e0eea495SMark { 153e0eea495SMark PetscFunctionBeginUser; 154e0eea495SMark PetscFunctionReturn(0); 155e0eea495SMark } 156e0eea495SMark 157e0eea495SMark /* */ 1588a6f2e61SMark Adams static PetscErrorCode testSpitzer(TS ts, Vec X, PetscInt stepi, PetscReal time, PetscBool islast, LandauCtx *ctx, REctx *rectx) 159e0eea495SMark { 160e0eea495SMark PetscErrorCode ierr; 161e0eea495SMark PetscInt ii; 162e0eea495SMark PetscDS prob; 163a587d139SMark static PetscReal old_ratio = 1e10; 164a587d139SMark TSConvergedReason reason; 165a587d139SMark PetscReal J,J_re,spit_eta,Te_kev=0,E,ratio,Z,n_e,v,v2; 1668a6f2e61SMark Adams PetscScalar user[2] = {0.,ctx->charges[0]}, q[LANDAU_MAX_SPECIES],tt[LANDAU_MAX_SPECIES],vz; 167930e68a5SMark Adams PetscReal dt; 1688a6f2e61SMark Adams DM pack, plexe = ctx->plex[0], plexi = (ctx->num_grids==1) ? NULL : ctx->plex[1]; 1698fdabdddSMark Adams Vec XsubArray[LANDAU_MAX_GRIDS*LANDAU_MAX_BATCH_SZ]; 170e0eea495SMark 171e0eea495SMark PetscFunctionBeginUser; 1728a6f2e61SMark Adams if (ctx->num_species!=2) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "ctx->num_species %D != 2",ctx->num_species); 1738a6f2e61SMark Adams ierr = VecGetDM(X, &pack);CHKERRQ(ierr); 1748a6f2e61SMark Adams if (!pack) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "no DM"); 1758fdabdddSMark Adams ierr = DMCompositeGetAccessArray(pack, X, ctx->num_grids*ctx->batch_sz, NULL, XsubArray);CHKERRQ(ierr); // read only 176930e68a5SMark Adams ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr); 1778a6f2e61SMark Adams /* get current for each grid */ 1788a6f2e61SMark Adams for (ii=0;ii<ctx->num_species;ii++) q[ii] = ctx->charges[ii]; 1798a6f2e61SMark Adams ierr = DMGetDS(plexe, &prob);CHKERRQ(ierr); 1808a6f2e61SMark Adams ierr = PetscDSSetConstants(prob, 2, &q[0]);CHKERRQ(ierr); 181e0eea495SMark ierr = PetscDSSetObjective(prob, 0, &f0_jz_sum);CHKERRQ(ierr); 1828fdabdddSMark Adams //ierr = DMPlexComputeIntegralFEM(plexe,XsubArray[ctx->batch_view_idx*ctx->num_grids],tt,NULL);CHKERRQ(ierr); 1838fdabdddSMark Adams ierr = DMPlexComputeIntegralFEM(plexe,XsubArray[ LAND_PACK_IDX(ctx->batch_view_idx,0) ],tt,NULL);CHKERRQ(ierr); 184e0eea495SMark J = -ctx->n_0*ctx->v_0*PetscRealPart(tt[0]); 1858fdabdddSMark Adams if (plexi) { // add first (only) ion 1868a6f2e61SMark Adams ierr = DMGetDS(plexi, &prob);CHKERRQ(ierr); 1878a6f2e61SMark Adams ierr = PetscDSSetConstants(prob, 1, &q[1]);CHKERRQ(ierr); 1888a6f2e61SMark Adams ierr = PetscDSSetObjective(prob, 0, &f0_jz_sum);CHKERRQ(ierr); 1898fdabdddSMark Adams ierr = DMPlexComputeIntegralFEM(plexi,XsubArray[LAND_PACK_IDX(ctx->batch_view_idx,1)],tt,NULL);CHKERRQ(ierr); 1908a6f2e61SMark Adams J += -ctx->n_0*ctx->v_0*PetscRealPart(tt[0]); 1918a6f2e61SMark Adams } 192a587d139SMark /* get N_e */ 1938a6f2e61SMark Adams ierr = DMGetDS(plexe, &prob);CHKERRQ(ierr); 1948a6f2e61SMark Adams ierr = PetscDSSetConstants(prob, 1, user);CHKERRQ(ierr); 195a587d139SMark ierr = PetscDSSetObjective(prob, 0, &f0_n);CHKERRQ(ierr); 1968fdabdddSMark Adams ierr = DMPlexComputeIntegralFEM(plexe,XsubArray[LAND_PACK_IDX(ctx->batch_view_idx,0)],tt,NULL);CHKERRQ(ierr); 197a587d139SMark n_e = PetscRealPart(tt[0])*ctx->n_0; 198a587d139SMark /* Z */ 199a587d139SMark Z = -ctx->charges[1]/ctx->charges[0]; 200a587d139SMark /* remove drift */ 2018a6f2e61SMark Adams if (0) { 2028a6f2e61SMark Adams user[0] = 0; // electrons 2038a6f2e61SMark Adams ierr = DMGetDS(plexe, &prob);CHKERRQ(ierr); 2048a6f2e61SMark Adams ierr = PetscDSSetConstants(prob, 1, user);CHKERRQ(ierr); 205a587d139SMark ierr = PetscDSSetObjective(prob, 0, &f0_vz);CHKERRQ(ierr); 2068fdabdddSMark Adams ierr = DMPlexComputeIntegralFEM(plexe,XsubArray[LAND_PACK_IDX(ctx->batch_view_idx,0)],tt,NULL);CHKERRQ(ierr); 207a587d139SMark vz = ctx->n_0*PetscRealPart(tt[0])/n_e; /* non-dimensional */ 208a587d139SMark } else vz = 0; 209a587d139SMark /* thermal velocity */ 2108a6f2e61SMark Adams ierr = DMGetDS(plexe, &prob);CHKERRQ(ierr); 211a587d139SMark ierr = PetscDSSetConstants(prob, 1, &vz);CHKERRQ(ierr); 212a587d139SMark ierr = PetscDSSetObjective(prob, 0, &f0_ve_shift);CHKERRQ(ierr); 2138fdabdddSMark Adams ierr = DMPlexComputeIntegralFEM(plexe,XsubArray[LAND_PACK_IDX(ctx->batch_view_idx,0)],tt,NULL);CHKERRQ(ierr); 214a587d139SMark v = ctx->n_0*ctx->v_0*PetscRealPart(tt[0])/n_e; /* remove number density to get velocity */ 215a587d139SMark v2 = PetscSqr(v); /* use real space: m^2 / s^2 */ 216a587d139SMark Te_kev = (v2*ctx->masses[0]*PETSC_PI/8)*kev_joul; /* temperature in kev */ 2178fdabdddSMark Adams //Te_kev = ctx->thermal_temps[0]/1.1604525e7; 218a587d139SMark spit_eta = Spitzer(ctx->masses[0],-ctx->charges[0],Z,ctx->epsilon0,ctx->lnLam,Te_kev/kev_joul); /* kev --> J (kT) */ 2198a6f2e61SMark Adams if (0) { 2208a6f2e61SMark Adams ierr = DMGetDS(plexe, &prob);CHKERRQ(ierr); 2218a6f2e61SMark Adams ierr = PetscDSSetConstants(prob, 1, q);CHKERRQ(ierr); 222e0eea495SMark ierr = PetscDSSetObjective(prob, 0, &f0_j_re);CHKERRQ(ierr); 2238fdabdddSMark Adams ierr = DMPlexComputeIntegralFEM(plexe,XsubArray[LAND_PACK_IDX(ctx->batch_view_idx,0)],tt,NULL);CHKERRQ(ierr); 224a587d139SMark } else tt[0] = 0; 225e0eea495SMark J_re = -ctx->n_0*ctx->v_0*PetscRealPart(tt[0]); 2268fdabdddSMark Adams ierr = DMCompositeRestoreAccessArray(pack, X, ctx->num_grids*ctx->batch_sz, NULL, XsubArray);CHKERRQ(ierr); // read only 227a587d139SMark 22890d163c9SMark Adams if (rectx->use_spitzer_eta) { 22990d163c9SMark Adams E = ctx->Ez = spit_eta*(rectx->j-J_re); 23090d163c9SMark Adams } else { 23190d163c9SMark Adams E = ctx->Ez; /* keep real E */ 23290d163c9SMark Adams rectx->j = J; /* cache */ 23390d163c9SMark Adams } 23490d163c9SMark Adams 235e0eea495SMark ratio = E/J/spit_eta; 2368a6f2e61SMark Adams if (stepi>10 && !rectx->use_spitzer_eta && ( 2378a6f2e61SMark Adams //(old_ratio-ratio < 1.e-3 && ratio > 0.99 && ratio < 1.01) || 2388a6f2e61SMark Adams //(old_ratio-ratio < 1.e-4 && ratio > 0.98 && ratio < 1.02) || 2398a6f2e61SMark Adams (old_ratio-ratio < 1.e-6))) { 24090d163c9SMark Adams rectx->pulse_start = time + 0.98*dt; 241a587d139SMark rectx->use_spitzer_eta = PETSC_TRUE; 242a587d139SMark } 243e0eea495SMark ierr = TSGetConvergedReason(ts,&reason);CHKERRQ(ierr); 244930e68a5SMark Adams ierr = TSGetConvergedReason(ts,&reason);CHKERRQ(ierr); 24590d163c9SMark Adams if ((rectx->plotting) || stepi == 0 || reason || rectx->pulse_start == time + 0.98*dt) { 2468fdabdddSMark 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 spit_eta=%g\n",stepi,time,n_e/ctx->n_0,ctx->Ez,J,J_re,100*J_re/J, Te_kev,Z,ratio,old_ratio-ratio, rectx->use_spitzer_eta ? "using Spitzer eta*J E" : "constant E",rectx->pulse_start != time + 0.98*dt ? "normal" : "transition",spit_eta);CHKERRQ(ierr); 2478fdabdddSMark Adams if (rectx->pulse_start == time + 0.98*dt) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG,"Spitzer complete ratio=%g",ratio); 248930e68a5SMark Adams } 249e0eea495SMark old_ratio = ratio; 250e0eea495SMark PetscFunctionReturn(0); 251e0eea495SMark } 252e0eea495SMark 253e0eea495SMark static const double ppp = 2; 254e0eea495SMark static void f0_0_diff_lp(PetscInt dim, PetscInt Nf, PetscInt NfAux, 255e0eea495SMark const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 256e0eea495SMark const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 257e0eea495SMark PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0) 258e0eea495SMark { 259e0eea495SMark LandauCtx *ctx = (LandauCtx *)constants; 260e0eea495SMark REctx *rectx = (REctx*)ctx->data; 261e0eea495SMark PetscInt ii = rectx->idx, i; 262e0eea495SMark const PetscReal kT_m = ctx->k*ctx->thermal_temps[ii]/ctx->masses[ii]; /* kT/m */ 263e0eea495SMark const PetscReal n = ctx->n[ii]; 264e0eea495SMark PetscReal diff, f_maxwell, v2 = 0, theta = 2*kT_m/(ctx->v_0*ctx->v_0); /* theta = 2kT/mc^2 */ 265e0eea495SMark for (i = 0; i < dim; ++i) v2 += x[i]*x[i]; 266e0eea495SMark f_maxwell = n*PetscPowReal(PETSC_PI*theta,-1.5)*(PetscExpReal(-v2/theta)); 267e0eea495SMark diff = 2.*PETSC_PI*x[0]*(PetscRealPart(u[ii]) - f_maxwell); 268e0eea495SMark f0[0] = PetscPowReal(diff,ppp); 269e0eea495SMark } 270e0eea495SMark static void f0_0_maxwellian_lp(PetscInt dim, PetscInt Nf, PetscInt NfAux, 271e0eea495SMark const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 272e0eea495SMark const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 273e0eea495SMark PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0) 274e0eea495SMark { 275e0eea495SMark LandauCtx *ctx = (LandauCtx *)constants; 276e0eea495SMark REctx *rectx = (REctx*)ctx->data; 277e0eea495SMark PetscInt ii = rectx->idx, i; 278e0eea495SMark const PetscReal kT_m = ctx->k*ctx->thermal_temps[ii]/ctx->masses[ii]; /* kT/m */ 279e0eea495SMark const PetscReal n = ctx->n[ii]; 280e0eea495SMark PetscReal f_maxwell, v2 = 0, theta = 2*kT_m/(ctx->v_0*ctx->v_0); /* theta = 2kT/mc^2 */ 281e0eea495SMark for (i = 0; i < dim; ++i) v2 += x[i]*x[i]; 282e0eea495SMark f_maxwell = 2.*PETSC_PI*x[0] * n*PetscPowReal(PETSC_PI*theta,-1.5)*(PetscExpReal(-v2/theta)); 283e0eea495SMark f0[0] = PetscPowReal(f_maxwell,ppp); 284e0eea495SMark } 285e0eea495SMark 286e0eea495SMark /* */ 2878a6f2e61SMark Adams static PetscErrorCode testStable(TS ts, Vec X, PetscInt stepi, PetscReal time, PetscBool islast, LandauCtx *ctx, REctx *rectx) 288e0eea495SMark { 289e0eea495SMark PetscErrorCode ierr; 290e0eea495SMark PetscDS prob; 291e0eea495SMark Vec X2; 292e0eea495SMark PetscReal ediff,idiff=0,lpm0,lpm1=1; 293e0eea495SMark PetscScalar tt[LANDAU_MAX_SPECIES]; 2948a6f2e61SMark Adams DM dm, plex = ctx->plex[0]; 295e0eea495SMark 296e0eea495SMark PetscFunctionBeginUser; 297e0eea495SMark ierr = VecGetDM(X, &dm);CHKERRQ(ierr); 298e0eea495SMark ierr = DMGetDS(plex, &prob);CHKERRQ(ierr); 299e0eea495SMark ierr = VecDuplicate(X,&X2);CHKERRQ(ierr); 300e0eea495SMark ierr = VecCopy(X,X2);CHKERRQ(ierr); 301e0eea495SMark if (!rectx->X_0) { 302e0eea495SMark ierr = VecDuplicate(X,&rectx->X_0);CHKERRQ(ierr); 303e0eea495SMark ierr = VecCopy(X,rectx->X_0);CHKERRQ(ierr); 304e0eea495SMark } 305e0eea495SMark ierr = VecAXPY(X,-1.0,rectx->X_0);CHKERRQ(ierr); 306e0eea495SMark ierr = PetscDSSetConstants(prob, sizeof(LandauCtx)/sizeof(PetscScalar), (PetscScalar*)ctx);CHKERRQ(ierr); 307e0eea495SMark rectx->idx = 0; 308e0eea495SMark ierr = PetscDSSetObjective(prob, 0, &f0_0_diff_lp);CHKERRQ(ierr); 309e0eea495SMark ierr = DMPlexComputeIntegralFEM(plex,X2,tt,NULL);CHKERRQ(ierr); 310e0eea495SMark ediff = PetscPowReal(PetscRealPart(tt[0]),1./ppp); 311e0eea495SMark ierr = PetscDSSetObjective(prob, 0, &f0_0_maxwellian_lp);CHKERRQ(ierr); 312e0eea495SMark ierr = DMPlexComputeIntegralFEM(plex,X2,tt,NULL);CHKERRQ(ierr); 313e0eea495SMark lpm0 = PetscPowReal(PetscRealPart(tt[0]),1./ppp); 314e0eea495SMark if (ctx->num_species>1) { 315e0eea495SMark rectx->idx = 1; 316e0eea495SMark ierr = PetscDSSetObjective(prob, 0, &f0_0_diff_lp);CHKERRQ(ierr); 317e0eea495SMark ierr = DMPlexComputeIntegralFEM(plex,X2,tt,NULL);CHKERRQ(ierr); 318e0eea495SMark idiff = PetscPowReal(PetscRealPart(tt[0]),1./ppp); 319e0eea495SMark ierr = PetscDSSetObjective(prob, 0, &f0_0_maxwellian_lp);CHKERRQ(ierr); 320e0eea495SMark ierr = DMPlexComputeIntegralFEM(plex,X2,tt,NULL);CHKERRQ(ierr); 321e0eea495SMark lpm1 = PetscPowReal(PetscRealPart(tt[0]),1./ppp); 322e0eea495SMark } 323e0eea495SMark 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); 324e0eea495SMark /* view */ 325e0eea495SMark ierr = VecCopy(X2,X);CHKERRQ(ierr); 326e0eea495SMark ierr = VecDestroy(&X2);CHKERRQ(ierr); 327e0eea495SMark if (islast) { 328e0eea495SMark ierr = VecDestroy(&rectx->X_0);CHKERRQ(ierr); 329e0eea495SMark rectx->X_0 = NULL; 330e0eea495SMark } 331e0eea495SMark PetscFunctionReturn(0); 332e0eea495SMark } 333e0eea495SMark 334e0eea495SMark static PetscErrorCode EInduction(Vec X, Vec X_t, PetscInt step, PetscReal time, LandauCtx *ctx, PetscReal *a_E) 335e0eea495SMark { 336e0eea495SMark REctx *rectx = (REctx*)ctx->data; 337e0eea495SMark PetscErrorCode ierr; 338e0eea495SMark PetscInt ii; 339e0eea495SMark DM dm,plex; 3408a6f2e61SMark Adams PetscScalar tt[LANDAU_MAX_SPECIES], qv0[LANDAU_MAX_SPECIES]; 341e0eea495SMark PetscReal dJ_dt; 342e0eea495SMark PetscDS prob; 343e0eea495SMark 344e0eea495SMark PetscFunctionBeginUser; 3458a6f2e61SMark Adams for (ii=0;ii<ctx->num_species;ii++) qv0[ii] = ctx->charges[ii]*ctx->v_0; 346e0eea495SMark ierr = VecGetDM(X, &dm);CHKERRQ(ierr); 347e0eea495SMark ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 348e0eea495SMark ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr); 349e0eea495SMark /* get d current / dt */ 3508a6f2e61SMark Adams ierr = PetscDSSetConstants(prob, ctx->num_species, qv0);CHKERRQ(ierr); 351e0eea495SMark ierr = PetscDSSetObjective(prob, 0, &f0_jz_sum);CHKERRQ(ierr); 352e0eea495SMark if (!X_t) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "X_t"); 353e0eea495SMark ierr = DMPlexComputeIntegralFEM(plex,X_t,tt,NULL);CHKERRQ(ierr); 3548a6f2e61SMark Adams dJ_dt = -ctx->n_0*PetscRealPart(tt[0])/ctx->t_0; 355e0eea495SMark /* E induction */ 356e0eea495SMark *a_E = -rectx->L*dJ_dt + rectx->Ez_initial; 357e0eea495SMark ierr = DMDestroy(&plex);CHKERRQ(ierr); 358e0eea495SMark PetscFunctionReturn(0); 359e0eea495SMark } 360e0eea495SMark 361e0eea495SMark static PetscErrorCode EConstant(Vec X, Vec X_t, PetscInt step, PetscReal time, LandauCtx *ctx, PetscReal *a_E) 362e0eea495SMark { 363e0eea495SMark PetscFunctionBeginUser; 364e0eea495SMark *a_E = ctx->Ez; 365e0eea495SMark PetscFunctionReturn(0); 366e0eea495SMark } 367e0eea495SMark 368e0eea495SMark static PetscErrorCode ENone(Vec X, Vec X_t, PetscInt step, PetscReal time, LandauCtx *ctx, PetscReal *a_E) 369e0eea495SMark { 370e0eea495SMark PetscFunctionBeginUser; 371e0eea495SMark *a_E = 0; 372e0eea495SMark PetscFunctionReturn(0); 373e0eea495SMark } 374e0eea495SMark 375e0eea495SMark /* ------------------------------------------------------------------- */ 376e0eea495SMark /* 377e0eea495SMark FormSource - Evaluates source terms F(t). 378e0eea495SMark 379e0eea495SMark Input Parameters: 380e0eea495SMark . ts - the TS context 381e0eea495SMark . time - 382e0eea495SMark . X_dummmy - input vector 383e0eea495SMark . dummy - optional user-defined context, as set by SNESSetFunction() 384e0eea495SMark 385e0eea495SMark Output Parameter: 386e0eea495SMark . F - function vector 387e0eea495SMark */ 388e0eea495SMark static PetscErrorCode FormSource(TS ts, PetscReal ftime, Vec X_dummmy, Vec F, void *dummy) 389e0eea495SMark { 390e0eea495SMark PetscReal new_imp_rate; 391e0eea495SMark LandauCtx *ctx; 3928a6f2e61SMark Adams DM pack; 393e0eea495SMark PetscErrorCode ierr; 394e0eea495SMark REctx *rectx; 395e0eea495SMark 396e0eea495SMark PetscFunctionBeginUser; 3978a6f2e61SMark Adams ierr = TSGetDM(ts,&pack);CHKERRQ(ierr); 3988a6f2e61SMark Adams ierr = DMGetApplicationContext(pack, &ctx);CHKERRQ(ierr); 399e0eea495SMark rectx = (REctx*)ctx->data; 400e0eea495SMark /* check for impurities */ 401e0eea495SMark ierr = rectx->impuritySrcRate(ftime,&new_imp_rate,ctx);CHKERRQ(ierr); 402e0eea495SMark if (new_imp_rate != 0) { 403e0eea495SMark if (new_imp_rate != rectx->current_rate) { 404e0eea495SMark PetscInt ii; 405e0eea495SMark PetscReal dne_dt,dni_dt,tilda_ns[LANDAU_MAX_SPECIES],temps[LANDAU_MAX_SPECIES]; 4068fdabdddSMark Adams Vec globFarray[LANDAU_MAX_GRIDS*LANDAU_MAX_BATCH_SZ]; 407e0eea495SMark rectx->current_rate = new_imp_rate; 408e0eea495SMark for (ii=1;ii<LANDAU_MAX_SPECIES;ii++) tilda_ns[ii] = 0; 409e0eea495SMark for (ii=1;ii<LANDAU_MAX_SPECIES;ii++) temps[ii] = 1; 4108a6f2e61SMark Adams dni_dt = new_imp_rate /* *ctx->t_0 */; /* fully ionized immediately, no normalize, stay in non-dim */ 4118a6f2e61SMark Adams dne_dt = new_imp_rate*rectx->Ne_ion /* *ctx->t_0 */; 412e0eea495SMark tilda_ns[0] = dne_dt; tilda_ns[rectx->imp_idx] = dni_dt; 413e0eea495SMark temps[0] = rectx->T_cold; temps[rectx->imp_idx] = rectx->T_cold; 4148a6f2e61SMark 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); 4158fdabdddSMark Adams ierr = DMCompositeGetAccessArray(pack, F, ctx->num_grids*ctx->batch_sz, NULL, globFarray);CHKERRQ(ierr); 4168a6f2e61SMark Adams for (PetscInt grid=0 ; grid<ctx->num_grids ; grid++) { 417e0eea495SMark /* add it */ 4188fdabdddSMark Adams ierr = LandauAddMaxwellians(ctx->plex[grid],globFarray[ LAND_PACK_IDX(0,grid) ],ftime,temps,tilda_ns,grid,0,ctx);CHKERRQ(ierr); 4198fdabdddSMark Adams ierr = VecViewFromOptions(globFarray[ LAND_PACK_IDX(0,grid) ],NULL,"-vec_view_sources");CHKERRQ(ierr); 420e0eea495SMark } 4218fdabdddSMark Adams // Does DMCompositeRestoreAccessArray copy the data back? (no) 4228fdabdddSMark Adams ierr = DMCompositeRestoreAccessArray(pack, F, ctx->num_grids*ctx->batch_sz, NULL, globFarray);CHKERRQ(ierr); 423e0eea495SMark } 424e0eea495SMark } else { 425a587d139SMark ierr = VecZeroEntries(F);CHKERRQ(ierr); 426e0eea495SMark rectx->current_rate = 0; 427e0eea495SMark } 428e0eea495SMark PetscFunctionReturn(0); 429e0eea495SMark } 430e0eea495SMark PetscErrorCode Monitor(TS ts, PetscInt stepi, PetscReal time, Vec X, void *actx) 431e0eea495SMark { 432e0eea495SMark LandauCtx *ctx = (LandauCtx*) actx; /* user-defined application context */ 433e0eea495SMark REctx *rectx = (REctx*)ctx->data; 4348a6f2e61SMark Adams DM pack; 4358fdabdddSMark Adams Vec globXArray[LANDAU_MAX_GRIDS*LANDAU_MAX_BATCH_SZ]; 436e0eea495SMark TSConvergedReason reason; 437e0eea495SMark PetscErrorCode ierr; 438e0eea495SMark PetscFunctionBeginUser; 4398a6f2e61SMark Adams ierr = VecGetDM(X, &pack);CHKERRQ(ierr); 4408fdabdddSMark Adams ierr = DMCompositeGetAccessArray(pack, X, ctx->num_grids*ctx->batch_sz, NULL, globXArray);CHKERRQ(ierr); 441e0eea495SMark if (stepi > rectx->plotStep && rectx->plotting) { 442e0eea495SMark rectx->plotting = PETSC_FALSE; /* was doing diagnostics, now done */ 443e0eea495SMark rectx->plotIdx++; 444e0eea495SMark } 445e0eea495SMark /* view */ 446e0eea495SMark ierr = TSGetConvergedReason(ts,&reason);CHKERRQ(ierr); 447e0eea495SMark if (time/rectx->plotDt >= (PetscReal)rectx->plotIdx || reason) { 4488a6f2e61SMark Adams if ((reason || stepi==0 || rectx->plotIdx%rectx->print_period==0) && ctx->verbose > 0) { 449e0eea495SMark /* print norms */ 450e0eea495SMark ierr = LandauPrintNorms(X, stepi);CHKERRQ(ierr); 451a587d139SMark } 452930e68a5SMark Adams if (!rectx->plotting) { /* first step of possible backtracks */ 453930e68a5SMark Adams rectx->plotting = PETSC_TRUE; 454930e68a5SMark Adams /* diagnostics + change E field with Sptizer (not just a monitor) */ 4558a6f2e61SMark Adams ierr = rectx->test(ts,X,stepi,time,reason ? PETSC_TRUE : PETSC_FALSE, ctx, rectx);CHKERRQ(ierr); 456930e68a5SMark Adams } else { 457930e68a5SMark Adams PetscPrintf(PETSC_COMM_WORLD, "\t\t ERROR SKIP test spit ------\n"); 458930e68a5SMark Adams rectx->plotting = PETSC_TRUE; 459930e68a5SMark Adams } 4608fdabdddSMark Adams ierr = PetscObjectSetName((PetscObject) globXArray[ LAND_PACK_IDX(ctx->batch_view_idx,rectx->grid_view_idx) ], rectx->grid_view_idx==0 ? "ue" : "ui");CHKERRQ(ierr); 461930e68a5SMark Adams /* view, overwrite step when back tracked */ 4628a6f2e61SMark Adams ierr = DMSetOutputSequenceNumber(pack, rectx->plotIdx, time*ctx->t_0);CHKERRQ(ierr); 4638fdabdddSMark Adams ierr = VecViewFromOptions(globXArray[ LAND_PACK_IDX(ctx->batch_view_idx, rectx->grid_view_idx) ],NULL,"-vec_view");CHKERRQ(ierr); 4648fdabdddSMark Adams 465e0eea495SMark rectx->plotStep = stepi; 466930e68a5SMark Adams } else { 467930e68a5SMark Adams if (rectx->plotting) PetscPrintf(PETSC_COMM_WORLD," ERROR rectx->plotting=%D step %D\n",rectx->plotting,stepi); 468930e68a5SMark Adams /* diagnostics + change E field with Sptizer (not just a monitor) - can we lag this? */ 4698a6f2e61SMark Adams ierr = rectx->test(ts,X,stepi,time,reason ? PETSC_TRUE : PETSC_FALSE, ctx, rectx);CHKERRQ(ierr); 470e0eea495SMark } 4718fdabdddSMark Adams ierr = DMCompositeRestoreAccessArray(pack, X, ctx->num_grids*ctx->batch_sz, NULL, globXArray);CHKERRQ(ierr); 472e0eea495SMark /* parallel check */ 473930e68a5SMark Adams if (reason && ctx->verbose > 0) { 474e0eea495SMark PetscReal val,rval; 475e0eea495SMark PetscMPIInt rank; 476ffc4695bSBarry Smith ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank);CHKERRMPI(ierr); 477e0eea495SMark ierr = TSGetSolution(ts, &X);CHKERRQ(ierr); 478e0eea495SMark ierr = VecNorm(X,NORM_2,&val);CHKERRQ(ierr); 479820f2d46SBarry Smith ierr = MPIU_Allreduce(&val,&rval,1,MPIU_REAL,MPIU_MAX,PETSC_COMM_WORLD);CHKERRMPI(ierr); 480e0eea495SMark if (rval != val) { 481a587d139SMark 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); 482e0eea495SMark } else { 483e0eea495SMark ierr = PetscPrintf(PETSC_COMM_WORLD, "[%D] parallel consistency check OK\n",rank);CHKERRQ(ierr); 484e0eea495SMark } 485e0eea495SMark } 486e0eea495SMark rectx->idx = 0; 487e0eea495SMark PetscFunctionReturn(0); 488e0eea495SMark } 489e0eea495SMark 490e0eea495SMark PetscErrorCode PreStep(TS ts) 491e0eea495SMark { 492e0eea495SMark PetscErrorCode ierr; 493e0eea495SMark LandauCtx *ctx; 494e0eea495SMark REctx *rectx; 495e0eea495SMark DM dm; 496e0eea495SMark PetscInt stepi; 497e0eea495SMark PetscReal time; 498e0eea495SMark Vec X; 499e0eea495SMark 500e0eea495SMark PetscFunctionBeginUser; 501930e68a5SMark Adams /* not used */ 502e0eea495SMark ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 503e0eea495SMark ierr = TSGetTime(ts,&time);CHKERRQ(ierr); 504e0eea495SMark ierr = TSGetSolution(ts,&X);CHKERRQ(ierr); 505e0eea495SMark ierr = DMGetApplicationContext(dm, &ctx);CHKERRQ(ierr); 506e0eea495SMark rectx = (REctx*)ctx->data; 507e0eea495SMark ierr = TSGetStepNumber(ts, &stepi);CHKERRQ(ierr); 508e0eea495SMark /* update E */ 509e0eea495SMark ierr = rectx->E(X, NULL, stepi, time, ctx, &ctx->Ez);CHKERRQ(ierr); 510e0eea495SMark PetscFunctionReturn(0); 511e0eea495SMark } 512e0eea495SMark 513e0eea495SMark /* 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) */ 514e0eea495SMark static PetscErrorCode stepSrc(PetscReal time, PetscReal *rho, LandauCtx *ctx) 515e0eea495SMark { 516e0eea495SMark REctx *rectx = (REctx*)ctx->data; 517e0eea495SMark 518e0eea495SMark PetscFunctionBeginUser; 519e0eea495SMark if (time >= rectx->pulse_start) *rho = rectx->pulse_rate; 520e0eea495SMark else *rho = 0.; 521e0eea495SMark PetscFunctionReturn(0); 522e0eea495SMark } 523e0eea495SMark static PetscErrorCode zeroSrc(PetscReal time, PetscReal *rho, LandauCtx *ctx) 524e0eea495SMark { 525e0eea495SMark PetscFunctionBeginUser; 526e0eea495SMark *rho = 0.; 527e0eea495SMark PetscFunctionReturn(0); 528e0eea495SMark } 529e0eea495SMark static PetscErrorCode pulseSrc(PetscReal time, PetscReal *rho, LandauCtx *ctx) 530e0eea495SMark { 531e0eea495SMark REctx *rectx = (REctx*)ctx->data; 532e0eea495SMark 533e0eea495SMark PetscFunctionBeginUser; 534a587d139SMark 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'"); 535e0eea495SMark if (time < rectx->pulse_start || time > rectx->pulse_start + 3*rectx->pulse_width) *rho = 0; 536a587d139SMark /* else if (0) { */ 537a587d139SMark /* 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); */ 538a587d139SMark /* *rho = rectx->pulse_rate * (cycle / (stop - start)) / (1 + PetscExpReal(steep*(PetscSinReal(2*PETSC_PI*((t - start)/cycle + xi)) - PetscSinReal(2*PETSC_PI*xi)))); */ 539a587d139SMark /* } else if (0) { */ 540a587d139SMark /* double x = 2*(time - rectx->pulse_start)/(3*rectx->pulse_width) - 1; */ 541a587d139SMark /* if (x==1 || x==-1) *rho = 0; */ 542a587d139SMark /* else *rho = rectx->pulse_rate * PetscExpReal(-1/(1-x*x)); */ 543a587d139SMark /* } */ 544a587d139SMark else { 545e0eea495SMark double x = PetscSinReal((time-rectx->pulse_start)/(3*rectx->pulse_width)*2*PETSC_PI - PETSC_PI/2) + 1; /* 0:2, integrates to 1.0 */ 546e0eea495SMark *rho = rectx->pulse_rate * x / (3*rectx->pulse_width); 547a587d139SMark if (!rectx->use_spitzer_eta) rectx->use_spitzer_eta = PETSC_TRUE; /* use it next time */ 548e0eea495SMark } 549e0eea495SMark PetscFunctionReturn(0); 550e0eea495SMark } 551e0eea495SMark 552e0eea495SMark #undef __FUNCT__ 553e0eea495SMark #define __FUNCT__ "ProcessREOptions" 554e0eea495SMark static PetscErrorCode ProcessREOptions(REctx *rectx, const LandauCtx *ctx, DM dm, const char prefix[]) 555e0eea495SMark { 556e0eea495SMark PetscErrorCode ierr; 557e0eea495SMark PetscFunctionList plist = NULL, testlist = NULL, elist = NULL; 558e0eea495SMark char pname[256],testname[256],ename[256]; 559930e68a5SMark Adams DM dm_dummy; 560a587d139SMark PetscBool Connor_E = PETSC_FALSE; 561e0eea495SMark 562e0eea495SMark PetscFunctionBeginUser; 563930e68a5SMark Adams ierr = DMCreate(PETSC_COMM_WORLD,&dm_dummy);CHKERRQ(ierr); 564e0eea495SMark rectx->Ne_ion = 1; /* number of electrons given up by impurity ion */ 565e0eea495SMark rectx->T_cold = .005; /* kev */ 566e0eea495SMark rectx->ion_potential = 15; /* ev */ 567e0eea495SMark rectx->L = 2; 568e0eea495SMark rectx->X_0 = NULL; 569e0eea495SMark rectx->imp_idx = ctx->num_species - 1; /* default ionized impurity as last one */ 570a587d139SMark rectx->pulse_start = PETSC_MAX_REAL; 571e0eea495SMark rectx->pulse_width = 1; 572e0eea495SMark rectx->plotStep = PETSC_MAX_INT; 573e0eea495SMark rectx->pulse_rate = 1.e-1; 574e0eea495SMark rectx->current_rate = 0; 575e0eea495SMark rectx->plotIdx = 0; 576e0eea495SMark rectx->j = 0; 577e0eea495SMark rectx->plotDt = 1.0; 578930e68a5SMark Adams rectx->plotting = PETSC_FALSE; 579e0eea495SMark rectx->use_spitzer_eta = PETSC_FALSE; 580e0eea495SMark rectx->idx = 0; 5818a6f2e61SMark Adams rectx->print_period = 10; 5828fdabdddSMark Adams rectx->grid_view_idx = 0; 583e0eea495SMark /* Register the available impurity sources */ 584e0eea495SMark ierr = PetscFunctionListAdd(&plist,"step",&stepSrc);CHKERRQ(ierr); 585e0eea495SMark ierr = PetscFunctionListAdd(&plist,"none",&zeroSrc);CHKERRQ(ierr); 586e0eea495SMark ierr = PetscFunctionListAdd(&plist,"pulse",&pulseSrc);CHKERRQ(ierr); 587e0eea495SMark ierr = PetscStrcpy(pname,"none");CHKERRQ(ierr); 588e0eea495SMark ierr = PetscFunctionListAdd(&testlist,"none",&testNone);CHKERRQ(ierr); 589e0eea495SMark ierr = PetscFunctionListAdd(&testlist,"spitzer",&testSpitzer);CHKERRQ(ierr); 590e0eea495SMark ierr = PetscFunctionListAdd(&testlist,"stable",&testStable);CHKERRQ(ierr); 591e0eea495SMark ierr = PetscStrcpy(testname,"none");CHKERRQ(ierr); 592e0eea495SMark ierr = PetscFunctionListAdd(&elist,"none",&ENone);CHKERRQ(ierr); 593e0eea495SMark ierr = PetscFunctionListAdd(&elist,"induction",&EInduction);CHKERRQ(ierr); 594e0eea495SMark ierr = PetscFunctionListAdd(&elist,"constant",&EConstant);CHKERRQ(ierr); 595e0eea495SMark ierr = PetscStrcpy(ename,"constant");CHKERRQ(ierr); 596e0eea495SMark 597e0eea495SMark ierr = PetscOptionsBegin(PETSC_COMM_SELF, prefix, "Options for Runaway/seed electron model", "none");CHKERRQ(ierr); 5988fdabdddSMark Adams ierr = PetscOptionsReal("-ex2_plot_dt", "Plotting interval", "ex2.c", rectx->plotDt, &rectx->plotDt, NULL);CHKERRQ(ierr); 599e0eea495SMark if (rectx->plotDt < 0) rectx->plotDt = 1e30; 600e0eea495SMark if (rectx->plotDt == 0) rectx->plotDt = 1e-30; 6018fdabdddSMark Adams ierr = PetscOptionsInt("-ex2_print_period", "Plotting interval", "ex2.c", rectx->print_period, &rectx->print_period, NULL);CHKERRQ(ierr); 6028fdabdddSMark Adams ierr = PetscOptionsInt("-ex2_grid_view_idx", "grid_view_idx", "ex2.c", rectx->grid_view_idx, &rectx->grid_view_idx, NULL);CHKERRQ(ierr); 6038fdabdddSMark Adams if (rectx->grid_view_idx >= ctx->num_grids) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG,"rectx->grid_view_idx (%D) >= ctx->num_grids (%D)",rectx->imp_idx,ctx->num_grids); 604e0eea495SMark ierr = PetscOptionsFList("-ex2_impurity_source_type","Name of impurity source to run","",plist,pname,pname,sizeof(pname),NULL);CHKERRQ(ierr); 605e0eea495SMark ierr = PetscOptionsFList("-ex2_test_type","Name of test to run","",testlist,testname,testname,sizeof(testname),NULL);CHKERRQ(ierr); 606e0eea495SMark ierr = PetscOptionsInt("-ex2_impurity_index", "index of sink for impurities", "none", rectx->imp_idx, &rectx->imp_idx, NULL);CHKERRQ(ierr); 607e0eea495SMark 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); 608e0eea495SMark ierr = PetscOptionsFList("-ex2_e_field_type","Electric field type","",elist,ename,ename,sizeof(ename),NULL);CHKERRQ(ierr); 609e0eea495SMark rectx->Ne_ion = -ctx->charges[rectx->imp_idx]/ctx->charges[0]; 610e0eea495SMark 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); 611e0eea495SMark ierr = PetscOptionsReal("-ex2_pulse_start_time","Time at which pulse happens for 'pulse' source","none",rectx->pulse_start,&rectx->pulse_start, NULL);CHKERRQ(ierr); 612e0eea495SMark ierr = PetscOptionsReal("-ex2_pulse_width_time","Width of pulse 'pulse' source","none",rectx->pulse_width,&rectx->pulse_width, NULL);CHKERRQ(ierr); 613e0eea495SMark ierr = PetscOptionsReal("-ex2_pulse_rate","Number density of pulse for 'pulse' source","none",rectx->pulse_rate,&rectx->pulse_rate, NULL);CHKERRQ(ierr); 614e0eea495SMark rectx->T_cold *= 1.16e7; /* convert to Kelvin */ 615e0eea495SMark ierr = PetscOptionsReal("-ex2_ion_potential","Potential to ionize impurity (should be array) in ev","none",rectx->ion_potential,&rectx->ion_potential, NULL);CHKERRQ(ierr); 616a587d139SMark ierr = PetscOptionsReal("-ex2_inductance","Inductance E feild","none",rectx->L,&rectx->L, NULL);CHKERRQ(ierr); 617a587d139SMark ierr = PetscOptionsBool("-ex2_connor_e_field_units","Scale Ex but Connor-Hastie E_c","none",Connor_E,&Connor_E, NULL);CHKERRQ(ierr); 618930e68a5SMark 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); 619e0eea495SMark ierr = PetscOptionsEnd();CHKERRQ(ierr); 620e0eea495SMark /* get impurity source rate function */ 621e0eea495SMark ierr = PetscFunctionListFind(plist,pname,&rectx->impuritySrcRate);CHKERRQ(ierr); 622e0eea495SMark if (!rectx->impuritySrcRate) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG,"No impurity source function found '%s'",pname); 623e0eea495SMark ierr = PetscFunctionListFind(testlist,testname,&rectx->test);CHKERRQ(ierr); 624e0eea495SMark if (!rectx->test) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG,"No test found '%s'",testname); 625e0eea495SMark ierr = PetscFunctionListFind(elist,ename,&rectx->E);CHKERRQ(ierr); 626e0eea495SMark if (!rectx->E) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG,"No E field function found '%s'",ename); 627e0eea495SMark ierr = PetscFunctionListDestroy(&plist);CHKERRQ(ierr); 628e0eea495SMark ierr = PetscFunctionListDestroy(&testlist);CHKERRQ(ierr); 629e0eea495SMark ierr = PetscFunctionListDestroy(&elist);CHKERRQ(ierr); 630930e68a5SMark Adams 631a587d139SMark /* convert E from Connor-Hastie E_c units to real if doing Spitzer E */ 632a587d139SMark if (Connor_E) { 633e0eea495SMark PetscReal E = ctx->Ez, Tev = ctx->thermal_temps[0]*8.621738e-5, n = ctx->n_0*ctx->n[0]; 634e0eea495SMark CalculateE(Tev, n, ctx->lnLam, ctx->epsilon0, &E); 635e0eea495SMark ((LandauCtx*)ctx)->Ez *= E; 636e0eea495SMark } 637930e68a5SMark Adams ierr = DMDestroy(&dm_dummy);CHKERRQ(ierr); 638e0eea495SMark PetscFunctionReturn(0); 639e0eea495SMark } 640e0eea495SMark 641e0eea495SMark int main(int argc, char **argv) 642e0eea495SMark { 6438a6f2e61SMark Adams DM pack; 6448fdabdddSMark Adams Vec X,XsubArray[LANDAU_MAX_GRIDS*LANDAU_MAX_BATCH_SZ]; 645e0eea495SMark PetscErrorCode ierr; 646e0eea495SMark PetscInt dim = 2; 647e0eea495SMark TS ts; 648e0eea495SMark Mat J; 649e0eea495SMark PetscDS prob; 650e0eea495SMark LandauCtx *ctx; 651e0eea495SMark REctx *rectx; 652a587d139SMark #if defined PETSC_USE_LOG 653a587d139SMark PetscLogStage stage; 654a587d139SMark #endif 655930e68a5SMark Adams PetscMPIInt rank; 6568fdabdddSMark Adams #if defined(PETSC_HAVE_THREADSAFETY) 6578fdabdddSMark Adams double starttime, endtime; 6588fdabdddSMark Adams #endif 659e0eea495SMark ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr; 66055b25c41SPierre Jolivet ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank);CHKERRMPI(ierr); 661930e68a5SMark Adams if (rank) { /* turn off output stuff for duplicate runs */ 6628fdabdddSMark Adams ierr = PetscOptionsClearValue(NULL,"-dm_view");CHKERRQ(ierr); 6638fdabdddSMark Adams ierr = PetscOptionsClearValue(NULL,"-vec_view");CHKERRQ(ierr); 6648fdabdddSMark Adams ierr = PetscOptionsClearValue(NULL,"-dm_view_diff");CHKERRQ(ierr); 6658fdabdddSMark Adams ierr = PetscOptionsClearValue(NULL,"-vec_view_diff");CHKERRQ(ierr); 6668fdabdddSMark Adams ierr = PetscOptionsClearValue(NULL,"-dm_view_sources");CHKERRQ(ierr); 6678fdabdddSMark Adams ierr = PetscOptionsClearValue(NULL,"-vec_view_0");CHKERRQ(ierr); 6688fdabdddSMark Adams ierr = PetscOptionsClearValue(NULL,"-dm_view_0");CHKERRQ(ierr); 6698fdabdddSMark Adams ierr = PetscOptionsClearValue(NULL,"-vec_view_sources");CHKERRQ(ierr); 670930e68a5SMark Adams ierr = PetscOptionsClearValue(NULL,"-info");CHKERRQ(ierr); /* this does not work */ 671930e68a5SMark Adams } 672e0eea495SMark ierr = PetscOptionsGetInt(NULL,NULL, "-dim", &dim, NULL);CHKERRQ(ierr); 673e0eea495SMark /* Create a mesh */ 6748a6f2e61SMark Adams ierr = LandauCreateVelocitySpace(PETSC_COMM_WORLD, dim, "", &X, &J, &pack);CHKERRQ(ierr); 675930e68a5SMark Adams ierr = PetscObjectSetName((PetscObject)J, "Jacobian");CHKERRQ(ierr); 6768a6f2e61SMark Adams ierr = PetscObjectSetName((PetscObject)X, "f");CHKERRQ(ierr); 6778a6f2e61SMark Adams ierr = LandauCreateMassMatrix(pack, NULL);CHKERRQ(ierr); 6788a6f2e61SMark Adams ierr = DMGetApplicationContext(pack, &ctx);CHKERRQ(ierr); 6798a6f2e61SMark Adams ierr = DMSetUp(pack);CHKERRQ(ierr); 680e0eea495SMark /* context */ 681e0eea495SMark ierr = PetscNew(&rectx);CHKERRQ(ierr); 682e0eea495SMark ctx->data = rectx; 6838a6f2e61SMark Adams ierr = ProcessREOptions(rectx,ctx,pack,"");CHKERRQ(ierr); 6848fdabdddSMark Adams ierr = DMGetDS(pack, &prob);CHKERRQ(ierr); 6858fdabdddSMark Adams ierr = DMCompositeGetAccessArray(pack, X, ctx->num_grids*ctx->batch_sz, NULL, XsubArray);CHKERRQ(ierr); // read only 6868fdabdddSMark Adams ierr = PetscObjectSetName((PetscObject) XsubArray[ LAND_PACK_IDX(ctx->batch_view_idx, rectx->grid_view_idx) ], rectx->grid_view_idx==0 ? "ue" : "ui");CHKERRQ(ierr); 6878fdabdddSMark Adams ierr = DMViewFromOptions(ctx->plex[rectx->grid_view_idx],NULL,"-dm_view");CHKERRQ(ierr); 6888fdabdddSMark Adams ierr = DMViewFromOptions(ctx->plex[rectx->grid_view_idx], NULL,"-dm_view_0");CHKERRQ(ierr); 6898fdabdddSMark Adams ierr = VecViewFromOptions(XsubArray[ LAND_PACK_IDX(ctx->batch_view_idx,rectx->grid_view_idx) ], NULL,"-vec_view_0");CHKERRQ(ierr); // initial condition (monitor plots after step) 6908fdabdddSMark Adams ierr = DMCompositeRestoreAccessArray(pack, X, ctx->num_grids*ctx->batch_sz, NULL, XsubArray);CHKERRQ(ierr); // read only 6918fdabdddSMark Adams ierr = VecViewFromOptions(X, NULL,"-vec_view_global");CHKERRQ(ierr); // initial condition (monitor plots after step) 6928a6f2e61SMark Adams ierr = DMSetOutputSequenceNumber(pack, 0, 0.0);CHKERRQ(ierr); 693e0eea495SMark /* Create timestepping solver context */ 694e0eea495SMark ierr = TSCreate(PETSC_COMM_SELF,&ts);CHKERRQ(ierr); 6958a6f2e61SMark Adams ierr = TSSetDM(ts,pack);CHKERRQ(ierr); 696e0eea495SMark ierr = TSSetIFunction(ts,NULL,LandauIFunction,NULL);CHKERRQ(ierr); 697e0eea495SMark ierr = TSSetIJacobian(ts,J,J,LandauIJacobian,NULL);CHKERRQ(ierr); 698e0eea495SMark ierr = TSSetRHSFunction(ts,NULL,FormSource,NULL);CHKERRQ(ierr); 699e0eea495SMark ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 700e0eea495SMark ierr = TSSetSolution(ts,X);CHKERRQ(ierr); 701e0eea495SMark ierr = TSSetApplicationContext(ts, ctx);CHKERRQ(ierr); 702e0eea495SMark ierr = TSMonitorSet(ts,Monitor,ctx,NULL);CHKERRQ(ierr); 703e0eea495SMark ierr = TSSetPreStep(ts,PreStep);CHKERRQ(ierr); 7048fdabdddSMark Adams 705e0eea495SMark rectx->Ez_initial = ctx->Ez; /* cache for induction caclulation - applied E field */ 706c8a6034eSMark if (1) { /* warm up an test just LandauIJacobian */ 707e0eea495SMark Vec vec; 708a587d139SMark PetscInt nsteps; 709a587d139SMark PetscReal dt; 710e0eea495SMark ierr = PetscLogStageRegister("Warmup", &stage);CHKERRQ(ierr); 711e0eea495SMark ierr = PetscLogStagePush(stage);CHKERRQ(ierr); 712a587d139SMark ierr = VecDuplicate(X,&vec);CHKERRQ(ierr); 713a587d139SMark ierr = VecCopy(X,vec);CHKERRQ(ierr); 7141e1ea65dSPierre Jolivet ierr = TSGetMaxSteps(ts,&nsteps);CHKERRQ(ierr); 715a587d139SMark ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr); 7161e1ea65dSPierre Jolivet ierr = TSSetMaxSteps(ts,1);CHKERRQ(ierr); 7178a6f2e61SMark Adams ierr = TSSolve(ts,X);CHKERRQ(ierr); 7181e1ea65dSPierre Jolivet ierr = TSSetMaxSteps(ts,nsteps);CHKERRQ(ierr); 719a587d139SMark ierr = TSSetStepNumber(ts,0);CHKERRQ(ierr); 720a587d139SMark ierr = TSSetTime(ts,0);CHKERRQ(ierr); 721a587d139SMark ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr); 722a587d139SMark rectx->plotIdx = 0; 723930e68a5SMark Adams rectx->plotting = PETSC_FALSE; 724e0eea495SMark ierr = PetscLogStagePop();CHKERRQ(ierr); 7258a6f2e61SMark Adams ierr = VecCopy(vec,X);CHKERRQ(ierr); 726e0eea495SMark ierr = VecDestroy(&vec);CHKERRQ(ierr); 727e8d2b73aSMark Adams ctx->aux_bool = PETSC_FALSE; // flag for not a clean Jacobian 728e0eea495SMark } 729e0eea495SMark /* go */ 730a587d139SMark ierr = PetscLogStageRegister("Solve", &stage);CHKERRQ(ierr); 731c751c0a2SMark Adams ierr = PetscLogStageRegister("Landau", &ctx->stage);CHKERRQ(ierr); 7328fdabdddSMark Adams #if defined(PETSC_HAVE_THREADSAFETY) 7338fdabdddSMark Adams ctx->stage = 1; // not set with thread safty 7348fdabdddSMark Adams #endif 735c751c0a2SMark Adams ierr = TSSetSolution(ts,X);CHKERRQ(ierr); 736a587d139SMark ierr = PetscLogStagePush(stage);CHKERRQ(ierr); 7378fdabdddSMark Adams #if defined(PETSC_HAVE_THREADSAFETY) 7388fdabdddSMark Adams starttime = MPI_Wtime(); 7398fdabdddSMark Adams #endif 740e0eea495SMark ierr = TSSolve(ts,X);CHKERRQ(ierr); 741a587d139SMark ierr = PetscLogStagePop();CHKERRQ(ierr); 7428fdabdddSMark Adams #if defined(PETSC_HAVE_THREADSAFETY) 7438fdabdddSMark Adams endtime = MPI_Wtime(); 7448fdabdddSMark Adams ctx->times[LANDAU_EX2_TSSOLVE] += (endtime - starttime); 7458fdabdddSMark Adams #endif 7468fdabdddSMark Adams ierr = VecViewFromOptions(X, NULL,"-vec_view_global");CHKERRQ(ierr); 747e0eea495SMark /* clean up */ 7488a6f2e61SMark Adams ierr = LandauDestroyVelocitySpace(&pack);CHKERRQ(ierr); 749e0eea495SMark ierr = TSDestroy(&ts);CHKERRQ(ierr); 750e0eea495SMark ierr = VecDestroy(&X);CHKERRQ(ierr); 751e0eea495SMark ierr = PetscFree(rectx);CHKERRQ(ierr); 752e0eea495SMark ierr = PetscFinalize(); 753e0eea495SMark return ierr; 754e0eea495SMark } 755e0eea495SMark 756e0eea495SMark /*TEST 757a587d139SMark 758*5dac466eSMark Adams testset: 759*5dac466eSMark Adams requires: p4est !complex double 760*5dac466eSMark Adams output_file: output/ex2_0.out 761*5dac466eSMark 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 -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_batch_size 2 762*5dac466eSMark Adams test: 763*5dac466eSMark Adams suffix: cpu 764*5dac466eSMark Adams args: -dm_landau_device_type cpu 765a587d139SMark test: 766a587d139SMark suffix: kokkos 767*5dac466eSMark Adams requires: kokkos_kernels 768*5dac466eSMark Adams args: -dm_landau_device_type kokkos -dm_mat_type aijkokkos -dm_vec_type kokkos 76990d163c9SMark Adams test: 77090d163c9SMark Adams suffix: cuda 771*5dac466eSMark Adams requires: cuda 772*5dac466eSMark Adams args: -dm_landau_device_type cuda -dm_mat_type aijcusparse -dm_vec_type cuda -mat_cusparse_use_cpu_solve 77390d163c9SMark Adams 774e0eea495SMark TEST*/ 775