Actual source code: plexland.c

  1: #include <../src/mat/impls/aij/seq/aij.h>
  2: #include <petsc/private/dmpleximpl.h>
  3: #include <petsclandau.h>
  4: #include <petscts.h>
  5: #include <petscdmforest.h>
  6: #include <petscdmcomposite.h>

  8: /* Landau collision operator */

 10: /* relativistic terms */
 11: #if defined(PETSC_USE_REAL_SINGLE)
 12:   #define SPEED_OF_LIGHT 2.99792458e8F
 13:   #define C_0(v0)        (SPEED_OF_LIGHT / v0) /* needed for relativistic tensor on all architectures */
 14: #else
 15:   #define SPEED_OF_LIGHT 2.99792458e8
 16:   #define C_0(v0)        (SPEED_OF_LIGHT / v0) /* needed for relativistic tensor on all architectures */
 17: #endif

 19: #include "land_tensors.h"

 21: #if defined(PETSC_HAVE_OPENMP)
 22:   #include <omp.h>
 23: #endif

 25: static PetscErrorCode LandauGPUMapsDestroy(void **ptr)
 26: {
 27:   P4estVertexMaps *maps = (P4estVertexMaps *)*ptr;

 29:   PetscFunctionBegin;
 30:   // free device data
 31:   if (maps[0].deviceType != LANDAU_CPU) {
 32: #if defined(PETSC_HAVE_KOKKOS)
 33:     if (maps[0].deviceType == LANDAU_KOKKOS) PetscCall(LandauKokkosDestroyMatMaps(maps, maps[0].numgrids)); // implies Kokkos does
 34: #endif
 35:   }
 36:   // free host data
 37:   for (PetscInt grid = 0; grid < maps[0].numgrids; grid++) {
 38:     PetscCall(PetscFree(maps[grid].c_maps));
 39:     PetscCall(PetscFree(maps[grid].gIdx));
 40:   }
 41:   PetscCall(PetscFree(maps));
 42:   PetscFunctionReturn(PETSC_SUCCESS);
 43: }
 44: static PetscErrorCode energy_f(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf_dummy, PetscScalar *u, void *actx)
 45: {
 46:   PetscReal v2 = 0;

 48:   PetscFunctionBegin;
 49:   /* compute v^2 / 2 */
 50:   for (PetscInt i = 0; i < dim; ++i) v2 += x[i] * x[i];
 51:   /* evaluate the Maxwellian */
 52:   u[0] = v2 / 2;
 53:   PetscFunctionReturn(PETSC_SUCCESS);
 54: }

 56: /* needs double */
 57: static PetscErrorCode gamma_m1_f(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf_dummy, PetscScalar *u, void *actx)
 58: {
 59:   PetscReal *c2_0_arr = ((PetscReal *)actx);
 60:   double     u2 = 0, c02 = (double)*c2_0_arr, xx;

 62:   PetscFunctionBegin;
 63:   /* compute u^2 / 2 */
 64:   for (PetscInt i = 0; i < dim; ++i) u2 += x[i] * x[i];
 65:   /* gamma - 1 = g_eps, for conditioning and we only take derivatives */
 66:   xx = u2 / c02;
 67: #if defined(PETSC_USE_DEBUG)
 68:   u[0] = PetscSqrtReal(1. + xx);
 69: #else
 70:   u[0] = xx / (PetscSqrtReal(1. + xx) + 1.) - 1.; // better conditioned. -1 might help condition and only used for derivative
 71: #endif
 72:   PetscFunctionReturn(PETSC_SUCCESS);
 73: }

 75: /*
 76:  LandauFormJacobian_Internal - Evaluates Jacobian matrix.

 78:  Input Parameters:
 79:  .  globX - input vector
 80:  .  actx - optional user-defined context
 81:  .  dim - dimension

 83:  Output Parameter:
 84:  .  J0acP - Jacobian matrix filled, not created
 85:  */
 86: static PetscErrorCode LandauFormJacobian_Internal(Vec a_X, Mat JacP, const PetscInt dim, PetscReal shift, void *a_ctx)
 87: {
 88:   LandauCtx         *ctx = (LandauCtx *)a_ctx;
 89:   PetscInt           numCells[LANDAU_MAX_GRIDS], Nq, Nb;
 90:   PetscQuadrature    quad;
 91:   PetscReal          Eq_m[LANDAU_MAX_SPECIES]; // could be static data w/o quench (ex2)
 92:   PetscScalar       *cellClosure = NULL;
 93:   const PetscScalar *xdata       = NULL;
 94:   PetscDS            prob;
 95:   PetscContainer     container;
 96:   P4estVertexMaps   *maps;
 97:   Mat                subJ[LANDAU_MAX_GRIDS * LANDAU_MAX_BATCH_SZ];

 99:   PetscFunctionBegin;
102:   PetscAssertPointer(ctx, 5);
103:   /* check for matrix container for GPU assembly. Support CPU assembly for debugging */
104:   PetscCheck(ctx->plex[0] != NULL, ctx->comm, PETSC_ERR_ARG_WRONG, "Plex not created");
105:   PetscCall(PetscLogEventBegin(ctx->events[10], 0, 0, 0, 0));
106:   PetscCall(DMGetDS(ctx->plex[0], &prob)); // same DS for all grids
107:   PetscCall(PetscObjectQuery((PetscObject)JacP, "assembly_maps", (PetscObject *)&container));
108:   if (container) {
109:     PetscCheck(ctx->gpu_assembly, ctx->comm, PETSC_ERR_ARG_WRONG, "maps but no GPU assembly");
110:     PetscCall(PetscContainerGetPointer(container, (void **)&maps));
111:     PetscCheck(maps, ctx->comm, PETSC_ERR_ARG_WRONG, "empty GPU matrix container");
112:     for (PetscInt i = 0; i < ctx->num_grids * ctx->batch_sz; i++) subJ[i] = NULL;
113:   } else {
114:     PetscCheck(!ctx->gpu_assembly, ctx->comm, PETSC_ERR_ARG_WRONG, "No maps but GPU assembly");
115:     for (PetscInt tid = 0; tid < ctx->batch_sz; tid++) {
116:       for (PetscInt grid = 0; grid < ctx->num_grids; grid++) PetscCall(DMCreateMatrix(ctx->plex[grid], &subJ[LAND_PACK_IDX(tid, grid)]));
117:     }
118:     maps = NULL;
119:   }
120:   // get dynamic data (Eq is odd, for quench and Spitzer test) for CPU assembly and raw data for Jacobian GPU assembly. Get host numCells[], Nq (yuck)
121:   PetscCall(PetscFEGetQuadrature(ctx->fe[0], &quad));
122:   PetscCall(PetscQuadratureGetData(quad, NULL, NULL, &Nq, NULL, NULL));
123:   PetscCall(PetscFEGetDimension(ctx->fe[0], &Nb));
124:   PetscCheck(Nq <= LANDAU_MAX_NQND, ctx->comm, PETSC_ERR_ARG_WRONG, "Order too high. Nq = %" PetscInt_FMT " > LANDAU_MAX_NQND (%d)", Nq, LANDAU_MAX_NQND);
125:   PetscCheck(Nb <= LANDAU_MAX_NQND, ctx->comm, PETSC_ERR_ARG_WRONG, "Order too high. Nb = %" PetscInt_FMT " > LANDAU_MAX_NQND (%d)", Nb, LANDAU_MAX_NQND);
126:   // get metadata for collecting dynamic data
127:   for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
128:     PetscInt cStart, cEnd;
129:     PetscCheck(ctx->plex[grid] != NULL, ctx->comm, PETSC_ERR_ARG_WRONG, "Plex not created");
130:     PetscCall(DMPlexGetHeightStratum(ctx->plex[grid], 0, &cStart, &cEnd));
131:     numCells[grid] = cEnd - cStart; // grids can have different topology
132:   }
133:   PetscCall(PetscLogEventEnd(ctx->events[10], 0, 0, 0, 0));
134:   if (shift == 0) { /* create dynamic point data: f_alpha for closure of each cell (cellClosure[nbatch,ngrids,ncells[g],f[Nb,ns[g]]]) or xdata */
135:     DM pack;
136:     PetscCall(VecGetDM(a_X, &pack));
137:     PetscCheck(pack, PETSC_COMM_SELF, PETSC_ERR_PLIB, "pack has no DM");
138:     PetscCall(PetscLogEventBegin(ctx->events[1], 0, 0, 0, 0));
139:     for (PetscInt fieldA = 0; fieldA < ctx->num_species; fieldA++) {
140:       Eq_m[fieldA] = ctx->Ez * ctx->t_0 * ctx->charges[fieldA] / (ctx->v_0 * ctx->masses[fieldA]); /* normalize dimensionless */
141:       if (dim == 2) Eq_m[fieldA] *= 2 * PETSC_PI;                                                  /* add the 2pi term that is not in Landau */
142:     }
143:     if (!ctx->gpu_assembly) {
144:       Vec         *locXArray, *globXArray;
145:       PetscScalar *cellClosure_it;
146:       PetscInt     cellClosure_sz = 0, nDMs, Nf[LANDAU_MAX_GRIDS];
147:       PetscSection section[LANDAU_MAX_GRIDS], globsection[LANDAU_MAX_GRIDS];
148:       for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
149:         PetscCall(DMGetLocalSection(ctx->plex[grid], &section[grid]));
150:         PetscCall(DMGetGlobalSection(ctx->plex[grid], &globsection[grid]));
151:         PetscCall(PetscSectionGetNumFields(section[grid], &Nf[grid]));
152:       }
153:       /* count cellClosure size */
154:       PetscCall(DMCompositeGetNumberDM(pack, &nDMs));
155:       for (PetscInt grid = 0; grid < ctx->num_grids; grid++) cellClosure_sz += Nb * Nf[grid] * numCells[grid];
156:       PetscCall(PetscMalloc1(cellClosure_sz * ctx->batch_sz, &cellClosure));
157:       cellClosure_it = cellClosure;
158:       PetscCall(PetscMalloc(sizeof(*locXArray) * nDMs, &locXArray));
159:       PetscCall(PetscMalloc(sizeof(*globXArray) * nDMs, &globXArray));
160:       PetscCall(DMCompositeGetLocalAccessArray(pack, a_X, nDMs, NULL, locXArray));
161:       PetscCall(DMCompositeGetAccessArray(pack, a_X, nDMs, NULL, globXArray));
162:       for (PetscInt b_id = 0; b_id < ctx->batch_sz; b_id++) { // OpenMP (once)
163:         for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
164:           Vec      locX = locXArray[LAND_PACK_IDX(b_id, grid)], globX = globXArray[LAND_PACK_IDX(b_id, grid)], locX2;
165:           PetscInt cStart, cEnd, ei;
166:           PetscCall(VecDuplicate(locX, &locX2));
167:           PetscCall(DMGlobalToLocalBegin(ctx->plex[grid], globX, INSERT_VALUES, locX2));
168:           PetscCall(DMGlobalToLocalEnd(ctx->plex[grid], globX, INSERT_VALUES, locX2));
169:           PetscCall(DMPlexGetHeightStratum(ctx->plex[grid], 0, &cStart, &cEnd));
170:           for (ei = cStart; ei < cEnd; ++ei) {
171:             PetscScalar *coef = NULL;
172:             PetscCall(DMPlexVecGetClosure(ctx->plex[grid], section[grid], locX2, ei, NULL, &coef));
173:             PetscCall(PetscMemcpy(cellClosure_it, coef, Nb * Nf[grid] * sizeof(*cellClosure_it))); /* change if LandauIPReal != PetscScalar */
174:             PetscCall(DMPlexVecRestoreClosure(ctx->plex[grid], section[grid], locX2, ei, NULL, &coef));
175:             cellClosure_it += Nb * Nf[grid];
176:           }
177:           PetscCall(VecDestroy(&locX2));
178:         }
179:       }
180:       PetscCheck(cellClosure_it - cellClosure == cellClosure_sz * ctx->batch_sz, PETSC_COMM_SELF, PETSC_ERR_PLIB, "iteration wrong %" PetscCount_FMT " != cellClosure_sz = %" PetscInt_FMT, cellClosure_it - cellClosure, cellClosure_sz * ctx->batch_sz);
181:       PetscCall(DMCompositeRestoreLocalAccessArray(pack, a_X, nDMs, NULL, locXArray));
182:       PetscCall(DMCompositeRestoreAccessArray(pack, a_X, nDMs, NULL, globXArray));
183:       PetscCall(PetscFree(locXArray));
184:       PetscCall(PetscFree(globXArray));
185:       xdata = NULL;
186:     } else {
187:       PetscMemType mtype;
188:       if (ctx->jacobian_field_major_order) { // get data in batch ordering
189:         PetscCall(VecScatterBegin(ctx->plex_batch, a_X, ctx->work_vec, INSERT_VALUES, SCATTER_FORWARD));
190:         PetscCall(VecScatterEnd(ctx->plex_batch, a_X, ctx->work_vec, INSERT_VALUES, SCATTER_FORWARD));
191:         PetscCall(VecGetArrayReadAndMemType(ctx->work_vec, &xdata, &mtype));
192:       } else {
193:         PetscCall(VecGetArrayReadAndMemType(a_X, &xdata, &mtype));
194:       }
195:       PetscCheck(mtype == PETSC_MEMTYPE_HOST || ctx->deviceType != LANDAU_CPU, ctx->comm, PETSC_ERR_ARG_WRONG, "CPU run with device data: use -mat_type aij");
196:       cellClosure = NULL;
197:     }
198:     PetscCall(PetscLogEventEnd(ctx->events[1], 0, 0, 0, 0));
199:   } else xdata = cellClosure = NULL;

201:   /* do it */
202:   if (ctx->deviceType == LANDAU_KOKKOS) {
203: #if defined(PETSC_HAVE_KOKKOS)
204:     PetscCall(LandauKokkosJacobian(ctx->plex, Nq, Nb, ctx->batch_sz, ctx->num_grids, numCells, Eq_m, cellClosure, xdata, &ctx->SData_d, shift, ctx->events, ctx->mat_offset, ctx->species_offset, subJ, JacP));
205: #else
206:     SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONG, "-landau_device_type %s not built", "kokkos");
207: #endif
208:   } else {               /* CPU version */
209:     PetscTabulation *Tf; // used for CPU and print info. Same on all grids and all species
210:     PetscInt         ip_offset[LANDAU_MAX_GRIDS + 1], ipf_offset[LANDAU_MAX_GRIDS + 1], elem_offset[LANDAU_MAX_GRIDS + 1], IPf_sz_glb, IPf_sz_tot, num_grids = ctx->num_grids, Nf[LANDAU_MAX_GRIDS];
211:     PetscReal       *ff, *dudx, *dudy, *dudz, *invJ_a = (PetscReal *)ctx->SData_d.invJ, *xx = (PetscReal *)ctx->SData_d.x, *yy = (PetscReal *)ctx->SData_d.y, *zz = (PetscReal *)ctx->SData_d.z, *ww = (PetscReal *)ctx->SData_d.w;
212:     PetscReal       *nu_alpha = (PetscReal *)ctx->SData_d.alpha, *nu_beta = (PetscReal *)ctx->SData_d.beta, *invMass = (PetscReal *)ctx->SData_d.invMass;
213:     PetscReal (*lambdas)[LANDAU_MAX_GRIDS][LANDAU_MAX_GRIDS] = (PetscReal (*)[LANDAU_MAX_GRIDS][LANDAU_MAX_GRIDS])ctx->SData_d.lambdas;
214:     PetscSection section[LANDAU_MAX_GRIDS], globsection[LANDAU_MAX_GRIDS];
215:     PetscScalar *coo_vals = NULL;
216:     for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
217:       PetscCall(DMGetLocalSection(ctx->plex[grid], &section[grid]));
218:       PetscCall(DMGetGlobalSection(ctx->plex[grid], &globsection[grid]));
219:       PetscCall(PetscSectionGetNumFields(section[grid], &Nf[grid]));
220:     }
221:     /* count IPf size, etc */
222:     PetscCall(PetscDSGetTabulation(prob, &Tf)); // Bf, &Df same for all grids
223:     const PetscReal *const BB = Tf[0]->T[0], *const DD = Tf[0]->T[1];
224:     ip_offset[0] = ipf_offset[0] = elem_offset[0] = 0;
225:     for (PetscInt grid = 0; grid < num_grids; grid++) {
226:       PetscInt nfloc        = ctx->species_offset[grid + 1] - ctx->species_offset[grid];
227:       elem_offset[grid + 1] = elem_offset[grid] + numCells[grid];
228:       ip_offset[grid + 1]   = ip_offset[grid] + numCells[grid] * Nq;
229:       ipf_offset[grid + 1]  = ipf_offset[grid] + Nq * nfloc * numCells[grid];
230:     }
231:     IPf_sz_glb = ipf_offset[num_grids];
232:     IPf_sz_tot = IPf_sz_glb * ctx->batch_sz;
233:     // prep COO
234:     PetscCall(PetscMalloc1(ctx->SData_d.coo_size, &coo_vals)); // allocate every time?
235:     if (shift == 0.0) {                                        /* compute dynamic data f and df and init data for Jacobian */
236: #if defined(PETSC_HAVE_THREADSAFETY)
237:       double starttime, endtime;
238:       starttime = MPI_Wtime();
239: #endif
240:       PetscCall(PetscLogEventBegin(ctx->events[8], 0, 0, 0, 0));
241:       PetscCall(PetscMalloc4(IPf_sz_tot, &ff, IPf_sz_tot, &dudx, IPf_sz_tot, &dudy, (dim == 3 ? IPf_sz_tot : 0), &dudz));
242:       // F df/dx
243:       for (PetscInt tid = 0; tid < ctx->batch_sz * elem_offset[num_grids]; tid++) {                        // for each element
244:         const PetscInt b_Nelem = elem_offset[num_grids], b_elem_idx = tid % b_Nelem, b_id = tid / b_Nelem; // b_id == OMP thd_id in batch
245:         // find my grid:
246:         PetscInt grid = 0;
247:         while (b_elem_idx >= elem_offset[grid + 1]) grid++; // yuck search for grid
248:         {
249:           const PetscInt loc_nip = numCells[grid] * Nq, loc_Nf = ctx->species_offset[grid + 1] - ctx->species_offset[grid], loc_elem = b_elem_idx - elem_offset[grid];
250:           const PetscInt moffset = LAND_MOFFSET(b_id, grid, ctx->batch_sz, ctx->num_grids, ctx->mat_offset); //b_id*b_N + ctx->mat_offset[grid];
251:           PetscScalar   *coef, coef_buff[LANDAU_MAX_SPECIES * LANDAU_MAX_NQND];
252:           PetscReal     *invJe = &invJ_a[(ip_offset[grid] + loc_elem * Nq) * dim * dim]; // ingJ is static data on batch 0
253:           PetscInt       b, f, q;
254:           if (cellClosure) {
255:             coef = &cellClosure[b_id * IPf_sz_glb + ipf_offset[grid] + loc_elem * Nb * loc_Nf]; // this is const
256:           } else {
257:             coef = coef_buff;
258:             for (f = 0; f < loc_Nf; ++f) {
259:               LandauIdx *const Idxs = &maps[grid].gIdx[loc_elem][f][0];
260:               for (b = 0; b < Nb; ++b) {
261:                 PetscInt idx = Idxs[b];
262:                 if (idx >= 0) {
263:                   coef[f * Nb + b] = xdata[idx + moffset];
264:                 } else {
265:                   idx              = -idx - 1;
266:                   coef[f * Nb + b] = 0;
267:                   for (q = 0; q < maps[grid].num_face; q++) {
268:                     PetscInt    id    = maps[grid].c_maps[idx][q].gid;
269:                     PetscScalar scale = maps[grid].c_maps[idx][q].scale;
270:                     coef[f * Nb + b] += scale * xdata[id + moffset];
271:                   }
272:                 }
273:               }
274:             }
275:           }
276:           /* get f and df */
277:           for (PetscInt qi = 0; qi < Nq; qi++) {
278:             const PetscReal *invJ = &invJe[qi * dim * dim];
279:             const PetscReal *Bq   = &BB[qi * Nb];
280:             const PetscReal *Dq   = &DD[qi * Nb * dim];
281:             PetscReal        u_x[LANDAU_DIM];
282:             /* get f & df */
283:             for (f = 0; f < loc_Nf; ++f) {
284:               const PetscInt idx = b_id * IPf_sz_glb + ipf_offset[grid] + f * loc_nip + loc_elem * Nq + qi;
285:               PetscInt       b, e;
286:               PetscReal      refSpaceDer[LANDAU_DIM];
287:               ff[idx] = 0.0;
288:               for (PetscInt d = 0; d < LANDAU_DIM; ++d) refSpaceDer[d] = 0.0;
289:               for (b = 0; b < Nb; ++b) {
290:                 const PetscInt cidx = b;
291:                 ff[idx] += Bq[cidx] * PetscRealPart(coef[f * Nb + cidx]);
292:                 for (PetscInt d = 0; d < dim; ++d) refSpaceDer[d] += Dq[cidx * dim + d] * PetscRealPart(coef[f * Nb + cidx]);
293:               }
294:               for (PetscInt d = 0; d < LANDAU_DIM; ++d) {
295:                 for (e = 0, u_x[d] = 0.0; e < LANDAU_DIM; ++e) u_x[d] += invJ[e * dim + d] * refSpaceDer[e];
296:               }
297:               dudx[idx] = u_x[0];
298:               dudy[idx] = u_x[1];
299: #if LANDAU_DIM == 3
300:               dudz[idx] = u_x[2];
301: #endif
302:             }
303:           } // q
304:         } // grid
305:       } // grid*batch
306:       PetscCall(PetscLogEventEnd(ctx->events[8], 0, 0, 0, 0));
307: #if defined(PETSC_HAVE_THREADSAFETY)
308:       endtime = MPI_Wtime();
309:       if (ctx->stage) ctx->times[LANDAU_F_DF] += (endtime - starttime);
310: #endif
311:     } // Jacobian setup
312:     // assemble Jacobian (or mass)
313:     for (PetscInt tid = 0; tid < ctx->batch_sz * elem_offset[num_grids]; tid++) { // for each element
314:       const PetscInt b_Nelem      = elem_offset[num_grids];
315:       const PetscInt glb_elem_idx = tid % b_Nelem, b_id = tid / b_Nelem;
316:       PetscInt       grid = 0;
317: #if defined(PETSC_HAVE_THREADSAFETY)
318:       double starttime, endtime;
319:       starttime = MPI_Wtime();
320: #endif
321:       while (glb_elem_idx >= elem_offset[grid + 1]) grid++;
322:       {
323:         const PetscInt   loc_Nf = ctx->species_offset[grid + 1] - ctx->species_offset[grid], loc_elem = glb_elem_idx - elem_offset[grid];
324:         const PetscInt   moffset = LAND_MOFFSET(b_id, grid, ctx->batch_sz, ctx->num_grids, ctx->mat_offset), totDim = loc_Nf * Nq, elemMatSize = totDim * totDim;
325:         PetscScalar     *elemMat;
326:         const PetscReal *invJe = &invJ_a[(ip_offset[grid] + loc_elem * Nq) * dim * dim];
327:         PetscCall(PetscMalloc1(elemMatSize, &elemMat));
328:         PetscCall(PetscMemzero(elemMat, elemMatSize * sizeof(*elemMat)));
329:         if (shift == 0.0) { // Jacobian
330:           PetscCall(PetscLogEventBegin(ctx->events[4], 0, 0, 0, 0));
331:         } else { // mass
332:           PetscCall(PetscLogEventBegin(ctx->events[16], 0, 0, 0, 0));
333:         }
334:         for (PetscInt qj = 0; qj < Nq; ++qj) {
335:           const PetscInt jpidx_glb = ip_offset[grid] + qj + loc_elem * Nq;
336:           PetscReal      g0[LANDAU_MAX_SPECIES], g2[LANDAU_MAX_SPECIES][LANDAU_DIM], g3[LANDAU_MAX_SPECIES][LANDAU_DIM][LANDAU_DIM]; // could make a LANDAU_MAX_SPECIES_GRID ~ number of ions - 1
337:           PetscInt       d, d2, dp, d3, IPf_idx;
338:           if (shift == 0.0) { // Jacobian
339:             const PetscReal *const invJj = &invJe[qj * dim * dim];
340:             PetscReal              gg2[LANDAU_MAX_SPECIES][LANDAU_DIM], gg3[LANDAU_MAX_SPECIES][LANDAU_DIM][LANDAU_DIM], gg2_temp[LANDAU_DIM], gg3_temp[LANDAU_DIM][LANDAU_DIM];
341:             const PetscReal        vj[3] = {xx[jpidx_glb], yy[jpidx_glb], zz ? zz[jpidx_glb] : 0}, wj = ww[jpidx_glb];
342:             // create g2 & g3
343:             for (d = 0; d < LANDAU_DIM; d++) { // clear accumulation data D & K
344:               gg2_temp[d] = 0;
345:               for (d2 = 0; d2 < LANDAU_DIM; d2++) gg3_temp[d][d2] = 0;
346:             }
347:             /* inner beta reduction */
348:             IPf_idx = 0;
349:             for (PetscInt grid_r = 0, f_off = 0, ipidx = 0; grid_r < ctx->num_grids; grid_r++, f_off = ctx->species_offset[grid_r]) { // IPf_idx += nip_loc_r*Nfloc_r
350:               PetscInt nip_loc_r = numCells[grid_r] * Nq, Nfloc_r = Nf[grid_r];
351:               for (PetscInt ei_r = 0; ei_r < numCells[grid_r]; ++ei_r) {
352:                 for (PetscInt qi = 0; qi < Nq; qi++, ipidx++) {
353:                   const PetscReal wi = ww[ipidx], x = xx[ipidx], y = yy[ipidx];
354:                   PetscReal       temp1[3] = {0, 0, 0}, temp2 = 0;
355: #if LANDAU_DIM == 2
356:                   PetscReal Ud[2][2], Uk[2][2], mask = (PetscAbs(vj[0] - x) < 100 * PETSC_SQRT_MACHINE_EPSILON && PetscAbs(vj[1] - y) < 100 * PETSC_SQRT_MACHINE_EPSILON) ? 0. : 1.;
357:                   LandauTensor2D(vj, x, y, Ud, Uk, mask);
358: #else
359:                   PetscReal U[3][3], z = zz[ipidx], mask = (PetscAbs(vj[0] - x) < 100 * PETSC_SQRT_MACHINE_EPSILON && PetscAbs(vj[1] - y) < 100 * PETSC_SQRT_MACHINE_EPSILON && PetscAbs(vj[2] - z) < 100 * PETSC_SQRT_MACHINE_EPSILON) ? 0. : 1.;
360:                   if (ctx->use_relativistic_corrections) {
361:                     LandauTensor3DRelativistic(vj, x, y, z, U, mask, C_0(ctx->v_0));
362:                   } else {
363:                     LandauTensor3D(vj, x, y, z, U, mask);
364:                   }
365: #endif
366:                   for (PetscInt f = 0; f < Nfloc_r; ++f) {
367:                     const PetscInt idx = b_id * IPf_sz_glb + ipf_offset[grid_r] + f * nip_loc_r + ei_r * Nq + qi;

369:                     temp1[0] += dudx[idx] * nu_beta[f + f_off] * invMass[f + f_off] * (*lambdas)[grid][grid_r];
370:                     temp1[1] += dudy[idx] * nu_beta[f + f_off] * invMass[f + f_off] * (*lambdas)[grid][grid_r];
371: #if LANDAU_DIM == 3
372:                     temp1[2] += dudz[idx] * nu_beta[f + f_off] * invMass[f + f_off] * (*lambdas)[grid][grid_r];
373: #endif
374:                     temp2 += ff[idx] * nu_beta[f + f_off] * (*lambdas)[grid][grid_r];
375:                   }
376:                   temp1[0] *= wi;
377:                   temp1[1] *= wi;
378: #if LANDAU_DIM == 3
379:                   temp1[2] *= wi;
380: #endif
381:                   temp2 *= wi;
382: #if LANDAU_DIM == 2
383:                   for (d2 = 0; d2 < 2; d2++) {
384:                     for (d3 = 0; d3 < 2; ++d3) {
385:                       /* K = U * grad(f): g2=e: i,A */
386:                       gg2_temp[d2] += Uk[d2][d3] * temp1[d3];
387:                       /* D = -U * (I \kron (fx)): g3=f: i,j,A */
388:                       gg3_temp[d2][d3] += Ud[d2][d3] * temp2;
389:                     }
390:                   }
391: #else
392:                   for (d2 = 0; d2 < 3; ++d2) {
393:                     for (d3 = 0; d3 < 3; ++d3) {
394:                       /* K = U * grad(f): g2 = e: i,A */
395:                       gg2_temp[d2] += U[d2][d3] * temp1[d3];
396:                       /* D = -U * (I \kron (fx)): g3 = f: i,j,A */
397:                       gg3_temp[d2][d3] += U[d2][d3] * temp2;
398:                     }
399:                   }
400: #endif
401:                 } // qi
402:               } // ei_r
403:               IPf_idx += nip_loc_r * Nfloc_r;
404:             } /* grid_r - IPs */
405:             PetscCheck(IPf_idx == IPf_sz_glb, PETSC_COMM_SELF, PETSC_ERR_PLIB, "IPf_idx != IPf_sz %" PetscInt_FMT " %" PetscInt_FMT, IPf_idx, IPf_sz_glb);
406:             // add alpha and put in gg2/3
407:             for (PetscInt fieldA = 0, f_off = ctx->species_offset[grid]; fieldA < loc_Nf; ++fieldA) {
408:               for (d2 = 0; d2 < LANDAU_DIM; d2++) {
409:                 gg2[fieldA][d2] = gg2_temp[d2] * nu_alpha[fieldA + f_off];
410:                 for (d3 = 0; d3 < LANDAU_DIM; d3++) gg3[fieldA][d2][d3] = -gg3_temp[d2][d3] * nu_alpha[fieldA + f_off] * invMass[fieldA + f_off];
411:               }
412:             }
413:             /* add electric field term once per IP */
414:             for (PetscInt fieldA = 0, f_off = ctx->species_offset[grid]; fieldA < loc_Nf; ++fieldA) gg2[fieldA][LANDAU_DIM - 1] += Eq_m[fieldA + f_off];
415:             /* Jacobian transform - g2, g3 */
416:             for (PetscInt fieldA = 0; fieldA < loc_Nf; ++fieldA) {
417:               for (d = 0; d < dim; ++d) {
418:                 g2[fieldA][d] = 0.0;
419:                 for (d2 = 0; d2 < dim; ++d2) {
420:                   g2[fieldA][d] += invJj[d * dim + d2] * gg2[fieldA][d2];
421:                   g3[fieldA][d][d2] = 0.0;
422:                   for (d3 = 0; d3 < dim; ++d3) {
423:                     for (dp = 0; dp < dim; ++dp) g3[fieldA][d][d2] += invJj[d * dim + d3] * gg3[fieldA][d3][dp] * invJj[d2 * dim + dp];
424:                   }
425:                   g3[fieldA][d][d2] *= wj;
426:                 }
427:                 g2[fieldA][d] *= wj;
428:               }
429:             }
430:           } else { // mass
431:             PetscReal wj = ww[jpidx_glb];
432:             /* Jacobian transform - g0 */
433:             for (PetscInt fieldA = 0; fieldA < loc_Nf; ++fieldA) {
434:               if (dim == 2) {
435:                 g0[fieldA] = wj * shift * 2. * PETSC_PI; // move this to below and remove g0
436:               } else {
437:                 g0[fieldA] = wj * shift; // move this to below and remove g0
438:               }
439:             }
440:           }
441:           /* FE matrix construction */
442:           {
443:             PetscInt         fieldA, d, f, d2, g;
444:             const PetscReal *BJq = &BB[qj * Nb], *DIq = &DD[qj * Nb * dim];
445:             /* assemble - on the diagonal (I,I) */
446:             for (fieldA = 0; fieldA < loc_Nf; fieldA++) {
447:               for (f = 0; f < Nb; f++) {
448:                 const PetscInt i = fieldA * Nb + f; /* Element matrix row */
449:                 for (g = 0; g < Nb; ++g) {
450:                   const PetscInt j    = fieldA * Nb + g; /* Element matrix column */
451:                   const PetscInt fOff = i * totDim + j;
452:                   if (shift == 0.0) {
453:                     for (d = 0; d < dim; ++d) {
454:                       elemMat[fOff] += DIq[f * dim + d] * g2[fieldA][d] * BJq[g];
455:                       for (d2 = 0; d2 < dim; ++d2) elemMat[fOff] += DIq[f * dim + d] * g3[fieldA][d][d2] * DIq[g * dim + d2];
456:                     }
457:                   } else { // mass
458:                     elemMat[fOff] += BJq[f] * g0[fieldA] * BJq[g];
459:                   }
460:                 }
461:               }
462:             }
463:           }
464:         } /* qj loop */
465:         if (shift == 0.0) { // Jacobian
466:           PetscCall(PetscLogEventEnd(ctx->events[4], 0, 0, 0, 0));
467:         } else {
468:           PetscCall(PetscLogEventEnd(ctx->events[16], 0, 0, 0, 0));
469:         }
470: #if defined(PETSC_HAVE_THREADSAFETY)
471:         endtime = MPI_Wtime();
472:         if (ctx->stage) ctx->times[LANDAU_KERNEL] += (endtime - starttime);
473: #endif
474:         /* assemble matrix */
475:         if (!container) {
476:           PetscInt cStart;
477:           PetscCall(PetscLogEventBegin(ctx->events[6], 0, 0, 0, 0));
478:           PetscCall(DMPlexGetHeightStratum(ctx->plex[grid], 0, &cStart, NULL));
479:           PetscCall(DMPlexMatSetClosure(ctx->plex[grid], section[grid], globsection[grid], subJ[LAND_PACK_IDX(b_id, grid)], loc_elem + cStart, elemMat, ADD_VALUES));
480:           PetscCall(PetscLogEventEnd(ctx->events[6], 0, 0, 0, 0));
481:         } else { // GPU like assembly for debugging
482:           PetscInt    fieldA, q, f, g, d, nr, nc, rows0[LANDAU_MAX_Q_FACE] = {0}, cols0[LANDAU_MAX_Q_FACE] = {0}, rows[LANDAU_MAX_Q_FACE], cols[LANDAU_MAX_Q_FACE];
483:           PetscScalar vals[LANDAU_MAX_Q_FACE * LANDAU_MAX_Q_FACE] = {0}, row_scale[LANDAU_MAX_Q_FACE] = {0}, col_scale[LANDAU_MAX_Q_FACE] = {0};
484:           LandauIdx *coo_elem_offsets = (LandauIdx *)ctx->SData_d.coo_elem_offsets, *coo_elem_fullNb = (LandauIdx *)ctx->SData_d.coo_elem_fullNb, (*coo_elem_point_offsets)[LANDAU_MAX_NQND + 1] = (LandauIdx(*)[LANDAU_MAX_NQND + 1]) ctx->SData_d.coo_elem_point_offsets;
485:           /* assemble - from the diagonal (I,I) in this format for DMPlexMatSetClosure */
486:           for (fieldA = 0; fieldA < loc_Nf; fieldA++) {
487:             LandauIdx *const Idxs = &maps[grid].gIdx[loc_elem][fieldA][0];
488:             for (f = 0; f < Nb; f++) {
489:               PetscInt idx = Idxs[f];
490:               if (idx >= 0) {
491:                 nr           = 1;
492:                 rows0[0]     = idx;
493:                 row_scale[0] = 1.;
494:               } else {
495:                 idx = -idx - 1;
496:                 for (q = 0, nr = 0; q < maps[grid].num_face; q++, nr++) {
497:                   if (maps[grid].c_maps[idx][q].gid < 0) break;
498:                   rows0[q]     = maps[grid].c_maps[idx][q].gid;
499:                   row_scale[q] = maps[grid].c_maps[idx][q].scale;
500:                 }
501:               }
502:               for (g = 0; g < Nb; ++g) {
503:                 idx = Idxs[g];
504:                 if (idx >= 0) {
505:                   nc           = 1;
506:                   cols0[0]     = idx;
507:                   col_scale[0] = 1.;
508:                 } else {
509:                   idx = -idx - 1;
510:                   nc  = maps[grid].num_face;
511:                   for (q = 0, nc = 0; q < maps[grid].num_face; q++, nc++) {
512:                     if (maps[grid].c_maps[idx][q].gid < 0) break;
513:                     cols0[q]     = maps[grid].c_maps[idx][q].gid;
514:                     col_scale[q] = maps[grid].c_maps[idx][q].scale;
515:                   }
516:                 }
517:                 const PetscInt    i   = fieldA * Nb + f; /* Element matrix row */
518:                 const PetscInt    j   = fieldA * Nb + g; /* Element matrix column */
519:                 const PetscScalar Aij = elemMat[i * totDim + j];
520:                 if (coo_vals) { // mirror (i,j) in CreateStaticGPUData
521:                   const PetscInt fullNb = coo_elem_fullNb[glb_elem_idx], fullNb2 = fullNb * fullNb;
522:                   const PetscInt idx0 = b_id * coo_elem_offsets[elem_offset[num_grids]] + coo_elem_offsets[glb_elem_idx] + fieldA * fullNb2 + fullNb * coo_elem_point_offsets[glb_elem_idx][f] + nr * coo_elem_point_offsets[glb_elem_idx][g];
523:                   for (PetscInt q = 0, idx2 = idx0; q < nr; q++) {
524:                     for (PetscInt d = 0; d < nc; d++, idx2++) coo_vals[idx2] = row_scale[q] * col_scale[d] * Aij;
525:                   }
526:                 } else {
527:                   for (q = 0; q < nr; q++) rows[q] = rows0[q] + moffset;
528:                   for (d = 0; d < nc; d++) cols[d] = cols0[d] + moffset;
529:                   for (q = 0; q < nr; q++) {
530:                     for (d = 0; d < nc; d++) vals[q * nc + d] = row_scale[q] * col_scale[d] * Aij;
531:                   }
532:                   PetscCall(MatSetValues(JacP, nr, rows, nc, cols, vals, ADD_VALUES));
533:                 }
534:               }
535:             }
536:           }
537:         }
538:         if (loc_elem == -1) {
539:           PetscCall(PetscPrintf(ctx->comm, "CPU Element matrix\n"));
540:           for (PetscInt d = 0; d < totDim; ++d) {
541:             for (PetscInt f = 0; f < totDim; ++f) PetscCall(PetscPrintf(ctx->comm, " %12.5e", (double)PetscRealPart(elemMat[d * totDim + f])));
542:             PetscCall(PetscPrintf(ctx->comm, "\n"));
543:           }
544:           exit(12);
545:         }
546:         PetscCall(PetscFree(elemMat));
547:       } /* grid */
548:     } /* outer element & batch loop */
549:     if (shift == 0.0) { // mass
550:       PetscCall(PetscFree4(ff, dudx, dudy, dudz));
551:     }
552:     if (!container) {                                         // 'CPU' assembly move nest matrix to global JacP
553:       for (PetscInt b_id = 0; b_id < ctx->batch_sz; b_id++) { // OpenMP
554:         for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
555:           const PetscInt     moffset = LAND_MOFFSET(b_id, grid, ctx->batch_sz, ctx->num_grids, ctx->mat_offset); // b_id*b_N + ctx->mat_offset[grid];
556:           PetscInt           nloc, nzl, colbuf[1024], row;
557:           const PetscInt    *cols;
558:           const PetscScalar *vals;
559:           Mat                B = subJ[LAND_PACK_IDX(b_id, grid)];
560:           PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
561:           PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
562:           PetscCall(MatGetSize(B, &nloc, NULL));
563:           for (PetscInt i = 0; i < nloc; i++) {
564:             PetscCall(MatGetRow(B, i, &nzl, &cols, &vals));
565:             PetscCheck(nzl <= 1024, PetscObjectComm((PetscObject)B), PETSC_ERR_PLIB, "Row too big: %" PetscInt_FMT, nzl);
566:             for (PetscInt j = 0; j < nzl; j++) colbuf[j] = moffset + cols[j];
567:             row = moffset + i;
568:             PetscCall(MatSetValues(JacP, 1, &row, nzl, colbuf, vals, ADD_VALUES));
569:             PetscCall(MatRestoreRow(B, i, &nzl, &cols, &vals));
570:           }
571:           PetscCall(MatDestroy(&B));
572:         }
573:       }
574:     }
575:     if (coo_vals) {
576:       PetscCall(MatSetValuesCOO(JacP, coo_vals, ADD_VALUES));
577:       PetscCall(PetscFree(coo_vals));
578:     }
579:   } /* CPU version */
580:   PetscCall(MatAssemblyBegin(JacP, MAT_FINAL_ASSEMBLY));
581:   PetscCall(MatAssemblyEnd(JacP, MAT_FINAL_ASSEMBLY));
582:   /* clean up */
583:   if (cellClosure) PetscCall(PetscFree(cellClosure));
584:   if (xdata) PetscCall(VecRestoreArrayReadAndMemType(a_X, &xdata));
585:   PetscFunctionReturn(PETSC_SUCCESS);
586: }

588: static PetscErrorCode GeometryDMLandau(DM base, PetscInt point, PetscInt dim, const PetscReal abc[], PetscReal xyz[], void *a_ctx)
589: {
590:   PetscReal  r = abc[0], z = abc[1];
591:   LandauCtx *ctx = (LandauCtx *)a_ctx;

593:   PetscFunctionBegin;
594:   if (ctx->sphere && dim == 3) { // make sphere: works for one AMR and Q2
595:     PetscInt nzero = 0, idx = 0;
596:     xyz[0] = r;
597:     xyz[1] = z;
598:     xyz[2] = abc[2];
599:     for (PetscInt i = 0; i < 3; i++) {
600:       if (PetscAbs(xyz[i]) < PETSC_SQRT_MACHINE_EPSILON) nzero++;
601:       else idx = i;
602:     }
603:     if (nzero == 2) xyz[idx] *= 1.732050807568877; // sqrt(3)
604:     else if (nzero == 1) {
605:       for (PetscInt i = 0; i < 3; i++) xyz[i] *= 1.224744871391589; // sqrt(3/2)
606:     }
607:   } else {
608:     xyz[0] = r;
609:     xyz[1] = z;
610:     if (dim == 3) xyz[2] = abc[2];
611:   }
612:   PetscFunctionReturn(PETSC_SUCCESS);
613: }

615: /* create DMComposite of meshes for each species group */
616: static PetscErrorCode LandauDMCreateVMeshes(MPI_Comm comm_self, const PetscInt dim, const char prefix[], LandauCtx *ctx, DM pack)
617: {
618:   PetscFunctionBegin;
619:   { /* p4est, quads */
620:     /* Create plex mesh of Landau domain */
621:     for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
622:       PetscReal par_radius = ctx->radius_par[grid], perp_radius = ctx->radius_perp[grid];
623:       if (!ctx->sphere && !ctx->simplex) { // 2 or 3D (only 3D option)
624:         PetscReal      lo[] = {-perp_radius, -par_radius, -par_radius}, hi[] = {perp_radius, par_radius, par_radius};
625:         DMBoundaryType periodicity[3] = {DM_BOUNDARY_NONE, dim == 2 ? DM_BOUNDARY_NONE : DM_BOUNDARY_NONE, DM_BOUNDARY_NONE};
626:         if (dim == 2) lo[0] = 0;
627:         else {
628:           lo[1] = -perp_radius;
629:           hi[1] = perp_radius; // 3D y is a perp
630:         }
631:         PetscCall(DMPlexCreateBoxMesh(comm_self, dim, PETSC_FALSE, ctx->cells0, lo, hi, periodicity, PETSC_TRUE, 0, PETSC_TRUE, &ctx->plex[grid])); // TODO: make composite and create dm[grid] here
632:         PetscCall(DMLocalizeCoordinates(ctx->plex[grid]));                                                                                          /* needed for periodic */
633:         if (dim == 3) PetscCall(PetscObjectSetName((PetscObject)ctx->plex[grid], "cube"));
634:         else PetscCall(PetscObjectSetName((PetscObject)ctx->plex[grid], "half-plane"));
635:       } else if (dim == 2) {
636:         size_t len;
637:         PetscCall(PetscStrlen(ctx->filename, &len));
638:         if (len) {
639:           Vec          coords;
640:           PetscScalar *x;
641:           PetscInt     N;
642:           char         str[] = "-dm_landau_view_file_0";
643:           str[21] += grid;
644:           PetscCall(DMPlexCreateFromFile(comm_self, ctx->filename, "plexland.c", PETSC_TRUE, &ctx->plex[grid]));
645:           PetscCall(DMPlexOrient(ctx->plex[grid]));
646:           PetscCall(DMGetCoordinatesLocal(ctx->plex[grid], &coords));
647:           PetscCall(VecGetSize(coords, &N));
648:           PetscCall(VecGetArray(coords, &x));
649:           /* scale by domain size */
650:           for (PetscInt i = 0; i < N; i += 2) {
651:             x[i + 0] *= ctx->radius_perp[grid];
652:             x[i + 1] *= ctx->radius_par[grid];
653:           }
654:           PetscCall(VecRestoreArray(coords, &x));
655:           PetscCall(PetscObjectSetName((PetscObject)ctx->plex[grid], ctx->filename));
656:           PetscCall(PetscInfo(ctx->plex[grid], "%" PetscInt_FMT ") Read %s mesh file (%s)\n", grid, ctx->filename, str));
657:           PetscCall(DMViewFromOptions(ctx->plex[grid], NULL, str));
658:         } else { // simplex forces a sphere
659:           PetscInt       numCells = ctx->simplex ? 12 : 6, cell_size = ctx->simplex ? 3 : 4, j;
660:           const PetscInt numVerts    = 11;
661:           PetscInt       cellsT[][4] = {
662:             {0,  1, 6, 5 },
663:             {1,  2, 7, 6 },
664:             {2,  3, 8, 7 },
665:             {3,  4, 9, 8 },
666:             {5,  6, 7, 10},
667:             {10, 7, 8, 9 }
668:           };
669:           PetscInt cellsS[][3] = {
670:             {0,  1, 6 },
671:             {1,  2, 6 },
672:             {6,  2, 7 },
673:             {7,  2, 8 },
674:             {8,  2, 3 },
675:             {8,  3, 4 },
676:             {0,  6, 5 },
677:             {5,  6, 7 },
678:             {5,  7, 10},
679:             {10, 7, 9 },
680:             {9,  7, 8 },
681:             {9,  8, 4 }
682:           };
683:           const PetscInt *pcell = (const PetscInt *)(ctx->simplex ? &cellsS[0][0] : &cellsT[0][0]);
684:           PetscReal       coords[11][2], *flatCoords = &coords[0][0];
685:           PetscReal       rad = ctx->radius[grid];
686:           for (j = 0; j < 5; j++) { // outside edge
687:             PetscReal z, r, theta = -PETSC_PI / 2 + (j % 5) * PETSC_PI / 4;
688:             r            = rad * PetscCosReal(theta);
689:             coords[j][0] = r;
690:             z            = rad * PetscSinReal(theta);
691:             coords[j][1] = z;
692:           }
693:           coords[j][0]   = 0;
694:           coords[j++][1] = -rad * ctx->sphere_inner_radius_90degree[grid];
695:           coords[j][0]   = rad * ctx->sphere_inner_radius_45degree[grid] * 0.707106781186548;
696:           coords[j++][1] = -rad * ctx->sphere_inner_radius_45degree[grid] * 0.707106781186548;
697:           coords[j][0]   = rad * ctx->sphere_inner_radius_90degree[grid];
698:           coords[j++][1] = 0;
699:           coords[j][0]   = rad * ctx->sphere_inner_radius_45degree[grid] * 0.707106781186548;
700:           coords[j++][1] = rad * ctx->sphere_inner_radius_45degree[grid] * 0.707106781186548;
701:           coords[j][0]   = 0;
702:           coords[j++][1] = rad * ctx->sphere_inner_radius_90degree[grid];
703:           coords[j][0]   = 0;
704:           coords[j++][1] = 0;
705:           PetscCall(DMPlexCreateFromCellListPetsc(comm_self, 2, numCells, numVerts, cell_size, ctx->interpolate, pcell, 2, flatCoords, &ctx->plex[grid]));
706:           PetscCall(PetscObjectSetName((PetscObject)ctx->plex[grid], "semi-circle"));
707:           PetscCall(PetscInfo(ctx->plex[grid], "\t%" PetscInt_FMT ") Make circle %s mesh\n", grid, ctx->simplex ? "simplex" : "tensor"));
708:         }
709:       } else {
710:         PetscCheck(dim == 3 && ctx->sphere && !ctx->simplex, ctx->comm, PETSC_ERR_ARG_WRONG, "not: dim == 3 && ctx->sphere && !ctx->simplex");
711:         PetscReal      rad = ctx->radius[grid] / 1.732050807568877, inner_rad = rad * ctx->sphere_inner_radius_45degree[grid], outer_rad = rad;
712:         const PetscInt numCells = 7, cell_size = 8, numVerts = 16;
713:         const PetscInt cells[][8] = {
714:           {0, 3, 2, 1, 4,  5,  6,  7 },
715:           {0, 4, 5, 1, 8,  9,  13, 12},
716:           {1, 5, 6, 2, 9,  10, 14, 13},
717:           {2, 6, 7, 3, 10, 11, 15, 14},
718:           {0, 3, 7, 4, 8,  12, 15, 11},
719:           {0, 1, 2, 3, 8,  11, 10, 9 },
720:           {4, 7, 6, 5, 12, 13, 14, 15}
721:         };
722:         PetscReal coords[16 /* numVerts */][3];
723:         for (PetscInt j = 0; j < 4; j++) { // inner edge, low
724:           coords[j][0] = inner_rad * (j == 0 || j == 3 ? 1 : -1);
725:           coords[j][1] = inner_rad * (j / 2 < 1 ? 1 : -1);
726:           coords[j][2] = inner_rad * -1;
727:         }
728:         for (PetscInt j = 0, jj = 4; j < 4; j++, jj++) { // inner edge, hi
729:           coords[jj][0] = inner_rad * (j == 0 || j == 3 ? 1 : -1);
730:           coords[jj][1] = inner_rad * (j / 2 < 1 ? 1 : -1);
731:           coords[jj][2] = inner_rad * 1;
732:         }
733:         for (PetscInt j = 0, jj = 8; j < 4; j++, jj++) { // outer edge, low
734:           coords[jj][0] = outer_rad * (j == 0 || j == 3 ? 1 : -1);
735:           coords[jj][1] = outer_rad * (j / 2 < 1 ? 1 : -1);
736:           coords[jj][2] = outer_rad * -1;
737:         }
738:         for (PetscInt j = 0, jj = 12; j < 4; j++, jj++) { // outer edge, hi
739:           coords[jj][0] = outer_rad * (j == 0 || j == 3 ? 1 : -1);
740:           coords[jj][1] = outer_rad * (j / 2 < 1 ? 1 : -1);
741:           coords[jj][2] = outer_rad * 1;
742:         }
743:         PetscCall(DMPlexCreateFromCellListPetsc(comm_self, 3, numCells, numVerts, cell_size, ctx->interpolate, (const PetscInt *)cells, 3, (const PetscReal *)coords, &ctx->plex[grid]));
744:         PetscCall(PetscObjectSetName((PetscObject)ctx->plex[grid], "cubed sphere"));
745:         PetscCall(PetscInfo(ctx->plex[grid], "\t%" PetscInt_FMT ") Make cubed sphere %s mesh\n", grid, ctx->simplex ? "simplex" : "tensor"));
746:       }
747:       PetscCall(DMSetFromOptions(ctx->plex[grid]));
748:     } // grid loop
749:     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)pack, prefix));
750:     { /* convert to p4est (or whatever), wait for discretization to create pack */
751:       char      convType[256];
752:       PetscBool flg;

754:       PetscOptionsBegin(ctx->comm, prefix, "Mesh conversion options", "DMPLEX");
755:       PetscCall(PetscOptionsFList("-dm_landau_type", "Convert DMPlex to another format (p4est)", "plexland.c", DMList, DMPLEX, convType, 256, &flg));
756:       PetscOptionsEnd();
757:       if (flg) {
758:         ctx->use_p4est = PETSC_TRUE; /* flag for Forest */
759:         for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
760:           DM        dmforest;
761:           PetscBool isForest;

763:           PetscCall(DMConvert(ctx->plex[grid], convType, &dmforest));
764:           PetscCheck(dmforest, ctx->comm, PETSC_ERR_PLIB, "Convert failed?");
765:           PetscCall(PetscObjectSetOptionsPrefix((PetscObject)dmforest, prefix));
766:           PetscCall(DMIsForest(dmforest, &isForest));
767:           PetscCheck(isForest, ctx->comm, PETSC_ERR_PLIB, "Converted to non Forest?");
768:           if (ctx->sphere) PetscCall(DMForestSetBaseCoordinateMapping(dmforest, GeometryDMLandau, ctx));
769:           PetscCall(DMDestroy(&ctx->plex[grid]));
770:           ctx->plex[grid] = dmforest; // Forest for adaptivity
771:         }
772:       } else ctx->use_p4est = PETSC_FALSE; /* flag for Forest */
773:     }
774:   } /* non-file */
775:   PetscCall(DMSetDimension(pack, dim));
776:   PetscCall(PetscObjectSetName((PetscObject)pack, "Mesh"));
777:   PetscCall(DMSetApplicationContext(pack, ctx));
778:   PetscFunctionReturn(PETSC_SUCCESS);
779: }

781: static PetscErrorCode SetupDS(DM pack, PetscInt dim, PetscInt grid, LandauCtx *ctx)
782: {
783:   PetscInt     ii, i0;
784:   char         buf[256];
785:   PetscSection section;

787:   PetscFunctionBegin;
788:   for (ii = ctx->species_offset[grid], i0 = 0; ii < ctx->species_offset[grid + 1]; ii++, i0++) {
789:     if (ii == 0) PetscCall(PetscSNPrintf(buf, sizeof(buf), "e"));
790:     else PetscCall(PetscSNPrintf(buf, sizeof(buf), "i%" PetscInt_FMT, ii));
791:     /* Setup Discretization - FEM */
792:     PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, ctx->simplex, NULL, PETSC_DECIDE, &ctx->fe[ii]));
793:     PetscCall(PetscObjectSetName((PetscObject)ctx->fe[ii], buf));
794:     PetscCall(DMSetField(ctx->plex[grid], i0, NULL, (PetscObject)ctx->fe[ii]));
795:   }
796:   PetscCall(DMCreateDS(ctx->plex[grid]));
797:   PetscCall(DMGetLocalSection(ctx->plex[grid], &section));
798:   for (PetscInt ii = ctx->species_offset[grid], i0 = 0; ii < ctx->species_offset[grid + 1]; ii++, i0++) {
799:     if (ii == 0) PetscCall(PetscSNPrintf(buf, sizeof(buf), "se"));
800:     else PetscCall(PetscSNPrintf(buf, sizeof(buf), "si%" PetscInt_FMT, ii));
801:     PetscCall(PetscSectionSetComponentName(section, i0, 0, buf));
802:   }
803:   PetscFunctionReturn(PETSC_SUCCESS);
804: }

806: /* Define a Maxwellian function for testing out the operator. */

808: /* Using cartesian velocity space coordinates, the particle */
809: /* density, [1/m^3], is defined according to */

811: /* $$ n=\int_{R^3} dv^3 \left(\frac{m}{2\pi T}\right)^{3/2}\exp [- mv^2/(2T)] $$ */

813: /* Using some constant, c, we normalize the velocity vector into a */
814: /* dimensionless variable according to v=c*x. Thus the density, $n$, becomes */

816: /* $$ n=\int_{R^3} dx^3 \left(\frac{mc^2}{2\pi T}\right)^{3/2}\exp [- mc^2/(2T)*x^2] $$ */

818: /* Defining $\theta=2T/mc^2$, we thus find that the probability density */
819: /* for finding the particle within the interval in a box dx^3 around x is */

821: /* f(x;\theta)=\left(\frac{1}{\pi\theta}\right)^{3/2} \exp [ -x^2/\theta ] */

823: typedef struct {
824:   PetscReal v_0;
825:   PetscReal kT_m;
826:   PetscReal n;
827:   PetscReal shift;
828: } MaxwellianCtx;

830: static PetscErrorCode maxwellian(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf_dummy, PetscScalar *u, void *actx)
831: {
832:   MaxwellianCtx *mctx = (MaxwellianCtx *)actx;
833:   PetscInt       i;
834:   PetscReal      v2 = 0, theta = 2 * mctx->kT_m / (mctx->v_0 * mctx->v_0), shift; /* theta = 2kT/mc^2 */

836:   PetscFunctionBegin;
837:   /* compute the exponents, v^2 */
838:   for (i = 0; i < dim; ++i) v2 += x[i] * x[i];
839:   /* evaluate the Maxwellian */
840:   if (mctx->shift < 0) shift = -mctx->shift;
841:   else {
842:     u[0]  = mctx->n * PetscPowReal(PETSC_PI * theta, -1.5) * (PetscExpReal(-v2 / theta));
843:     shift = mctx->shift;
844:   }
845:   if (shift != 0.) {
846:     v2 = 0;
847:     for (i = 0; i < dim - 1; ++i) v2 += x[i] * x[i];
848:     v2 += (x[dim - 1] - shift) * (x[dim - 1] - shift);
849:     /* evaluate the shifted Maxwellian */
850:     u[0] += mctx->n * PetscPowReal(PETSC_PI * theta, -1.5) * (PetscExpReal(-v2 / theta));
851:   }
852:   PetscFunctionReturn(PETSC_SUCCESS);
853: }

855: /*@
856:   DMPlexLandauAddMaxwellians - Add a Maxwellian distribution to a state

858:   Collective

860:   Input Parameters:
861: + dm      - The mesh (local)
862: . time    - Current time
863: . temps   - Temperatures of each species (global)
864: . ns      - Number density of each species (global)
865: . grid    - index into current grid - just used for offset into `temp` and `ns`
866: . b_id    - batch index
867: . n_batch - number of batches
868: - actx    - Landau context

870:   Output Parameter:
871: . X - The state (local to this grid)

873:   Level: beginner

875: .seealso: `DMPlexLandauCreateVelocitySpace()`
876:  @*/
877: PetscErrorCode DMPlexLandauAddMaxwellians(DM dm, Vec X, PetscReal time, PetscReal temps[], PetscReal ns[], PetscInt grid, PetscInt b_id, PetscInt n_batch, void *actx)
878: {
879:   LandauCtx *ctx = (LandauCtx *)actx;
880:   PetscErrorCode (*initu[LANDAU_MAX_SPECIES])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *);
881:   PetscInt       dim;
882:   MaxwellianCtx *mctxs[LANDAU_MAX_SPECIES], data[LANDAU_MAX_SPECIES];

884:   PetscFunctionBegin;
885:   PetscCall(DMGetDimension(dm, &dim));
886:   if (!ctx) PetscCall(DMGetApplicationContext(dm, &ctx));
887:   for (PetscInt ii = ctx->species_offset[grid], i0 = 0; ii < ctx->species_offset[grid + 1]; ii++, i0++) {
888:     mctxs[i0]      = &data[i0];
889:     data[i0].v_0   = ctx->v_0;                             // v_0 same for all grids
890:     data[i0].kT_m  = ctx->k * temps[ii] / ctx->masses[ii]; /* kT/m */
891:     data[i0].n     = ns[ii];
892:     initu[i0]      = maxwellian;
893:     data[i0].shift = 0;
894:   }
895:   data[0].shift = ctx->electronShift;
896:   /* need to make ADD_ALL_VALUES work - TODO */
897:   PetscCall(DMProjectFunction(dm, time, initu, (void **)mctxs, INSERT_ALL_VALUES, X));
898:   PetscFunctionReturn(PETSC_SUCCESS);
899: }

901: /*
902:  LandauSetInitialCondition - Adds Maxwellians with context

904:  Collective

906:  Input Parameters:
907:  .   dm - The mesh
908:  -   grid - index into current grid - just used for offset into temp and ns
909:  .   b_id - batch index
910:  -   n_batch - number of batches
911:  +   actx - Landau context with T and n

913:  Output Parameter:
914:  .   X  - The state

916:  Level: beginner

918: .seealso: `DMPlexLandauCreateVelocitySpace()`, `DMPlexLandauAddMaxwellians()`
919:  */
920: static PetscErrorCode LandauSetInitialCondition(DM dm, Vec X, PetscInt grid, PetscInt b_id, PetscInt n_batch, void *actx)
921: {
922:   LandauCtx *ctx = (LandauCtx *)actx;

924:   PetscFunctionBegin;
925:   if (!ctx) PetscCall(DMGetApplicationContext(dm, &ctx));
926:   PetscCall(VecZeroEntries(X));
927:   PetscCall(DMPlexLandauAddMaxwellians(dm, X, 0.0, ctx->thermal_temps, ctx->n, grid, b_id, n_batch, ctx));
928:   PetscFunctionReturn(PETSC_SUCCESS);
929: }

931: // adapt a level once. Forest in/out
932: #if defined(PETSC_USE_INFO)
933: static const char *s_refine_names[] = {"RE", "Z1", "Origin", "Z2", "Uniform"};
934: #endif
935: static PetscErrorCode adaptToleranceFEM(PetscFE fem, Vec sol, PetscInt type, PetscInt grid, LandauCtx *ctx, DM *newForest)
936: {
937:   DM              forest, plex, adaptedDM = NULL;
938:   PetscDS         prob;
939:   PetscBool       isForest;
940:   PetscQuadrature quad;
941:   PetscInt        Nq, Nb, *Nb2, cStart, cEnd, c, dim, qj, k;
942:   DMLabel         adaptLabel = NULL;

944:   PetscFunctionBegin;
945:   forest = ctx->plex[grid];
946:   PetscCall(DMCreateDS(forest));
947:   PetscCall(DMGetDS(forest, &prob));
948:   PetscCall(DMGetDimension(forest, &dim));
949:   PetscCall(DMIsForest(forest, &isForest));
950:   PetscCheck(isForest, ctx->comm, PETSC_ERR_ARG_WRONG, "! Forest");
951:   PetscCall(DMConvert(forest, DMPLEX, &plex));
952:   PetscCall(DMPlexGetHeightStratum(plex, 0, &cStart, &cEnd));
953:   PetscCall(DMLabelCreate(PETSC_COMM_SELF, "adapt", &adaptLabel));
954:   PetscCall(PetscFEGetQuadrature(fem, &quad));
955:   PetscCall(PetscQuadratureGetData(quad, NULL, NULL, &Nq, NULL, NULL));
956:   PetscCheck(Nq <= LANDAU_MAX_NQND, ctx->comm, PETSC_ERR_ARG_WRONG, "Order too high. Nq = %" PetscInt_FMT " > LANDAU_MAX_NQND (%d)", Nq, LANDAU_MAX_NQND);
957:   PetscCall(PetscFEGetDimension(ctx->fe[0], &Nb));
958:   PetscCall(PetscDSGetDimensions(prob, &Nb2));
959:   PetscCheck(Nb2[0] == Nb, ctx->comm, PETSC_ERR_ARG_WRONG, " Nb = %" PetscInt_FMT " != Nb (%" PetscInt_FMT ")", Nb, Nb2[0]);
960:   PetscCheck(Nb <= LANDAU_MAX_NQND, ctx->comm, PETSC_ERR_ARG_WRONG, "Order too high. Nb = %" PetscInt_FMT " > LANDAU_MAX_NQND (%d)", Nb, LANDAU_MAX_NQND);
961:   PetscCall(PetscInfo(sol, "%" PetscInt_FMT ") Refine phase: %s\n", grid, s_refine_names[type]));
962:   if (type == 4) {
963:     for (c = cStart; c < cEnd; c++) PetscCall(DMLabelSetValue(adaptLabel, c, DM_ADAPT_REFINE));
964:   } else if (type == 2) {
965:     PetscInt  rCellIdx[8], nr = 0, nrmax = (dim == 3) ? 8 : 2;
966:     PetscReal minRad = PETSC_INFINITY, r;
967:     for (c = cStart; c < cEnd; c++) {
968:       PetscReal tt, v0[LANDAU_MAX_NQND * 3], J[LANDAU_MAX_NQND * 9], invJ[LANDAU_MAX_NQND * 9], detJ[LANDAU_MAX_NQND];
969:       PetscCall(DMPlexComputeCellGeometryFEM(plex, c, quad, v0, J, invJ, detJ));
970:       (void)J;
971:       (void)invJ;
972:       for (qj = 0; qj < Nq; ++qj) {
973:         tt = PetscSqr(v0[dim * qj + 0]) + PetscSqr(v0[dim * qj + 1]) + PetscSqr((dim == 3) ? v0[dim * qj + 2] : 0);
974:         r  = PetscSqrtReal(tt);
975:         if (r < minRad - PETSC_SQRT_MACHINE_EPSILON * 10.) {
976:           minRad         = r;
977:           nr             = 0;
978:           rCellIdx[nr++] = c;
979:           PetscCall(PetscInfo(sol, "\t\t%" PetscInt_FMT ") Found first inner r=%e, cell %" PetscInt_FMT ", qp %" PetscInt_FMT "/%" PetscInt_FMT "\n", grid, (double)r, c, qj + 1, Nq));
980:         } else if ((r - minRad) < PETSC_SQRT_MACHINE_EPSILON * 100. && nr < nrmax) {
981:           for (k = 0; k < nr; k++)
982:             if (c == rCellIdx[k]) break;
983:           if (k == nr) {
984:             rCellIdx[nr++] = c;
985:             PetscCall(PetscInfo(sol, "\t\t\t%" PetscInt_FMT ") Found another inner r=%e, cell %" PetscInt_FMT ", qp %" PetscInt_FMT "/%" PetscInt_FMT ", d=%e\n", grid, (double)r, c, qj + 1, Nq, (double)(r - minRad)));
986:           }
987:         }
988:       }
989:     }
990:     for (k = 0; k < nr; k++) PetscCall(DMLabelSetValue(adaptLabel, rCellIdx[k], DM_ADAPT_REFINE));
991:     PetscCall(PetscInfo(sol, "\t\t\t%" PetscInt_FMT ") Refined %" PetscInt_FMT " origin cells %" PetscInt_FMT ",%" PetscInt_FMT " r=%g\n", grid, nr, rCellIdx[0], rCellIdx[1], (double)minRad));
992:   } else if (type == 0 || type == 1 || type == 3) { /* refine along r=0 axis */
993:     PetscScalar *coef = NULL;
994:     Vec          coords;
995:     PetscInt     csize, Nv, d, nz, nrefined = 0;
996:     DM           cdm;
997:     PetscSection cs;
998:     PetscCall(DMGetCoordinatesLocal(forest, &coords));
999:     PetscCall(DMGetCoordinateDM(forest, &cdm));
1000:     PetscCall(DMGetLocalSection(cdm, &cs));
1001:     for (c = cStart; c < cEnd; c++) {
1002:       PetscInt doit = 0, outside = 0;
1003:       PetscCall(DMPlexVecGetClosure(cdm, cs, coords, c, &csize, &coef));
1004:       Nv = csize / dim;
1005:       for (nz = d = 0; d < Nv; d++) {
1006:         PetscReal z = PetscRealPart(coef[d * dim + (dim - 1)]), x = PetscSqr(PetscRealPart(coef[d * dim + 0])) + ((dim == 3) ? PetscSqr(PetscRealPart(coef[d * dim + 1])) : 0);
1007:         x = PetscSqrtReal(x);
1008:         if (type == 0) {
1009:           if (ctx->re_radius > PETSC_SQRT_MACHINE_EPSILON && (z < -PETSC_MACHINE_EPSILON * 10. || z > ctx->re_radius + PETSC_MACHINE_EPSILON * 10.)) outside++; /* first pass don't refine bottom */
1010:         } else if (type == 1 && (z > ctx->vperp0_radius1 || z < -ctx->vperp0_radius1)) {
1011:           outside++; /* don't refine outside electron refine radius */
1012:           PetscCall(PetscInfo(sol, "\t%" PetscInt_FMT ") (debug) found %s cells\n", grid, s_refine_names[type]));
1013:         } else if (type == 3 && (z > ctx->vperp0_radius2 || z < -ctx->vperp0_radius2)) {
1014:           outside++; /* refine r=0 cells on refinement front */
1015:           PetscCall(PetscInfo(sol, "\t%" PetscInt_FMT ") (debug) found %s cells\n", grid, s_refine_names[type]));
1016:         }
1017:         if (x < PETSC_MACHINE_EPSILON * 10. && (type != 0 || ctx->re_radius > PETSC_SQRT_MACHINE_EPSILON)) nz++;
1018:       }
1019:       PetscCall(DMPlexVecRestoreClosure(cdm, cs, coords, c, &csize, &coef));
1020:       if (doit || (outside < Nv && nz)) {
1021:         PetscCall(DMLabelSetValue(adaptLabel, c, DM_ADAPT_REFINE));
1022:         nrefined++;
1023:       }
1024:     }
1025:     PetscCall(PetscInfo(sol, "\t%" PetscInt_FMT ") Refined %" PetscInt_FMT " cells\n", grid, nrefined));
1026:   }
1027:   PetscCall(DMDestroy(&plex));
1028:   PetscCall(DMAdaptLabel(forest, adaptLabel, &adaptedDM));
1029:   PetscCall(DMLabelDestroy(&adaptLabel));
1030:   *newForest = adaptedDM;
1031:   if (adaptedDM) {
1032:     if (isForest) PetscCall(DMForestSetAdaptivityForest(adaptedDM, NULL)); // ????
1033:     PetscCall(DMConvert(adaptedDM, DMPLEX, &plex));
1034:     PetscCall(DMPlexGetHeightStratum(plex, 0, &cStart, &cEnd));
1035:     PetscCall(PetscInfo(sol, "\t\t\t\t%" PetscInt_FMT ") %" PetscInt_FMT " cells, %" PetscInt_FMT " total quadrature points\n", grid, cEnd - cStart, Nq * (cEnd - cStart)));
1036:     PetscCall(DMDestroy(&plex));
1037:   } else *newForest = NULL;
1038:   PetscFunctionReturn(PETSC_SUCCESS);
1039: }

1041: // forest goes in (ctx->plex[grid]), plex comes out
1042: static PetscErrorCode adapt(PetscInt grid, LandauCtx *ctx, Vec *uu)
1043: {
1044:   PetscInt adaptIter;

1046:   PetscFunctionBegin;
1047:   PetscInt type, limits[5] = {(grid == 0) ? ctx->numRERefine : 0, (grid == 0) ? ctx->nZRefine1 : 0, ctx->numAMRRefine[grid], (grid == 0) ? ctx->nZRefine2 : 0, ctx->postAMRRefine[grid]};
1048:   for (type = 0; type < 5; type++) {
1049:     for (adaptIter = 0; adaptIter < limits[type]; adaptIter++) {
1050:       DM newForest = NULL;
1051:       PetscCall(adaptToleranceFEM(ctx->fe[0], *uu, type, grid, ctx, &newForest));
1052:       if (newForest) {
1053:         PetscCall(DMDestroy(&ctx->plex[grid]));
1054:         PetscCall(VecDestroy(uu));
1055:         PetscCall(DMCreateGlobalVector(newForest, uu));
1056:         PetscCall(LandauSetInitialCondition(newForest, *uu, grid, 0, 1, ctx));
1057:         ctx->plex[grid] = newForest;
1058:       } else {
1059:         PetscCall(PetscInfo(*uu, "No refinement\n"));
1060:       }
1061:     }
1062:   }
1063:   PetscCall(PetscObjectSetName((PetscObject)*uu, "uAMR"));
1064:   PetscFunctionReturn(PETSC_SUCCESS);
1065: }

1067: // make log(Lambdas) from NRL Plasma formulary
1068: static PetscErrorCode makeLambdas(LandauCtx *ctx)
1069: {
1070:   PetscFunctionBegin;
1071:   for (PetscInt gridi = 0; gridi < ctx->num_grids; gridi++) {
1072:     PetscInt  iii   = ctx->species_offset[gridi];
1073:     PetscReal Ti_ev = (ctx->thermal_temps[iii] / 1.1604525e7) * 1000; // convert (back) to eV
1074:     PetscReal ni    = ctx->n[iii] * ctx->n_0;
1075:     for (PetscInt gridj = gridi; gridj < ctx->num_grids; gridj++) {
1076:       PetscInt  jjj = ctx->species_offset[gridj];
1077:       PetscReal Zj  = ctx->charges[jjj] / 1.6022e-19;
1078:       if (gridi == 0) {
1079:         if (gridj == 0) { // lam_ee
1080:           ctx->lambdas[gridi][gridj] = 23.5 - PetscLogReal(PetscSqrtReal(ni) * PetscPowReal(Ti_ev, -1.25)) - PetscSqrtReal(1e-5 + PetscSqr(PetscLogReal(Ti_ev) - 2) / 16);
1081:         } else { // lam_ei == lam_ie
1082:           if (10 * Zj * Zj > Ti_ev) {
1083:             ctx->lambdas[gridi][gridj] = ctx->lambdas[gridj][gridi] = 23 - PetscLogReal(PetscSqrtReal(ni) * Zj * PetscPowReal(Ti_ev, -1.5));
1084:           } else {
1085:             ctx->lambdas[gridi][gridj] = ctx->lambdas[gridj][gridi] = 24 - PetscLogReal(PetscSqrtReal(ni) / Ti_ev);
1086:           }
1087:         }
1088:       } else { // lam_ii'
1089:         PetscReal mui = ctx->masses[iii] / 1.6720e-27, Zi = ctx->charges[iii] / 1.6022e-19;
1090:         PetscReal Tj_ev            = (ctx->thermal_temps[jjj] / 1.1604525e7) * 1000; // convert (back) to eV
1091:         PetscReal muj              = ctx->masses[jjj] / 1.6720e-27;
1092:         PetscReal nj               = ctx->n[jjj] * ctx->n_0;
1093:         ctx->lambdas[gridi][gridj] = ctx->lambdas[gridj][gridi] = 23 - PetscLogReal(Zi * Zj * (mui + muj) / (mui * Tj_ev + muj * Ti_ev) * PetscSqrtReal(ni * Zi * Zi / Ti_ev + nj * Zj * Zj / Tj_ev));
1094:       }
1095:     }
1096:   }
1097:   //PetscReal v0 = PetscSqrtReal(ctx->k * ctx->thermal_temps[iii] / ctx->masses[iii]); /* arbitrary units for non-dimensionalization: plasma formulary def */
1098:   PetscFunctionReturn(PETSC_SUCCESS);
1099: }

1101: static PetscErrorCode ProcessOptions(LandauCtx *ctx, const char prefix[])
1102: {
1103:   PetscBool flg, fileflg;
1104:   PetscInt  ii, nt, nm, nc, num_species_grid[LANDAU_MAX_GRIDS], non_dim_grid;
1105:   PetscReal lnLam = 10;
1106:   DM        dummy;

1108:   PetscFunctionBegin;
1109:   PetscCall(DMCreate(ctx->comm, &dummy));
1110:   /* get options - initialize context */
1111:   ctx->verbose        = 1; // should be 0 for silent compliance
1112:   ctx->batch_sz       = 1;
1113:   ctx->batch_view_idx = 0;
1114:   ctx->interpolate    = PETSC_TRUE;
1115:   ctx->gpu_assembly   = PETSC_TRUE;
1116:   ctx->norm_state     = 0;
1117:   ctx->electronShift  = 0;
1118:   ctx->M              = NULL;
1119:   ctx->J              = NULL;
1120:   /* geometry and grids */
1121:   ctx->sphere    = PETSC_FALSE;
1122:   ctx->use_p4est = PETSC_FALSE;
1123:   ctx->simplex   = PETSC_FALSE;
1124:   for (PetscInt grid = 0; grid < LANDAU_MAX_GRIDS; grid++) {
1125:     ctx->radius[grid]             = 5.; /* thermal radius (velocity) */
1126:     ctx->radius_perp[grid]        = 5.; /* thermal radius (velocity) */
1127:     ctx->radius_par[grid]         = 5.; /* thermal radius (velocity) */
1128:     ctx->numAMRRefine[grid]       = 0;
1129:     ctx->postAMRRefine[grid]      = 0;
1130:     ctx->species_offset[grid + 1] = 1; // one species default
1131:     num_species_grid[grid]        = 0;
1132:     ctx->plex[grid]               = NULL; /* cache as expensive to Convert */
1133:   }
1134:   ctx->species_offset[0] = 0;
1135:   ctx->re_radius         = 0.;
1136:   ctx->vperp0_radius1    = 0;
1137:   ctx->vperp0_radius2    = 0;
1138:   ctx->nZRefine1         = 0;
1139:   ctx->nZRefine2         = 0;
1140:   ctx->numRERefine       = 0;
1141:   num_species_grid[0]    = 1; // one species default
1142:   /* species - [0] electrons, [1] one ion species eg, duetarium, [2] heavy impurity ion, ... */
1143:   ctx->charges[0]       = -1;                       /* electron charge (MKS) */
1144:   ctx->masses[0]        = 1 / 1835.469965278441013; /* temporary value in proton mass */
1145:   ctx->n[0]             = 1;
1146:   ctx->v_0              = 1; /* thermal velocity, we could start with a scale != 1 */
1147:   ctx->thermal_temps[0] = 1;
1148:   /* constants, etc. */
1149:   ctx->epsilon0 = 8.8542e-12;     /* permittivity of free space (MKS) F/m */
1150:   ctx->k        = 1.38064852e-23; /* Boltzmann constant (MKS) J/K */
1151:   ctx->n_0      = 1.e20;          /* typical plasma n, but could set it to 1 */
1152:   ctx->Ez       = 0;
1153:   for (PetscInt grid = 0; grid < LANDAU_NUM_TIMERS; grid++) ctx->times[grid] = 0;
1154:   for (PetscInt ii = 0; ii < LANDAU_DIM; ii++) ctx->cells0[ii] = 2;
1155:   if (LANDAU_DIM == 2) ctx->cells0[0] = 1;
1156:   ctx->use_matrix_mass                = PETSC_FALSE;
1157:   ctx->use_relativistic_corrections   = PETSC_FALSE;
1158:   ctx->use_energy_tensor_trick        = PETSC_FALSE; /* Use Eero's trick for energy conservation v --> grad(v^2/2) */
1159:   ctx->SData_d.w                      = NULL;
1160:   ctx->SData_d.x                      = NULL;
1161:   ctx->SData_d.y                      = NULL;
1162:   ctx->SData_d.z                      = NULL;
1163:   ctx->SData_d.invJ                   = NULL;
1164:   ctx->jacobian_field_major_order     = PETSC_FALSE;
1165:   ctx->SData_d.coo_elem_offsets       = NULL;
1166:   ctx->SData_d.coo_elem_point_offsets = NULL;
1167:   ctx->SData_d.coo_elem_fullNb        = NULL;
1168:   ctx->SData_d.coo_size               = 0;
1169:   PetscOptionsBegin(ctx->comm, prefix, "Options for Fokker-Plank-Landau collision operator", "none");
1170: #if defined(PETSC_HAVE_KOKKOS)
1171:   ctx->deviceType = LANDAU_KOKKOS;
1172:   PetscCall(PetscStrncpy(ctx->filename, "kokkos", sizeof(ctx->filename)));
1173: #else
1174:   ctx->deviceType = LANDAU_CPU;
1175:   PetscCall(PetscStrncpy(ctx->filename, "cpu", sizeof(ctx->filename)));
1176: #endif
1177:   PetscCall(PetscOptionsString("-dm_landau_device_type", "Use kernels on 'cpu' 'kokkos'", "plexland.c", ctx->filename, ctx->filename, sizeof(ctx->filename), NULL));
1178:   PetscCall(PetscStrcmp("cpu", ctx->filename, &flg));
1179:   if (flg) {
1180:     ctx->deviceType = LANDAU_CPU;
1181:   } else {
1182:     PetscCall(PetscStrcmp("kokkos", ctx->filename, &flg));
1183:     PetscCheck(flg, ctx->comm, PETSC_ERR_ARG_WRONG, "-dm_landau_device_type %s", ctx->filename);
1184:     ctx->deviceType = LANDAU_KOKKOS;
1185:   }
1186:   ctx->filename[0] = '\0';
1187:   PetscCall(PetscOptionsString("-dm_landau_filename", "file to read mesh from", "plexland.c", ctx->filename, ctx->filename, sizeof(ctx->filename), &fileflg));
1188:   PetscCall(PetscOptionsReal("-dm_landau_electron_shift", "Shift in thermal velocity of electrons", "none", ctx->electronShift, &ctx->electronShift, NULL));
1189:   PetscCall(PetscOptionsInt("-dm_landau_verbose", "Level of verbosity output", "plexland.c", ctx->verbose, &ctx->verbose, NULL));
1190:   PetscCall(PetscOptionsInt("-dm_landau_batch_size", "Number of 'vertices' to batch", "ex2.c", ctx->batch_sz, &ctx->batch_sz, NULL));
1191:   PetscCheck(LANDAU_MAX_BATCH_SZ >= ctx->batch_sz, ctx->comm, PETSC_ERR_ARG_WRONG, "LANDAU_MAX_BATCH_SZ %d < ctx->batch_sz %" PetscInt_FMT, LANDAU_MAX_BATCH_SZ, ctx->batch_sz);
1192:   PetscCall(PetscOptionsInt("-dm_landau_batch_view_idx", "Index of batch for diagnostics like plotting", "ex2.c", ctx->batch_view_idx, &ctx->batch_view_idx, NULL));
1193:   PetscCheck(ctx->batch_view_idx < ctx->batch_sz, ctx->comm, PETSC_ERR_ARG_WRONG, "-ctx->batch_view_idx %" PetscInt_FMT " > ctx->batch_sz %" PetscInt_FMT, ctx->batch_view_idx, ctx->batch_sz);
1194:   PetscCall(PetscOptionsReal("-dm_landau_Ez", "Initial parallel electric field in unites of Conner-Hastie critical field", "plexland.c", ctx->Ez, &ctx->Ez, NULL));
1195:   PetscCall(PetscOptionsReal("-dm_landau_n_0", "Normalization constant for number density", "plexland.c", ctx->n_0, &ctx->n_0, NULL));
1196:   PetscCall(PetscOptionsBool("-dm_landau_use_mataxpy_mass", "Use fast but slightly fragile MATAXPY to add mass term", "plexland.c", ctx->use_matrix_mass, &ctx->use_matrix_mass, NULL));
1197:   PetscCall(PetscOptionsBool("-dm_landau_use_relativistic_corrections", "Use relativistic corrections", "plexland.c", ctx->use_relativistic_corrections, &ctx->use_relativistic_corrections, NULL));
1198:   PetscCall(PetscOptionsBool("-dm_landau_simplex", "Use simplex elements", "plexland.c", ctx->simplex, &ctx->simplex, NULL));
1199:   PetscCall(PetscOptionsBool("-dm_landau_sphere", "use sphere/semi-circle domain instead of rectangle", "plexland.c", ctx->sphere, &ctx->sphere, NULL));
1200:   if (LANDAU_DIM == 2 && ctx->use_relativistic_corrections) ctx->use_relativistic_corrections = PETSC_FALSE; // should warn
1201:   PetscCall(PetscOptionsBool("-dm_landau_use_energy_tensor_trick", "Use Eero's trick of using grad(v^2/2) instead of v as args to Landau tensor to conserve energy with relativistic corrections and Q1 elements", "plexland.c", ctx->use_energy_tensor_trick,
1202:                              &ctx->use_energy_tensor_trick, NULL));

1204:   /* get num species with temperature, set defaults */
1205:   for (ii = 1; ii < LANDAU_MAX_SPECIES; ii++) {
1206:     ctx->thermal_temps[ii] = 1;
1207:     ctx->charges[ii]       = 1;
1208:     ctx->masses[ii]        = 1;
1209:     ctx->n[ii]             = 1;
1210:   }
1211:   nt = LANDAU_MAX_SPECIES;
1212:   PetscCall(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));
1213:   PetscCheck(flg, ctx->comm, PETSC_ERR_ARG_WRONG, "-dm_landau_thermal_temps ,t1,t2,.. must be provided to set the number of species");
1214:   PetscCall(PetscInfo(dummy, "num_species set to number of thermal temps provided (%" PetscInt_FMT ")\n", nt));
1215:   ctx->num_species = nt;
1216:   for (ii = 0; ii < ctx->num_species; ii++) ctx->thermal_temps[ii] *= 1.1604525e7; /* convert to Kelvin */
1217:   nm = LANDAU_MAX_SPECIES - 1;
1218:   PetscCall(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));
1219:   PetscCheck(!flg || nm == ctx->num_species - 1, ctx->comm, PETSC_ERR_ARG_WRONG, "num ion masses %" PetscInt_FMT " != num species %" PetscInt_FMT, nm, ctx->num_species - 1);
1220:   nm = LANDAU_MAX_SPECIES;
1221:   PetscCall(PetscOptionsRealArray("-dm_landau_n", "Number density of each species = n_s * n_0", "plexland.c", ctx->n, &nm, &flg));
1222:   PetscCheck(!flg || nm == ctx->num_species, ctx->comm, PETSC_ERR_ARG_WRONG, "wrong num n: %" PetscInt_FMT " != num species %" PetscInt_FMT, nm, ctx->num_species);
1223:   for (ii = 0; ii < LANDAU_MAX_SPECIES; ii++) ctx->masses[ii] *= 1.6720e-27; /* scale by proton mass kg */
1224:   ctx->masses[0] = 9.10938356e-31;                                           /* electron mass kg (should be about right already) */
1225:   nc             = LANDAU_MAX_SPECIES - 1;
1226:   PetscCall(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));
1227:   if (flg) PetscCheck(nc == ctx->num_species - 1, ctx->comm, PETSC_ERR_ARG_WRONG, "num charges %" PetscInt_FMT " != num species %" PetscInt_FMT, nc, ctx->num_species - 1);
1228:   for (ii = 0; ii < LANDAU_MAX_SPECIES; ii++) ctx->charges[ii] *= 1.6022e-19; /* electron/proton charge (MKS) */
1229:   /* geometry and grids */
1230:   nt = LANDAU_MAX_GRIDS;
1231:   PetscCall(PetscOptionsIntArray("-dm_landau_num_species_grid", "Number of species on each grid: [ 1, ....] or [S, 0 ....] for single grid", "plexland.c", num_species_grid, &nt, &flg));
1232:   if (flg) {
1233:     ctx->num_grids = nt;
1234:     for (ii = nt = 0; ii < ctx->num_grids; ii++) nt += num_species_grid[ii];
1235:     PetscCheck(ctx->num_species == nt, ctx->comm, PETSC_ERR_ARG_WRONG, "-dm_landau_num_species_grid: sum %" PetscInt_FMT " != num_species = %" PetscInt_FMT ". %" PetscInt_FMT " grids (check that number of grids <= LANDAU_MAX_GRIDS = %d)", nt, ctx->num_species,
1236:                ctx->num_grids, LANDAU_MAX_GRIDS);
1237:   } else {
1238:     if (ctx->num_species > LANDAU_MAX_GRIDS) {
1239:       num_species_grid[0] = 1;
1240:       num_species_grid[1] = ctx->num_species - 1;
1241:       ctx->num_grids      = 2;
1242:     } else {
1243:       ctx->num_grids = ctx->num_species;
1244:       for (ii = 0; ii < ctx->num_grids; ii++) num_species_grid[ii] = 1;
1245:     }
1246:   }
1247:   for (ctx->species_offset[0] = ii = 0; ii < ctx->num_grids; ii++) ctx->species_offset[ii + 1] = ctx->species_offset[ii] + num_species_grid[ii];
1248:   PetscCheck(ctx->species_offset[ctx->num_grids] == ctx->num_species, ctx->comm, PETSC_ERR_ARG_WRONG, "ctx->species_offset[ctx->num_grids] %" PetscInt_FMT " != ctx->num_species = %" PetscInt_FMT " ???????????", ctx->species_offset[ctx->num_grids],
1249:              ctx->num_species);
1250:   for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
1251:     PetscInt iii             = ctx->species_offset[grid];                                          // normalize with first (arbitrary) species on grid
1252:     ctx->thermal_speed[grid] = PetscSqrtReal(ctx->k * ctx->thermal_temps[iii] / ctx->masses[iii]); /* arbitrary units for non-dimensionalization: plasma formulary def */
1253:   }
1254:   // get lambdas here because we need them for t_0 etc
1255:   PetscCall(PetscOptionsReal("-dm_landau_ln_lambda", "Universal cross section parameter. Default uses NRL formulas", "plexland.c", lnLam, &lnLam, &flg));
1256:   if (flg) {
1257:     for (PetscInt grid = 0; grid < LANDAU_MAX_GRIDS; grid++) {
1258:       for (PetscInt gridj = 0; gridj < LANDAU_MAX_GRIDS; gridj++) ctx->lambdas[gridj][grid] = lnLam; /* cross section ratio large - small angle collisions */
1259:     }
1260:   } else {
1261:     PetscCall(makeLambdas(ctx));
1262:   }
1263:   non_dim_grid = 0;
1264:   PetscCall(PetscOptionsInt("-dm_landau_normalization_grid", "Index of grid to use for setting v_0, m_0, t_0. (Not recommended)", "plexland.c", non_dim_grid, &non_dim_grid, &flg));
1265:   if (non_dim_grid != 0) PetscCall(PetscInfo(dummy, "Normalization grid set to %" PetscInt_FMT ", but non-default not well verified\n", non_dim_grid));
1266:   PetscCheck(non_dim_grid >= 0 && non_dim_grid < ctx->num_species, ctx->comm, PETSC_ERR_ARG_WRONG, "Normalization grid wrong: %" PetscInt_FMT, non_dim_grid);
1267:   ctx->v_0 = ctx->thermal_speed[non_dim_grid]; /* arbitrary units for non dimensionalization: global mean velocity in 1D of electrons */
1268:   ctx->m_0 = ctx->masses[non_dim_grid];        /* arbitrary reference mass, electrons */
1269:   ctx->t_0 = 8 * PETSC_PI * PetscSqr(ctx->epsilon0 * ctx->m_0 / PetscSqr(ctx->charges[non_dim_grid])) / ctx->lambdas[non_dim_grid][non_dim_grid] / ctx->n_0 * PetscPowReal(ctx->v_0, 3); /* note, this t_0 makes nu[non_dim_grid,non_dim_grid]=1 */
1270:   /* domain */
1271:   nt = LANDAU_MAX_GRIDS;
1272:   PetscCall(PetscOptionsRealArray("-dm_landau_domain_radius", "Phase space size in units of thermal velocity of grid", "plexland.c", ctx->radius, &nt, &flg));
1273:   if (flg) {
1274:     PetscCheck(nt >= ctx->num_grids, ctx->comm, PETSC_ERR_ARG_WRONG, "-dm_landau_domain_radius: given %" PetscInt_FMT " radius != number grids %" PetscInt_FMT, nt, ctx->num_grids);
1275:     while (nt--) ctx->radius_par[nt] = ctx->radius_perp[nt] = ctx->radius[nt];
1276:   } else {
1277:     nt = LANDAU_MAX_GRIDS;
1278:     PetscCall(PetscOptionsRealArray("-dm_landau_domain_max_par", "Parallel velocity domain size in units of thermal velocity of grid", "plexland.c", ctx->radius_par, &nt, &flg));
1279:     if (flg) PetscCheck(nt >= ctx->num_grids, ctx->comm, PETSC_ERR_ARG_WRONG, "-dm_landau_domain_max_par: given %" PetscInt_FMT " radius != number grids %" PetscInt_FMT, nt, ctx->num_grids);
1280:     PetscCall(PetscOptionsRealArray("-dm_landau_domain_max_perp", "Perpendicular velocity domain size in units of thermal velocity of grid", "plexland.c", ctx->radius_perp, &nt, &flg));
1281:     if (flg) PetscCheck(nt >= ctx->num_grids, ctx->comm, PETSC_ERR_ARG_WRONG, "-dm_landau_domain_max_perp: given %" PetscInt_FMT " radius != number grids %" PetscInt_FMT, nt, ctx->num_grids);
1282:   }
1283:   for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
1284:     if (flg && ctx->radius[grid] <= 0) { /* negative is ratio of c - need to set par and perp with this -- todo */
1285:       if (ctx->radius[grid] == 0) ctx->radius[grid] = 0.75;
1286:       else ctx->radius[grid] = -ctx->radius[grid];
1287:       ctx->radius[grid] = ctx->radius[grid] * SPEED_OF_LIGHT / ctx->v_0; // use any species on grid to normalize (v_0 same for all on grid)
1288:       PetscCall(PetscInfo(dummy, "Change domain radius to %g for grid %" PetscInt_FMT "\n", (double)ctx->radius[grid], grid));
1289:     }
1290:     ctx->radius[grid] *= ctx->thermal_speed[grid] / ctx->v_0;      // scale domain by thermal radius relative to v_0
1291:     ctx->radius_perp[grid] *= ctx->thermal_speed[grid] / ctx->v_0; // scale domain by thermal radius relative to v_0
1292:     ctx->radius_par[grid] *= ctx->thermal_speed[grid] / ctx->v_0;  // scale domain by thermal radius relative to v_0
1293:   }
1294:   /* amr parameters */
1295:   if (!fileflg) {
1296:     nt = LANDAU_MAX_GRIDS;
1297:     PetscCall(PetscOptionsIntArray("-dm_landau_amr_levels_max", "Number of AMR levels of refinement around origin, after (RE) refinements along z", "plexland.c", ctx->numAMRRefine, &nt, &flg));
1298:     PetscCheck(!flg || nt >= ctx->num_grids, ctx->comm, PETSC_ERR_ARG_WRONG, "-dm_landau_amr_levels_max: given %" PetscInt_FMT " != number grids %" PetscInt_FMT, nt, ctx->num_grids);
1299:     nt = LANDAU_MAX_GRIDS;
1300:     PetscCall(PetscOptionsIntArray("-dm_landau_amr_post_refine", "Number of levels to uniformly refine after AMR", "plexland.c", ctx->postAMRRefine, &nt, &flg));
1301:     for (ii = 1; ii < ctx->num_grids; ii++) ctx->postAMRRefine[ii] = ctx->postAMRRefine[0]; // all grids the same now
1302:     PetscCall(PetscOptionsInt("-dm_landau_amr_re_levels", "Number of levels to refine along v_perp=0, z>0", "plexland.c", ctx->numRERefine, &ctx->numRERefine, &flg));
1303:     PetscCall(PetscOptionsInt("-dm_landau_amr_z_refine_pre", "Number of levels to refine along v_perp=0 before origin refine", "plexland.c", ctx->nZRefine1, &ctx->nZRefine1, &flg));
1304:     PetscCall(PetscOptionsInt("-dm_landau_amr_z_refine_post", "Number of levels to refine along v_perp=0 after origin refine", "plexland.c", ctx->nZRefine2, &ctx->nZRefine2, &flg));
1305:     PetscCall(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));
1306:     PetscCall(PetscOptionsReal("-dm_landau_z_radius_pre", "velocity range to refine r=0 axis (for electrons)", "plexland.c", ctx->vperp0_radius1, &ctx->vperp0_radius1, &flg));
1307:     PetscCall(PetscOptionsReal("-dm_landau_z_radius_post", "velocity range to refine r=0 axis (for electrons) after origin AMR", "plexland.c", ctx->vperp0_radius2, &ctx->vperp0_radius2, &flg));
1308:     /* spherical domain */
1309:     if (ctx->sphere || ctx->simplex) {
1310:       ctx->sphere_uniform_normal = PETSC_FALSE;
1311:       PetscCall(PetscOptionsBool("-dm_landau_sphere_uniform_normal", "Scaling of circle radius to get uniform particles per cell with Maxwellians (not used)", "plexland.c", ctx->sphere_uniform_normal, &ctx->sphere_uniform_normal, NULL));
1312:       if (!ctx->sphere_uniform_normal) { // true
1313:         nt = LANDAU_MAX_GRIDS;
1314:         PetscCall(PetscOptionsRealArray("-dm_landau_sphere_inner_radius_90degree_scale", "Scaling of radius for inner circle on 90 degree grid", "plexland.c", ctx->sphere_inner_radius_90degree, &nt, &flg));
1315:         if (flg && nt < ctx->num_grids) {
1316:           for (PetscInt grid = nt; grid < ctx->num_grids; grid++) ctx->sphere_inner_radius_90degree[grid] = ctx->sphere_inner_radius_90degree[0];
1317:         } else if (!flg || nt == 0) {
1318:           if (LANDAU_DIM == 2) {
1319:             for (PetscInt grid = 0; grid < ctx->num_grids; grid++) ctx->sphere_inner_radius_90degree[grid] = 0.4; // optimized for R=5, Q4, AMR=0
1320:           } else {
1321:             for (PetscInt grid = 0; grid < ctx->num_grids; grid++) ctx->sphere_inner_radius_90degree[grid] = 0.577 * 0.40;
1322:           }
1323:         }
1324:         nt = LANDAU_MAX_GRIDS;
1325:         PetscCall(PetscOptionsRealArray("-dm_landau_sphere_inner_radius_45degree_scale", "Scaling of radius for inner circle on 45 degree grid", "plexland.c", ctx->sphere_inner_radius_45degree, &nt, &flg));
1326:         if (flg && nt < ctx->num_grids) {
1327:           for (PetscInt grid = nt; grid < ctx->num_grids; grid++) ctx->sphere_inner_radius_45degree[grid] = ctx->sphere_inner_radius_45degree[0];
1328:         } else if (!flg || nt == 0) {
1329:           if (LANDAU_DIM == 2) {
1330:             for (PetscInt grid = 0; grid < ctx->num_grids; grid++) ctx->sphere_inner_radius_45degree[grid] = 0.45; // optimized for R=5, Q4, AMR=0
1331:           } else {
1332:             for (PetscInt grid = 0; grid < ctx->num_grids; grid++) ctx->sphere_inner_radius_45degree[grid] = 0.4; // 3D sphere
1333:           }
1334:         }
1335:         if (ctx->sphere) PetscCall(PetscInfo(ctx->plex[0], "sphere : , 45 degree scaling = %g; 90 degree scaling = %g\n", (double)ctx->sphere_inner_radius_45degree[0], (double)ctx->sphere_inner_radius_90degree[0]));
1336:       } else {
1337:         for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
1338:           switch (ctx->numAMRRefine[grid]) {
1339:           case 0:
1340:           case 1:
1341:           case 2:
1342:           case 3:
1343:           default:
1344:             if (LANDAU_DIM == 2) {
1345:               ctx->sphere_inner_radius_90degree[grid] = 0.40;
1346:               ctx->sphere_inner_radius_45degree[grid] = 0.45;
1347:             } else {
1348:               ctx->sphere_inner_radius_45degree[grid] = 0.25;
1349:             }
1350:           }
1351:         }
1352:       }
1353:     } else {
1354:       nt = LANDAU_DIM;
1355:       PetscCall(PetscOptionsIntArray("-dm_landau_num_cells", "Number of cells in each dimension of base grid", "plexland.c", ctx->cells0, &nt, &flg));
1356:     }
1357:   }
1358:   /* processing options */
1359:   PetscCall(PetscOptionsBool("-dm_landau_gpu_assembly", "Assemble Jacobian on GPU", "plexland.c", ctx->gpu_assembly, &ctx->gpu_assembly, NULL));
1360:   PetscCall(PetscOptionsBool("-dm_landau_jacobian_field_major_order", "Reorder Jacobian for GPU assembly with field major, or block diagonal, ordering (DEPRECATED)", "plexland.c", ctx->jacobian_field_major_order, &ctx->jacobian_field_major_order, NULL));
1361:   if (ctx->jacobian_field_major_order) PetscCheck(ctx->gpu_assembly, ctx->comm, PETSC_ERR_ARG_WRONG, "-dm_landau_jacobian_field_major_order requires -dm_landau_gpu_assembly");
1362:   PetscCheck(!ctx->jacobian_field_major_order, ctx->comm, PETSC_ERR_ARG_WRONG, "-dm_landau_jacobian_field_major_order DEPRECATED");
1363:   PetscOptionsEnd();

1365:   for (ii = ctx->num_species; ii < LANDAU_MAX_SPECIES; ii++) ctx->masses[ii] = ctx->thermal_temps[ii] = ctx->charges[ii] = 0;
1366:   if (ctx->verbose != 0) {
1367:     PetscReal pmassunit = PetscRealConstant(1.6720e-27);

1369:     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "masses:        e=%10.3e; ions in proton mass units:   %10.3e %10.3e ...\n", (double)ctx->masses[0], (double)(ctx->masses[1] / pmassunit), (double)(ctx->num_species > 2 ? ctx->masses[2] / pmassunit : 0)));
1370:     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "charges:       e=%10.3e; charges in elementary units: %10.3e %10.3e\n", (double)ctx->charges[0], (double)(-ctx->charges[1] / ctx->charges[0]), (double)(ctx->num_species > 2 ? -ctx->charges[2] / ctx->charges[0] : 0)));
1371:     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "n:             e: %10.3e                           i: %10.3e %10.3e\n", (double)ctx->n[0], (double)ctx->n[1], (double)(ctx->num_species > 2 ? ctx->n[2] : 0)));
1372:     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "thermal T (K): e=%10.3e i=%10.3e %10.3e. Normalization grid %" PetscInt_FMT ": v_0=%10.3e (%10.3ec) n_0=%10.3e t_0=%10.3e %" PetscInt_FMT " batched, view batch %" PetscInt_FMT "\n", (double)ctx->thermal_temps[0],
1373:                           (double)ctx->thermal_temps[1], (double)((ctx->num_species > 2) ? ctx->thermal_temps[2] : 0), non_dim_grid, (double)ctx->v_0, (double)(ctx->v_0 / SPEED_OF_LIGHT), (double)ctx->n_0, (double)ctx->t_0, ctx->batch_sz, ctx->batch_view_idx));
1374:     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Domain radius (AMR levels) grid %d: par=%10.3e perp=%10.3e (%" PetscInt_FMT ") ", 0, (double)ctx->radius_par[0], (double)ctx->radius_perp[0], ctx->numAMRRefine[0]));
1375:     for (ii = 1; ii < ctx->num_grids; ii++) PetscCall(PetscPrintf(PETSC_COMM_WORLD, ", %" PetscInt_FMT ": par=%10.3e perp=%10.3e (%" PetscInt_FMT ") ", ii, (double)ctx->radius_par[ii], (double)ctx->radius_perp[ii], ctx->numAMRRefine[ii]));
1376:     if (ctx->use_relativistic_corrections) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nUse relativistic corrections\n"));
1377:     else PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n"));
1378:   }
1379:   PetscCall(DMDestroy(&dummy));
1380:   {
1381:     PetscMPIInt rank;
1382:     PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
1383:     ctx->stage = 0;
1384:     PetscCall(PetscLogEventRegister("Landau Create", DM_CLASSID, &ctx->events[13]));   /* 13 */
1385:     PetscCall(PetscLogEventRegister(" GPU ass. setup", DM_CLASSID, &ctx->events[2]));  /* 2 */
1386:     PetscCall(PetscLogEventRegister(" Build matrix", DM_CLASSID, &ctx->events[12]));   /* 12 */
1387:     PetscCall(PetscLogEventRegister(" Assembly maps", DM_CLASSID, &ctx->events[15]));  /* 15 */
1388:     PetscCall(PetscLogEventRegister("Landau Mass mat", DM_CLASSID, &ctx->events[14])); /* 14 */
1389:     PetscCall(PetscLogEventRegister("Landau Operator", DM_CLASSID, &ctx->events[11])); /* 11 */
1390:     PetscCall(PetscLogEventRegister("Landau Jacobian", DM_CLASSID, &ctx->events[0]));  /* 0 */
1391:     PetscCall(PetscLogEventRegister("Landau Mass", DM_CLASSID, &ctx->events[9]));      /* 9 */
1392:     PetscCall(PetscLogEventRegister(" Preamble", DM_CLASSID, &ctx->events[10]));       /* 10 */
1393:     PetscCall(PetscLogEventRegister(" static IP Data", DM_CLASSID, &ctx->events[7]));  /* 7 */
1394:     PetscCall(PetscLogEventRegister(" dynamic IP-Jac", DM_CLASSID, &ctx->events[1]));  /* 1 */
1395:     PetscCall(PetscLogEventRegister(" Kernel-init", DM_CLASSID, &ctx->events[3]));     /* 3 */
1396:     PetscCall(PetscLogEventRegister(" Jac-f-df (GPU)", DM_CLASSID, &ctx->events[8]));  /* 8 */
1397:     PetscCall(PetscLogEventRegister(" J Kernel (GPU)", DM_CLASSID, &ctx->events[4]));  /* 4 */
1398:     PetscCall(PetscLogEventRegister(" M Kernel (GPU)", DM_CLASSID, &ctx->events[16])); /* 16 */
1399:     PetscCall(PetscLogEventRegister(" Copy to CPU", DM_CLASSID, &ctx->events[5]));     /* 5 */
1400:     PetscCall(PetscLogEventRegister(" CPU assemble", DM_CLASSID, &ctx->events[6]));    /* 6 */

1402:     if (rank) { /* turn off output stuff for duplicate runs - do we need to add the prefix to all this? */
1403:       PetscCall(PetscOptionsClearValue(NULL, "-snes_converged_reason"));
1404:       PetscCall(PetscOptionsClearValue(NULL, "-ksp_converged_reason"));
1405:       PetscCall(PetscOptionsClearValue(NULL, "-snes_monitor"));
1406:       PetscCall(PetscOptionsClearValue(NULL, "-ksp_monitor"));
1407:       PetscCall(PetscOptionsClearValue(NULL, "-ts_monitor"));
1408:       PetscCall(PetscOptionsClearValue(NULL, "-ts_view"));
1409:       PetscCall(PetscOptionsClearValue(NULL, "-ts_adapt_monitor"));
1410:       PetscCall(PetscOptionsClearValue(NULL, "-dm_landau_amr_dm_view"));
1411:       PetscCall(PetscOptionsClearValue(NULL, "-dm_landau_amr_vec_view"));
1412:       PetscCall(PetscOptionsClearValue(NULL, "-dm_landau_mass_dm_view"));
1413:       PetscCall(PetscOptionsClearValue(NULL, "-dm_landau_mass_view"));
1414:       PetscCall(PetscOptionsClearValue(NULL, "-dm_landau_jacobian_view"));
1415:       PetscCall(PetscOptionsClearValue(NULL, "-dm_landau_mat_view"));
1416:       PetscCall(PetscOptionsClearValue(NULL, "-pc_bjkokkos_ksp_converged_reason"));
1417:       PetscCall(PetscOptionsClearValue(NULL, "-pc_bjkokkos_ksp_monitor"));
1418:       PetscCall(PetscOptionsClearValue(NULL, "-"));
1419:       PetscCall(PetscOptionsClearValue(NULL, "-info"));
1420:     }
1421:   }
1422:   PetscFunctionReturn(PETSC_SUCCESS);
1423: }

1425: static PetscErrorCode CreateStaticData(PetscInt dim, IS grid_batch_is_inv[], LandauCtx *ctx)
1426: {
1427:   PetscSection     section[LANDAU_MAX_GRIDS], globsection[LANDAU_MAX_GRIDS];
1428:   PetscQuadrature  quad;
1429:   const PetscReal *quadWeights;
1430:   PetscReal        invMass[LANDAU_MAX_SPECIES], nu_alpha[LANDAU_MAX_SPECIES], nu_beta[LANDAU_MAX_SPECIES];
1431:   PetscInt         numCells[LANDAU_MAX_GRIDS], Nq, Nb, Nf[LANDAU_MAX_GRIDS], ncellsTot = 0, MAP_BF_SIZE = 64 * LANDAU_DIM * LANDAU_DIM * LANDAU_MAX_Q_FACE * LANDAU_MAX_SPECIES;
1432:   PetscTabulation *Tf;
1433:   PetscDS          prob;

1435:   PetscFunctionBegin;
1436:   PetscCall(PetscFEGetDimension(ctx->fe[0], &Nb));
1437:   PetscCheck(Nb <= LANDAU_MAX_NQND, ctx->comm, PETSC_ERR_ARG_WRONG, "Order too high. Nb = %" PetscInt_FMT " > LANDAU_MAX_NQND (%d)", Nb, LANDAU_MAX_NQND);
1438:   for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
1439:     for (PetscInt ii = ctx->species_offset[grid]; ii < ctx->species_offset[grid + 1]; ii++) {
1440:       invMass[ii]  = ctx->m_0 / ctx->masses[ii];
1441:       nu_alpha[ii] = PetscSqr(ctx->charges[ii] / ctx->m_0) * ctx->m_0 / ctx->masses[ii];
1442:       nu_beta[ii]  = PetscSqr(ctx->charges[ii] / ctx->epsilon0) / (8 * PETSC_PI) * ctx->t_0 * ctx->n_0 / PetscPowReal(ctx->v_0, 3);
1443:     }
1444:   }
1445:   if (ctx->verbose == 4) {
1446:     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "nu_alpha: "));
1447:     for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
1448:       PetscInt iii = ctx->species_offset[grid];
1449:       for (PetscInt ii = iii; ii < ctx->species_offset[grid + 1]; ii++) PetscCall(PetscPrintf(PETSC_COMM_WORLD, " %e", (double)nu_alpha[ii]));
1450:     }
1451:     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nnu_beta: "));
1452:     for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
1453:       PetscInt iii = ctx->species_offset[grid];
1454:       for (PetscInt ii = iii; ii < ctx->species_offset[grid + 1]; ii++) PetscCall(PetscPrintf(PETSC_COMM_WORLD, " %e", (double)nu_beta[ii]));
1455:     }
1456:     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nnu_alpha[i]*nu_beta[j]*lambda[i][j]:\n"));
1457:     for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
1458:       PetscInt iii = ctx->species_offset[grid];
1459:       for (PetscInt ii = iii; ii < ctx->species_offset[grid + 1]; ii++) {
1460:         for (PetscInt gridj = 0; gridj < ctx->num_grids; gridj++) {
1461:           PetscInt jjj = ctx->species_offset[gridj];
1462:           for (PetscInt jj = jjj; jj < ctx->species_offset[gridj + 1]; jj++) PetscCall(PetscPrintf(PETSC_COMM_WORLD, " %14.9e", (double)(nu_alpha[ii] * nu_beta[jj] * ctx->lambdas[grid][gridj])));
1463:         }
1464:         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n"));
1465:       }
1466:     }
1467:     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "lambda[i][j]:\n"));
1468:     for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
1469:       PetscInt iii = ctx->species_offset[grid];
1470:       for (PetscInt ii = iii; ii < ctx->species_offset[grid + 1]; ii++) {
1471:         for (PetscInt gridj = 0; gridj < ctx->num_grids; gridj++) {
1472:           PetscInt jjj = ctx->species_offset[gridj];
1473:           for (PetscInt jj = jjj; jj < ctx->species_offset[gridj + 1]; jj++) PetscCall(PetscPrintf(PETSC_COMM_WORLD, " %14.9e", (double)ctx->lambdas[grid][gridj]));
1474:         }
1475:         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n"));
1476:       }
1477:     }
1478:   }
1479:   PetscCall(DMGetDS(ctx->plex[0], &prob));    // same DS for all grids
1480:   PetscCall(PetscDSGetTabulation(prob, &Tf)); // Bf, &Df same for all grids
1481:   /* DS, Tab and quad is same on all grids */
1482:   PetscCheck(ctx->plex[0], ctx->comm, PETSC_ERR_ARG_WRONG, "Plex not created");
1483:   PetscCall(PetscFEGetQuadrature(ctx->fe[0], &quad));
1484:   PetscCall(PetscQuadratureGetData(quad, NULL, NULL, &Nq, NULL, &quadWeights));
1485:   PetscCheck(Nq <= LANDAU_MAX_NQND, ctx->comm, PETSC_ERR_ARG_WRONG, "Order too high. Nq = %" PetscInt_FMT " > LANDAU_MAX_NQND (%d)", Nq, LANDAU_MAX_NQND);
1486:   /* setup each grid */
1487:   for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
1488:     PetscInt cStart, cEnd;
1489:     PetscCheck(ctx->plex[grid] != NULL, ctx->comm, PETSC_ERR_ARG_WRONG, "Plex not created");
1490:     PetscCall(DMPlexGetHeightStratum(ctx->plex[grid], 0, &cStart, &cEnd));
1491:     numCells[grid] = cEnd - cStart; // grids can have different topology
1492:     PetscCall(DMGetLocalSection(ctx->plex[grid], &section[grid]));
1493:     PetscCall(DMGetGlobalSection(ctx->plex[grid], &globsection[grid]));
1494:     PetscCall(PetscSectionGetNumFields(section[grid], &Nf[grid]));
1495:     ncellsTot += numCells[grid];
1496:   }
1497:   /* create GPU assembly data */
1498:   if (ctx->gpu_assembly) { /* we need GPU object with GPU assembly */
1499:     PetscContainer container;
1500:     PetscScalar   *elemMatrix, *elMat;
1501:     pointInterpolationP4est(*pointMaps)[LANDAU_MAX_Q_FACE];
1502:     P4estVertexMaps *maps;
1503:     const PetscInt  *plex_batch = NULL, elMatSz = Nb * Nb * ctx->num_species * ctx->num_species;
1504:     LandauIdx       *coo_elem_offsets = NULL, *coo_elem_fullNb = NULL, (*coo_elem_point_offsets)[LANDAU_MAX_NQND + 1] = NULL;
1505:     /* create GPU assembly data */
1506:     PetscCall(PetscInfo(ctx->plex[0], "Make GPU maps %d\n", 1));
1507:     PetscCall(PetscLogEventBegin(ctx->events[2], 0, 0, 0, 0));
1508:     PetscCall(PetscMalloc(sizeof(*maps) * ctx->num_grids, &maps));
1509:     PetscCall(PetscMalloc(sizeof(*pointMaps) * MAP_BF_SIZE, &pointMaps));
1510:     PetscCall(PetscMalloc(sizeof(*elemMatrix) * elMatSz, &elemMatrix));

1512:     {                                                                                                                             // setup COO assembly -- put COO metadata directly in ctx->SData_d
1513:       PetscCall(PetscMalloc3(ncellsTot + 1, &coo_elem_offsets, ncellsTot, &coo_elem_fullNb, ncellsTot, &coo_elem_point_offsets)); // array of integer pointers
1514:       coo_elem_offsets[0] = 0;                                                                                                    // finish later
1515:       PetscCall(PetscInfo(ctx->plex[0], "COO initialization, %" PetscInt_FMT " cells\n", ncellsTot));
1516:       ctx->SData_d.coo_n_cellsTot         = ncellsTot;
1517:       ctx->SData_d.coo_elem_offsets       = (void *)coo_elem_offsets;
1518:       ctx->SData_d.coo_elem_fullNb        = (void *)coo_elem_fullNb;
1519:       ctx->SData_d.coo_elem_point_offsets = (void *)coo_elem_point_offsets;
1520:     }

1522:     ctx->SData_d.coo_max_fullnb = 0;
1523:     for (PetscInt grid = 0, glb_elem_idx = 0; grid < ctx->num_grids; grid++) {
1524:       PetscInt cStart, cEnd, Nfloc = Nf[grid], totDim = Nfloc * Nb;
1525:       if (grid_batch_is_inv[grid]) PetscCall(ISGetIndices(grid_batch_is_inv[grid], &plex_batch));
1526:       PetscCheck(!plex_batch, ctx->comm, PETSC_ERR_ARG_WRONG, "-dm_landau_jacobian_field_major_order DEPRECATED");
1527:       PetscCall(DMPlexGetHeightStratum(ctx->plex[grid], 0, &cStart, &cEnd));
1528:       // make maps
1529:       maps[grid].d_self       = NULL;
1530:       maps[grid].num_elements = numCells[grid];
1531:       maps[grid].num_face     = (PetscInt)(pow(Nq, 1. / ((double)dim)) + .001);                 // Q
1532:       maps[grid].num_face     = (PetscInt)(pow(maps[grid].num_face, (double)(dim - 1)) + .001); // Q^2
1533:       maps[grid].num_reduced  = 0;
1534:       maps[grid].deviceType   = ctx->deviceType;
1535:       maps[grid].numgrids     = ctx->num_grids;
1536:       // count reduced and get
1537:       PetscCall(PetscMalloc(maps[grid].num_elements * sizeof(*maps[grid].gIdx), &maps[grid].gIdx));
1538:       for (PetscInt ej = cStart, eidx = 0; ej < cEnd; ++ej, ++eidx, glb_elem_idx++) {
1539:         if (coo_elem_offsets) coo_elem_offsets[glb_elem_idx + 1] = coo_elem_offsets[glb_elem_idx]; // start with last one, then add
1540:         for (PetscInt fieldA = 0; fieldA < Nf[grid]; fieldA++) {
1541:           PetscInt fullNb = 0;
1542:           for (PetscInt q = 0; q < Nb; ++q) {
1543:             PetscInt     numindices, *indices;
1544:             PetscScalar *valuesOrig = elMat = elemMatrix;
1545:             PetscCall(PetscArrayzero(elMat, totDim * totDim));
1546:             elMat[(fieldA * Nb + q) * totDim + fieldA * Nb + q] = 1;
1547:             PetscCall(DMPlexGetClosureIndices(ctx->plex[grid], section[grid], globsection[grid], ej, PETSC_TRUE, &numindices, &indices, NULL, &elMat));
1548:             if (ctx->simplex) {
1549:               PetscCheck(numindices == Nb, ctx->comm, PETSC_ERR_ARG_WRONG, "numindices != Nb numindices=%" PetscInt_FMT " Nb=%" PetscInt_FMT, numindices, Nb);
1550:               for (PetscInt q = 0; q < numindices; ++q) maps[grid].gIdx[eidx][fieldA][q] = indices[q];
1551:               fullNb++;
1552:             } else {
1553:               for (PetscInt f = 0; f < numindices; ++f) { // look for a non-zero on the diagonal (is this too complicated for simplices?)
1554:                 if (PetscAbs(PetscRealPart(elMat[f * numindices + f])) > PETSC_MACHINE_EPSILON) {
1555:                   // found it
1556:                   if (PetscAbs(PetscRealPart(elMat[f * numindices + f] - 1.)) < PETSC_MACHINE_EPSILON) { // normal vertex 1.0
1557:                     if (plex_batch) {
1558:                       maps[grid].gIdx[eidx][fieldA][q] = plex_batch[indices[f]];
1559:                     } else {
1560:                       maps[grid].gIdx[eidx][fieldA][q] = indices[f];
1561:                     }
1562:                     fullNb++;
1563:                   } else { //found a constraint
1564:                     PetscInt       jj                = 0;
1565:                     PetscReal      sum               = 0;
1566:                     const PetscInt ff                = f;
1567:                     maps[grid].gIdx[eidx][fieldA][q] = -maps[grid].num_reduced - 1; // store (-)index: id = -(idx+1): idx = -id - 1
1568:                     PetscCheck(!ctx->simplex, ctx->comm, PETSC_ERR_ARG_WRONG, "No constraints with simplex");
1569:                     do {                                                                                              // constraints are continuous in Plex - exploit that here
1570:                       PetscInt ii;                                                                                    // get 'scale'
1571:                       for (ii = 0, pointMaps[maps[grid].num_reduced][jj].scale = 0; ii < maps[grid].num_face; ii++) { // sum row of outer product to recover vector value
1572:                         if (ff + ii < numindices) {                                                                   // 3D has Q and Q^2 interps so might run off end. We could test that elMat[f*numindices + ff + ii] > 0, and break if not
1573:                           pointMaps[maps[grid].num_reduced][jj].scale += PetscRealPart(elMat[f * numindices + ff + ii]);
1574:                         }
1575:                       }
1576:                       sum += pointMaps[maps[grid].num_reduced][jj].scale; // diagnostic
1577:                       // get 'gid'
1578:                       if (pointMaps[maps[grid].num_reduced][jj].scale == 0) pointMaps[maps[grid].num_reduced][jj].gid = -1; // 3D has Q and Q^2 interps
1579:                       else {
1580:                         if (plex_batch) {
1581:                           pointMaps[maps[grid].num_reduced][jj].gid = plex_batch[indices[f]];
1582:                         } else {
1583:                           pointMaps[maps[grid].num_reduced][jj].gid = indices[f];
1584:                         }
1585:                         fullNb++;
1586:                       }
1587:                     } while (++jj < maps[grid].num_face && ++f < numindices); // jj is incremented if we hit the end
1588:                     while (jj < maps[grid].num_face) {
1589:                       pointMaps[maps[grid].num_reduced][jj].scale = 0;
1590:                       pointMaps[maps[grid].num_reduced][jj].gid   = -1;
1591:                       jj++;
1592:                     }
1593:                     if (PetscAbs(sum - 1.0) > 10 * PETSC_MACHINE_EPSILON) { // debug
1594:                       PetscInt  d, f;
1595:                       PetscReal tmp = 0;
1596:                       PetscCall(
1597:                         PetscPrintf(PETSC_COMM_SELF, "\t\t%" PetscInt_FMT ".%" PetscInt_FMT ".%" PetscInt_FMT ") ERROR total I = %22.16e (LANDAU_MAX_Q_FACE=%d, #face=%" PetscInt_FMT ")\n", eidx, q, fieldA, (double)sum, LANDAU_MAX_Q_FACE, maps[grid].num_face));
1598:                       for (d = 0, tmp = 0; d < numindices; ++d) {
1599:                         if (tmp != 0 && PetscAbs(tmp - 1.0) > 10 * PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%3" PetscInt_FMT ") %3" PetscInt_FMT ": ", d, indices[d]));
1600:                         for (f = 0; f < numindices; ++f) tmp += PetscRealPart(elMat[d * numindices + f]);
1601:                         if (tmp != 0) PetscCall(PetscPrintf(ctx->comm, " | %22.16e\n", (double)tmp));
1602:                       }
1603:                     }
1604:                     maps[grid].num_reduced++;
1605:                     PetscCheck(maps[grid].num_reduced < MAP_BF_SIZE, PETSC_COMM_SELF, PETSC_ERR_PLIB, "maps[grid].num_reduced %" PetscInt_FMT " > %" PetscInt_FMT, maps[grid].num_reduced, MAP_BF_SIZE);
1606:                   }
1607:                   break;
1608:                 }
1609:               }
1610:             } // !simplex
1611:             // cleanup
1612:             PetscCall(DMPlexRestoreClosureIndices(ctx->plex[grid], section[grid], globsection[grid], ej, PETSC_TRUE, &numindices, &indices, NULL, &elMat));
1613:             if (elMat != valuesOrig) PetscCall(DMRestoreWorkArray(ctx->plex[grid], numindices * numindices, MPIU_SCALAR, &elMat));
1614:           }
1615:           {                                                        // setup COO assembly
1616:             coo_elem_offsets[glb_elem_idx + 1] += fullNb * fullNb; // one species block, adds a block for each species, on this element in this grid
1617:             if (fieldA == 0) {                                     // cache full Nb for this element, on this grid per species
1618:               coo_elem_fullNb[glb_elem_idx] = fullNb;
1619:               if (fullNb > ctx->SData_d.coo_max_fullnb) ctx->SData_d.coo_max_fullnb = fullNb;
1620:             } else PetscCheck(coo_elem_fullNb[glb_elem_idx] == fullNb, PETSC_COMM_SELF, PETSC_ERR_PLIB, "full element size change with species %" PetscInt_FMT " %" PetscInt_FMT, coo_elem_fullNb[glb_elem_idx], fullNb);
1621:           }
1622:         } // field
1623:       } // cell
1624:       // allocate and copy point data maps[grid].gIdx[eidx][field][q]
1625:       PetscCall(PetscMalloc(maps[grid].num_reduced * sizeof(*maps[grid].c_maps), &maps[grid].c_maps));
1626:       for (PetscInt ej = 0; ej < maps[grid].num_reduced; ++ej) {
1627:         for (PetscInt q = 0; q < maps[grid].num_face; ++q) {
1628:           maps[grid].c_maps[ej][q].scale = pointMaps[ej][q].scale;
1629:           maps[grid].c_maps[ej][q].gid   = pointMaps[ej][q].gid;
1630:         }
1631:       }
1632: #if defined(PETSC_HAVE_KOKKOS)
1633:       if (ctx->deviceType == LANDAU_KOKKOS) PetscCall(LandauKokkosCreateMatMaps(maps, pointMaps, Nf, grid)); // implies Kokkos does
1634: #endif
1635:       if (plex_batch) {
1636:         PetscCall(ISRestoreIndices(grid_batch_is_inv[grid], &plex_batch));
1637:         PetscCall(ISDestroy(&grid_batch_is_inv[grid])); // we are done with this
1638:       }
1639:     } /* grids */
1640:     // finish COO
1641:     { // setup COO assembly
1642:       PetscInt *oor, *ooc;
1643:       ctx->SData_d.coo_size = coo_elem_offsets[ncellsTot] * ctx->batch_sz;
1644:       PetscCall(PetscMalloc2(ctx->SData_d.coo_size, &oor, ctx->SData_d.coo_size, &ooc));
1645:       for (PetscInt i = 0; i < ctx->SData_d.coo_size; i++) oor[i] = ooc[i] = -1;
1646:       // get
1647:       for (PetscInt grid = 0, glb_elem_idx = 0; grid < ctx->num_grids; grid++) {
1648:         for (PetscInt ej = 0; ej < numCells[grid]; ++ej, glb_elem_idx++) {
1649:           const PetscInt         fullNb           = coo_elem_fullNb[glb_elem_idx];
1650:           const LandauIdx *const Idxs             = &maps[grid].gIdx[ej][0][0]; // just use field-0 maps, They should be the same but this is just for COO storage
1651:           coo_elem_point_offsets[glb_elem_idx][0] = 0;
1652:           for (PetscInt f = 0, cnt2 = 0; f < Nb; f++) {
1653:             PetscInt idx                                = Idxs[f];
1654:             coo_elem_point_offsets[glb_elem_idx][f + 1] = coo_elem_point_offsets[glb_elem_idx][f]; // start at last
1655:             if (idx >= 0) {
1656:               cnt2++;
1657:               coo_elem_point_offsets[glb_elem_idx][f + 1]++; // inc
1658:             } else {
1659:               idx = -idx - 1;
1660:               for (PetscInt q = 0; q < maps[grid].num_face; q++) {
1661:                 if (maps[grid].c_maps[idx][q].gid < 0) break;
1662:                 cnt2++;
1663:                 coo_elem_point_offsets[glb_elem_idx][f + 1]++; // inc
1664:               }
1665:             }
1666:             PetscCheck(cnt2 <= fullNb, PETSC_COMM_SELF, PETSC_ERR_PLIB, "wrong count %" PetscInt_FMT " < %" PetscInt_FMT, fullNb, cnt2);
1667:           }
1668:           PetscCheck(coo_elem_point_offsets[glb_elem_idx][Nb] == fullNb, PETSC_COMM_SELF, PETSC_ERR_PLIB, "coo_elem_point_offsets size %" PetscInt_FMT " != fullNb=%" PetscInt_FMT, coo_elem_point_offsets[glb_elem_idx][Nb], fullNb);
1669:         }
1670:       }
1671:       // set
1672:       for (PetscInt b_id = 0; b_id < ctx->batch_sz; b_id++) {
1673:         for (PetscInt grid = 0, glb_elem_idx = 0; grid < ctx->num_grids; grid++) {
1674:           const PetscInt moffset = LAND_MOFFSET(b_id, grid, ctx->batch_sz, ctx->num_grids, ctx->mat_offset);
1675:           for (PetscInt ej = 0; ej < numCells[grid]; ++ej, glb_elem_idx++) {
1676:             const PetscInt fullNb = coo_elem_fullNb[glb_elem_idx], fullNb2 = fullNb * fullNb;
1677:             // set (i,j)
1678:             for (PetscInt fieldA = 0; fieldA < Nf[grid]; fieldA++) {
1679:               const LandauIdx *const Idxs = &maps[grid].gIdx[ej][fieldA][0];
1680:               PetscInt               rows[LANDAU_MAX_Q_FACE], cols[LANDAU_MAX_Q_FACE];
1681:               for (PetscInt f = 0; f < Nb; ++f) {
1682:                 const PetscInt nr = coo_elem_point_offsets[glb_elem_idx][f + 1] - coo_elem_point_offsets[glb_elem_idx][f];
1683:                 if (nr == 1) rows[0] = Idxs[f];
1684:                 else {
1685:                   const PetscInt idx = -Idxs[f] - 1;
1686:                   for (PetscInt q = 0; q < nr; q++) rows[q] = maps[grid].c_maps[idx][q].gid;
1687:                 }
1688:                 for (PetscInt g = 0; g < Nb; ++g) {
1689:                   const PetscInt nc = coo_elem_point_offsets[glb_elem_idx][g + 1] - coo_elem_point_offsets[glb_elem_idx][g];
1690:                   if (nc == 1) cols[0] = Idxs[g];
1691:                   else {
1692:                     const PetscInt idx = -Idxs[g] - 1;
1693:                     for (PetscInt q = 0; q < nc; q++) cols[q] = maps[grid].c_maps[idx][q].gid;
1694:                   }
1695:                   const PetscInt idx0 = b_id * coo_elem_offsets[ncellsTot] + coo_elem_offsets[glb_elem_idx] + fieldA * fullNb2 + fullNb * coo_elem_point_offsets[glb_elem_idx][f] + nr * coo_elem_point_offsets[glb_elem_idx][g];
1696:                   for (PetscInt q = 0, idx = idx0; q < nr; q++) {
1697:                     for (PetscInt d = 0; d < nc; d++, idx++) {
1698:                       oor[idx] = rows[q] + moffset;
1699:                       ooc[idx] = cols[d] + moffset;
1700:                     }
1701:                   }
1702:                 }
1703:               }
1704:             }
1705:           } // cell
1706:         } // grid
1707:       } // batch
1708:       PetscCall(MatSetPreallocationCOO(ctx->J, ctx->SData_d.coo_size, oor, ooc));
1709:       PetscCall(PetscFree2(oor, ooc));
1710:     }
1711:     PetscCall(PetscFree(pointMaps));
1712:     PetscCall(PetscFree(elemMatrix));
1713:     PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &container));
1714:     PetscCall(PetscContainerSetPointer(container, (void *)maps));
1715:     PetscCall(PetscContainerSetCtxDestroy(container, LandauGPUMapsDestroy));
1716:     PetscCall(PetscObjectCompose((PetscObject)ctx->J, "assembly_maps", (PetscObject)container));
1717:     PetscCall(PetscContainerDestroy(&container));
1718:     PetscCall(PetscLogEventEnd(ctx->events[2], 0, 0, 0, 0));
1719:   } // end GPU assembly
1720:   { /* create static point data, Jacobian called first, only one vertex copy */
1721:     PetscReal *invJe, *ww, *xx, *yy, *zz = NULL, *invJ_a;
1722:     PetscInt   outer_ipidx, outer_ej, grid, nip_glb = 0;
1723:     PetscFE    fe;
1724:     PetscCall(PetscLogEventBegin(ctx->events[7], 0, 0, 0, 0));
1725:     PetscCall(PetscInfo(ctx->plex[0], "Initialize static data\n"));
1726:     for (PetscInt grid = 0; grid < ctx->num_grids; grid++) nip_glb += Nq * numCells[grid];
1727:     /* collect f data, first time is for Jacobian, but make mass now */
1728:     if (ctx->verbose != 0) {
1729:       PetscInt ncells = 0, N;
1730:       MatInfo  info;
1731:       PetscCall(MatGetInfo(ctx->J, MAT_LOCAL, &info));
1732:       PetscCall(MatGetSize(ctx->J, &N, NULL));
1733:       for (PetscInt grid = 0; grid < ctx->num_grids; grid++) ncells += numCells[grid];
1734:       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%d) %s %" PetscInt_FMT " IPs, %" PetscInt_FMT " cells total, Nb=%" PetscInt_FMT ", Nq=%" PetscInt_FMT ", dim=%" PetscInt_FMT ", Tab: Nb=%" PetscInt_FMT " Nf=%" PetscInt_FMT " Np=%" PetscInt_FMT " cdim=%" PetscInt_FMT " N=%" PetscInt_FMT " nnz= %" PetscInt_FMT "\n", 0, "FormLandau", nip_glb, ncells, Nb, Nq, dim, Nb,
1735:                             ctx->num_species, Nb, dim, N, (PetscInt)info.nz_used));
1736:     }
1737:     PetscCall(PetscMalloc4(nip_glb, &ww, nip_glb, &xx, nip_glb, &yy, nip_glb * dim * dim, &invJ_a));
1738:     if (dim == 3) PetscCall(PetscMalloc1(nip_glb, &zz));
1739:     if (ctx->use_energy_tensor_trick) {
1740:       PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, ctx->simplex, NULL, PETSC_DECIDE, &fe));
1741:       PetscCall(PetscObjectSetName((PetscObject)fe, "energy"));
1742:     }
1743:     /* init each grids static data - no batch */
1744:     for (grid = 0, outer_ipidx = 0, outer_ej = 0; grid < ctx->num_grids; grid++) { // OpenMP (once)
1745:       Vec          v2_2 = NULL;                                                    // projected function: v^2/2 for non-relativistic, gamma... for relativistic
1746:       PetscSection e_section;
1747:       DM           dmEnergy;
1748:       PetscInt     cStart, cEnd, ej;

1750:       PetscCall(DMPlexGetHeightStratum(ctx->plex[grid], 0, &cStart, &cEnd));
1751:       // prep energy trick, get v^2 / 2 vector
1752:       if (ctx->use_energy_tensor_trick) {
1753:         PetscErrorCode (*energyf[1])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *) = {ctx->use_relativistic_corrections ? gamma_m1_f : energy_f};
1754:         Vec        glob_v2;
1755:         PetscReal *c2_0[1], data[1] = {PetscSqr(C_0(ctx->v_0))};

1757:         PetscCall(DMClone(ctx->plex[grid], &dmEnergy));
1758:         PetscCall(PetscObjectSetName((PetscObject)dmEnergy, "energy"));
1759:         PetscCall(DMSetField(dmEnergy, 0, NULL, (PetscObject)fe));
1760:         PetscCall(DMCreateDS(dmEnergy));
1761:         PetscCall(DMGetLocalSection(dmEnergy, &e_section));
1762:         PetscCall(DMGetGlobalVector(dmEnergy, &glob_v2));
1763:         PetscCall(PetscObjectSetName((PetscObject)glob_v2, "trick"));
1764:         c2_0[0] = &data[0];
1765:         PetscCall(DMProjectFunction(dmEnergy, 0., energyf, (void **)c2_0, INSERT_ALL_VALUES, glob_v2));
1766:         PetscCall(DMGetLocalVector(dmEnergy, &v2_2));
1767:         PetscCall(VecZeroEntries(v2_2)); /* zero BCs so don't set */
1768:         PetscCall(DMGlobalToLocalBegin(dmEnergy, glob_v2, INSERT_VALUES, v2_2));
1769:         PetscCall(DMGlobalToLocalEnd(dmEnergy, glob_v2, INSERT_VALUES, v2_2));
1770:         PetscCall(DMViewFromOptions(dmEnergy, NULL, "-energy_dm_view"));
1771:         PetscCall(VecViewFromOptions(glob_v2, NULL, "-energy_vec_view"));
1772:         PetscCall(DMRestoreGlobalVector(dmEnergy, &glob_v2));
1773:       }
1774:       /* append part of the IP data for each grid */
1775:       for (ej = 0; ej < numCells[grid]; ++ej, ++outer_ej) {
1776:         PetscScalar *coefs = NULL;
1777:         PetscReal    vj[LANDAU_MAX_NQND * LANDAU_DIM], detJj[LANDAU_MAX_NQND], Jdummy[LANDAU_MAX_NQND * LANDAU_DIM * LANDAU_DIM], c0 = C_0(ctx->v_0), c02 = PetscSqr(c0);
1778:         invJe = invJ_a + outer_ej * Nq * dim * dim;
1779:         PetscCall(DMPlexComputeCellGeometryFEM(ctx->plex[grid], ej + cStart, quad, vj, Jdummy, invJe, detJj));
1780:         if (ctx->use_energy_tensor_trick) PetscCall(DMPlexVecGetClosure(dmEnergy, e_section, v2_2, ej + cStart, NULL, &coefs));
1781:         /* create static point data */
1782:         for (PetscInt qj = 0; qj < Nq; qj++, outer_ipidx++) {
1783:           const PetscInt   gidx = outer_ipidx;
1784:           const PetscReal *invJ = &invJe[qj * dim * dim];
1785:           ww[gidx]              = detJj[qj] * quadWeights[qj];
1786:           if (dim == 2) ww[gidx] *= vj[qj * dim + 0]; /* cylindrical coordinate, w/o 2pi */
1787:           // get xx, yy, zz
1788:           if (ctx->use_energy_tensor_trick) {
1789:             double                 refSpaceDer[3], eGradPhi[3];
1790:             const PetscReal *const DD = Tf[0]->T[1];
1791:             const PetscReal       *Dq = &DD[qj * Nb * dim];
1792:             for (PetscInt d = 0; d < 3; ++d) refSpaceDer[d] = eGradPhi[d] = 0.0;
1793:             for (PetscInt b = 0; b < Nb; ++b) {
1794:               for (PetscInt d = 0; d < dim; ++d) refSpaceDer[d] += Dq[b * dim + d] * PetscRealPart(coefs[b]);
1795:             }
1796:             xx[gidx] = 1e10;
1797:             if (ctx->use_relativistic_corrections) {
1798:               double dg2_c2 = 0;
1799:               //for (PetscInt d = 0; d < dim; ++d) refSpaceDer[d] *= c02;
1800:               for (PetscInt d = 0; d < dim; ++d) dg2_c2 += PetscSqr(refSpaceDer[d]);
1801:               dg2_c2 *= (double)c02;
1802:               if (dg2_c2 >= .999) {
1803:                 xx[gidx] = vj[qj * dim + 0]; /* coordinate */
1804:                 yy[gidx] = vj[qj * dim + 1];
1805:                 if (dim == 3) zz[gidx] = vj[qj * dim + 2];
1806:                 PetscCall(PetscPrintf(ctx->comm, "Error: %12.5e %" PetscInt_FMT ".%" PetscInt_FMT ") dg2/c02 = %12.5e x= %12.5e %12.5e %12.5e\n", (double)PetscSqrtReal(xx[gidx] * xx[gidx] + yy[gidx] * yy[gidx] + zz[gidx] * zz[gidx]), ej, qj, dg2_c2, (double)xx[gidx], (double)yy[gidx], (double)zz[gidx]));
1807:               } else {
1808:                 PetscReal fact = c02 / PetscSqrtReal(1. - dg2_c2);
1809:                 for (PetscInt d = 0; d < dim; ++d) refSpaceDer[d] *= fact;
1810:                 // could test with other point u' that (grad - grad') * U (refSpaceDer, refSpaceDer') == 0
1811:               }
1812:             }
1813:             if (xx[gidx] == 1e10) {
1814:               for (PetscInt d = 0; d < dim; ++d) {
1815:                 for (PetscInt e = 0; e < dim; ++e) eGradPhi[d] += invJ[e * dim + d] * refSpaceDer[e];
1816:               }
1817:               xx[gidx] = eGradPhi[0];
1818:               yy[gidx] = eGradPhi[1];
1819:               if (dim == 3) zz[gidx] = eGradPhi[2];
1820:             }
1821:           } else {
1822:             xx[gidx] = vj[qj * dim + 0]; /* coordinate */
1823:             yy[gidx] = vj[qj * dim + 1];
1824:             if (dim == 3) zz[gidx] = vj[qj * dim + 2];
1825:           }
1826:         } /* q */
1827:         if (ctx->use_energy_tensor_trick) PetscCall(DMPlexVecRestoreClosure(dmEnergy, e_section, v2_2, ej + cStart, NULL, &coefs));
1828:       } /* ej */
1829:       if (ctx->use_energy_tensor_trick) {
1830:         PetscCall(DMRestoreLocalVector(dmEnergy, &v2_2));
1831:         PetscCall(DMDestroy(&dmEnergy));
1832:       }
1833:     } /* grid */
1834:     if (ctx->use_energy_tensor_trick) PetscCall(PetscFEDestroy(&fe));
1835:     /* cache static data */
1836:     if (ctx->deviceType == LANDAU_KOKKOS) {
1837: #if defined(PETSC_HAVE_KOKKOS)
1838:       PetscCall(LandauKokkosStaticDataSet(ctx->plex[0], Nq, Nb, ctx->batch_sz, ctx->num_grids, numCells, ctx->species_offset, ctx->mat_offset, nu_alpha, nu_beta, invMass, (PetscReal *)ctx->lambdas, invJ_a, xx, yy, zz, ww, &ctx->SData_d));
1839:       /* free */
1840:       PetscCall(PetscFree4(ww, xx, yy, invJ_a));
1841:       if (dim == 3) PetscCall(PetscFree(zz));
1842: #else
1843:       SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONG, "-landau_device_type kokkos not built");
1844: #endif
1845:     } else {                                                                                                                                                                   /* CPU version, just copy in, only use part */
1846:       PetscReal *nu_alpha_p = (PetscReal *)ctx->SData_d.alpha, *nu_beta_p = (PetscReal *)ctx->SData_d.beta, *invMass_p = (PetscReal *)ctx->SData_d.invMass, *lambdas_p = NULL; // why set these ?
1847:       ctx->SData_d.w    = (void *)ww;
1848:       ctx->SData_d.x    = (void *)xx;
1849:       ctx->SData_d.y    = (void *)yy;
1850:       ctx->SData_d.z    = (void *)zz;
1851:       ctx->SData_d.invJ = (void *)invJ_a;
1852:       PetscCall(PetscMalloc4(ctx->num_species, &nu_alpha_p, ctx->num_species, &nu_beta_p, ctx->num_species, &invMass_p, LANDAU_MAX_GRIDS * LANDAU_MAX_GRIDS, &lambdas_p));
1853:       for (PetscInt ii = 0; ii < ctx->num_species; ii++) {
1854:         nu_alpha_p[ii] = nu_alpha[ii];
1855:         nu_beta_p[ii]  = nu_beta[ii];
1856:         invMass_p[ii]  = invMass[ii];
1857:       }
1858:       ctx->SData_d.alpha   = (void *)nu_alpha_p;
1859:       ctx->SData_d.beta    = (void *)nu_beta_p;
1860:       ctx->SData_d.invMass = (void *)invMass_p;
1861:       ctx->SData_d.lambdas = (void *)lambdas_p;
1862:       for (PetscInt grid = 0; grid < LANDAU_MAX_GRIDS; grid++) {
1863:         PetscReal (*lambdas)[LANDAU_MAX_GRIDS][LANDAU_MAX_GRIDS] = (PetscReal (*)[LANDAU_MAX_GRIDS][LANDAU_MAX_GRIDS])ctx->SData_d.lambdas;
1864:         for (PetscInt gridj = 0; gridj < LANDAU_MAX_GRIDS; gridj++) (*lambdas)[grid][gridj] = ctx->lambdas[grid][gridj];
1865:       }
1866:     }
1867:     PetscCall(PetscLogEventEnd(ctx->events[7], 0, 0, 0, 0));
1868:   } // initialize
1869:   PetscFunctionReturn(PETSC_SUCCESS);
1870: }

1872: /* < v, u > */
1873: static void g0_1(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
1874: {
1875:   g0[0] = 1.;
1876: }

1878: /* < v, u > */
1879: static void g0_fake(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
1880: {
1881:   static double ttt = 1e-12;
1882:   g0[0]             = ttt++;
1883: }

1885: /* < v, u > */
1886: static void g0_r(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
1887: {
1888:   g0[0] = 2. * PETSC_PI * x[0];
1889: }

1891: /*
1892:  LandauCreateJacobianMatrix - creates ctx->J with without real data. Hard to keep sparse.
1893:   - Like DMPlexLandauCreateMassMatrix. Should remove one and combine
1894:   - has old support for field major ordering
1895:  */
1896: static PetscErrorCode LandauCreateJacobianMatrix(MPI_Comm comm, Vec X, IS grid_batch_is_inv[LANDAU_MAX_GRIDS], LandauCtx *ctx)
1897: {
1898:   PetscInt *idxs = NULL;
1899:   Mat       subM[LANDAU_MAX_GRIDS];

1901:   PetscFunctionBegin;
1902:   if (!ctx->gpu_assembly) { /* we need GPU object with GPU assembly */
1903:     PetscFunctionReturn(PETSC_SUCCESS);
1904:   }
1905:   // get the RCM for this grid to separate out species into blocks -- create 'idxs' & 'ctx->batch_is' -- not used
1906:   if (ctx->gpu_assembly && ctx->jacobian_field_major_order) PetscCall(PetscMalloc1(ctx->mat_offset[ctx->num_grids] * ctx->batch_sz, &idxs));
1907:   for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
1908:     const PetscInt *values, n = ctx->mat_offset[grid + 1] - ctx->mat_offset[grid];
1909:     Mat             gMat;
1910:     DM              massDM;
1911:     PetscDS         prob;
1912:     Vec             tvec;
1913:     // get "mass" matrix for reordering
1914:     PetscCall(DMClone(ctx->plex[grid], &massDM));
1915:     PetscCall(DMCopyFields(ctx->plex[grid], PETSC_DETERMINE, PETSC_DETERMINE, massDM));
1916:     PetscCall(DMCreateDS(massDM));
1917:     PetscCall(DMGetDS(massDM, &prob));
1918:     for (PetscInt ix = 0, ii = ctx->species_offset[grid]; ii < ctx->species_offset[grid + 1]; ii++, ix++) PetscCall(PetscDSSetJacobian(prob, ix, ix, g0_fake, NULL, NULL, NULL));
1919:     PetscCall(PetscOptionsInsertString(NULL, "-dm_preallocate_only")); // this trick is need to both sparsify the matrix and avoid runtime error
1920:     PetscCall(DMCreateMatrix(massDM, &gMat));
1921:     PetscCall(PetscOptionsInsertString(NULL, "-dm_preallocate_only false"));
1922:     PetscCall(MatSetOption(gMat, MAT_STRUCTURALLY_SYMMETRIC, PETSC_TRUE));
1923:     PetscCall(MatSetOption(gMat, MAT_IGNORE_ZERO_ENTRIES, PETSC_TRUE));
1924:     PetscCall(DMCreateLocalVector(ctx->plex[grid], &tvec));
1925:     PetscCall(DMPlexSNESComputeJacobianFEM(massDM, tvec, gMat, gMat, ctx));
1926:     PetscCall(MatViewFromOptions(gMat, NULL, "-dm_landau_reorder_mat_view"));
1927:     PetscCall(DMDestroy(&massDM));
1928:     PetscCall(VecDestroy(&tvec));
1929:     subM[grid] = gMat;
1930:     if (ctx->gpu_assembly && ctx->jacobian_field_major_order) {
1931:       MatOrderingType rtype = MATORDERINGRCM;
1932:       IS              isrow, isicol;
1933:       PetscCall(MatGetOrdering(gMat, rtype, &isrow, &isicol));
1934:       PetscCall(ISInvertPermutation(isrow, PETSC_DECIDE, &grid_batch_is_inv[grid]));
1935:       PetscCall(ISGetIndices(isrow, &values));
1936:       for (PetscInt b_id = 0; b_id < ctx->batch_sz; b_id++) { // add batch size DMs for this species grid
1937: #if !defined(LANDAU_SPECIES_MAJOR)
1938:         PetscInt N = ctx->mat_offset[ctx->num_grids], n0 = ctx->mat_offset[grid] + b_id * N;
1939:         for (PetscInt ii = 0; ii < n; ++ii) idxs[n0 + ii] = values[ii] + n0;
1940: #else
1941:         PetscInt n0 = ctx->mat_offset[grid] * ctx->batch_sz + b_id * n;
1942:         for (PetscInt ii = 0; ii < n; ++ii) idxs[n0 + ii] = values[ii] + n0;
1943: #endif
1944:       }
1945:       PetscCall(ISRestoreIndices(isrow, &values));
1946:       PetscCall(ISDestroy(&isrow));
1947:       PetscCall(ISDestroy(&isicol));
1948:     }
1949:   }
1950:   if (ctx->gpu_assembly && ctx->jacobian_field_major_order) PetscCall(ISCreateGeneral(comm, ctx->mat_offset[ctx->num_grids] * ctx->batch_sz, idxs, PETSC_OWN_POINTER, &ctx->batch_is));
1951:   // get a block matrix
1952:   for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
1953:     Mat      B = subM[grid];
1954:     PetscInt nloc, nzl, *colbuf, row, COL_BF_SIZE = 1024;
1955:     PetscCall(PetscMalloc(sizeof(*colbuf) * COL_BF_SIZE, &colbuf));
1956:     PetscCall(MatGetSize(B, &nloc, NULL));
1957:     for (PetscInt b_id = 0; b_id < ctx->batch_sz; b_id++) {
1958:       const PetscInt     moffset = LAND_MOFFSET(b_id, grid, ctx->batch_sz, ctx->num_grids, ctx->mat_offset);
1959:       const PetscInt    *cols;
1960:       const PetscScalar *vals;
1961:       for (PetscInt i = 0; i < nloc; i++) {
1962:         PetscCall(MatGetRow(B, i, &nzl, NULL, NULL));
1963:         if (nzl > COL_BF_SIZE) {
1964:           PetscCall(PetscFree(colbuf));
1965:           PetscCall(PetscInfo(ctx->plex[grid], "Realloc buffer %" PetscInt_FMT " to %" PetscInt_FMT " (row size %" PetscInt_FMT ") \n", COL_BF_SIZE, 2 * COL_BF_SIZE, nzl));
1966:           COL_BF_SIZE = nzl;
1967:           PetscCall(PetscMalloc(sizeof(*colbuf) * COL_BF_SIZE, &colbuf));
1968:         }
1969:         PetscCall(MatGetRow(B, i, &nzl, &cols, &vals));
1970:         for (PetscInt j = 0; j < nzl; j++) colbuf[j] = cols[j] + moffset;
1971:         row = i + moffset;
1972:         PetscCall(MatSetValues(ctx->J, 1, &row, nzl, colbuf, vals, INSERT_VALUES));
1973:         PetscCall(MatRestoreRow(B, i, &nzl, &cols, &vals));
1974:       }
1975:     }
1976:     PetscCall(PetscFree(colbuf));
1977:   }
1978:   for (PetscInt grid = 0; grid < ctx->num_grids; grid++) PetscCall(MatDestroy(&subM[grid]));
1979:   PetscCall(MatAssemblyBegin(ctx->J, MAT_FINAL_ASSEMBLY));
1980:   PetscCall(MatAssemblyEnd(ctx->J, MAT_FINAL_ASSEMBLY));

1982:   // debug
1983:   PetscCall(MatViewFromOptions(ctx->J, NULL, "-dm_landau_mat_view"));
1984:   if (ctx->gpu_assembly && ctx->jacobian_field_major_order) {
1985:     Mat mat_block_order;
1986:     PetscCall(MatCreateSubMatrix(ctx->J, ctx->batch_is, ctx->batch_is, MAT_INITIAL_MATRIX, &mat_block_order)); // use MatPermute
1987:     PetscCall(MatViewFromOptions(mat_block_order, NULL, "-dm_landau_mat_view"));
1988:     PetscCall(MatDestroy(&mat_block_order));
1989:     PetscCall(VecScatterCreate(X, ctx->batch_is, X, NULL, &ctx->plex_batch));
1990:     PetscCall(VecDuplicate(X, &ctx->work_vec));
1991:   }
1992:   PetscFunctionReturn(PETSC_SUCCESS);
1993: }

1995: PetscErrorCode DMPlexLandauCreateMassMatrix(DM pack, Mat *Amat);
1996: /*@C
1997:   DMPlexLandauCreateVelocitySpace - Create a `DMPLEX` velocity space mesh

1999:   Collective

2001:   Input Parameters:
2002: + comm   - The MPI communicator
2003: . dim    - velocity space dimension (2 for axisymmetric, 3 for full 3X + 3V solver)
2004: - prefix - prefix for options (not tested)

2006:   Output Parameters:
2007: + pack - The `DM` object representing the mesh
2008: . X    - A vector (user destroys)
2009: - J    - Optional matrix (object destroys)

2011:   Level: beginner

2013: .seealso: `DMPlexCreate()`, `DMPlexLandauDestroyVelocitySpace()`
2014:  @*/
2015: PetscErrorCode DMPlexLandauCreateVelocitySpace(MPI_Comm comm, PetscInt dim, const char prefix[], Vec *X, Mat *J, DM *pack)
2016: {
2017:   LandauCtx *ctx;
2018:   Vec        Xsub[LANDAU_MAX_GRIDS];
2019:   IS         grid_batch_is_inv[LANDAU_MAX_GRIDS];

2021:   PetscFunctionBegin;
2022:   PetscCheck(dim == 2 || dim == 3, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Only 2D and 3D supported");
2023:   PetscCheck(LANDAU_DIM == dim, PETSC_COMM_SELF, PETSC_ERR_PLIB, "dim %" PetscInt_FMT " != LANDAU_DIM %d", dim, LANDAU_DIM);
2024:   PetscCall(PetscNew(&ctx));
2025:   ctx->comm = comm; /* used for diagnostics and global errors */
2026:   /* process options */
2027:   PetscCall(ProcessOptions(ctx, prefix));
2028:   if (dim == 2) ctx->use_relativistic_corrections = PETSC_FALSE;
2029:   /* Create Mesh */
2030:   PetscCall(DMCompositeCreate(PETSC_COMM_SELF, pack));
2031:   PetscCall(PetscLogEventBegin(ctx->events[13], 0, 0, 0, 0));
2032:   PetscCall(PetscLogEventBegin(ctx->events[15], 0, 0, 0, 0));
2033:   PetscCall(LandauDMCreateVMeshes(PETSC_COMM_SELF, dim, prefix, ctx, *pack)); // creates grids (Forest of AMR)
2034:   for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
2035:     /* create FEM */
2036:     PetscCall(SetupDS(ctx->plex[grid], dim, grid, ctx));
2037:     /* set initial state */
2038:     PetscCall(DMCreateGlobalVector(ctx->plex[grid], &Xsub[grid]));
2039:     PetscCall(PetscObjectSetName((PetscObject)Xsub[grid], "u_orig"));
2040:     /* initial static refinement, no solve */
2041:     PetscCall(LandauSetInitialCondition(ctx->plex[grid], Xsub[grid], grid, 0, 1, ctx));
2042:     /* forest refinement - forest goes in (if forest), plex comes out */
2043:     if (ctx->use_p4est) {
2044:       DM plex;
2045:       PetscCall(adapt(grid, ctx, &Xsub[grid])); // forest goes in, plex comes out
2046:       if (grid == 0) {
2047:         PetscCall(DMViewFromOptions(ctx->plex[grid], NULL, "-dm_landau_amr_dm_view")); // need to differentiate - todo
2048:         PetscCall(VecViewFromOptions(Xsub[grid], NULL, "-dm_landau_amr_vec_view"));
2049:       }
2050:       // convert to plex, all done with this level
2051:       PetscCall(DMConvert(ctx->plex[grid], DMPLEX, &plex));
2052:       PetscCall(DMDestroy(&ctx->plex[grid]));
2053:       ctx->plex[grid] = plex;
2054:     }
2055: #if !defined(LANDAU_SPECIES_MAJOR)
2056:     PetscCall(DMCompositeAddDM(*pack, ctx->plex[grid]));
2057: #else
2058:     for (PetscInt b_id = 0; b_id < ctx->batch_sz; b_id++) { // add batch size DMs for this species grid
2059:       PetscCall(DMCompositeAddDM(*pack, ctx->plex[grid]));
2060:     }
2061: #endif
2062:     PetscCall(DMSetApplicationContext(ctx->plex[grid], ctx));
2063:   }
2064: #if !defined(LANDAU_SPECIES_MAJOR)
2065:   // stack the batched DMs, could do it all here!!! b_id=0
2066:   for (PetscInt b_id = 1; b_id < ctx->batch_sz; b_id++) {
2067:     for (PetscInt grid = 0; grid < ctx->num_grids; grid++) PetscCall(DMCompositeAddDM(*pack, ctx->plex[grid]));
2068:   }
2069: #endif
2070:   // create ctx->mat_offset
2071:   ctx->mat_offset[0] = 0;
2072:   for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
2073:     PetscInt n;
2074:     PetscCall(VecGetLocalSize(Xsub[grid], &n));
2075:     ctx->mat_offset[grid + 1] = ctx->mat_offset[grid] + n;
2076:   }
2077:   // creat DM & Jac
2078:   PetscCall(DMSetApplicationContext(*pack, ctx));
2079:   PetscCall(PetscOptionsInsertString(NULL, "-dm_preallocate_only"));
2080:   PetscCall(DMCreateMatrix(*pack, &ctx->J));
2081:   PetscCall(PetscOptionsInsertString(NULL, "-dm_preallocate_only false"));
2082:   PetscCall(MatSetOption(ctx->J, MAT_STRUCTURALLY_SYMMETRIC, PETSC_TRUE));
2083:   PetscCall(MatSetOption(ctx->J, MAT_IGNORE_ZERO_ENTRIES, PETSC_TRUE));
2084:   PetscCall(PetscObjectSetName((PetscObject)ctx->J, "Jac"));
2085:   // construct initial conditions in X
2086:   PetscCall(DMCreateGlobalVector(*pack, X));
2087:   for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
2088:     PetscInt n;
2089:     PetscCall(VecGetLocalSize(Xsub[grid], &n));
2090:     for (PetscInt b_id = 0; b_id < ctx->batch_sz; b_id++) {
2091:       PetscScalar const *values;
2092:       const PetscInt     moffset = LAND_MOFFSET(b_id, grid, ctx->batch_sz, ctx->num_grids, ctx->mat_offset);
2093:       PetscCall(LandauSetInitialCondition(ctx->plex[grid], Xsub[grid], grid, b_id, ctx->batch_sz, ctx));
2094:       PetscCall(VecGetArrayRead(Xsub[grid], &values)); // Drop whole grid in Plex ordering
2095:       for (PetscInt i = 0, idx = moffset; i < n; i++, idx++) PetscCall(VecSetValue(*X, idx, values[i], INSERT_VALUES));
2096:       PetscCall(VecRestoreArrayRead(Xsub[grid], &values));
2097:     }
2098:   }
2099:   // cleanup
2100:   for (PetscInt grid = 0; grid < ctx->num_grids; grid++) PetscCall(VecDestroy(&Xsub[grid]));
2101:   /* check for correct matrix type */
2102:   if (ctx->gpu_assembly) { /* we need GPU object with GPU assembly */
2103:     PetscBool flg;
2104:     if (ctx->deviceType == LANDAU_KOKKOS) {
2105:       PetscCall(PetscObjectTypeCompareAny((PetscObject)ctx->J, &flg, MATSEQAIJKOKKOS, MATMPIAIJKOKKOS, MATAIJKOKKOS, ""));
2106: #if defined(PETSC_HAVE_KOKKOS)
2107:       PetscCheck(flg, ctx->comm, PETSC_ERR_ARG_WRONG, "must use '-dm_mat_type aijkokkos -dm_vec_type kokkos' for GPU assembly and Kokkos or use '-dm_landau_device_type cpu'");
2108: #else
2109:       PetscCheck(flg, ctx->comm, PETSC_ERR_ARG_WRONG, "must configure with '--download-kokkos-kernels' for GPU assembly and Kokkos or use '-dm_landau_device_type cpu'");
2110: #endif
2111:     }
2112:   }
2113:   PetscCall(PetscLogEventEnd(ctx->events[15], 0, 0, 0, 0));

2115:   // create field major ordering
2116:   ctx->work_vec   = NULL;
2117:   ctx->plex_batch = NULL;
2118:   ctx->batch_is   = NULL;
2119:   for (PetscInt i = 0; i < LANDAU_MAX_GRIDS; i++) grid_batch_is_inv[i] = NULL;
2120:   PetscCall(PetscLogEventBegin(ctx->events[12], 0, 0, 0, 0));
2121:   PetscCall(LandauCreateJacobianMatrix(comm, *X, grid_batch_is_inv, ctx));
2122:   PetscCall(PetscLogEventEnd(ctx->events[12], 0, 0, 0, 0));

2124:   // create AMR GPU assembly maps and static GPU data
2125:   PetscCall(CreateStaticData(dim, grid_batch_is_inv, ctx));

2127:   PetscCall(PetscLogEventEnd(ctx->events[13], 0, 0, 0, 0));

2129:   // create mass matrix
2130:   PetscCall(DMPlexLandauCreateMassMatrix(*pack, NULL));

2132:   if (J) *J = ctx->J;

2134:   if (ctx->gpu_assembly && ctx->jacobian_field_major_order) {
2135:     PetscContainer container;
2136:     // cache ctx for KSP with batch/field major Jacobian ordering -ksp_type gmres/etc -dm_landau_jacobian_field_major_order
2137:     PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &container));
2138:     PetscCall(PetscContainerSetPointer(container, (void *)ctx));
2139:     PetscCall(PetscObjectCompose((PetscObject)ctx->J, "LandauCtx", (PetscObject)container));
2140:     PetscCall(PetscContainerDestroy(&container));
2141:     // batch solvers need to map -- can batch solvers work
2142:     PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &container));
2143:     PetscCall(PetscContainerSetPointer(container, (void *)ctx->plex_batch));
2144:     PetscCall(PetscObjectCompose((PetscObject)ctx->J, "plex_batch_is", (PetscObject)container));
2145:     PetscCall(PetscContainerDestroy(&container));
2146:   }
2147:   // for batch solvers
2148:   {
2149:     PetscContainer container;
2150:     PetscInt      *pNf;
2151:     PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &container));
2152:     PetscCall(PetscMalloc1(sizeof(*pNf), &pNf));
2153:     *pNf = ctx->batch_sz;
2154:     PetscCall(PetscContainerSetPointer(container, (void *)pNf));
2155:     PetscCall(PetscContainerSetCtxDestroy(container, PetscCtxDestroyDefault));
2156:     PetscCall(PetscObjectCompose((PetscObject)ctx->J, "batch size", (PetscObject)container));
2157:     PetscCall(PetscContainerDestroy(&container));
2158:   }
2159:   PetscFunctionReturn(PETSC_SUCCESS);
2160: }

2162: /*@C
2163:   DMPlexLandauAccess - Access to the distribution function with user callback

2165:   Collective

2167:   Input Parameters:
2168: + pack     - the `DMCOMPOSITE`
2169: . func     - call back function
2170: - user_ctx - user context

2172:   Input/Output Parameter:
2173: . X - Vector to data to

2175:   Level: advanced

2177: .seealso: `DMPlexLandauCreateVelocitySpace()`
2178:  @*/
2179: PetscErrorCode DMPlexLandauAccess(DM pack, Vec X, PetscErrorCode (*func)(DM, Vec, PetscInt, PetscInt, PetscInt, void *), void *user_ctx)
2180: {
2181:   LandauCtx *ctx;

2183:   PetscFunctionBegin;
2184:   PetscCall(DMGetApplicationContext(pack, &ctx)); // uses ctx->num_grids; ctx->plex[grid]; ctx->batch_sz; ctx->mat_offset
2185:   for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
2186:     PetscInt dim, n;
2187:     PetscCall(DMGetDimension(pack, &dim));
2188:     for (PetscInt sp = ctx->species_offset[grid], i0 = 0; sp < ctx->species_offset[grid + 1]; sp++, i0++) {
2189:       Vec      vec;
2190:       PetscInt vf[1] = {i0};
2191:       IS       vis;
2192:       DM       vdm;
2193:       PetscCall(DMCreateSubDM(ctx->plex[grid], 1, vf, &vis, &vdm));
2194:       PetscCall(DMSetApplicationContext(vdm, ctx)); // the user might want this
2195:       PetscCall(DMCreateGlobalVector(vdm, &vec));
2196:       PetscCall(VecGetSize(vec, &n));
2197:       for (PetscInt b_id = 0; b_id < ctx->batch_sz; b_id++) {
2198:         const PetscInt moffset = LAND_MOFFSET(b_id, grid, ctx->batch_sz, ctx->num_grids, ctx->mat_offset);
2199:         PetscCall(VecZeroEntries(vec));
2200:         /* Add your data with 'dm' for species 'sp' to 'vec' */
2201:         PetscCall(func(vdm, vec, i0, grid, b_id, user_ctx));
2202:         /* add to global */
2203:         PetscScalar const *values;
2204:         const PetscInt    *offsets;
2205:         PetscCall(VecGetArrayRead(vec, &values));
2206:         PetscCall(ISGetIndices(vis, &offsets));
2207:         for (PetscInt i = 0; i < n; i++) PetscCall(VecSetValue(X, moffset + offsets[i], values[i], ADD_VALUES));
2208:         PetscCall(VecRestoreArrayRead(vec, &values));
2209:         PetscCall(ISRestoreIndices(vis, &offsets));
2210:       } // batch
2211:       PetscCall(VecDestroy(&vec));
2212:       PetscCall(ISDestroy(&vis));
2213:       PetscCall(DMDestroy(&vdm));
2214:     }
2215:   } // grid
2216:   PetscFunctionReturn(PETSC_SUCCESS);
2217: }

2219: /*@
2220:   DMPlexLandauDestroyVelocitySpace - Destroy a `DMPLEX` velocity space mesh

2222:   Collective

2224:   Input/Output Parameters:
2225: . dm - the `DM` to destroy

2227:   Level: beginner

2229: .seealso: `DMPlexLandauCreateVelocitySpace()`
2230:  @*/
2231: PetscErrorCode DMPlexLandauDestroyVelocitySpace(DM *dm)
2232: {
2233:   LandauCtx *ctx;

2235:   PetscFunctionBegin;
2236:   PetscCall(DMGetApplicationContext(*dm, &ctx));
2237:   PetscCall(MatDestroy(&ctx->M));
2238:   PetscCall(MatDestroy(&ctx->J));
2239:   for (PetscInt ii = 0; ii < ctx->num_species; ii++) PetscCall(PetscFEDestroy(&ctx->fe[ii]));
2240:   PetscCall(ISDestroy(&ctx->batch_is));
2241:   PetscCall(VecDestroy(&ctx->work_vec));
2242:   PetscCall(VecScatterDestroy(&ctx->plex_batch));
2243:   if (ctx->deviceType == LANDAU_KOKKOS) {
2244: #if defined(PETSC_HAVE_KOKKOS)
2245:     PetscCall(LandauKokkosStaticDataClear(&ctx->SData_d));
2246: #else
2247:     SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONG, "-landau_device_type %s not built", "kokkos");
2248: #endif
2249:   } else {
2250:     if (ctx->SData_d.x) { /* in a CPU run */
2251:       PetscReal *invJ = (PetscReal *)ctx->SData_d.invJ, *xx = (PetscReal *)ctx->SData_d.x, *yy = (PetscReal *)ctx->SData_d.y, *zz = (PetscReal *)ctx->SData_d.z, *ww = (PetscReal *)ctx->SData_d.w;
2252:       LandauIdx *coo_elem_offsets = (LandauIdx *)ctx->SData_d.coo_elem_offsets, *coo_elem_fullNb = (LandauIdx *)ctx->SData_d.coo_elem_fullNb, (*coo_elem_point_offsets)[LANDAU_MAX_NQND + 1] = (LandauIdx(*)[LANDAU_MAX_NQND + 1]) ctx->SData_d.coo_elem_point_offsets;
2253:       PetscCall(PetscFree4(ww, xx, yy, invJ));
2254:       if (zz) PetscCall(PetscFree(zz));
2255:       if (coo_elem_offsets) PetscCall(PetscFree3(coo_elem_offsets, coo_elem_fullNb, coo_elem_point_offsets)); // could be NULL
2256:       PetscCall(PetscFree4(ctx->SData_d.alpha, ctx->SData_d.beta, ctx->SData_d.invMass, ctx->SData_d.lambdas));
2257:     }
2258:   }

2260:   if (ctx->times[LANDAU_MATRIX_TOTAL] > 0) { // OMP timings
2261:     PetscCall(PetscPrintf(ctx->comm, "TSStep               N  1.0 %10.3e\n", ctx->times[LANDAU_EX2_TSSOLVE]));
2262:     PetscCall(PetscPrintf(ctx->comm, "2:           Solve:  %10.3e with %" PetscInt_FMT " threads\n", ctx->times[LANDAU_EX2_TSSOLVE] - ctx->times[LANDAU_MATRIX_TOTAL], ctx->batch_sz));
2263:     PetscCall(PetscPrintf(ctx->comm, "3:          Landau:  %10.3e\n", ctx->times[LANDAU_MATRIX_TOTAL]));
2264:     PetscCall(PetscPrintf(ctx->comm, "Landau Jacobian       %" PetscInt_FMT " 1.0 %10.3e\n", (PetscInt)ctx->times[LANDAU_JACOBIAN_COUNT], ctx->times[LANDAU_JACOBIAN]));
2265:     PetscCall(PetscPrintf(ctx->comm, "Landau Operator       N 1.0  %10.3e\n", ctx->times[LANDAU_OPERATOR]));
2266:     PetscCall(PetscPrintf(ctx->comm, "Landau Mass           N 1.0  %10.3e\n", ctx->times[LANDAU_MASS]));
2267:     PetscCall(PetscPrintf(ctx->comm, " Jac-f-df (GPU)       N 1.0  %10.3e\n", ctx->times[LANDAU_F_DF]));
2268:     PetscCall(PetscPrintf(ctx->comm, " Kernel (GPU)         N 1.0  %10.3e\n", ctx->times[LANDAU_KERNEL]));
2269:     PetscCall(PetscPrintf(ctx->comm, "MatLUFactorNum        X 1.0 %10.3e\n", ctx->times[KSP_FACTOR]));
2270:     PetscCall(PetscPrintf(ctx->comm, "MatSolve              X 1.0 %10.3e\n", ctx->times[KSP_SOLVE]));
2271:   }
2272:   for (PetscInt grid = 0; grid < ctx->num_grids; grid++) PetscCall(DMDestroy(&ctx->plex[grid]));
2273:   PetscCall(PetscFree(ctx));
2274:   PetscCall(DMDestroy(dm));
2275:   PetscFunctionReturn(PETSC_SUCCESS);
2276: }

2278: /* < v, ru > */
2279: static void f0_s_den(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0)
2280: {
2281:   PetscInt ii = (PetscInt)PetscRealPart(constants[0]);
2282:   f0[0]       = u[ii];
2283: }

2285: /* < v, ru > */
2286: static void f0_s_mom(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0)
2287: {
2288:   PetscInt ii = (PetscInt)PetscRealPart(constants[0]), jj = (PetscInt)PetscRealPart(constants[1]);
2289:   f0[0] = x[jj] * u[ii]; /* x momentum */
2290: }

2292: static void f0_s_v2(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0)
2293: {
2294:   PetscInt i, ii = (PetscInt)PetscRealPart(constants[0]);
2295:   double   tmp1 = 0.;
2296:   for (i = 0; i < dim; ++i) tmp1 += x[i] * x[i];
2297:   f0[0] = tmp1 * u[ii];
2298: }

2300: static PetscErrorCode gamma_n_f(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *actx)
2301: {
2302:   const PetscReal *c2_0_arr = ((PetscReal *)actx);
2303:   const PetscReal  c02      = c2_0_arr[0];

2305:   PetscFunctionBegin;
2306:   for (PetscInt s = 0; s < Nf; s++) {
2307:     PetscReal tmp1 = 0.;
2308:     for (PetscInt i = 0; i < dim; ++i) tmp1 += x[i] * x[i];
2309: #if defined(PETSC_USE_DEBUG)
2310:     u[s] = PetscSqrtReal(1. + tmp1 / c02); //  u[0] = PetscSqrtReal(1. + xx);
2311: #else
2312:     {
2313:       PetscReal xx = tmp1 / c02;
2314:       u[s]         = xx / (PetscSqrtReal(1. + xx) + 1.); // better conditioned = xx/(PetscSqrtReal(1. + xx) + 1.)
2315:     }
2316: #endif
2317:   }
2318:   PetscFunctionReturn(PETSC_SUCCESS);
2319: }

2321: /* < v, ru > */
2322: static void f0_s_rden(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0)
2323: {
2324:   PetscInt ii = (PetscInt)PetscRealPart(constants[0]);
2325:   f0[0]       = 2. * PETSC_PI * x[0] * u[ii];
2326: }

2328: /* < v, ru > */
2329: static void f0_s_rmom(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0)
2330: {
2331:   PetscInt ii = (PetscInt)PetscRealPart(constants[0]);
2332:   f0[0]       = 2. * PETSC_PI * x[0] * x[1] * u[ii];
2333: }

2335: static void f0_s_rv2(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar *f0)
2336: {
2337:   PetscInt ii = (PetscInt)PetscRealPart(constants[0]);
2338:   f0[0]       = 2. * PETSC_PI * x[0] * (x[0] * x[0] + x[1] * x[1]) * u[ii];
2339: }

2341: /*@
2342:   DMPlexLandauPrintNorms - collects moments and prints them

2344:   Collective

2346:   Input Parameters:
2347: + X     - the state
2348: - stepi - current step to print

2350:   Level: beginner

2352: .seealso: `DMPlexLandauCreateVelocitySpace()`
2353:  @*/
2354: PetscErrorCode DMPlexLandauPrintNorms(Vec X, PetscInt stepi)
2355: {
2356:   LandauCtx  *ctx;
2357:   PetscDS     prob;
2358:   DM          pack;
2359:   PetscInt    cStart, cEnd, dim, ii, i0, nDMs;
2360:   PetscScalar xmomentumtot = 0, ymomentumtot = 0, zmomentumtot = 0, energytot = 0, densitytot = 0, tt[LANDAU_MAX_SPECIES];
2361:   PetscScalar xmomentum[LANDAU_MAX_SPECIES], ymomentum[LANDAU_MAX_SPECIES], zmomentum[LANDAU_MAX_SPECIES], energy[LANDAU_MAX_SPECIES], density[LANDAU_MAX_SPECIES];
2362:   Vec        *globXArray;

2364:   PetscFunctionBegin;
2365:   PetscCall(VecGetDM(X, &pack));
2366:   PetscCheck(pack, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Vector has no DM");
2367:   PetscCall(DMGetDimension(pack, &dim));
2368:   PetscCheck(dim == 2 || dim == 3, PETSC_COMM_SELF, PETSC_ERR_PLIB, "dim %" PetscInt_FMT " not in [2,3]", dim);
2369:   PetscCall(DMGetApplicationContext(pack, &ctx));
2370:   PetscCheck(ctx, PETSC_COMM_SELF, PETSC_ERR_PLIB, "no context");
2371:   /* print momentum and energy */
2372:   PetscCall(DMCompositeGetNumberDM(pack, &nDMs));
2373:   PetscCheck(nDMs == ctx->num_grids * ctx->batch_sz, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "#DM wrong %" PetscInt_FMT " %" PetscInt_FMT, nDMs, ctx->num_grids * ctx->batch_sz);
2374:   PetscCall(PetscMalloc(sizeof(*globXArray) * nDMs, &globXArray));
2375:   PetscCall(DMCompositeGetAccessArray(pack, X, nDMs, NULL, globXArray));
2376:   for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
2377:     Vec Xloc = globXArray[LAND_PACK_IDX(ctx->batch_view_idx, grid)];
2378:     PetscCall(DMGetDS(ctx->plex[grid], &prob));
2379:     for (ii = ctx->species_offset[grid], i0 = 0; ii < ctx->species_offset[grid + 1]; ii++, i0++) {
2380:       PetscScalar user[2] = {(PetscScalar)i0, ctx->charges[ii]};
2381:       PetscCall(PetscDSSetConstants(prob, 2, user));
2382:       if (dim == 2) { /* 2/3X + 3V (cylindrical coordinates) */
2383:         PetscCall(PetscDSSetObjective(prob, 0, &f0_s_rden));
2384:         PetscCall(DMPlexComputeIntegralFEM(ctx->plex[grid], Xloc, tt, ctx));
2385:         density[ii] = tt[0] * ctx->n_0 * ctx->charges[ii];
2386:         PetscCall(PetscDSSetObjective(prob, 0, &f0_s_rmom));
2387:         PetscCall(DMPlexComputeIntegralFEM(ctx->plex[grid], Xloc, tt, ctx));
2388:         zmomentum[ii] = tt[0] * ctx->n_0 * ctx->v_0 * ctx->masses[ii];
2389:         PetscCall(PetscDSSetObjective(prob, 0, &f0_s_rv2));
2390:         PetscCall(DMPlexComputeIntegralFEM(ctx->plex[grid], Xloc, tt, ctx));
2391:         energy[ii] = tt[0] * 0.5 * ctx->n_0 * ctx->v_0 * ctx->v_0 * ctx->masses[ii];
2392:         zmomentumtot += zmomentum[ii];
2393:         energytot += energy[ii];
2394:         densitytot += density[ii];
2395:         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%3" PetscInt_FMT ") species-%" PetscInt_FMT ": charge density= %20.13e z-momentum= %20.13e energy= %20.13e", stepi, ii, (double)PetscRealPart(density[ii]), (double)PetscRealPart(zmomentum[ii]), (double)PetscRealPart(energy[ii])));
2396:       } else { /* 2/3Xloc + 3V */
2397:         PetscCall(PetscDSSetObjective(prob, 0, &f0_s_den));
2398:         PetscCall(DMPlexComputeIntegralFEM(ctx->plex[grid], Xloc, tt, ctx));
2399:         density[ii] = tt[0] * ctx->n_0 * ctx->charges[ii];
2400:         PetscCall(PetscDSSetObjective(prob, 0, &f0_s_mom));
2401:         user[1] = 0;
2402:         PetscCall(PetscDSSetConstants(prob, 2, user));
2403:         PetscCall(DMPlexComputeIntegralFEM(ctx->plex[grid], Xloc, tt, ctx));
2404:         xmomentum[ii] = tt[0] * ctx->n_0 * ctx->v_0 * ctx->masses[ii];
2405:         user[1]       = 1;
2406:         PetscCall(PetscDSSetConstants(prob, 2, user));
2407:         PetscCall(DMPlexComputeIntegralFEM(ctx->plex[grid], Xloc, tt, ctx));
2408:         ymomentum[ii] = tt[0] * ctx->n_0 * ctx->v_0 * ctx->masses[ii];
2409:         user[1]       = 2;
2410:         PetscCall(PetscDSSetConstants(prob, 2, user));
2411:         PetscCall(DMPlexComputeIntegralFEM(ctx->plex[grid], Xloc, tt, ctx));
2412:         zmomentum[ii] = tt[0] * ctx->n_0 * ctx->v_0 * ctx->masses[ii];
2413:         if (ctx->use_relativistic_corrections) {
2414:           /* gamma * M * f */
2415:           if (ii == 0 && grid == 0) { // do all at once
2416:             Vec Mf, globGamma, *globMfArray, *globGammaArray;
2417:             PetscErrorCode (*gammaf[1])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *) = {gamma_n_f};
2418:             PetscReal *c2_0[1], data[1];

2420:             PetscCall(VecDuplicate(X, &globGamma));
2421:             PetscCall(VecDuplicate(X, &Mf));
2422:             PetscCall(PetscMalloc(sizeof(*globMfArray) * nDMs, &globMfArray));
2423:             PetscCall(PetscMalloc(sizeof(*globMfArray) * nDMs, &globGammaArray));
2424:             /* M * f */
2425:             PetscCall(MatMult(ctx->M, X, Mf));
2426:             /* gamma */
2427:             PetscCall(DMCompositeGetAccessArray(pack, globGamma, nDMs, NULL, globGammaArray));
2428:             for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // yes a grid loop in a grid loop to print nice, need to fix for batching
2429:               Vec v1  = globGammaArray[LAND_PACK_IDX(ctx->batch_view_idx, grid)];
2430:               data[0] = PetscSqr(C_0(ctx->v_0));
2431:               c2_0[0] = &data[0];
2432:               PetscCall(DMProjectFunction(ctx->plex[grid], 0., gammaf, (void **)c2_0, INSERT_ALL_VALUES, v1));
2433:             }
2434:             PetscCall(DMCompositeRestoreAccessArray(pack, globGamma, nDMs, NULL, globGammaArray));
2435:             /* gamma * Mf */
2436:             PetscCall(DMCompositeGetAccessArray(pack, globGamma, nDMs, NULL, globGammaArray));
2437:             PetscCall(DMCompositeGetAccessArray(pack, Mf, nDMs, NULL, globMfArray));
2438:             for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // yes a grid loop in a grid loop to print nice
2439:               PetscInt Nf    = ctx->species_offset[grid + 1] - ctx->species_offset[grid], N, bs;
2440:               Vec      Mfsub = globMfArray[LAND_PACK_IDX(ctx->batch_view_idx, grid)], Gsub = globGammaArray[LAND_PACK_IDX(ctx->batch_view_idx, grid)], v1, v2;
2441:               // get each component
2442:               PetscCall(VecGetSize(Mfsub, &N));
2443:               PetscCall(VecCreate(ctx->comm, &v1));
2444:               PetscCall(VecSetSizes(v1, PETSC_DECIDE, N / Nf));
2445:               PetscCall(VecCreate(ctx->comm, &v2));
2446:               PetscCall(VecSetSizes(v2, PETSC_DECIDE, N / Nf));
2447:               PetscCall(VecSetFromOptions(v1)); // ???
2448:               PetscCall(VecSetFromOptions(v2));
2449:               // get each component
2450:               PetscCall(VecGetBlockSize(Gsub, &bs));
2451:               PetscCheck(bs == Nf, PETSC_COMM_SELF, PETSC_ERR_PLIB, "bs %" PetscInt_FMT " != num_species %" PetscInt_FMT " in Gsub", bs, Nf);
2452:               PetscCall(VecGetBlockSize(Mfsub, &bs));
2453:               PetscCheck(bs == Nf, PETSC_COMM_SELF, PETSC_ERR_PLIB, "bs %" PetscInt_FMT " != num_species %" PetscInt_FMT, bs, Nf);
2454:               for (PetscInt i = 0, ix = ctx->species_offset[grid]; i < Nf; i++, ix++) {
2455:                 PetscScalar val;
2456:                 PetscCall(VecStrideGather(Gsub, i, v1, INSERT_VALUES)); // this is not right -- TODO
2457:                 PetscCall(VecStrideGather(Mfsub, i, v2, INSERT_VALUES));
2458:                 PetscCall(VecDot(v1, v2, &val));
2459:                 energy[ix] = PetscRealPart(val) * ctx->n_0 * ctx->v_0 * ctx->v_0 * ctx->masses[ix];
2460:               }
2461:               PetscCall(VecDestroy(&v1));
2462:               PetscCall(VecDestroy(&v2));
2463:             } /* grids */
2464:             PetscCall(DMCompositeRestoreAccessArray(pack, globGamma, nDMs, NULL, globGammaArray));
2465:             PetscCall(DMCompositeRestoreAccessArray(pack, Mf, nDMs, NULL, globMfArray));
2466:             PetscCall(PetscFree(globGammaArray));
2467:             PetscCall(PetscFree(globMfArray));
2468:             PetscCall(VecDestroy(&globGamma));
2469:             PetscCall(VecDestroy(&Mf));
2470:           }
2471:         } else {
2472:           PetscCall(PetscDSSetObjective(prob, 0, &f0_s_v2));
2473:           PetscCall(DMPlexComputeIntegralFEM(ctx->plex[grid], Xloc, tt, ctx));
2474:           energy[ii] = 0.5 * tt[0] * ctx->n_0 * ctx->v_0 * ctx->v_0 * ctx->masses[ii];
2475:         }
2476:         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%3" PetscInt_FMT ") species %" PetscInt_FMT ": density=%20.13e, x-momentum=%20.13e, y-momentum=%20.13e, z-momentum=%20.13e, energy=%21.13e", stepi, ii, (double)PetscRealPart(density[ii]), (double)PetscRealPart(xmomentum[ii]), (double)PetscRealPart(ymomentum[ii]), (double)PetscRealPart(zmomentum[ii]), (double)PetscRealPart(energy[ii])));
2477:         xmomentumtot += xmomentum[ii];
2478:         ymomentumtot += ymomentum[ii];
2479:         zmomentumtot += zmomentum[ii];
2480:         energytot += energy[ii];
2481:         densitytot += density[ii];
2482:       }
2483:       if (ctx->num_species > 1) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n"));
2484:     }
2485:   }
2486:   PetscCall(DMCompositeRestoreAccessArray(pack, X, nDMs, NULL, globXArray));
2487:   PetscCall(PetscFree(globXArray));
2488:   /* totals */
2489:   PetscCall(DMPlexGetHeightStratum(ctx->plex[0], 0, &cStart, &cEnd));
2490:   if (ctx->num_species > 1) {
2491:     if (dim == 2) {
2492:       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\t%3" PetscInt_FMT ") Total: charge density=%21.13e, momentum=%21.13e, energy=%21.13e (m_i[0]/m_e = %g, %" PetscInt_FMT " cells on electron grid)", stepi, (double)PetscRealPart(densitytot), (double)PetscRealPart(zmomentumtot), (double)PetscRealPart(energytot),
2493:                             (double)(ctx->masses[1] / ctx->masses[0]), cEnd - cStart));
2494:     } else {
2495:       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\t%3" PetscInt_FMT ") 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, %" PetscInt_FMT " cells)", stepi, (double)PetscRealPart(densitytot), (double)PetscRealPart(xmomentumtot), (double)PetscRealPart(ymomentumtot), (double)PetscRealPart(zmomentumtot), (double)PetscRealPart(energytot),
2496:                             (double)(ctx->masses[1] / ctx->masses[0]), cEnd - cStart));
2497:     }
2498:   } else PetscCall(PetscPrintf(PETSC_COMM_WORLD, " -- %" PetscInt_FMT " cells", cEnd - cStart));
2499:   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n"));
2500:   PetscFunctionReturn(PETSC_SUCCESS);
2501: }

2503: /*@
2504:   DMPlexLandauCreateMassMatrix - Create mass matrix for Landau in Plex space (not field major order of Jacobian)
2505:   - puts mass matrix into ctx->M

2507:   Collective

2509:   Input Parameter:
2510: . pack - the `DM` object. Puts matrix in Landau context M field

2512:   Output Parameter:
2513: . Amat - The mass matrix (optional), mass matrix is added to the `DM` context

2515:   Level: beginner

2517: .seealso: `DMPlexLandauCreateVelocitySpace()`
2518:  @*/
2519: PetscErrorCode DMPlexLandauCreateMassMatrix(DM pack, Mat *Amat)
2520: {
2521:   DM         mass_pack, massDM[LANDAU_MAX_GRIDS];
2522:   PetscDS    prob;
2523:   PetscInt   ii, dim, N1 = 1, N2;
2524:   LandauCtx *ctx;
2525:   Mat        packM, subM[LANDAU_MAX_GRIDS];

2527:   PetscFunctionBegin;
2529:   if (Amat) PetscAssertPointer(Amat, 2);
2530:   PetscCall(DMGetApplicationContext(pack, &ctx));
2531:   PetscCheck(ctx, PETSC_COMM_SELF, PETSC_ERR_PLIB, "no context");
2532:   PetscCall(PetscLogEventBegin(ctx->events[14], 0, 0, 0, 0));
2533:   PetscCall(DMGetDimension(pack, &dim));
2534:   PetscCall(DMCompositeCreate(PetscObjectComm((PetscObject)pack), &mass_pack));
2535:   /* create pack mass matrix */
2536:   for (PetscInt grid = 0, ix = 0; grid < ctx->num_grids; grid++) {
2537:     PetscCall(DMClone(ctx->plex[grid], &massDM[grid]));
2538:     PetscCall(DMCopyFields(ctx->plex[grid], PETSC_DETERMINE, PETSC_DETERMINE, massDM[grid]));
2539:     PetscCall(DMCreateDS(massDM[grid]));
2540:     PetscCall(DMGetDS(massDM[grid], &prob));
2541:     for (ix = 0, ii = ctx->species_offset[grid]; ii < ctx->species_offset[grid + 1]; ii++, ix++) {
2542:       if (dim == 3) PetscCall(PetscDSSetJacobian(prob, ix, ix, g0_1, NULL, NULL, NULL));
2543:       else PetscCall(PetscDSSetJacobian(prob, ix, ix, g0_r, NULL, NULL, NULL));
2544:     }
2545: #if !defined(LANDAU_SPECIES_MAJOR)
2546:     PetscCall(DMCompositeAddDM(mass_pack, massDM[grid]));
2547: #else
2548:     for (PetscInt b_id = 0; b_id < ctx->batch_sz; b_id++) { // add batch size DMs for this species grid
2549:       PetscCall(DMCompositeAddDM(mass_pack, massDM[grid]));
2550:     }
2551: #endif
2552:     PetscCall(DMCreateMatrix(massDM[grid], &subM[grid]));
2553:   }
2554: #if !defined(LANDAU_SPECIES_MAJOR)
2555:   // stack the batched DMs
2556:   for (PetscInt b_id = 1; b_id < ctx->batch_sz; b_id++) {
2557:     for (PetscInt grid = 0; grid < ctx->num_grids; grid++) PetscCall(DMCompositeAddDM(mass_pack, massDM[grid]));
2558:   }
2559: #endif
2560:   PetscCall(PetscOptionsInsertString(NULL, "-dm_preallocate_only"));
2561:   PetscCall(DMCreateMatrix(mass_pack, &packM));
2562:   PetscCall(PetscOptionsInsertString(NULL, "-dm_preallocate_only false"));
2563:   PetscCall(MatSetOption(packM, MAT_STRUCTURALLY_SYMMETRIC, PETSC_TRUE));
2564:   PetscCall(MatSetOption(packM, MAT_IGNORE_ZERO_ENTRIES, PETSC_TRUE));
2565:   PetscCall(DMDestroy(&mass_pack));
2566:   /* make mass matrix for each block */
2567:   for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
2568:     Vec locX;
2569:     DM  plex = massDM[grid];
2570:     PetscCall(DMGetLocalVector(plex, &locX));
2571:     /* Mass matrix is independent of the input, so no need to fill locX */
2572:     PetscCall(DMPlexSNESComputeJacobianFEM(plex, locX, subM[grid], subM[grid], ctx));
2573:     PetscCall(DMRestoreLocalVector(plex, &locX));
2574:     PetscCall(DMDestroy(&massDM[grid]));
2575:   }
2576:   PetscCall(MatGetSize(ctx->J, &N1, NULL));
2577:   PetscCall(MatGetSize(packM, &N2, NULL));
2578:   PetscCheck(N1 == N2, PetscObjectComm((PetscObject)pack), PETSC_ERR_PLIB, "Incorrect matrix sizes: |Jacobian| = %" PetscInt_FMT ", |Mass|=%" PetscInt_FMT, N1, N2);
2579:   /* assemble block diagonals */
2580:   for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
2581:     Mat      B = subM[grid];
2582:     PetscInt nloc, nzl, *colbuf, COL_BF_SIZE = 1024, row;
2583:     PetscCall(PetscMalloc(sizeof(*colbuf) * COL_BF_SIZE, &colbuf));
2584:     PetscCall(MatGetSize(B, &nloc, NULL));
2585:     for (PetscInt b_id = 0; b_id < ctx->batch_sz; b_id++) {
2586:       const PetscInt     moffset = LAND_MOFFSET(b_id, grid, ctx->batch_sz, ctx->num_grids, ctx->mat_offset);
2587:       const PetscInt    *cols;
2588:       const PetscScalar *vals;
2589:       for (PetscInt i = 0; i < nloc; i++) {
2590:         PetscCall(MatGetRow(B, i, &nzl, NULL, NULL));
2591:         if (nzl > COL_BF_SIZE) {
2592:           PetscCall(PetscFree(colbuf));
2593:           PetscCall(PetscInfo(pack, "Realloc buffer %" PetscInt_FMT " to %" PetscInt_FMT " (row size %" PetscInt_FMT ") \n", COL_BF_SIZE, 2 * COL_BF_SIZE, nzl));
2594:           COL_BF_SIZE = nzl;
2595:           PetscCall(PetscMalloc(sizeof(*colbuf) * COL_BF_SIZE, &colbuf));
2596:         }
2597:         PetscCall(MatGetRow(B, i, &nzl, &cols, &vals));
2598:         for (PetscInt j = 0; j < nzl; j++) colbuf[j] = cols[j] + moffset;
2599:         row = i + moffset;
2600:         PetscCall(MatSetValues(packM, 1, &row, nzl, colbuf, vals, INSERT_VALUES));
2601:         PetscCall(MatRestoreRow(B, i, &nzl, &cols, &vals));
2602:       }
2603:     }
2604:     PetscCall(PetscFree(colbuf));
2605:   }
2606:   // cleanup
2607:   for (PetscInt grid = 0; grid < ctx->num_grids; grid++) PetscCall(MatDestroy(&subM[grid]));
2608:   PetscCall(MatAssemblyBegin(packM, MAT_FINAL_ASSEMBLY));
2609:   PetscCall(MatAssemblyEnd(packM, MAT_FINAL_ASSEMBLY));
2610:   PetscCall(PetscObjectSetName((PetscObject)packM, "mass"));
2611:   PetscCall(MatViewFromOptions(packM, NULL, "-dm_landau_mass_view"));
2612:   ctx->M = packM;
2613:   if (Amat) *Amat = packM;
2614:   PetscCall(PetscLogEventEnd(ctx->events[14], 0, 0, 0, 0));
2615:   PetscFunctionReturn(PETSC_SUCCESS);
2616: }

2618: /*@
2619:   DMPlexLandauIFunction - `TS` residual calculation, confusingly this computes the Jacobian w/o mass

2621:   Collective

2623:   Input Parameters:
2624: + ts         - The time stepping context
2625: . time_dummy - current time (not used)
2626: . X          - Current state
2627: . X_t        - Time derivative of current state
2628: - actx       - Landau context

2630:   Output Parameter:
2631: . F - The residual

2633:   Level: beginner

2635: .seealso: `DMPlexLandauCreateVelocitySpace()`, `DMPlexLandauIJacobian()`
2636:  @*/
2637: PetscErrorCode DMPlexLandauIFunction(TS ts, PetscReal time_dummy, Vec X, Vec X_t, Vec F, void *actx)
2638: {
2639:   LandauCtx *ctx = (LandauCtx *)actx;
2640:   PetscInt   dim;
2641:   DM         pack;
2642: #if defined(PETSC_HAVE_THREADSAFETY)
2643:   double starttime, endtime;
2644: #endif
2645:   PetscObjectState state;

2647:   PetscFunctionBegin;
2648:   PetscCall(TSGetDM(ts, &pack));
2649:   PetscCall(DMGetApplicationContext(pack, &ctx));
2650:   PetscCheck(ctx, PETSC_COMM_SELF, PETSC_ERR_PLIB, "no context");
2651:   if (ctx->stage) PetscCall(PetscLogStagePush(ctx->stage));
2652:   PetscCall(PetscLogEventBegin(ctx->events[11], 0, 0, 0, 0));
2653:   PetscCall(PetscLogEventBegin(ctx->events[0], 0, 0, 0, 0));
2654: #if defined(PETSC_HAVE_THREADSAFETY)
2655:   starttime = MPI_Wtime();
2656: #endif
2657:   PetscCall(DMGetDimension(pack, &dim));
2658:   PetscCall(PetscObjectStateGet((PetscObject)ctx->J, &state));
2659:   if (state != ctx->norm_state) {
2660:     PetscCall(MatZeroEntries(ctx->J));
2661:     PetscCall(LandauFormJacobian_Internal(X, ctx->J, dim, 0.0, (void *)ctx));
2662:     PetscCall(MatViewFromOptions(ctx->J, NULL, "-dm_landau_jacobian_view"));
2663:     PetscCall(PetscObjectStateGet((PetscObject)ctx->J, &state));
2664:     ctx->norm_state = state;
2665:   } else {
2666:     PetscCall(PetscInfo(ts, "WARNING Skip forming Jacobian, has not changed %" PetscInt64_FMT "\n", state));
2667:   }
2668:   /* mat vec for op */
2669:   PetscCall(MatMult(ctx->J, X, F)); /* C*f */
2670:   /* add time term */
2671:   if (X_t) PetscCall(MatMultAdd(ctx->M, X_t, F, F));
2672: #if defined(PETSC_HAVE_THREADSAFETY)
2673:   if (ctx->stage) {
2674:     endtime = MPI_Wtime();
2675:     ctx->times[LANDAU_OPERATOR] += (endtime - starttime);
2676:     ctx->times[LANDAU_JACOBIAN] += (endtime - starttime);
2677:     ctx->times[LANDAU_MATRIX_TOTAL] += (endtime - starttime);
2678:     ctx->times[LANDAU_JACOBIAN_COUNT] += 1;
2679:   }
2680: #endif
2681:   PetscCall(PetscLogEventEnd(ctx->events[0], 0, 0, 0, 0));
2682:   PetscCall(PetscLogEventEnd(ctx->events[11], 0, 0, 0, 0));
2683:   if (ctx->stage) PetscCall(PetscLogStagePop());
2684:   PetscFunctionReturn(PETSC_SUCCESS);
2685: }

2687: /*@
2688:   DMPlexLandauIJacobian - `TS` Jacobian construction, confusingly this adds mass

2690:   Collective

2692:   Input Parameters:
2693: + ts         - The time stepping context
2694: . time_dummy - current time (not used)
2695: . X          - Current state
2696: . U_tdummy   - Time derivative of current state (not used)
2697: . shift      - shift for du/dt term
2698: - actx       - Landau context

2700:   Output Parameters:
2701: + Amat - Jacobian
2702: - Pmat - same as Amat

2704:   Level: beginner

2706: .seealso: `DMPlexLandauCreateVelocitySpace()`, `DMPlexLandauIFunction()`
2707:  @*/
2708: PetscErrorCode DMPlexLandauIJacobian(TS ts, PetscReal time_dummy, Vec X, Vec U_tdummy, PetscReal shift, Mat Amat, Mat Pmat, void *actx)
2709: {
2710:   LandauCtx *ctx = NULL;
2711:   PetscInt   dim;
2712:   DM         pack;
2713: #if defined(PETSC_HAVE_THREADSAFETY)
2714:   double starttime, endtime;
2715: #endif
2716:   PetscObjectState state;

2718:   PetscFunctionBegin;
2719:   PetscCall(TSGetDM(ts, &pack));
2720:   PetscCall(DMGetApplicationContext(pack, &ctx));
2721:   PetscCheck(ctx, PETSC_COMM_SELF, PETSC_ERR_PLIB, "no context");
2722:   PetscCheck(Amat == Pmat && Amat == ctx->J, ctx->comm, PETSC_ERR_PLIB, "Amat!=Pmat || Amat!=ctx->J");
2723:   PetscCall(DMGetDimension(pack, &dim));
2724:   /* get collision Jacobian into A */
2725:   if (ctx->stage) PetscCall(PetscLogStagePush(ctx->stage));
2726:   PetscCall(PetscLogEventBegin(ctx->events[11], 0, 0, 0, 0));
2727:   PetscCall(PetscLogEventBegin(ctx->events[9], 0, 0, 0, 0));
2728: #if defined(PETSC_HAVE_THREADSAFETY)
2729:   starttime = MPI_Wtime();
2730: #endif
2731:   PetscCheck(shift != 0.0, ctx->comm, PETSC_ERR_PLIB, "zero shift");
2732:   PetscCall(PetscObjectStateGet((PetscObject)ctx->J, &state));
2733:   PetscCheck(state == ctx->norm_state, ctx->comm, PETSC_ERR_PLIB, "wrong state, %" PetscInt64_FMT " %" PetscInt64_FMT, ctx->norm_state, state);
2734:   if (!ctx->use_matrix_mass) {
2735:     PetscCall(LandauFormJacobian_Internal(X, ctx->J, dim, shift, (void *)ctx));
2736:   } else { /* add mass */
2737:     PetscCall(MatAXPY(Pmat, shift, ctx->M, SAME_NONZERO_PATTERN));
2738:   }
2739: #if defined(PETSC_HAVE_THREADSAFETY)
2740:   if (ctx->stage) {
2741:     endtime = MPI_Wtime();
2742:     ctx->times[LANDAU_OPERATOR] += (endtime - starttime);
2743:     ctx->times[LANDAU_MASS] += (endtime - starttime);
2744:     ctx->times[LANDAU_MATRIX_TOTAL] += (endtime - starttime);
2745:   }
2746: #endif
2747:   PetscCall(PetscLogEventEnd(ctx->events[9], 0, 0, 0, 0));
2748:   PetscCall(PetscLogEventEnd(ctx->events[11], 0, 0, 0, 0));
2749:   if (ctx->stage) PetscCall(PetscLogStagePop());
2750:   PetscFunctionReturn(PETSC_SUCCESS);
2751: }