xref: /petsc/src/ts/utils/dmplexlandau/plexland.c (revision e0eea49528ddcce32a775a862994a11c74fcac87)
1*e0eea495SMark #include <petsc/private/dmpleximpl.h>   /*I "petscdmplex.h" I*/
2*e0eea495SMark #include <petsclandau.h>                /*I "petsclandau.h"   I*/
3*e0eea495SMark #include <petscts.h>
4*e0eea495SMark #include <petscdmforest.h>
5*e0eea495SMark 
6*e0eea495SMark /* Landau collision operator */
7*e0eea495SMark #define PETSC_THREAD_SYNC
8*e0eea495SMark #define PETSC_DEVICE_FUNC_DECL static
9*e0eea495SMark #include "land_kernel.h"
10*e0eea495SMark 
11*e0eea495SMark #define LANDAU_VL  1
12*e0eea495SMark static PetscErrorCode LandauPointDataCreate(PetscReal **IPData, PetscInt dim, PetscInt nip, PetscInt Ns)
13*e0eea495SMark {
14*e0eea495SMark   PetscErrorCode  ierr, d, s, jj, nip_pad = LANDAU_VL*(nip/LANDAU_VL + !!(nip%LANDAU_VL)), pnt_sz = (dim + Ns*(1+dim));
15*e0eea495SMark   PetscReal       *pdata;
16*e0eea495SMark   PetscFunctionBegin;
17*e0eea495SMark   ierr = PetscMalloc(nip_pad*pnt_sz*sizeof(PetscReal),IPData);CHKERRQ(ierr);
18*e0eea495SMark   /* pad with zeros in case we vectorize into this */
19*e0eea495SMark   for (jj=nip, pdata = *IPData + nip*pnt_sz; jj < nip_pad; jj++, pdata += pnt_sz){
20*e0eea495SMark     LandauPointData *fplpt = (LandauPointData*)pdata; /* [dim + NS*(1+dim)] */
21*e0eea495SMark     for (d=0;d<dim;d++) fplpt->crd[d] = -1;
22*e0eea495SMark     for (s=0;s<Ns;s++) {
23*e0eea495SMark       fplpt->fdf[s].f = 0;
24*e0eea495SMark       for (d=0;d<dim;d++) fplpt->fdf[s].df[d] = 0;
25*e0eea495SMark     }
26*e0eea495SMark   }
27*e0eea495SMark   PetscFunctionReturn(0);
28*e0eea495SMark }
29*e0eea495SMark 
30*e0eea495SMark static PetscErrorCode LandauPointDataDestroy(PetscReal *IPData)
31*e0eea495SMark {
32*e0eea495SMark   PetscErrorCode   ierr;
33*e0eea495SMark   PetscFunctionBegin;
34*e0eea495SMark   ierr = PetscFree(IPData);CHKERRQ(ierr);
35*e0eea495SMark   PetscFunctionReturn(0);
36*e0eea495SMark }
37*e0eea495SMark /* ------------------------------------------------------------------- */
38*e0eea495SMark /*
39*e0eea495SMark   LandauFormJacobian_Internal - Evaluates Jacobian matrix.
40*e0eea495SMark 
41*e0eea495SMark   Input Parameters:
42*e0eea495SMark   .  globX - input vector
43*e0eea495SMark   .  actx - optional user-defined context
44*e0eea495SMark   .  dim - dimension
45*e0eea495SMark 
46*e0eea495SMark   Output Parameters:
47*e0eea495SMark   .  J0acP - Jacobian matrix filled, not created
48*e0eea495SMark */
49*e0eea495SMark PetscErrorCode LandauFormJacobian_Internal(Vec a_X, Mat JacP, const PetscInt dim, void *a_ctx)
50*e0eea495SMark {
51*e0eea495SMark   LandauCtx         *ctx = (LandauCtx*)a_ctx;
52*e0eea495SMark   PetscErrorCode    ierr;
53*e0eea495SMark   PetscInt          cStart, cEnd, elemMatSize;
54*e0eea495SMark   DM                plex = NULL;
55*e0eea495SMark   PetscDS           prob;
56*e0eea495SMark   PetscSection      section,globsection;
57*e0eea495SMark   PetscScalar       *elemMat;
58*e0eea495SMark   PetscInt          numCells,totDim,ej,Nq,*Nbf,*Ncf,Nb,Ncx,Nf,d,f,fieldA,Nip,nip_pad,ipdata_sz;
59*e0eea495SMark   PetscQuadrature   quad;
60*e0eea495SMark   PetscTabulation   *Tf;
61*e0eea495SMark   PetscReal         *wiGlob, nu_alpha[LANDAU_MAX_SPECIES], nu_beta[LANDAU_MAX_SPECIES];
62*e0eea495SMark   const PetscReal   *quadWeights;
63*e0eea495SMark   PetscReal         *IPData,*invJ,*invJ_a;
64*e0eea495SMark   PetscReal         invMass[LANDAU_MAX_SPECIES],Eq_m[LANDAU_MAX_SPECIES],m_0=ctx->m_0; /* normalize mass -- not needed! */
65*e0eea495SMark   PetscLogDouble    flops;
66*e0eea495SMark   Vec               locX;
67*e0eea495SMark 
68*e0eea495SMark   PetscFunctionBegin;
69*e0eea495SMark   PetscValidHeaderSpecific(a_X,VEC_CLASSID,1);
70*e0eea495SMark   PetscValidHeaderSpecific(JacP,MAT_CLASSID,2);
71*e0eea495SMark   PetscValidPointer(ctx,4);
72*e0eea495SMark 
73*e0eea495SMark   ierr = PetscLogEventBegin(ctx->events[1],0,0,0,0);CHKERRQ(ierr);
74*e0eea495SMark   ierr = DMConvert(ctx->dmv, DMPLEX, &plex);CHKERRQ(ierr);
75*e0eea495SMark   ierr = DMCreateLocalVector(plex, &locX);CHKERRQ(ierr);
76*e0eea495SMark   ierr = VecZeroEntries(locX);CHKERRQ(ierr); /* zero BCs so don't set */
77*e0eea495SMark   ierr = DMGlobalToLocalBegin(plex, a_X, INSERT_VALUES, locX);CHKERRQ(ierr);
78*e0eea495SMark   ierr = DMGlobalToLocalEnd  (plex, a_X, INSERT_VALUES, locX);CHKERRQ(ierr);
79*e0eea495SMark   ierr = DMPlexGetHeightStratum(plex, 0, &cStart, &cEnd);CHKERRQ(ierr);
80*e0eea495SMark   ierr = DMGetLocalSection(plex, &section);CHKERRQ(ierr);
81*e0eea495SMark   ierr = DMGetGlobalSection(plex, &globsection);CHKERRQ(ierr);
82*e0eea495SMark   ierr = DMGetDS(plex, &prob);CHKERRQ(ierr);
83*e0eea495SMark   ierr = PetscDSGetTabulation(prob, &Tf);CHKERRQ(ierr); // Bf, &Df
84*e0eea495SMark   ierr = PetscDSGetDimensions(prob, &Nbf);CHKERRQ(ierr); Nb = Nbf[0]; /* number of vertices*S */
85*e0eea495SMark   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);         if (Nf!=ctx->num_species) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Nf %D != S",Nf);
86*e0eea495SMark   ierr = PetscDSGetComponents(prob, &Ncf);CHKERRQ(ierr); Ncx = Ncf[0]; if (Ncx!=1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Nc %D != 1",Ncx);
87*e0eea495SMark   for (fieldA=0;fieldA<Nf;fieldA++) {
88*e0eea495SMark     invMass[fieldA] = m_0/ctx->masses[fieldA];
89*e0eea495SMark     Eq_m[fieldA] = -ctx->Ez * ctx->t_0 * ctx->charges[fieldA] / (ctx->v_0 * ctx->masses[fieldA]); /* normalize dimensionless */
90*e0eea495SMark     if (dim==2) Eq_m[fieldA] *=  2 * PETSC_PI; /* add the 2pi term that is not in Landau */
91*e0eea495SMark     nu_alpha[fieldA] = PetscSqr(ctx->charges[fieldA]/m_0)*m_0/ctx->masses[fieldA];
92*e0eea495SMark     nu_beta[fieldA] = PetscSqr(ctx->charges[fieldA]/ctx->epsilon0)*ctx->lnLam / (8*PETSC_PI) * ctx->t_0*ctx->n_0/PetscPowReal(ctx->v_0,3);
93*e0eea495SMark   }
94*e0eea495SMark   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
95*e0eea495SMark   numCells = cEnd - cStart;
96*e0eea495SMark   ierr = PetscFEGetQuadrature(ctx->fe[0], &quad);CHKERRQ(ierr);
97*e0eea495SMark   ierr = PetscQuadratureGetData(quad, NULL, NULL, &Nq, NULL, &quadWeights);CHKERRQ(ierr);
98*e0eea495SMark   if (Nb!=Nq) SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Nb!=Nq %D %D over integration or simplices? Tf[0]->Nb=%D dim=%D",Nb,Nq,Tf[0]->Nb,dim);
99*e0eea495SMark   if (Nq >LANDAU_MAX_NQ) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG,"Order too high. Nq = %D > LANDAU_MAX_NQ (%D)",Nq,LANDAU_MAX_NQ);
100*e0eea495SMark   Nip = numCells*Nq;
101*e0eea495SMark   nip_pad = LANDAU_VL*(Nip/LANDAU_VL + !!(Nip%LANDAU_VL));
102*e0eea495SMark   flops = (PetscLogDouble)numCells*(PetscLogDouble)Nq*(PetscLogDouble)(5*dim*dim*Nf*Nf + 165);
103*e0eea495SMark   ierr = MatZeroEntries(JacP);CHKERRQ(ierr);
104*e0eea495SMark   elemMatSize = totDim*totDim;
105*e0eea495SMark   {
106*e0eea495SMark     static int         cc = 0;
107*e0eea495SMark     PetscScalar        uu[LANDAU_MAX_SPECIES],u_x[LANDAU_MAX_SPECIES*LANDAU_DIM];
108*e0eea495SMark     /* collect f data */
109*e0eea495SMark     if (ctx->verbose > 2 || (ctx->verbose > 0 && cc++ == 0)) {
110*e0eea495SMark       PetscInt N,Nloc;
111*e0eea495SMark       ierr = MatGetSize(JacP,&N,NULL);CHKERRQ(ierr);
112*e0eea495SMark       ierr = VecGetSize(locX,&Nloc);CHKERRQ(ierr);
113*e0eea495SMark       ierr = PetscPrintf(PETSC_COMM_WORLD,"[%D]%s: %D IPs, %D cells, %s elements, totDim=%D, Nb=%D, Nq=%D, elemMatSize=%D, dim=%D, Tab: Nb=%D Nf=%D Np=%D cdim=%D N=%D N_loc=%D\n",
114*e0eea495SMark                          0,"FormLandau",Nq*numCells,numCells,ctx->simplex ? "SIMPLEX" : "TENSOR", totDim, Nb, Nq, elemMatSize, dim, Tf[0]->Nb, Nf, Tf[0]->Np, Tf[0]->cdim, N, Nloc);CHKERRQ(ierr);
115*e0eea495SMark     }
116*e0eea495SMark     ierr = LandauPointDataCreate(&IPData, dim, Nq*numCells, Nf);CHKERRQ(ierr);
117*e0eea495SMark     ipdata_sz = (dim + Nf*(1+dim));
118*e0eea495SMark     ierr = PetscMalloc3(elemMatSize,&elemMat,nip_pad,&wiGlob,nip_pad*dim*dim,&invJ_a);CHKERRQ(ierr);
119*e0eea495SMark     /* cache geometry and x, f and df/dx at IPs */
120*e0eea495SMark     for (ej = 0, invJ = invJ_a ; ej < numCells; ++ej, invJ += Nq*dim*dim) {
121*e0eea495SMark       PetscReal    vj[LANDAU_MAX_NQ*LANDAU_DIM],detJj[LANDAU_MAX_NQ], Jdummy[LANDAU_MAX_NQ*LANDAU_DIM*LANDAU_DIM];
122*e0eea495SMark       PetscInt     qj,f;
123*e0eea495SMark       PetscScalar *coef = NULL;
124*e0eea495SMark       ierr = DMPlexComputeCellGeometryFEM(plex, cStart+ej, quad, vj, Jdummy, invJ, detJj);CHKERRQ(ierr);
125*e0eea495SMark       ierr = DMPlexVecGetClosure(plex, section, locX, cStart+ej, NULL, &coef);CHKERRQ(ierr);
126*e0eea495SMark       /* create point data for cell i for Landau tensor: x, f(x), grad f(x) */
127*e0eea495SMark       for (qj = 0; qj < Nq; ++qj) {
128*e0eea495SMark         PetscInt       gidx = (ej*Nq + qj);
129*e0eea495SMark         LandauPointData  *pnt_data = (LandauPointData*)(IPData + gidx*ipdata_sz);
130*e0eea495SMark         PetscScalar    refSpaceDer[LANDAU_DIM];
131*e0eea495SMark         PetscInt       dOffset = 0, fOffset = 0;
132*e0eea495SMark         for (d = 0; d < dim; ++d) pnt_data->crd[d] = vj[qj * dim + d]; /* coordinate */
133*e0eea495SMark         wiGlob[gidx] = detJj[qj] * quadWeights[qj];
134*e0eea495SMark         if (dim==2) wiGlob[gidx] *= pnt_data->crd[0];  /* cylindrical coordinate, w/o 2pi */
135*e0eea495SMark         /* get u & du (EvaluateFieldJets) */
136*e0eea495SMark         for (f = 0; f < Nf; ++f) {
137*e0eea495SMark           const PetscReal *Bq = &Tf[f]->T[0][qj*Nb];
138*e0eea495SMark           const PetscReal *Dq = &Tf[f]->T[1][qj*Nb*dim];
139*e0eea495SMark           PetscInt         b, e;
140*e0eea495SMark           uu[fOffset] = 0.0;
141*e0eea495SMark           for (d = 0; d < dim; ++d) refSpaceDer[d] = 0.0;
142*e0eea495SMark           for (b = 0; b < Nb; ++b) {
143*e0eea495SMark             const PetscInt    cidx = b;
144*e0eea495SMark             uu[fOffset] += Bq[cidx]*coef[dOffset+cidx];
145*e0eea495SMark             for (d = 0; d < dim; ++d) refSpaceDer[d] += Dq[cidx*dim+d]*coef[dOffset+cidx];
146*e0eea495SMark           }
147*e0eea495SMark           for (d = 0; d < dim; ++d) {
148*e0eea495SMark             for (e = 0, u_x[fOffset*dim+d] = 0.0; e < dim; ++e) { // should add directly to point data here!!!
149*e0eea495SMark               u_x[fOffset*dim+d] += invJ[qj * dim * dim + e*dim+d]*refSpaceDer[e];
150*e0eea495SMark             }
151*e0eea495SMark           }
152*e0eea495SMark           fOffset += 1;
153*e0eea495SMark           dOffset += Nb;
154*e0eea495SMark         }
155*e0eea495SMark         /* copy to IPDataLocal */
156*e0eea495SMark         for (f=0;f<Nf;f++) {
157*e0eea495SMark           pnt_data->fdf[f].f = PetscRealPart(uu[f]);
158*e0eea495SMark           for (d = 0; d < dim; ++d) pnt_data->fdf[f].df[d] = PetscRealPart(u_x[f*dim+d]);
159*e0eea495SMark         }
160*e0eea495SMark       } /* q */
161*e0eea495SMark       ierr = DMPlexVecRestoreClosure(plex, section, locX, cStart+ej, NULL, &coef);CHKERRQ(ierr);
162*e0eea495SMark     } /* e */
163*e0eea495SMark   }
164*e0eea495SMark   ierr = DMRestoreLocalVector(plex, &locX);CHKERRQ(ierr);
165*e0eea495SMark   ierr = PetscLogEventEnd(ctx->events[1],0,0,0,0);CHKERRQ(ierr);
166*e0eea495SMark 
167*e0eea495SMark   /* outer element loop j is like a regular assembly loop */
168*e0eea495SMark   if (ctx->deviceType == LANDAU_CUDA) {
169*e0eea495SMark #if defined(PETSC_HAVE_CUDA)
170*e0eea495SMark     ierr = LandauCUDAJacobian(plex,Nq,nu_alpha,nu_beta,invMass,Eq_m,IPData,wiGlob,invJ_a,ctx->subThreadBlockSize,ctx->events,ctx->quarter3DDomain,JacP);CHKERRQ(ierr);
171*e0eea495SMark #else
172*e0eea495SMark     SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG,"-landau_device_type %s not built","cuda");
173*e0eea495SMark #endif
174*e0eea495SMark   } else if (ctx->deviceType == LANDAU_KOKKOS) {
175*e0eea495SMark #if defined(PETSC_HAVE_KOKKOS)
176*e0eea495SMark     ierr = LandauKokkosJacobian(plex,Nq,nu_alpha,nu_beta,invMass,Eq_m,IPData,wiGlob,invJ_a,ctx->subThreadBlockSize,ctx->events,ctx->quarter3DDomain,JacP);CHKERRQ(ierr);
177*e0eea495SMark #else
178*e0eea495SMark     SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG,"-landau_device_type %s not built","kokkos");
179*e0eea495SMark #endif
180*e0eea495SMark   } else { /* CPU version */
181*e0eea495SMark     for (ej = cStart, invJ = invJ_a; ej < cEnd; ++ej, invJ += Nq*dim*dim) {
182*e0eea495SMark       PetscInt     qj;
183*e0eea495SMark       ierr = PetscLogEventBegin(ctx->events[3],0,0,0,0);CHKERRQ(ierr);
184*e0eea495SMark       ierr = PetscMemzero(elemMat, totDim *totDim * sizeof(PetscScalar));CHKERRQ(ierr);
185*e0eea495SMark       ierr = PetscLogEventEnd(ctx->events[3],0,0,0,0);CHKERRQ(ierr);
186*e0eea495SMark       for (qj = 0; qj < Nq; ++qj) {
187*e0eea495SMark         PetscReal       g2[1][LANDAU_MAX_SUB_THREAD_BLOCKS][LANDAU_MAX_SPECIES][LANDAU_DIM], g3[1][LANDAU_MAX_SUB_THREAD_BLOCKS][LANDAU_MAX_SPECIES][LANDAU_DIM][LANDAU_DIM];
188*e0eea495SMark         const PetscInt  nip = numCells*Nq, jpidx = Nq*(ej-cStart) + qj, one = 1, zero = 0; /* length of inner global interation, outer integration point */
189*e0eea495SMark 
190*e0eea495SMark         ierr = PetscLogEventBegin(ctx->events[4],0,0,0,0);CHKERRQ(ierr);
191*e0eea495SMark         ierr = PetscLogFlops(flops);CHKERRQ(ierr);
192*e0eea495SMark         landau_inner_integral(zero, one, zero, one, zero, nip, 1, jpidx, Nf, dim, IPData, wiGlob, &invJ[qj*dim*dim], nu_alpha, nu_beta, invMass, Eq_m, ctx->quarter3DDomain, Nq, Nb, qj, qj+1, Tf[0]->T[0], Tf[0]->T[1], elemMat, g2, g3, ej);
193*e0eea495SMark         ierr = PetscLogEventEnd(ctx->events[4],0,0,0,0);CHKERRQ(ierr);
194*e0eea495SMark       } /* qj loop */
195*e0eea495SMark       /* assemble matrix */
196*e0eea495SMark       ierr = PetscLogEventBegin(ctx->events[6],0,0,0,0);CHKERRQ(ierr);
197*e0eea495SMark       ierr = DMPlexMatSetClosure(plex, section, globsection, JacP, ej, elemMat, ADD_VALUES);CHKERRQ(ierr);
198*e0eea495SMark       ierr = PetscLogEventEnd(ctx->events[6],0,0,0,0);CHKERRQ(ierr);
199*e0eea495SMark 
200*e0eea495SMark       if (ej==-1) {
201*e0eea495SMark         ierr = PetscPrintf(PETSC_COMM_SELF, "CPU Element matrix\n");CHKERRQ(ierr);
202*e0eea495SMark         for (d = 0; d < totDim; ++d){
203*e0eea495SMark           for (f = 0; f < totDim; ++f) PetscPrintf(PETSC_COMM_SELF," %17.9e",  PetscRealPart(elemMat[d*totDim + f]));
204*e0eea495SMark            PetscPrintf(PETSC_COMM_SELF,"\n");CHKERRQ(ierr);
205*e0eea495SMark         }
206*e0eea495SMark       }
207*e0eea495SMark     } /* ej cells loop, not cuda */
208*e0eea495SMark   }
209*e0eea495SMark   // PetscSleep(2); exit(13);
210*e0eea495SMark 
211*e0eea495SMark   /* assemble matrix or vector */
212*e0eea495SMark   ierr = PetscLogEventBegin(ctx->events[7],0,0,0,0);CHKERRQ(ierr);
213*e0eea495SMark   ierr = MatAssemblyBegin(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
214*e0eea495SMark   ierr = MatAssemblyEnd(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
215*e0eea495SMark   ierr = MatScale(JacP, -1.0);CHKERRQ(ierr); /* The code reflect the papers: du/dt = C, whereas PETSc use the form G(u) = du/dt - C(u) = 0 */
216*e0eea495SMark   ierr = PetscLogEventEnd(ctx->events[7],0,0,0,0);CHKERRQ(ierr);
217*e0eea495SMark 
218*e0eea495SMark   /* clean up */
219*e0eea495SMark   ierr = PetscFree3(elemMat,wiGlob,invJ_a);CHKERRQ(ierr);
220*e0eea495SMark   ierr = DMDestroy(&plex);CHKERRQ(ierr);
221*e0eea495SMark   ierr = LandauPointDataDestroy(IPData);CHKERRQ(ierr);
222*e0eea495SMark   PetscFunctionReturn(0);
223*e0eea495SMark }
224*e0eea495SMark 
225*e0eea495SMark #if defined(LANDAU_ADD_BCS)
226*e0eea495SMark static void zero_bc(PetscInt dim, PetscInt Nf, PetscInt NfAux,
227*e0eea495SMark                     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
228*e0eea495SMark                     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
229*e0eea495SMark                     PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar uexact[])
230*e0eea495SMark {
231*e0eea495SMark   uexact[0] = 0;
232*e0eea495SMark }
233*e0eea495SMark #endif
234*e0eea495SMark 
235*e0eea495SMark #define MATVEC2(__a,__x,__p) {int i,j; for (i=0.; i<2; i++) {__p[i] = 0; for (j=0.; j<2; j++) __p[i] += __a[i][j]*__x[j]; }}
236*e0eea495SMark static void CircleInflate(PetscReal r1, PetscReal r2, PetscReal r0, PetscInt num_sections, PetscReal x, PetscReal y,
237*e0eea495SMark                           PetscReal *outX, PetscReal *outY)
238*e0eea495SMark {
239*e0eea495SMark   PetscReal rr = PetscSqrtReal(x*x + y*y), outfact, efact;
240*e0eea495SMark   if (rr < r1 + 1.e-8) {
241*e0eea495SMark     *outX = x; *outY = y;
242*e0eea495SMark   } else {
243*e0eea495SMark     const PetscReal xy[2] = {x,y}, sinphi=y/rr, cosphi=x/rr;
244*e0eea495SMark     PetscReal       cth,sth,xyprime[2],Rth[2][2],rotcos,newrr;
245*e0eea495SMark     if (num_sections==2) {
246*e0eea495SMark       rotcos = 0.70710678118654;
247*e0eea495SMark       outfact = 1.5; efact = 2.5;
248*e0eea495SMark       /* rotate normalized vector into [-pi/4,pi/4) */
249*e0eea495SMark       if (sinphi >= 0.) {         /* top cell, -pi/2 */
250*e0eea495SMark         cth = 0.707106781186548; sth = -0.707106781186548;
251*e0eea495SMark       } else {                    /* bottom cell -pi/8 */
252*e0eea495SMark         cth = 0.707106781186548; sth = .707106781186548;
253*e0eea495SMark       }
254*e0eea495SMark     } else if (num_sections==3) {
255*e0eea495SMark       rotcos = 0.86602540378443;
256*e0eea495SMark       outfact = 1.5; efact = 2.5;
257*e0eea495SMark       /* rotate normalized vector into [-pi/6,pi/6) */
258*e0eea495SMark       if (sinphi >= 0.5) {         /* top cell, -pi/3 */
259*e0eea495SMark         cth = 0.5; sth = -0.866025403784439;
260*e0eea495SMark       } else if (sinphi >= -.5) {  /* mid cell 0 */
261*e0eea495SMark         cth = 1.; sth = .0;
262*e0eea495SMark       } else { /* bottom cell +pi/3 */
263*e0eea495SMark         cth = 0.5; sth = 0.866025403784439;
264*e0eea495SMark       }
265*e0eea495SMark     } else if (num_sections==4) {
266*e0eea495SMark       rotcos = 0.9238795325112;
267*e0eea495SMark       outfact = 1.5; efact = 3;
268*e0eea495SMark       /* rotate normalized vector into [-pi/8,pi/8) */
269*e0eea495SMark       if (sinphi >= 0.707106781186548) {         /* top cell, -3pi/8 */
270*e0eea495SMark         cth = 0.38268343236509; sth = -0.923879532511287;
271*e0eea495SMark       } else if (sinphi >= 0.) {                 /* mid top cell -pi/8 */
272*e0eea495SMark         cth = 0.923879532511287; sth = -.38268343236509;
273*e0eea495SMark       } else if (sinphi >= -0.707106781186548) { /* mid bottom cell + pi/8 */
274*e0eea495SMark         cth = 0.923879532511287; sth = 0.38268343236509;
275*e0eea495SMark       } else {                                   /* bottom cell + 3pi/8 */
276*e0eea495SMark         cth = 0.38268343236509; sth = .923879532511287;
277*e0eea495SMark       }
278*e0eea495SMark     } else {
279*e0eea495SMark       cth = 0.; sth = 0.; rotcos = 0; efact = 0;
280*e0eea495SMark     }
281*e0eea495SMark     Rth[0][0] = cth; Rth[0][1] =-sth;
282*e0eea495SMark     Rth[1][0] = sth; Rth[1][1] = cth;
283*e0eea495SMark     MATVEC2(Rth,xy,xyprime);
284*e0eea495SMark     if (num_sections==2) {
285*e0eea495SMark       newrr = xyprime[0]/rotcos;
286*e0eea495SMark     } else {
287*e0eea495SMark       PetscReal newcosphi=xyprime[0]/rr, rin = r1, rout = rr - rin;
288*e0eea495SMark       PetscReal routmax = r0*rotcos/newcosphi - rin, nroutmax = r0 - rin, routfrac = rout/routmax;
289*e0eea495SMark       newrr = rin + routfrac*nroutmax;
290*e0eea495SMark     }
291*e0eea495SMark     *outX = cosphi*newrr; *outY = sinphi*newrr;
292*e0eea495SMark     /* grade */
293*e0eea495SMark     PetscReal fact,tt,rs,re, rr = PetscSqrtReal(PetscSqr(*outX) + PetscSqr(*outY));
294*e0eea495SMark     if (rr > r2) { rs = r2; re = r0; fact = outfact;} /* outer zone */
295*e0eea495SMark     else {         rs = r1; re = r2; fact = efact;} /* electron zone */
296*e0eea495SMark     tt = (rs + PetscPowReal((rr - rs)/(re - rs),fact) * (re-rs)) / rr;
297*e0eea495SMark     *outX *= tt;
298*e0eea495SMark     *outY *= tt;
299*e0eea495SMark   }
300*e0eea495SMark }
301*e0eea495SMark 
302*e0eea495SMark static PetscErrorCode GeometryDMLandau(DM base, PetscInt point, PetscInt dim, const PetscReal abc[], PetscReal xyz[], void *a_ctx)
303*e0eea495SMark {
304*e0eea495SMark   LandauCtx   *ctx = (LandauCtx*)a_ctx;
305*e0eea495SMark   PetscReal   r = abc[0], z = abc[1];
306*e0eea495SMark   if (ctx->inflate) {
307*e0eea495SMark     PetscReal absR, absZ;
308*e0eea495SMark     absR = PetscAbsReal(r);
309*e0eea495SMark     absZ = PetscAbsReal(z);
310*e0eea495SMark     CircleInflate(ctx->i_radius,ctx->e_radius,ctx->radius,ctx->num_sections,absR,absZ,&absR,&absZ);
311*e0eea495SMark     r = (r > 0) ? absR : -absR;
312*e0eea495SMark     z = (z > 0) ? absZ : -absZ;
313*e0eea495SMark   }
314*e0eea495SMark   xyz[0] = r;
315*e0eea495SMark   xyz[1] = z;
316*e0eea495SMark   if (dim==3) xyz[2] = abc[2];
317*e0eea495SMark 
318*e0eea495SMark   PetscFunctionReturn(0);
319*e0eea495SMark }
320*e0eea495SMark 
321*e0eea495SMark static PetscErrorCode ErrorIndicator_Simple(PetscInt dim, PetscReal volume, PetscReal x[], PetscInt Nc, const PetscInt Nf[], const PetscScalar u[], const PetscScalar u_x[], PetscReal *error, void *actx)
322*e0eea495SMark {
323*e0eea495SMark   PetscReal err = 0.0;
324*e0eea495SMark   PetscInt  f = *(PetscInt*)actx, j;
325*e0eea495SMark   PetscFunctionBegin;
326*e0eea495SMark   for (j = 0; j < dim; ++j) {
327*e0eea495SMark     err += PetscSqr(PetscRealPart(u_x[f*dim+j]));
328*e0eea495SMark   }
329*e0eea495SMark   err = PetscRealPart(u[f]); /* just use rho */
330*e0eea495SMark   *error = volume * err; /* * (ctx->axisymmetric ? 2.*PETSC_PI * r : 1); */
331*e0eea495SMark   PetscFunctionReturn(0);
332*e0eea495SMark }
333*e0eea495SMark 
334*e0eea495SMark static PetscErrorCode LandauDMCreateVMesh(MPI_Comm comm, const PetscInt dim, const char prefix[], LandauCtx *ctx, DM *dm)
335*e0eea495SMark {
336*e0eea495SMark   PetscErrorCode ierr;
337*e0eea495SMark   PetscReal      radius = ctx->radius;
338*e0eea495SMark   size_t         len;
339*e0eea495SMark   char           fname[128] = ""; /* we can add a file if we want */
340*e0eea495SMark 
341*e0eea495SMark   PetscFunctionBegin;
342*e0eea495SMark   /* create DM */
343*e0eea495SMark   ierr = PetscStrlen(fname, &len);CHKERRQ(ierr);
344*e0eea495SMark   if (len) {
345*e0eea495SMark     PetscInt dim2;
346*e0eea495SMark     ierr = DMPlexCreateFromFile(comm, fname, ctx->interpolate, dm);CHKERRQ(ierr);
347*e0eea495SMark     ierr = DMGetDimension(*dm, &dim2);CHKERRQ(ierr);
348*e0eea495SMark   } else {    /* p4est, quads */
349*e0eea495SMark     /* Create plex mesh of Landau domain */
350*e0eea495SMark     if (!ctx->sphere) {
351*e0eea495SMark       PetscInt       cells[] = {2,2,2};
352*e0eea495SMark       PetscReal      lo[] = {-radius,-radius,-radius}, hi[] = {radius,radius,radius};
353*e0eea495SMark       DMBoundaryType periodicity[3] = {DM_BOUNDARY_NONE, dim==2 ? DM_BOUNDARY_NONE : DM_BOUNDARY_NONE, DM_BOUNDARY_NONE};
354*e0eea495SMark       if (dim==2) { lo[0] = 0; cells[0] = 1; }
355*e0eea495SMark       else if (ctx->quarter3DDomain) { lo[0] = lo[1] = 0; cells[0] = cells[1] = 2; }
356*e0eea495SMark       ierr = DMPlexCreateBoxMesh(comm, dim, PETSC_FALSE, cells, lo, hi, periodicity, PETSC_TRUE, dm);CHKERRQ(ierr);
357*e0eea495SMark       ierr = DMLocalizeCoordinates(*dm);CHKERRQ(ierr); /* needed for periodic */
358*e0eea495SMark       if (dim==3) {ierr = PetscObjectSetName((PetscObject) *dm, "cube");CHKERRQ(ierr);}
359*e0eea495SMark       else {ierr = PetscObjectSetName((PetscObject) *dm, "half-plane");CHKERRQ(ierr);}
360*e0eea495SMark     } else if (dim==2) {
361*e0eea495SMark       PetscInt       numCells,cells[16][4],i,j;
362*e0eea495SMark       PetscInt       numVerts;
363*e0eea495SMark       PetscReal      inner_radius1 = ctx->i_radius, inner_radius2 = ctx->e_radius;
364*e0eea495SMark       PetscReal      *flatCoords = NULL;
365*e0eea495SMark       PetscInt       *flatCells = NULL, *pcell;
366*e0eea495SMark       if (ctx->num_sections==2) {
367*e0eea495SMark #if 1
368*e0eea495SMark         numCells = 5;
369*e0eea495SMark         numVerts = 10;
370*e0eea495SMark         int cells2[][4] = { {0,1,4,3},
371*e0eea495SMark                             {1,2,5,4},
372*e0eea495SMark                             {3,4,7,6},
373*e0eea495SMark                             {4,5,8,7},
374*e0eea495SMark                             {6,7,8,9} };
375*e0eea495SMark         for (i = 0; i < numCells; i++) for (j = 0; j < 4; j++) cells[i][j] = cells2[i][j];
376*e0eea495SMark         ierr = PetscMalloc2(numVerts * 2, &flatCoords, numCells * 4, &flatCells);CHKERRQ(ierr);
377*e0eea495SMark         {
378*e0eea495SMark           PetscReal (*coords)[2] = (PetscReal (*) [2]) flatCoords;
379*e0eea495SMark           for (j = 0; j < numVerts-1; j++) {
380*e0eea495SMark             PetscReal z, r, theta = -PETSC_PI/2 + (j%3) * PETSC_PI/2;
381*e0eea495SMark             PetscReal rad = (j >= 6) ? inner_radius1 : (j >= 3) ? inner_radius2 : ctx->radius;
382*e0eea495SMark             z = rad * PetscSinReal(theta);
383*e0eea495SMark             coords[j][1] = z;
384*e0eea495SMark             r = rad * PetscCosReal(theta);
385*e0eea495SMark             coords[j][0] = r;
386*e0eea495SMark           }
387*e0eea495SMark           coords[numVerts-1][0] = coords[numVerts-1][1] = 0;
388*e0eea495SMark         }
389*e0eea495SMark #else
390*e0eea495SMark         numCells = 4;
391*e0eea495SMark         numVerts = 8;
392*e0eea495SMark         static int     cells2[][4] = {{0,1,2,3},
393*e0eea495SMark                                      {4,5,1,0},
394*e0eea495SMark                                      {5,6,2,1},
395*e0eea495SMark                                      {6,7,3,2}};
396*e0eea495SMark         for (i = 0; i < numCells; i++) for (j = 0; j < 4; j++) cells[i][j] = cells2[i][j];
397*e0eea495SMark         ierr = PetscMalloc2(numVerts * 2, &flatCoords, numCells * 4, &flatCells);CHKERRQ(ierr);
398*e0eea495SMark         {
399*e0eea495SMark           PetscReal (*coords)[2] = (PetscReal (*) [2]) flatCoords;
400*e0eea495SMark           PetscInt j;
401*e0eea495SMark           for (j = 0; j < 8; j++) {
402*e0eea495SMark             PetscReal z, r;
403*e0eea495SMark             PetscReal theta = -PETSC_PI/2 + (j%4) * PETSC_PI/3.;
404*e0eea495SMark             PetscReal rad = ctx->radius * ((j < 4) ? 0.5 : 1.0);
405*e0eea495SMark             z = rad * PetscSinReal(theta);
406*e0eea495SMark             coords[j][1] = z;
407*e0eea495SMark             r = rad * PetscCosReal(theta);
408*e0eea495SMark             coords[j][0] = r;
409*e0eea495SMark           }
410*e0eea495SMark         }
411*e0eea495SMark #endif
412*e0eea495SMark       } else if (ctx->num_sections==3) {
413*e0eea495SMark         numCells = 7;
414*e0eea495SMark         numVerts = 12;
415*e0eea495SMark         int cells2[][4] = { {0,1,5,4},
416*e0eea495SMark                             {1,2,6,5},
417*e0eea495SMark                             {2,3,7,6},
418*e0eea495SMark                             {4,5,9,8},
419*e0eea495SMark                             {5,6,10,9},
420*e0eea495SMark                             {6,7,11,10},
421*e0eea495SMark                             {8,9,10,11} };
422*e0eea495SMark         for (i = 0; i < numCells; i++) for (j = 0; j < 4; j++) cells[i][j] = cells2[i][j];
423*e0eea495SMark         ierr = PetscMalloc2(numVerts * 2, &flatCoords, numCells * 4, &flatCells);CHKERRQ(ierr);
424*e0eea495SMark         {
425*e0eea495SMark           PetscReal (*coords)[2] = (PetscReal (*) [2]) flatCoords;
426*e0eea495SMark           for (j = 0; j < numVerts; j++) {
427*e0eea495SMark             PetscReal z, r, theta = -PETSC_PI/2 + (j%4) * PETSC_PI/3;
428*e0eea495SMark             PetscReal rad = (j >= 8) ? inner_radius1 : (j >= 4) ? inner_radius2 : ctx->radius;
429*e0eea495SMark             z = rad * PetscSinReal(theta);
430*e0eea495SMark             coords[j][1] = z;
431*e0eea495SMark             r = rad * PetscCosReal(theta);
432*e0eea495SMark             coords[j][0] = r;
433*e0eea495SMark           }
434*e0eea495SMark         }
435*e0eea495SMark       } else if (ctx->num_sections==4) {
436*e0eea495SMark         numCells = 10;
437*e0eea495SMark         numVerts = 16;
438*e0eea495SMark         int cells2[][4] = { {0,1,6,5},
439*e0eea495SMark                             {1,2,7,6},
440*e0eea495SMark                             {2,3,8,7},
441*e0eea495SMark                             {3,4,9,8},
442*e0eea495SMark                             {5,6,11,10},
443*e0eea495SMark                             {6,7,12,11},
444*e0eea495SMark                             {7,8,13,12},
445*e0eea495SMark                             {8,9,14,13},
446*e0eea495SMark                             {10,11,12,15},
447*e0eea495SMark                             {12,13,14,15}};
448*e0eea495SMark         for (i = 0; i < numCells; i++) for (j = 0; j < 4; j++) cells[i][j] = cells2[i][j];
449*e0eea495SMark         ierr = PetscMalloc2(numVerts * 2, &flatCoords, numCells * 4, &flatCells);CHKERRQ(ierr);
450*e0eea495SMark         {
451*e0eea495SMark           PetscReal (*coords)[2] = (PetscReal (*) [2]) flatCoords;
452*e0eea495SMark           for (j = 0; j < numVerts-1; j++) {
453*e0eea495SMark             PetscReal z, r, theta = -PETSC_PI/2 + (j%5) * PETSC_PI/4;
454*e0eea495SMark             PetscReal rad = (j >= 10) ? inner_radius1 : (j >= 5) ? inner_radius2 : ctx->radius;
455*e0eea495SMark             z = rad * PetscSinReal(theta);
456*e0eea495SMark             coords[j][1] = z;
457*e0eea495SMark             r = rad * PetscCosReal(theta);
458*e0eea495SMark             coords[j][0] = r;
459*e0eea495SMark           }
460*e0eea495SMark           coords[numVerts-1][0] = coords[numVerts-1][1] = 0;
461*e0eea495SMark         }
462*e0eea495SMark       } else {
463*e0eea495SMark         numCells = 0;
464*e0eea495SMark         numVerts = 0;
465*e0eea495SMark       }
466*e0eea495SMark       for (j = 0, pcell = flatCells; j < numCells; j++, pcell += 4) {
467*e0eea495SMark         pcell[0] = cells[j][0]; pcell[1] = cells[j][1];
468*e0eea495SMark         pcell[2] = cells[j][2]; pcell[3] = cells[j][3];
469*e0eea495SMark       }
470*e0eea495SMark       ierr = DMPlexCreateFromCellListPetsc(comm,2,numCells,numVerts,4,ctx->interpolate,flatCells,2,flatCoords,dm);CHKERRQ(ierr);
471*e0eea495SMark       ierr = PetscFree2(flatCoords,flatCells);CHKERRQ(ierr);
472*e0eea495SMark       ierr = PetscObjectSetName((PetscObject) *dm, "semi-circle");CHKERRQ(ierr);
473*e0eea495SMark     } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Velocity space meshes does not support cubed sphere");
474*e0eea495SMark   }
475*e0eea495SMark   ierr = PetscObjectSetOptionsPrefix((PetscObject)*dm,prefix);CHKERRQ(ierr);
476*e0eea495SMark 
477*e0eea495SMark   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); /* Plex refine */
478*e0eea495SMark 
479*e0eea495SMark   { /* p4est? */
480*e0eea495SMark     char      convType[256];
481*e0eea495SMark     PetscBool flg;
482*e0eea495SMark     ierr = PetscOptionsBegin(PETSC_COMM_WORLD, prefix, "Mesh conversion options", "DMPLEX");CHKERRQ(ierr);
483*e0eea495SMark     ierr = PetscOptionsFList("-dm_landau_type","Convert DMPlex to another format (should not be Plex!)","ex6f.c",DMList,DMPLEX,convType,256,&flg);CHKERRQ(ierr);
484*e0eea495SMark     ierr = PetscOptionsEnd();
485*e0eea495SMark     if (flg) {
486*e0eea495SMark       DM dmforest;
487*e0eea495SMark       ierr = DMConvert(*dm,convType,&dmforest);CHKERRQ(ierr);
488*e0eea495SMark       if (dmforest) {
489*e0eea495SMark         PetscBool isForest;
490*e0eea495SMark         ierr = PetscObjectSetOptionsPrefix((PetscObject)dmforest,prefix);CHKERRQ(ierr);
491*e0eea495SMark         ierr = DMIsForest(dmforest,&isForest);CHKERRQ(ierr);
492*e0eea495SMark         if (isForest) {
493*e0eea495SMark           if (ctx->sphere && ctx->inflate) {
494*e0eea495SMark             ierr = DMForestSetBaseCoordinateMapping(dmforest,GeometryDMLandau,ctx);CHKERRQ(ierr);
495*e0eea495SMark           }
496*e0eea495SMark           ierr = DMDestroy(dm);CHKERRQ(ierr);
497*e0eea495SMark           *dm = dmforest;
498*e0eea495SMark           ctx->errorIndicator = ErrorIndicator_Simple; /* flag for Forest */
499*e0eea495SMark         } else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_USER, "Converted to non Forest?");
500*e0eea495SMark       } else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_USER, "Convert failed?");
501*e0eea495SMark     }
502*e0eea495SMark   }
503*e0eea495SMark   ierr = PetscObjectSetName((PetscObject) *dm, "Mesh");CHKERRQ(ierr);
504*e0eea495SMark   PetscFunctionReturn(0);
505*e0eea495SMark }
506*e0eea495SMark 
507*e0eea495SMark static PetscErrorCode SetupDS(DM dm, PetscInt dim, LandauCtx *ctx)
508*e0eea495SMark {
509*e0eea495SMark   PetscErrorCode  ierr;
510*e0eea495SMark   PetscInt        ii;
511*e0eea495SMark   PetscFunctionBegin;
512*e0eea495SMark   for (ii=0;ii<ctx->num_species;ii++) {
513*e0eea495SMark     char     buf[256];
514*e0eea495SMark     if (ii==0) ierr = PetscSNPrintf(buf, 256, "e");
515*e0eea495SMark     else {ierr = PetscSNPrintf(buf, 256, "i%D", ii);CHKERRQ(ierr);}
516*e0eea495SMark     /* Setup Discretization - FEM */
517*e0eea495SMark     ierr = PetscFECreateDefault(PetscObjectComm((PetscObject) dm), dim, 1, ctx->simplex, NULL, PETSC_DECIDE, &ctx->fe[ii]);CHKERRQ(ierr);
518*e0eea495SMark     ierr = PetscObjectSetName((PetscObject) ctx->fe[ii], buf);CHKERRQ(ierr);
519*e0eea495SMark     ierr = DMSetField(dm, ii, NULL, (PetscObject) ctx->fe[ii]);CHKERRQ(ierr);
520*e0eea495SMark   }
521*e0eea495SMark   ierr = DMCreateDS(dm);CHKERRQ(ierr);
522*e0eea495SMark   if (1) {
523*e0eea495SMark     PetscInt        ii;
524*e0eea495SMark     PetscSection    section;
525*e0eea495SMark     ierr = DMGetSection(dm, &section);CHKERRQ(ierr);
526*e0eea495SMark     for (ii=0;ii<ctx->num_species;ii++){
527*e0eea495SMark       char buf[256];
528*e0eea495SMark       if (ii==0) ierr = PetscSNPrintf(buf, 256, "se");
529*e0eea495SMark       else ierr = PetscSNPrintf(buf, 256, "si%D", ii);
530*e0eea495SMark       ierr = PetscSectionSetComponentName(section, ii, 0, buf);CHKERRQ(ierr);
531*e0eea495SMark     }
532*e0eea495SMark   }
533*e0eea495SMark   PetscFunctionReturn(0);
534*e0eea495SMark }
535*e0eea495SMark 
536*e0eea495SMark /* Define a Maxwellian function for testing out the operator. */
537*e0eea495SMark 
538*e0eea495SMark  /* Using cartesian velocity space coordinates, the particle */
539*e0eea495SMark  /* density, [1/m^3], is defined according to */
540*e0eea495SMark 
541*e0eea495SMark  /* $$ n=\int_{R^3} dv^3 \left(\frac{m}{2\pi T}\right)^{3/2}\exp [- mv^2/(2T)] $$ */
542*e0eea495SMark 
543*e0eea495SMark  /* Using some constant, c, we normalize the velocity vector into a */
544*e0eea495SMark  /* dimensionless variable according to v=c*x. Thus the density, $n$, becomes */
545*e0eea495SMark 
546*e0eea495SMark  /* $$ n=\int_{R^3} dx^3 \left(\frac{mc^2}{2\pi T}\right)^{3/2}\exp [- mc^2/(2T)*x^2] $$ */
547*e0eea495SMark 
548*e0eea495SMark  /* Defining $\theta=2T/mc^2$, we thus find that the probability density */
549*e0eea495SMark  /* for finding the particle within the interval in a box dx^3 around x is */
550*e0eea495SMark 
551*e0eea495SMark  /* f(x;\theta)=\left(\frac{1}{\pi\theta}\right)^{3/2} \exp [ -x^2/\theta ] */
552*e0eea495SMark 
553*e0eea495SMark typedef struct {
554*e0eea495SMark   LandauCtx   *ctx;
555*e0eea495SMark   PetscReal kT_m;
556*e0eea495SMark   PetscReal n;
557*e0eea495SMark   PetscReal shift;
558*e0eea495SMark } MaxwellianCtx;
559*e0eea495SMark 
560*e0eea495SMark static PetscErrorCode maxwellian(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf_dummy, PetscScalar *u, void *actx)
561*e0eea495SMark {
562*e0eea495SMark   MaxwellianCtx *mctx = (MaxwellianCtx*)actx;
563*e0eea495SMark   LandauCtx     *ctx = mctx->ctx;
564*e0eea495SMark   PetscInt      i;
565*e0eea495SMark   PetscReal     v2 = 0, theta = 2*mctx->kT_m/(ctx->v_0*ctx->v_0); /* theta = 2kT/mc^2 */
566*e0eea495SMark   PetscFunctionBegin;
567*e0eea495SMark   /* compute the exponents, v^2 */
568*e0eea495SMark   for (i = 0; i < dim; ++i) v2 += x[i]*x[i];
569*e0eea495SMark   /* evaluate the Maxwellian */
570*e0eea495SMark   u[0] = mctx->n*PetscPowReal(PETSC_PI*theta,-1.5)*(PetscExpReal(-v2/theta));
571*e0eea495SMark   if (mctx->shift!=0.) {
572*e0eea495SMark     v2 = 0;
573*e0eea495SMark     for (i = 0; i < dim-1; ++i) v2 += x[i]*x[i];
574*e0eea495SMark     v2 += (x[dim-1]-mctx->shift)*(x[dim-1]-mctx->shift);
575*e0eea495SMark     /* evaluate the shifted Maxwellian */
576*e0eea495SMark     u[0] += mctx->n*PetscPowReal(PETSC_PI*theta,-1.5)*(PetscExpReal(-v2/theta));
577*e0eea495SMark   }
578*e0eea495SMark   PetscFunctionReturn(0);
579*e0eea495SMark }
580*e0eea495SMark 
581*e0eea495SMark /*@
582*e0eea495SMark  LandauAddMaxwellians - Add a Maxwellian distribution to a state
583*e0eea495SMark 
584*e0eea495SMark  Collective on X
585*e0eea495SMark 
586*e0eea495SMark  Input Parameters:
587*e0eea495SMark  .   dm - The mesh
588*e0eea495SMark  +   time - Current time
589*e0eea495SMark  -   temps - Temperatures of each species
590*e0eea495SMark  .   ns - Number density of each species
591*e0eea495SMark  +   actx - Landau context
592*e0eea495SMark 
593*e0eea495SMark  Output Parameter:
594*e0eea495SMark  .   X  - The state
595*e0eea495SMark 
596*e0eea495SMark .keywords: mesh
597*e0eea495SMark .seealso: LandauCreateVelocitySpace()
598*e0eea495SMark @*/
599*e0eea495SMark PetscErrorCode LandauAddMaxwellians(DM dm, Vec X, PetscReal time, PetscReal temps[], PetscReal ns[], void *actx)
600*e0eea495SMark {
601*e0eea495SMark   LandauCtx      *ctx = (LandauCtx*)actx;
602*e0eea495SMark   PetscErrorCode (*initu[LANDAU_MAX_SPECIES])(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar [], void *);
603*e0eea495SMark   PetscErrorCode ierr,ii;
604*e0eea495SMark   PetscInt       dim;
605*e0eea495SMark   MaxwellianCtx  *mctxs[LANDAU_MAX_SPECIES], data[LANDAU_MAX_SPECIES];
606*e0eea495SMark 
607*e0eea495SMark   PetscFunctionBegin;
608*e0eea495SMark   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
609*e0eea495SMark   if (!ctx) { ierr = DMGetApplicationContext(dm, &ctx);CHKERRQ(ierr); }
610*e0eea495SMark   for (ii=0;ii<ctx->num_species;ii++) {
611*e0eea495SMark     mctxs[ii] = &data[ii];
612*e0eea495SMark     data[ii].ctx = ctx;
613*e0eea495SMark     data[ii].kT_m = ctx->k*temps[ii]/ctx->masses[ii]; /* kT/m */
614*e0eea495SMark     data[ii].n = ns[ii];
615*e0eea495SMark     initu[ii] = maxwellian;
616*e0eea495SMark     data[ii].shift = 0;
617*e0eea495SMark   }
618*e0eea495SMark   data[0].shift = ctx->electronShift;
619*e0eea495SMark   /* need to make ADD_ALL_VALUES work - TODO */
620*e0eea495SMark   ierr = DMProjectFunction(dm, time, initu, (void**)mctxs, INSERT_ALL_VALUES, X);CHKERRQ(ierr);
621*e0eea495SMark   PetscFunctionReturn(0);
622*e0eea495SMark }
623*e0eea495SMark 
624*e0eea495SMark /*
625*e0eea495SMark  LandauSetInitialCondition - Addes Maxwellians with context
626*e0eea495SMark 
627*e0eea495SMark Collective on X
628*e0eea495SMark 
629*e0eea495SMark  Input Parameters:
630*e0eea495SMark  .   dm - The mesh
631*e0eea495SMark  +   actx - Landau context with T and n
632*e0eea495SMark 
633*e0eea495SMark  Output Parameter:
634*e0eea495SMark  .   X  - The state
635*e0eea495SMark 
636*e0eea495SMark  Level: beginner
637*e0eea495SMark 
638*e0eea495SMark .keywords: mesh
639*e0eea495SMark .seealso: LandauCreateVelocitySpace(), LandauAddMaxwellians()
640*e0eea495SMark */
641*e0eea495SMark static PetscErrorCode LandauSetInitialCondition(DM dm, Vec X, void *actx)
642*e0eea495SMark {
643*e0eea495SMark   LandauCtx        *ctx = (LandauCtx*)actx;
644*e0eea495SMark   PetscErrorCode ierr;
645*e0eea495SMark   PetscFunctionBegin;
646*e0eea495SMark   if (!ctx) { ierr = DMGetApplicationContext(dm, &ctx);CHKERRQ(ierr); }
647*e0eea495SMark   ierr = VecZeroEntries(X);CHKERRQ(ierr);
648*e0eea495SMark   ierr = LandauAddMaxwellians(dm, X, 0.0, ctx->thermal_temps, ctx->n, ctx);CHKERRQ(ierr);
649*e0eea495SMark   PetscFunctionReturn(0);
650*e0eea495SMark }
651*e0eea495SMark 
652*e0eea495SMark static PetscErrorCode adaptToleranceFEM(PetscFE fem, Vec sol, PetscReal refineTol[], PetscReal coarsenTol[], PetscInt type, LandauCtx *ctx, DM *newDM)
653*e0eea495SMark {
654*e0eea495SMark   DM               dm, plex, adaptedDM = NULL;
655*e0eea495SMark   PetscDS          prob;
656*e0eea495SMark   PetscBool        isForest;
657*e0eea495SMark   PetscQuadrature  quad;
658*e0eea495SMark   PetscInt         Nq, *Nb, cStart, cEnd, c, dim, qj, k;
659*e0eea495SMark   DMLabel          adaptLabel = NULL;
660*e0eea495SMark   PetscErrorCode   ierr;
661*e0eea495SMark 
662*e0eea495SMark   PetscFunctionBegin;
663*e0eea495SMark   ierr = VecGetDM(sol, &dm);CHKERRQ(ierr);
664*e0eea495SMark   ierr = DMCreateDS(dm);CHKERRQ(ierr);
665*e0eea495SMark   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
666*e0eea495SMark   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
667*e0eea495SMark   ierr = DMIsForest(dm, &isForest);CHKERRQ(ierr);
668*e0eea495SMark   ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr);
669*e0eea495SMark   ierr = DMPlexGetHeightStratum(plex,0,&cStart,&cEnd);CHKERRQ(ierr);
670*e0eea495SMark   ierr = DMLabelCreate(PETSC_COMM_SELF,"adapt",&adaptLabel);CHKERRQ(ierr);
671*e0eea495SMark   ierr = PetscFEGetQuadrature(fem, &quad);CHKERRQ(ierr);
672*e0eea495SMark   ierr = PetscQuadratureGetData(quad, NULL, NULL, &Nq, NULL, NULL);CHKERRQ(ierr);
673*e0eea495SMark   if (Nq >LANDAU_MAX_NQ) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG,"Order too high. Nq = %D > LANDAU_MAX_NQ (%D)",Nq,LANDAU_MAX_NQ);
674*e0eea495SMark   ierr = PetscDSGetDimensions(prob, &Nb);CHKERRQ(ierr);
675*e0eea495SMark   if (type==4) {
676*e0eea495SMark     for (c = cStart; c < cEnd; c++) {
677*e0eea495SMark       ierr = DMLabelSetValue(adaptLabel, c, DM_ADAPT_REFINE);CHKERRQ(ierr);
678*e0eea495SMark     }
679*e0eea495SMark     ierr = PetscInfo1(sol, "Phase:%s: Uniform refinement\n","adaptToleranceFEM");
680*e0eea495SMark   } else if (type==2) {
681*e0eea495SMark     PetscInt  rCellIdx[8], eCellIdx[64], iCellIdx[64], eMaxIdx = -1, iMaxIdx = -1, nr = 0, nrmax = (dim==3 && !ctx->quarter3DDomain) ? 8 : 2;
682*e0eea495SMark     PetscReal minRad = 1.e100, r, eMinRad = 1.e100, iMinRad = 1.e100;
683*e0eea495SMark     for (c = 0; c < 64; c++) { eCellIdx[c] = iCellIdx[c] = -1; }
684*e0eea495SMark     for (c = cStart; c < cEnd; c++) {
685*e0eea495SMark       PetscReal    tt, v0[LANDAU_MAX_NQ*3], detJ[LANDAU_MAX_NQ];
686*e0eea495SMark       ierr = DMPlexComputeCellGeometryFEM(plex, c, quad, v0, NULL, NULL, detJ);CHKERRQ(ierr);
687*e0eea495SMark       for (qj = 0; qj < Nq; ++qj) {
688*e0eea495SMark         tt = PetscSqr(v0[dim*qj+0]) + PetscSqr(v0[dim*qj+1]) + PetscSqr(((dim==3) ? v0[dim*qj+2] : 0));
689*e0eea495SMark         r = PetscSqrtReal(tt);
690*e0eea495SMark         if (r < minRad - 1.e-6) {
691*e0eea495SMark           minRad = r;
692*e0eea495SMark           nr = 0;
693*e0eea495SMark           rCellIdx[nr++]= c;
694*e0eea495SMark           ierr = PetscInfo4(sol, "\t\tPhase: adaptToleranceFEM Found first inner r=%e, cell %D, qp %D/%D\n", r, c, qj+1, Nq);CHKERRQ(ierr);
695*e0eea495SMark         } else if ((r-minRad) < 1.e-8 && nr < nrmax) {
696*e0eea495SMark           for (k=0;k<nr;k++) if (c == rCellIdx[k]) break;
697*e0eea495SMark           if (k==nr) {
698*e0eea495SMark             rCellIdx[nr++]= c;
699*e0eea495SMark             ierr = PetscInfo5(sol, "\t\t\tPhase: adaptToleranceFEM Found another inner r=%e, cell %D, qp %D/%D, d=%e\n", r, c, qj+1, Nq, r-minRad);CHKERRQ(ierr);
700*e0eea495SMark           }
701*e0eea495SMark         }
702*e0eea495SMark         if (ctx->sphere) {
703*e0eea495SMark           if ((tt=r-ctx->e_radius) > 0) {
704*e0eea495SMark             PetscInfo2(sol, "\t\t\t %D cell r=%g\n",c,tt);
705*e0eea495SMark             if (tt < eMinRad - 1.e-5) {
706*e0eea495SMark               eMinRad = tt;
707*e0eea495SMark               eMaxIdx = 0;
708*e0eea495SMark               eCellIdx[eMaxIdx++] = c;
709*e0eea495SMark             } else if (eMaxIdx > 0 && (tt-eMinRad) <= 1.e-5 && c != eCellIdx[eMaxIdx-1]) {
710*e0eea495SMark               eCellIdx[eMaxIdx++] = c;
711*e0eea495SMark             }
712*e0eea495SMark           }
713*e0eea495SMark           if ((tt=r-ctx->i_radius) > 0) {
714*e0eea495SMark             if (tt < iMinRad - 1.e-5) {
715*e0eea495SMark               iMinRad = tt;
716*e0eea495SMark               iMaxIdx = 0;
717*e0eea495SMark               iCellIdx[iMaxIdx++] = c;
718*e0eea495SMark             } else if (iMaxIdx > 0 && (tt-iMinRad) <= 1.e-5  && c != iCellIdx[iMaxIdx-1]) {
719*e0eea495SMark               iCellIdx[iMaxIdx++] = c;
720*e0eea495SMark             }
721*e0eea495SMark           }
722*e0eea495SMark         }
723*e0eea495SMark       }
724*e0eea495SMark     }
725*e0eea495SMark     for (k=0;k<nr;k++) {
726*e0eea495SMark       ierr = DMLabelSetValue(adaptLabel, rCellIdx[k], DM_ADAPT_REFINE);CHKERRQ(ierr);
727*e0eea495SMark     }
728*e0eea495SMark     if (ctx->sphere) {
729*e0eea495SMark       for (c = 0; c < eMaxIdx; c++) {
730*e0eea495SMark         ierr = DMLabelSetValue(adaptLabel, eCellIdx[c], DM_ADAPT_REFINE);CHKERRQ(ierr);
731*e0eea495SMark         ierr = PetscInfo3(sol, "\t\tPhase:%s: refine sphere e cell %D r=%g\n","adaptToleranceFEM",eCellIdx[c],eMinRad);
732*e0eea495SMark       }
733*e0eea495SMark       for (c = 0; c < iMaxIdx; c++) {
734*e0eea495SMark         ierr = DMLabelSetValue(adaptLabel, iCellIdx[c], DM_ADAPT_REFINE);CHKERRQ(ierr);
735*e0eea495SMark         ierr = PetscInfo3(sol, "\t\tPhase:%s: refine sphere i cell %D r=%g\n","adaptToleranceFEM",iCellIdx[c],iMinRad);
736*e0eea495SMark       }
737*e0eea495SMark     }
738*e0eea495SMark     ierr = PetscInfo4(sol, "Phase:%s: Adaptive refine origin cells %D,%D r=%g\n","adaptToleranceFEM",rCellIdx[0],rCellIdx[1],minRad);
739*e0eea495SMark   } else if (type==0 || type==1 || type==3) { /* refine along r=0 axis */
740*e0eea495SMark     PetscScalar  *coef = NULL;
741*e0eea495SMark     Vec          coords;
742*e0eea495SMark     PetscInt     csize,Nv,d,nz;
743*e0eea495SMark     DM           cdm;
744*e0eea495SMark     PetscSection cs;
745*e0eea495SMark     ierr = DMGetCoordinatesLocal(dm, &coords);CHKERRQ(ierr);
746*e0eea495SMark     ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
747*e0eea495SMark     ierr = DMGetLocalSection(cdm, &cs);CHKERRQ(ierr);
748*e0eea495SMark     for (c = cStart; c < cEnd; c++) {
749*e0eea495SMark       PetscInt doit = 0, outside = 0;
750*e0eea495SMark       ierr = DMPlexVecGetClosure(cdm, cs, coords, c, &csize, &coef);CHKERRQ(ierr);
751*e0eea495SMark       Nv = csize/dim;
752*e0eea495SMark       for (nz = d = 0; d < Nv; d++) {
753*e0eea495SMark         PetscReal z = PetscRealPart(coef[d*dim + (dim-1)]), x = PetscSqr(PetscRealPart(coef[d*dim + 0])) + ((dim==3) ? PetscSqr(PetscRealPart(coef[d*dim + 1])) : 0);
754*e0eea495SMark         x = PetscSqrtReal(x);
755*e0eea495SMark         if (x < 1e-12 && PetscAbsReal(z)<1e-12) doit = 1;             /* refine origin */
756*e0eea495SMark         else if (type==0 && (z < -1e-12 || z > ctx->re_radius+1e-12)) outside++;   /* first pass don't refine bottom */
757*e0eea495SMark         else if (type==1 && (z > ctx->vperp0_radius1 || z < -ctx->vperp0_radius1)) outside++; /* don't refine outside electron refine radius */
758*e0eea495SMark         else if (type==3 && (z > ctx->vperp0_radius2 || z < -ctx->vperp0_radius2)) outside++; /* don't refine outside ion refine radius */
759*e0eea495SMark         if (x < 1e-12) nz++;
760*e0eea495SMark       }
761*e0eea495SMark       ierr = DMPlexVecRestoreClosure(cdm, cs, coords, c, &csize, &coef);CHKERRQ(ierr);
762*e0eea495SMark       if (doit || (outside<Nv && nz)) {
763*e0eea495SMark         ierr = DMLabelSetValue(adaptLabel, c, DM_ADAPT_REFINE);CHKERRQ(ierr);
764*e0eea495SMark       }
765*e0eea495SMark     }
766*e0eea495SMark     ierr = PetscInfo1(sol, "Phase:%s: RE refinement\n","adaptToleranceFEM");
767*e0eea495SMark   }
768*e0eea495SMark   /* ierr = VecDestroy(&locX);CHKERRQ(ierr); */
769*e0eea495SMark   ierr = DMDestroy(&plex);CHKERRQ(ierr);
770*e0eea495SMark   ierr = DMAdaptLabel(dm, adaptLabel, &adaptedDM);CHKERRQ(ierr);
771*e0eea495SMark   ierr = DMLabelDestroy(&adaptLabel);CHKERRQ(ierr);
772*e0eea495SMark   *newDM = adaptedDM;
773*e0eea495SMark   if (adaptedDM) {
774*e0eea495SMark     if (isForest) {
775*e0eea495SMark       ierr = DMForestSetAdaptivityForest(adaptedDM,NULL);CHKERRQ(ierr);
776*e0eea495SMark     }
777*e0eea495SMark     ierr = DMConvert(adaptedDM, DMPLEX, &plex);CHKERRQ(ierr);
778*e0eea495SMark     ierr = DMPlexGetHeightStratum(plex,0,&cStart,&cEnd);CHKERRQ(ierr);
779*e0eea495SMark     ierr = PetscInfo2(sol, "\tPhase: adaptToleranceFEM: %D cells, %d total quadrature points\n",cEnd-cStart,Nq*(cEnd-cStart));CHKERRQ(ierr);
780*e0eea495SMark     ierr = DMDestroy(&plex);CHKERRQ(ierr);
781*e0eea495SMark   }
782*e0eea495SMark   PetscFunctionReturn(0);
783*e0eea495SMark }
784*e0eea495SMark 
785*e0eea495SMark static PetscErrorCode adapt(DM *dm, LandauCtx *ctx, Vec *uu)
786*e0eea495SMark {
787*e0eea495SMark   PetscErrorCode  ierr;
788*e0eea495SMark   PetscInt        type, limits[5] = {ctx->numRERefine,ctx->nZRefine1,ctx->maxRefIts,ctx->nZRefine2,ctx->postAMRRefine};
789*e0eea495SMark   PetscInt        adaptIter;
790*e0eea495SMark 
791*e0eea495SMark   PetscFunctionBegin;
792*e0eea495SMark   for (type=0;type<5;type++) {
793*e0eea495SMark     for (adaptIter = 0; adaptIter<limits[type];adaptIter++) {
794*e0eea495SMark       DM  dmNew = NULL;
795*e0eea495SMark       ierr = adaptToleranceFEM(ctx->fe[0], *uu, ctx->refineTol, ctx->coarsenTol, type, ctx, &dmNew);CHKERRQ(ierr);
796*e0eea495SMark       if (!dmNew) {
797*e0eea495SMark         exit(113);
798*e0eea495SMark         break;
799*e0eea495SMark       } else {
800*e0eea495SMark         ierr = DMDestroy(dm);CHKERRQ(ierr);
801*e0eea495SMark         ierr = VecDestroy(uu);CHKERRQ(ierr);
802*e0eea495SMark         ierr = DMCreateGlobalVector(dmNew,uu);CHKERRQ(ierr);
803*e0eea495SMark         ierr = PetscObjectSetName((PetscObject) *uu, "u");CHKERRQ(ierr);
804*e0eea495SMark         ierr = LandauSetInitialCondition(dmNew, *uu, ctx);CHKERRQ(ierr);
805*e0eea495SMark         *dm = dmNew;
806*e0eea495SMark       }
807*e0eea495SMark     }
808*e0eea495SMark   }
809*e0eea495SMark   PetscFunctionReturn(0);
810*e0eea495SMark }
811*e0eea495SMark 
812*e0eea495SMark static PetscErrorCode ProcessOptions(LandauCtx *ctx, const char prefix[])
813*e0eea495SMark {
814*e0eea495SMark   PetscErrorCode    ierr;
815*e0eea495SMark   PetscBool         flg, sph_flg;
816*e0eea495SMark   PetscInt          ii,nt,nm,nc;
817*e0eea495SMark   DM                dummy;
818*e0eea495SMark 
819*e0eea495SMark   PetscFunctionBegin;
820*e0eea495SMark   ierr = DMCreate(PETSC_COMM_WORLD,&dummy);CHKERRQ(ierr);
821*e0eea495SMark   /* get options - initialize context */
822*e0eea495SMark   ctx->normJ = 0;
823*e0eea495SMark   ctx->verbose = 1;
824*e0eea495SMark   ctx->interpolate = PETSC_TRUE;
825*e0eea495SMark   ctx->simplex = PETSC_FALSE;
826*e0eea495SMark   ctx->sphere = PETSC_FALSE;
827*e0eea495SMark   ctx->inflate = PETSC_FALSE;
828*e0eea495SMark   ctx->electronShift = 0;
829*e0eea495SMark   ctx->errorIndicator = NULL;
830*e0eea495SMark   ctx->radius = 5.; /* electron thermal radius (velocity) */
831*e0eea495SMark   ctx->re_radius = 0.;
832*e0eea495SMark   ctx->vperp0_radius1 = 0;
833*e0eea495SMark   ctx->vperp0_radius2 = 0;
834*e0eea495SMark   ctx->e_radius = .1;
835*e0eea495SMark   ctx->i_radius = .01;
836*e0eea495SMark   ctx->maxRefIts = 5;
837*e0eea495SMark   ctx->postAMRRefine = 0;
838*e0eea495SMark   ctx->nZRefine1 = 0;
839*e0eea495SMark   ctx->nZRefine2 = 0;
840*e0eea495SMark   ctx->numRERefine = 0;
841*e0eea495SMark   ctx->num_sections = 3; /* 2, 3 or 4 */
842*e0eea495SMark   /* species - [0] electrons, [1] one ion species eg, duetarium, [2] heavy impurity ion, ... */
843*e0eea495SMark   ctx->charges[0] = -1;  /* electron charge (MKS) */
844*e0eea495SMark   ctx->masses[0] = 1/1835.5; /* temporary value in proton mass */
845*e0eea495SMark   ctx->n[0] = 1;
846*e0eea495SMark   ctx->thermal_temps[0] = 1;
847*e0eea495SMark   /* constants, etc. */
848*e0eea495SMark   ctx->epsilon0 = 8.8542e-12; /* permittivity of free space (MKS) F/m */
849*e0eea495SMark   ctx->k = 1.38064852e-23; /* Boltzmann constant (MKS) J/K */
850*e0eea495SMark   ctx->lnLam = 10;         /* cross section ratio large - small angle collisions */
851*e0eea495SMark   ctx->n_0 = 1.e20;        /* typical plasma n, but could set it to 1 */
852*e0eea495SMark   ctx->Ez = 0;
853*e0eea495SMark   ctx->v_0 = 1; /* in electron thermal velocity */
854*e0eea495SMark   ctx->subThreadBlockSize = 1; /* for device and maybe OMP */
855*e0eea495SMark   ctx->quarter3DDomain = PETSC_FALSE;
856*e0eea495SMark   ierr = PetscOptionsBegin(PETSC_COMM_WORLD, prefix, "Options for Fokker-Plank-Landau collision operator", "none");CHKERRQ(ierr);
857*e0eea495SMark   {
858*e0eea495SMark     char opstring[256];
859*e0eea495SMark #if defined(PETSC_HAVE_KOKKOS)
860*e0eea495SMark     ctx->deviceType = LANDAU_KOKKOS;
861*e0eea495SMark     ierr = PetscStrcpy(opstring,"kokkos");CHKERRQ(ierr);
862*e0eea495SMark #if !defined(PETSC_HAVE_CUDA)
863*e0eea495SMark     ctx->subThreadBlockSize = 0;
864*e0eea495SMark #endif
865*e0eea495SMark #elif defined(PETSC_HAVE_CUDA)
866*e0eea495SMark     ctx->deviceType = LANDAU_CUDA;
867*e0eea495SMark     ierr = PetscStrcpy(opstring,"cuda");CHKERRQ(ierr);
868*e0eea495SMark #else
869*e0eea495SMark     ctx->deviceType = LANDAU_CPU;
870*e0eea495SMark     ierr = PetscStrcpy(opstring,"cpu");CHKERRQ(ierr);
871*e0eea495SMark     ctx->subThreadBlockSize = 0;
872*e0eea495SMark #endif
873*e0eea495SMark     ierr = PetscOptionsString("-dm_landau_device_type","Use kernels on 'cpu', 'cuda', or 'kokkos'","plexland.c",opstring,opstring,256,NULL);CHKERRQ(ierr);
874*e0eea495SMark     ierr = PetscStrcmp("cpu",opstring,&flg);CHKERRQ(ierr);
875*e0eea495SMark     if (flg) {
876*e0eea495SMark       ctx->deviceType = LANDAU_CPU;
877*e0eea495SMark       ctx->subThreadBlockSize = 0;
878*e0eea495SMark     } else {
879*e0eea495SMark       ierr = PetscStrcmp("cuda",opstring,&flg);CHKERRQ(ierr);
880*e0eea495SMark       if (flg) ctx->deviceType = LANDAU_CUDA;
881*e0eea495SMark       else {
882*e0eea495SMark         ierr = PetscStrcmp("kokkos",opstring,&flg);CHKERRQ(ierr);
883*e0eea495SMark         if (flg) ctx->deviceType = LANDAU_KOKKOS;
884*e0eea495SMark         else SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG,"-dm_landau_device_type %s",opstring);
885*e0eea495SMark       }
886*e0eea495SMark     }
887*e0eea495SMark   }
888*e0eea495SMark   ierr = PetscOptionsReal("-dm_landau_electron_shift","Shift in thermal velocity of electrons","none",ctx->electronShift,&ctx->electronShift, NULL);CHKERRQ(ierr);
889*e0eea495SMark   ierr = PetscOptionsBool("-dm_landau_sphere", "use sphere/semi-circle domain instead of rectangle", "plexland.c", ctx->sphere, &ctx->sphere, &sph_flg);CHKERRQ(ierr);
890*e0eea495SMark   ierr = PetscOptionsBool("-dm_landau_inflate", "With sphere, inflate for curved edges (no AMR)", "plexland.c", ctx->inflate, &ctx->inflate, NULL);CHKERRQ(ierr);
891*e0eea495SMark   /* ierr = PetscOptionsBool("-dm_landau_quarter_3d_domain", "Use symmetry in 3D to model 1/4 of domain", "plexland.c", ctx->quarter3DDomain, &ctx->quarter3DDomain, NULL);CHKERRQ(ierr); */
892*e0eea495SMark   ierr = PetscOptionsInt("-dm_landau_amr_re_levels", "Number of levels to refine along v_perp=0, z>0", "plexland.c", ctx->numRERefine, &ctx->numRERefine, NULL);CHKERRQ(ierr);
893*e0eea495SMark   ierr = PetscOptionsInt("-dm_landau_amr_z_refine1",  "Number of levels to refine along v_perp=0", "plexland.c", ctx->nZRefine1, &ctx->nZRefine1, NULL);CHKERRQ(ierr);
894*e0eea495SMark   ierr = PetscOptionsInt("-dm_landau_amr_z_refine2",  "Number of levels to refine along v_perp=0", "plexland.c", ctx->nZRefine2, &ctx->nZRefine2, NULL);CHKERRQ(ierr);
895*e0eea495SMark   ierr = PetscOptionsInt("-dm_landau_amr_levels_max", "Number of AMR levels of refinement around origin after r=0 refinements", "plexland.c", ctx->maxRefIts, &ctx->maxRefIts, NULL);CHKERRQ(ierr);
896*e0eea495SMark   ierr = PetscOptionsInt("-dm_landau_amr_post_refine", "Number of levels to uniformly refine after AMR", "plexland.c", ctx->postAMRRefine, &ctx->postAMRRefine, NULL);CHKERRQ(ierr);
897*e0eea495SMark   ierr = PetscOptionsInt("-dm_landau_verbose", "", "plexland.c", ctx->verbose, &ctx->verbose, NULL);CHKERRQ(ierr);
898*e0eea495SMark   ierr = PetscOptionsReal("-dm_landau_re_radius","velocity range to refine on positive (z>0) r=0 axis for runaways","plexland.c",ctx->re_radius,&ctx->re_radius, &flg);CHKERRQ(ierr);
899*e0eea495SMark   ierr = PetscOptionsReal("-dm_landau_z_radius1","velocity range to refine r=0 axis (for electrons)","plexland.c",ctx->vperp0_radius1,&ctx->vperp0_radius1, &flg);CHKERRQ(ierr);
900*e0eea495SMark   ierr = PetscOptionsReal("-dm_landau_z_radius2","velocity range to refine r=0 axis (for ions) after origin AMR","plexland.c",ctx->vperp0_radius2,&ctx->vperp0_radius2, &flg);CHKERRQ(ierr);
901*e0eea495SMark   ierr = PetscOptionsReal("-dm_landau_Ez","Initial parallel electric field in unites of Conner-Hastie criticle field","plexland.c",ctx->Ez,&ctx->Ez, NULL);CHKERRQ(ierr);
902*e0eea495SMark   ierr = PetscOptionsReal("-dm_landau_n_0","Normalization constant for number density","plexland.c",ctx->n_0,&ctx->n_0, NULL);CHKERRQ(ierr);
903*e0eea495SMark   ierr = PetscOptionsReal("-dm_landau_ln_lambda","Cross section parameter","plexland.c",ctx->lnLam,&ctx->lnLam, NULL);CHKERRQ(ierr);
904*e0eea495SMark   ierr = PetscOptionsInt("-dm_landau_num_sections", "Number of tangential section in (2D) grid, 2, 3, of 4", "plexland.c", ctx->num_sections, &ctx->num_sections, NULL);CHKERRQ(ierr);
905*e0eea495SMark   ctx->simplex = PETSC_FALSE;
906*e0eea495SMark   /* get num species with tempurature*/
907*e0eea495SMark   {
908*e0eea495SMark     PetscReal arr[100];
909*e0eea495SMark     nt = 100;
910*e0eea495SMark     ierr = PetscOptionsRealArray("-dm_landau_thermal_temps", "Temperature of each species [e,i_0,i_1,...] in keV", "plexland.c", arr, &nt, &flg);CHKERRQ(ierr);
911*e0eea495SMark     if (flg && nt > LANDAU_MAX_SPECIES) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG,"-thermal_temps ,t1,t2,.. number of species %D > MAX %D",nt,LANDAU_MAX_SPECIES);
912*e0eea495SMark   }
913*e0eea495SMark   nt = LANDAU_MAX_SPECIES;
914*e0eea495SMark   for (ii=1;ii<LANDAU_MAX_SPECIES;ii++) {
915*e0eea495SMark     ctx->thermal_temps[ii] = 1.;
916*e0eea495SMark     ctx->charges[ii] = 1;
917*e0eea495SMark     ctx->masses[ii] = 1;
918*e0eea495SMark     ctx->n[ii] = (ii==1) ? 1 : 0;
919*e0eea495SMark   }
920*e0eea495SMark   ierr = PetscOptionsRealArray("-dm_landau_thermal_temps", "Temperature of each species [e,i_0,i_1,...] in keV (must be set to set number of species)", "plexland.c", ctx->thermal_temps, &nt, &flg);CHKERRQ(ierr);
921*e0eea495SMark   if (flg) {
922*e0eea495SMark     PetscInfo1(dummy, "num_species set to number of thermal temps provided (%D)\n",nt);
923*e0eea495SMark     ctx->num_species = nt;
924*e0eea495SMark   } else SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG,"-dm_landau_thermal_temps ,t1,t2,.. must be provided to set the number of species");
925*e0eea495SMark   for (ii=0;ii<ctx->num_species;ii++) ctx->thermal_temps[ii] *= 1.1604525e7; /* convert to Kelvin */
926*e0eea495SMark   nm = LANDAU_MAX_SPECIES-1;
927*e0eea495SMark   ierr = PetscOptionsRealArray("-dm_landau_ion_masses", "Mass of each species in units of proton mass [i_0=2,i_1=40...]", "plexland.c", &ctx->masses[1], &nm, &flg);CHKERRQ(ierr);
928*e0eea495SMark   if (flg && nm != ctx->num_species-1) {
929*e0eea495SMark     SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG,"num ion masses %D != num species %D",nm,ctx->num_species-1);
930*e0eea495SMark   }
931*e0eea495SMark   nm = LANDAU_MAX_SPECIES;
932*e0eea495SMark   ierr = PetscOptionsRealArray("-dm_landau_n", "Normalized (by -n_0) number density of each species", "plexland.c", ctx->n, &nm, &flg);CHKERRQ(ierr);
933*e0eea495SMark   if (flg && nm != ctx->num_species) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG,"wrong num n: %D != num species %D",nm,ctx->num_species);
934*e0eea495SMark   ctx->n_0 *= ctx->n[0]; /* normalized number density */
935*e0eea495SMark   for (ii=1;ii<ctx->num_species;ii++) ctx->n[ii] = ctx->n[ii]/ctx->n[0];
936*e0eea495SMark   ctx->n[0] = 1;
937*e0eea495SMark   for (ii=0;ii<LANDAU_MAX_SPECIES;ii++) ctx->masses[ii] *= 1.6720e-27; /* scale by proton mass kg */
938*e0eea495SMark   ctx->masses[0] = 9.10938356e-31; /* electron mass kg (should be about right already) */
939*e0eea495SMark   ctx->m_0 = ctx->masses[0]; /* arbitrary reference mass, electrons */
940*e0eea495SMark   ierr = PetscOptionsReal("-dm_landau_v_0","Velocity to normalize with in units of initial electrons thermal velocity (not recommended to change default)","plexland.c",ctx->v_0,&ctx->v_0, NULL);CHKERRQ(ierr);
941*e0eea495SMark   ctx->v_0 *= PetscSqrtReal(ctx->k*ctx->thermal_temps[0]/(ctx->masses[0])); /* electron mean velocity in 1D (need 3D form in computing T from FE integral) */
942*e0eea495SMark   nc = LANDAU_MAX_SPECIES-1;
943*e0eea495SMark   ierr = PetscOptionsRealArray("-dm_landau_ion_charges", "Charge of each species in units of proton charge [i_0=2,i_1=18,...]", "plexland.c", &ctx->charges[1], &nc, &flg);CHKERRQ(ierr);
944*e0eea495SMark   if (flg && nc != ctx->num_species-1) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG,"num charges %D != num species %D",nc,ctx->num_species-1);
945*e0eea495SMark   for (ii=0;ii<LANDAU_MAX_SPECIES;ii++) ctx->charges[ii] *= 1.6022e-19; /* electron/proton charge (MKS) */
946*e0eea495SMark   ctx->t_0 = 8*PETSC_PI*PetscSqr(ctx->epsilon0*ctx->m_0/PetscSqr(ctx->charges[0]))/ctx->lnLam/ctx->n_0*PetscPowReal(ctx->v_0,3); /* note, this t_0 makes nu[0,0]=1 */
947*e0eea495SMark   /* geometry */
948*e0eea495SMark   for (ii=0;ii<ctx->num_species;ii++) ctx->refineTol[ii]  = PETSC_MAX_REAL;
949*e0eea495SMark   for (ii=0;ii<ctx->num_species;ii++) ctx->coarsenTol[ii] = 0.;
950*e0eea495SMark   ii = LANDAU_MAX_SPECIES;
951*e0eea495SMark   ierr = PetscOptionsRealArray("-dm_landau_refine_tol","tolerance for refining cells in AMR","plexland.c",ctx->refineTol, &ii, &flg);CHKERRQ(ierr);
952*e0eea495SMark   if (flg && ii != ctx->num_species) ierr = PetscInfo2(dummy, "Phase: Warning, #refine_tol %D != num_species %D\n",ii,ctx->num_species);CHKERRQ(ierr);
953*e0eea495SMark   ii = LANDAU_MAX_SPECIES;
954*e0eea495SMark   ierr = PetscOptionsRealArray("-dm_landau_coarsen_tol","tolerance for coarsening cells in AMR","plexland.c",ctx->coarsenTol, &ii, &flg);CHKERRQ(ierr);
955*e0eea495SMark   if (flg && ii != ctx->num_species) ierr = PetscInfo2(dummy, "Phase: Warning, #coarsen_tol %D != num_species %D\n",ii,ctx->num_species);CHKERRQ(ierr);
956*e0eea495SMark   ierr = PetscOptionsReal("-dm_landau_domain_radius","Phase space size in units of electron thermal velocity","plexland.c",ctx->radius,&ctx->radius, &flg);CHKERRQ(ierr);
957*e0eea495SMark   if (flg && ctx->radius <= 0) { /* negative is ratio of c */
958*e0eea495SMark     if (ctx->radius == 0) ctx->radius = 0.75;
959*e0eea495SMark     else ctx->radius = -ctx->radius;
960*e0eea495SMark     ctx->radius = ctx->radius*299792458/ctx->v_0;
961*e0eea495SMark     ierr = PetscInfo1(dummy, "Change domain radius to %e\n",ctx->radius);CHKERRQ(ierr);
962*e0eea495SMark   }
963*e0eea495SMark   ierr = PetscOptionsReal("-dm_landau_i_radius","Ion thermal velocity, used for circular meshes","plexland.c",ctx->i_radius,&ctx->i_radius, &flg);CHKERRQ(ierr);
964*e0eea495SMark   if (flg && !sph_flg) ctx->sphere = PETSC_TRUE; /* you gave me an ion radius but did not set sphere, user error really */
965*e0eea495SMark   if (!flg) {
966*e0eea495SMark     ctx->i_radius = 1.5*PetscSqrtReal(8*ctx->k*ctx->thermal_temps[1]/ctx->masses[1]/PETSC_PI)/ctx->v_0; /* normalized radius with thermal velocity of first ion */
967*e0eea495SMark   }
968*e0eea495SMark   ierr = PetscOptionsReal("-dm_landau_e_radius","Electron thermal velocity, used for circular meshes","plexland.c",ctx->e_radius,&ctx->e_radius, &flg);CHKERRQ(ierr);
969*e0eea495SMark   if (flg && !sph_flg) ctx->sphere = PETSC_TRUE; /* you gave me an e radius but did not set sphere, user error really */
970*e0eea495SMark   if (!flg) {
971*e0eea495SMark     ctx->e_radius = 1.5*PetscSqrtReal(8*ctx->k*ctx->thermal_temps[0]/ctx->masses[0]/PETSC_PI)/ctx->v_0; /* normalized radius with thermal velocity of electrons */
972*e0eea495SMark   }
973*e0eea495SMark   if (ctx->sphere && (ctx->e_radius <= ctx->i_radius || ctx->radius <= ctx->e_radius)) SETERRQ3(PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG,"bad radii: %g < %g < %g",ctx->i_radius,ctx->e_radius,ctx->radius);
974*e0eea495SMark   ierr = PetscOptionsInt("-dm_landau_sub_thread_block_size", "Number of threads in CUDA integration point subblock", "plexland.c", ctx->subThreadBlockSize, &ctx->subThreadBlockSize, NULL);CHKERRQ(ierr);
975*e0eea495SMark   if (ctx->subThreadBlockSize > LANDAU_MAX_SUB_THREAD_BLOCKS) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG,"num sub threads %D > MAX %D",ctx->subThreadBlockSize,LANDAU_MAX_SUB_THREAD_BLOCKS);
976*e0eea495SMark   ierr = PetscOptionsEnd();CHKERRQ(ierr);
977*e0eea495SMark   for (ii=ctx->num_species;ii<LANDAU_MAX_SPECIES;ii++) ctx->masses[ii] = ctx->thermal_temps[ii]  = ctx->charges[ii] = 0;
978*e0eea495SMark   if (ctx->verbose > 0) {
979*e0eea495SMark     ierr = PetscPrintf(PETSC_COMM_WORLD, "masses:        e=%10.3e; ions in proton mass units:   %10.3e %10.3e ...\n",ctx->masses[0],ctx->masses[1]/1.6720e-27,ctx->num_species>2 ? ctx->masses[2]/1.6720e-27 : 0);CHKERRQ(ierr);
980*e0eea495SMark     ierr = PetscPrintf(PETSC_COMM_WORLD, "charges:       e=%10.3e; charges in elementary units: %10.3e %10.3e\n", ctx->charges[0],-ctx->charges[1]/ctx->charges[0],ctx->num_species>2 ? -ctx->charges[2]/ctx->charges[0] : 0);CHKERRQ(ierr);
981*e0eea495SMark     ierr = PetscPrintf(PETSC_COMM_WORLD, "thermal T (K): e=%10.3e i=%10.3e imp=%10.3e. v_0=%10.3e n_0=%10.3e t_0=%10.3e domain=%10.3e\n",ctx->thermal_temps[0],ctx->thermal_temps[1],ctx->num_species>2 ? ctx->thermal_temps[2] : 0,ctx->v_0,ctx->n_0,ctx->t_0,ctx->radius);CHKERRQ(ierr);
982*e0eea495SMark   }
983*e0eea495SMark   ierr = DMDestroy(&dummy);CHKERRQ(ierr);
984*e0eea495SMark   {
985*e0eea495SMark     PetscMPIInt    rank;
986*e0eea495SMark     ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank);CHKERRQ(ierr);
987*e0eea495SMark     /* PetscLogStage  setup_stage; */
988*e0eea495SMark     ierr = PetscLogEventRegister("Landau Operator", DM_CLASSID, &ctx->events[0]);CHKERRQ(ierr); /* 0 */
989*e0eea495SMark     ierr = PetscLogEventRegister(" Jac-vector", DM_CLASSID, &ctx->events[1]);CHKERRQ(ierr); /* 1 */
990*e0eea495SMark     ierr = PetscLogEventRegister(" Jac-kern-init", DM_CLASSID, &ctx->events[3]);CHKERRQ(ierr); /* 3 */
991*e0eea495SMark     ierr = PetscLogEventRegister(" Jac-kernel", DM_CLASSID, &ctx->events[4]);CHKERRQ(ierr); /* 4 */
992*e0eea495SMark     ierr = PetscLogEventRegister(" Jac-kernel-post", DM_CLASSID, &ctx->events[5]);CHKERRQ(ierr); /* 5 */
993*e0eea495SMark     ierr = PetscLogEventRegister(" Jac-assemble", DM_CLASSID, &ctx->events[6]);CHKERRQ(ierr); /* 6 */
994*e0eea495SMark     ierr = PetscLogEventRegister(" Jac-end", DM_CLASSID, &ctx->events[7]);CHKERRQ(ierr); /* 7 */
995*e0eea495SMark     ierr = PetscLogEventRegister("  Jac-geo-color", DM_CLASSID, &ctx->events[8]);CHKERRQ(ierr); /* 8 */
996*e0eea495SMark     ierr = PetscLogEventRegister("  Jac-cuda-sum", DM_CLASSID, &ctx->events[2]);CHKERRQ(ierr); /* 2 */
997*e0eea495SMark     ierr = PetscLogEventRegister("Landau Jacobian", DM_CLASSID, &ctx->events[9]);CHKERRQ(ierr); /* 9 */
998*e0eea495SMark     if (rank) { /* turn off output stuff for duplicate runs - do we need to add the prefix to all this? */
999*e0eea495SMark       ierr = PetscOptionsClearValue(NULL,"-snes_converged_reason");CHKERRQ(ierr);
1000*e0eea495SMark       ierr = PetscOptionsClearValue(NULL,"-ksp_converged_reason");CHKERRQ(ierr);
1001*e0eea495SMark       ierr = PetscOptionsClearValue(NULL,"-snes_monitor");CHKERRQ(ierr);
1002*e0eea495SMark       ierr = PetscOptionsClearValue(NULL,"-ksp_monitor");CHKERRQ(ierr);
1003*e0eea495SMark       ierr = PetscOptionsClearValue(NULL,"-ts_monitor");CHKERRQ(ierr);
1004*e0eea495SMark       ierr = PetscOptionsClearValue(NULL,"-ts_adapt_monitor");CHKERRQ(ierr);
1005*e0eea495SMark       ierr = PetscOptionsClearValue(NULL,"-dm_landau_amr_dm_view");CHKERRQ(ierr);
1006*e0eea495SMark       ierr = PetscOptionsClearValue(NULL,"-dm_landau_amr_vec_view");CHKERRQ(ierr);
1007*e0eea495SMark       ierr = PetscOptionsClearValue(NULL,"-dm_landau_pre_dm_view");CHKERRQ(ierr);
1008*e0eea495SMark       ierr = PetscOptionsClearValue(NULL,"-dm_landau_pre_vec_view");CHKERRQ(ierr);
1009*e0eea495SMark       ierr = PetscOptionsClearValue(NULL,"-info");CHKERRQ(ierr);
1010*e0eea495SMark     }
1011*e0eea495SMark   }
1012*e0eea495SMark   PetscFunctionReturn(0);
1013*e0eea495SMark }
1014*e0eea495SMark 
1015*e0eea495SMark /*@C
1016*e0eea495SMark   LandauCreateVelocitySpace - Create a DMPlex velocity space mesh
1017*e0eea495SMark 
1018*e0eea495SMark   Collective on comm
1019*e0eea495SMark 
1020*e0eea495SMark   Input Parameters:
1021*e0eea495SMark +   comm  - The MPI communicator
1022*e0eea495SMark .   dim - velocity space dimension (2 for axisymmetric, 3 for full 3X + 3V solver)
1023*e0eea495SMark -   prefix - prefix for options
1024*e0eea495SMark 
1025*e0eea495SMark   Output Parameter:
1026*e0eea495SMark .   dm  - The DM object representing the mesh
1027*e0eea495SMark +   X - A vector (user destroys)
1028*e0eea495SMark -   J - Optional matrix (object destroys)
1029*e0eea495SMark 
1030*e0eea495SMark   Level: beginner
1031*e0eea495SMark 
1032*e0eea495SMark .keywords: mesh
1033*e0eea495SMark .seealso: DMPlexCreate(), LandauDestroyVelocitySpace()
1034*e0eea495SMark @*/
1035*e0eea495SMark PetscErrorCode LandauCreateVelocitySpace(MPI_Comm comm, PetscInt dim, const char prefix[], Vec *X, Mat *J, DM *dm)
1036*e0eea495SMark {
1037*e0eea495SMark   PetscMPIInt    size;
1038*e0eea495SMark   PetscErrorCode ierr;
1039*e0eea495SMark   LandauCtx      *ctx;
1040*e0eea495SMark 
1041*e0eea495SMark   PetscFunctionBegin;
1042*e0eea495SMark   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
1043*e0eea495SMark   if (size!=1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Velocity space meshes should be serial (but should work in parallel)");
1044*e0eea495SMark   if (dim!=2 && dim!=3) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Only 2D and 3D supported");
1045*e0eea495SMark   ctx = (LandauCtx*)malloc(sizeof(LandauCtx));
1046*e0eea495SMark   /* process options */
1047*e0eea495SMark   ierr = ProcessOptions(ctx,prefix);CHKERRQ(ierr);
1048*e0eea495SMark   /* Create Mesh */
1049*e0eea495SMark   ierr = LandauDMCreateVMesh(comm, dim, prefix, ctx, dm);CHKERRQ(ierr);
1050*e0eea495SMark   ierr = DMViewFromOptions(*dm,NULL,"-dm_landau_pre_dm_view");CHKERRQ(ierr);
1051*e0eea495SMark   ierr = DMSetApplicationContext(*dm, ctx);CHKERRQ(ierr);
1052*e0eea495SMark   /* create FEM */
1053*e0eea495SMark   ierr = SetupDS(*dm,dim,ctx);CHKERRQ(ierr);
1054*e0eea495SMark   /* set initial state */
1055*e0eea495SMark   ierr = DMCreateGlobalVector(*dm,X);CHKERRQ(ierr);
1056*e0eea495SMark   ierr = PetscObjectSetName((PetscObject) *X, "u");CHKERRQ(ierr);
1057*e0eea495SMark   /* initial static refinement, no solve */
1058*e0eea495SMark   ierr = LandauSetInitialCondition(*dm, *X, ctx);CHKERRQ(ierr);
1059*e0eea495SMark   ierr = VecViewFromOptions(*X, NULL, "-dm_landau_pre_vec_view");CHKERRQ(ierr);
1060*e0eea495SMark   /* forest refinement */
1061*e0eea495SMark   if (ctx->errorIndicator) {
1062*e0eea495SMark     /* AMR */
1063*e0eea495SMark     ierr = adapt(dm,ctx,X);CHKERRQ(ierr);
1064*e0eea495SMark     ierr = DMViewFromOptions(*dm,NULL,"-dm_landau_amr_dm_view");CHKERRQ(ierr);
1065*e0eea495SMark     ierr = VecViewFromOptions(*X, NULL, "-dm_landau_amr_vec_view");CHKERRQ(ierr);
1066*e0eea495SMark   }
1067*e0eea495SMark   ierr = DMSetApplicationContext(*dm, ctx);CHKERRQ(ierr);
1068*e0eea495SMark   ctx->dmv = *dm;
1069*e0eea495SMark   ierr = DMCreateMatrix(ctx->dmv, &ctx->J);CHKERRQ(ierr);
1070*e0eea495SMark   if (J) *J = ctx->J;
1071*e0eea495SMark   PetscFunctionReturn(0);
1072*e0eea495SMark }
1073*e0eea495SMark 
1074*e0eea495SMark /*@
1075*e0eea495SMark   LandauDestroyVelocitySpace - Destroy a DMPlex velocity space mesh
1076*e0eea495SMark 
1077*e0eea495SMark   Collective on dm
1078*e0eea495SMark 
1079*e0eea495SMark   Input/Output Parameters:
1080*e0eea495SMark   .   dm - the dm to destroy
1081*e0eea495SMark 
1082*e0eea495SMark   Level: beginner
1083*e0eea495SMark 
1084*e0eea495SMark .keywords: mesh
1085*e0eea495SMark .seealso: LandauCreateVelocitySpace()
1086*e0eea495SMark @*/
1087*e0eea495SMark PetscErrorCode LandauDestroyVelocitySpace(DM *dm)
1088*e0eea495SMark {
1089*e0eea495SMark   PetscErrorCode ierr,ii;
1090*e0eea495SMark   LandauCtx        *ctx;
1091*e0eea495SMark   PetscContainer container = NULL;
1092*e0eea495SMark   PetscFunctionBegin;
1093*e0eea495SMark   ierr = DMGetApplicationContext(*dm, &ctx);CHKERRQ(ierr);
1094*e0eea495SMark   ierr = PetscObjectQuery((PetscObject)ctx->J,"coloring", (PetscObject*)&container);CHKERRQ(ierr);
1095*e0eea495SMark   if (container) {
1096*e0eea495SMark     ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
1097*e0eea495SMark   }
1098*e0eea495SMark   ierr = MatDestroy(&ctx->M);CHKERRQ(ierr);
1099*e0eea495SMark   ierr = MatDestroy(&ctx->J);CHKERRQ(ierr);
1100*e0eea495SMark   for (ii=0;ii<ctx->num_species;ii++) {
1101*e0eea495SMark     ierr = PetscFEDestroy(&ctx->fe[ii]);CHKERRQ(ierr);
1102*e0eea495SMark   }
1103*e0eea495SMark   free(ctx);
1104*e0eea495SMark   ierr = DMDestroy(dm);CHKERRQ(ierr);
1105*e0eea495SMark   PetscFunctionReturn(0);
1106*e0eea495SMark }
1107*e0eea495SMark 
1108*e0eea495SMark /* < v, ru > */
1109*e0eea495SMark static void f0_s_den(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1110*e0eea495SMark                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1111*e0eea495SMark                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1112*e0eea495SMark                      PetscReal t, const PetscReal x[],  PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0)
1113*e0eea495SMark {
1114*e0eea495SMark   PetscInt ii = (PetscInt)PetscRealPart(constants[0]);
1115*e0eea495SMark   f0[0] = u[ii];
1116*e0eea495SMark }
1117*e0eea495SMark 
1118*e0eea495SMark /* < v, ru > */
1119*e0eea495SMark static void f0_s_mom(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1120*e0eea495SMark                     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1121*e0eea495SMark                     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1122*e0eea495SMark                     PetscReal t, const PetscReal x[],  PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0)
1123*e0eea495SMark {
1124*e0eea495SMark   PetscInt ii = (PetscInt)PetscRealPart(constants[0]), jj = (PetscInt)PetscRealPart(constants[1]);
1125*e0eea495SMark   f0[0] = x[jj]*u[ii]; /* x momentum */
1126*e0eea495SMark }
1127*e0eea495SMark 
1128*e0eea495SMark static void f0_s_v2(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1129*e0eea495SMark                     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1130*e0eea495SMark                     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1131*e0eea495SMark                     PetscReal t, const PetscReal x[],  PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0)
1132*e0eea495SMark {
1133*e0eea495SMark   PetscInt i, ii = (PetscInt)PetscRealPart(constants[0]);
1134*e0eea495SMark   double tmp1 = 0.;
1135*e0eea495SMark   for (i = 0; i < dim; ++i) tmp1 += x[i]*x[i];
1136*e0eea495SMark   f0[0] = tmp1*u[ii];
1137*e0eea495SMark }
1138*e0eea495SMark 
1139*e0eea495SMark /* < v, ru > */
1140*e0eea495SMark static void f0_s_rden(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1141*e0eea495SMark                       const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1142*e0eea495SMark                       const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1143*e0eea495SMark                       PetscReal t, const PetscReal x[],  PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0)
1144*e0eea495SMark {
1145*e0eea495SMark   PetscInt ii = (PetscInt)PetscRealPart(constants[0]);
1146*e0eea495SMark   f0[0] = 2.*PETSC_PI*x[0]*u[ii];
1147*e0eea495SMark }
1148*e0eea495SMark 
1149*e0eea495SMark /* < v, ru > */
1150*e0eea495SMark static void f0_s_rmom(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1151*e0eea495SMark                       const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1152*e0eea495SMark                       const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1153*e0eea495SMark                       PetscReal t, const PetscReal x[],  PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0)
1154*e0eea495SMark {
1155*e0eea495SMark   PetscInt ii = (PetscInt)PetscRealPart(constants[0]);
1156*e0eea495SMark   f0[0] = 2.*PETSC_PI*x[0]*x[1]*u[ii];
1157*e0eea495SMark }
1158*e0eea495SMark 
1159*e0eea495SMark static void f0_s_rv2(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1160*e0eea495SMark                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1161*e0eea495SMark                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1162*e0eea495SMark                      PetscReal t, const PetscReal x[],  PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0)
1163*e0eea495SMark {
1164*e0eea495SMark   PetscInt ii = (PetscInt)PetscRealPart(constants[0]);
1165*e0eea495SMark   f0[0] =  2.*PETSC_PI*x[0]*(x[0]*x[0] + x[1]*x[1])*u[ii];
1166*e0eea495SMark }
1167*e0eea495SMark 
1168*e0eea495SMark /*@
1169*e0eea495SMark   LandauPrintNorms - collects moments and prints them
1170*e0eea495SMark 
1171*e0eea495SMark   Collective on dm
1172*e0eea495SMark 
1173*e0eea495SMark   Input Parameters:
1174*e0eea495SMark +   X  - the state
1175*e0eea495SMark .   stepi - current step to print
1176*e0eea495SMark 
1177*e0eea495SMark   Level: beginner
1178*e0eea495SMark 
1179*e0eea495SMark .keywords: mesh
1180*e0eea495SMark .seealso: LandauCreateVelocitySpace()
1181*e0eea495SMark @*/
1182*e0eea495SMark PetscErrorCode LandauPrintNorms(Vec X, PetscInt stepi)
1183*e0eea495SMark {
1184*e0eea495SMark   PetscErrorCode ierr;
1185*e0eea495SMark   LandauCtx      *ctx;
1186*e0eea495SMark   PetscDS        prob;
1187*e0eea495SMark   DM             plex,dm;
1188*e0eea495SMark   PetscInt       cStart, cEnd, dim, ii;
1189*e0eea495SMark   PetscScalar    xmomentumtot=0, ymomentumtot=0, zmomentumtot=0, energytot=0, densitytot=0, tt[LANDAU_MAX_SPECIES];
1190*e0eea495SMark   PetscScalar    xmomentum[LANDAU_MAX_SPECIES],  ymomentum[LANDAU_MAX_SPECIES],  zmomentum[LANDAU_MAX_SPECIES], energy[LANDAU_MAX_SPECIES], density[LANDAU_MAX_SPECIES];
1191*e0eea495SMark 
1192*e0eea495SMark   PetscFunctionBegin;
1193*e0eea495SMark   ierr = VecGetDM(X, &dm);CHKERRQ(ierr);
1194*e0eea495SMark   if (!dm) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "no DM");
1195*e0eea495SMark   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1196*e0eea495SMark   ierr = DMGetApplicationContext(dm, &ctx);CHKERRQ(ierr);
1197*e0eea495SMark   if (!ctx) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "no context");
1198*e0eea495SMark   ierr = DMConvert(ctx->dmv, DMPLEX, &plex);CHKERRQ(ierr);
1199*e0eea495SMark   ierr = DMCreateDS(plex);CHKERRQ(ierr);
1200*e0eea495SMark   ierr = DMGetDS(plex, &prob);CHKERRQ(ierr);
1201*e0eea495SMark   /* print momentum and energy */
1202*e0eea495SMark   for (ii=0;ii<ctx->num_species;ii++) {
1203*e0eea495SMark     PetscScalar user[2] = { (PetscScalar)ii, (PetscScalar)ctx->charges[ii]};
1204*e0eea495SMark     ierr = PetscDSSetConstants(prob, 2, user);CHKERRQ(ierr);
1205*e0eea495SMark     if (dim==2) { /* 2/3X + 3V (cylindrical coordinates) */
1206*e0eea495SMark       ierr = PetscDSSetObjective(prob, 0, &f0_s_rden);CHKERRQ(ierr);
1207*e0eea495SMark       ierr = DMPlexComputeIntegralFEM(plex,X,tt,ctx);CHKERRQ(ierr);
1208*e0eea495SMark       density[ii] = tt[0]*ctx->n_0*ctx->charges[ii];
1209*e0eea495SMark       ierr = PetscDSSetObjective(prob, 0, &f0_s_rmom);CHKERRQ(ierr);
1210*e0eea495SMark       ierr = DMPlexComputeIntegralFEM(plex,X,tt,ctx);CHKERRQ(ierr);
1211*e0eea495SMark       zmomentum[ii] = tt[0]*ctx->n_0*ctx->v_0*ctx->masses[ii];
1212*e0eea495SMark       ierr = PetscDSSetObjective(prob, 0, &f0_s_rv2);CHKERRQ(ierr);
1213*e0eea495SMark       ierr = DMPlexComputeIntegralFEM(plex,X,tt,ctx);CHKERRQ(ierr);
1214*e0eea495SMark       energy[ii] = tt[0]*0.5*ctx->n_0*ctx->v_0*ctx->v_0*ctx->masses[ii];
1215*e0eea495SMark       zmomentumtot += zmomentum[ii];
1216*e0eea495SMark       energytot  += energy[ii];
1217*e0eea495SMark       densitytot += density[ii];
1218*e0eea495SMark       ierr = PetscPrintf(PETSC_COMM_WORLD, "%3D) species-%D: charge density= %20.13e z-momentum= %20.13e energy= %20.13e",stepi,ii,density[ii],zmomentum[ii],energy[ii]);CHKERRQ(ierr);
1219*e0eea495SMark     } else { /* 2/3X + 3V */
1220*e0eea495SMark       ierr = PetscDSSetObjective(prob, 0, &f0_s_den);CHKERRQ(ierr);
1221*e0eea495SMark       ierr = DMPlexComputeIntegralFEM(plex,X,tt,ctx);CHKERRQ(ierr);
1222*e0eea495SMark       density[ii] = tt[0]*ctx->n_0*ctx->charges[ii];
1223*e0eea495SMark       ierr = PetscDSSetObjective(prob, 0, &f0_s_mom);CHKERRQ(ierr);
1224*e0eea495SMark       user[1] = 0;
1225*e0eea495SMark       ierr = DMPlexComputeIntegralFEM(plex,X,tt,ctx);CHKERRQ(ierr);
1226*e0eea495SMark       xmomentum[ii]  = tt[0]*ctx->n_0*ctx->v_0*ctx->masses[ii];
1227*e0eea495SMark       user[1] = 1;
1228*e0eea495SMark       ierr = DMPlexComputeIntegralFEM(plex,X,tt,ctx);CHKERRQ(ierr);
1229*e0eea495SMark       ymomentum[ii] = tt[0]*ctx->n_0*ctx->v_0*ctx->masses[ii];
1230*e0eea495SMark       user[1] = 2;
1231*e0eea495SMark       ierr = DMPlexComputeIntegralFEM(plex,X,tt,ctx);CHKERRQ(ierr);
1232*e0eea495SMark       zmomentum[ii] = tt[0]*ctx->n_0*ctx->v_0*ctx->masses[ii];
1233*e0eea495SMark       ierr = PetscDSSetObjective(prob, 0, &f0_s_v2);CHKERRQ(ierr);
1234*e0eea495SMark       ierr = DMPlexComputeIntegralFEM(plex,X,tt,ctx);CHKERRQ(ierr);
1235*e0eea495SMark       energy[ii]    = 0.5*tt[0]*ctx->n_0*ctx->v_0*ctx->v_0*ctx->masses[ii];
1236*e0eea495SMark       ierr = PetscPrintf(PETSC_COMM_WORLD, "%3D) species %D: density=%20.13e, x-momentum=%20.13e, y-momentum=%20.13e, z-momentum=%20.13e, energy=%21.13e",
1237*e0eea495SMark                          stepi,ii,density[ii],xmomentum[ii],ymomentum[ii],zmomentum[ii],energy[ii]);CHKERRQ(ierr);
1238*e0eea495SMark       xmomentumtot += xmomentum[ii];
1239*e0eea495SMark       ymomentumtot += ymomentum[ii];
1240*e0eea495SMark       zmomentumtot += zmomentum[ii];
1241*e0eea495SMark       energytot  += energy[ii];
1242*e0eea495SMark       densitytot += density[ii];
1243*e0eea495SMark     }
1244*e0eea495SMark     if (ctx->num_species>1) PetscPrintf(PETSC_COMM_WORLD, "\n");
1245*e0eea495SMark   }
1246*e0eea495SMark   /* totals */
1247*e0eea495SMark   ierr = DMPlexGetHeightStratum(plex,0,&cStart,&cEnd);CHKERRQ(ierr);
1248*e0eea495SMark   ierr = DMDestroy(&plex);CHKERRQ(ierr);
1249*e0eea495SMark   if (ctx->num_species>1) {
1250*e0eea495SMark     if (dim==2) {
1251*e0eea495SMark       ierr = PetscPrintf(PETSC_COMM_WORLD, "\t%3D) Total: charge density=%21.13e, momentum=%21.13e, energy=%21.13e (m_i[0]/m_e = %g, %D cells)",
1252*e0eea495SMark                   stepi,densitytot,zmomentumtot,energytot,ctx->masses[1]/ctx->masses[0],cEnd-cStart);CHKERRQ(ierr);
1253*e0eea495SMark     } else {
1254*e0eea495SMark       ierr = PetscPrintf(PETSC_COMM_WORLD, "\t%3D) Total: charge density=%21.13e, x-momentum=%21.13e, y-momentum=%21.13e, z-momentum=%21.13e, energy=%21.13e (m_i[0]/m_e = %g, %D cells)",
1255*e0eea495SMark                   stepi,densitytot,xmomentumtot,ymomentumtot,zmomentumtot,energytot,ctx->masses[1]/ctx->masses[0],cEnd-cStart);CHKERRQ(ierr);
1256*e0eea495SMark     }
1257*e0eea495SMark   } else {
1258*e0eea495SMark     ierr = PetscPrintf(PETSC_COMM_WORLD, " -- %D cells",cEnd-cStart);CHKERRQ(ierr);
1259*e0eea495SMark   }
1260*e0eea495SMark   if (ctx->verbose > 1) {ierr = PetscPrintf(PETSC_COMM_WORLD,", %D sub (vector) threads\n",ctx->subThreadBlockSize);CHKERRQ(ierr);}
1261*e0eea495SMark   else {ierr = PetscPrintf(PETSC_COMM_WORLD,"\n");CHKERRQ(ierr);}
1262*e0eea495SMark   PetscFunctionReturn(0);
1263*e0eea495SMark }
1264*e0eea495SMark 
1265*e0eea495SMark static PetscErrorCode destroy_coloring (void *is)
1266*e0eea495SMark {
1267*e0eea495SMark   ISColoring tmp = (ISColoring)is;
1268*e0eea495SMark   return ISColoringDestroy(&tmp);
1269*e0eea495SMark }
1270*e0eea495SMark 
1271*e0eea495SMark /*@
1272*e0eea495SMark   LandauCreateColoring - create a coloring and add to matrix (Landau context used just for 'print' flag, should be in DMPlex)
1273*e0eea495SMark 
1274*e0eea495SMark   Collective on JacP
1275*e0eea495SMark 
1276*e0eea495SMark   Input Parameters:
1277*e0eea495SMark +   JacP  - matrix to add coloring to
1278*e0eea495SMark .   plex - The DM
1279*e0eea495SMark 
1280*e0eea495SMark   Output Parameter:
1281*e0eea495SMark .   container  - Container with coloring
1282*e0eea495SMark 
1283*e0eea495SMark   Level: beginner
1284*e0eea495SMark 
1285*e0eea495SMark .keywords: mesh
1286*e0eea495SMark .seealso: LandauCreateVelocitySpace()
1287*e0eea495SMark @*/
1288*e0eea495SMark PetscErrorCode LandauCreateColoring(Mat JacP, DM plex, PetscContainer *container)
1289*e0eea495SMark {
1290*e0eea495SMark   PetscErrorCode  ierr;
1291*e0eea495SMark   PetscInt        dim,cell,i,ej,nc,Nv,totDim,numGCells,cStart,cEnd;
1292*e0eea495SMark   ISColoring      iscoloring = NULL;
1293*e0eea495SMark   Mat             G,Q;
1294*e0eea495SMark   PetscScalar     ones[128];
1295*e0eea495SMark   MatColoring     mc;
1296*e0eea495SMark   IS             *is;
1297*e0eea495SMark   PetscInt        csize,colour,j,k;
1298*e0eea495SMark   const PetscInt *indices;
1299*e0eea495SMark   PetscInt       numComp[1];
1300*e0eea495SMark   PetscInt       numDof[4];
1301*e0eea495SMark   PetscFE        fe;
1302*e0eea495SMark   DM             colordm;
1303*e0eea495SMark   PetscSection   csection, section, globalSection;
1304*e0eea495SMark   PetscDS        prob;
1305*e0eea495SMark   LandauCtx      *ctx;
1306*e0eea495SMark 
1307*e0eea495SMark   PetscFunctionBegin;
1308*e0eea495SMark   ierr = DMGetApplicationContext(plex, &ctx);CHKERRQ(ierr);
1309*e0eea495SMark   ierr = DMGetLocalSection(plex, &section);CHKERRQ(ierr);
1310*e0eea495SMark   ierr = DMGetGlobalSection(plex, &globalSection);CHKERRQ(ierr);
1311*e0eea495SMark   ierr = DMGetDimension(plex, &dim);CHKERRQ(ierr);
1312*e0eea495SMark   ierr = DMGetDS(plex, &prob);CHKERRQ(ierr);
1313*e0eea495SMark   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
1314*e0eea495SMark   ierr = DMPlexGetHeightStratum(plex,0,&cStart,&cEnd);CHKERRQ(ierr);
1315*e0eea495SMark   numGCells = cEnd - cStart;
1316*e0eea495SMark   /* create cell centered DM */
1317*e0eea495SMark   ierr = DMClone(plex, &colordm);CHKERRQ(ierr);
1318*e0eea495SMark   ierr = PetscFECreateDefault(PetscObjectComm((PetscObject) plex), dim, 1, PETSC_FALSE, "color_", PETSC_DECIDE, &fe);CHKERRQ(ierr);
1319*e0eea495SMark   ierr = PetscObjectSetName((PetscObject) fe, "color");CHKERRQ(ierr);
1320*e0eea495SMark   ierr = DMSetField(colordm, 0, NULL, (PetscObject)fe);CHKERRQ(ierr);
1321*e0eea495SMark   ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
1322*e0eea495SMark   for (i = 0; i < (dim+1); ++i) numDof[i] = 0;
1323*e0eea495SMark   numDof[dim] = 1;
1324*e0eea495SMark   numComp[0] = 1;
1325*e0eea495SMark   ierr = DMPlexCreateSection(colordm, NULL, numComp, numDof, 0, NULL, NULL, NULL, NULL, &csection);CHKERRQ(ierr);
1326*e0eea495SMark   ierr = PetscSectionSetFieldName(csection, 0, "color");CHKERRQ(ierr);
1327*e0eea495SMark   ierr = DMSetLocalSection(colordm, csection);CHKERRQ(ierr);
1328*e0eea495SMark   ierr = DMViewFromOptions(colordm,NULL,"-color_dm_view");CHKERRQ(ierr);
1329*e0eea495SMark   /* get vertex to element map Q and colroing graph G */
1330*e0eea495SMark   ierr = MatGetSize(JacP,NULL,&Nv);CHKERRQ(ierr);
1331*e0eea495SMark   ierr = MatCreateAIJ(PETSC_COMM_SELF,PETSC_DECIDE,PETSC_DECIDE,numGCells,Nv,totDim,NULL,0,NULL,&Q);CHKERRQ(ierr);
1332*e0eea495SMark   for (i=0;i<128;i++) ones[i] = 1.0;
1333*e0eea495SMark   for (cell = cStart, ej = 0 ; cell < cEnd; ++cell, ++ej) {
1334*e0eea495SMark     PetscInt numindices,*indices;
1335*e0eea495SMark     ierr = DMPlexGetClosureIndices(plex, section, globalSection, cell, PETSC_TRUE, &numindices, &indices, NULL, NULL);CHKERRQ(ierr);
1336*e0eea495SMark     if (numindices>128) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "too many indices. %D > %D",numindices,128);
1337*e0eea495SMark     ierr = MatSetValues(Q,1,&ej,numindices,indices,ones,ADD_VALUES);CHKERRQ(ierr);
1338*e0eea495SMark     ierr = DMPlexRestoreClosureIndices(plex, section, globalSection, cell, PETSC_TRUE, &numindices, &indices, NULL, NULL);CHKERRQ(ierr);
1339*e0eea495SMark   }
1340*e0eea495SMark   ierr = MatAssemblyBegin(Q, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1341*e0eea495SMark   ierr = MatAssemblyEnd(Q, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1342*e0eea495SMark   ierr = MatMatTransposeMult(Q,Q,MAT_INITIAL_MATRIX,4.0,&G);CHKERRQ(ierr);
1343*e0eea495SMark   ierr = PetscObjectSetName((PetscObject) Q, "Q");CHKERRQ(ierr);
1344*e0eea495SMark   ierr = PetscObjectSetName((PetscObject) G, "coloring graph");CHKERRQ(ierr);
1345*e0eea495SMark   ierr = MatViewFromOptions(G,NULL,"-coloring_mat_view");CHKERRQ(ierr);
1346*e0eea495SMark   ierr = MatViewFromOptions(Q,NULL,"-coloring_mat_view");CHKERRQ(ierr);
1347*e0eea495SMark   ierr = MatDestroy(&Q);CHKERRQ(ierr);
1348*e0eea495SMark   /* coloring */
1349*e0eea495SMark   ierr = MatColoringCreate(G,&mc);CHKERRQ(ierr);
1350*e0eea495SMark   ierr = MatColoringSetDistance(mc,1);CHKERRQ(ierr);
1351*e0eea495SMark   ierr = MatColoringSetType(mc,MATCOLORINGJP);CHKERRQ(ierr);
1352*e0eea495SMark   ierr = MatColoringSetFromOptions(mc);CHKERRQ(ierr);
1353*e0eea495SMark   ierr = MatColoringApply(mc,&iscoloring);CHKERRQ(ierr);
1354*e0eea495SMark   ierr = MatColoringDestroy(&mc);CHKERRQ(ierr);
1355*e0eea495SMark   /* view */
1356*e0eea495SMark   ierr = ISColoringViewFromOptions(iscoloring,NULL,"-coloring_is_view");CHKERRQ(ierr);
1357*e0eea495SMark   ierr = ISColoringGetIS(iscoloring,PETSC_USE_POINTER,&nc,&is);CHKERRQ(ierr);
1358*e0eea495SMark   if (ctx && ctx->verbose > 5) {
1359*e0eea495SMark     PetscViewer    viewer;
1360*e0eea495SMark     Vec            color_vec, eidx_vec;
1361*e0eea495SMark     ierr = DMGetGlobalVector(colordm, &color_vec);CHKERRQ(ierr);
1362*e0eea495SMark     ierr = DMGetGlobalVector(colordm, &eidx_vec);CHKERRQ(ierr);
1363*e0eea495SMark     for (colour=0; colour<nc; colour++) {
1364*e0eea495SMark       ierr = ISGetLocalSize(is[colour],&csize);CHKERRQ(ierr);
1365*e0eea495SMark       ierr = ISGetIndices(is[colour],&indices);CHKERRQ(ierr);
1366*e0eea495SMark       for (j=0; j<csize; j++) {
1367*e0eea495SMark         PetscScalar v = (PetscScalar)colour;
1368*e0eea495SMark         k = indices[j];
1369*e0eea495SMark         ierr = VecSetValues(color_vec,1,&k,&v,INSERT_VALUES);
1370*e0eea495SMark         v = (PetscScalar)k;
1371*e0eea495SMark         ierr = VecSetValues(eidx_vec,1,&k,&v,INSERT_VALUES);
1372*e0eea495SMark       }
1373*e0eea495SMark       ierr = ISRestoreIndices(is[colour],&indices);CHKERRQ(ierr);
1374*e0eea495SMark     }
1375*e0eea495SMark     /* view */
1376*e0eea495SMark     ierr = PetscViewerCreate(PETSC_COMM_WORLD, &viewer);CHKERRQ(ierr);
1377*e0eea495SMark     ierr = PetscViewerSetType(viewer, PETSCVIEWERVTK);CHKERRQ(ierr);
1378*e0eea495SMark     ierr = PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_VTK);CHKERRQ(ierr);
1379*e0eea495SMark     ierr = PetscViewerFileSetName(viewer, "color.vtk");CHKERRQ(ierr);
1380*e0eea495SMark     ierr = PetscObjectSetName((PetscObject) color_vec, "color");CHKERRQ(ierr);
1381*e0eea495SMark     ierr = VecView(color_vec, viewer);CHKERRQ(ierr);
1382*e0eea495SMark     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1383*e0eea495SMark     ierr = PetscViewerCreate(PETSC_COMM_WORLD, &viewer);CHKERRQ(ierr);
1384*e0eea495SMark     ierr = PetscViewerSetType(viewer, PETSCVIEWERVTK);CHKERRQ(ierr);
1385*e0eea495SMark     ierr = PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_VTK);CHKERRQ(ierr);
1386*e0eea495SMark     ierr = PetscViewerFileSetName(viewer, "eidx.vtk");CHKERRQ(ierr);
1387*e0eea495SMark     ierr = PetscObjectSetName((PetscObject) eidx_vec, "element-idx");CHKERRQ(ierr);
1388*e0eea495SMark     ierr = VecView(eidx_vec, viewer);CHKERRQ(ierr);
1389*e0eea495SMark     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1390*e0eea495SMark     ierr = DMRestoreGlobalVector(colordm, &color_vec);CHKERRQ(ierr);
1391*e0eea495SMark     ierr = DMRestoreGlobalVector(colordm, &eidx_vec);CHKERRQ(ierr);
1392*e0eea495SMark   }
1393*e0eea495SMark   ierr = PetscSectionDestroy(&csection);CHKERRQ(ierr);
1394*e0eea495SMark   ierr = DMDestroy(&colordm);CHKERRQ(ierr);
1395*e0eea495SMark   ierr = ISColoringRestoreIS(iscoloring,PETSC_USE_POINTER,&is);CHKERRQ(ierr);
1396*e0eea495SMark   ierr = MatDestroy(&G);CHKERRQ(ierr);
1397*e0eea495SMark   /* stash coloring */
1398*e0eea495SMark   ierr = PetscContainerCreate(PETSC_COMM_SELF, container);CHKERRQ(ierr);
1399*e0eea495SMark   ierr = PetscContainerSetPointer(*container,(void*)iscoloring);CHKERRQ(ierr);
1400*e0eea495SMark   ierr = PetscContainerSetUserDestroy(*container, destroy_coloring);CHKERRQ(ierr);
1401*e0eea495SMark   ierr = PetscObjectCompose((PetscObject)JacP,"coloring",(PetscObject)*container);CHKERRQ(ierr);
1402*e0eea495SMark   if (ctx && ctx->verbose > 0) {
1403*e0eea495SMark     ierr = PetscPrintf(PETSC_COMM_WORLD, "Made coloring with %D colors\n", nc);CHKERRQ(ierr);
1404*e0eea495SMark   }
1405*e0eea495SMark   PetscFunctionReturn(0);
1406*e0eea495SMark   }
1407*e0eea495SMark 
1408*e0eea495SMark PetscErrorCode LandauAssembleOpenMP(PetscInt cStart, PetscInt cEnd, PetscInt totDim, DM plex, PetscSection section, PetscSection globalSection, Mat JacP, PetscScalar elemMats[], PetscContainer container)
1409*e0eea495SMark {
1410*e0eea495SMark   PetscErrorCode  ierr;
1411*e0eea495SMark   IS             *is;
1412*e0eea495SMark   PetscInt        nc,colour,j;
1413*e0eea495SMark   const PetscInt *clr_idxs;
1414*e0eea495SMark   ISColoring      iscoloring;
1415*e0eea495SMark   PetscFunctionBegin;
1416*e0eea495SMark   ierr = PetscContainerGetPointer(container,(void**)&iscoloring);CHKERRQ(ierr);
1417*e0eea495SMark   ierr = ISColoringGetIS(iscoloring,PETSC_USE_POINTER,&nc,&is);CHKERRQ(ierr);
1418*e0eea495SMark   for (colour=0; colour<nc; colour++) {
1419*e0eea495SMark     PetscInt    *idx_arr[1024]; /* need to make dynamic for general use */
1420*e0eea495SMark     PetscScalar *new_el_mats[1024];
1421*e0eea495SMark     PetscInt     idx_size[1024],csize;
1422*e0eea495SMark     ierr = ISGetLocalSize(is[colour],&csize);CHKERRQ(ierr);
1423*e0eea495SMark     if (csize>1024) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "too many elements in color. %D > %D",csize,1024);
1424*e0eea495SMark     ierr = ISGetIndices(is[colour],&clr_idxs);CHKERRQ(ierr);
1425*e0eea495SMark     /* get indices and mats */
1426*e0eea495SMark     for (j=0; j<csize; j++) {
1427*e0eea495SMark       PetscInt    cell = cStart + clr_idxs[j];
1428*e0eea495SMark       PetscInt    numindices,*indices;
1429*e0eea495SMark       PetscScalar *elMat = &elemMats[clr_idxs[j]*totDim*totDim];
1430*e0eea495SMark       PetscScalar *valuesOrig = elMat;
1431*e0eea495SMark       ierr = DMPlexGetClosureIndices(plex, section, globalSection, cell, PETSC_TRUE, &numindices, &indices, NULL, (PetscScalar **) &elMat);CHKERRQ(ierr);
1432*e0eea495SMark       idx_size[j] = numindices;
1433*e0eea495SMark       ierr = PetscMalloc2(numindices,&idx_arr[j],numindices*numindices,&new_el_mats[j]);CHKERRQ(ierr);
1434*e0eea495SMark       ierr = PetscMemcpy(idx_arr[j],indices,numindices*sizeof(PetscInt));CHKERRQ(ierr);
1435*e0eea495SMark       ierr = PetscMemcpy(new_el_mats[j],elMat,numindices*numindices*sizeof(PetscScalar));CHKERRQ(ierr);
1436*e0eea495SMark       ierr = DMPlexRestoreClosureIndices(plex, section, globalSection, cell, PETSC_TRUE, &numindices, &indices, NULL, (PetscScalar **) &elMat);CHKERRQ(ierr);
1437*e0eea495SMark       if (elMat != valuesOrig) {ierr = DMRestoreWorkArray(plex, numindices*numindices, MPIU_SCALAR, &elMat);}
1438*e0eea495SMark     }
1439*e0eea495SMark     /* assemble matrix - pragmas break CI ? */
1440*e0eea495SMark     //#pragma omp parallel default(JacP,idx_size,idx_arr,new_el_mats,colour,clr_idxs)  private(j)
1441*e0eea495SMark     //#pragma omp parallel for private(j)
1442*e0eea495SMark     for (j=0; j<csize; j++) {
1443*e0eea495SMark       PetscInt    numindices = idx_size[j], *indices = idx_arr[j];
1444*e0eea495SMark       PetscScalar *elMat = new_el_mats[j];
1445*e0eea495SMark       MatSetValues(JacP,numindices,indices,numindices,indices,elMat,ADD_VALUES);
1446*e0eea495SMark     }
1447*e0eea495SMark     /* free */
1448*e0eea495SMark     ierr = ISRestoreIndices(is[colour],&clr_idxs);CHKERRQ(ierr);
1449*e0eea495SMark     for (j=0; j<csize; j++) {
1450*e0eea495SMark       ierr = PetscFree2(idx_arr[j],new_el_mats[j]);CHKERRQ(ierr);
1451*e0eea495SMark     }
1452*e0eea495SMark   }
1453*e0eea495SMark   ierr = ISColoringRestoreIS(iscoloring,PETSC_USE_POINTER,&is);CHKERRQ(ierr);
1454*e0eea495SMark   PetscFunctionReturn(0);
1455*e0eea495SMark }
1456*e0eea495SMark 
1457*e0eea495SMark /* < v, u > */
1458*e0eea495SMark static void g0_1(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1459*e0eea495SMark                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1460*e0eea495SMark                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1461*e0eea495SMark                   PetscReal t, PetscReal u_tShift, const PetscReal x[],  PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
1462*e0eea495SMark {
1463*e0eea495SMark   g0[0] = 1.;
1464*e0eea495SMark }
1465*e0eea495SMark 
1466*e0eea495SMark /* < v, u > */
1467*e0eea495SMark static void g0_r(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1468*e0eea495SMark                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1469*e0eea495SMark                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1470*e0eea495SMark                   PetscReal t, PetscReal u_tShift, const PetscReal x[],  PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
1471*e0eea495SMark {
1472*e0eea495SMark   g0[0] = 2.*PETSC_PI*x[0];
1473*e0eea495SMark }
1474*e0eea495SMark 
1475*e0eea495SMark /*@
1476*e0eea495SMark   LandauCreateMassMatrix - Create mass matrix for Landau
1477*e0eea495SMark 
1478*e0eea495SMark   Collective on dm
1479*e0eea495SMark 
1480*e0eea495SMark   Input Parameters:
1481*e0eea495SMark + dm     - the DM object
1482*e0eea495SMark 
1483*e0eea495SMark  Output Parameters:
1484*e0eea495SMark + Amat - The mass matrix (optional), mass matrix is added to the DM context
1485*e0eea495SMark 
1486*e0eea495SMark   Level: beginner
1487*e0eea495SMark 
1488*e0eea495SMark .keywords: mesh
1489*e0eea495SMark .seealso: LandauCreateVelocitySpace()
1490*e0eea495SMark @*/
1491*e0eea495SMark PetscErrorCode LandauCreateMassMatrix(DM dm, Mat *Amat)
1492*e0eea495SMark {
1493*e0eea495SMark   DM             massDM;
1494*e0eea495SMark   PetscDS        prob;
1495*e0eea495SMark   PetscInt       ii,dim,N1=1,N2;
1496*e0eea495SMark   PetscErrorCode ierr;
1497*e0eea495SMark   LandauCtx      *ctx;
1498*e0eea495SMark   Mat            M;
1499*e0eea495SMark 
1500*e0eea495SMark   PetscFunctionBegin;
1501*e0eea495SMark   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1502*e0eea495SMark   if (Amat) PetscValidPointer(Amat,3);
1503*e0eea495SMark   ierr = DMGetApplicationContext(dm, &ctx);CHKERRQ(ierr);
1504*e0eea495SMark   if (!ctx) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "no context");
1505*e0eea495SMark   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1506*e0eea495SMark   ierr = DMClone(dm, &massDM);CHKERRQ(ierr);
1507*e0eea495SMark   ierr = DMCopyFields(dm, massDM);CHKERRQ(ierr);
1508*e0eea495SMark   ierr = DMCreateDS(massDM);CHKERRQ(ierr);
1509*e0eea495SMark   ierr = DMGetDS(massDM, &prob);CHKERRQ(ierr);
1510*e0eea495SMark   for (ii=0;ii<ctx->num_species;ii++) {
1511*e0eea495SMark     if (dim==3) {ierr = PetscDSSetJacobian(prob, ii, ii, g0_1, NULL, NULL, NULL);CHKERRQ(ierr);}
1512*e0eea495SMark     else        {ierr = PetscDSSetJacobian(prob, ii, ii, g0_r, NULL, NULL, NULL);CHKERRQ(ierr);}
1513*e0eea495SMark   }
1514*e0eea495SMark   ierr = DMViewFromOptions(massDM,NULL,"-dm_landau_mass_dm_view");CHKERRQ(ierr);
1515*e0eea495SMark   ierr = DMCreateMatrix(massDM, &M);CHKERRQ(ierr);
1516*e0eea495SMark   {
1517*e0eea495SMark     Vec locX;
1518*e0eea495SMark     DM  plex;
1519*e0eea495SMark     ierr = DMConvert(massDM, DMPLEX, &plex);CHKERRQ(ierr);
1520*e0eea495SMark     ierr = DMGetLocalVector(massDM, &locX);CHKERRQ(ierr);
1521*e0eea495SMark     /* Mass matrix is independent of the input, so no need to fill locX */
1522*e0eea495SMark     ierr = DMPlexSNESComputeJacobianFEM(plex, locX, M, M, ctx);CHKERRQ(ierr);
1523*e0eea495SMark     ierr = DMRestoreLocalVector(massDM, &locX);CHKERRQ(ierr);
1524*e0eea495SMark     ierr = DMDestroy(&plex);CHKERRQ(ierr);
1525*e0eea495SMark   }
1526*e0eea495SMark   ierr = DMDestroy(&massDM);CHKERRQ(ierr);
1527*e0eea495SMark   ierr = MatGetSize(ctx->J, &N1, NULL);CHKERRQ(ierr);
1528*e0eea495SMark   ierr = MatGetSize(M, &N2, NULL);CHKERRQ(ierr);
1529*e0eea495SMark   if (N1 != N2) SETERRQ2(PetscObjectComm((PetscObject) dm), PETSC_ERR_PLIB, "Incorrect matrix sizes: |Jacobian| = %D, |Mass|=%D",N1,N2);
1530*e0eea495SMark   ierr = MatViewFromOptions(M,NULL,"-dm_landau_mass_mat_view");CHKERRQ(ierr);
1531*e0eea495SMark   ctx->M = M; /* this could be a noop, a = a */
1532*e0eea495SMark   if (Amat) *Amat = M;
1533*e0eea495SMark   PetscFunctionReturn(0);
1534*e0eea495SMark }
1535*e0eea495SMark 
1536*e0eea495SMark /*@
1537*e0eea495SMark   LandauIFunction - TS residual calculation
1538*e0eea495SMark 
1539*e0eea495SMark   Collective on ts
1540*e0eea495SMark 
1541*e0eea495SMark   Input Parameters:
1542*e0eea495SMark +   TS  - The time stepping context
1543*e0eea495SMark .   time_dummy - current time (not used)
1544*e0eea495SMark -   X - Current state
1545*e0eea495SMark +   X_t - Time derivative of current state
1546*e0eea495SMark .   actx - Landau context
1547*e0eea495SMark 
1548*e0eea495SMark   Output Parameter:
1549*e0eea495SMark .   F  - The residual
1550*e0eea495SMark 
1551*e0eea495SMark   Level: beginner
1552*e0eea495SMark 
1553*e0eea495SMark .keywords: mesh
1554*e0eea495SMark .seealso: LandauCreateVelocitySpace(), LandauIJacobian()
1555*e0eea495SMark @*/
1556*e0eea495SMark PetscErrorCode LandauIFunction(TS ts,PetscReal time_dummy,Vec X,Vec X_t,Vec F,void *actx)
1557*e0eea495SMark {
1558*e0eea495SMark   PetscErrorCode ierr;
1559*e0eea495SMark   LandauCtx      *ctx=(LandauCtx*)actx;
1560*e0eea495SMark   PetscReal      unorm;
1561*e0eea495SMark   PetscInt       dim;
1562*e0eea495SMark   DM dm;
1563*e0eea495SMark 
1564*e0eea495SMark   PetscFunctionBegin;
1565*e0eea495SMark   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
1566*e0eea495SMark   ierr = DMGetApplicationContext(dm, &ctx);CHKERRQ(ierr);
1567*e0eea495SMark   if (!ctx) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "no context");
1568*e0eea495SMark   ierr = VecNorm(X,NORM_2,&unorm);CHKERRQ(ierr);
1569*e0eea495SMark   ierr = MPI_Barrier(PETSC_COMM_WORLD);CHKERRQ(ierr); // remove in real application
1570*e0eea495SMark   ierr = PetscLogEventBegin(ctx->events[0],0,0,0,0);CHKERRQ(ierr);
1571*e0eea495SMark   ierr = DMGetDimension(ctx->dmv, &dim);CHKERRQ(ierr);
1572*e0eea495SMark   if (ctx->normJ!=unorm) {
1573*e0eea495SMark     ctx->normJ = unorm;
1574*e0eea495SMark     ierr = LandauFormJacobian_Internal(X,ctx->J,dim,(void*)ctx);CHKERRQ(ierr);
1575*e0eea495SMark     ctx->aux_bool = PETSC_TRUE; /* debug: set flag that we made a new Jacobian */
1576*e0eea495SMark   } else ctx->aux_bool = PETSC_FALSE;
1577*e0eea495SMark   /* mat vec for op */
1578*e0eea495SMark   ierr = MatMult(ctx->J,X,F);CHKERRQ(ierr);CHKERRQ(ierr); /* C*f */
1579*e0eea495SMark   /* add time term */
1580*e0eea495SMark   if (X_t) {
1581*e0eea495SMark     ierr = MatMultAdd(ctx->M,X_t,F,F);CHKERRQ(ierr);
1582*e0eea495SMark   }
1583*e0eea495SMark   ierr = PetscLogEventEnd(ctx->events[0],0,0,0,0);CHKERRQ(ierr);
1584*e0eea495SMark   PetscFunctionReturn(0);
1585*e0eea495SMark }
1586*e0eea495SMark 
1587*e0eea495SMark /*@
1588*e0eea495SMark   LandauIJacobian - TS Jacobian construction
1589*e0eea495SMark 
1590*e0eea495SMark   Collective on ts
1591*e0eea495SMark 
1592*e0eea495SMark   Input Parameters:
1593*e0eea495SMark +   TS  - The time stepping context
1594*e0eea495SMark .   time_dummy - current time (not used)
1595*e0eea495SMark -   X - Current state
1596*e0eea495SMark +   U_tdummy - Time derivative of current state (not used)
1597*e0eea495SMark .   shift - shift for du/dt term
1598*e0eea495SMark -   actx - Landau context
1599*e0eea495SMark 
1600*e0eea495SMark   Output Parameter:
1601*e0eea495SMark .   Amat  - Jacobian
1602*e0eea495SMark +   Pmat  - same as Amat
1603*e0eea495SMark 
1604*e0eea495SMark   Level: beginner
1605*e0eea495SMark 
1606*e0eea495SMark .keywords: mesh
1607*e0eea495SMark .seealso: LandauCreateVelocitySpace(), LandauIFunction()
1608*e0eea495SMark @*/
1609*e0eea495SMark PetscErrorCode LandauIJacobian(TS ts,PetscReal time_dummy,Vec X,Vec U_tdummy,PetscReal shift,Mat Amat,Mat Pmat,void *actx)
1610*e0eea495SMark {
1611*e0eea495SMark   PetscErrorCode ierr;
1612*e0eea495SMark   LandauCtx      *ctx=(LandauCtx*)actx;
1613*e0eea495SMark   PetscReal      unorm;
1614*e0eea495SMark   PetscInt       dim;
1615*e0eea495SMark   DM dm;
1616*e0eea495SMark 
1617*e0eea495SMark   PetscFunctionBegin;
1618*e0eea495SMark   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
1619*e0eea495SMark   ierr = DMGetApplicationContext(dm, &ctx);CHKERRQ(ierr);
1620*e0eea495SMark   if (!ctx) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "no context");
1621*e0eea495SMark   if (Amat!=Pmat || Amat!=ctx->J) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Amat!=Pmat || Amat!=ctx->J");
1622*e0eea495SMark   ierr = DMGetDimension(ctx->dmv, &dim);CHKERRQ(ierr);
1623*e0eea495SMark   /* get collision Jacobian into A */
1624*e0eea495SMark   ierr = PetscLogEventBegin(ctx->events[9],0,0,0,0);CHKERRQ(ierr);
1625*e0eea495SMark   ierr = VecNorm(X,NORM_2,&unorm);CHKERRQ(ierr);
1626*e0eea495SMark   if (ctx->normJ!=unorm) {
1627*e0eea495SMark     ierr = LandauFormJacobian_Internal(X,ctx->J,dim,(void*)ctx);CHKERRQ(ierr);
1628*e0eea495SMark     ctx->normJ = unorm;
1629*e0eea495SMark     ctx->aux_bool = PETSC_TRUE; /* debug: set flag that we made a new Jacobian */
1630*e0eea495SMark   } else ctx->aux_bool = PETSC_FALSE;
1631*e0eea495SMark   /* add C */
1632*e0eea495SMark   ierr = MatCopy(ctx->J,Pmat,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
1633*e0eea495SMark   /* add mass */
1634*e0eea495SMark   ierr = MatAXPY(Pmat,shift,ctx->M,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
1635*e0eea495SMark   ierr = PetscLogEventEnd(ctx->events[9],0,0,0,0);CHKERRQ(ierr);
1636*e0eea495SMark   PetscFunctionReturn(0);
1637*e0eea495SMark }
1638