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