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