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