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