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