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