xref: /petsc/src/ts/utils/dmplexlandau/plexland.c (revision ebcb266d3a667e450268eb330045d5d1ecf6fdff)
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_tensors.h"
10 
11 /* vector padding not supported */
12 #define LANDAU_VL  1
13 
14 int LandauGetIPDataSize(const LandauIPData *const d) {
15   return d->nip_*(1 + d->dim_ + d->ns_); /* assumes Nq == Nd */
16 }
17 
18 static PetscErrorCode LandauPointDataCreate(LandauIPData *IPData, PetscInt dim, PetscInt nip, PetscInt Ns)
19 {
20   PetscErrorCode  ierr;
21   PetscInt        sz, nip_pad = nip ; /* LANDAU_VL*(nip/LANDAU_VL + !!(nip%LANDAU_VL)); */
22   LandauIPReal    *pdata;
23   PetscFunctionBegin;
24   IPData->dim_ = dim;
25   IPData->nip_ = nip_pad;
26   IPData->ns_  = Ns;
27   sz = LandauGetIPDataSize(IPData);
28   ierr = PetscMalloc(sizeof(LandauIPReal)*sz,&pdata);CHKERRQ(ierr);
29   /* pack data */
30   IPData->w    = pdata + 0; /* w */
31   IPData->x    = pdata + 1*nip_pad;
32   IPData->y    = pdata + 2*nip_pad;
33   IPData->z    = pdata + 3*nip_pad;
34   IPData->coefs= pdata + (dim+1)*nip_pad;
35   PetscFunctionReturn(0);
36 }
37 static PetscErrorCode LandauGPUDataDestroy(void *ptr)
38 {
39   P4estVertexMaps *maps = (P4estVertexMaps *)ptr;
40   PetscErrorCode  ierr;
41   PetscFunctionBegin;
42   if (maps->deviceType != LANDAU_CPU) {
43 #if defined(PETSC_HAVE_KOKKOS_KERNELS)
44     if (maps->deviceType == LANDAU_KOKKOS) {
45       ierr = LandauKokkosDestroyMatMaps(maps);CHKERRQ(ierr); // imples Kokkos does
46     } // else could be CUDA
47 #elif defined(PETSC_HAVE_CUDA)
48     if (maps->deviceType == LANDAU_CUDA){
49       ierr = LandauCUDADestroyMatMaps(maps);CHKERRQ(ierr);
50     } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "maps->deviceType %D ?????",maps->deviceType);
51 #endif
52   }
53   ierr = PetscFree(maps->c_maps);CHKERRQ(ierr);
54   ierr = PetscFree(maps->gIdx);CHKERRQ(ierr);
55   ierr = PetscFree(maps);CHKERRQ(ierr);
56   PetscFunctionReturn(0);
57 }
58 static PetscErrorCode LandauPointDataDestroy(LandauIPData *IPData)
59 {
60   PetscErrorCode   ierr;
61   PetscFunctionBegin;
62   ierr = PetscFree(IPData->w);CHKERRQ(ierr);
63   PetscFunctionReturn(0);
64 }
65 /* ------------------------------------------------------------------- */
66 /*
67  LandauFormJacobian_Internal - Evaluates Jacobian matrix.
68 
69  Input Parameters:
70  .  globX - input vector
71  .  actx - optional user-defined context
72  .  dim - dimension
73 
74  Output Parameters:
75  .  J0acP - Jacobian matrix filled, not created
76  */
77 static PetscErrorCode LandauFormJacobian_Internal(Vec a_X, Mat JacP, const PetscInt dim, PetscReal shift, void *a_ctx)
78 {
79   LandauCtx         *ctx = (LandauCtx*)a_ctx;
80   PetscErrorCode    ierr;
81   PetscInt          cStart, cEnd, elemMatSize;
82   DM                plex = NULL;
83   PetscDS           prob;
84   PetscSection      section,globsection;
85   PetscInt          numCells,totDim,ej,Nq,*Nbf,*Ncf,Nb,Ncx,Nf,d,f,fieldA,qj;
86   PetscQuadrature   quad;
87   PetscTabulation   *Tf;
88   PetscReal         nu_alpha[LANDAU_MAX_SPECIES], nu_beta[LANDAU_MAX_SPECIES];
89   const PetscReal   *quadWeights;
90   PetscReal         *invJ,*invJ_a=NULL,*mass_w=NULL;
91   PetscReal         invMass[LANDAU_MAX_SPECIES],Eq_m[LANDAU_MAX_SPECIES],m_0=ctx->m_0; /* normalize mass -- not needed! */
92   PetscLogDouble    flops;
93   Vec               locX;
94   LandauIPData      IPData;
95   PetscContainer    container;
96   P4estVertexMaps   *maps=NULL;
97 
98   PetscFunctionBegin;
99   PetscValidHeaderSpecific(a_X,VEC_CLASSID,1);
100   PetscValidHeaderSpecific(JacP,MAT_CLASSID,2);
101   PetscValidPointer(ctx,4);
102 
103   /* check for matrix container for GPU assembly */
104   ierr = PetscLogEventBegin(ctx->events[10],0,0,0,0);CHKERRQ(ierr);
105   ierr = PetscObjectQuery((PetscObject) JacP, "assembly_maps", (PetscObject *) &container);CHKERRQ(ierr);
106   if (container /* && ctx->deviceType != LANDAU_CPU */) {
107     if (!ctx->gpu_assembly) SETERRQ(ctx->comm,PETSC_ERR_ARG_WRONG,"GPU matrix container but no GPU assembly");
108     ierr = PetscContainerGetPointer(container, (void **) &maps);CHKERRQ(ierr);
109     if (!maps) SETERRQ(ctx->comm,PETSC_ERR_ARG_WRONG,"empty GPU matrix container");
110   }
111   ierr = DMConvert(ctx->dmv, DMPLEX, &plex);CHKERRQ(ierr);
112   ierr = DMCreateLocalVector(plex, &locX);CHKERRQ(ierr);
113   ierr = VecZeroEntries(locX);CHKERRQ(ierr); /* zero BCs so don't set */
114   ierr = DMGlobalToLocalBegin(plex, a_X, INSERT_VALUES, locX);CHKERRQ(ierr);
115   ierr = DMGlobalToLocalEnd  (plex, a_X, INSERT_VALUES, locX);CHKERRQ(ierr);
116   ierr = DMPlexGetHeightStratum(plex, 0, &cStart, &cEnd);CHKERRQ(ierr);
117   ierr = DMGetLocalSection(plex, &section);CHKERRQ(ierr);
118   ierr = DMGetGlobalSection(plex, &globsection);CHKERRQ(ierr);
119   ierr = DMGetDS(plex, &prob);CHKERRQ(ierr);
120   ierr = PetscDSGetTabulation(prob, &Tf);CHKERRQ(ierr); // Bf, &Df
121   ierr = PetscDSGetDimensions(prob, &Nbf);CHKERRQ(ierr); Nb = Nbf[0]; /* number of vertices*S */
122   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);         if (Nf!=ctx->num_species) SETERRQ1(ctx->comm, PETSC_ERR_PLIB, "Nf %D != S",Nf);
123   ierr = PetscDSGetComponents(prob, &Ncf);CHKERRQ(ierr); Ncx = Ncf[0]; if (Ncx!=1) SETERRQ1(ctx->comm, PETSC_ERR_PLIB, "Nc %D != 1",Ncx);
124   if (shift==0.0) {
125     for (fieldA=0;fieldA<Nf;fieldA++) {
126       invMass[fieldA] = m_0/ctx->masses[fieldA];
127       Eq_m[fieldA] = ctx->Ez * ctx->t_0 * ctx->charges[fieldA] / (ctx->v_0 * ctx->masses[fieldA]); /* normalize dimensionless */
128       if (dim==2) Eq_m[fieldA] *=  2 * PETSC_PI; /* add the 2pi term that is not in Landau */
129       nu_alpha[fieldA] = PetscSqr(ctx->charges[fieldA]/m_0)*m_0/ctx->masses[fieldA];
130       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);
131     }
132   }
133   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
134   numCells = cEnd - cStart;
135   ierr = PetscFEGetQuadrature(ctx->fe[0], &quad);CHKERRQ(ierr);
136   ierr = PetscQuadratureGetData(quad, NULL, NULL, &Nq, NULL, &quadWeights);CHKERRQ(ierr);
137   if (Nb!=Nq) SETERRQ4(ctx->comm, PETSC_ERR_PLIB, "Nb!=Nq %D %D over integration or simplices? Tf[0]->Nb=%D dim=%D",Nb,Nq,Tf[0]->Nb,dim);
138   if (Nq >LANDAU_MAX_NQ) SETERRQ2(ctx->comm,PETSC_ERR_ARG_WRONG,"Order too high. Nq = %D > LANDAU_MAX_NQ (%D)",Nq,LANDAU_MAX_NQ);
139   if (LANDAU_DIM != dim) SETERRQ2(ctx->comm, PETSC_ERR_PLIB, "dim %D != LANDAU_DIM %d",dim,LANDAU_DIM);
140   if (shift==0.0) {
141     ierr = MatZeroEntries(JacP);CHKERRQ(ierr);
142     flops = (PetscLogDouble)numCells*(PetscLogDouble)Nq*(PetscLogDouble)(5*dim*dim*Nf*Nf + 165);
143   } else {
144   flops = (PetscLogDouble)numCells*(PetscLogDouble)Nq*(PetscLogDouble)(5*dim*dim*Nf*Nf);
145   }
146   elemMatSize = totDim*totDim;
147   ierr = PetscLogEventEnd(ctx->events[10],0,0,0,0);CHKERRQ(ierr);
148   {
149     static int         cc = 0;
150     /* collect f data */
151     if (ctx->verbose > 1 || (ctx->verbose > 0 && cc++ == 0)) {
152       PetscInt N,Nloc;
153       ierr = MatGetSize(JacP,&N,NULL);CHKERRQ(ierr);
154       ierr = VecGetSize(locX,&Nloc);CHKERRQ(ierr);
155       ierr = PetscPrintf(ctx->comm,"[%D]%s: %D IPs, %D cells, totDim=%D, Nb=%D, Nq=%D, elemMatSize=%D, dim=%D, Tab: Nb=%D Nf=%D Np=%D cdim=%D N=%D N+hang=%D, shift=%g\n",
156                          0,"FormLandau",Nq*numCells,numCells, totDim, Nb, Nq, elemMatSize, dim, Tf[0]->Nb, Nf, Tf[0]->Np, Tf[0]->cdim, N, Nloc, shift);CHKERRQ(ierr);
157     }
158     if (shift==0.0) {
159       ierr = LandauPointDataCreate(&IPData, dim, Nq*numCells, Nf);CHKERRQ(ierr);
160       ierr = PetscMalloc1(numCells*Nq*dim*dim,&invJ_a);CHKERRQ(ierr);
161       ierr = PetscLogEventBegin(ctx->events[7],0,0,0,0);CHKERRQ(ierr);
162     } else { // mass
163       ierr = PetscMalloc1(numCells*Nq,&mass_w);CHKERRQ(ierr);
164       IPData.w = NULL;
165       ierr = PetscLogEventBegin(ctx->events[1],0,0,0,0);CHKERRQ(ierr);
166     }
167     /* cache geometry and x, f and df/dx at IPs */
168     for (ej = 0 ; ej < numCells; ++ej) {
169       PetscReal    vj[LANDAU_MAX_NQ*LANDAU_DIM],detJj[LANDAU_MAX_NQ], Jdummy[LANDAU_MAX_NQ*LANDAU_DIM*LANDAU_DIM];
170       PetscScalar *coef = NULL;
171       invJ = invJ_a ? invJ_a + ej * Nq*dim*dim : NULL;
172       ierr = DMPlexComputeCellGeometryFEM(plex, cStart+ej, quad, vj, Jdummy, invJ, detJj);CHKERRQ(ierr);
173       if (shift!=0.0) { // mass
174         for (qj = 0; qj < Nq; ++qj) {
175           PetscInt         gidx = (ej*Nq + qj);
176           mass_w[gidx] = detJj[qj] * quadWeights[qj];
177           if (dim==2) mass_w[gidx] *=  2.*PETSC_PI*vj[qj * dim + 0]; /* cylindrical coordinate, w/o 2pi */
178         }
179       } else {
180         ierr = DMPlexVecGetClosure(plex, section, locX, cStart+ej, NULL, &coef);CHKERRQ(ierr);
181         ierr = PetscMemcpy(&IPData.coefs[ej*Nb*Nf],coef,Nb*Nf*sizeof(PetscScalar));CHKERRQ(ierr); /* change if LandauIPReal != PetscScalar */
182         /* create point data for cell i for Landau tensor: x, f(x), grad f(x) */
183         for (qj = 0; qj < Nq; ++qj) {
184           PetscInt         gidx = (ej*Nq + qj);
185           IPData.x[gidx] = vj[qj * dim + 0]; /* coordinate */
186           IPData.y[gidx] = vj[qj * dim + 1];
187           if (dim==3) IPData.z[gidx] = vj[qj * dim + 2];
188           IPData.w[gidx] = detJj[qj] * quadWeights[qj];
189           if (dim==2) IPData.w[gidx] *= IPData.x[gidx];  /* cylindrical coordinate, w/o 2pi */
190         } /* q */
191         ierr = DMPlexVecRestoreClosure(plex, section, locX, cStart+ej, NULL, &coef);CHKERRQ(ierr);
192       }
193     } /* ej */
194     if (shift==0.0) {
195       ierr = PetscLogEventEnd(ctx->events[7],0,0,0,0);CHKERRQ(ierr);
196     } else { // mass
197       ierr = PetscLogEventEnd(ctx->events[1],0,0,0,0);CHKERRQ(ierr);
198     }
199   }
200   ierr = DMRestoreLocalVector(plex, &locX);CHKERRQ(ierr);
201 
202   /* do it */
203   if (ctx->deviceType == LANDAU_CUDA || ctx->deviceType == LANDAU_KOKKOS) {
204     if (ctx->deviceType == LANDAU_CUDA) {
205 #if defined(PETSC_HAVE_CUDA)
206       ierr = LandauCUDAJacobian(plex,Nq,nu_alpha,nu_beta,invMass,Eq_m,&IPData,invJ_a,mass_w,shift,ctx->events,JacP);CHKERRQ(ierr);
207 #else
208       SETERRQ1(ctx->comm,PETSC_ERR_ARG_WRONG,"-landau_device_type %s not built","cuda");
209 #endif
210     } else if (ctx->deviceType == LANDAU_KOKKOS) {
211 #if defined(PETSC_HAVE_KOKKOS)
212       ierr = LandauKokkosJacobian(plex,Nq,nu_alpha,nu_beta,invMass,Eq_m,&IPData,invJ_a,ctx->subThreadBlockSize,mass_w,shift,ctx->events,JacP);CHKERRQ(ierr);
213 #else
214       SETERRQ1(ctx->comm,PETSC_ERR_ARG_WRONG,"-landau_device_type %s not built","kokkos");
215 #endif
216     }
217   } else { /* CPU version */
218     PetscInt                ei, qi;
219     PetscScalar             *elemMat;
220     PetscReal               *ff, *dudx, *dudy, *dudz;
221     const PetscReal * const BB = Tf[0]->T[0], * const DD = Tf[0]->T[1];
222     if (shift!=0.0) { // mass
223       ierr = PetscMalloc1(elemMatSize, &elemMat);CHKERRQ(ierr);
224     } else {
225       ierr = PetscLogEventBegin(ctx->events[8],0,0,0,0);CHKERRQ(ierr);
226       ierr = PetscMalloc5(elemMatSize, &elemMat, IPData.nip_*Nf, &ff,IPData.nip_*Nf, &dudx, IPData.nip_*Nf, &dudy, dim==3 ? IPData.nip_*Nf : 0, &dudz);CHKERRQ(ierr);
227       /* compute f and df */
228       for (ei = cStart, invJ = invJ_a; ei < cEnd; ++ei, invJ += Nq*dim*dim) {
229         LandauIPReal  *coef = &IPData.coefs[ei*Nb*Nf];
230         PetscScalar   u_x[LANDAU_MAX_SPECIES][LANDAU_DIM];
231         /* get f and df */
232         for (qi = 0; qi < Nq; ++qi) {
233           const PetscReal  *Bq = &BB[qi*Nb];
234           const PetscReal  *Dq = &DD[qi*Nb*dim];
235           const PetscInt   gidx = ei*Nq + qi;
236           /* get f & df */
237           for (f = 0; f < Nf; ++f) {
238             PetscInt    b, e;
239             PetscScalar refSpaceDer[LANDAU_DIM];
240             ff[gidx + f*IPData.nip_] = 0.0;
241             for (d = 0; d < LANDAU_DIM; ++d) refSpaceDer[d] = 0.0;
242             for (b = 0; b < Nb; ++b) {
243               const PetscInt    cidx = b;
244               ff[gidx + f*IPData.nip_] += Bq[cidx]*coef[f*Nb+cidx];
245               for (d = 0; d < dim; ++d) refSpaceDer[d] += Dq[cidx*dim+d]*coef[f*Nb+cidx];
246             }
247             for (d = 0; d < dim; ++d) {
248               for (e = 0, u_x[f][d] = 0.0; e < dim; ++e) {
249                 u_x[f][d] += invJ[qi * dim * dim + e*dim+d]*refSpaceDer[e];
250               }
251             }
252           }
253           for (f=0;f<Nf;f++) {
254             dudx[gidx + f*IPData.nip_] = PetscRealPart(u_x[f][0]);
255             dudy[gidx + f*IPData.nip_] = PetscRealPart(u_x[f][1]);
256 #if LANDAU_DIM==3
257             dudz[gidx + f*IPData.nip_] = PetscRealPart(u_x[f][2]);
258 #endif
259           }
260         }
261       }
262       ierr = PetscLogEventEnd(ctx->events[8],0,0,0,0);CHKERRQ(ierr);
263     }
264     for (ej = cStart; ej < cEnd; ++ej) {
265       ierr = PetscLogEventBegin(ctx->events[3],0,0,0,0);CHKERRQ(ierr);
266       ierr = PetscMemzero(elemMat, totDim *totDim * sizeof(PetscScalar));CHKERRQ(ierr);
267       ierr = PetscLogEventEnd(ctx->events[3],0,0,0,0);CHKERRQ(ierr);
268       ierr = PetscLogEventBegin(ctx->events[4],0,0,0,0);CHKERRQ(ierr);
269       ierr = PetscLogFlops((PetscLogDouble)Nq*flops);CHKERRQ(ierr);
270       invJ = invJ_a ? invJ_a + ej * Nq*dim*dim : NULL;
271       for (qj = 0; qj < Nq; ++qj) {
272         const PetscReal * const BB = Tf[0]->T[0], * const DD = Tf[0]->T[1];
273         PetscReal               g0[LANDAU_MAX_SPECIES], g2[LANDAU_MAX_SPECIES][LANDAU_DIM], g3[LANDAU_MAX_SPECIES][LANDAU_DIM][LANDAU_DIM];
274         PetscInt                d,d2,dp,d3,ipidx,fieldA;
275         const PetscInt          jpidx = Nq*(ej-cStart) + qj;
276         if (shift==0.0) {
277           const PetscReal * const invJj = &invJ[qj*dim*dim];
278           PetscReal               gg2[LANDAU_MAX_SPECIES][LANDAU_DIM],gg3[LANDAU_MAX_SPECIES][LANDAU_DIM][LANDAU_DIM], gg2_temp[LANDAU_DIM], gg3_temp[LANDAU_DIM][LANDAU_DIM];
279           const PetscReal         vj[3] = {IPData.x[jpidx], IPData.y[jpidx], IPData.z ? IPData.z[jpidx] : 0}, wj = IPData.w[jpidx];
280 
281           // create g2 & g3
282           for (d=0;d<dim;d++) { // clear accumulation data D & K
283             gg2_temp[d] = 0;
284             for (d2=0;d2<dim;d2++) gg3_temp[d][d2] = 0;
285           }
286           for (ipidx = 0; ipidx < IPData.nip_; ipidx++) {
287             const PetscReal wi = IPData.w[ipidx], x = IPData.x[ipidx], y = IPData.y[ipidx];
288             PetscReal       temp1[3] = {0, 0, 0}, temp2 = 0;
289 #if LANDAU_DIM==2
290             PetscReal       Ud[2][2], Uk[2][2];
291             LandauTensor2D(vj, x, y, Ud, Uk, (ipidx==jpidx) ? 0. : 1.);
292 #else
293             PetscReal U[3][3], z = IPData.z[ipidx];
294             LandauTensor3D(vj, x, y, z, U, (ipidx==jpidx) ? 0. : 1.);
295 #endif
296             for (fieldA = 0; fieldA < Nf; ++fieldA) {
297               temp1[0] += dudx[ipidx + fieldA*IPData.nip_]*nu_beta[fieldA]*invMass[fieldA];
298               temp1[1] += dudy[ipidx + fieldA*IPData.nip_]*nu_beta[fieldA]*invMass[fieldA];
299 #if LANDAU_DIM==3
300               temp1[2] += dudz[ipidx + fieldA*IPData.nip_]*nu_beta[fieldA]*invMass[fieldA];
301 #endif
302               temp2    += ff[ipidx + fieldA*IPData.nip_]*nu_beta[fieldA];
303             }
304             temp1[0] *= wi;
305             temp1[1] *= wi;
306 #if LANDAU_DIM==3
307             temp1[2] *= wi;
308 #endif
309             temp2    *= wi;
310 #if LANDAU_DIM==2
311             for (d2 = 0; d2 < 2; d2++) {
312               for (d3 = 0; d3 < 2; ++d3) {
313                 /* K = U * grad(f): g2=e: i,A */
314                 gg2_temp[d2] += Uk[d2][d3]*temp1[d3];
315                 /* D = -U * (I \kron (fx)): g3=f: i,j,A */
316                 gg3_temp[d2][d3] += Ud[d2][d3]*temp2;
317               }
318             }
319 #else
320             for (d2 = 0; d2 < 3; ++d2) {
321               for (d3 = 0; d3 < 3; ++d3) {
322                 /* K = U * grad(f): g2 = e: i,A */
323                 gg2_temp[d2] += U[d2][d3]*temp1[d3];
324                 /* D = -U * (I \kron (fx)): g3 = f: i,j,A */
325                 gg3_temp[d2][d3] += U[d2][d3]*temp2;
326               }
327             }
328 #endif
329           } /* IPs */
330           //if (ej==0) printf("\t:%d.%d) temp gg3=%e %e %e %e\n",ej,qj,gg3_temp[0][0],gg3_temp[1][0],gg3_temp[0][1],gg3_temp[1][1]);
331           // add alpha and put in gg2/3
332           for (fieldA = 0; fieldA < Nf; ++fieldA) {
333             for (d2 = 0; d2 < dim; d2++) {
334               gg2[fieldA][d2] = gg2_temp[d2]*nu_alpha[fieldA];
335               for (d3 = 0; d3 < dim; d3++) {
336                 gg3[fieldA][d2][d3] = -gg3_temp[d2][d3]*nu_alpha[fieldA]*invMass[fieldA];
337               }
338             }
339           }
340           /* add electric field term once per IP */
341           for (fieldA = 0; fieldA < Nf; ++fieldA) {
342             gg2[fieldA][dim-1] += Eq_m[fieldA];
343           }
344           /* Jacobian transform - g2, g3 */
345           for (fieldA = 0; fieldA < Nf; ++fieldA) {
346             for (d = 0; d < dim; ++d) {
347               g2[fieldA][d] = 0.0;
348               for (d2 = 0; d2 < dim; ++d2) {
349                 g2[fieldA][d] += invJj[d*dim+d2]*gg2[fieldA][d2];
350                 g3[fieldA][d][d2] = 0.0;
351                 for (d3 = 0; d3 < dim; ++d3) {
352                   for (dp = 0; dp < dim; ++dp) {
353                     g3[fieldA][d][d2] += invJj[d*dim + d3]*gg3[fieldA][d3][dp]*invJj[d2*dim + dp];
354                   }
355                 }
356                 g3[fieldA][d][d2] *= wj;
357               }
358               g2[fieldA][d] *= wj;
359             }
360           }
361         } else { // mass
362           /* Jacobian transform - g0 */
363           for (fieldA = 0; fieldA < Nf; ++fieldA) {
364             g0[fieldA] = mass_w[jpidx] * shift; // move this to below and remove g0
365           }
366         }
367         /* FE matrix construction */
368         {
369           PetscInt  fieldA,d,f,d2,g;
370           const PetscReal *BJq = &BB[qj*Nb], *DIq = &DD[qj*Nb*dim];
371           /* assemble - on the diagonal (I,I) */
372           for (fieldA = 0; fieldA < Nf ; fieldA++) {
373             for (f = 0; f < Nb ; f++) {
374               const PetscInt i = fieldA*Nb + f; /* Element matrix row */
375               for (g = 0; g < Nb; ++g) {
376                 const PetscInt j    = fieldA*Nb + g; /* Element matrix column */
377                 const PetscInt fOff = i*totDim + j;
378                 if (shift==0.0) {
379                   for (d = 0; d < dim; ++d) {
380                     elemMat[fOff] += DIq[f*dim+d]*g2[fieldA][d]*BJq[g];
381                     for (d2 = 0; d2 < dim; ++d2) {
382                       elemMat[fOff] += DIq[f*dim + d]*g3[fieldA][d][d2]*DIq[g*dim + d2];
383                     }
384                   }
385                 } else { // mass
386                   elemMat[fOff] += BJq[f]*g0[fieldA]*BJq[g];
387                 }
388               }
389             }
390           }
391         }
392       } /* qj loop */
393       ierr = PetscLogEventEnd(ctx->events[4],0,0,0,0);CHKERRQ(ierr);
394       /* assemble matrix */
395       ierr = PetscLogEventBegin(ctx->events[6],0,0,0,0);CHKERRQ(ierr);
396       if (!maps) {
397         ierr = DMPlexMatSetClosure(plex, section, globsection, JacP, ej, elemMat, ADD_VALUES);CHKERRQ(ierr);
398       } else {  // GPU like assembly for debugging
399         PetscInt      fieldA,idx,q,f,g,d,nr,nc,rows0[LANDAU_MAX_Q_FACE],cols0[LANDAU_MAX_Q_FACE],rows[LANDAU_MAX_Q_FACE],cols[LANDAU_MAX_Q_FACE];
400         PetscScalar   vals[LANDAU_MAX_Q_FACE*LANDAU_MAX_Q_FACE],row_scale[LANDAU_MAX_Q_FACE],col_scale[LANDAU_MAX_Q_FACE];
401         for (g=0;g<LANDAU_MAX_Q_FACE;g++) { col_scale[g]=0.; cols0[g]=0.; }
402         /* assemble - from the diagonal (I,I) in this format for DMPlexMatSetClosure */
403         for (fieldA = 0; fieldA < Nf ; fieldA++) {
404           LandauIdx *const Idxs = &maps->gIdx[ej-cStart][fieldA][0];
405           for (f = 0; f < Nb ; f++) {
406             idx = Idxs[f];
407             if (idx >= 0) {
408               nr = 1;
409               rows0[0] = idx;
410               row_scale[0] = 1.;
411             } else {
412               idx = -idx - 1;
413               nr = maps->num_face;
414               for (q = 0; q < maps->num_face; q++) {
415                 rows0[q]     = maps->c_maps[idx][q].gid;
416                 row_scale[q] = maps->c_maps[idx][q].scale;
417               }
418             }
419             for (g = 0; g < Nb; ++g) {
420               idx = Idxs[g];
421               if (idx >= 0) {
422                 nc = 1;
423                 cols0[0] = idx;
424                 col_scale[0] = 1.;
425               } else {
426                 idx = -idx - 1;
427                 nc = maps->num_face;
428                 for (q = 0; q < maps->num_face; q++) {
429                   cols0[q]     = maps->c_maps[idx][q].gid;
430                   col_scale[q] = maps->c_maps[idx][q].scale;
431                 }
432               }
433               const PetscInt    i = fieldA*Nb + f; /* Element matrix row */
434               const PetscInt    j = fieldA*Nb + g; /* Element matrix column */
435               const PetscScalar Aij = elemMat[i*totDim + j];
436               for (q = 0; q < nr; q++) rows[q] = rows0[q];
437               for (q = 0; q < nc; q++) cols[q] = cols0[q];
438               for (q = 0; q < nr; q++) {
439                 for (d = 0; d < nc; d++) {
440                   vals[q*nc + d] = row_scale[q]*col_scale[d]*Aij;
441                 }
442               }
443               ierr = MatSetValues(JacP,nr,rows,nc,cols,vals,ADD_VALUES);CHKERRQ(ierr);
444             }
445           }
446         }
447       }
448       if (ej==-3) {
449         PetscErrorCode    ierr2;
450         ierr2 = PetscPrintf(ctx->comm,"CPU Element matrix\n");CHKERRQ(ierr2);
451         for (d = 0; d < totDim; ++d){
452           for (f = 0; f < totDim; ++f) {ierr2 = PetscPrintf(ctx->comm," %12.5e",  PetscRealPart(elemMat[d*totDim + f]));CHKERRQ(ierr2);}
453           ierr2 = PetscPrintf(ctx->comm,"\n");CHKERRQ(ierr2);
454         }
455         exit(12);
456       }
457       ierr = PetscLogEventEnd(ctx->events[6],0,0,0,0);CHKERRQ(ierr);
458     } /* ej cells loop, not cuda */
459     if (shift!=0.0) { // mass
460       ierr = PetscFree(elemMat);CHKERRQ(ierr);
461     } else {
462       ierr = PetscFree5(elemMat, ff, dudx, dudy, dudz);CHKERRQ(ierr);
463     }
464   } /* CPU version */
465   /* assemble matrix or vector */
466   ierr = MatAssemblyBegin(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
467   ierr = MatAssemblyEnd(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
468 #define MAP_BF_SIZE (128*LANDAU_DIM*LANDAU_MAX_Q_FACE*LANDAU_MAX_SPECIES)
469   if (ctx->gpu_assembly && !container) {
470     PetscScalar             elemMatrix[LANDAU_MAX_NQ*LANDAU_MAX_NQ*LANDAU_MAX_SPECIES*LANDAU_MAX_SPECIES], *elMat;
471     pointInterpolationP4est pointMaps[MAP_BF_SIZE][LANDAU_MAX_Q_FACE];
472     PetscInt                q,eidx,fieldA;
473     MatType                 type;
474     ierr = PetscInfo1(JacP, "Make GPU maps %D\n",1);CHKERRQ(ierr);
475     ierr = MatGetType(JacP,&type);CHKERRQ(ierr);
476     ierr = PetscLogEventBegin(ctx->events[2],0,0,0,0);CHKERRQ(ierr);
477     ierr = PetscMalloc(sizeof(P4estVertexMaps), &maps);CHKERRQ(ierr);
478     ierr = PetscContainerCreate(PETSC_COMM_SELF, &container);CHKERRQ(ierr);
479     ierr = PetscContainerSetPointer(container, (void *)maps);CHKERRQ(ierr);
480     ierr = PetscContainerSetUserDestroy(container, LandauGPUDataDestroy);CHKERRQ(ierr);
481     ierr = PetscObjectCompose((PetscObject) JacP, "assembly_maps", (PetscObject) container);CHKERRQ(ierr);
482     ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
483     // make maps
484     maps->data = NULL;
485     maps->num_elements = numCells;
486     maps->num_face = (PetscInt)(pow(Nq,1./((double)dim))+.001); // Q
487     maps->num_face = (PetscInt)(pow(maps->num_face,(double)(dim-1))+.001); // Q^2
488     maps->num_reduced = 0;
489     maps->deviceType = ctx->deviceType;
490     // count reduced and get
491     ierr = PetscMalloc(maps->num_elements * sizeof *maps->gIdx, &maps->gIdx);CHKERRQ(ierr);
492     for (fieldA=0;fieldA<Nf;fieldA++) {
493       for (ej = cStart, eidx = 0 ; ej < cEnd; ++ej, ++eidx) {
494         for (q = 0; q < Nb; ++q) {
495           PetscInt    numindices,*indices;
496           PetscScalar *valuesOrig = elMat = elemMatrix;
497           ierr = PetscMemzero(elMat, totDim*totDim*sizeof(PetscScalar));CHKERRQ(ierr);
498           elMat[ (fieldA*Nb + q)*totDim + fieldA*Nb + q] = 1;
499           ierr = DMPlexGetClosureIndices(plex, section, globsection, ej, PETSC_TRUE, &numindices, &indices, NULL, (PetscScalar **) &elMat);CHKERRQ(ierr);
500           for (f = 0 ; f < numindices ; ++f) { // look for a non-zero on the diagonal
501             if (PetscAbs(PetscRealPart(elMat[f*numindices + f])) > PETSC_MACHINE_EPSILON) {
502               // found it
503               if (PetscAbs(PetscRealPart(elMat[f*numindices + f] - 1.)) < PETSC_MACHINE_EPSILON) {
504                 maps->gIdx[eidx][fieldA][q] = (LandauIdx)indices[f]; // normal vertex 1.0
505               } else { //found a constraint
506                 int       jj = 0;
507                 PetscReal sum = 0;
508                 const PetscInt ff = f;
509                 maps->gIdx[eidx][fieldA][q] = -maps->num_reduced - 1; // gid = -(idx+1): idx = -gid - 1
510                 do {  // constraints are continous in Plex - exploit that here
511                   int ii;
512                   for (ii = 0, pointMaps[maps->num_reduced][jj].scale = 0; ii < maps->num_face; ii++) { // DMPlex puts them all together
513                     if (ff + ii < numindices) {
514                       pointMaps[maps->num_reduced][jj].scale += PetscRealPart(elMat[f*numindices + ff + ii]);
515                     }
516                   }
517                   sum += pointMaps[maps->num_reduced][jj].scale;
518                   if (pointMaps[maps->num_reduced][jj].scale == 0) pointMaps[maps->num_reduced][jj].gid = -1; // 3D has Q and Q^2 interps -- all contiguous???
519                   else                                             pointMaps[maps->num_reduced][jj].gid = indices[f];
520                 } while (++jj < maps->num_face && ++f < numindices); // jj is incremented if we hit the end
521                 while (jj++ < maps->num_face) {
522                   pointMaps[maps->num_reduced][jj].scale = 0;
523                   pointMaps[maps->num_reduced][jj].gid = -1;
524                 }
525                 if (PetscAbs(sum-1.0)>PETSC_MACHINE_EPSILON*2.0) { // debug
526                   int       d,f;
527                   PetscReal tmp = 0;
528                   PetscPrintf(PETSC_COMM_SELF,"\t\t%D.%D.%D) ERROR total I = %22.16e (LANDAU_MAX_Q_FACE=%d, #face=%D)\n",eidx,q,fieldA,tmp,LANDAU_MAX_Q_FACE,maps->num_face);
529                   for (d = 0, tmp = 0; d < numindices; ++d){
530                     if (tmp!=0 && PetscAbs(tmp-1.0)>2*PETSC_MACHINE_EPSILON) ierr = PetscPrintf(PETSC_COMM_WORLD,"%3D) %3D: ",d,indices[d]);CHKERRQ(ierr);
531                     for (f = 0; f < numindices; ++f) {
532                       tmp += PetscRealPart(elMat[d*numindices + f]);
533                     }
534                     if (tmp!=0) ierr = PetscPrintf(ctx->comm," | %22.16e\n",tmp);CHKERRQ(ierr);
535                   }
536                 }
537                 maps->num_reduced++;
538                 if (maps->num_reduced>=MAP_BF_SIZE) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "maps->num_reduced %d > %d",maps->num_reduced,MAP_BF_SIZE);
539               }
540               break;
541             }
542           }
543           // cleanup
544           ierr = DMPlexRestoreClosureIndices(plex, section, globsection, ej, PETSC_TRUE, &numindices, &indices, NULL, (PetscScalar **) &elMat);CHKERRQ(ierr);
545           if (elMat != valuesOrig) {ierr = DMRestoreWorkArray(plex, numindices*numindices, MPIU_SCALAR, &elMat);}
546         }
547       }
548     }
549     // allocate and copy point datamaps->gIdx[eidx][field][q] -- for CPU version of this code, for debugging
550     ierr = PetscMalloc(maps->num_reduced * sizeof *maps->c_maps, &maps->c_maps);CHKERRQ(ierr);
551     for (ej = 0; ej < maps->num_reduced; ++ej) {
552       for (q = 0; q < maps->num_face; ++q) {
553         maps->c_maps[ej][q].scale = pointMaps[ej][q].scale;
554         maps->c_maps[ej][q].gid = pointMaps[ej][q].gid;
555       }
556     }
557 #if defined(PETSC_HAVE_KOKKOS_KERNELS)
558     if (ctx->deviceType == LANDAU_KOKKOS) {
559       ierr = LandauKokkosCreateMatMaps(maps, pointMaps,Nf,Nq);CHKERRQ(ierr); // imples Kokkos does
560     } // else could be CUDA
561 #endif
562 #if defined(PETSC_HAVE_CUDA)
563     if (ctx->deviceType == LANDAU_CUDA){
564       ierr = LandauCUDACreateMatMaps(maps, pointMaps,Nf,Nq);CHKERRQ(ierr);
565     }
566 #endif
567     ierr = PetscLogEventEnd(ctx->events[2],0,0,0,0);CHKERRQ(ierr);
568   }
569   /* clean up */
570   ierr = DMDestroy(&plex);CHKERRQ(ierr);
571   if (shift==0.0) {
572     ierr = LandauPointDataDestroy(&IPData);CHKERRQ(ierr);
573     ierr = PetscFree(invJ_a);CHKERRQ(ierr);
574   } else {
575     ierr = PetscFree(mass_w);CHKERRQ(ierr);
576   }
577   PetscFunctionReturn(0);
578 }
579 
580 #if defined(LANDAU_ADD_BCS)
581 static void zero_bc(PetscInt dim, PetscInt Nf, PetscInt NfAux,
582                     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
583                     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
584                     PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar uexact[])
585 {
586   uexact[0] = 0;
587 }
588 #endif
589 
590 #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]; }}
591 static void CircleInflate(PetscReal r1, PetscReal r2, PetscReal r0, PetscInt num_sections, PetscReal x, PetscReal y,
592                           PetscReal *outX, PetscReal *outY)
593 {
594   PetscReal rr = PetscSqrtReal(x*x + y*y), outfact, efact;
595   if (rr < r1 + PETSC_SQRT_MACHINE_EPSILON) {
596     *outX = x; *outY = y;
597   } else {
598     const PetscReal xy[2] = {x,y}, sinphi=y/rr, cosphi=x/rr;
599     PetscReal       cth,sth,xyprime[2],Rth[2][2],rotcos,newrr;
600     if (num_sections==2) {
601       rotcos = 0.70710678118654;
602       outfact = 1.5; efact = 2.5;
603       /* rotate normalized vector into [-pi/4,pi/4) */
604       if (sinphi >= 0.) {         /* top cell, -pi/2 */
605         cth = 0.707106781186548; sth = -0.707106781186548;
606       } else {                    /* bottom cell -pi/8 */
607         cth = 0.707106781186548; sth = .707106781186548;
608       }
609     } else if (num_sections==3) {
610       rotcos = 0.86602540378443;
611       outfact = 1.5; efact = 2.5;
612       /* rotate normalized vector into [-pi/6,pi/6) */
613       if (sinphi >= 0.5) {         /* top cell, -pi/3 */
614         cth = 0.5; sth = -0.866025403784439;
615       } else if (sinphi >= -.5) {  /* mid cell 0 */
616         cth = 1.; sth = .0;
617       } else { /* bottom cell +pi/3 */
618         cth = 0.5; sth = 0.866025403784439;
619       }
620     } else if (num_sections==4) {
621       rotcos = 0.9238795325112;
622       outfact = 1.5; efact = 3;
623       /* rotate normalized vector into [-pi/8,pi/8) */
624       if (sinphi >= 0.707106781186548) {         /* top cell, -3pi/8 */
625         cth = 0.38268343236509; sth = -0.923879532511287;
626       } else if (sinphi >= 0.) {                 /* mid top cell -pi/8 */
627         cth = 0.923879532511287; sth = -.38268343236509;
628       } else if (sinphi >= -0.707106781186548) { /* mid bottom cell + pi/8 */
629         cth = 0.923879532511287; sth = 0.38268343236509;
630       } else {                                   /* bottom cell + 3pi/8 */
631         cth = 0.38268343236509; sth = .923879532511287;
632       }
633     } else {
634       cth = 0.; sth = 0.; rotcos = 0; efact = 0;
635     }
636     Rth[0][0] = cth; Rth[0][1] =-sth;
637     Rth[1][0] = sth; Rth[1][1] = cth;
638     MATVEC2(Rth,xy,xyprime);
639     if (num_sections==2) {
640       newrr = xyprime[0]/rotcos;
641     } else {
642       PetscReal newcosphi=xyprime[0]/rr, rin = r1, rout = rr - rin;
643       PetscReal routmax = r0*rotcos/newcosphi - rin, nroutmax = r0 - rin, routfrac = rout/routmax;
644       newrr = rin + routfrac*nroutmax;
645     }
646     *outX = cosphi*newrr; *outY = sinphi*newrr;
647     /* grade */
648     PetscReal fact,tt,rs,re, rr = PetscSqrtReal(PetscSqr(*outX) + PetscSqr(*outY));
649     if (rr > r2) { rs = r2; re = r0; fact = outfact;} /* outer zone */
650     else {         rs = r1; re = r2; fact = efact;} /* electron zone */
651     tt = (rs + PetscPowReal((rr - rs)/(re - rs),fact) * (re-rs)) / rr;
652     *outX *= tt;
653     *outY *= tt;
654   }
655 }
656 
657 static PetscErrorCode GeometryDMLandau(DM base, PetscInt point, PetscInt dim, const PetscReal abc[], PetscReal xyz[], void *a_ctx)
658 {
659   LandauCtx   *ctx = (LandauCtx*)a_ctx;
660   PetscReal   r = abc[0], z = abc[1];
661   if (ctx->inflate) {
662     PetscReal absR, absZ;
663     absR = PetscAbs(r);
664     absZ = PetscAbs(z);
665     CircleInflate(ctx->i_radius,ctx->e_radius,ctx->radius,ctx->num_sections,absR,absZ,&absR,&absZ);
666     r = (r > 0) ? absR : -absR;
667     z = (z > 0) ? absZ : -absZ;
668   }
669   xyz[0] = r;
670   xyz[1] = z;
671   if (dim==3) xyz[2] = abc[2];
672 
673   PetscFunctionReturn(0);
674 }
675 
676 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)
677 {
678   PetscReal err = 0.0;
679   PetscInt  f = *(PetscInt*)actx, j;
680   PetscFunctionBegin;
681   for (j = 0; j < dim; ++j) {
682     err += PetscSqr(PetscRealPart(u_x[f*dim+j]));
683   }
684   err = PetscRealPart(u[f]); /* just use rho */
685   *error = volume * err; /* * (ctx->axisymmetric ? 2.*PETSC_PI * r : 1); */
686   PetscFunctionReturn(0);
687 }
688 
689 static PetscErrorCode LandauDMCreateVMesh(MPI_Comm comm, const PetscInt dim, const char prefix[], LandauCtx *ctx, DM *dm)
690 {
691   PetscErrorCode ierr;
692   PetscReal      radius = ctx->radius;
693   size_t         len;
694   char           fname[128] = ""; /* we can add a file if we want */
695 
696   PetscFunctionBegin;
697   /* create DM */
698   ierr = PetscStrlen(fname, &len);CHKERRQ(ierr);
699   if (len) {
700     PetscInt dim2;
701     ierr = DMPlexCreateFromFile(comm, fname, ctx->interpolate, dm);CHKERRQ(ierr);
702     ierr = DMGetDimension(*dm, &dim2);CHKERRQ(ierr);
703     if (LANDAU_DIM != dim2) SETERRQ2(comm, PETSC_ERR_PLIB, "dim %D != LANDAU_DIM %d",dim2,LANDAU_DIM);
704   } else {    /* p4est, quads */
705     /* Create plex mesh of Landau domain */
706     if (!ctx->sphere) {
707       PetscInt       cells[] = {2,2,2};
708       PetscReal      lo[] = {-radius,-radius,-radius}, hi[] = {radius,radius,radius};
709       DMBoundaryType periodicity[3] = {DM_BOUNDARY_NONE, dim==2 ? DM_BOUNDARY_NONE : DM_BOUNDARY_NONE, DM_BOUNDARY_NONE};
710       if (dim==2) { lo[0] = 0; cells[0] = 1; }
711       ierr = DMPlexCreateBoxMesh(comm, dim, PETSC_FALSE, cells, lo, hi, periodicity, PETSC_TRUE, dm);CHKERRQ(ierr);
712       ierr = DMLocalizeCoordinates(*dm);CHKERRQ(ierr); /* needed for periodic */
713       if (dim==3) {ierr = PetscObjectSetName((PetscObject) *dm, "cube");CHKERRQ(ierr);}
714       else {ierr = PetscObjectSetName((PetscObject) *dm, "half-plane");CHKERRQ(ierr);}
715     } else if (dim==2) {
716       PetscInt       numCells,cells[16][4],i,j;
717       PetscInt       numVerts;
718       PetscReal      inner_radius1 = ctx->i_radius, inner_radius2 = ctx->e_radius;
719       PetscReal      *flatCoords = NULL;
720       PetscInt       *flatCells = NULL, *pcell;
721       if (ctx->num_sections==2) {
722 #if 1
723         numCells = 5;
724         numVerts = 10;
725         int cells2[][4] = { {0,1,4,3},
726                             {1,2,5,4},
727                             {3,4,7,6},
728                             {4,5,8,7},
729                             {6,7,8,9} };
730         for (i = 0; i < numCells; i++) for (j = 0; j < 4; j++) cells[i][j] = cells2[i][j];
731         ierr = PetscMalloc2(numVerts * 2, &flatCoords, numCells * 4, &flatCells);CHKERRQ(ierr);
732         {
733           PetscReal (*coords)[2] = (PetscReal (*) [2]) flatCoords;
734           for (j = 0; j < numVerts-1; j++) {
735             PetscReal z, r, theta = -PETSC_PI/2 + (j%3) * PETSC_PI/2;
736             PetscReal rad = (j >= 6) ? inner_radius1 : (j >= 3) ? inner_radius2 : ctx->radius;
737             z = rad * PetscSinReal(theta);
738             coords[j][1] = z;
739             r = rad * PetscCosReal(theta);
740             coords[j][0] = r;
741           }
742           coords[numVerts-1][0] = coords[numVerts-1][1] = 0;
743         }
744 #else
745         numCells = 4;
746         numVerts = 8;
747         static int     cells2[][4] = {{0,1,2,3},
748                                       {4,5,1,0},
749                                       {5,6,2,1},
750                                       {6,7,3,2}};
751         for (i = 0; i < numCells; i++) for (j = 0; j < 4; j++) cells[i][j] = cells2[i][j];
752         ierr = loc2(numVerts * 2, &flatCoords, numCells * 4, &flatCells);CHKERRQ(ierr);
753         {
754           PetscReal (*coords)[2] = (PetscReal (*) [2]) flatCoords;
755           PetscInt j;
756           for (j = 0; j < 8; j++) {
757             PetscReal z, r;
758             PetscReal theta = -PETSC_PI/2 + (j%4) * PETSC_PI/3.;
759             PetscReal rad = ctx->radius * ((j < 4) ? 0.5 : 1.0);
760             z = rad * PetscSinReal(theta);
761             coords[j][1] = z;
762             r = rad * PetscCosReal(theta);
763             coords[j][0] = r;
764           }
765         }
766 #endif
767       } else if (ctx->num_sections==3) {
768         numCells = 7;
769         numVerts = 12;
770         int cells2[][4] = { {0,1,5,4},
771                             {1,2,6,5},
772                             {2,3,7,6},
773                             {4,5,9,8},
774                             {5,6,10,9},
775                             {6,7,11,10},
776                             {8,9,10,11} };
777         for (i = 0; i < numCells; i++) for (j = 0; j < 4; j++) cells[i][j] = cells2[i][j];
778         ierr = PetscMalloc2(numVerts * 2, &flatCoords, numCells * 4, &flatCells);CHKERRQ(ierr);
779         {
780           PetscReal (*coords)[2] = (PetscReal (*) [2]) flatCoords;
781           for (j = 0; j < numVerts; j++) {
782             PetscReal z, r, theta = -PETSC_PI/2 + (j%4) * PETSC_PI/3;
783             PetscReal rad = (j >= 8) ? inner_radius1 : (j >= 4) ? inner_radius2 : ctx->radius;
784             z = rad * PetscSinReal(theta);
785             coords[j][1] = z;
786             r = rad * PetscCosReal(theta);
787             coords[j][0] = r;
788           }
789         }
790       } else if (ctx->num_sections==4) {
791         numCells = 10;
792         numVerts = 16;
793         int cells2[][4] = { {0,1,6,5},
794                             {1,2,7,6},
795                             {2,3,8,7},
796                             {3,4,9,8},
797                             {5,6,11,10},
798                             {6,7,12,11},
799                             {7,8,13,12},
800                             {8,9,14,13},
801                             {10,11,12,15},
802                             {12,13,14,15}};
803         for (i = 0; i < numCells; i++) for (j = 0; j < 4; j++) cells[i][j] = cells2[i][j];
804         ierr = PetscMalloc2(numVerts * 2, &flatCoords, numCells * 4, &flatCells);CHKERRQ(ierr);
805         {
806           PetscReal (*coords)[2] = (PetscReal (*) [2]) flatCoords;
807           for (j = 0; j < numVerts-1; j++) {
808             PetscReal z, r, theta = -PETSC_PI/2 + (j%5) * PETSC_PI/4;
809             PetscReal rad = (j >= 10) ? inner_radius1 : (j >= 5) ? inner_radius2 : ctx->radius;
810             z = rad * PetscSinReal(theta);
811             coords[j][1] = z;
812             r = rad * PetscCosReal(theta);
813             coords[j][0] = r;
814           }
815           coords[numVerts-1][0] = coords[numVerts-1][1] = 0;
816         }
817       } else {
818         numCells = 0;
819         numVerts = 0;
820       }
821       for (j = 0, pcell = flatCells; j < numCells; j++, pcell += 4) {
822         pcell[0] = cells[j][0]; pcell[1] = cells[j][1];
823         pcell[2] = cells[j][2]; pcell[3] = cells[j][3];
824       }
825       ierr = DMPlexCreateFromCellListPetsc(comm,2,numCells,numVerts,4,ctx->interpolate,flatCells,2,flatCoords,dm);CHKERRQ(ierr);
826       ierr = PetscFree2(flatCoords,flatCells);CHKERRQ(ierr);
827       ierr = PetscObjectSetName((PetscObject) *dm, "semi-circle");CHKERRQ(ierr);
828     } else SETERRQ(ctx->comm, PETSC_ERR_PLIB, "Velocity space meshes does not support cubed sphere");
829   }
830   ierr = PetscObjectSetOptionsPrefix((PetscObject)*dm,prefix);CHKERRQ(ierr);
831 
832   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); /* Plex refine */
833 
834   { /* p4est? */
835     char      convType[256];
836     PetscBool flg;
837     ierr = PetscOptionsBegin(ctx->comm, prefix, "Mesh conversion options", "DMPLEX");CHKERRQ(ierr);
838     ierr = PetscOptionsFList("-dm_landau_type","Convert DMPlex to another format (should not be Plex!)","plexland.c",DMList,DMPLEX,convType,256,&flg);CHKERRQ(ierr);
839     ierr = PetscOptionsEnd();
840     if (flg) {
841       DM dmforest;
842       ierr = DMConvert(*dm,convType,&dmforest);CHKERRQ(ierr);
843       if (dmforest) {
844         PetscBool isForest;
845         if (dmforest->prealloc_only != (*dm)->prealloc_only) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"plex->prealloc_only != dm->prealloc_only");
846         ierr = PetscObjectSetOptionsPrefix((PetscObject)dmforest,prefix);CHKERRQ(ierr);
847         ierr = DMIsForest(dmforest,&isForest);CHKERRQ(ierr);
848         if (isForest) {
849           if (ctx->sphere && ctx->inflate) {
850             ierr = DMForestSetBaseCoordinateMapping(dmforest,GeometryDMLandau,ctx);CHKERRQ(ierr);
851           }
852           if (dmforest->prealloc_only != (*dm)->prealloc_only) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"plex->prealloc_only != dm->prealloc_only");
853           ierr = DMDestroy(dm);CHKERRQ(ierr);
854           *dm = dmforest;
855           ctx->errorIndicator = ErrorIndicator_Simple; /* flag for Forest */
856         } else SETERRQ(ctx->comm, PETSC_ERR_USER, "Converted to non Forest?");
857       } else SETERRQ(ctx->comm, PETSC_ERR_USER, "Convert failed?");
858     }
859   }
860   ierr = PetscObjectSetName((PetscObject) *dm, "Mesh");CHKERRQ(ierr);
861   PetscFunctionReturn(0);
862 }
863 
864 static PetscErrorCode SetupDS(DM dm, PetscInt dim, LandauCtx *ctx)
865 {
866   PetscErrorCode  ierr;
867   PetscInt        ii;
868   PetscFunctionBegin;
869   for (ii=0;ii<ctx->num_species;ii++) {
870     char     buf[256];
871     if (ii==0) ierr = PetscSNPrintf(buf, 256, "e");
872     else {ierr = PetscSNPrintf(buf, 256, "i%D", ii);CHKERRQ(ierr);}
873     /* Setup Discretization - FEM */
874     ierr = PetscFECreateDefault(PetscObjectComm((PetscObject) dm), dim, 1, PETSC_FALSE, NULL, PETSC_DECIDE, &ctx->fe[ii]);CHKERRQ(ierr);
875     ierr = PetscObjectSetName((PetscObject) ctx->fe[ii], buf);CHKERRQ(ierr);
876     ierr = DMSetField(dm, ii, NULL, (PetscObject) ctx->fe[ii]);CHKERRQ(ierr);
877   }
878   ierr = DMCreateDS(dm);CHKERRQ(ierr);
879   if (1) {
880     PetscInt        ii;
881     PetscSection    section;
882     ierr = DMGetSection(dm, &section);CHKERRQ(ierr);
883     for (ii=0;ii<ctx->num_species;ii++){
884       char buf[256];
885       if (ii==0) ierr = PetscSNPrintf(buf, 256, "se");
886       else ierr = PetscSNPrintf(buf, 256, "si%D", ii);
887       ierr = PetscSectionSetComponentName(section, ii, 0, buf);CHKERRQ(ierr);
888     }
889   }
890   PetscFunctionReturn(0);
891 }
892 
893 /* Define a Maxwellian function for testing out the operator. */
894 
895 /* Using cartesian velocity space coordinates, the particle */
896 /* density, [1/m^3], is defined according to */
897 
898 /* $$ n=\int_{R^3} dv^3 \left(\frac{m}{2\pi T}\right)^{3/2}\exp [- mv^2/(2T)] $$ */
899 
900 /* Using some constant, c, we normalize the velocity vector into a */
901 /* dimensionless variable according to v=c*x. Thus the density, $n$, becomes */
902 
903 /* $$ n=\int_{R^3} dx^3 \left(\frac{mc^2}{2\pi T}\right)^{3/2}\exp [- mc^2/(2T)*x^2] $$ */
904 
905 /* Defining $\theta=2T/mc^2$, we thus find that the probability density */
906 /* for finding the particle within the interval in a box dx^3 around x is */
907 
908 /* f(x;\theta)=\left(\frac{1}{\pi\theta}\right)^{3/2} \exp [ -x^2/\theta ] */
909 
910 typedef struct {
911   LandauCtx   *ctx;
912   PetscReal kT_m;
913   PetscReal n;
914   PetscReal shift;
915 } MaxwellianCtx;
916 
917 static PetscErrorCode maxwellian(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf_dummy, PetscScalar *u, void *actx)
918 {
919   MaxwellianCtx *mctx = (MaxwellianCtx*)actx;
920   LandauCtx     *ctx = mctx->ctx;
921   PetscInt      i;
922   PetscReal     v2 = 0, theta = 2*mctx->kT_m/(ctx->v_0*ctx->v_0); /* theta = 2kT/mc^2 */
923   PetscFunctionBegin;
924   /* compute the exponents, v^2 */
925   for (i = 0; i < dim; ++i) v2 += x[i]*x[i];
926   /* evaluate the Maxwellian */
927   u[0] = mctx->n*PetscPowReal(PETSC_PI*theta,-1.5)*(PetscExpReal(-v2/theta));
928   if (mctx->shift!=0.) {
929     v2 = 0;
930     for (i = 0; i < dim-1; ++i) v2 += x[i]*x[i];
931     v2 += (x[dim-1]-mctx->shift)*(x[dim-1]-mctx->shift);
932     /* evaluate the shifted Maxwellian */
933     u[0] += mctx->n*PetscPowReal(PETSC_PI*theta,-1.5)*(PetscExpReal(-v2/theta));
934   }
935   PetscFunctionReturn(0);
936 }
937 
938 /*@
939  LandauAddMaxwellians - Add a Maxwellian distribution to a state
940 
941  Collective on X
942 
943  Input Parameters:
944  .   dm - The mesh
945  +   time - Current time
946  -   temps - Temperatures of each species
947  .   ns - Number density of each species
948  +   actx - Landau context
949 
950  Output Parameter:
951  .   X  - The state
952 
953  Level: beginner
954 
955  .keywords: mesh
956  .seealso: LandauCreateVelocitySpace()
957  @*/
958 PetscErrorCode LandauAddMaxwellians(DM dm, Vec X, PetscReal time, PetscReal temps[], PetscReal ns[], void *actx)
959 {
960   LandauCtx      *ctx = (LandauCtx*)actx;
961   PetscErrorCode (*initu[LANDAU_MAX_SPECIES])(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar [], void *);
962   PetscErrorCode ierr,ii;
963   PetscInt       dim;
964   MaxwellianCtx  *mctxs[LANDAU_MAX_SPECIES], data[LANDAU_MAX_SPECIES];
965 
966   PetscFunctionBegin;
967   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
968   if (!ctx) { ierr = DMGetApplicationContext(dm, &ctx);CHKERRQ(ierr); }
969   for (ii=0;ii<ctx->num_species;ii++) {
970     mctxs[ii] = &data[ii];
971     data[ii].ctx = ctx;
972     data[ii].kT_m = ctx->k*temps[ii]/ctx->masses[ii]; /* kT/m */
973     data[ii].n = ns[ii];
974     initu[ii] = maxwellian;
975     data[ii].shift = 0;
976   }
977   data[0].shift = ctx->electronShift;
978   /* need to make ADD_ALL_VALUES work - TODO */
979   ierr = DMProjectFunction(dm, time, initu, (void**)mctxs, INSERT_ALL_VALUES, X);CHKERRQ(ierr);
980   PetscFunctionReturn(0);
981 }
982 
983 /*
984  LandauSetInitialCondition - Addes Maxwellians with context
985 
986  Collective on X
987 
988  Input Parameters:
989  .   dm - The mesh
990  +   actx - Landau context with T and n
991 
992  Output Parameter:
993  .   X  - The state
994 
995  Level: beginner
996 
997  .keywords: mesh
998  .seealso: LandauCreateVelocitySpace(), LandauAddMaxwellians()
999  */
1000 static PetscErrorCode LandauSetInitialCondition(DM dm, Vec X, void *actx)
1001 {
1002   LandauCtx        *ctx = (LandauCtx*)actx;
1003   PetscErrorCode ierr;
1004   PetscFunctionBegin;
1005   if (!ctx) { ierr = DMGetApplicationContext(dm, &ctx);CHKERRQ(ierr); }
1006   ierr = VecZeroEntries(X);CHKERRQ(ierr);
1007   ierr = LandauAddMaxwellians(dm, X, 0.0, ctx->thermal_temps, ctx->n, ctx);CHKERRQ(ierr);
1008   PetscFunctionReturn(0);
1009 }
1010 
1011 static PetscErrorCode adaptToleranceFEM(PetscFE fem, Vec sol, PetscReal refineTol[], PetscReal coarsenTol[], PetscInt type, LandauCtx *ctx, DM *newDM)
1012 {
1013   DM               dm, plex, adaptedDM = NULL;
1014   PetscDS          prob;
1015   PetscBool        isForest;
1016   PetscQuadrature  quad;
1017   PetscInt         Nq, *Nb, cStart, cEnd, c, dim, qj, k;
1018   DMLabel          adaptLabel = NULL;
1019   PetscErrorCode   ierr;
1020 
1021   PetscFunctionBegin;
1022   ierr = VecGetDM(sol, &dm);CHKERRQ(ierr);
1023   ierr = DMCreateDS(dm);CHKERRQ(ierr);
1024   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
1025   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1026   ierr = DMIsForest(dm, &isForest);CHKERRQ(ierr);
1027   ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr);
1028   ierr = DMPlexGetHeightStratum(plex,0,&cStart,&cEnd);CHKERRQ(ierr);
1029   ierr = DMLabelCreate(PETSC_COMM_SELF,"adapt",&adaptLabel);CHKERRQ(ierr);
1030   ierr = PetscFEGetQuadrature(fem, &quad);CHKERRQ(ierr);
1031   ierr = PetscQuadratureGetData(quad, NULL, NULL, &Nq, NULL, NULL);CHKERRQ(ierr);
1032   if (Nq >LANDAU_MAX_NQ) SETERRQ2(ctx->comm,PETSC_ERR_ARG_WRONG,"Order too high. Nq = %D > LANDAU_MAX_NQ (%D)",Nq,LANDAU_MAX_NQ);
1033   ierr = PetscDSGetDimensions(prob, &Nb);CHKERRQ(ierr);
1034   if (type==4) {
1035     for (c = cStart; c < cEnd; c++) {
1036       ierr = DMLabelSetValue(adaptLabel, c, DM_ADAPT_REFINE);CHKERRQ(ierr);
1037     }
1038     ierr = PetscInfo1(sol, "Phase:%s: Uniform refinement\n","adaptToleranceFEM");CHKERRQ(ierr);
1039   } else if (type==2) {
1040     PetscInt  rCellIdx[8], eCellIdx[64], iCellIdx[64], eMaxIdx = -1, iMaxIdx = -1, nr = 0, nrmax = (dim==3) ? 8 : 2;
1041     PetscReal minRad = PETSC_INFINITY, r, eMinRad = PETSC_INFINITY, iMinRad = PETSC_INFINITY;
1042     for (c = 0; c < 64; c++) { eCellIdx[c] = iCellIdx[c] = -1; }
1043     for (c = cStart; c < cEnd; c++) {
1044       PetscReal    tt, v0[LANDAU_MAX_NQ*3], detJ[LANDAU_MAX_NQ];
1045       ierr = DMPlexComputeCellGeometryFEM(plex, c, quad, v0, NULL, NULL, detJ);CHKERRQ(ierr);
1046       for (qj = 0; qj < Nq; ++qj) {
1047         tt = PetscSqr(v0[dim*qj+0]) + PetscSqr(v0[dim*qj+1]) + PetscSqr(((dim==3) ? v0[dim*qj+2] : 0));
1048         r = PetscSqrtReal(tt);
1049         if (r < minRad - PETSC_SQRT_MACHINE_EPSILON*10.) {
1050           minRad = r;
1051           nr = 0;
1052           rCellIdx[nr++]= c;
1053           ierr = PetscInfo4(sol, "\t\tPhase: adaptToleranceFEM Found first inner r=%e, cell %D, qp %D/%D\n", r, c, qj+1, Nq);CHKERRQ(ierr);
1054         } else if ((r-minRad) < PETSC_SQRT_MACHINE_EPSILON*100. && nr < nrmax) {
1055           for (k=0;k<nr;k++) if (c == rCellIdx[k]) break;
1056           if (k==nr) {
1057             rCellIdx[nr++]= c;
1058             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);
1059           }
1060         }
1061         if (ctx->sphere) {
1062           if ((tt=r-ctx->e_radius) > 0) {
1063             PetscInfo2(sol, "\t\t\t %D cell r=%g\n",c,tt);
1064             if (tt < eMinRad - PETSC_SQRT_MACHINE_EPSILON*100.) {
1065               eMinRad = tt;
1066               eMaxIdx = 0;
1067               eCellIdx[eMaxIdx++] = c;
1068             } else if (eMaxIdx > 0 && (tt-eMinRad) <= PETSC_SQRT_MACHINE_EPSILON && c != eCellIdx[eMaxIdx-1]) {
1069               eCellIdx[eMaxIdx++] = c;
1070             }
1071           }
1072           if ((tt=r-ctx->i_radius) > 0) {
1073             if (tt < iMinRad - 1.e-5) {
1074               iMinRad = tt;
1075               iMaxIdx = 0;
1076               iCellIdx[iMaxIdx++] = c;
1077             } else if (iMaxIdx > 0 && (tt-iMinRad) <= PETSC_SQRT_MACHINE_EPSILON && c != iCellIdx[iMaxIdx-1]) {
1078               iCellIdx[iMaxIdx++] = c;
1079             }
1080           }
1081         }
1082       }
1083     }
1084     for (k=0;k<nr;k++) {
1085       ierr = DMLabelSetValue(adaptLabel, rCellIdx[k], DM_ADAPT_REFINE);CHKERRQ(ierr);
1086     }
1087     if (ctx->sphere) {
1088       for (c = 0; c < eMaxIdx; c++) {
1089         ierr = DMLabelSetValue(adaptLabel, eCellIdx[c], DM_ADAPT_REFINE);CHKERRQ(ierr);
1090         ierr = PetscInfo3(sol, "\t\tPhase:%s: refine sphere e cell %D r=%g\n","adaptToleranceFEM",eCellIdx[c],eMinRad);
1091       }
1092       for (c = 0; c < iMaxIdx; c++) {
1093         ierr = DMLabelSetValue(adaptLabel, iCellIdx[c], DM_ADAPT_REFINE);CHKERRQ(ierr);
1094         ierr = PetscInfo3(sol, "\t\tPhase:%s: refine sphere i cell %D r=%g\n","adaptToleranceFEM",iCellIdx[c],iMinRad);
1095       }
1096     }
1097     ierr = PetscInfo4(sol, "Phase:%s: Adaptive refine origin cells %D,%D r=%g\n","adaptToleranceFEM",rCellIdx[0],rCellIdx[1],minRad);
1098   } else if (type==0 || type==1 || type==3) { /* refine along r=0 axis */
1099     PetscScalar  *coef = NULL;
1100     Vec          coords;
1101     PetscInt     csize,Nv,d,nz;
1102     DM           cdm;
1103     PetscSection cs;
1104     ierr = DMGetCoordinatesLocal(dm, &coords);CHKERRQ(ierr);
1105     ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
1106     ierr = DMGetLocalSection(cdm, &cs);CHKERRQ(ierr);
1107     for (c = cStart; c < cEnd; c++) {
1108       PetscInt doit = 0, outside = 0;
1109       ierr = DMPlexVecGetClosure(cdm, cs, coords, c, &csize, &coef);CHKERRQ(ierr);
1110       Nv = csize/dim;
1111       for (nz = d = 0; d < Nv; d++) {
1112         PetscReal z = PetscRealPart(coef[d*dim + (dim-1)]), x = PetscSqr(PetscRealPart(coef[d*dim + 0])) + ((dim==3) ? PetscSqr(PetscRealPart(coef[d*dim + 1])) : 0);
1113         x = PetscSqrtReal(x);
1114         if (x < PETSC_MACHINE_EPSILON*10. && PetscAbs(z)<PETSC_MACHINE_EPSILON*10.) doit = 1;             /* refine origin */
1115         else if (type==0 && (z < -PETSC_MACHINE_EPSILON*10. || z > ctx->re_radius+PETSC_MACHINE_EPSILON*10.)) outside++;   /* first pass don't refine bottom */
1116         else if (type==1 && (z > ctx->vperp0_radius1 || z < -ctx->vperp0_radius1)) outside++; /* don't refine outside electron refine radius */
1117         else if (type==3 && (z > ctx->vperp0_radius2 || z < -ctx->vperp0_radius2)) outside++; /* don't refine outside ion refine radius */
1118         if (x < PETSC_MACHINE_EPSILON*10.) nz++;
1119       }
1120       ierr = DMPlexVecRestoreClosure(cdm, cs, coords, c, &csize, &coef);CHKERRQ(ierr);
1121       if (doit || (outside<Nv && nz)) {
1122         ierr = DMLabelSetValue(adaptLabel, c, DM_ADAPT_REFINE);CHKERRQ(ierr);
1123       }
1124     }
1125     ierr = PetscInfo1(sol, "Phase:%s: RE refinement\n","adaptToleranceFEM");
1126   }
1127   /* ierr = VecDestroy(&locX);CHKERRQ(ierr); */
1128   ierr = DMDestroy(&plex);CHKERRQ(ierr);
1129   ierr = DMAdaptLabel(dm, adaptLabel, &adaptedDM);CHKERRQ(ierr);
1130   ierr = DMLabelDestroy(&adaptLabel);CHKERRQ(ierr);
1131   *newDM = adaptedDM;
1132   if (adaptedDM) {
1133     if (isForest) {
1134       ierr = DMForestSetAdaptivityForest(adaptedDM,NULL);CHKERRQ(ierr);
1135     }
1136     ierr = DMConvert(adaptedDM, DMPLEX, &plex);CHKERRQ(ierr);
1137     ierr = DMPlexGetHeightStratum(plex,0,&cStart,&cEnd);CHKERRQ(ierr);
1138     ierr = PetscInfo2(sol, "\tPhase: adaptToleranceFEM: %D cells, %d total quadrature points\n",cEnd-cStart,Nq*(cEnd-cStart));CHKERRQ(ierr);
1139     ierr = DMDestroy(&plex);CHKERRQ(ierr);
1140   }
1141   PetscFunctionReturn(0);
1142 }
1143 
1144 static PetscErrorCode adapt(DM *dm, LandauCtx *ctx, Vec *uu)
1145 {
1146   PetscErrorCode  ierr;
1147   PetscInt        type, limits[5] = {ctx->numRERefine,ctx->nZRefine1,ctx->maxRefIts,ctx->nZRefine2,ctx->postAMRRefine};
1148   PetscInt        adaptIter;
1149 
1150   PetscFunctionBegin;
1151   for (type=0;type<5;type++) {
1152     for (adaptIter = 0; adaptIter<limits[type];adaptIter++) {
1153       DM  dmNew = NULL;
1154       ierr = adaptToleranceFEM(ctx->fe[0], *uu, ctx->refineTol, ctx->coarsenTol, type, ctx, &dmNew);CHKERRQ(ierr);
1155       if (!dmNew) SETERRQ(ctx->comm,PETSC_ERR_ARG_WRONG,"should not happen");
1156       else {
1157         ierr = DMDestroy(dm);CHKERRQ(ierr);
1158         ierr = VecDestroy(uu);CHKERRQ(ierr);
1159         ierr = DMCreateGlobalVector(dmNew,uu);CHKERRQ(ierr);
1160         ierr = PetscObjectSetName((PetscObject) *uu, "u");CHKERRQ(ierr);
1161         ierr = LandauSetInitialCondition(dmNew, *uu, ctx);CHKERRQ(ierr);
1162         *dm = dmNew;
1163       }
1164     }
1165   }
1166   PetscFunctionReturn(0);
1167 }
1168 
1169 static PetscErrorCode ProcessOptions(LandauCtx *ctx, const char prefix[])
1170 {
1171   PetscErrorCode    ierr;
1172   PetscBool         flg, sph_flg;
1173   PetscInt          ii,nt,nm,nc;
1174   DM                dummy;
1175 
1176   PetscFunctionBegin;
1177   ierr = DMCreate(ctx->comm,&dummy);CHKERRQ(ierr);
1178   /* get options - initialize context */
1179   ctx->normJ = 0;
1180   ctx->verbose = 1;
1181   ctx->interpolate = PETSC_TRUE;
1182   ctx->gpu_assembly = PETSC_TRUE;
1183   ctx->sphere = PETSC_FALSE;
1184   ctx->inflate = PETSC_FALSE;
1185   ctx->electronShift = 0;
1186   ctx->errorIndicator = NULL;
1187   ctx->radius = 5.; /* electron thermal radius (velocity) */
1188   ctx->re_radius = 0.;
1189   ctx->vperp0_radius1 = 0;
1190   ctx->vperp0_radius2 = 0;
1191   ctx->e_radius = .1;
1192   ctx->i_radius = .01;
1193   ctx->maxRefIts = 5;
1194   ctx->postAMRRefine = 0;
1195   ctx->nZRefine1 = 0;
1196   ctx->nZRefine2 = 0;
1197   ctx->numRERefine = 0;
1198   ctx->aux_bool = PETSC_FALSE;
1199   ctx->num_sections = 3; /* 2, 3 or 4 */
1200   /* species - [0] electrons, [1] one ion species eg, duetarium, [2] heavy impurity ion, ... */
1201   ctx->charges[0] = -1;  /* electron charge (MKS) */
1202   ctx->masses[0] = 1/1835.5; /* temporary value in proton mass */
1203   ctx->n[0] = 1;
1204   ctx->thermal_temps[0] = 1;
1205   /* constants, etc. */
1206   ctx->epsilon0 = 8.8542e-12; /* permittivity of free space (MKS) F/m */
1207   ctx->k = 1.38064852e-23; /* Boltzmann constant (MKS) J/K */
1208   ctx->lnLam = 10;         /* cross section ratio large - small angle collisions */
1209   ctx->n_0 = 1.e20;        /* typical plasma n, but could set it to 1 */
1210   ctx->Ez = 0;
1211   ctx->v_0 = 1; /* in electron thermal velocity */
1212   ctx->subThreadBlockSize = 1; /* for device and maybe OMP */
1213   ctx->numConcurrency = 1; /* for device */
1214   ierr = PetscOptionsBegin(ctx->comm, prefix, "Options for Fokker-Plank-Landau collision operator", "none");CHKERRQ(ierr);
1215   {
1216     char opstring[256];
1217 #if defined(PETSC_HAVE_KOKKOS)
1218     ctx->deviceType = LANDAU_KOKKOS;
1219     ierr = PetscStrcpy(opstring,"kokkos");CHKERRQ(ierr);
1220 #if defined(PETSC_HAVE_CUDA)
1221     ctx->subThreadBlockSize = 8;
1222 #endif
1223 #elif defined(PETSC_HAVE_CUDA)
1224     ctx->deviceType = LANDAU_CUDA;
1225     ierr = PetscStrcpy(opstring,"cuda");CHKERRQ(ierr);
1226 #else
1227     ctx->deviceType = LANDAU_CPU;
1228     ierr = PetscStrcpy(opstring,"cpu");CHKERRQ(ierr);
1229     ctx->subThreadBlockSize = 0;
1230 #endif
1231     ierr = PetscOptionsString("-dm_landau_device_type","Use kernels on 'cpu', 'cuda', or 'kokkos'","plexland.c",opstring,opstring,256,NULL);CHKERRQ(ierr);
1232     ierr = PetscStrcmp("cpu",opstring,&flg);CHKERRQ(ierr);
1233     if (flg) {
1234       ctx->deviceType = LANDAU_CPU;
1235       ctx->subThreadBlockSize = 0;
1236     } else {
1237       ierr = PetscStrcmp("cuda",opstring,&flg);CHKERRQ(ierr);
1238       if (flg) {
1239         ctx->deviceType = LANDAU_CUDA;
1240         ctx->subThreadBlockSize = 0;
1241       } else {
1242         ierr = PetscStrcmp("kokkos",opstring,&flg);CHKERRQ(ierr);
1243         if (flg) ctx->deviceType = LANDAU_KOKKOS;
1244         else SETERRQ1(ctx->comm,PETSC_ERR_ARG_WRONG,"-dm_landau_device_type %s",opstring);
1245       }
1246     }
1247   }
1248   ierr = PetscOptionsBool("-dm_landau_gpu_assembly", "Assemble Jacobian on GPU", "plexland.c", ctx->gpu_assembly, &ctx->gpu_assembly, NULL);CHKERRQ(ierr);
1249   ierr = PetscOptionsReal("-dm_landau_electron_shift","Shift in thermal velocity of electrons","none",ctx->electronShift,&ctx->electronShift, NULL);CHKERRQ(ierr);
1250   ierr = PetscOptionsBool("-dm_landau_sphere", "use sphere/semi-circle domain instead of rectangle", "plexland.c", ctx->sphere, &ctx->sphere, &sph_flg);CHKERRQ(ierr);
1251   ierr = PetscOptionsBool("-dm_landau_inflate", "With sphere, inflate for curved edges (no AMR)", "plexland.c", ctx->inflate, &ctx->inflate, NULL);CHKERRQ(ierr);
1252   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);
1253   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);
1254   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);
1255   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);
1256   ierr = PetscOptionsInt("-dm_landau_amr_post_refine", "Number of levels to uniformly refine after AMR", "plexland.c", ctx->postAMRRefine, &ctx->postAMRRefine, NULL);CHKERRQ(ierr);
1257   ierr = PetscOptionsInt("-dm_landau_verbose", "", "plexland.c", ctx->verbose, &ctx->verbose, NULL);CHKERRQ(ierr);
1258   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);
1259   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);
1260   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);
1261   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);
1262   ierr = PetscOptionsReal("-dm_landau_n_0","Normalization constant for number density","plexland.c",ctx->n_0,&ctx->n_0, NULL);CHKERRQ(ierr);
1263   ierr = PetscOptionsReal("-dm_landau_ln_lambda","Cross section parameter","plexland.c",ctx->lnLam,&ctx->lnLam, NULL);CHKERRQ(ierr);
1264   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);
1265   ierr = PetscOptionsInt("-dm_landau_num_thread_teams", "The number of other concurrent runs to make room for", "plexland.c", ctx->numConcurrency, &ctx->numConcurrency, NULL);CHKERRQ(ierr);
1266 
1267   /* get num species with tempurature*/
1268   {
1269     PetscReal arr[100];
1270     nt = 100;
1271     ierr = PetscOptionsRealArray("-dm_landau_thermal_temps", "Temperature of each species [e,i_0,i_1,...] in keV", "plexland.c", arr, &nt, &flg);CHKERRQ(ierr);
1272     if (flg && nt > LANDAU_MAX_SPECIES) SETERRQ2(ctx->comm,PETSC_ERR_ARG_WRONG,"-thermal_temps ,t1,t2,.. number of species %D > MAX %D",nt,LANDAU_MAX_SPECIES);
1273   }
1274   nt = LANDAU_MAX_SPECIES;
1275   for (ii=1;ii<LANDAU_MAX_SPECIES;ii++) {
1276     ctx->thermal_temps[ii] = 1.;
1277     ctx->charges[ii] = 1;
1278     ctx->masses[ii] = 1;
1279     ctx->n[ii] = (ii==1) ? 1 : 0;
1280   }
1281   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);
1282   if (flg) {
1283     PetscInfo1(dummy, "num_species set to number of thermal temps provided (%D)\n",nt);
1284     ctx->num_species = nt;
1285   } else SETERRQ(ctx->comm,PETSC_ERR_ARG_WRONG,"-dm_landau_thermal_temps ,t1,t2,.. must be provided to set the number of species");
1286   for (ii=0;ii<ctx->num_species;ii++) ctx->thermal_temps[ii] *= 1.1604525e7; /* convert to Kelvin */
1287   nm = LANDAU_MAX_SPECIES-1;
1288   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);
1289   if (flg && nm != ctx->num_species-1) {
1290     SETERRQ2(ctx->comm,PETSC_ERR_ARG_WRONG,"num ion masses %D != num species %D",nm,ctx->num_species-1);
1291   }
1292   nm = LANDAU_MAX_SPECIES;
1293   ierr = PetscOptionsRealArray("-dm_landau_n", "Normalized (by -n_0) number density of each species", "plexland.c", ctx->n, &nm, &flg);CHKERRQ(ierr);
1294   if (flg && nm != ctx->num_species) SETERRQ2(ctx->comm,PETSC_ERR_ARG_WRONG,"wrong num n: %D != num species %D",nm,ctx->num_species);
1295   ctx->n_0 *= ctx->n[0]; /* normalized number density */
1296   for (ii=1;ii<ctx->num_species;ii++) ctx->n[ii] = ctx->n[ii]/ctx->n[0];
1297   ctx->n[0] = 1;
1298   for (ii=0;ii<LANDAU_MAX_SPECIES;ii++) ctx->masses[ii] *= 1.6720e-27; /* scale by proton mass kg */
1299   ctx->masses[0] = 9.10938356e-31; /* electron mass kg (should be about right already) */
1300   ctx->m_0 = ctx->masses[0]; /* arbitrary reference mass, electrons */
1301   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);
1302   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) */
1303   nc = LANDAU_MAX_SPECIES-1;
1304   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);
1305   if (flg && nc != ctx->num_species-1) SETERRQ2(ctx->comm,PETSC_ERR_ARG_WRONG,"num charges %D != num species %D",nc,ctx->num_species-1);
1306   for (ii=0;ii<LANDAU_MAX_SPECIES;ii++) ctx->charges[ii] *= 1.6022e-19; /* electron/proton charge (MKS) */
1307   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 */
1308   /* geometry */
1309   for (ii=0;ii<ctx->num_species;ii++) ctx->refineTol[ii]  = PETSC_MAX_REAL;
1310   for (ii=0;ii<ctx->num_species;ii++) ctx->coarsenTol[ii] = 0.;
1311   ii = LANDAU_MAX_SPECIES;
1312   ierr = PetscOptionsRealArray("-dm_landau_refine_tol","tolerance for refining cells in AMR","plexland.c",ctx->refineTol, &ii, &flg);CHKERRQ(ierr);
1313   if (flg && ii != ctx->num_species) ierr = PetscInfo2(dummy, "Phase: Warning, #refine_tol %D != num_species %D\n",ii,ctx->num_species);CHKERRQ(ierr);
1314   ii = LANDAU_MAX_SPECIES;
1315   ierr = PetscOptionsRealArray("-dm_landau_coarsen_tol","tolerance for coarsening cells in AMR","plexland.c",ctx->coarsenTol, &ii, &flg);CHKERRQ(ierr);
1316   if (flg && ii != ctx->num_species) ierr = PetscInfo2(dummy, "Phase: Warning, #coarsen_tol %D != num_species %D\n",ii,ctx->num_species);CHKERRQ(ierr);
1317   ierr = PetscOptionsReal("-dm_landau_domain_radius","Phase space size in units of electron thermal velocity","plexland.c",ctx->radius,&ctx->radius, &flg);CHKERRQ(ierr);
1318   if (flg && ctx->radius <= 0) { /* negative is ratio of c */
1319     if (ctx->radius == 0) ctx->radius = 0.75;
1320     else ctx->radius = -ctx->radius;
1321     ctx->radius = ctx->radius*299792458.0/ctx->v_0;
1322     ierr = PetscInfo1(dummy, "Change domain radius to %e\n",ctx->radius);CHKERRQ(ierr);
1323   }
1324   ierr = PetscOptionsReal("-dm_landau_i_radius","Ion thermal velocity, used for circular meshes","plexland.c",ctx->i_radius,&ctx->i_radius, &flg);CHKERRQ(ierr);
1325   if (flg && !sph_flg) ctx->sphere = PETSC_TRUE; /* you gave me an ion radius but did not set sphere, user error really */
1326   if (!flg) {
1327     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 */
1328   }
1329   ierr = PetscOptionsReal("-dm_landau_e_radius","Electron thermal velocity, used for circular meshes","plexland.c",ctx->e_radius,&ctx->e_radius, &flg);CHKERRQ(ierr);
1330   if (flg && !sph_flg) ctx->sphere = PETSC_TRUE; /* you gave me an e radius but did not set sphere, user error really */
1331   if (!flg) {
1332     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 */
1333   }
1334   if (ctx->sphere && (ctx->e_radius <= ctx->i_radius || ctx->radius <= ctx->e_radius)) SETERRQ3(ctx->comm,PETSC_ERR_ARG_WRONG,"bad radii: %g < %g < %g",ctx->i_radius,ctx->e_radius,ctx->radius);
1335   ierr = PetscOptionsInt("-dm_landau_sub_thread_block_size", "Number of threads in Kokkos integration point subblock", "plexland.c", ctx->subThreadBlockSize, &ctx->subThreadBlockSize, NULL);CHKERRQ(ierr);
1336   ierr = PetscOptionsEnd();CHKERRQ(ierr);
1337   for (ii=ctx->num_species;ii<LANDAU_MAX_SPECIES;ii++) ctx->masses[ii] = ctx->thermal_temps[ii]  = ctx->charges[ii] = 0;
1338   if (ctx->verbose > 0) {
1339     ierr = PetscPrintf(ctx->comm, "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);
1340     ierr = PetscPrintf(ctx->comm, "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);
1341     ierr = PetscPrintf(ctx->comm, "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);
1342   }
1343   ierr = DMDestroy(&dummy);CHKERRQ(ierr);
1344   {
1345     PetscMPIInt    rank;
1346     ierr = MPI_Comm_rank(ctx->comm, &rank);CHKERRMPI(ierr);
1347     /* PetscLogStage  setup_stage; */
1348     ierr = PetscLogEventRegister("Landau Jacobian", DM_CLASSID, &ctx->events[0]);CHKERRQ(ierr); /* 0 */
1349     ierr = PetscLogEventRegister(" Initialize", DM_CLASSID, &ctx->events[10]);CHKERRQ(ierr); /* 10 */
1350     ierr = PetscLogEventRegister("  IP Data-jac", DM_CLASSID, &ctx->events[7]);CHKERRQ(ierr); /* 7 */
1351     ierr = PetscLogEventRegister(" Kernal-init", DM_CLASSID, &ctx->events[3]);CHKERRQ(ierr); /* 3 */
1352     ierr = PetscLogEventRegister(" GPU Kernel", DM_CLASSID, &ctx->events[4]);CHKERRQ(ierr); /* 4 */
1353     ierr = PetscLogEventRegister(" Copy to CPU", DM_CLASSID, &ctx->events[5]);CHKERRQ(ierr); /* 5 */
1354     ierr = PetscLogEventRegister(" Jac-assemble", DM_CLASSID, &ctx->events[6]);CHKERRQ(ierr); /* 6 */
1355     ierr = PetscLogEventRegister(" Jac-f-df", DM_CLASSID, &ctx->events[8]);CHKERRQ(ierr); /* 8 */
1356     ierr = PetscLogEventRegister(" Jac asmbl setup", DM_CLASSID, &ctx->events[2]);CHKERRQ(ierr); /* 2 */
1357     ierr = PetscLogEventRegister("Mass Operator", DM_CLASSID, &ctx->events[9]);CHKERRQ(ierr); /* 9 */
1358     ierr = PetscLogEventRegister("  IP Data-mass", DM_CLASSID, &ctx->events[1]);CHKERRQ(ierr); /* 1 */
1359 
1360     if (rank) { /* turn off output stuff for duplicate runs - do we need to add the prefix to all this? */
1361       ierr = PetscOptionsClearValue(NULL,"-snes_converged_reason");CHKERRQ(ierr);
1362       ierr = PetscOptionsClearValue(NULL,"-ksp_converged_reason");CHKERRQ(ierr);
1363       ierr = PetscOptionsClearValue(NULL,"-snes_monitor");CHKERRQ(ierr);
1364       ierr = PetscOptionsClearValue(NULL,"-ksp_monitor");CHKERRQ(ierr);
1365       ierr = PetscOptionsClearValue(NULL,"-ts_monitor");CHKERRQ(ierr);
1366       ierr = PetscOptionsClearValue(NULL,"-ts_adapt_monitor");CHKERRQ(ierr);
1367       ierr = PetscOptionsClearValue(NULL,"-dm_landau_amr_dm_view");CHKERRQ(ierr);
1368       ierr = PetscOptionsClearValue(NULL,"-dm_landau_amr_vec_view");CHKERRQ(ierr);
1369       ierr = PetscOptionsClearValue(NULL,"-dm_landau_pre_dm_view");CHKERRQ(ierr);
1370       ierr = PetscOptionsClearValue(NULL,"-dm_landau_pre_vec_view");CHKERRQ(ierr);
1371       ierr = PetscOptionsClearValue(NULL,"-info");CHKERRQ(ierr);
1372     }
1373   }
1374   PetscFunctionReturn(0);
1375 }
1376 
1377 /*@C
1378  LandauCreateVelocitySpace - Create a DMPlex velocity space mesh
1379 
1380  Collective on comm
1381 
1382  Input Parameters:
1383  +   comm  - The MPI communicator
1384  .   dim - velocity space dimension (2 for axisymmetric, 3 for full 3X + 3V solver)
1385  -   prefix - prefix for options
1386 
1387  Output Parameter:
1388  .   dm  - The DM object representing the mesh
1389  +   X - A vector (user destroys)
1390  -   J - Optional matrix (object destroys)
1391 
1392  Level: beginner
1393 
1394  .keywords: mesh
1395  .seealso: DMPlexCreate(), LandauDestroyVelocitySpace()
1396  @*/
1397 PetscErrorCode LandauCreateVelocitySpace(MPI_Comm comm, PetscInt dim, const char prefix[], Vec *X, Mat *J, DM *dm)
1398 {
1399   PetscErrorCode ierr;
1400   LandauCtx      *ctx;
1401   PetscBool      prealloc_only,flg;
1402 
1403   PetscFunctionBegin;
1404   if (dim!=2 && dim!=3) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Only 2D and 3D supported");
1405   ctx = (LandauCtx*)malloc(sizeof(LandauCtx));
1406   ctx->comm = comm; /* used for diagnostics and global errors */
1407   /* process options */
1408   ierr = ProcessOptions(ctx,prefix);CHKERRQ(ierr);
1409   /* Create Mesh */
1410   ierr = LandauDMCreateVMesh(PETSC_COMM_SELF, dim, prefix, ctx, dm);CHKERRQ(ierr);
1411   prealloc_only = (*dm)->prealloc_only;
1412   ierr = DMViewFromOptions(*dm,NULL,"-dm_landau_pre_dm_view");CHKERRQ(ierr);
1413   ierr = DMSetApplicationContext(*dm, ctx);CHKERRQ(ierr);
1414   /* create FEM */
1415   ierr = SetupDS(*dm,dim,ctx);CHKERRQ(ierr);
1416   /* set initial state */
1417   ierr = DMCreateGlobalVector(*dm,X);CHKERRQ(ierr);
1418   ierr = PetscObjectSetName((PetscObject) *X, "u");CHKERRQ(ierr);
1419   /* initial static refinement, no solve */
1420   ierr = LandauSetInitialCondition(*dm, *X, ctx);CHKERRQ(ierr);
1421   ierr = VecViewFromOptions(*X, NULL, "-dm_landau_pre_vec_view");CHKERRQ(ierr);
1422   /* forest refinement */
1423   if (ctx->errorIndicator) {
1424     /* AMR */
1425     ierr = adapt(dm,ctx,X);CHKERRQ(ierr);
1426     if ((*dm)->prealloc_only != prealloc_only) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"(*dm)->prealloc_only != prealloc_only");
1427     ierr = DMViewFromOptions(*dm,NULL,"-dm_landau_amr_dm_view");CHKERRQ(ierr);
1428     ierr = VecViewFromOptions(*X, NULL, "-dm_landau_amr_vec_view");CHKERRQ(ierr);
1429   }
1430   ierr = DMSetApplicationContext(*dm, ctx);CHKERRQ(ierr);
1431   ctx->dmv = *dm;
1432   if (ctx->dmv->prealloc_only != prealloc_only) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"ctx->dmv->prealloc_only != prealloc_only");
1433   ierr = DMCreateMatrix(ctx->dmv, &ctx->J);CHKERRQ(ierr);
1434   ierr = MatSetOption(ctx->J, MAT_IGNORE_ZERO_ENTRIES, PETSC_TRUE);CHKERRQ(ierr);
1435   ierr = MatSetOption(ctx->J, MAT_STRUCTURALLY_SYMMETRIC, PETSC_TRUE);CHKERRQ(ierr);
1436   if (J) *J = ctx->J;
1437   /* check for types that we need */
1438 #if defined(PETSC_HAVE_KOKKOS)
1439   if (ctx->deviceType == LANDAU_CPU) {
1440     ierr = PetscObjectTypeCompareAny((PetscObject)ctx->J,&flg,MATSEQAIJKOKKOS,MATMPIAIJKOKKOS,MATAIJKOKKOS,"");CHKERRQ(ierr);
1441     //if (flg) SETERRQ(ctx->comm,PETSC_ERR_ARG_WRONG,"with device=cpu must not use '-dm_mat_type aijkokkos -dm_vec_type kokkos' for GPU assembly and Kokkos");
1442   }
1443 #elif defined(PETSC_HAVE_CUDA)
1444   if (ctx->deviceType == LANDAU_CPU) {
1445     ierr = PetscObjectTypeCompareAny((PetscObject)ctx->J,&flg,MATSEQAIJCUSPARSE,MATMPIAIJCUSPARSE,MATAIJCUSPARSE,"");CHKERRQ(ierr);
1446     //if (flg) SETERRQ(ctx->comm,PETSC_ERR_ARG_WRONG,"with device=cpu must not use '-dm_mat_type aijcusparse -dm_vec_type cuda' for GPU assembly and Cuda");
1447   }
1448 #endif
1449   if (ctx->gpu_assembly) { /* we need GPU object with GPU assembly */
1450     if (ctx->deviceType == LANDAU_CUDA) {
1451       ierr = PetscObjectTypeCompareAny((PetscObject)ctx->J,&flg,MATSEQAIJCUSPARSE,MATMPIAIJCUSPARSE,MATAIJCUSPARSE,"");CHKERRQ(ierr);
1452       if (!flg) SETERRQ(ctx->comm,PETSC_ERR_ARG_WRONG,"must use '-dm_mat_type aijcusparse -dm_vec_type cuda' for GPU assembly and Cuda");
1453     } else if (ctx->deviceType == LANDAU_KOKKOS) {
1454       ierr = PetscObjectTypeCompareAny((PetscObject)ctx->J,&flg,MATSEQAIJKOKKOS,MATMPIAIJKOKKOS,MATAIJKOKKOS,"");CHKERRQ(ierr);
1455 #if defined(PETSC_HAVE_KOKKOS_KERNELS)
1456       if (!flg) SETERRQ(ctx->comm,PETSC_ERR_ARG_WRONG,"must use '-dm_mat_type aijkokkos -dm_vec_type kokkos' for GPU assembly and Kokkos");
1457 #else
1458       if (!flg) SETERRQ(ctx->comm,PETSC_ERR_ARG_WRONG,"must configure with '--download-kokkos-kernels=1' for GPU assembly and Kokkos");
1459 #endif
1460     }
1461   }
1462   PetscFunctionReturn(0);
1463 }
1464 
1465 /*@
1466  LandauDestroyVelocitySpace - Destroy a DMPlex velocity space mesh
1467 
1468  Collective on dm
1469 
1470  Input/Output Parameters:
1471  .   dm - the dm to destroy
1472 
1473  Level: beginner
1474 
1475  .keywords: mesh
1476  .seealso: LandauCreateVelocitySpace()
1477  @*/
1478 PetscErrorCode LandauDestroyVelocitySpace(DM *dm)
1479 {
1480   PetscErrorCode ierr,ii;
1481   LandauCtx      *ctx;
1482   PetscContainer container = NULL;
1483   PetscFunctionBegin;
1484   ierr = DMGetApplicationContext(*dm, &ctx);CHKERRQ(ierr);
1485   ierr = PetscObjectQuery((PetscObject)ctx->J,"coloring", (PetscObject*)&container);CHKERRQ(ierr);
1486   if (container) {
1487     ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
1488   }
1489   ierr = MatDestroy(&ctx->M);CHKERRQ(ierr);
1490   ierr = MatDestroy(&ctx->J);CHKERRQ(ierr);
1491   for (ii=0;ii<ctx->num_species;ii++) {
1492     ierr = PetscFEDestroy(&ctx->fe[ii]);CHKERRQ(ierr);
1493   }
1494   free(ctx);
1495   ierr = DMDestroy(dm);CHKERRQ(ierr);
1496   PetscFunctionReturn(0);
1497 }
1498 
1499 /* < v, ru > */
1500 static void f0_s_den(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1501                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1502                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1503                      PetscReal t, const PetscReal x[],  PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0)
1504 {
1505   PetscInt ii = (PetscInt)PetscRealPart(constants[0]);
1506   f0[0] = u[ii];
1507 }
1508 
1509 /* < v, ru > */
1510 static void f0_s_mom(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1511                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1512                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1513                      PetscReal t, const PetscReal x[],  PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0)
1514 {
1515   PetscInt ii = (PetscInt)PetscRealPart(constants[0]), jj = (PetscInt)PetscRealPart(constants[1]);
1516   f0[0] = x[jj]*u[ii]; /* x momentum */
1517 }
1518 
1519 static void f0_s_v2(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1520                     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1521                     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1522                     PetscReal t, const PetscReal x[],  PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0)
1523 {
1524   PetscInt i, ii = (PetscInt)PetscRealPart(constants[0]);
1525   double tmp1 = 0.;
1526   for (i = 0; i < dim; ++i) tmp1 += x[i]*x[i];
1527   f0[0] = tmp1*u[ii];
1528 }
1529 
1530 /* < v, ru > */
1531 static void f0_s_rden(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1532                       const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1533                       const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1534                       PetscReal t, const PetscReal x[],  PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0)
1535 {
1536   PetscInt ii = (PetscInt)PetscRealPart(constants[0]);
1537   f0[0] = 2.*PETSC_PI*x[0]*u[ii];
1538 }
1539 
1540 /* < v, ru > */
1541 static void f0_s_rmom(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1542                       const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1543                       const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1544                       PetscReal t, const PetscReal x[],  PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0)
1545 {
1546   PetscInt ii = (PetscInt)PetscRealPart(constants[0]);
1547   f0[0] = 2.*PETSC_PI*x[0]*x[1]*u[ii];
1548 }
1549 
1550 static void f0_s_rv2(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1551                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1552                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1553                      PetscReal t, const PetscReal x[],  PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0)
1554 {
1555   PetscInt ii = (PetscInt)PetscRealPart(constants[0]);
1556   f0[0] =  2.*PETSC_PI*x[0]*(x[0]*x[0] + x[1]*x[1])*u[ii];
1557 }
1558 
1559 /*@
1560  LandauPrintNorms - collects moments and prints them
1561 
1562  Collective on dm
1563 
1564  Input Parameters:
1565  +   X  - the state
1566  -   stepi - current step to print
1567 
1568  Level: beginner
1569 
1570  .keywords: mesh
1571  .seealso: LandauCreateVelocitySpace()
1572  @*/
1573 PetscErrorCode LandauPrintNorms(Vec X, PetscInt stepi)
1574 {
1575   PetscErrorCode ierr;
1576   LandauCtx      *ctx;
1577   PetscDS        prob;
1578   DM             plex,dm;
1579   PetscInt       cStart, cEnd, dim, ii;
1580   PetscScalar    xmomentumtot=0, ymomentumtot=0, zmomentumtot=0, energytot=0, densitytot=0, tt[LANDAU_MAX_SPECIES];
1581   PetscScalar    xmomentum[LANDAU_MAX_SPECIES],  ymomentum[LANDAU_MAX_SPECIES],  zmomentum[LANDAU_MAX_SPECIES], energy[LANDAU_MAX_SPECIES], density[LANDAU_MAX_SPECIES];
1582 
1583   PetscFunctionBegin;
1584   ierr = VecGetDM(X, &dm);CHKERRQ(ierr);
1585   if (!dm) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "no DM");
1586   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1587   ierr = DMGetApplicationContext(dm, &ctx);CHKERRQ(ierr);
1588   if (!ctx) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "no context");
1589   ierr = DMConvert(ctx->dmv, DMPLEX, &plex);CHKERRQ(ierr);
1590   ierr = DMCreateDS(plex);CHKERRQ(ierr);
1591   ierr = DMGetDS(plex, &prob);CHKERRQ(ierr);
1592   /* print momentum and energy */
1593   for (ii=0;ii<ctx->num_species;ii++) {
1594     PetscScalar user[2] = { (PetscScalar)ii, (PetscScalar)ctx->charges[ii]};
1595     ierr = PetscDSSetConstants(prob, 2, user);CHKERRQ(ierr);
1596     if (dim==2) { /* 2/3X + 3V (cylindrical coordinates) */
1597       ierr = PetscDSSetObjective(prob, 0, &f0_s_rden);CHKERRQ(ierr);
1598       ierr = DMPlexComputeIntegralFEM(plex,X,tt,ctx);CHKERRQ(ierr);
1599       density[ii] = tt[0]*ctx->n_0*ctx->charges[ii];
1600       ierr = PetscDSSetObjective(prob, 0, &f0_s_rmom);CHKERRQ(ierr);
1601       ierr = DMPlexComputeIntegralFEM(plex,X,tt,ctx);CHKERRQ(ierr);
1602       zmomentum[ii] = tt[0]*ctx->n_0*ctx->v_0*ctx->masses[ii];
1603       ierr = PetscDSSetObjective(prob, 0, &f0_s_rv2);CHKERRQ(ierr);
1604       ierr = DMPlexComputeIntegralFEM(plex,X,tt,ctx);CHKERRQ(ierr);
1605       energy[ii] = tt[0]*0.5*ctx->n_0*ctx->v_0*ctx->v_0*ctx->masses[ii];
1606       zmomentumtot += zmomentum[ii];
1607       energytot  += energy[ii];
1608       densitytot += density[ii];
1609       ierr = PetscPrintf(ctx->comm, "%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);
1610     } else { /* 2/3X + 3V */
1611       ierr = PetscDSSetObjective(prob, 0, &f0_s_den);CHKERRQ(ierr);
1612       ierr = DMPlexComputeIntegralFEM(plex,X,tt,ctx);CHKERRQ(ierr);
1613       density[ii] = tt[0]*ctx->n_0*ctx->charges[ii];
1614       ierr = PetscDSSetObjective(prob, 0, &f0_s_mom);CHKERRQ(ierr);
1615       user[1] = 0;
1616       ierr = DMPlexComputeIntegralFEM(plex,X,tt,ctx);CHKERRQ(ierr);
1617       xmomentum[ii]  = tt[0]*ctx->n_0*ctx->v_0*ctx->masses[ii];
1618       user[1] = 1;
1619       ierr = DMPlexComputeIntegralFEM(plex,X,tt,ctx);CHKERRQ(ierr);
1620       ymomentum[ii] = tt[0]*ctx->n_0*ctx->v_0*ctx->masses[ii];
1621       user[1] = 2;
1622       ierr = DMPlexComputeIntegralFEM(plex,X,tt,ctx);CHKERRQ(ierr);
1623       zmomentum[ii] = tt[0]*ctx->n_0*ctx->v_0*ctx->masses[ii];
1624       ierr = PetscDSSetObjective(prob, 0, &f0_s_v2);CHKERRQ(ierr);
1625       ierr = DMPlexComputeIntegralFEM(plex,X,tt,ctx);CHKERRQ(ierr);
1626       energy[ii]    = 0.5*tt[0]*ctx->n_0*ctx->v_0*ctx->v_0*ctx->masses[ii];
1627       ierr = PetscPrintf(ctx->comm, "%3D) species %D: density=%20.13e, x-momentum=%20.13e, y-momentum=%20.13e, z-momentum=%20.13e, energy=%21.13e",
1628                          stepi,ii,PetscRealPart(density[ii]),PetscRealPart(xmomentum[ii]),PetscRealPart(ymomentum[ii]),PetscRealPart(zmomentum[ii]),PetscRealPart(energy[ii]));CHKERRQ(ierr);
1629       xmomentumtot += xmomentum[ii];
1630       ymomentumtot += ymomentum[ii];
1631       zmomentumtot += zmomentum[ii];
1632       energytot  += energy[ii];
1633       densitytot += density[ii];
1634     }
1635     if (ctx->num_species>1) PetscPrintf(ctx->comm, "\n");
1636   }
1637   /* totals */
1638   ierr = DMPlexGetHeightStratum(plex,0,&cStart,&cEnd);CHKERRQ(ierr);
1639   ierr = DMDestroy(&plex);CHKERRQ(ierr);
1640   if (ctx->num_species>1) {
1641     if (dim==2) {
1642       ierr = PetscPrintf(ctx->comm, "\t%3D) Total: charge density=%21.13e, momentum=%21.13e, energy=%21.13e (m_i[0]/m_e = %g, %D cells)",
1643                          stepi,(double)PetscRealPart(densitytot),(double)PetscRealPart(zmomentumtot),(double)PetscRealPart(energytot),(double)(ctx->masses[1]/ctx->masses[0]),cEnd-cStart);CHKERRQ(ierr);
1644     } else {
1645       ierr = PetscPrintf(ctx->comm, "\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)",
1646                          stepi,(double)PetscRealPart(densitytot),(double)PetscRealPart(xmomentumtot),(double)PetscRealPart(ymomentumtot),(double)PetscRealPart(zmomentumtot),(double)PetscRealPart(energytot),(double)(ctx->masses[1]/ctx->masses[0]),cEnd-cStart);CHKERRQ(ierr);
1647     }
1648   } else {
1649     ierr = PetscPrintf(ctx->comm, " -- %D cells",cEnd-cStart);CHKERRQ(ierr);
1650   }
1651   if (ctx->verbose > 1) {ierr = PetscPrintf(ctx->comm,", %D sub (vector) threads\n",ctx->subThreadBlockSize);CHKERRQ(ierr);}
1652   else {ierr = PetscPrintf(ctx->comm,"\n");CHKERRQ(ierr);}
1653   PetscFunctionReturn(0);
1654 }
1655 
1656 static PetscErrorCode destroy_coloring (void *is)
1657 {
1658   ISColoring tmp = (ISColoring)is;
1659   return ISColoringDestroy(&tmp);
1660 }
1661 
1662 /*@
1663  LandauCreateColoring - create a coloring and add to matrix (Landau context used just for 'print' flag, should be in DMPlex)
1664 
1665  Collective on JacP
1666 
1667  Input Parameters:
1668  +   JacP  - matrix to add coloring to
1669  -   plex - The DM
1670 
1671  Output Parameter:
1672  .   container  - Container with coloring
1673 
1674  Level: beginner
1675 
1676  .keywords: mesh
1677  .seealso: LandauCreateVelocitySpace()
1678  @*/
1679 PetscErrorCode LandauCreateColoring(Mat JacP, DM plex, PetscContainer *container)
1680 {
1681   PetscErrorCode  ierr;
1682   PetscInt        dim,cell,i,ej,nc,Nv,totDim,numGCells,cStart,cEnd;
1683   ISColoring      iscoloring = NULL;
1684   Mat             G,Q;
1685   PetscScalar     ones[128];
1686   MatColoring     mc;
1687   IS             *is;
1688   PetscInt        csize,colour,j,k;
1689   const PetscInt *indices;
1690   PetscInt       numComp[1];
1691   PetscInt       numDof[4];
1692   PetscFE        fe;
1693   DM             colordm;
1694   PetscSection   csection, section, globalSection;
1695   PetscDS        prob;
1696   LandauCtx      *ctx;
1697 
1698   PetscFunctionBegin;
1699   ierr = DMGetApplicationContext(plex, &ctx);CHKERRQ(ierr);
1700   ierr = DMGetLocalSection(plex, &section);CHKERRQ(ierr);
1701   ierr = DMGetGlobalSection(plex, &globalSection);CHKERRQ(ierr);
1702   ierr = DMGetDimension(plex, &dim);CHKERRQ(ierr);
1703   ierr = DMGetDS(plex, &prob);CHKERRQ(ierr);
1704   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
1705   ierr = DMPlexGetHeightStratum(plex,0,&cStart,&cEnd);CHKERRQ(ierr);
1706   numGCells = cEnd - cStart;
1707   /* create cell centered DM */
1708   ierr = DMClone(plex, &colordm);CHKERRQ(ierr);
1709   ierr = PetscFECreateDefault(PetscObjectComm((PetscObject) plex), dim, 1, PETSC_FALSE, "color_", PETSC_DECIDE, &fe);CHKERRQ(ierr);
1710   ierr = PetscObjectSetName((PetscObject) fe, "color");CHKERRQ(ierr);
1711   ierr = DMSetField(colordm, 0, NULL, (PetscObject)fe);CHKERRQ(ierr);
1712   ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
1713   for (i = 0; i < (dim+1); ++i) numDof[i] = 0;
1714   numDof[dim] = 1;
1715   numComp[0] = 1;
1716   ierr = DMPlexCreateSection(colordm, NULL, numComp, numDof, 0, NULL, NULL, NULL, NULL, &csection);CHKERRQ(ierr);
1717   ierr = PetscSectionSetFieldName(csection, 0, "color");CHKERRQ(ierr);
1718   ierr = DMSetLocalSection(colordm, csection);CHKERRQ(ierr);
1719   ierr = DMViewFromOptions(colordm,NULL,"-color_dm_view");CHKERRQ(ierr);
1720   /* get vertex to element map Q and colroing graph G */
1721   ierr = MatGetSize(JacP,NULL,&Nv);CHKERRQ(ierr);
1722   ierr = MatCreateAIJ(PETSC_COMM_SELF,PETSC_DECIDE,PETSC_DECIDE,numGCells,Nv,totDim,NULL,0,NULL,&Q);CHKERRQ(ierr);
1723   for (i=0;i<128;i++) ones[i] = 1.0;
1724   for (cell = cStart, ej = 0 ; cell < cEnd; ++cell, ++ej) {
1725     PetscInt numindices,*indices;
1726     ierr = DMPlexGetClosureIndices(plex, section, globalSection, cell, PETSC_TRUE, &numindices, &indices, NULL, NULL);CHKERRQ(ierr);
1727     if (numindices>128) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "too many indices. %D > %D",numindices,128);
1728     ierr = MatSetValues(Q,1,&ej,numindices,indices,ones,ADD_VALUES);CHKERRQ(ierr);
1729     ierr = DMPlexRestoreClosureIndices(plex, section, globalSection, cell, PETSC_TRUE, &numindices, &indices, NULL, NULL);CHKERRQ(ierr);
1730   }
1731   ierr = MatAssemblyBegin(Q, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1732   ierr = MatAssemblyEnd(Q, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1733   ierr = MatMatTransposeMult(Q,Q,MAT_INITIAL_MATRIX,4.0,&G);CHKERRQ(ierr);
1734   ierr = PetscObjectSetName((PetscObject) Q, "Q");CHKERRQ(ierr);
1735   ierr = PetscObjectSetName((PetscObject) G, "coloring graph");CHKERRQ(ierr);
1736   ierr = MatViewFromOptions(G,NULL,"-coloring_mat_view");CHKERRQ(ierr);
1737   ierr = MatViewFromOptions(Q,NULL,"-coloring_mat_view");CHKERRQ(ierr);
1738   ierr = MatDestroy(&Q);CHKERRQ(ierr);
1739   /* coloring */
1740   ierr = MatColoringCreate(G,&mc);CHKERRQ(ierr);
1741   ierr = MatColoringSetDistance(mc,1);CHKERRQ(ierr);
1742   ierr = MatColoringSetType(mc,MATCOLORINGJP);CHKERRQ(ierr);
1743   ierr = MatColoringSetFromOptions(mc);CHKERRQ(ierr);
1744   ierr = MatColoringApply(mc,&iscoloring);CHKERRQ(ierr);
1745   ierr = MatColoringDestroy(&mc);CHKERRQ(ierr);
1746   /* view */
1747   ierr = ISColoringViewFromOptions(iscoloring,NULL,"-coloring_is_view");CHKERRQ(ierr);
1748   ierr = ISColoringGetIS(iscoloring,PETSC_USE_POINTER,&nc,&is);CHKERRQ(ierr);
1749   if (ctx && ctx->verbose > 2) {
1750     PetscViewer    viewer;
1751     Vec            color_vec, eidx_vec;
1752     ierr = DMGetGlobalVector(colordm, &color_vec);CHKERRQ(ierr);
1753     ierr = DMGetGlobalVector(colordm, &eidx_vec);CHKERRQ(ierr);
1754     for (colour=0; colour<nc; colour++) {
1755       ierr = ISGetLocalSize(is[colour],&csize);CHKERRQ(ierr);
1756       ierr = ISGetIndices(is[colour],&indices);CHKERRQ(ierr);
1757       for (j=0; j<csize; j++) {
1758         PetscScalar v = (PetscScalar)colour;
1759         k = indices[j];
1760         ierr = VecSetValues(color_vec,1,&k,&v,INSERT_VALUES);
1761         v = (PetscScalar)k;
1762         ierr = VecSetValues(eidx_vec,1,&k,&v,INSERT_VALUES);
1763       }
1764       ierr = ISRestoreIndices(is[colour],&indices);CHKERRQ(ierr);
1765     }
1766     /* view */
1767     ierr = PetscViewerVTKOpen(ctx->comm, "color.vtu", FILE_MODE_WRITE, &viewer);CHKERRQ(ierr);
1768     ierr = PetscObjectSetName((PetscObject) color_vec, "color");CHKERRQ(ierr);
1769     ierr = VecView(color_vec, viewer);CHKERRQ(ierr);
1770     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1771     ierr = PetscViewerVTKOpen(ctx->comm, "eidx.vtu", FILE_MODE_WRITE, &viewer);CHKERRQ(ierr);
1772     ierr = PetscObjectSetName((PetscObject) eidx_vec, "element-idx");CHKERRQ(ierr);
1773     ierr = VecView(eidx_vec, viewer);CHKERRQ(ierr);
1774     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1775     ierr = DMRestoreGlobalVector(colordm, &color_vec);CHKERRQ(ierr);
1776     ierr = DMRestoreGlobalVector(colordm, &eidx_vec);CHKERRQ(ierr);
1777   }
1778   ierr = PetscSectionDestroy(&csection);CHKERRQ(ierr);
1779   ierr = DMDestroy(&colordm);CHKERRQ(ierr);
1780   ierr = ISColoringRestoreIS(iscoloring,PETSC_USE_POINTER,&is);CHKERRQ(ierr);
1781   ierr = MatDestroy(&G);CHKERRQ(ierr);
1782   /* stash coloring */
1783   ierr = PetscContainerCreate(PETSC_COMM_SELF, container);CHKERRQ(ierr);
1784   ierr = PetscContainerSetPointer(*container,(void*)iscoloring);CHKERRQ(ierr);
1785   ierr = PetscContainerSetUserDestroy(*container, destroy_coloring);CHKERRQ(ierr);
1786   ierr = PetscObjectCompose((PetscObject)JacP,"coloring",(PetscObject)*container);CHKERRQ(ierr);
1787   if (ctx && ctx->verbose > 0) {
1788     ierr = PetscPrintf(ctx->comm, "Made coloring with %D colors\n", nc);CHKERRQ(ierr);
1789   }
1790   PetscFunctionReturn(0);
1791 }
1792 
1793 PetscErrorCode LandauAssembleOpenMP(PetscInt cStart, PetscInt cEnd, PetscInt totDim, DM plex, PetscSection section, PetscSection globalSection, Mat JacP, PetscScalar elemMats[], PetscContainer container)
1794 {
1795   PetscErrorCode  ierr;
1796   IS             *is;
1797   PetscInt        nc,colour,j;
1798   const PetscInt *clr_idxs;
1799   ISColoring      iscoloring;
1800   PetscFunctionBegin;
1801   ierr = PetscContainerGetPointer(container,(void**)&iscoloring);CHKERRQ(ierr);
1802   ierr = ISColoringGetIS(iscoloring,PETSC_USE_POINTER,&nc,&is);CHKERRQ(ierr);
1803   for (colour=0; colour<nc; colour++) {
1804     PetscInt    *idx_arr[1024]; /* need to make dynamic for general use */
1805     PetscScalar *new_el_mats[1024];
1806     PetscInt     idx_size[1024],csize;
1807     ierr = ISGetLocalSize(is[colour],&csize);CHKERRQ(ierr);
1808     if (csize>1024) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "too many elements in color. %D > %D",csize,1024);
1809     ierr = ISGetIndices(is[colour],&clr_idxs);CHKERRQ(ierr);
1810     /* get indices and mats */
1811     for (j=0; j<csize; j++) {
1812       PetscInt    cell = cStart + clr_idxs[j];
1813       PetscInt    numindices,*indices;
1814       PetscScalar *elMat = &elemMats[clr_idxs[j]*totDim*totDim];
1815       PetscScalar *valuesOrig = elMat;
1816       ierr = DMPlexGetClosureIndices(plex, section, globalSection, cell, PETSC_TRUE, &numindices, &indices, NULL, (PetscScalar **) &elMat);CHKERRQ(ierr);
1817       idx_size[j] = numindices;
1818       ierr = PetscMalloc2(numindices,&idx_arr[j],numindices*numindices,&new_el_mats[j]);CHKERRQ(ierr);
1819       ierr = PetscMemcpy(idx_arr[j],indices,numindices*sizeof(PetscInt));CHKERRQ(ierr);
1820       ierr = PetscMemcpy(new_el_mats[j],elMat,numindices*numindices*sizeof(PetscScalar));CHKERRQ(ierr);
1821       ierr = DMPlexRestoreClosureIndices(plex, section, globalSection, cell, PETSC_TRUE, &numindices, &indices, NULL, (PetscScalar **) &elMat);CHKERRQ(ierr);
1822       if (elMat != valuesOrig) {ierr = DMRestoreWorkArray(plex, numindices*numindices, MPIU_SCALAR, &elMat);}
1823     }
1824     /* assemble matrix - pragmas break CI ? */
1825     //#pragma omp parallel default(JacP,idx_size,idx_arr,new_el_mats,colour,clr_idxs)  private(j)
1826     //#pragma omp parallel for private(j)
1827     for (j=0; j<csize; j++) {
1828       PetscInt    numindices = idx_size[j], *indices = idx_arr[j];
1829       PetscScalar *elMat = new_el_mats[j];
1830       MatSetValues(JacP,numindices,indices,numindices,indices,elMat,ADD_VALUES);
1831     }
1832     /* free */
1833     ierr = ISRestoreIndices(is[colour],&clr_idxs);CHKERRQ(ierr);
1834     for (j=0; j<csize; j++) {
1835       ierr = PetscFree2(idx_arr[j],new_el_mats[j]);CHKERRQ(ierr);
1836     }
1837   }
1838   ierr = ISColoringRestoreIS(iscoloring,PETSC_USE_POINTER,&is);CHKERRQ(ierr);
1839   PetscFunctionReturn(0);
1840 }
1841 
1842 /* < v, u > */
1843 static void g0_1(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1844                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1845                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1846                  PetscReal t, PetscReal u_tShift, const PetscReal x[],  PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
1847 {
1848   g0[0] = 1.;
1849 }
1850 
1851 /* < v, u > */
1852 static void g0_r(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1853                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1854                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1855                  PetscReal t, PetscReal u_tShift, const PetscReal x[],  PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
1856 {
1857   g0[0] = 2.*PETSC_PI*x[0];
1858 }
1859 
1860 /*@
1861  LandauCreateMassMatrix - Create mass matrix for Landau
1862 
1863  Collective on dm
1864 
1865  Input Parameters:
1866  . dm     - the DM object
1867 
1868  Output Parameters:
1869  . Amat - The mass matrix (optional), mass matrix is added to the DM context
1870 
1871  Level: beginner
1872 
1873  .keywords: mesh
1874  .seealso: LandauCreateVelocitySpace()
1875  @*/
1876 PetscErrorCode LandauCreateMassMatrix(DM dm, Mat *Amat)
1877 {
1878   DM             massDM;
1879   PetscDS        prob;
1880   PetscInt       ii,dim,N1=1,N2;
1881   PetscErrorCode ierr;
1882   LandauCtx      *ctx;
1883   Mat            M;
1884 
1885   PetscFunctionBegin;
1886   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1887   if (Amat) PetscValidPointer(Amat,3);
1888   ierr = DMGetApplicationContext(dm, &ctx);CHKERRQ(ierr);
1889   if (!ctx) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "no context");
1890   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1891   ierr = DMClone(dm, &massDM);CHKERRQ(ierr);
1892   ierr = DMCopyFields(dm, massDM);CHKERRQ(ierr);
1893   ierr = DMCreateDS(massDM);CHKERRQ(ierr);
1894   ierr = DMGetDS(massDM, &prob);CHKERRQ(ierr);
1895   for (ii=0;ii<ctx->num_species;ii++) {
1896     if (dim==3) {ierr = PetscDSSetJacobian(prob, ii, ii, g0_1, NULL, NULL, NULL);CHKERRQ(ierr);}
1897     else        {ierr = PetscDSSetJacobian(prob, ii, ii, g0_r, NULL, NULL, NULL);CHKERRQ(ierr);}
1898   }
1899   ierr = DMViewFromOptions(massDM,NULL,"-dm_landau_mass_dm_view");CHKERRQ(ierr);
1900   ierr = DMCreateMatrix(massDM, &M);CHKERRQ(ierr);
1901   ierr = MatSetOption(M, MAT_IGNORE_ZERO_ENTRIES, PETSC_TRUE);CHKERRQ(ierr);
1902   {
1903     Vec locX;
1904     DM  plex;
1905     ierr = DMConvert(massDM, DMPLEX, &plex);CHKERRQ(ierr);
1906     ierr = DMGetLocalVector(massDM, &locX);CHKERRQ(ierr);
1907     /* Mass matrix is independent of the input, so no need to fill locX */
1908     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);
1909     ierr = DMPlexSNESComputeJacobianFEM(plex, locX, M, M, ctx);CHKERRQ(ierr);
1910     ierr = DMRestoreLocalVector(massDM, &locX);CHKERRQ(ierr);
1911     ierr = DMDestroy(&plex);CHKERRQ(ierr);
1912   }
1913   ierr = DMDestroy(&massDM);CHKERRQ(ierr);
1914   ierr = MatGetSize(ctx->J, &N1, NULL);CHKERRQ(ierr);
1915   ierr = MatGetSize(M, &N2, NULL);CHKERRQ(ierr);
1916   if (N1 != N2) SETERRQ2(PetscObjectComm((PetscObject) dm), PETSC_ERR_PLIB, "Incorrect matrix sizes: |Jacobian| = %D, |Mass|=%D",N1,N2);
1917   ierr = PetscObjectSetName((PetscObject)M, "mass");CHKERRQ(ierr);
1918   ierr = MatViewFromOptions(M,NULL,"-dm_landau_mass_mat_view");CHKERRQ(ierr);
1919   ctx->M = M; /* this could be a noop, a = a */
1920   if (Amat) *Amat = M;
1921   PetscFunctionReturn(0);
1922 }
1923 
1924 /*@
1925  LandauIFunction - TS residual calculation
1926 
1927  Collective on ts
1928 
1929  Input Parameters:
1930  +   TS  - The time stepping context
1931  .   time_dummy - current time (not used)
1932  -   X - Current state
1933  +   X_t - Time derivative of current state
1934  .   actx - Landau context
1935 
1936  Output Parameter:
1937  .   F  - The residual
1938 
1939  Level: beginner
1940 
1941  .keywords: mesh
1942  .seealso: LandauCreateVelocitySpace(), LandauIJacobian()
1943  @*/
1944 PetscErrorCode LandauIFunction(TS ts, PetscReal time_dummy, Vec X, Vec X_t, Vec F, void *actx)
1945 {
1946   PetscErrorCode ierr;
1947   LandauCtx      *ctx=(LandauCtx*)actx;
1948   PetscInt       dim;
1949   DM             dm;
1950 
1951   PetscFunctionBegin;
1952   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
1953   ierr = DMGetApplicationContext(dm, &ctx);CHKERRQ(ierr);
1954   if (!ctx) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "no context");
1955   ierr = PetscLogEventBegin(ctx->events[0],0,0,0,0);CHKERRQ(ierr);
1956   ierr = DMGetDimension(ctx->dmv, &dim);CHKERRQ(ierr);
1957   ierr = PetscInfo3(ts, "Create Landau Jacobian t=%g X'=%p %s\n",time_dummy,X_t,ctx->aux_bool ? " -- seems to be in line search" : "");CHKERRQ(ierr);
1958   ierr = LandauFormJacobian_Internal(X,ctx->J,dim,0.0,(void*)ctx);CHKERRQ(ierr);
1959   ctx->aux_bool = PETSC_TRUE;
1960   ierr = MatViewFromOptions(ctx->J,NULL,"-landau_jacobian_mat_view");CHKERRQ(ierr);
1961   /* mat vec for op */
1962   ierr = MatMult(ctx->J,X,F);CHKERRQ(ierr);CHKERRQ(ierr); /* C*f */
1963   /* add time term */
1964   if (X_t) {
1965     ierr = MatMultAdd(ctx->M,X_t,F,F);CHKERRQ(ierr);
1966   }
1967   ierr = PetscLogEventEnd(ctx->events[0],0,0,0,0);CHKERRQ(ierr);
1968   PetscFunctionReturn(0);
1969 }
1970 static PetscErrorCode MatrixNfDestroy(void *ptr)
1971 {
1972   PetscInt *nf = (PetscInt *)ptr;
1973   PetscErrorCode  ierr;
1974   PetscFunctionBegin;
1975   ierr = PetscFree(nf);CHKERRQ(ierr);
1976   PetscFunctionReturn(0);
1977 }
1978 /*@
1979  LandauIJacobian - TS Jacobian construction
1980 
1981  Collective on ts
1982 
1983  Input Parameters:
1984  +   TS  - The time stepping context
1985  .   time_dummy - current time (not used)
1986  -   X - Current state
1987  +   U_tdummy - Time derivative of current state (not used)
1988  .   shift - shift for du/dt term
1989  -   actx - Landau context
1990 
1991  Output Parameter:
1992  .   Amat  - Jacobian
1993  +   Pmat  - same as Amat
1994 
1995  Level: beginner
1996 
1997  .keywords: mesh
1998  .seealso: LandauCreateVelocitySpace(), LandauIFunction()
1999  @*/
2000 PetscErrorCode LandauIJacobian(TS ts, PetscReal time_dummy, Vec X, Vec U_tdummy, PetscReal shift, Mat Amat, Mat Pmat, void *actx)
2001 {
2002   PetscErrorCode ierr;
2003   LandauCtx      *ctx=(LandauCtx*)actx;
2004   PetscInt       dim;
2005   DM             dm;
2006   PetscContainer container;
2007   PetscFunctionBegin;
2008   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
2009   ierr = DMGetApplicationContext(dm, &ctx);CHKERRQ(ierr);
2010   if (!ctx) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "no context");
2011   if (Amat!=Pmat || Amat!=ctx->J) SETERRQ(ctx->comm, PETSC_ERR_PLIB, "Amat!=Pmat || Amat!=ctx->J");
2012   ierr = DMGetDimension(ctx->dmv, &dim);CHKERRQ(ierr);
2013   /* get collision Jacobian into A */
2014   ierr = PetscLogEventBegin(ctx->events[9],0,0,0,0);CHKERRQ(ierr);
2015   ierr = PetscInfo2(ts, "Adding mass to Jacobian t=%g, shift=%g\n",(double)time_dummy,(double)shift);CHKERRQ(ierr);
2016   if (shift==0.0) SETERRQ(ctx->comm, PETSC_ERR_PLIB, "zero shift");
2017   if (!ctx->aux_bool) SETERRQ(ctx->comm, PETSC_ERR_PLIB, "wrong state");
2018   ierr = LandauFormJacobian_Internal(X,ctx->J,dim,shift,(void*)ctx);CHKERRQ(ierr);
2019   ctx->aux_bool = PETSC_FALSE;
2020   ierr = MatViewFromOptions(Pmat,NULL,"-landau_mat_view");CHKERRQ(ierr);
2021   ierr = PetscLogEventEnd(ctx->events[9],0,0,0,0);CHKERRQ(ierr);
2022   /* set number species in Jacobian */
2023   ierr = PetscObjectQuery((PetscObject) ctx->J, "Nf", (PetscObject *) &container);CHKERRQ(ierr);
2024   if (!container) {
2025     PetscInt *pNf;
2026     ierr = PetscContainerCreate(PETSC_COMM_SELF, &container);CHKERRQ(ierr);
2027     ierr = PetscMalloc(sizeof(PetscInt), &pNf);CHKERRQ(ierr);
2028     *pNf = ctx->num_species + 1000*ctx->numConcurrency;
2029     ierr = PetscContainerSetPointer(container, (void *)pNf);CHKERRQ(ierr);
2030     ierr = PetscContainerSetUserDestroy(container, MatrixNfDestroy);CHKERRQ(ierr);
2031     ierr = PetscObjectCompose((PetscObject)ctx->J, "Nf", (PetscObject) container);CHKERRQ(ierr);
2032     ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
2033   }
2034 
2035   PetscFunctionReturn(0);
2036 }
2037