xref: /petsc/src/ts/tests/ex30.c (revision bbdffb7f254ff885b1eb49b202c4ae026f152a3b)
1f3237bfeSMark Adams static char help[] = "Grid based Landau collision operator with PIC interface with OpenMP setup. (one species per grid)\n";
2f3237bfeSMark Adams 
3f3237bfeSMark Adams /*
4f3237bfeSMark Adams    Support 2.5V with axisymmetric coordinates
5f3237bfeSMark Adams      - r,z coordinates
6f3237bfeSMark Adams      - Domain and species data input by Landau operator
7f3237bfeSMark Adams      - "radius" for each grid, normalized with electron thermal velocity
8f3237bfeSMark Adams      - Domain: (0,radius) x (-radius,radius), thus first coordinate x[0] is perpendicular velocity and 2pi*x[0] term is added for axisymmetric
9f3237bfeSMark Adams    Supports full 3V
10f3237bfeSMark Adams 
11f3237bfeSMark Adams  */
12f3237bfeSMark Adams 
1346233b44SBarry Smith #include <petscdmplex.h>
1446233b44SBarry Smith #include <petscds.h>
1546233b44SBarry Smith #include <petscdmswarm.h>
1646233b44SBarry Smith #include <petscksp.h>
17f3237bfeSMark Adams #include <petsc/private/petscimpl.h>
18f3237bfeSMark Adams #if defined(PETSC_HAVE_OPENMP) && defined(PETSC_HAVE_THREADSAFETY)
19f3237bfeSMark Adams   #include <omp.h>
20f3237bfeSMark Adams #endif
21f3237bfeSMark Adams #include <petsclandau.h>
22f3237bfeSMark Adams #include <petscdmcomposite.h>
23f3237bfeSMark Adams 
24f3237bfeSMark Adams typedef struct {
25f3237bfeSMark Adams   Mat MpTrans;
26f3237bfeSMark Adams   Mat Mp;
27f3237bfeSMark Adams   Vec ff;
28f3237bfeSMark Adams   Vec uu;
29f3237bfeSMark Adams } MatShellCtx;
30f3237bfeSMark Adams 
313cdbf90fSMark Adams typedef struct {
323cdbf90fSMark Adams   PetscInt   v_target;
331b431c67SMark Adams   PetscInt   g_target;
341b431c67SMark Adams   PetscInt   global_vertex_id_0;
353cdbf90fSMark Adams   DM        *globSwarmArray;
363cdbf90fSMark Adams   LandauCtx *ctx;
373cdbf90fSMark Adams   DM        *grid_dm;
383cdbf90fSMark Adams   Mat       *g_Mass;
393cdbf90fSMark Adams   Mat       *globMpArray;
403cdbf90fSMark Adams   Vec       *globXArray;
413cdbf90fSMark Adams   PetscBool  print;
421b431c67SMark Adams   PetscBool  print_entropy;
433cdbf90fSMark Adams } PrintCtx;
443cdbf90fSMark Adams 
45d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultMtM_SeqAIJ(Mat MtM, Vec xx, Vec yy)
46d71ae5a4SJacob Faibussowitsch {
47f3237bfeSMark Adams   MatShellCtx *matshellctx;
48f3237bfeSMark Adams 
49f3237bfeSMark Adams   PetscFunctionBeginUser;
50f3237bfeSMark Adams   PetscCall(MatShellGetContext(MtM, &matshellctx));
51f3237bfeSMark Adams   PetscCheck(matshellctx, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "No context");
52f3237bfeSMark Adams   PetscCall(MatMult(matshellctx->Mp, xx, matshellctx->ff));
53f3237bfeSMark Adams   PetscCall(MatMult(matshellctx->MpTrans, matshellctx->ff, yy));
543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
55f3237bfeSMark Adams }
56f3237bfeSMark Adams 
57d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultAddMtM_SeqAIJ(Mat MtM, Vec xx, Vec yy, Vec zz)
58d71ae5a4SJacob Faibussowitsch {
59f3237bfeSMark Adams   MatShellCtx *matshellctx;
60f3237bfeSMark Adams 
61f3237bfeSMark Adams   PetscFunctionBeginUser;
62f3237bfeSMark Adams   PetscCall(MatShellGetContext(MtM, &matshellctx));
63f3237bfeSMark Adams   PetscCheck(matshellctx, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "No context");
64f3237bfeSMark Adams   PetscCall(MatMult(matshellctx->Mp, xx, matshellctx->ff));
65f3237bfeSMark Adams   PetscCall(MatMultAdd(matshellctx->MpTrans, matshellctx->ff, yy, zz));
663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
67f3237bfeSMark Adams }
68f3237bfeSMark Adams 
69d71ae5a4SJacob Faibussowitsch PetscErrorCode createSwarm(const DM dm, PetscInt dim, DM *sw)
70d71ae5a4SJacob Faibussowitsch {
71f3237bfeSMark Adams   PetscInt Nc = 1;
72f3237bfeSMark Adams 
73f3237bfeSMark Adams   PetscFunctionBeginUser;
74f3237bfeSMark Adams   PetscCall(DMCreate(PETSC_COMM_SELF, sw));
75f3237bfeSMark Adams   PetscCall(DMSetType(*sw, DMSWARM));
76f3237bfeSMark Adams   PetscCall(DMSetDimension(*sw, dim));
77f3237bfeSMark Adams   PetscCall(DMSwarmSetType(*sw, DMSWARM_PIC));
78f3237bfeSMark Adams   PetscCall(DMSwarmSetCellDM(*sw, dm));
793cdbf90fSMark Adams   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "w_q", Nc, PETSC_REAL));
80f3237bfeSMark Adams   PetscCall(DMSwarmFinalizeFieldRegister(*sw));
81f3237bfeSMark Adams   PetscCall(DMSetFromOptions(*sw));
823cdbf90fSMark Adams   PetscCall(PetscObjectSetName((PetscObject)*sw, "Particle Grid"));
833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
84f3237bfeSMark Adams }
85f3237bfeSMark Adams 
863cdbf90fSMark Adams static PetscErrorCode makeSwarm(DM sw, const PetscInt dim, const PetscInt Np, const PetscReal xx[], const PetscReal yy[], const PetscReal zz[])
873cdbf90fSMark Adams {
883cdbf90fSMark Adams   PetscReal    *coords;
893cdbf90fSMark Adams   PetscDataType dtype;
903cdbf90fSMark Adams   PetscInt      bs, p, zero = 0;
913cdbf90fSMark Adams 
924d86920dSPierre Jolivet   PetscFunctionBeginUser;
933cdbf90fSMark Adams   PetscCall(DMSwarmSetLocalSizes(sw, Np, zero));
943cdbf90fSMark Adams   PetscCall(DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
953cdbf90fSMark Adams   for (p = 0; p < Np; p++) {
963cdbf90fSMark Adams     coords[p * dim + 0] = xx[p];
973cdbf90fSMark Adams     coords[p * dim + 1] = yy[p];
983cdbf90fSMark Adams     if (dim == 3) coords[p * dim + 2] = zz[p];
993cdbf90fSMark Adams   }
1003cdbf90fSMark Adams   PetscCall(DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
1013cdbf90fSMark Adams   PetscFunctionReturn(PETSC_SUCCESS);
1023cdbf90fSMark Adams }
1033cdbf90fSMark Adams 
1043cdbf90fSMark Adams static PetscErrorCode createMp(const DM dm, DM sw, Mat *Mp_out)
1053cdbf90fSMark Adams {
1063cdbf90fSMark Adams   PetscBool removePoints = PETSC_TRUE;
1073cdbf90fSMark Adams   Mat       M_p;
1084d86920dSPierre Jolivet 
1093cdbf90fSMark Adams   PetscFunctionBeginUser;
1103cdbf90fSMark Adams   // migrate after coords are set
1113cdbf90fSMark Adams   PetscCall(DMSwarmMigrate(sw, removePoints));
1123cdbf90fSMark Adams   PetscCall(PetscObjectSetName((PetscObject)sw, "Particle Grid"));
1133cdbf90fSMark Adams   /* This gives M f = \int_\Omega \phi f, which looks like a rhs for a PDE */
1143cdbf90fSMark Adams   PetscCall(DMCreateMassMatrix(sw, dm, &M_p));
1153cdbf90fSMark Adams   PetscCall(DMViewFromOptions(sw, NULL, "-ex30_sw_view"));
1163cdbf90fSMark Adams   // output
1173cdbf90fSMark Adams   *Mp_out = M_p;
1183cdbf90fSMark Adams   PetscFunctionReturn(PETSC_SUCCESS);
1193cdbf90fSMark Adams }
1203cdbf90fSMark Adams 
1211b431c67SMark Adams static PetscErrorCode particlesToGrid(const DM dm, DM sw, const PetscInt a_tid, const PetscInt dim, const PetscReal a_wp[], Vec rho, Mat M_p)
1223cdbf90fSMark Adams {
1233cdbf90fSMark Adams   PetscReal    *wq;
1243cdbf90fSMark Adams   PetscDataType dtype;
1253cdbf90fSMark Adams   Vec           ff;
1261b431c67SMark Adams   PetscInt      bs, p, Np;
1273cdbf90fSMark Adams 
1283cdbf90fSMark Adams   PetscFunctionBeginUser;
1293cdbf90fSMark Adams   PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void **)&wq));
1301b431c67SMark Adams   PetscCall(DMSwarmGetLocalSize(sw, &Np));
1313cdbf90fSMark Adams   for (p = 0; p < Np; p++) wq[p] = a_wp[p];
1323cdbf90fSMark Adams   PetscCall(DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void **)&wq));
1333cdbf90fSMark Adams   PetscCall(PetscObjectSetName((PetscObject)rho, "rho"));
1343cdbf90fSMark Adams   PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &ff));
1353cdbf90fSMark Adams   PetscCall(PetscObjectSetName((PetscObject)ff, "weights"));
1363cdbf90fSMark Adams   PetscCall(MatMultTranspose(M_p, ff, rho));
1373cdbf90fSMark Adams   PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &ff));
1383cdbf90fSMark Adams   PetscFunctionReturn(PETSC_SUCCESS);
1393cdbf90fSMark Adams }
1403cdbf90fSMark Adams 
1413cdbf90fSMark Adams //
1423cdbf90fSMark Adams // add grid to arg 'sw.w_q'
1433cdbf90fSMark Adams //
1443cdbf90fSMark Adams PetscErrorCode gridToParticles(const DM dm, DM sw, const Vec rhs, Vec work, Mat M_p, Mat Mass)
145d71ae5a4SJacob Faibussowitsch {
146f3237bfeSMark Adams   PetscBool    is_lsqr;
147f3237bfeSMark Adams   KSP          ksp;
148f3237bfeSMark Adams   Mat          PM_p = NULL, MtM, D;
149f3237bfeSMark Adams   Vec          ff;
150f3237bfeSMark Adams   PetscInt     N, M, nzl;
151f3237bfeSMark Adams   MatShellCtx *matshellctx;
1521b431c67SMark Adams   PC           pc;
153f3237bfeSMark Adams 
154f3237bfeSMark Adams   PetscFunctionBeginUser;
1551b431c67SMark Adams   // (Mp Mp)^-1 M
156f3237bfeSMark Adams   PetscCall(MatMult(Mass, rhs, work));
157f3237bfeSMark Adams   // pseudo-inverse
158f3237bfeSMark Adams   PetscCall(KSPCreate(PETSC_COMM_SELF, &ksp));
1591b431c67SMark Adams   PetscCall(KSPSetType(ksp, KSPCG));
1601b431c67SMark Adams   PetscCall(KSPGetPC(ksp, &pc));
1611b431c67SMark Adams   PetscCall(PCSetType(pc, PCJACOBI));
162f3237bfeSMark Adams   PetscCall(KSPSetOptionsPrefix(ksp, "ftop_"));
163f3237bfeSMark Adams   PetscCall(KSPSetFromOptions(ksp));
164f3237bfeSMark Adams   PetscCall(PetscObjectTypeCompare((PetscObject)ksp, KSPLSQR, &is_lsqr));
165f3237bfeSMark Adams   if (!is_lsqr) {
166f3237bfeSMark Adams     PetscCall(MatGetLocalSize(M_p, &M, &N));
167f3237bfeSMark Adams     if (N > M) {
1681b431c67SMark Adams       PetscCall(PetscInfo(ksp, " M (%" PetscInt_FMT ") < M (%" PetscInt_FMT ") more vertices than particles: revert to lsqr\n", M, N));
169f3237bfeSMark Adams       is_lsqr = PETSC_TRUE;
170f3237bfeSMark Adams       PetscCall(KSPSetType(ksp, KSPLSQR));
1711b431c67SMark Adams       PetscCall(PCSetType(pc, PCNONE)); // should not happen, but could solve stable (Mp Mp^T), move projection Mp before solve
172f3237bfeSMark Adams     } else {
173f3237bfeSMark Adams       PetscCall(PetscNew(&matshellctx));
174f3237bfeSMark Adams       PetscCall(MatCreateShell(PetscObjectComm((PetscObject)dm), N, N, PETSC_DECIDE, PETSC_DECIDE, matshellctx, &MtM));
175f3237bfeSMark Adams       PetscCall(MatTranspose(M_p, MAT_INITIAL_MATRIX, &matshellctx->MpTrans));
176f3237bfeSMark Adams       matshellctx->Mp = M_p;
177f3237bfeSMark Adams       PetscCall(MatShellSetOperation(MtM, MATOP_MULT, (void (*)(void))MatMultMtM_SeqAIJ));
178f3237bfeSMark Adams       PetscCall(MatShellSetOperation(MtM, MATOP_MULT_ADD, (void (*)(void))MatMultAddMtM_SeqAIJ));
179f3237bfeSMark Adams       PetscCall(MatCreateVecs(M_p, &matshellctx->uu, &matshellctx->ff));
180f3237bfeSMark Adams       PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, N, N, 1, NULL, &D));
181f3237bfeSMark Adams       PetscCall(MatViewFromOptions(matshellctx->MpTrans, NULL, "-ftop2_Mp_mat_view"));
182f3237bfeSMark Adams       for (int i = 0; i < N; i++) {
183f3237bfeSMark Adams         const PetscScalar *vals;
184f3237bfeSMark Adams         const PetscInt    *cols;
185f3237bfeSMark Adams         PetscScalar        dot = 0;
186f3237bfeSMark Adams         PetscCall(MatGetRow(matshellctx->MpTrans, i, &nzl, &cols, &vals));
187f3237bfeSMark Adams         for (int ii = 0; ii < nzl; ii++) dot += PetscSqr(vals[ii]);
1881b431c67SMark Adams         if (dot == 0.0) dot = 1; // empty rows
189f3237bfeSMark Adams         PetscCall(MatSetValue(D, i, i, dot, INSERT_VALUES));
190f3237bfeSMark Adams       }
191f3237bfeSMark Adams       PetscCall(MatAssemblyBegin(D, MAT_FINAL_ASSEMBLY));
192f3237bfeSMark Adams       PetscCall(MatAssemblyEnd(D, MAT_FINAL_ASSEMBLY));
193f3237bfeSMark Adams       PetscCall(PetscInfo(M_p, "createMtMKSP Have %" PetscInt_FMT " eqs, nzl = %" PetscInt_FMT "\n", N, nzl));
194f3237bfeSMark Adams       PetscCall(KSPSetOperators(ksp, MtM, D));
195f3237bfeSMark Adams       PetscCall(MatViewFromOptions(D, NULL, "-ftop2_D_mat_view"));
196f3237bfeSMark Adams       PetscCall(MatViewFromOptions(M_p, NULL, "-ftop2_Mp_mat_view"));
197f3237bfeSMark Adams       PetscCall(MatViewFromOptions(matshellctx->MpTrans, NULL, "-ftop2_MpTranspose_mat_view"));
1986b664d00Smarkadams4       PetscCall(MatViewFromOptions(MtM, NULL, "-ftop2_MtM_mat_view"));
199f3237bfeSMark Adams     }
200f3237bfeSMark Adams   }
201f3237bfeSMark Adams   if (is_lsqr) {
202f3237bfeSMark Adams     PC        pc;
203f3237bfeSMark Adams     PetscBool is_bjac;
204f3237bfeSMark Adams     PetscCall(KSPGetPC(ksp, &pc));
205f3237bfeSMark Adams     PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCBJACOBI, &is_bjac));
206f3237bfeSMark Adams     if (is_bjac) {
207f3237bfeSMark Adams       PetscCall(DMSwarmCreateMassMatrixSquare(sw, dm, &PM_p));
208f3237bfeSMark Adams       PetscCall(KSPSetOperators(ksp, M_p, PM_p));
209f3237bfeSMark Adams     } else {
210f3237bfeSMark Adams       PetscCall(KSPSetOperators(ksp, M_p, M_p));
211f3237bfeSMark Adams     }
212f3237bfeSMark Adams   }
213f3237bfeSMark Adams   PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &ff)); // this grabs access
214f3237bfeSMark Adams   if (!is_lsqr) {
2153cdbf90fSMark Adams     PetscCall(KSPSolve(ksp, work, matshellctx->uu));
216f3237bfeSMark Adams     PetscCall(MatMult(M_p, matshellctx->uu, ff));
217f3237bfeSMark Adams     PetscCall(MatDestroy(&matshellctx->MpTrans));
218f3237bfeSMark Adams     PetscCall(VecDestroy(&matshellctx->ff));
219f3237bfeSMark Adams     PetscCall(VecDestroy(&matshellctx->uu));
220f3237bfeSMark Adams     PetscCall(MatDestroy(&D));
221f3237bfeSMark Adams     PetscCall(MatDestroy(&MtM));
222f3237bfeSMark Adams     PetscCall(PetscFree(matshellctx));
223f3237bfeSMark Adams   } else {
2243cdbf90fSMark Adams     PetscCall(KSPSolveTranspose(ksp, work, ff));
225f3237bfeSMark Adams   }
226f3237bfeSMark Adams   PetscCall(KSPDestroy(&ksp));
227f3237bfeSMark Adams   PetscCall(MatDestroy(&PM_p));
228f3237bfeSMark Adams   PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &ff));
2293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
230f3237bfeSMark Adams }
231f3237bfeSMark Adams 
2323cdbf90fSMark Adams #define EX30_MAX_NUM_THRDS 12
2333cdbf90fSMark Adams #define EX30_MAX_BATCH_SZ  1024
2343cdbf90fSMark Adams //
2353cdbf90fSMark Adams // add grid to arg 'globSwarmArray[].w_q'
2363cdbf90fSMark Adams //
2373cdbf90fSMark Adams PetscErrorCode gridToParticles_private(DM grid_dm[], DM globSwarmArray[], const PetscInt dim, const PetscInt v_target, const PetscInt numthreads, const PetscInt num_vertices, const PetscInt global_vertex_id, Mat globMpArray[], Mat g_Mass[], Vec t_fhat[][EX30_MAX_NUM_THRDS], PetscReal moments[], Vec globXArray[], LandauCtx *ctx)
238d71ae5a4SJacob Faibussowitsch {
2393cdbf90fSMark Adams   PetscErrorCode ierr = (PetscErrorCode)0; // used for inside thread loops
240f3237bfeSMark Adams 
241f3237bfeSMark Adams   PetscFunctionBeginUser;
2423cdbf90fSMark Adams   // map back to particles
2433cdbf90fSMark Adams   for (PetscInt v_id_0 = 0; v_id_0 < ctx->batch_sz; v_id_0 += numthreads) {
2443cdbf90fSMark Adams     PetscCall(PetscInfo(grid_dm[0], "g2p: global batch %" PetscInt_FMT " of %" PetscInt_FMT ", Landau batch %" PetscInt_FMT " of %" PetscInt_FMT ": map back to particles\n", global_vertex_id + 1, num_vertices, v_id_0 + 1, ctx->batch_sz));
2453cdbf90fSMark Adams     //PetscPragmaOMP(parallel for)
2463cdbf90fSMark Adams     for (int tid = 0; tid < numthreads; tid++) {
2473cdbf90fSMark Adams       const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id + v_id;
2483cdbf90fSMark Adams       if (glb_v_id < num_vertices) {
2493cdbf90fSMark Adams         for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids
2503cdbf90fSMark Adams           PetscErrorCode ierr_t;
2513cdbf90fSMark Adams           ierr_t = PetscInfo(grid_dm[0], "gridToParticles: global batch %" PetscInt_FMT ", local batch b=%" PetscInt_FMT ", grid g=%" PetscInt_FMT ", index(b,g) %" PetscInt_FMT "\n", global_vertex_id, v_id, grid, LAND_PACK_IDX(v_id, grid));
2523cdbf90fSMark Adams           ierr_t = gridToParticles(grid_dm[grid], globSwarmArray[LAND_PACK_IDX(v_id, grid)], globXArray[LAND_PACK_IDX(v_id, grid)], t_fhat[grid][tid], globMpArray[LAND_PACK_IDX(v_id, grid)], g_Mass[grid]);
2533cdbf90fSMark Adams           if (ierr_t) ierr = ierr_t;
254f3237bfeSMark Adams         }
2553cdbf90fSMark Adams       }
2563cdbf90fSMark Adams     }
2573cdbf90fSMark Adams     PetscCheck(!ierr, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Error in OMP loop. ierr = %d", (int)ierr);
2583cdbf90fSMark Adams     /* Get moments */
2593cdbf90fSMark Adams     PetscCall(PetscInfo(grid_dm[0], "Cleanup batches %" PetscInt_FMT " to %" PetscInt_FMT "\n", v_id_0, v_id_0 + numthreads));
2603cdbf90fSMark Adams     for (int tid = 0; tid < numthreads; tid++) {
2613cdbf90fSMark Adams       const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id + v_id;
2623cdbf90fSMark Adams       if (glb_v_id == v_target) {
2633cdbf90fSMark Adams         for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
2643cdbf90fSMark Adams           PetscDataType dtype;
2653cdbf90fSMark Adams           PetscReal    *wp, *coords;
2663cdbf90fSMark Adams           DM            sw = globSwarmArray[LAND_PACK_IDX(v_id, grid)];
2673cdbf90fSMark Adams           PetscInt      npoints, bs = 1;
2683cdbf90fSMark Adams           PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void **)&wp)); // take data out here
2693cdbf90fSMark Adams           PetscCall(DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
2703cdbf90fSMark Adams           PetscCall(DMSwarmGetLocalSize(sw, &npoints));
2713cdbf90fSMark Adams           for (int p = 0; p < npoints; p++) {
2723cdbf90fSMark Adams             PetscReal v2 = 0, fact = (dim == 2) ? 2.0 * PETSC_PI * coords[p * dim + 0] : 1, w = fact * wp[p] * ctx->n_0 * ctx->masses[ctx->species_offset[grid]];
2733cdbf90fSMark Adams             for (int i = 0; i < dim; ++i) v2 += PetscSqr(coords[p * dim + i]);
2743cdbf90fSMark Adams             moments[0] += w;
2753cdbf90fSMark Adams             moments[1] += w * ctx->v_0 * coords[p * dim + 1]; // z-momentum
2761b431c67SMark Adams             moments[2] += w * 0.5 * ctx->v_0 * ctx->v_0 * v2;
2773cdbf90fSMark Adams           }
2783cdbf90fSMark Adams           PetscCall(DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void **)&wp));
279f3237bfeSMark Adams           PetscCall(DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
2803cdbf90fSMark Adams         }
2813cdbf90fSMark Adams         const PetscReal N_inv = 1 / moments[0];
2821b431c67SMark Adams         PetscCall(PetscInfo(grid_dm[0], "gridToParticles_private [%" PetscInt_FMT "], n = %g\n", v_id, (double)moments[0]));
2833cdbf90fSMark Adams         for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
2843cdbf90fSMark Adams           PetscDataType dtype;
2853cdbf90fSMark Adams           PetscReal    *wp, *coords;
2863cdbf90fSMark Adams           DM            sw = globSwarmArray[LAND_PACK_IDX(v_id, grid)];
2873cdbf90fSMark Adams           PetscInt      npoints, bs = 1;
2883cdbf90fSMark Adams           PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void **)&wp)); // take data out here
2893cdbf90fSMark Adams           PetscCall(DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
2903cdbf90fSMark Adams           PetscCall(DMSwarmGetLocalSize(sw, &npoints));
2913cdbf90fSMark Adams           for (int p = 0; p < npoints; p++) {
2923cdbf90fSMark Adams             const PetscReal fact = dim == 2 ? 2.0 * PETSC_PI * coords[p * dim + 0] : 1, w = fact * wp[p] * ctx->n_0 * ctx->masses[ctx->species_offset[grid]], ww = w * N_inv;
2931b431c67SMark Adams             if (w > PETSC_REAL_MIN) {
2943cdbf90fSMark Adams               moments[3] -= ww * PetscLogReal(ww);
2953cdbf90fSMark Adams               PetscCheck(ww < 1 - PETSC_MACHINE_EPSILON, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "ww (%g) > 1", (double)ww);
2961b431c67SMark Adams             } else moments[4] -= w;
2973cdbf90fSMark Adams           }
2983cdbf90fSMark Adams           PetscCall(DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void **)&wp));
2993cdbf90fSMark Adams           PetscCall(DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
3003cdbf90fSMark Adams         }
3013cdbf90fSMark Adams       }
3023cdbf90fSMark Adams     } // thread batch
3033cdbf90fSMark Adams   } // batch
3043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
305f3237bfeSMark Adams }
3063cdbf90fSMark Adams 
3073cdbf90fSMark Adams static void maxwellian(PetscInt dim, const PetscReal x[], PetscReal kt_m, PetscReal n, PetscReal shift, PetscScalar *u)
308d71ae5a4SJacob Faibussowitsch {
309f3237bfeSMark Adams   PetscInt  i;
310f3237bfeSMark Adams   PetscReal v2 = 0, theta = 2.0 * kt_m; /* theta = 2kT/mc^2 */
311f3237bfeSMark Adams 
3123cdbf90fSMark Adams   if (shift != 0.) {
3133cdbf90fSMark Adams     v2 = 0;
3143cdbf90fSMark Adams     for (i = 0; i < dim - 1; ++i) v2 += x[i] * x[i];
3153cdbf90fSMark Adams     v2 += (x[dim - 1] - shift) * (x[dim - 1] - shift);
3163cdbf90fSMark Adams     /* evaluate the shifted Maxwellian */
3173cdbf90fSMark Adams     u[0] += n * PetscPowReal(PETSC_PI * theta, -1.5) * (PetscExpReal(-v2 / theta));
3181b431c67SMark Adams   } else {
3191b431c67SMark Adams     /* compute the exponents, v^2 */
3201b431c67SMark Adams     for (i = 0; i < dim; ++i) v2 += x[i] * x[i];
3211b431c67SMark Adams     /* evaluate the Maxwellian */
3221b431c67SMark Adams     u[0] += n * PetscPowReal(PETSC_PI * theta, -1.5) * (PetscExpReal(-v2 / theta));
3233cdbf90fSMark Adams   }
324f3237bfeSMark Adams }
325f3237bfeSMark Adams 
3263cdbf90fSMark Adams static PetscErrorCode PostStep(TS ts)
3273cdbf90fSMark Adams {
3281b431c67SMark Adams   PetscInt   n, dim, nDMs, v_id;
3293cdbf90fSMark Adams   PetscReal  t;
3303cdbf90fSMark Adams   LandauCtx *ctx;
3313cdbf90fSMark Adams   Vec        X;
3323cdbf90fSMark Adams   PrintCtx  *printCtx;
3333cdbf90fSMark Adams   DM         pack;
3341b431c67SMark Adams   PetscReal  moments[5], e_grid[LANDAU_MAX_GRIDS];
3353cdbf90fSMark Adams 
3363cdbf90fSMark Adams   PetscFunctionBeginUser;
3373cdbf90fSMark Adams   PetscCall(TSGetApplicationContext(ts, &printCtx));
3381b431c67SMark Adams   if (!printCtx->print && !printCtx->print_entropy) PetscFunctionReturn(PETSC_SUCCESS);
3393cdbf90fSMark Adams   ctx = printCtx->ctx;
3401b431c67SMark Adams   if (printCtx->v_target < printCtx->global_vertex_id_0 || printCtx->v_target >= printCtx->global_vertex_id_0 + ctx->batch_sz) PetscFunctionReturn(PETSC_SUCCESS);
3411b431c67SMark Adams   for (int i = 0; i < 5; i++) moments[i] = 0;
3421b431c67SMark Adams   for (int i = 0; i < LANDAU_MAX_GRIDS; i++) e_grid[i] = 0;
3431b431c67SMark Adams   v_id = printCtx->v_target % ctx->batch_sz;
3443cdbf90fSMark Adams   PetscCall(TSGetDM(ts, &pack));
3453cdbf90fSMark Adams   PetscCall(DMGetDimension(pack, &dim));
3463cdbf90fSMark Adams   PetscCall(DMCompositeGetNumberDM(pack, &nDMs)); // number of vertices * number of grids
3473cdbf90fSMark Adams   PetscCall(TSGetSolution(ts, &X));
3483cdbf90fSMark Adams   PetscCall(TSGetStepNumber(ts, &n));
3493cdbf90fSMark Adams   PetscCall(TSGetTime(ts, &t));
3503cdbf90fSMark Adams   PetscCall(DMCompositeGetAccessArray(pack, X, nDMs, NULL, printCtx->globXArray));
3511b431c67SMark Adams   if (printCtx->print_entropy && printCtx->v_target >= 0) {
3523cdbf90fSMark Adams     for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
3533cdbf90fSMark Adams       PetscDataType dtype;
3543cdbf90fSMark Adams       PetscReal    *wp, *coords;
3553cdbf90fSMark Adams       DM            sw = printCtx->globSwarmArray[LAND_PACK_IDX(v_id, grid)];
3563cdbf90fSMark Adams       Vec           work, subX = printCtx->globXArray[LAND_PACK_IDX(v_id, grid)];
3571b431c67SMark Adams       PetscInt      bs, NN;
3583cdbf90fSMark Adams       // C-G moments
3593cdbf90fSMark Adams       PetscCall(VecDuplicate(subX, &work));
3603cdbf90fSMark Adams       PetscCall(gridToParticles(printCtx->grid_dm[grid], sw, subX, work, printCtx->globMpArray[LAND_PACK_IDX(v_id, grid)], printCtx->g_Mass[grid]));
3613cdbf90fSMark Adams       PetscCall(VecDestroy(&work));
3623cdbf90fSMark Adams       // moments
3633cdbf90fSMark Adams       PetscCall(DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
3641b431c67SMark Adams       PetscCall(DMSwarmGetLocalSize(sw, &NN));
3651b431c67SMark Adams       PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void **)&wp));
3663cdbf90fSMark Adams       for (int pp = 0; pp < NN; pp++) {
3671b431c67SMark Adams         PetscReal v2 = 0, fact = (dim == 2) ? 2.0 * PETSC_PI * coords[pp * dim + 0] : 1, w = fact * wp[pp] * ctx->n_0 * ctx->masses[ctx->species_offset[grid]];
3683cdbf90fSMark Adams         for (int i = 0; i < dim; ++i) v2 += PetscSqr(coords[pp * dim + i]);
3693cdbf90fSMark Adams         moments[0] += w;
3703cdbf90fSMark Adams         moments[1] += w * ctx->v_0 * coords[pp * dim + 1]; // z-momentum
3711b431c67SMark Adams         moments[2] += w * 0.5 * ctx->v_0 * ctx->v_0 * v2;
3721b431c67SMark Adams         e_grid[grid] += w * 0.5 * ctx->v_0 * ctx->v_0 * v2;
3733cdbf90fSMark Adams       }
3743cdbf90fSMark Adams       PetscCall(DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
3753cdbf90fSMark Adams       PetscCall(DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void **)&wp));
3763cdbf90fSMark Adams     }
3771b431c67SMark Adams     // entropy
3781b431c67SMark Adams     const PetscReal N_inv = 1 / moments[0];
3791b431c67SMark Adams     for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
3801b431c67SMark Adams       PetscDataType dtype;
3811b431c67SMark Adams       PetscReal    *wp, *coords;
3821b431c67SMark Adams       DM            sw = printCtx->globSwarmArray[LAND_PACK_IDX(v_id, grid)];
3831b431c67SMark Adams       PetscInt      bs, NN;
3841b431c67SMark Adams       PetscCall(DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
3851b431c67SMark Adams       PetscCall(DMSwarmGetLocalSize(sw, &NN));
3861b431c67SMark Adams       PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void **)&wp));
3871b431c67SMark Adams       for (int pp = 0; pp < NN; pp++) {
3881b431c67SMark Adams         PetscReal fact = (dim == 2) ? 2.0 * PETSC_PI * coords[pp * dim + 0] : 1, w = fact * wp[pp] * ctx->n_0 * ctx->masses[ctx->species_offset[grid]], ww = w * N_inv;
3891b431c67SMark Adams         if (w > PETSC_REAL_MIN) {
3901b431c67SMark Adams           moments[3] -= ww * PetscLogReal(ww);
3911b431c67SMark Adams           PetscCheck(ww < 1 - PETSC_MACHINE_EPSILON, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "ww (%g) > 1", (double)ww);
3921b431c67SMark Adams         } else moments[4] -= w;
3931b431c67SMark Adams       }
3941b431c67SMark Adams       PetscCall(DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
3951b431c67SMark Adams       PetscCall(DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void **)&wp));
3961b431c67SMark Adams     }
3971b431c67SMark Adams     PetscCall(PetscInfo(X, "%4d) time %e, Landau particle moments: 0: %18.12e 1: %19.12e 2: %18.12e entropy: %e loss %e. energy = %e + %e + %e\n", (int)n, (double)t, (double)moments[0], (double)moments[1], (double)moments[2], (double)moments[3], (double)moments[4], (double)e_grid[0], (double)e_grid[1], (double)e_grid[2]));
3981b431c67SMark Adams   }
3991b431c67SMark Adams   if (printCtx->print && printCtx->g_target >= 0) {
4001b431c67SMark Adams     PetscInt         grid   = printCtx->g_target, id;
4011b431c67SMark Adams     static PetscReal last_t = -100000, period = .5;
4021b431c67SMark Adams     if (last_t == -100000) last_t = -period + t;
4031b431c67SMark Adams     if (t >= last_t + period) {
4041b431c67SMark Adams       last_t = t;
4051b431c67SMark Adams       PetscCall(DMGetOutputSequenceNumber(ctx->plex[grid], &id, NULL));
4061b431c67SMark Adams       PetscCall(DMSetOutputSequenceNumber(ctx->plex[grid], id + 1, t));
4071b431c67SMark Adams       PetscCall(VecViewFromOptions(printCtx->globXArray[LAND_PACK_IDX(v_id % ctx->batch_sz, grid)], NULL, "-ex30_vec_view"));
4081b431c67SMark Adams       if (ctx->num_grids > grid + 1) {
4091b431c67SMark Adams         PetscCall(DMSetOutputSequenceNumber(ctx->plex[grid + 1], id + 1, t));
4101b431c67SMark Adams         PetscCall(VecViewFromOptions(printCtx->globXArray[LAND_PACK_IDX(v_id % ctx->batch_sz, grid + 1)], NULL, "-ex30_vec_view2"));
4111b431c67SMark Adams       }
4121b431c67SMark Adams       PetscCall(PetscInfo(X, "%4d) time %e View\n", (int)n, (double)t));
4131b431c67SMark Adams     }
4141b431c67SMark Adams   }
4153cdbf90fSMark Adams   PetscCall(DMCompositeRestoreAccessArray(pack, X, nDMs, NULL, printCtx->globXArray));
4163cdbf90fSMark Adams   PetscFunctionReturn(PETSC_SUCCESS);
4173cdbf90fSMark Adams }
4183cdbf90fSMark Adams 
4191b431c67SMark Adams PetscErrorCode go(TS ts, Vec X, const PetscInt num_vertices, const PetscInt a_Np, const PetscInt dim, const PetscInt v_target, const PetscInt g_target, PetscReal shift, PetscBool use_uniform_particle_grid)
420d71ae5a4SJacob Faibussowitsch {
421f3237bfeSMark Adams   DM             pack, *globSwarmArray, grid_dm[LANDAU_MAX_GRIDS];
422f3237bfeSMark Adams   Mat           *globMpArray, g_Mass[LANDAU_MAX_GRIDS];
4233cdbf90fSMark Adams   KSP            t_ksp[LANDAU_MAX_GRIDS][EX30_MAX_NUM_THRDS];
4243cdbf90fSMark Adams   Vec            t_fhat[LANDAU_MAX_GRIDS][EX30_MAX_NUM_THRDS];
4251b431c67SMark Adams   PetscInt       nDMs;
4263cdbf90fSMark Adams   PetscErrorCode ierr = (PetscErrorCode)0; // used for inside thread loops
427f3237bfeSMark Adams #if defined(PETSC_HAVE_OPENMP) && defined(PETSC_HAVE_THREADSAFETY)
428f3237bfeSMark Adams   PetscInt numthreads = PetscNumOMPThreads;
429f3237bfeSMark Adams #else
430f3237bfeSMark Adams   PetscInt numthreads = 1;
431f3237bfeSMark Adams #endif
432f3237bfeSMark Adams   LandauCtx *ctx;
433f3237bfeSMark Adams   Vec       *globXArray;
4341b431c67SMark Adams   PetscReal  moments_0[5], moments_1a[5], moments_1b[5], dt_init;
4353cdbf90fSMark Adams   PrintCtx  *printCtx;
436f3237bfeSMark Adams 
437f3237bfeSMark Adams   PetscFunctionBeginUser;
4383cdbf90fSMark Adams   PetscCheck(numthreads <= EX30_MAX_NUM_THRDS, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Too many threads %" PetscInt_FMT " > %d", numthreads, EX30_MAX_NUM_THRDS);
4393cdbf90fSMark Adams   PetscCheck(numthreads > 0, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Number threads %" PetscInt_FMT " > %d", numthreads, EX30_MAX_NUM_THRDS);
440f3237bfeSMark Adams   PetscCall(TSGetDM(ts, &pack));
441f3237bfeSMark Adams   PetscCall(DMGetApplicationContext(pack, &ctx));
442f3237bfeSMark Adams   PetscCheck(ctx->batch_sz % numthreads == 0, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "batch size (-dm_landau_batch_size) %" PetscInt_FMT "  mod #threads %" PetscInt_FMT " must equal zero", ctx->batch_sz, numthreads);
4433cdbf90fSMark Adams   PetscCall(DMCompositeGetNumberDM(pack, &nDMs)); // number of vertices * number of grids
4441b431c67SMark Adams   PetscCall(PetscInfo(pack, "Have %" PetscInt_FMT " total grids, with %" PetscInt_FMT " Landau local batched and %" PetscInt_FMT " global items (vertices) %d DMs\n", ctx->num_grids, ctx->batch_sz, num_vertices, (int)nDMs));
445f3237bfeSMark Adams   PetscCall(PetscMalloc(sizeof(*globXArray) * nDMs, &globXArray));
446f3237bfeSMark Adams   PetscCall(PetscMalloc(sizeof(*globMpArray) * nDMs, &globMpArray));
447f3237bfeSMark Adams   PetscCall(PetscMalloc(sizeof(*globSwarmArray) * nDMs, &globSwarmArray));
4483cdbf90fSMark Adams   // print ctx
4493cdbf90fSMark Adams   PetscCall(PetscNew(&printCtx));
4503cdbf90fSMark Adams   PetscCall(TSSetApplicationContext(ts, printCtx));
4513cdbf90fSMark Adams   printCtx->v_target       = v_target;
4521b431c67SMark Adams   printCtx->g_target       = g_target;
4533cdbf90fSMark Adams   printCtx->ctx            = ctx;
4543cdbf90fSMark Adams   printCtx->globSwarmArray = globSwarmArray;
4553cdbf90fSMark Adams   printCtx->grid_dm        = grid_dm;
4563cdbf90fSMark Adams   printCtx->globMpArray    = globMpArray;
4573cdbf90fSMark Adams   printCtx->g_Mass         = g_Mass;
4583cdbf90fSMark Adams   printCtx->globXArray     = globXArray;
4591b431c67SMark Adams   printCtx->print_entropy  = PETSC_FALSE;
46015229ffcSPierre Jolivet   PetscOptionsBegin(PETSC_COMM_SELF, "", "Print Options", "DMPLEX");
4611b431c67SMark Adams   PetscCall(PetscOptionsBool("-print_entropy", "Print entropy and moments at each time step", "ex30.c", printCtx->print_entropy, &printCtx->print_entropy, NULL));
4621b431c67SMark Adams   PetscOptionsEnd();
4633cdbf90fSMark Adams   // view
464f3237bfeSMark Adams   PetscCall(DMViewFromOptions(ctx->plex[g_target], NULL, "-ex30_dm_view"));
4651b431c67SMark Adams   if (ctx->num_grids > g_target + 1) { PetscCall(DMViewFromOptions(ctx->plex[g_target + 1], NULL, "-ex30_dm_view2")); }
4663cdbf90fSMark Adams   // create mesh mass matrices
4676b664d00Smarkadams4   PetscCall(VecZeroEntries(X));
468f3237bfeSMark Adams   PetscCall(DMCompositeGetAccessArray(pack, X, nDMs, NULL, globXArray)); // just to duplicate
469f3237bfeSMark Adams   for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {               // add same particels for all grids
470f3237bfeSMark Adams     Vec          subX = globXArray[LAND_PACK_IDX(0, grid)];
471f3237bfeSMark Adams     DM           dm   = ctx->plex[grid];
472f3237bfeSMark Adams     PetscSection s;
473f3237bfeSMark Adams     grid_dm[grid] = dm;
474f3237bfeSMark Adams     PetscCall(DMCreateMassMatrix(dm, dm, &g_Mass[grid]));
475f3237bfeSMark Adams     //
476f3237bfeSMark Adams     PetscCall(DMGetLocalSection(dm, &s));
477f3237bfeSMark Adams     PetscCall(DMPlexCreateClosureIndex(dm, s));
478f3237bfeSMark Adams     for (int tid = 0; tid < numthreads; tid++) {
4791b431c67SMark Adams       PC pc;
480f3237bfeSMark Adams       PetscCall(VecDuplicate(subX, &t_fhat[grid][tid]));
481f3237bfeSMark Adams       PetscCall(KSPCreate(PETSC_COMM_SELF, &t_ksp[grid][tid]));
4821b431c67SMark Adams       PetscCall(KSPSetType(t_ksp[grid][tid], KSPCG));
4831b431c67SMark Adams       PetscCall(KSPGetPC(t_ksp[grid][tid], &pc));
4841b431c67SMark Adams       PetscCall(PCSetType(pc, PCJACOBI));
485f3237bfeSMark Adams       PetscCall(KSPSetOptionsPrefix(t_ksp[grid][tid], "ptof_"));
486f3237bfeSMark Adams       PetscCall(KSPSetOperators(t_ksp[grid][tid], g_Mass[grid], g_Mass[grid]));
487f3237bfeSMark Adams       PetscCall(KSPSetFromOptions(t_ksp[grid][tid]));
488f3237bfeSMark Adams     }
489f3237bfeSMark Adams   }
490f3237bfeSMark Adams   PetscCall(DMCompositeRestoreAccessArray(pack, X, nDMs, NULL, globXArray));
491f3237bfeSMark Adams   PetscCall(TSGetTimeStep(ts, &dt_init)); // we could have an adaptive time stepper
4923cdbf90fSMark Adams   // loop over all vertices in chucks that are batched for TSSolve
4931b431c67SMark Adams   for (int i = 0; i < 5; i++) moments_0[i] = moments_1a[i] = moments_1b[i] = 0;
4943cdbf90fSMark Adams   for (PetscInt global_vertex_id_0 = 0; global_vertex_id_0 < num_vertices; global_vertex_id_0 += ctx->batch_sz, shift /= 2) { // outer vertex loop
4956b664d00Smarkadams4     PetscCall(TSSetTime(ts, 0));
4966b664d00Smarkadams4     PetscCall(TSSetStepNumber(ts, 0));
4976b664d00Smarkadams4     PetscCall(TSSetTimeStep(ts, dt_init));
498f3237bfeSMark Adams     PetscCall(DMCompositeGetAccessArray(pack, X, nDMs, NULL, globXArray));
4991b431c67SMark Adams     printCtx->global_vertex_id_0 = global_vertex_id_0;
5003cdbf90fSMark Adams     if (v_target >= global_vertex_id_0 && v_target < global_vertex_id_0 + ctx->batch_sz) {
5013cdbf90fSMark Adams       PetscCall(PetscObjectSetName((PetscObject)globXArray[LAND_PACK_IDX(v_target % ctx->batch_sz, g_target)], "rho"));
5023cdbf90fSMark Adams       printCtx->print = PETSC_TRUE;
5033cdbf90fSMark Adams     } else printCtx->print = PETSC_FALSE;
5043cdbf90fSMark Adams     // create fake particles in batches with threads
5053cdbf90fSMark Adams     for (PetscInt v_id_0 = 0; v_id_0 < ctx->batch_sz; v_id_0 += numthreads) {
5063cdbf90fSMark Adams       PetscReal *xx_t[LANDAU_MAX_GRIDS][EX30_MAX_NUM_THRDS], *yy_t[LANDAU_MAX_GRIDS][EX30_MAX_NUM_THRDS], *zz_t[LANDAU_MAX_GRIDS][EX30_MAX_NUM_THRDS], *wp_t[LANDAU_MAX_GRIDS][EX30_MAX_NUM_THRDS];
5073cdbf90fSMark Adams       PetscInt   Np_t[LANDAU_MAX_GRIDS][EX30_MAX_NUM_THRDS];
508f3237bfeSMark Adams       // make particles
509f3237bfeSMark Adams       for (int tid = 0; tid < numthreads; tid++) {
5103cdbf90fSMark Adams         const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id;
5113cdbf90fSMark Adams         if (glb_v_id < num_vertices) {                                                                                                                                            // the ragged edge (in last batch)
5121b431c67SMark Adams           PetscInt Npp0 = a_Np + (glb_v_id % (a_Np / 10 + 1)), nTargetP[LANDAU_MAX_GRIDS];                                                                                        // n of particels in each dim with load imbalance
513f3237bfeSMark Adams           for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {                                                                                                                // add same particels for all grids
5143cdbf90fSMark Adams             const PetscReal kT_m  = ctx->k * ctx->thermal_temps[ctx->species_offset[grid]] / ctx->masses[ctx->species_offset[grid]] / (ctx->v_0 * ctx->v_0);                      /* theta = 2kT/mc^2 per species */
515f3237bfeSMark Adams             PetscReal       lo[3] = {-ctx->radius[grid], -ctx->radius[grid], -ctx->radius[grid]}, hi[3] = {ctx->radius[grid], ctx->radius[grid], ctx->radius[grid]}, hp[3], vole; // would be nice to get box from DM
516f3237bfeSMark Adams             PetscInt        Npi = Npp0, Npj = 2 * Npp0, Npk = 1;
5171b431c67SMark Adams             PetscRandom     rand;
5181b431c67SMark Adams             PetscReal       sigma = ctx->thermal_speed[grid] / ctx->thermal_speed[0], p2_shift = grid == 0 ? shift : -shift;
5191b431c67SMark Adams             PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rand));
5201b431c67SMark Adams             PetscCall(PetscRandomSetInterval(rand, 0., 1.));
5211b431c67SMark Adams             PetscCall(PetscRandomSetFromOptions(rand));
522f3237bfeSMark Adams             if (dim == 2) lo[0] = 0; // Landau coordinate (r,z)
523f3237bfeSMark Adams             else Npi = Npj = Npk = Npp0;
5243cdbf90fSMark Adams             // User: use glb_v_id to index into your data
5251b431c67SMark Adams             const PetscInt NN = Npi * Npj * Npk; // make a regular grid of particles Npp x Npp
526f3237bfeSMark Adams             Np_t[grid][tid]   = NN;
5273cdbf90fSMark Adams             if (glb_v_id == v_target) nTargetP[grid] = NN;
528f3237bfeSMark Adams             PetscCall(PetscMalloc4(NN, &xx_t[grid][tid], NN, &yy_t[grid][tid], NN, &wp_t[grid][tid], dim == 2 ? 1 : NN, &zz_t[grid][tid]));
529f3237bfeSMark Adams             hp[0] = (hi[0] - lo[0]) / Npi;
530f3237bfeSMark Adams             hp[1] = (hi[1] - lo[1]) / Npj;
531f3237bfeSMark Adams             hp[2] = (hi[2] - lo[2]) / Npk;
532f3237bfeSMark Adams             if (dim == 2) hp[2] = 1;
53363a3b9bcSJacob Faibussowitsch             PetscCall(PetscInfo(pack, " lo = %14.7e, hi = %14.7e; hp = %14.7e, %14.7e; kT_m = %g; \n", (double)lo[1], (double)hi[1], (double)hp[0], (double)hp[1], (double)kT_m)); // temp
534f3237bfeSMark Adams             vole = hp[0] * hp[1] * hp[2] * ctx->n[grid];                                                                                                                           // fix for multi-species
5353cdbf90fSMark Adams             PetscCall(PetscInfo(pack, "Vertex %" PetscInt_FMT ", grid %" PetscInt_FMT " with %" PetscInt_FMT " particles (diagnostic target = %" PetscInt_FMT ")\n", glb_v_id, grid, NN, v_target));
536f3237bfeSMark Adams             for (int pj = 0, pp = 0; pj < Npj; pj++) {
537f3237bfeSMark Adams               for (int pk = 0; pk < Npk; pk++) {
538f3237bfeSMark Adams                 for (int pi = 0; pi < Npi; pi++, pp++) {
5391b431c67SMark Adams                   PetscReal p_shift   = p2_shift;
5401b431c67SMark Adams                   wp_t[grid][tid][pp] = 0;
5411b431c67SMark Adams                   if (use_uniform_particle_grid) {
542f3237bfeSMark Adams                     xx_t[grid][tid][pp] = lo[0] + hp[0] / 2.0 + pi * hp[0];
543f3237bfeSMark Adams                     yy_t[grid][tid][pp] = lo[1] + hp[1] / 2.0 + pj * hp[1];
544f3237bfeSMark Adams                     if (dim == 3) zz_t[grid][tid][pp] = lo[2] + hp[2] / 2.0 + pk * hp[2];
545f3237bfeSMark Adams                     PetscReal x[] = {xx_t[grid][tid][pp], yy_t[grid][tid][pp], dim == 2 ? 0 : zz_t[grid][tid][pp]};
5461b431c67SMark Adams                     p_shift *= ctx->thermal_speed[grid] / ctx->v_0;
5471b431c67SMark Adams                     maxwellian(dim, x, kT_m, vole, p_shift, &wp_t[grid][tid][pp]);
5481b431c67SMark Adams                     if (ctx->num_grids == 1 && shift != 0) { // bi-maxwellian, electron plasma
5491b431c67SMark Adams                       maxwellian(dim, x, kT_m, vole, -p_shift, &wp_t[grid][tid][pp]);
5501b431c67SMark Adams                     }
5511b431c67SMark Adams                   } else {
5521b431c67SMark Adams                     PetscReal u1, u2;
5531b431c67SMark Adams                     do {
5541b431c67SMark Adams                       do {
5551b431c67SMark Adams                         PetscCall(PetscRandomGetValueReal(rand, &u1));
5561b431c67SMark Adams                       } while (u1 == 0);
5571b431c67SMark Adams                       PetscCall(PetscRandomGetValueReal(rand, &u2));
5581b431c67SMark Adams                       //compute z0 and z1
5591b431c67SMark Adams                       PetscReal mag       = sigma * PetscSqrtReal(-2.0 * PetscLogReal(u1));
5601b431c67SMark Adams                       xx_t[grid][tid][pp] = mag * PetscCosReal(2.0 * PETSC_PI * u2);
5611b431c67SMark Adams                       yy_t[grid][tid][pp] = mag * PetscSinReal(2.0 * PETSC_PI * u2);
5621b431c67SMark Adams                       if (dim == 2 && xx_t[grid][tid][pp] < lo[0]) xx_t[grid][tid][pp] = -xx_t[grid][tid][pp];
5631b431c67SMark Adams                       if (dim == 3) zz_t[grid][tid][pp] = lo[2] + hp[2] / 2.0 + pk * hp[2];
5641b431c67SMark Adams                       if (ctx->num_grids == 1 && pp % 2 == 0) p_shift = 0; // one species, split bi-max
5651b431c67SMark Adams                       p_shift *= ctx->thermal_speed[grid] / ctx->v_0;
5661b431c67SMark Adams                       if (dim == 3) zz_t[grid][tid][pp] += p_shift;
5671b431c67SMark Adams                       else yy_t[grid][tid][pp] += p_shift;
5681b431c67SMark Adams                       wp_t[grid][tid][pp] += ctx->n[grid] / NN * PetscSqrtReal(ctx->masses[ctx->species_offset[grid]] / ctx->masses[0]);
5691b431c67SMark Adams                       if (p_shift <= 0) break; // add bi-max for electron plasma only
5701b431c67SMark Adams                       p_shift = -p_shift;
5711b431c67SMark Adams                     } while (ctx->num_grids == 1); // add bi-max for electron plasma only
5721b431c67SMark Adams                   }
5731b431c67SMark Adams                   {
5743cdbf90fSMark Adams                     if (glb_v_id == v_target) {
5751b431c67SMark Adams                       PetscReal x[] = {xx_t[grid][tid][pp], yy_t[grid][tid][pp], dim == 2 ? 0 : zz_t[grid][tid][pp]};
5763cdbf90fSMark Adams                       PetscReal v2 = 0, fact = dim == 2 ? 2.0 * PETSC_PI * x[0] : 1, w = fact * wp_t[grid][tid][pp] * ctx->n_0 * ctx->masses[ctx->species_offset[grid]];
577f3237bfeSMark Adams                       for (int i = 0; i < dim; ++i) v2 += PetscSqr(x[i]);
5783cdbf90fSMark Adams                       moments_0[0] += w;                   // not thread safe
5793cdbf90fSMark Adams                       moments_0[1] += w * ctx->v_0 * x[1]; // z-momentum
5801b431c67SMark Adams                       moments_0[2] += w * 0.5 * ctx->v_0 * ctx->v_0 * v2;
581f3237bfeSMark Adams                     }
582f3237bfeSMark Adams                   }
583f3237bfeSMark Adams                 }
584f3237bfeSMark Adams               }
585f3237bfeSMark Adams             }
5861b431c67SMark Adams             PetscCall(PetscRandomDestroy(&rand));
5873cdbf90fSMark Adams           }
5881b431c67SMark Adams           // entropy init, need global n
5893cdbf90fSMark Adams           if (glb_v_id == v_target) {
5901b431c67SMark Adams             const PetscReal N_inv = 1 / moments_0[0];
5913cdbf90fSMark Adams             PetscCall(PetscInfo(pack, "Target %" PetscInt_FMT " with %" PetscInt_FMT " particels\n", glb_v_id, nTargetP[0]));
5923cdbf90fSMark Adams             for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
5931b431c67SMark Adams               const PetscInt NN = nTargetP[grid];
5943cdbf90fSMark Adams               for (int pp = 0; pp < NN; pp++) {
5951b431c67SMark Adams                 const PetscReal fact = dim == 2 ? 2.0 * PETSC_PI * xx_t[grid][tid][pp] : 1, w = fact * ctx->n_0 * ctx->masses[ctx->species_offset[grid]] * wp_t[grid][tid][pp], ww = w * N_inv;
5961b431c67SMark Adams                 if (w > PETSC_REAL_MIN) {
5973cdbf90fSMark Adams                   moments_0[3] -= ww * PetscLogReal(ww);
5983cdbf90fSMark Adams                   PetscCheck(ww < 1 - PETSC_MACHINE_EPSILON, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "ww (%g) > 1", (double)ww);
5991b431c67SMark Adams                 } else moments_0[4] -= w;
6003cdbf90fSMark Adams               }
601f3237bfeSMark Adams             } // grid
6021b431c67SMark Adams           } // target
603f3237bfeSMark Adams         } // active
6043cdbf90fSMark Adams       } // threads
605f3237bfeSMark Adams       /* Create particle swarm */
606f37bacd1SJacob Faibussowitsch       for (int tid = 0; tid < numthreads; tid++) {
6073cdbf90fSMark Adams         const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id;
6083cdbf90fSMark Adams         if (glb_v_id < num_vertices) {                             // the ragged edge of the last batch
609f3237bfeSMark Adams           for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids
610f3237bfeSMark Adams             PetscErrorCode ierr_t;
611f3237bfeSMark Adams             PetscSection   section;
612f3237bfeSMark Adams             PetscInt       Nf;
613f3237bfeSMark Adams             DM             dm = grid_dm[grid];
614f3237bfeSMark Adams             ierr_t            = DMGetLocalSection(dm, &section);
615f3237bfeSMark Adams             ierr_t            = PetscSectionGetNumFields(section, &Nf);
6163cdbf90fSMark Adams             if (Nf != 1) ierr_t = (PetscErrorCode)9999;
617f3237bfeSMark Adams             else {
618f3237bfeSMark Adams               ierr_t = DMViewFromOptions(dm, NULL, "-dm_view");
6191b431c67SMark Adams               ierr_t = PetscInfo(pack, "call createSwarm [%" PetscInt_FMT ".%" PetscInt_FMT "] local block index %" PetscInt_FMT "\n", v_id, grid, LAND_PACK_IDX(v_id, grid));
6203cdbf90fSMark Adams               ierr_t = createSwarm(dm, dim, &globSwarmArray[LAND_PACK_IDX(v_id, grid)]);
621f3237bfeSMark Adams             }
622f3237bfeSMark Adams             if (ierr_t) ierr = ierr_t;
623f3237bfeSMark Adams           }
624f3237bfeSMark Adams         } // active
6253cdbf90fSMark Adams       } // threads
626cad9d221SBarry Smith       PetscCheck(ierr != 9999, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Only support one species per grid");
6276b664d00Smarkadams4       PetscCheck(!ierr, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Error in OMP loop. ierr = %d", (int)ierr);
6283cdbf90fSMark Adams       // make globMpArray
629f3237bfeSMark Adams       PetscPragmaOMP(parallel for)
630f37bacd1SJacob Faibussowitsch       for (int tid = 0; tid < numthreads; tid++) {
6313cdbf90fSMark Adams         const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id;
6323cdbf90fSMark Adams         if (glb_v_id < num_vertices) {
6333cdbf90fSMark Adams           for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids
6343cdbf90fSMark Adams             PetscErrorCode ierr_t;
6353cdbf90fSMark Adams             DM             sw = globSwarmArray[LAND_PACK_IDX(v_id, grid)];
6361b431c67SMark Adams             ierr_t            = PetscInfo(pack, "makeSwarm %" PetscInt_FMT ".%" PetscInt_FMT ") for block %" PetscInt_FMT "\n", v_id, grid, LAND_PACK_IDX(v_id, grid));
6373cdbf90fSMark Adams             ierr_t            = makeSwarm(sw, dim, Np_t[grid][tid], xx_t[grid][tid], yy_t[grid][tid], zz_t[grid][tid]);
6383cdbf90fSMark Adams             if (ierr_t) ierr = ierr_t;
6393cdbf90fSMark Adams           }
6403cdbf90fSMark Adams         }
6413cdbf90fSMark Adams       }
6423cdbf90fSMark Adams       //PetscPragmaOMP(parallel for)
6433cdbf90fSMark Adams       for (int tid = 0; tid < numthreads; tid++) {
6443cdbf90fSMark Adams         const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id;
6453cdbf90fSMark Adams         if (glb_v_id < num_vertices) {
646f3237bfeSMark Adams           for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids
647f3237bfeSMark Adams             PetscErrorCode ierr_t;
648f3237bfeSMark Adams             DM             dm = grid_dm[grid];
6493cdbf90fSMark Adams             DM             sw = globSwarmArray[LAND_PACK_IDX(v_id, grid)];
6501b431c67SMark Adams             ierr_t            = PetscInfo(pack, "createMp %" PetscInt_FMT ".%" PetscInt_FMT ") for block %" PetscInt_FMT "\n", v_id, grid, LAND_PACK_IDX(v_id, grid));
6513cdbf90fSMark Adams             ierr_t            = createMp(dm, sw, &globMpArray[LAND_PACK_IDX(v_id, grid)]);
6523cdbf90fSMark Adams             if (ierr_t) ierr = ierr_t;
6533cdbf90fSMark Adams           }
6543cdbf90fSMark Adams         }
6553cdbf90fSMark Adams       }
6563cdbf90fSMark Adams       PetscCheck(!ierr, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Error in OMP loop. ierr = %d", (int)ierr);
6573cdbf90fSMark Adams       // p --> g: set X
6583cdbf90fSMark Adams       // PetscPragmaOMP(parallel for)
6593cdbf90fSMark Adams       for (int tid = 0; tid < numthreads; tid++) {
6603cdbf90fSMark Adams         const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id;
6613cdbf90fSMark Adams         if (glb_v_id < num_vertices) {
6623cdbf90fSMark Adams           for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids
6633cdbf90fSMark Adams             PetscErrorCode ierr_t;
6643cdbf90fSMark Adams             DM             dm   = grid_dm[grid];
6653cdbf90fSMark Adams             DM             sw   = globSwarmArray[LAND_PACK_IDX(v_id, grid)];
6663cdbf90fSMark Adams             Vec            subX = globXArray[LAND_PACK_IDX(v_id, grid)], work = t_fhat[grid][tid];
6671b431c67SMark Adams             ierr_t = PetscInfo(pack, "particlesToGrid %" PetscInt_FMT ".%" PetscInt_FMT ") for block %" PetscInt_FMT "\n", v_id, grid, LAND_PACK_IDX(v_id, grid));
6681b431c67SMark Adams             ierr_t = particlesToGrid(dm, sw, tid, dim, wp_t[grid][tid], subX, globMpArray[LAND_PACK_IDX(v_id, grid)]);
669f3237bfeSMark Adams             if (ierr_t) ierr = ierr_t;
670f3237bfeSMark Adams             // u = M^_1 f_w
671f3237bfeSMark Adams             ierr_t = VecCopy(subX, work);
672f3237bfeSMark Adams             ierr_t = KSPSolve(t_ksp[grid][tid], work, subX);
673f3237bfeSMark Adams             if (ierr_t) ierr = ierr_t;
674f3237bfeSMark Adams           }
675f3237bfeSMark Adams         }
676f3237bfeSMark Adams       }
6776b664d00Smarkadams4       PetscCheck(!ierr, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Error in OMP loop. ierr = %d", (int)ierr);
678f3237bfeSMark Adams       /* Cleanup */
679f3237bfeSMark Adams       for (int tid = 0; tid < numthreads; tid++) {
6803cdbf90fSMark Adams         const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id;
6813cdbf90fSMark Adams         if (glb_v_id < num_vertices) {
682f3237bfeSMark Adams           for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids
683f3237bfeSMark Adams             PetscCall(PetscFree4(xx_t[grid][tid], yy_t[grid][tid], wp_t[grid][tid], zz_t[grid][tid]));
684f3237bfeSMark Adams           }
685f3237bfeSMark Adams         } // active
6863cdbf90fSMark Adams       } // threads
6873cdbf90fSMark Adams     } // (fake) particle loop
6883cdbf90fSMark Adams     // standard view of initial conditions
6893cdbf90fSMark Adams     if (v_target >= global_vertex_id_0 && v_target < global_vertex_id_0 + ctx->batch_sz) {
6903cdbf90fSMark Adams       PetscCall(DMSetOutputSequenceNumber(ctx->plex[g_target], 0, 0.0));
6913cdbf90fSMark Adams       PetscCall(VecViewFromOptions(globXArray[LAND_PACK_IDX(v_target % ctx->batch_sz, g_target)], NULL, "-ex30_vec_view"));
6921b431c67SMark Adams       if (ctx->num_grids > g_target + 1) {
6931b431c67SMark Adams         PetscCall(DMSetOutputSequenceNumber(ctx->plex[g_target + 1], 0, 0.0));
6941b431c67SMark Adams         PetscCall(VecViewFromOptions(globXArray[LAND_PACK_IDX(v_target % ctx->batch_sz, g_target + 1)], NULL, "-ex30_vec_view2"));
6951b431c67SMark Adams       }
696f3237bfeSMark Adams     }
6973cdbf90fSMark Adams     // coarse graining moments
6981b431c67SMark Adams     if (v_target >= global_vertex_id_0 && v_target < global_vertex_id_0 + ctx->batch_sz && printCtx->print_entropy) {
6993cdbf90fSMark Adams       const PetscInt v_id = v_target % ctx->batch_sz;
700f3237bfeSMark Adams       for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
701f3237bfeSMark Adams         PetscDataType dtype;
702f3237bfeSMark Adams         PetscReal    *wp, *coords;
7033cdbf90fSMark Adams         DM            sw = globSwarmArray[LAND_PACK_IDX(v_id, grid)];
7043cdbf90fSMark Adams         Vec           work, subX = globXArray[LAND_PACK_IDX(v_id, grid)];
7051b431c67SMark Adams         PetscInt      bs, NN;
7063cdbf90fSMark Adams         // C-G moments
7073cdbf90fSMark Adams         PetscCall(VecDuplicate(subX, &work));
7083cdbf90fSMark Adams         PetscCall(gridToParticles(grid_dm[grid], sw, subX, work, globMpArray[LAND_PACK_IDX(v_id, grid)], g_Mass[grid]));
7093cdbf90fSMark Adams         PetscCall(VecDestroy(&work));
7103cdbf90fSMark Adams         // moments
711f3237bfeSMark Adams         PetscCall(DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
7121b431c67SMark Adams         PetscCall(DMSwarmGetLocalSize(sw, &NN));
7131b431c67SMark Adams         PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void **)&wp));
7143cdbf90fSMark Adams         for (int pp = 0; pp < NN; pp++) {
7151b431c67SMark Adams           PetscReal v2 = 0, fact = (dim == 2) ? 2.0 * PETSC_PI * coords[pp * dim + 0] : 1, w = fact * wp[pp] * ctx->n_0 * ctx->masses[ctx->species_offset[grid]];
7163cdbf90fSMark Adams           for (int i = 0; i < dim; ++i) v2 += PetscSqr(coords[pp * dim + i]);
7173cdbf90fSMark Adams           moments_1a[0] += w;
7183cdbf90fSMark Adams           moments_1a[1] += w * ctx->v_0 * coords[pp * dim + 1]; // z-momentum
7191b431c67SMark Adams           moments_1a[2] += w * 0.5 * ctx->v_0 * ctx->v_0 * v2;
720f3237bfeSMark Adams         }
721f3237bfeSMark Adams         PetscCall(DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
722f3237bfeSMark Adams         PetscCall(DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void **)&wp));
723f3237bfeSMark Adams       }
7241b431c67SMark Adams       // entropy
7251b431c67SMark Adams       const PetscReal N_inv = 1 / moments_1a[0];
7261b431c67SMark Adams       PetscCall(PetscInfo(pack, "Entropy batch %" PetscInt_FMT " of %" PetscInt_FMT ", n = %g\n", v_target, num_vertices, (double)(1 / N_inv)));
7271b431c67SMark Adams       for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
7281b431c67SMark Adams         PetscDataType dtype;
7291b431c67SMark Adams         PetscReal    *wp, *coords;
7301b431c67SMark Adams         DM            sw = globSwarmArray[LAND_PACK_IDX(v_id, grid)];
7311b431c67SMark Adams         PetscInt      bs, NN;
7321b431c67SMark Adams         PetscCall(DMSwarmGetLocalSize(sw, &NN));
7331b431c67SMark Adams         PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void **)&wp));
7341b431c67SMark Adams         PetscCall(DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
7351b431c67SMark Adams         for (int pp = 0; pp < NN; pp++) {
7361b431c67SMark Adams           PetscReal fact = (dim == 2) ? 2.0 * PETSC_PI * coords[pp * dim + 0] : 1, w = fact * wp[pp] * ctx->n_0 * ctx->masses[ctx->species_offset[grid]], ww = w * N_inv;
7371b431c67SMark Adams           if (w > PETSC_REAL_MIN) {
7381b431c67SMark Adams             moments_1a[3] -= ww * PetscLogReal(ww);
7391b431c67SMark Adams             PetscCheck(ww < 1 - PETSC_MACHINE_EPSILON, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "ww (%g) > 1", (double)ww);
7401b431c67SMark Adams           } else moments_1a[4] -= w;
7411b431c67SMark Adams         }
7421b431c67SMark Adams         PetscCall(DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void **)&wp));
7431b431c67SMark Adams         PetscCall(DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
7441b431c67SMark Adams       }
745f3237bfeSMark Adams     }
7463cdbf90fSMark Adams     // restore vector
747f3237bfeSMark Adams     PetscCall(DMCompositeRestoreAccessArray(pack, X, nDMs, NULL, globXArray));
7483cdbf90fSMark Adams     // view
7491b431c67SMark Adams     if (v_target >= global_vertex_id_0 && v_target < global_vertex_id_0 + ctx->batch_sz) { PetscCall(DMPlexLandauPrintNorms(X, 0)); }
7503cdbf90fSMark Adams     // advance
7513cdbf90fSMark Adams     PetscCall(TSSetSolution(ts, X));
7521b431c67SMark Adams     PetscCall(PetscInfo(pack, "Advance vertex %" PetscInt_FMT " to %" PetscInt_FMT "\n", global_vertex_id_0, global_vertex_id_0 + ctx->batch_sz));
7533cdbf90fSMark Adams     PetscCall(TSSetPostStep(ts, PostStep));
7543cdbf90fSMark Adams     PetscCall(PostStep(ts));
7553cdbf90fSMark Adams     PetscCall(TSSolve(ts, X));
7563cdbf90fSMark Adams     // view
7573cdbf90fSMark Adams     PetscCall(DMCompositeGetAccessArray(pack, X, nDMs, NULL, globXArray));
7583cdbf90fSMark Adams     if (v_target >= global_vertex_id_0 && v_target < global_vertex_id_0 + ctx->batch_sz) {
7591b431c67SMark Adams       DM        sw = globSwarmArray[LAND_PACK_IDX(v_target % ctx->batch_sz, g_target)];
7601b431c67SMark Adams       PetscInt  id;
7611b431c67SMark Adams       PetscReal t;
7621b431c67SMark Adams       PetscCall(DMPlexLandauPrintNorms(X, 1));
7631b431c67SMark Adams       PetscCall(TSGetTime(ts, &t));
7641b431c67SMark Adams       PetscCall(DMGetOutputSequenceNumber(ctx->plex[g_target], &id, NULL));
7651b431c67SMark Adams       PetscCall(DMSetOutputSequenceNumber(ctx->plex[g_target], id + 1, t));
7663cdbf90fSMark Adams       PetscCall(VecViewFromOptions(globXArray[LAND_PACK_IDX(v_target % ctx->batch_sz, g_target)], NULL, "-ex30_vec_view"));
7671b431c67SMark Adams       if (ctx->num_grids > g_target + 1) {
7681b431c67SMark Adams         PetscCall(DMSetOutputSequenceNumber(ctx->plex[g_target + 1], id + 1, t));
7691b431c67SMark Adams         PetscCall(VecViewFromOptions(globXArray[LAND_PACK_IDX(v_target % ctx->batch_sz, g_target + 1)], NULL, "-ex30_vec_view2"));
7701b431c67SMark Adams       }
7711b431c67SMark Adams       /* Visualize particle field */
7721b431c67SMark Adams       Vec f;
7731b431c67SMark Adams       PetscCall(DMSetOutputSequenceNumber(sw, 0, 0.0));
7741b431c67SMark Adams       PetscCall(DMViewFromOptions(grid_dm[g_target], NULL, "-weights_dm_view"));
7751b431c67SMark Adams       PetscCall(DMViewFromOptions(sw, NULL, "-weights_sw_view"));
7761b431c67SMark Adams       PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &f));
7771b431c67SMark Adams       PetscCall(PetscObjectSetName((PetscObject)f, "weights"));
7781b431c67SMark Adams       PetscCall(VecViewFromOptions(f, NULL, "-weights_vec_view"));
7791b431c67SMark Adams       PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &f));
7803cdbf90fSMark Adams     }
7813cdbf90fSMark Adams     // particles to grid, compute moments and entropy
7821b431c67SMark Adams     if (v_target >= global_vertex_id_0 && v_target < global_vertex_id_0 + ctx->batch_sz && printCtx->print_entropy) {
7833cdbf90fSMark Adams       PetscCall(gridToParticles_private(grid_dm, globSwarmArray, dim, v_target, numthreads, num_vertices, global_vertex_id_0, globMpArray, g_Mass, t_fhat, moments_1b, globXArray, ctx));
7841b431c67SMark Adams       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Particle Moments:\t number density      momentum (par)     energy             entropy      loss : # OMP threads %g\n", (double)numthreads));
7851b431c67SMark Adams       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\tInitial:         %18.12e %19.12e %18.12e %e %e\n", (double)moments_0[0], (double)moments_0[1], (double)moments_0[2], (double)moments_0[3], (double)moments_0[4]));
7861b431c67SMark Adams       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\tCoarse-graining: %18.12e %19.12e %18.12e %e %e\n", (double)moments_1a[0], (double)moments_1a[1], (double)moments_1a[2], (double)moments_1a[3], (double)moments_1a[4]));
7871b431c67SMark Adams       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\tLandau:          %18.12e %19.12e %18.12e %e %e\n", (double)moments_1b[0], (double)moments_1b[1], (double)moments_1b[2], (double)moments_1b[3], (double)moments_1b[4]));
7881b431c67SMark Adams       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Coarse-graining entropy generation = %e ; Landau entropy generation = %e\n", (double)(moments_1a[3] - moments_0[3]), (double)(moments_1b[3] - moments_0[3])));
7891b431c67SMark Adams       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "(relative) energy conservation: Coarse-graining = %e ; Landau = %e\n", (double)(moments_1a[2] - moments_0[2]) / (double)moments_0[2], (double)(moments_1b[2] - moments_0[2]) / (double)moments_0[2]));
7901b431c67SMark Adams     }
7913cdbf90fSMark Adams     // restore vector
7923cdbf90fSMark Adams     PetscCall(DMCompositeRestoreAccessArray(pack, X, nDMs, NULL, globXArray));
7933cdbf90fSMark Adams     // cleanup
7943cdbf90fSMark Adams     for (PetscInt v_id_0 = 0; v_id_0 < ctx->batch_sz; v_id_0 += numthreads) {
7953cdbf90fSMark Adams       for (int tid = 0; tid < numthreads; tid++) {
7963cdbf90fSMark Adams         const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id;
7973cdbf90fSMark Adams         if (glb_v_id < num_vertices) {
7983cdbf90fSMark Adams           for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
7993cdbf90fSMark Adams             PetscCall(DMDestroy(&globSwarmArray[LAND_PACK_IDX(v_id, grid)]));
8003cdbf90fSMark Adams             PetscCall(MatDestroy(&globMpArray[LAND_PACK_IDX(v_id, grid)]));
8013cdbf90fSMark Adams           }
8023cdbf90fSMark Adams         }
8033cdbf90fSMark Adams       }
8043cdbf90fSMark Adams     }
8051b431c67SMark Adams   } // user batch, not used
806f3237bfeSMark Adams   /* Cleanup */
807f3237bfeSMark Adams   PetscCall(PetscFree(globXArray));
808f3237bfeSMark Adams   PetscCall(PetscFree(globSwarmArray));
809f3237bfeSMark Adams   PetscCall(PetscFree(globMpArray));
8103cdbf90fSMark Adams   PetscCall(PetscFree(printCtx));
811f3237bfeSMark Adams   // clean up mass matrices
812f3237bfeSMark Adams   for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids
813f3237bfeSMark Adams     PetscCall(MatDestroy(&g_Mass[grid]));
814f3237bfeSMark Adams     for (int tid = 0; tid < numthreads; tid++) {
815f3237bfeSMark Adams       PetscCall(VecDestroy(&t_fhat[grid][tid]));
816f3237bfeSMark Adams       PetscCall(KSPDestroy(&t_ksp[grid][tid]));
817f3237bfeSMark Adams     }
818f3237bfeSMark Adams   }
8193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
820f3237bfeSMark Adams }
821f3237bfeSMark Adams 
822d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
823d71ae5a4SJacob Faibussowitsch {
824f3237bfeSMark Adams   DM         pack;
825f3237bfeSMark Adams   Vec        X;
8261b431c67SMark Adams   PetscInt   dim = 2, num_vertices = 1, Np = 10, v_target = 0, g_target = 0;
827f3237bfeSMark Adams   TS         ts;
828f3237bfeSMark Adams   Mat        J;
829f3237bfeSMark Adams   LandauCtx *ctx;
8303cdbf90fSMark Adams   PetscReal  shift                     = 0;
8311b431c67SMark Adams   PetscBool  use_uniform_particle_grid = PETSC_TRUE;
832f3237bfeSMark Adams 
833327415f7SBarry Smith   PetscFunctionBeginUser;
834f3237bfeSMark Adams   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
835f3237bfeSMark Adams   // process args
836d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_SELF, "", "Collision Options", "DMPLEX");
837f3237bfeSMark Adams   PetscCall(PetscOptionsInt("-dim", "Velocity space dimension", "ex30.c", dim, &dim, NULL));
8381b431c67SMark Adams   PetscCall(PetscOptionsInt("-number_spatial_vertices", "Number of user spatial vertices to be batched for Landau", "ex30.c", num_vertices, &num_vertices, NULL));
839f3237bfeSMark Adams   PetscCall(PetscOptionsInt("-number_particles_per_dimension", "Number of particles per grid, with slight modification per spatial vertex, in each dimension of base Cartesian grid", "ex30.c", Np, &Np, NULL));
8401b431c67SMark Adams   PetscCall(PetscOptionsBool("-use_uniform_particle_grid", "Use uniform particle grid", "ex30.c", use_uniform_particle_grid, &use_uniform_particle_grid, NULL));
841f45b553cSPierre Jolivet   PetscCall(PetscOptionsRangeInt("-vertex_view_target", "Global vertex for diagnostics", "ex30.c", v_target, &v_target, NULL, 0, num_vertices - 1));
8423cdbf90fSMark Adams   PetscCall(PetscOptionsReal("-e_shift", "Bim-Maxwellian shift", "ex30.c", shift, &shift, NULL));
8431b431c67SMark Adams   PetscCall(PetscOptionsInt("-grid_view_target", "Grid to view with diagnostics", "ex30.c", g_target, &g_target, NULL));
844d0609cedSBarry Smith   PetscOptionsEnd();
845f3237bfeSMark Adams   /* Create a mesh */
846f3237bfeSMark Adams   PetscCall(DMPlexLandauCreateVelocitySpace(PETSC_COMM_SELF, dim, "", &X, &J, &pack));
847f3237bfeSMark Adams   PetscCall(DMSetUp(pack));
848f3237bfeSMark Adams   PetscCall(DMSetOutputSequenceNumber(pack, 0, 0.0));
849f3237bfeSMark Adams   PetscCall(DMGetApplicationContext(pack, &ctx));
8501b431c67SMark Adams   PetscCheck(g_target < ctx->num_grids, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Grid to view %" PetscInt_FMT " should be < number of grids %" PetscInt_FMT, g_target, ctx->num_grids);
8511b431c67SMark Adams   PetscCheck(ctx->batch_view_idx == v_target % ctx->batch_sz, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Global view index %" PetscInt_FMT " mode batch size %" PetscInt_FMT " != ctx->batch_view_idx %" PetscInt_FMT, v_target, ctx->batch_sz, ctx->batch_view_idx);
852f3237bfeSMark Adams   /* Create timestepping solver context */
853f3237bfeSMark Adams   PetscCall(TSCreate(PETSC_COMM_SELF, &ts));
854f3237bfeSMark Adams   PetscCall(TSSetDM(ts, pack));
855f3237bfeSMark Adams   PetscCall(TSSetIFunction(ts, NULL, DMPlexLandauIFunction, NULL));
856f3237bfeSMark Adams   PetscCall(TSSetIJacobian(ts, J, J, DMPlexLandauIJacobian, NULL));
857f3237bfeSMark Adams   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
858f3237bfeSMark Adams   PetscCall(TSSetFromOptions(ts));
859f3237bfeSMark Adams   PetscCall(PetscObjectSetName((PetscObject)X, "X"));
8606b664d00Smarkadams4   // do particle advance
8611b431c67SMark Adams   PetscCall(go(ts, X, num_vertices, Np, dim, v_target, g_target, shift, use_uniform_particle_grid));
862984ed092SMark Adams   PetscCall(MatZeroEntries(J)); // need to zero out so as to not reuse it in Landau's logic
863f3237bfeSMark Adams   /* clean up */
864f3237bfeSMark Adams   PetscCall(DMPlexLandauDestroyVelocitySpace(&pack));
865f3237bfeSMark Adams   PetscCall(TSDestroy(&ts));
866f3237bfeSMark Adams   PetscCall(VecDestroy(&X));
867f3237bfeSMark Adams   PetscCall(PetscFinalize());
868f3237bfeSMark Adams   return 0;
869f3237bfeSMark Adams }
870f3237bfeSMark Adams 
871f3237bfeSMark Adams /*TEST
872f3237bfeSMark Adams 
873f3237bfeSMark Adams   build:
8743cdbf90fSMark Adams     requires: !complex
875f3237bfeSMark Adams 
876f3237bfeSMark Adams   testset:
877984ed092SMark Adams     requires: double defined(PETSC_USE_DMLANDAU_2D)
878f3237bfeSMark Adams     output_file: output/ex30_0.out
8793cdbf90fSMark Adams     args: -dim 2 -petscspace_degree 3 -dm_landau_num_species_grid 1,1,1 -dm_refine 1 -number_particles_per_dimension 10 -dm_plex_hash_location \
8801b431c67SMark Adams           -dm_landau_batch_size 4 -number_spatial_vertices 6 -vertex_view_target 5 -grid_view_target 1 -dm_landau_batch_view_idx 1 \
881f3237bfeSMark Adams           -dm_landau_n 1.000018,1,1e-6 -dm_landau_thermal_temps 2,1,1 -dm_landau_ion_masses 2,180 -dm_landau_ion_charges 1,18 \
8823cdbf90fSMark Adams           -ftop_ksp_rtol 1e-10 -ftop_ksp_type lsqr -ftop_pc_type bjacobi -ftop_sub_pc_factor_shift_type nonzero -ftop_sub_pc_type lu \
8831b431c67SMark Adams           -ksp_type preonly -pc_type lu -dm_landau_verbose 4 -print_entropy \
8843cdbf90fSMark Adams           -ptof_ksp_type cg -ptof_pc_type jacobi -ptof_ksp_rtol 1e-12\
885f3237bfeSMark Adams           -snes_converged_reason -snes_monitor -snes_rtol 1e-14 -snes_stol 1e-14 \
8863cdbf90fSMark Adams           -ts_dt 0.01 -ts_rtol 1e-1 -ts_exact_final_time stepover -ts_max_snes_failures -1 -ts_max_steps 1 -ts_monitor -ts_type beuler
887f3237bfeSMark Adams 
888f3237bfeSMark Adams     test:
889f3237bfeSMark Adams       suffix: cpu
890f3237bfeSMark Adams       args: -dm_landau_device_type cpu
891f3237bfeSMark Adams     test:
892f3237bfeSMark Adams       suffix: kokkos
893*bbdffb7fSJunchao Zhang       # failed on Sunspot@ALCF with sycl
894*bbdffb7fSJunchao Zhang       requires: kokkos_kernels !openmp !sycl
8953cdbf90fSMark Adams       args: -dm_landau_device_type kokkos -dm_mat_type aijkokkos -dm_vec_type kokkos -pc_type bjkokkos -pc_bjkokkos_ksp_type tfqmr -pc_bjkokkos_pc_type jacobi
896f3237bfeSMark Adams 
897984ed092SMark Adams   testset:
898984ed092SMark Adams     requires: double !defined(PETSC_USE_DMLANDAU_2D)
899984ed092SMark Adams     output_file: output/ex30_3d.out
9001b431c67SMark Adams     args: -dim 3 -petscspace_degree 2 -dm_landau_num_species_grid 1,1,1 -dm_refine 0 -number_particles_per_dimension 10 -dm_plex_hash_location \
9011b431c67SMark Adams           -dm_landau_batch_size 1 -number_spatial_vertices 1 -vertex_view_target 0 -grid_view_target 0 -dm_landau_batch_view_idx 0 \
902984ed092SMark Adams           -dm_landau_n 1.000018,1,1e-6 -dm_landau_thermal_temps 2,1,1 -dm_landau_ion_masses 2,180 -dm_landau_ion_charges 1,18 \
9031b431c67SMark Adams           -ftop_ksp_rtol 1e-12 -ksp_type preonly -pc_type lu \
9043cdbf90fSMark Adams           -ptof_ksp_type cg -ptof_pc_type jacobi -ptof_ksp_rtol 1e-12\
9053cdbf90fSMark Adams           -snes_converged_reason -snes_monitor -snes_rtol 1e-12 -snes_stol 1e-12\
9061b431c67SMark Adams           -ts_dt 0.1 -ts_exact_final_time stepover -ts_max_snes_failures -1 -ts_max_steps 1 -ts_monitor -ts_type beuler -print_entropy
907984ed092SMark Adams 
908984ed092SMark Adams     test:
909984ed092SMark Adams       suffix: cpu_3d
910984ed092SMark Adams       args: -dm_landau_device_type cpu
911984ed092SMark Adams     test:
912984ed092SMark Adams       suffix: kokkos_3d
9133cdbf90fSMark Adams       requires: kokkos_kernels !openmp
9143cdbf90fSMark Adams       args: -dm_landau_device_type kokkos -dm_mat_type aijkokkos -dm_vec_type kokkos -pc_type bjkokkos -pc_bjkokkos_ksp_type tfqmr -pc_bjkokkos_pc_type jacobi
9153cdbf90fSMark Adams 
9163cdbf90fSMark Adams   testset:
9173cdbf90fSMark Adams     requires: !complex double defined(PETSC_USE_DMLANDAU_2D) !cuda
9181b431c67SMark Adams     args: -dm_refine 1 -dm_landau_num_species_grid 1 -dm_landau_thermal_temps 1 -petscspace_degree 3 -snes_converged_reason -ts_type beuler -ts_dt .01 \
9191b431c67SMark Adams           -ts_max_steps 1 -ksp_type preonly -pc_type lu -snes_rtol 1e-12 -snes_stol 1e-12 -dm_landau_device_type cpu -number_particles_per_dimension 40 \
9201b431c67SMark Adams           -ptof_ksp_rtol 1e-12 -dm_landau_batch_size 4 -number_spatial_vertices 4 -grid_view_target 0 \
9211b431c67SMark Adams           -vertex_view_target 3 -dm_landau_batch_view_idx 3
9223cdbf90fSMark Adams     test:
9233cdbf90fSMark Adams       suffix: simple
9243cdbf90fSMark Adams       args: -ex30_dm_view
9253cdbf90fSMark Adams     test:
9261b431c67SMark Adams       requires: cgns
9271b431c67SMark Adams       suffix: cgns
9281b431c67SMark Adams       args: -ex30_vec_view cgns:cgnsDi.cgns
9291b431c67SMark Adams     test:
9301b431c67SMark Adams       suffix: normal
9311b431c67SMark Adams       args: -ex30_dm_view -use_uniform_particle_grid false
932984ed092SMark Adams 
933f3237bfeSMark Adams TEST*/
934