xref: /petsc/src/ts/tests/ex30.c (revision bd1587442888ecfa2176fd832616cd1e28a092f3)
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));
112d043ef4cSMark Adams   //
1133cdbf90fSMark Adams   PetscCall(PetscObjectSetName((PetscObject)sw, "Particle Grid"));
114d043ef4cSMark Adams 
115d043ef4cSMark Adams   /* PetscInt  N,*count,nmin=10000,nmax=0,ntot=0; */
116d043ef4cSMark Adams   /* // count */
117d043ef4cSMark Adams   /* PetscCall(DMSwarmCreatePointPerCellCount(sw, &N, &count)); */
118d043ef4cSMark Adams   /* for (int i=0, n; i< N ; i++) { */
119d043ef4cSMark Adams   /*   if ((n=count[i]) > nmax) nmax = n; */
120d043ef4cSMark Adams   /*   if (n < nmin) nmin = n; */
121d043ef4cSMark Adams   /*   PetscCall(PetscInfo(dm, " %d) %d particles\n", i, n)); */
122d043ef4cSMark Adams   /*   ntot += n; */
123d043ef4cSMark Adams   /* } */
124d043ef4cSMark Adams   /* PetscCall(PetscFree(count)); */
125d043ef4cSMark Adams   /* PetscCall(PetscInfo(dm, " %" PetscInt_FMT " max particle / cell, and %" PetscInt_FMT " min, ratio = %g,  %" PetscInt_FMT " total\n", nmax, nmin, (double)nmax/(double)nmin,ntot)); */
126d043ef4cSMark Adams 
1273cdbf90fSMark Adams   /* This gives M f = \int_\Omega \phi f, which looks like a rhs for a PDE */
1283cdbf90fSMark Adams   PetscCall(DMCreateMassMatrix(sw, dm, &M_p));
1293cdbf90fSMark Adams   PetscCall(DMViewFromOptions(sw, NULL, "-ex30_sw_view"));
1303cdbf90fSMark Adams   // output
1313cdbf90fSMark Adams   *Mp_out = M_p;
1323cdbf90fSMark Adams   PetscFunctionReturn(PETSC_SUCCESS);
1333cdbf90fSMark Adams }
1343cdbf90fSMark Adams 
1351b431c67SMark 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)
1363cdbf90fSMark Adams {
1373cdbf90fSMark Adams   PetscReal    *wq;
1383cdbf90fSMark Adams   PetscDataType dtype;
1393cdbf90fSMark Adams   Vec           ff;
1401b431c67SMark Adams   PetscInt      bs, p, Np;
1413cdbf90fSMark Adams 
1423cdbf90fSMark Adams   PetscFunctionBeginUser;
1433cdbf90fSMark Adams   PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void **)&wq));
1441b431c67SMark Adams   PetscCall(DMSwarmGetLocalSize(sw, &Np));
1453cdbf90fSMark Adams   for (p = 0; p < Np; p++) wq[p] = a_wp[p];
1463cdbf90fSMark Adams   PetscCall(DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void **)&wq));
1473cdbf90fSMark Adams   PetscCall(PetscObjectSetName((PetscObject)rho, "rho"));
1483cdbf90fSMark Adams   PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &ff));
1493cdbf90fSMark Adams   PetscCall(PetscObjectSetName((PetscObject)ff, "weights"));
1503cdbf90fSMark Adams   PetscCall(MatMultTranspose(M_p, ff, rho));
1513cdbf90fSMark Adams   PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &ff));
1523cdbf90fSMark Adams   PetscFunctionReturn(PETSC_SUCCESS);
1533cdbf90fSMark Adams }
1543cdbf90fSMark Adams 
1553cdbf90fSMark Adams //
1563cdbf90fSMark Adams // add grid to arg 'sw.w_q'
1573cdbf90fSMark Adams //
158d043ef4cSMark Adams PetscErrorCode gridToParticles(const DM dm, DM sw, const Vec rhs, Vec work_ferhs, Mat M_p, Mat Mass)
159d71ae5a4SJacob Faibussowitsch {
160f3237bfeSMark Adams   PetscBool    is_lsqr;
161f3237bfeSMark Adams   KSP          ksp;
162d043ef4cSMark Adams   Mat          PM_p = NULL, MtM, D = NULL;
163f3237bfeSMark Adams   Vec          ff;
164f3237bfeSMark Adams   PetscInt     N, M, nzl;
165d043ef4cSMark Adams   MatShellCtx *matshellctx = NULL;
1661b431c67SMark Adams   PC           pc;
167f3237bfeSMark Adams 
168f3237bfeSMark Adams   PetscFunctionBeginUser;
1693643261fSMark Adams   // 1) apply M in, for Moore-Penrose with mass: Mp (Mp' Mp)^-1 M
170d043ef4cSMark Adams   PetscCall(MatMult(Mass, rhs, work_ferhs));
1713643261fSMark Adams   // 2) pseudo-inverse, first part: (Mp' Mp)^-1
172f3237bfeSMark Adams   PetscCall(KSPCreate(PETSC_COMM_SELF, &ksp));
1731b431c67SMark Adams   PetscCall(KSPSetType(ksp, KSPCG));
1741b431c67SMark Adams   PetscCall(KSPGetPC(ksp, &pc));
1751b431c67SMark Adams   PetscCall(PCSetType(pc, PCJACOBI));
176f3237bfeSMark Adams   PetscCall(KSPSetOptionsPrefix(ksp, "ftop_"));
177f3237bfeSMark Adams   PetscCall(KSPSetFromOptions(ksp));
178f3237bfeSMark Adams   PetscCall(PetscObjectTypeCompare((PetscObject)ksp, KSPLSQR, &is_lsqr));
179f3237bfeSMark Adams   if (!is_lsqr) {
180f3237bfeSMark Adams     PetscCall(MatGetLocalSize(M_p, &M, &N));
181f3237bfeSMark Adams     if (N > M) {
1821b431c67SMark Adams       PetscCall(PetscInfo(ksp, " M (%" PetscInt_FMT ") < M (%" PetscInt_FMT ") more vertices than particles: revert to lsqr\n", M, N));
183f3237bfeSMark Adams       is_lsqr = PETSC_TRUE;
184f3237bfeSMark Adams       PetscCall(KSPSetType(ksp, KSPLSQR));
185d043ef4cSMark Adams       PetscCall(PCSetType(pc, PCNONE)); // should not happen, but could solve stable (Mp^T Mp), move projection Mp before solve
186f3237bfeSMark Adams     } else {
187f3237bfeSMark Adams       PetscCall(PetscNew(&matshellctx));
188d043ef4cSMark Adams       PetscCall(MatCreateVecs(M_p, &matshellctx->uu, &matshellctx->ff));
189d043ef4cSMark Adams       if (0) {
190d043ef4cSMark Adams         PetscCall(MatTransposeMatMult(M_p, M_p, MAT_INITIAL_MATRIX, 4, &MtM));
191d043ef4cSMark Adams         PetscCall(KSPSetOperators(ksp, MtM, MtM));
192d043ef4cSMark Adams         PetscCall(PetscInfo(M_p, "createMtM KSP with explicit Mp'Mp\n"));
193d043ef4cSMark Adams         PetscCall(MatViewFromOptions(MtM, NULL, "-ftop2_MtM_mat_view"));
194d043ef4cSMark Adams       } else {
195f3237bfeSMark Adams         PetscCall(MatCreateShell(PetscObjectComm((PetscObject)dm), N, N, PETSC_DECIDE, PETSC_DECIDE, matshellctx, &MtM));
196f3237bfeSMark Adams         PetscCall(MatTranspose(M_p, MAT_INITIAL_MATRIX, &matshellctx->MpTrans));
197f3237bfeSMark Adams         matshellctx->Mp = M_p;
198f3237bfeSMark Adams         PetscCall(MatShellSetOperation(MtM, MATOP_MULT, (void (*)(void))MatMultMtM_SeqAIJ));
199f3237bfeSMark Adams         PetscCall(MatShellSetOperation(MtM, MATOP_MULT_ADD, (void (*)(void))MatMultAddMtM_SeqAIJ));
200f3237bfeSMark Adams         PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, N, N, 1, NULL, &D));
201d043ef4cSMark Adams         PetscCall(MatViewFromOptions(matshellctx->MpTrans, NULL, "-ftop2_MpT_mat_view"));
2026497c311SBarry Smith         for (PetscInt i = 0; i < N; i++) {
203f3237bfeSMark Adams           const PetscScalar *vals;
204f3237bfeSMark Adams           const PetscInt    *cols;
205f3237bfeSMark Adams           PetscScalar        dot = 0;
206f3237bfeSMark Adams           PetscCall(MatGetRow(matshellctx->MpTrans, i, &nzl, &cols, &vals));
2076497c311SBarry Smith           for (PetscInt ii = 0; ii < nzl; ii++) dot += PetscSqr(vals[ii]);
208d043ef4cSMark Adams           if (dot < PETSC_MACHINE_EPSILON) {
2096497c311SBarry Smith             PetscCall(PetscInfo(ksp, "empty row in pseudo-inverse %d\n", (int)i));
210d043ef4cSMark Adams             is_lsqr = PETSC_TRUE; // empty rows
211d043ef4cSMark Adams             PetscCall(KSPSetType(ksp, KSPLSQR));
212d043ef4cSMark Adams             PetscCall(PCSetType(pc, PCNONE)); // should not happen, but could solve stable (Mp Mp^T), move projection Mp before solve
213d043ef4cSMark Adams             // clean up
214d043ef4cSMark Adams             PetscCall(MatDestroy(&matshellctx->MpTrans));
215d043ef4cSMark Adams             PetscCall(VecDestroy(&matshellctx->ff));
216d043ef4cSMark Adams             PetscCall(VecDestroy(&matshellctx->uu));
217d043ef4cSMark Adams             PetscCall(MatDestroy(&D));
218d043ef4cSMark Adams             PetscCall(MatDestroy(&MtM));
219d043ef4cSMark Adams             PetscCall(PetscFree(matshellctx));
220d043ef4cSMark Adams             D = NULL;
221d043ef4cSMark Adams             break;
222d043ef4cSMark Adams           }
223f3237bfeSMark Adams           PetscCall(MatSetValue(D, i, i, dot, INSERT_VALUES));
224f3237bfeSMark Adams         }
225d043ef4cSMark Adams         if (D) {
226f3237bfeSMark Adams           PetscCall(MatAssemblyBegin(D, MAT_FINAL_ASSEMBLY));
227f3237bfeSMark Adams           PetscCall(MatAssemblyEnd(D, MAT_FINAL_ASSEMBLY));
228f3237bfeSMark Adams           PetscCall(PetscInfo(M_p, "createMtMKSP Have %" PetscInt_FMT " eqs, nzl = %" PetscInt_FMT "\n", N, nzl));
229f3237bfeSMark Adams           PetscCall(KSPSetOperators(ksp, MtM, D));
230f3237bfeSMark Adams           PetscCall(MatViewFromOptions(D, NULL, "-ftop2_D_mat_view"));
231f3237bfeSMark Adams           PetscCall(MatViewFromOptions(M_p, NULL, "-ftop2_Mp_mat_view"));
232f3237bfeSMark Adams           PetscCall(MatViewFromOptions(matshellctx->MpTrans, NULL, "-ftop2_MpTranspose_mat_view"));
2336b664d00Smarkadams4           PetscCall(MatViewFromOptions(MtM, NULL, "-ftop2_MtM_mat_view"));
234f3237bfeSMark Adams         }
235f3237bfeSMark Adams       }
236d043ef4cSMark Adams     }
237d043ef4cSMark Adams   }
238f3237bfeSMark Adams   if (is_lsqr) {
239d043ef4cSMark Adams     PC        pc2;
240f3237bfeSMark Adams     PetscBool is_bjac;
241d043ef4cSMark Adams     PetscCall(KSPGetPC(ksp, &pc2));
242d043ef4cSMark Adams     PetscCall(PetscObjectTypeCompare((PetscObject)pc2, PCBJACOBI, &is_bjac));
243f3237bfeSMark Adams     if (is_bjac) {
244f3237bfeSMark Adams       PetscCall(DMSwarmCreateMassMatrixSquare(sw, dm, &PM_p));
245f3237bfeSMark Adams       PetscCall(KSPSetOperators(ksp, M_p, PM_p));
246f3237bfeSMark Adams     } else {
247f3237bfeSMark Adams       PetscCall(KSPSetOperators(ksp, M_p, M_p));
248f3237bfeSMark Adams     }
249d043ef4cSMark Adams     PetscCall(MatViewFromOptions(M_p, NULL, "-ftop2_Mp_mat_view"));
250f3237bfeSMark Adams   }
251f3237bfeSMark Adams   PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &ff)); // this grabs access
252f3237bfeSMark Adams   if (!is_lsqr) {
253d043ef4cSMark Adams     PetscErrorCode ierr;
254d043ef4cSMark Adams     ierr = KSPSolve(ksp, work_ferhs, matshellctx->uu);
255d043ef4cSMark Adams     if (!ierr) {
2563643261fSMark Adams       // 3) with Moore-Penrose apply Mp: M_p (Mp' Mp)^-1 M
257f3237bfeSMark Adams       PetscCall(MatMult(M_p, matshellctx->uu, ff));
258d043ef4cSMark Adams     } else { // failed
259d043ef4cSMark Adams       PC        pc2;
260d043ef4cSMark Adams       PetscBool is_bjac;
261d043ef4cSMark Adams       PetscCall(PetscInfo(ksp, "Solver failed, probably singular, try lsqr\n"));
262d043ef4cSMark Adams       PetscCall(KSPReset(ksp));
263d043ef4cSMark Adams       PetscCall(KSPSetType(ksp, KSPLSQR));
264d043ef4cSMark Adams       PetscCall(KSPGetPC(ksp, &pc2));
265d043ef4cSMark Adams       PetscCall(PCSetType(pc2, PCNONE)); // should not happen, but could solve stable (Mp Mp^T), move projection Mp before solve
266d043ef4cSMark Adams       PetscCall(KSPSetOptionsPrefix(ksp, "ftop_"));
267d043ef4cSMark Adams       PetscCall(KSPSetFromOptions(ksp));
268d043ef4cSMark Adams       PetscCall(PetscObjectTypeCompare((PetscObject)pc2, PCBJACOBI, &is_bjac));
269d043ef4cSMark Adams       if (is_bjac) {
270d043ef4cSMark Adams         PetscCall(DMSwarmCreateMassMatrixSquare(sw, dm, &PM_p));
271d043ef4cSMark Adams         PetscCall(KSPSetOperators(ksp, M_p, PM_p));
272d043ef4cSMark Adams       } else {
273d043ef4cSMark Adams         PetscCall(KSPSetOperators(ksp, M_p, M_p));
274d043ef4cSMark Adams       }
275d043ef4cSMark Adams       ierr = KSPSolveTranspose(ksp, work_ferhs, ff);
276d043ef4cSMark Adams       if (ierr) { PetscCheck(!ierr, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "backup LSQR solver failed - need to add N_v > N_p Moore-Penrose pseudo-inverse"); }
277d043ef4cSMark Adams     }
278d043ef4cSMark Adams     if (D) PetscCall(MatDestroy(&D));
279d043ef4cSMark Adams     PetscCall(MatDestroy(&MtM));
280d043ef4cSMark Adams     if (matshellctx->MpTrans) PetscCall(MatDestroy(&matshellctx->MpTrans));
281f3237bfeSMark Adams     PetscCall(VecDestroy(&matshellctx->ff));
282f3237bfeSMark Adams     PetscCall(VecDestroy(&matshellctx->uu));
283f3237bfeSMark Adams     PetscCall(PetscFree(matshellctx));
284f3237bfeSMark Adams   } else {
285d043ef4cSMark Adams     PetscErrorCode ierr;
2863643261fSMark Adams     // finally with LSQR apply M_p^\dagger
287d043ef4cSMark Adams     ierr = KSPSolveTranspose(ksp, work_ferhs, ff);
288d043ef4cSMark Adams     if (ierr) { PetscCheck(!ierr, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "backup LSQR solver failed - need to add N_v > N_p Moore-Penrose pseudo-inverse"); }
289f3237bfeSMark Adams   }
290f3237bfeSMark Adams   PetscCall(KSPDestroy(&ksp));
291f3237bfeSMark Adams   PetscCall(MatDestroy(&PM_p));
292f3237bfeSMark Adams   PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &ff));
2933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
294f3237bfeSMark Adams }
295f3237bfeSMark Adams 
2963cdbf90fSMark Adams #define EX30_MAX_NUM_THRDS 12
2973cdbf90fSMark Adams #define EX30_MAX_BATCH_SZ  1024
2983cdbf90fSMark Adams //
2993cdbf90fSMark Adams // add grid to arg 'globSwarmArray[].w_q'
3003cdbf90fSMark Adams //
3013cdbf90fSMark 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)
302d71ae5a4SJacob Faibussowitsch {
3033cdbf90fSMark Adams   PetscErrorCode ierr = (PetscErrorCode)0; // used for inside thread loops
304f3237bfeSMark Adams 
305f3237bfeSMark Adams   PetscFunctionBeginUser;
3063cdbf90fSMark Adams   // map back to particles
3073cdbf90fSMark Adams   for (PetscInt v_id_0 = 0; v_id_0 < ctx->batch_sz; v_id_0 += numthreads) {
3083cdbf90fSMark 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));
3093cdbf90fSMark Adams     //PetscPragmaOMP(parallel for)
3106497c311SBarry Smith     for (PetscInt tid = 0; tid < numthreads; tid++) {
3113cdbf90fSMark Adams       const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id + v_id;
3123cdbf90fSMark Adams       if (glb_v_id < num_vertices) {
3133cdbf90fSMark Adams         for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids
3143cdbf90fSMark Adams           PetscErrorCode ierr_t;
3153cdbf90fSMark 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));
3163cdbf90fSMark 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]);
3173cdbf90fSMark Adams           if (ierr_t) ierr = ierr_t;
318f3237bfeSMark Adams         }
3193cdbf90fSMark Adams       }
3203cdbf90fSMark Adams     }
321d043ef4cSMark Adams     PetscCheck(!ierr, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Error in OMP loop. ierr = %d", (int)ierr);
3223cdbf90fSMark Adams     /* Get moments */
3233cdbf90fSMark Adams     PetscCall(PetscInfo(grid_dm[0], "Cleanup batches %" PetscInt_FMT " to %" PetscInt_FMT "\n", v_id_0, v_id_0 + numthreads));
3246497c311SBarry Smith     for (PetscInt tid = 0; tid < numthreads; tid++) {
3253cdbf90fSMark Adams       const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id + v_id;
3263cdbf90fSMark Adams       if (glb_v_id == v_target) {
3273cdbf90fSMark Adams         for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
3283cdbf90fSMark Adams           PetscDataType dtype;
3293cdbf90fSMark Adams           PetscReal    *wp, *coords;
3303cdbf90fSMark Adams           DM            sw = globSwarmArray[LAND_PACK_IDX(v_id, grid)];
3313cdbf90fSMark Adams           PetscInt      npoints, bs = 1;
3323cdbf90fSMark Adams           PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void **)&wp)); // take data out here
3333cdbf90fSMark Adams           PetscCall(DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
3343cdbf90fSMark Adams           PetscCall(DMSwarmGetLocalSize(sw, &npoints));
3356497c311SBarry Smith           for (PetscInt p = 0; p < npoints; p++) {
3363cdbf90fSMark 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]];
3376497c311SBarry Smith             for (PetscInt i = 0; i < dim; ++i) v2 += PetscSqr(coords[p * dim + i]);
3383cdbf90fSMark Adams             moments[0] += w;
3393cdbf90fSMark Adams             moments[1] += w * ctx->v_0 * coords[p * dim + 1]; // z-momentum
3401b431c67SMark Adams             moments[2] += w * 0.5 * ctx->v_0 * ctx->v_0 * v2;
3413cdbf90fSMark Adams           }
3423cdbf90fSMark Adams           PetscCall(DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void **)&wp));
343f3237bfeSMark Adams           PetscCall(DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
3443cdbf90fSMark Adams         }
3453cdbf90fSMark Adams         const PetscReal N_inv = 1 / moments[0];
3461b431c67SMark Adams         PetscCall(PetscInfo(grid_dm[0], "gridToParticles_private [%" PetscInt_FMT "], n = %g\n", v_id, (double)moments[0]));
3473cdbf90fSMark Adams         for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
3483cdbf90fSMark Adams           PetscDataType dtype;
3493cdbf90fSMark Adams           PetscReal    *wp, *coords;
3503cdbf90fSMark Adams           DM            sw = globSwarmArray[LAND_PACK_IDX(v_id, grid)];
3513cdbf90fSMark Adams           PetscInt      npoints, bs = 1;
3523cdbf90fSMark Adams           PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void **)&wp)); // take data out here
3533cdbf90fSMark Adams           PetscCall(DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
3543cdbf90fSMark Adams           PetscCall(DMSwarmGetLocalSize(sw, &npoints));
3556497c311SBarry Smith           for (PetscInt p = 0; p < npoints; p++) {
3563cdbf90fSMark 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;
3571b431c67SMark Adams             if (w > PETSC_REAL_MIN) {
3583cdbf90fSMark Adams               moments[3] -= ww * PetscLogReal(ww);
359d043ef4cSMark Adams               PetscCheck(ww < 1 - PETSC_MACHINE_EPSILON, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "ww (%g) > 1", (double)ww);
360d043ef4cSMark Adams             } else moments[4] -= w; // keep track of density that is lost
3613cdbf90fSMark Adams           }
3623cdbf90fSMark Adams           PetscCall(DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void **)&wp));
3633cdbf90fSMark Adams           PetscCall(DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
3643cdbf90fSMark Adams         }
3653cdbf90fSMark Adams       }
3663cdbf90fSMark Adams     } // thread batch
3673cdbf90fSMark Adams   } // batch
3683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
369f3237bfeSMark Adams }
3703cdbf90fSMark Adams 
3713cdbf90fSMark Adams static void maxwellian(PetscInt dim, const PetscReal x[], PetscReal kt_m, PetscReal n, PetscReal shift, PetscScalar *u)
372d71ae5a4SJacob Faibussowitsch {
373f3237bfeSMark Adams   PetscInt  i;
374f3237bfeSMark Adams   PetscReal v2 = 0, theta = 2.0 * kt_m; /* theta = 2kT/mc^2 */
375f3237bfeSMark Adams 
3763cdbf90fSMark Adams   if (shift != 0.) {
3773cdbf90fSMark Adams     v2 = 0;
3783cdbf90fSMark Adams     for (i = 0; i < dim - 1; ++i) v2 += x[i] * x[i];
3793cdbf90fSMark Adams     v2 += (x[dim - 1] - shift) * (x[dim - 1] - shift);
3803cdbf90fSMark Adams     /* evaluate the shifted Maxwellian */
3813cdbf90fSMark Adams     u[0] += n * PetscPowReal(PETSC_PI * theta, -1.5) * (PetscExpReal(-v2 / theta));
3821b431c67SMark Adams   } else {
3831b431c67SMark Adams     /* compute the exponents, v^2 */
3841b431c67SMark Adams     for (i = 0; i < dim; ++i) v2 += x[i] * x[i];
3851b431c67SMark Adams     /* evaluate the Maxwellian */
3861b431c67SMark Adams     u[0] += n * PetscPowReal(PETSC_PI * theta, -1.5) * (PetscExpReal(-v2 / theta));
3873cdbf90fSMark Adams   }
388f3237bfeSMark Adams }
389f3237bfeSMark Adams 
3903cdbf90fSMark Adams static PetscErrorCode PostStep(TS ts)
3913cdbf90fSMark Adams {
3921b431c67SMark Adams   PetscInt   n, dim, nDMs, v_id;
3933cdbf90fSMark Adams   PetscReal  t;
3943cdbf90fSMark Adams   LandauCtx *ctx;
3953cdbf90fSMark Adams   Vec        X;
3963cdbf90fSMark Adams   PrintCtx  *printCtx;
3973cdbf90fSMark Adams   DM         pack;
3981b431c67SMark Adams   PetscReal  moments[5], e_grid[LANDAU_MAX_GRIDS];
3993cdbf90fSMark Adams 
4003cdbf90fSMark Adams   PetscFunctionBeginUser;
4013cdbf90fSMark Adams   PetscCall(TSGetApplicationContext(ts, &printCtx));
4021b431c67SMark Adams   if (!printCtx->print && !printCtx->print_entropy) PetscFunctionReturn(PETSC_SUCCESS);
4033cdbf90fSMark Adams   ctx = printCtx->ctx;
4041b431c67SMark Adams   if (printCtx->v_target < printCtx->global_vertex_id_0 || printCtx->v_target >= printCtx->global_vertex_id_0 + ctx->batch_sz) PetscFunctionReturn(PETSC_SUCCESS);
4056497c311SBarry Smith   for (PetscInt i = 0; i < 5; i++) moments[i] = 0;
4066497c311SBarry Smith   for (PetscInt i = 0; i < LANDAU_MAX_GRIDS; i++) e_grid[i] = 0;
4071b431c67SMark Adams   v_id = printCtx->v_target % ctx->batch_sz;
4083cdbf90fSMark Adams   PetscCall(TSGetDM(ts, &pack));
4093cdbf90fSMark Adams   PetscCall(DMGetDimension(pack, &dim));
4103cdbf90fSMark Adams   PetscCall(DMCompositeGetNumberDM(pack, &nDMs)); // number of vertices * number of grids
4113cdbf90fSMark Adams   PetscCall(TSGetSolution(ts, &X));
4123cdbf90fSMark Adams   PetscCall(TSGetStepNumber(ts, &n));
4133cdbf90fSMark Adams   PetscCall(TSGetTime(ts, &t));
4143cdbf90fSMark Adams   PetscCall(DMCompositeGetAccessArray(pack, X, nDMs, NULL, printCtx->globXArray));
415d043ef4cSMark Adams   if (printCtx->print_entropy && printCtx->v_target >= 0 && 0) {
4163cdbf90fSMark Adams     for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
4173cdbf90fSMark Adams       PetscDataType dtype;
4183cdbf90fSMark Adams       PetscReal    *wp, *coords;
4193cdbf90fSMark Adams       DM            sw = printCtx->globSwarmArray[LAND_PACK_IDX(v_id, grid)];
4203cdbf90fSMark Adams       Vec           work, subX = printCtx->globXArray[LAND_PACK_IDX(v_id, grid)];
4211b431c67SMark Adams       PetscInt      bs, NN;
4223cdbf90fSMark Adams       // C-G moments
4233cdbf90fSMark Adams       PetscCall(VecDuplicate(subX, &work));
4243cdbf90fSMark Adams       PetscCall(gridToParticles(printCtx->grid_dm[grid], sw, subX, work, printCtx->globMpArray[LAND_PACK_IDX(v_id, grid)], printCtx->g_Mass[grid]));
4253cdbf90fSMark Adams       PetscCall(VecDestroy(&work));
4263cdbf90fSMark Adams       // moments
4273cdbf90fSMark Adams       PetscCall(DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
4281b431c67SMark Adams       PetscCall(DMSwarmGetLocalSize(sw, &NN));
4291b431c67SMark Adams       PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void **)&wp));
4306497c311SBarry Smith       for (PetscInt pp = 0; pp < NN; pp++) {
4311b431c67SMark 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]];
4326497c311SBarry Smith         for (PetscInt i = 0; i < dim; ++i) v2 += PetscSqr(coords[pp * dim + i]);
4333cdbf90fSMark Adams         moments[0] += w;
4343cdbf90fSMark Adams         moments[1] += w * ctx->v_0 * coords[pp * dim + 1]; // z-momentum
4351b431c67SMark Adams         moments[2] += w * 0.5 * ctx->v_0 * ctx->v_0 * v2;
4361b431c67SMark Adams         e_grid[grid] += w * 0.5 * ctx->v_0 * ctx->v_0 * v2;
4373cdbf90fSMark Adams       }
4383cdbf90fSMark Adams       PetscCall(DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
4393cdbf90fSMark Adams       PetscCall(DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void **)&wp));
4403cdbf90fSMark Adams     }
4411b431c67SMark Adams     // entropy
4421b431c67SMark Adams     const PetscReal N_inv = 1 / moments[0];
4431b431c67SMark Adams     for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
4441b431c67SMark Adams       PetscDataType dtype;
4451b431c67SMark Adams       PetscReal    *wp, *coords;
4461b431c67SMark Adams       DM            sw = printCtx->globSwarmArray[LAND_PACK_IDX(v_id, grid)];
4471b431c67SMark Adams       PetscInt      bs, NN;
4481b431c67SMark Adams       PetscCall(DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
4491b431c67SMark Adams       PetscCall(DMSwarmGetLocalSize(sw, &NN));
4501b431c67SMark Adams       PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void **)&wp));
4516497c311SBarry Smith       for (PetscInt pp = 0; pp < NN; pp++) {
4521b431c67SMark 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;
4531b431c67SMark Adams         if (w > PETSC_REAL_MIN) {
4541b431c67SMark Adams           moments[3] -= ww * PetscLogReal(ww);
4551b431c67SMark Adams         } else moments[4] -= w;
4561b431c67SMark Adams       }
4571b431c67SMark Adams       PetscCall(DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
4581b431c67SMark Adams       PetscCall(DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void **)&wp));
4591b431c67SMark Adams     }
460d043ef4cSMark 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] / moments[0]), (double)e_grid[0], (double)e_grid[1], (double)e_grid[2]));
4611b431c67SMark Adams   }
4621b431c67SMark Adams   if (printCtx->print && printCtx->g_target >= 0) {
4631b431c67SMark Adams     PetscInt         grid   = printCtx->g_target, id;
4641b431c67SMark Adams     static PetscReal last_t = -100000, period = .5;
4651b431c67SMark Adams     if (last_t == -100000) last_t = -period + t;
4661b431c67SMark Adams     if (t >= last_t + period) {
4671b431c67SMark Adams       last_t = t;
4681b431c67SMark Adams       PetscCall(DMGetOutputSequenceNumber(ctx->plex[grid], &id, NULL));
4691b431c67SMark Adams       PetscCall(DMSetOutputSequenceNumber(ctx->plex[grid], id + 1, t));
4701b431c67SMark Adams       PetscCall(VecViewFromOptions(printCtx->globXArray[LAND_PACK_IDX(v_id % ctx->batch_sz, grid)], NULL, "-ex30_vec_view"));
4711b431c67SMark Adams       if (ctx->num_grids > grid + 1) {
4721b431c67SMark Adams         PetscCall(DMSetOutputSequenceNumber(ctx->plex[grid + 1], id + 1, t));
4731b431c67SMark Adams         PetscCall(VecViewFromOptions(printCtx->globXArray[LAND_PACK_IDX(v_id % ctx->batch_sz, grid + 1)], NULL, "-ex30_vec_view2"));
4741b431c67SMark Adams       }
4751b431c67SMark Adams       PetscCall(PetscInfo(X, "%4d) time %e View\n", (int)n, (double)t));
4761b431c67SMark Adams     }
4771b431c67SMark Adams   }
4783cdbf90fSMark Adams   PetscCall(DMCompositeRestoreAccessArray(pack, X, nDMs, NULL, printCtx->globXArray));
4793cdbf90fSMark Adams   PetscFunctionReturn(PETSC_SUCCESS);
4803cdbf90fSMark Adams }
4813cdbf90fSMark Adams 
4821b431c67SMark 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)
483d71ae5a4SJacob Faibussowitsch {
484f3237bfeSMark Adams   DM             pack, *globSwarmArray, grid_dm[LANDAU_MAX_GRIDS];
485f3237bfeSMark Adams   Mat           *globMpArray, g_Mass[LANDAU_MAX_GRIDS];
4863cdbf90fSMark Adams   KSP            t_ksp[LANDAU_MAX_GRIDS][EX30_MAX_NUM_THRDS];
4873cdbf90fSMark Adams   Vec            t_fhat[LANDAU_MAX_GRIDS][EX30_MAX_NUM_THRDS];
4881b431c67SMark Adams   PetscInt       nDMs;
4893cdbf90fSMark Adams   PetscErrorCode ierr = (PetscErrorCode)0; // used for inside thread loops
490f3237bfeSMark Adams #if defined(PETSC_HAVE_OPENMP) && defined(PETSC_HAVE_THREADSAFETY)
491f3237bfeSMark Adams   PetscInt numthreads = PetscNumOMPThreads;
492f3237bfeSMark Adams #else
493f3237bfeSMark Adams   PetscInt numthreads = 1;
494f3237bfeSMark Adams #endif
495f3237bfeSMark Adams   LandauCtx *ctx;
496f3237bfeSMark Adams   Vec       *globXArray;
4971b431c67SMark Adams   PetscReal  moments_0[5], moments_1a[5], moments_1b[5], dt_init;
4983cdbf90fSMark Adams   PrintCtx  *printCtx;
499f3237bfeSMark Adams 
500f3237bfeSMark Adams   PetscFunctionBeginUser;
5013cdbf90fSMark 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);
5023cdbf90fSMark Adams   PetscCheck(numthreads > 0, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Number threads %" PetscInt_FMT " > %d", numthreads, EX30_MAX_NUM_THRDS);
503f3237bfeSMark Adams   PetscCall(TSGetDM(ts, &pack));
504f3237bfeSMark Adams   PetscCall(DMGetApplicationContext(pack, &ctx));
505f3237bfeSMark 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);
5063cdbf90fSMark Adams   PetscCall(DMCompositeGetNumberDM(pack, &nDMs)); // number of vertices * number of grids
5071b431c67SMark 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));
508f3237bfeSMark Adams   PetscCall(PetscMalloc(sizeof(*globXArray) * nDMs, &globXArray));
509f3237bfeSMark Adams   PetscCall(PetscMalloc(sizeof(*globMpArray) * nDMs, &globMpArray));
510f3237bfeSMark Adams   PetscCall(PetscMalloc(sizeof(*globSwarmArray) * nDMs, &globSwarmArray));
5113cdbf90fSMark Adams   // print ctx
5123cdbf90fSMark Adams   PetscCall(PetscNew(&printCtx));
5133cdbf90fSMark Adams   PetscCall(TSSetApplicationContext(ts, printCtx));
5143cdbf90fSMark Adams   printCtx->v_target       = v_target;
5151b431c67SMark Adams   printCtx->g_target       = g_target;
5163cdbf90fSMark Adams   printCtx->ctx            = ctx;
5173cdbf90fSMark Adams   printCtx->globSwarmArray = globSwarmArray;
5183cdbf90fSMark Adams   printCtx->grid_dm        = grid_dm;
5193cdbf90fSMark Adams   printCtx->globMpArray    = globMpArray;
5203cdbf90fSMark Adams   printCtx->g_Mass         = g_Mass;
5213cdbf90fSMark Adams   printCtx->globXArray     = globXArray;
5221b431c67SMark Adams   printCtx->print_entropy  = PETSC_FALSE;
52315229ffcSPierre Jolivet   PetscOptionsBegin(PETSC_COMM_SELF, "", "Print Options", "DMPLEX");
5241b431c67SMark Adams   PetscCall(PetscOptionsBool("-print_entropy", "Print entropy and moments at each time step", "ex30.c", printCtx->print_entropy, &printCtx->print_entropy, NULL));
5251b431c67SMark Adams   PetscOptionsEnd();
5263cdbf90fSMark Adams   // view
527f3237bfeSMark Adams   PetscCall(DMViewFromOptions(ctx->plex[g_target], NULL, "-ex30_dm_view"));
5281b431c67SMark Adams   if (ctx->num_grids > g_target + 1) { PetscCall(DMViewFromOptions(ctx->plex[g_target + 1], NULL, "-ex30_dm_view2")); }
5293cdbf90fSMark Adams   // create mesh mass matrices
5306b664d00Smarkadams4   PetscCall(VecZeroEntries(X));
531f3237bfeSMark Adams   PetscCall(DMCompositeGetAccessArray(pack, X, nDMs, NULL, globXArray)); // just to duplicate
532f3237bfeSMark Adams   for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {               // add same particels for all grids
533f3237bfeSMark Adams     Vec          subX = globXArray[LAND_PACK_IDX(0, grid)];
534f3237bfeSMark Adams     DM           dm   = ctx->plex[grid];
535f3237bfeSMark Adams     PetscSection s;
536f3237bfeSMark Adams     grid_dm[grid] = dm;
537f3237bfeSMark Adams     PetscCall(DMCreateMassMatrix(dm, dm, &g_Mass[grid]));
538f3237bfeSMark Adams     //
539f3237bfeSMark Adams     PetscCall(DMGetLocalSection(dm, &s));
540f3237bfeSMark Adams     PetscCall(DMPlexCreateClosureIndex(dm, s));
5416497c311SBarry Smith     for (PetscInt tid = 0; tid < numthreads; tid++) {
5421b431c67SMark Adams       PC pc;
543f3237bfeSMark Adams       PetscCall(VecDuplicate(subX, &t_fhat[grid][tid]));
544f3237bfeSMark Adams       PetscCall(KSPCreate(PETSC_COMM_SELF, &t_ksp[grid][tid]));
5451b431c67SMark Adams       PetscCall(KSPSetType(t_ksp[grid][tid], KSPCG));
5461b431c67SMark Adams       PetscCall(KSPGetPC(t_ksp[grid][tid], &pc));
5471b431c67SMark Adams       PetscCall(PCSetType(pc, PCJACOBI));
548f3237bfeSMark Adams       PetscCall(KSPSetOptionsPrefix(t_ksp[grid][tid], "ptof_"));
549f3237bfeSMark Adams       PetscCall(KSPSetOperators(t_ksp[grid][tid], g_Mass[grid], g_Mass[grid]));
550f3237bfeSMark Adams       PetscCall(KSPSetFromOptions(t_ksp[grid][tid]));
551f3237bfeSMark Adams     }
552f3237bfeSMark Adams   }
553f3237bfeSMark Adams   PetscCall(DMCompositeRestoreAccessArray(pack, X, nDMs, NULL, globXArray));
554f3237bfeSMark Adams   PetscCall(TSGetTimeStep(ts, &dt_init)); // we could have an adaptive time stepper
5553cdbf90fSMark Adams   // loop over all vertices in chucks that are batched for TSSolve
5566497c311SBarry Smith   for (PetscInt i = 0; i < 5; i++) moments_0[i] = moments_1a[i] = moments_1b[i] = 0;
5573cdbf90fSMark 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
5586b664d00Smarkadams4     PetscCall(TSSetTime(ts, 0));
5596b664d00Smarkadams4     PetscCall(TSSetStepNumber(ts, 0));
5606b664d00Smarkadams4     PetscCall(TSSetTimeStep(ts, dt_init));
561f3237bfeSMark Adams     PetscCall(DMCompositeGetAccessArray(pack, X, nDMs, NULL, globXArray));
5621b431c67SMark Adams     printCtx->global_vertex_id_0 = global_vertex_id_0;
5633cdbf90fSMark Adams     if (v_target >= global_vertex_id_0 && v_target < global_vertex_id_0 + ctx->batch_sz) {
5643cdbf90fSMark Adams       PetscCall(PetscObjectSetName((PetscObject)globXArray[LAND_PACK_IDX(v_target % ctx->batch_sz, g_target)], "rho"));
5653cdbf90fSMark Adams       printCtx->print = PETSC_TRUE;
5663cdbf90fSMark Adams     } else printCtx->print = PETSC_FALSE;
5673cdbf90fSMark Adams     // create fake particles in batches with threads
5683cdbf90fSMark Adams     for (PetscInt v_id_0 = 0; v_id_0 < ctx->batch_sz; v_id_0 += numthreads) {
569d043ef4cSMark 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] /* , radiuses[80000] */;
5703cdbf90fSMark Adams       PetscInt   Np_t[LANDAU_MAX_GRIDS][EX30_MAX_NUM_THRDS];
571f3237bfeSMark Adams       // make particles
5726497c311SBarry Smith       for (PetscInt tid = 0; tid < numthreads; tid++) {
5733cdbf90fSMark Adams         const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id;
5743cdbf90fSMark Adams         if (glb_v_id < num_vertices) {                                                     // the ragged edge (in last batch)
5751b431c67SMark 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
576f3237bfeSMark Adams           for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {                         // add same particels for all grids
577d043ef4cSMark Adams             // for (PetscInt sp = ctx->species_offset[grid], i0 = 0; sp < ctx->species_offset[grid + 1]; sp++, i0++) {
5783cdbf90fSMark 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 */
579f3237bfeSMark 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
580f3237bfeSMark Adams             PetscInt        Npi = Npp0, Npj = 2 * Npp0, Npk = 1;
5811b431c67SMark Adams             PetscRandom     rand;
582d043ef4cSMark Adams             PetscReal       sigma = ctx->thermal_speed[grid] / ctx->thermal_speed[0], p2_shift = grid == 0 ? shift : -shift; // symmetric shift of e vs ions
5831b431c67SMark Adams             PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rand));
5841b431c67SMark Adams             PetscCall(PetscRandomSetInterval(rand, 0., 1.));
5851b431c67SMark Adams             PetscCall(PetscRandomSetFromOptions(rand));
586f3237bfeSMark Adams             if (dim == 2) lo[0] = 0; // Landau coordinate (r,z)
587f3237bfeSMark Adams             else Npi = Npj = Npk = Npp0;
5883cdbf90fSMark Adams             // User: use glb_v_id to index into your data
589d043ef4cSMark Adams             const PetscInt NNreal = Npi * Npj * Npk, NN = NNreal + (dim == 2 ? 3 : 6); // make room for bounding box
590f3237bfeSMark Adams             Np_t[grid][tid] = NN;
5913cdbf90fSMark Adams             if (glb_v_id == v_target) nTargetP[grid] = NN;
592f3237bfeSMark 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]));
593f3237bfeSMark Adams             hp[0] = (hi[0] - lo[0]) / Npi;
594f3237bfeSMark Adams             hp[1] = (hi[1] - lo[1]) / Npj;
595f3237bfeSMark Adams             hp[2] = (hi[2] - lo[2]) / Npk;
596f3237bfeSMark Adams             if (dim == 2) hp[2] = 1;
59763a3b9bcSJacob 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
598f3237bfeSMark Adams             vole = hp[0] * hp[1] * hp[2] * ctx->n[grid];                                                                                                                           // fix for multi-species
5993cdbf90fSMark 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));
6006497c311SBarry Smith             for (PetscInt pj = 0, pp = 0; pj < Npj; pj++) {
6016497c311SBarry Smith               for (PetscInt pk = 0; pk < Npk; pk++) {
6026497c311SBarry Smith                 for (PetscInt pi = 0; pi < Npi; pi++, pp++) {
6031b431c67SMark Adams                   PetscReal p_shift   = p2_shift;
6041b431c67SMark Adams                   wp_t[grid][tid][pp] = 0;
6051b431c67SMark Adams                   if (use_uniform_particle_grid) {
606f3237bfeSMark Adams                     xx_t[grid][tid][pp] = lo[0] + hp[0] / 2.0 + pi * hp[0];
607f3237bfeSMark Adams                     yy_t[grid][tid][pp] = lo[1] + hp[1] / 2.0 + pj * hp[1];
608f3237bfeSMark Adams                     if (dim == 3) zz_t[grid][tid][pp] = lo[2] + hp[2] / 2.0 + pk * hp[2];
609f3237bfeSMark Adams                     PetscReal x[] = {xx_t[grid][tid][pp], yy_t[grid][tid][pp], dim == 2 ? 0 : zz_t[grid][tid][pp]};
6101b431c67SMark Adams                     p_shift *= ctx->thermal_speed[grid] / ctx->v_0;
611d043ef4cSMark Adams                     if (ctx->sphere && PetscSqrtReal(PetscSqr(xx_t[grid][tid][pp]) + PetscSqr(yy_t[grid][tid][pp])) > 0.92 * hi[0]) {
612d043ef4cSMark Adams                       wp_t[grid][tid][pp] = 0;
613d043ef4cSMark Adams                     } else {
6141b431c67SMark Adams                       maxwellian(dim, x, kT_m, vole, p_shift, &wp_t[grid][tid][pp]);
6151b431c67SMark Adams                       if (ctx->num_grids == 1 && shift != 0) {                          // bi-maxwellian, electron plasma
616d043ef4cSMark Adams                         maxwellian(dim, x, kT_m, vole, -p_shift, &wp_t[grid][tid][pp]); // symmetric shift of electron plasma
617d043ef4cSMark Adams                       }
6181b431c67SMark Adams                     }
6191b431c67SMark Adams                   } else {
6201b431c67SMark Adams                     PetscReal u1, u2;
6211b431c67SMark Adams                     do {
6221b431c67SMark Adams                       do {
6231b431c67SMark Adams                         PetscCall(PetscRandomGetValueReal(rand, &u1));
6241b431c67SMark Adams                       } while (u1 == 0);
6251b431c67SMark Adams                       PetscCall(PetscRandomGetValueReal(rand, &u2));
6261b431c67SMark Adams                       //compute z0 and z1
627d043ef4cSMark Adams                       PetscReal mag       = sigma * PetscSqrtReal(-2.0 * PetscLogReal(u1)); // is this the same scale grid Maxwellian? t_therm = sigma
6281b431c67SMark Adams                       xx_t[grid][tid][pp] = mag * PetscCosReal(2.0 * PETSC_PI * u2);
6291b431c67SMark Adams                       yy_t[grid][tid][pp] = mag * PetscSinReal(2.0 * PETSC_PI * u2);
6301b431c67SMark Adams                       if (dim == 2 && xx_t[grid][tid][pp] < lo[0]) xx_t[grid][tid][pp] = -xx_t[grid][tid][pp];
6311b431c67SMark Adams                       if (dim == 3) zz_t[grid][tid][pp] = lo[2] + hp[2] / 2.0 + pk * hp[2];
632d043ef4cSMark Adams                       if (!ctx->sphere) {
633d043ef4cSMark Adams                         if (dim == 2 && xx_t[grid][tid][pp] < 0) xx_t[grid][tid][pp] = -xx_t[grid][tid][pp]; // ???
634d043ef4cSMark Adams                         else if (dim == 3) {
635d043ef4cSMark Adams                           while (zz_t[grid][tid][pp] >= hi[2] || zz_t[grid][tid][pp] <= lo[2]) zz_t[grid][tid][pp] *= .9;
636d043ef4cSMark Adams                         }
637d043ef4cSMark Adams                         while (xx_t[grid][tid][pp] >= hi[0] || xx_t[grid][tid][pp] <= lo[0]) xx_t[grid][tid][pp] *= .9;
638d043ef4cSMark Adams                         while (yy_t[grid][tid][pp] >= hi[1] || yy_t[grid][tid][pp] <= lo[1]) yy_t[grid][tid][pp] *= .9;
639d043ef4cSMark Adams                       } else { // 2D
640d043ef4cSMark Adams                         //if (glb_v_id == v_target && pp < 80000) radiuses[pp] = PetscSqrtReal(PetscSqr(xx_t[grid][tid][pp]) + PetscSqr(yy_t[grid][tid][pp]));
641d043ef4cSMark Adams                         while (PetscSqrtReal(PetscSqr(xx_t[grid][tid][pp]) + PetscSqr(yy_t[grid][tid][pp])) > 0.92 * hi[0]) { // safety factor for facets of sphere
642d043ef4cSMark Adams                           xx_t[grid][tid][pp] *= .9;
643d043ef4cSMark Adams                           yy_t[grid][tid][pp] *= .9;
644d043ef4cSMark Adams                         }
645d043ef4cSMark Adams                       }
6461b431c67SMark Adams                       if (ctx->num_grids == 1 && pp % 2 == 0) p_shift = 0; // one species, split bi-max
6471b431c67SMark Adams                       p_shift *= ctx->thermal_speed[grid] / ctx->v_0;
6481b431c67SMark Adams                       if (dim == 3) zz_t[grid][tid][pp] += p_shift;
6491b431c67SMark Adams                       else yy_t[grid][tid][pp] += p_shift;
650d043ef4cSMark Adams                       wp_t[grid][tid][pp] += ctx->n[grid] / NNreal * PetscSqrtReal(ctx->masses[ctx->species_offset[grid]] / ctx->masses[0]);
6511b431c67SMark Adams                       if (p_shift <= 0) break; // add bi-max for electron plasma only
6521b431c67SMark Adams                       p_shift = -p_shift;
6531b431c67SMark Adams                     } while (ctx->num_grids == 1); // add bi-max for electron plasma only
6541b431c67SMark Adams                   }
6551b431c67SMark Adams                   {
6563cdbf90fSMark Adams                     if (glb_v_id == v_target) {
6571b431c67SMark Adams                       PetscReal x[] = {xx_t[grid][tid][pp], yy_t[grid][tid][pp], dim == 2 ? 0 : zz_t[grid][tid][pp]};
6583cdbf90fSMark 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]];
6596497c311SBarry Smith                       for (PetscInt i = 0; i < dim; ++i) v2 += PetscSqr(x[i]);
6603cdbf90fSMark Adams                       moments_0[0] += w;                   // not thread safe
6613cdbf90fSMark Adams                       moments_0[1] += w * ctx->v_0 * x[1]; // z-momentum
6621b431c67SMark Adams                       moments_0[2] += w * 0.5 * ctx->v_0 * ctx->v_0 * v2;
663f3237bfeSMark Adams                     }
664f3237bfeSMark Adams                   }
665f3237bfeSMark Adams                 }
666f3237bfeSMark Adams               }
667f3237bfeSMark Adams             }
668d043ef4cSMark Adams             if (dim == 2) { // fix bounding box
6696497c311SBarry Smith               PetscInt pp           = NNreal;
670d043ef4cSMark Adams               wp_t[grid][tid][pp]   = 0;
671d043ef4cSMark Adams               xx_t[grid][tid][pp]   = 1.e-7;
672d043ef4cSMark Adams               yy_t[grid][tid][pp++] = hi[1] - 5.e-7;
673d043ef4cSMark Adams               wp_t[grid][tid][pp]   = 0;
674d043ef4cSMark Adams               xx_t[grid][tid][pp]   = hi[0] - 5.e-7;
675d043ef4cSMark Adams               yy_t[grid][tid][pp++] = 0;
676d043ef4cSMark Adams               wp_t[grid][tid][pp]   = 0;
677d043ef4cSMark Adams               xx_t[grid][tid][pp]   = 1.e-7;
678d043ef4cSMark Adams               yy_t[grid][tid][pp++] = lo[1] + 5.e-7;
679d043ef4cSMark Adams             } else {
6806497c311SBarry Smith               const PetscInt p0 = NNreal;
6816497c311SBarry Smith               for (PetscInt pj = 0; pj < 6; pj++) { xx_t[grid][tid][p0 + pj] = yy_t[grid][tid][p0 + pj] = zz_t[grid][tid][p0 + pj] = wp_t[grid][tid][p0 + pj] = 0; }
682d043ef4cSMark Adams               xx_t[grid][tid][p0 + 0] = lo[0];
683d043ef4cSMark Adams               xx_t[grid][tid][p0 + 1] = hi[0];
684d043ef4cSMark Adams               yy_t[grid][tid][p0 + 2] = lo[1];
685d043ef4cSMark Adams               yy_t[grid][tid][p0 + 3] = hi[1];
686d043ef4cSMark Adams               zz_t[grid][tid][p0 + 4] = lo[2];
687d043ef4cSMark Adams               zz_t[grid][tid][p0 + 5] = hi[2];
688d043ef4cSMark Adams             }
6891b431c67SMark Adams             PetscCall(PetscRandomDestroy(&rand));
6903cdbf90fSMark Adams           }
6911b431c67SMark Adams           // entropy init, need global n
6923cdbf90fSMark Adams           if (glb_v_id == v_target) {
6931b431c67SMark Adams             const PetscReal N_inv = 1 / moments_0[0];
6943cdbf90fSMark Adams             PetscCall(PetscInfo(pack, "Target %" PetscInt_FMT " with %" PetscInt_FMT " particels\n", glb_v_id, nTargetP[0]));
6953cdbf90fSMark Adams             for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
6961b431c67SMark Adams               const PetscInt NN = nTargetP[grid];
6976497c311SBarry Smith               for (PetscInt pp = 0; pp < NN; pp++) {
6981b431c67SMark 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;
6991b431c67SMark Adams                 if (w > PETSC_REAL_MIN) {
7003cdbf90fSMark Adams                   moments_0[3] -= ww * PetscLogReal(ww);
7013cdbf90fSMark Adams                   PetscCheck(ww < 1 - PETSC_MACHINE_EPSILON, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "ww (%g) > 1", (double)ww);
7021b431c67SMark Adams                 } else moments_0[4] -= w;
7033cdbf90fSMark Adams               }
704f3237bfeSMark Adams             } // grid
7051b431c67SMark Adams           } // target
706f3237bfeSMark Adams         } // active
7073cdbf90fSMark Adams       } // threads
708f3237bfeSMark Adams       /* Create particle swarm */
7096497c311SBarry Smith       for (PetscInt tid = 0; tid < numthreads; tid++) {
7103cdbf90fSMark Adams         const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id;
7113cdbf90fSMark Adams         if (glb_v_id < num_vertices) {                             // the ragged edge of the last batch
712f3237bfeSMark Adams           for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids
713f3237bfeSMark Adams             PetscSection section;
714f3237bfeSMark Adams             PetscInt     Nf;
715f3237bfeSMark Adams             DM           dm = grid_dm[grid];
716d043ef4cSMark Adams             PetscCall(DMGetLocalSection(dm, &section));
717d043ef4cSMark Adams             PetscCall(PetscSectionGetNumFields(section, &Nf));
718*bd158744SPierre Jolivet             PetscCheck(Nf == 1, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Only one species per grid supported -- todo");
719d043ef4cSMark Adams             PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
720d043ef4cSMark Adams             PetscCall(PetscInfo(pack, "call createSwarm [%" PetscInt_FMT ".%" PetscInt_FMT "] local block index %" PetscInt_FMT "\n", v_id, grid, LAND_PACK_IDX(v_id, grid)));
721d043ef4cSMark Adams             PetscCall(createSwarm(dm, dim, &globSwarmArray[LAND_PACK_IDX(v_id, grid)]));
722f3237bfeSMark Adams           }
723f3237bfeSMark Adams         } // active
7243cdbf90fSMark Adams       } // threads
725d043ef4cSMark Adams       PetscCheck(ierr != 9999, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Only support one species per grid");
7263cdbf90fSMark Adams       // make globMpArray
727f3237bfeSMark Adams       PetscPragmaOMP(parallel for)
7286497c311SBarry Smith       for (PetscInt tid = 0; tid < numthreads; tid++) {
7293cdbf90fSMark Adams         const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id;
7303cdbf90fSMark Adams         if (glb_v_id < num_vertices) {
7313cdbf90fSMark Adams           for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids
732d043ef4cSMark Adams             // for (PetscInt sp = ctx->species_offset[grid], i0 = 0; sp < ctx->species_offset[grid + 1]; sp++, i0++) -- loop over species for Nf > 1 -- TODO
7333cdbf90fSMark Adams             PetscErrorCode ierr_t;
7343cdbf90fSMark Adams             DM             sw = globSwarmArray[LAND_PACK_IDX(v_id, grid)];
7351b431c67SMark Adams             ierr_t            = PetscInfo(pack, "makeSwarm %" PetscInt_FMT ".%" PetscInt_FMT ") for block %" PetscInt_FMT "\n", v_id, grid, LAND_PACK_IDX(v_id, grid));
7363cdbf90fSMark Adams             ierr_t            = makeSwarm(sw, dim, Np_t[grid][tid], xx_t[grid][tid], yy_t[grid][tid], zz_t[grid][tid]);
7373cdbf90fSMark Adams             if (ierr_t) ierr = ierr_t;
7383cdbf90fSMark Adams           }
7393cdbf90fSMark Adams         }
7403cdbf90fSMark Adams       }
7416497c311SBarry Smith       for (PetscInt tid = 0; tid < numthreads; tid++) {
7423cdbf90fSMark Adams         const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id;
7433cdbf90fSMark Adams         if (glb_v_id < num_vertices) {
744f3237bfeSMark Adams           for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids
745f3237bfeSMark Adams             DM dm = grid_dm[grid];
7463cdbf90fSMark Adams             DM sw = globSwarmArray[LAND_PACK_IDX(v_id, grid)];
747d043ef4cSMark Adams             PetscCall(PetscInfo(pack, "createMp %" PetscInt_FMT ".%" PetscInt_FMT ") for block %" PetscInt_FMT "\n", v_id, grid, LAND_PACK_IDX(v_id, grid)));
748d043ef4cSMark Adams             PetscCall(createMp(dm, sw, &globMpArray[LAND_PACK_IDX(v_id, grid)]));
749d043ef4cSMark Adams             PetscCall(MatViewFromOptions(globMpArray[LAND_PACK_IDX(v_id, grid)], NULL, "-mp_mat_view"));
7503cdbf90fSMark Adams           }
7513cdbf90fSMark Adams         }
7523cdbf90fSMark Adams       }
7533cdbf90fSMark Adams       // p --> g: set X
7543cdbf90fSMark Adams       // PetscPragmaOMP(parallel for)
7556497c311SBarry Smith       for (PetscInt tid = 0; tid < numthreads; tid++) {
7563cdbf90fSMark Adams         const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id;
7573cdbf90fSMark Adams         if (glb_v_id < num_vertices) {
7583cdbf90fSMark Adams           for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids
7593cdbf90fSMark Adams             PetscErrorCode ierr_t;
7603cdbf90fSMark Adams             DM             dm   = grid_dm[grid];
7613cdbf90fSMark Adams             DM             sw   = globSwarmArray[LAND_PACK_IDX(v_id, grid)];
7623cdbf90fSMark Adams             Vec            subX = globXArray[LAND_PACK_IDX(v_id, grid)], work = t_fhat[grid][tid];
7631b431c67SMark Adams             ierr_t = PetscInfo(pack, "particlesToGrid %" PetscInt_FMT ".%" PetscInt_FMT ") for block %" PetscInt_FMT "\n", v_id, grid, LAND_PACK_IDX(v_id, grid));
7641b431c67SMark Adams             ierr_t = particlesToGrid(dm, sw, tid, dim, wp_t[grid][tid], subX, globMpArray[LAND_PACK_IDX(v_id, grid)]);
765f3237bfeSMark Adams             if (ierr_t) ierr = ierr_t;
766f3237bfeSMark Adams             // u = M^_1 f_w
767f3237bfeSMark Adams             ierr_t = VecCopy(subX, work);
768f3237bfeSMark Adams             ierr_t = KSPSolve(t_ksp[grid][tid], work, subX);
769f3237bfeSMark Adams             if (ierr_t) ierr = ierr_t;
770f3237bfeSMark Adams           }
771f3237bfeSMark Adams         }
772f3237bfeSMark Adams       }
773d043ef4cSMark Adams       PetscCheck(!ierr, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Error in OMP loop. ierr = %d", (int)ierr);
774f3237bfeSMark Adams       /* Cleanup */
7756497c311SBarry Smith       for (PetscInt tid = 0; tid < numthreads; tid++) {
7763cdbf90fSMark Adams         const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id;
7773cdbf90fSMark Adams         if (glb_v_id < num_vertices) {
778f3237bfeSMark Adams           for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids
779f3237bfeSMark Adams             PetscCall(PetscFree4(xx_t[grid][tid], yy_t[grid][tid], wp_t[grid][tid], zz_t[grid][tid]));
780f3237bfeSMark Adams           }
781f3237bfeSMark Adams         } // active
7823cdbf90fSMark Adams       } // threads
7833cdbf90fSMark Adams     } // (fake) particle loop
7843cdbf90fSMark Adams     // standard view of initial conditions
7853cdbf90fSMark Adams     if (v_target >= global_vertex_id_0 && v_target < global_vertex_id_0 + ctx->batch_sz) {
7863cdbf90fSMark Adams       PetscCall(DMSetOutputSequenceNumber(ctx->plex[g_target], 0, 0.0));
7873cdbf90fSMark Adams       PetscCall(VecViewFromOptions(globXArray[LAND_PACK_IDX(v_target % ctx->batch_sz, g_target)], NULL, "-ex30_vec_view"));
7881b431c67SMark Adams       if (ctx->num_grids > g_target + 1) {
7891b431c67SMark Adams         PetscCall(DMSetOutputSequenceNumber(ctx->plex[g_target + 1], 0, 0.0));
7901b431c67SMark Adams         PetscCall(VecViewFromOptions(globXArray[LAND_PACK_IDX(v_target % ctx->batch_sz, g_target + 1)], NULL, "-ex30_vec_view2"));
7911b431c67SMark Adams       }
792d043ef4cSMark Adams       PetscCall(MatViewFromOptions(globMpArray[LAND_PACK_IDX(v_target % ctx->batch_sz, g_target)], NULL, "-ex30_mass_mat_view"));
793d043ef4cSMark Adams       PetscCall(DMViewFromOptions(globSwarmArray[LAND_PACK_IDX(v_target % ctx->batch_sz, g_target)], NULL, "-ex30_sw_view"));
794d043ef4cSMark Adams       PetscCall(DMSwarmViewXDMF(globSwarmArray[LAND_PACK_IDX(v_target % ctx->batch_sz, g_target)], "initial_swarm.xmf")); // writes a file by default!!!
795f3237bfeSMark Adams     }
796d043ef4cSMark Adams     // coarse graining moments_1a, bring f back from grid before advance
7971b431c67SMark Adams     if (v_target >= global_vertex_id_0 && v_target < global_vertex_id_0 + ctx->batch_sz && printCtx->print_entropy) {
7983cdbf90fSMark Adams       const PetscInt v_id = v_target % ctx->batch_sz;
799f3237bfeSMark Adams       for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
800f3237bfeSMark Adams         PetscDataType dtype;
801f3237bfeSMark Adams         PetscReal    *wp, *coords;
8023cdbf90fSMark Adams         DM            sw = globSwarmArray[LAND_PACK_IDX(v_id, grid)];
8033cdbf90fSMark Adams         Vec           work, subX = globXArray[LAND_PACK_IDX(v_id, grid)];
8041b431c67SMark Adams         PetscInt      bs, NN;
8053cdbf90fSMark Adams         // C-G moments
8063cdbf90fSMark Adams         PetscCall(VecDuplicate(subX, &work));
8073cdbf90fSMark Adams         PetscCall(gridToParticles(grid_dm[grid], sw, subX, work, globMpArray[LAND_PACK_IDX(v_id, grid)], g_Mass[grid]));
8083cdbf90fSMark Adams         PetscCall(VecDestroy(&work));
8093cdbf90fSMark Adams         // moments
810f3237bfeSMark Adams         PetscCall(DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
8111b431c67SMark Adams         PetscCall(DMSwarmGetLocalSize(sw, &NN));
8121b431c67SMark Adams         PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void **)&wp));
8136497c311SBarry Smith         for (PetscInt pp = 0; pp < NN; pp++) {
8141b431c67SMark 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]];
8156497c311SBarry Smith           for (PetscInt i = 0; i < dim; ++i) v2 += PetscSqr(coords[pp * dim + i]);
8163cdbf90fSMark Adams           moments_1a[0] += w;
8173cdbf90fSMark Adams           moments_1a[1] += w * ctx->v_0 * coords[pp * dim + 1]; // z-momentum
8181b431c67SMark Adams           moments_1a[2] += w * 0.5 * ctx->v_0 * ctx->v_0 * v2;
819f3237bfeSMark Adams         }
820f3237bfeSMark Adams         PetscCall(DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
821f3237bfeSMark Adams         PetscCall(DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void **)&wp));
822f3237bfeSMark Adams       }
8231b431c67SMark Adams       // entropy
8241b431c67SMark Adams       const PetscReal N_inv = 1 / moments_1a[0];
8251b431c67SMark Adams       PetscCall(PetscInfo(pack, "Entropy batch %" PetscInt_FMT " of %" PetscInt_FMT ", n = %g\n", v_target, num_vertices, (double)(1 / N_inv)));
8261b431c67SMark Adams       for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
8271b431c67SMark Adams         PetscDataType dtype;
8281b431c67SMark Adams         PetscReal    *wp, *coords;
8291b431c67SMark Adams         DM            sw = globSwarmArray[LAND_PACK_IDX(v_id, grid)];
8301b431c67SMark Adams         PetscInt      bs, NN;
8311b431c67SMark Adams         PetscCall(DMSwarmGetLocalSize(sw, &NN));
8321b431c67SMark Adams         PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void **)&wp));
8331b431c67SMark Adams         PetscCall(DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
8346497c311SBarry Smith         for (PetscInt pp = 0; pp < NN; pp++) {
8351b431c67SMark 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;
8361b431c67SMark Adams           if (w > PETSC_REAL_MIN) {
8371b431c67SMark Adams             moments_1a[3] -= ww * PetscLogReal(ww);
8381b431c67SMark Adams             PetscCheck(ww < 1 - PETSC_MACHINE_EPSILON, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "ww (%g) > 1", (double)ww);
8391b431c67SMark Adams           } else moments_1a[4] -= w;
8401b431c67SMark Adams         }
8411b431c67SMark Adams         PetscCall(DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void **)&wp));
8421b431c67SMark Adams         PetscCall(DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
8431b431c67SMark Adams       }
844f3237bfeSMark Adams     }
8453cdbf90fSMark Adams     // restore vector
846f3237bfeSMark Adams     PetscCall(DMCompositeRestoreAccessArray(pack, X, nDMs, NULL, globXArray));
847d043ef4cSMark Adams     // view initial grid
8481b431c67SMark Adams     if (v_target >= global_vertex_id_0 && v_target < global_vertex_id_0 + ctx->batch_sz) { PetscCall(DMPlexLandauPrintNorms(X, 0)); }
8493cdbf90fSMark Adams     // advance
8503cdbf90fSMark Adams     PetscCall(TSSetSolution(ts, X));
8511b431c67SMark Adams     PetscCall(PetscInfo(pack, "Advance vertex %" PetscInt_FMT " to %" PetscInt_FMT "\n", global_vertex_id_0, global_vertex_id_0 + ctx->batch_sz));
8523cdbf90fSMark Adams     PetscCall(TSSetPostStep(ts, PostStep));
8533cdbf90fSMark Adams     PetscCall(PostStep(ts));
8543cdbf90fSMark Adams     PetscCall(TSSolve(ts, X));
8553cdbf90fSMark Adams     // view
8563cdbf90fSMark Adams     PetscCall(DMCompositeGetAccessArray(pack, X, nDMs, NULL, globXArray));
8573cdbf90fSMark Adams     if (v_target >= global_vertex_id_0 && v_target < global_vertex_id_0 + ctx->batch_sz) {
858d043ef4cSMark Adams       /* Visualize original particle field */
8591b431c67SMark Adams       DM  sw = globSwarmArray[LAND_PACK_IDX(v_target % ctx->batch_sz, g_target)];
8601b431c67SMark Adams       Vec f;
8611b431c67SMark Adams       PetscCall(DMSetOutputSequenceNumber(sw, 0, 0.0));
8621b431c67SMark Adams       PetscCall(DMViewFromOptions(grid_dm[g_target], NULL, "-weights_dm_view"));
8631b431c67SMark Adams       PetscCall(DMViewFromOptions(sw, NULL, "-weights_sw_view"));
8641b431c67SMark Adams       PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &f));
8651b431c67SMark Adams       PetscCall(PetscObjectSetName((PetscObject)f, "weights"));
8661b431c67SMark Adams       PetscCall(VecViewFromOptions(f, NULL, "-weights_vec_view"));
8671b431c67SMark Adams       PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &f));
868d043ef4cSMark Adams       //
869d043ef4cSMark Adams       PetscCall(DMPlexLandauPrintNorms(X, 1));
8703cdbf90fSMark Adams     }
871d043ef4cSMark Adams     if (!use_uniform_particle_grid) { // resample to uniform grid
872d043ef4cSMark Adams       for (PetscInt v_id_0 = 0; v_id_0 < ctx->batch_sz; v_id_0 += numthreads) {
873d043ef4cSMark 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];
874d043ef4cSMark Adams         PetscInt   Np_t[LANDAU_MAX_GRIDS][EX30_MAX_NUM_THRDS];
8756497c311SBarry Smith         for (PetscInt tid = 0; tid < numthreads; tid++) {
876d043ef4cSMark Adams           const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id;
877d043ef4cSMark Adams           if (glb_v_id < num_vertices) {
878d043ef4cSMark Adams             // create uniform grid w/o weights & smaller
879d043ef4cSMark Adams             PetscInt Npp0 = (a_Np + (glb_v_id % (a_Np / 10 + 1))) / 2, Nv; // 1/2 of uniform particle grid size
880d043ef4cSMark Adams             for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
881d043ef4cSMark Adams               // for (PetscInt sp = ctx->species_offset[grid], i0 = 0; sp < ctx->species_offset[grid + 1]; sp++, i0++)
882d043ef4cSMark 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];
883d043ef4cSMark Adams               PetscInt  Npi = Npp0, Npj = 2 * Npp0, Npk = 1, NN;
884d043ef4cSMark Adams               // delete old particles and particle mass matrix
885d043ef4cSMark Adams               PetscCall(DMDestroy(&globSwarmArray[LAND_PACK_IDX(v_id, grid)]));
886d043ef4cSMark Adams               PetscCall(MatDestroy(&globMpArray[LAND_PACK_IDX(v_id, grid)]));
887d043ef4cSMark Adams               // create fake particles in batches with threads
888d043ef4cSMark Adams               PetscCall(MatGetLocalSize(g_Mass[grid], &Nv, NULL));
889d043ef4cSMark Adams               if (dim == 2) lo[0] = 0;
890d043ef4cSMark Adams               else Npi = Npj = Npk = Npp0;
891d043ef4cSMark Adams               NN = Npi * Npj * Npk + (dim == 2 ? 3 : 6); // make a regular grid of particles Npp x Npp
892d043ef4cSMark Adams               while (Npi * Npj * Npk < Nv) {             // make stable - no LS
893d043ef4cSMark Adams                 Npi++;
894d043ef4cSMark Adams                 Npj++;
895d043ef4cSMark Adams                 Npk++;
896d043ef4cSMark Adams                 NN = Npi * Npj * Npk + (dim == 2 ? 3 : 6);
897d043ef4cSMark Adams               }
898d043ef4cSMark Adams               Np_t[grid][tid] = NN;
899d043ef4cSMark 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]));
900d043ef4cSMark Adams               hp[0] = (hi[0] - lo[0]) / Npi;
901d043ef4cSMark Adams               hp[1] = (hi[1] - lo[1]) / Npj;
902d043ef4cSMark Adams               hp[2] = (hi[2] - lo[2]) / Npk;
903d043ef4cSMark Adams               if (dim == 2) hp[2] = 1;
904d043ef4cSMark Adams               PetscCall(PetscInfo(pack, "Resampling %d particles, %d vertices\n", (int)NN, (int)Nv)); // temp
9056497c311SBarry Smith               for (PetscInt pj = 0, pp = 0; pj < Npj; pj++) {
9066497c311SBarry Smith                 for (PetscInt pk = 0; pk < Npk; pk++) {
9076497c311SBarry Smith                   for (PetscInt pi = 0; pi < Npi; pi++, pp++) {
908d043ef4cSMark Adams                     wp_t[grid][tid][pp] = 0;
909d043ef4cSMark Adams                     xx_t[grid][tid][pp] = lo[0] + hp[0] / 2.0 + pi * hp[0];
910d043ef4cSMark Adams                     yy_t[grid][tid][pp] = lo[1] + hp[1] / 2.0 + pj * hp[1];
911d043ef4cSMark Adams                     if (dim == 3) zz_t[grid][tid][pp] = lo[2] + hp[2] / 2.0 + pk * hp[2];
912d043ef4cSMark Adams                   }
913d043ef4cSMark Adams                 }
914d043ef4cSMark Adams               }
915d043ef4cSMark Adams               if (dim == 2) { // fix bounding box
9166497c311SBarry Smith                 PetscInt pp           = NN - 3;
917d043ef4cSMark Adams                 wp_t[grid][tid][pp]   = 0;
918d043ef4cSMark Adams                 xx_t[grid][tid][pp]   = 1.e-7;
919d043ef4cSMark Adams                 yy_t[grid][tid][pp++] = hi[1] - 5.e-7;
920d043ef4cSMark Adams                 wp_t[grid][tid][pp]   = 0;
921d043ef4cSMark Adams                 xx_t[grid][tid][pp]   = hi[0] - 5.e-7;
922d043ef4cSMark Adams                 yy_t[grid][tid][pp++] = 0;
923d043ef4cSMark Adams                 wp_t[grid][tid][pp]   = 0;
924d043ef4cSMark Adams                 xx_t[grid][tid][pp]   = 1.e-7;
925d043ef4cSMark Adams                 yy_t[grid][tid][pp++] = lo[1] + 5.e-7;
926d043ef4cSMark Adams               } else {
9276497c311SBarry Smith                 const PetscInt p0 = NN - 6;
9286497c311SBarry Smith                 for (PetscInt pj = 0; pj < 6; pj++) { xx_t[grid][tid][p0 + pj] = yy_t[grid][tid][p0 + pj] = zz_t[grid][tid][p0 + pj] = wp_t[grid][tid][p0 + pj] = 0; }
929d043ef4cSMark Adams                 xx_t[grid][tid][p0 + 0] = lo[0];
930d043ef4cSMark Adams                 xx_t[grid][tid][p0 + 1] = hi[0];
931d043ef4cSMark Adams                 yy_t[grid][tid][p0 + 2] = lo[1];
932d043ef4cSMark Adams                 yy_t[grid][tid][p0 + 3] = hi[1];
933d043ef4cSMark Adams                 zz_t[grid][tid][p0 + 4] = lo[2];
934d043ef4cSMark Adams                 zz_t[grid][tid][p0 + 5] = hi[2];
935d043ef4cSMark Adams               }
936d043ef4cSMark Adams             }
937d043ef4cSMark Adams           } // active
938d043ef4cSMark Adams         } // threads
939d043ef4cSMark Adams         /* Create particle swarm */
9406497c311SBarry Smith         for (PetscInt tid = 0; tid < numthreads; tid++) {
941d043ef4cSMark Adams           const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id;
942d043ef4cSMark Adams           if (glb_v_id < num_vertices) {                             // the ragged edge of the last batch
943d043ef4cSMark Adams             for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids
944d043ef4cSMark Adams               // for (PetscInt sp = ctx->species_offset[grid], i0 = 0; sp < ctx->species_offset[grid + 1]; sp++, i0++) -- loop over species for Nf > 1 -- TODO
945d043ef4cSMark Adams               PetscErrorCode ierr_t;
946d043ef4cSMark Adams               PetscSection   section;
947d043ef4cSMark Adams               PetscInt       Nf;
948d043ef4cSMark Adams               DM             dm = grid_dm[grid];
949d043ef4cSMark Adams               ierr_t            = DMGetLocalSection(dm, &section);
950d043ef4cSMark Adams               ierr_t            = PetscSectionGetNumFields(section, &Nf);
951d043ef4cSMark Adams               if (Nf != 1) ierr_t = (PetscErrorCode)9999;
952d043ef4cSMark Adams               else {
953d043ef4cSMark Adams                 ierr_t = DMViewFromOptions(dm, NULL, "-dm_view");
954d043ef4cSMark 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));
955d043ef4cSMark Adams                 ierr_t = createSwarm(dm, dim, &globSwarmArray[LAND_PACK_IDX(v_id, grid)]);
956d043ef4cSMark Adams               }
957d043ef4cSMark Adams               if (ierr_t) ierr = ierr_t;
958d043ef4cSMark Adams             }
959d043ef4cSMark Adams           } // active
960d043ef4cSMark Adams         } // threads
961d043ef4cSMark Adams         PetscCheck(ierr != 9999, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Only support one species per grid");
962d043ef4cSMark Adams         PetscCheck(!ierr, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Error in OMP loop. ierr = %d", (int)ierr);
963d043ef4cSMark Adams         // make globMpArray
964d043ef4cSMark Adams         PetscPragmaOMP(parallel for)
9656497c311SBarry Smith         for (PetscInt tid = 0; tid < numthreads; tid++) {
966d043ef4cSMark Adams           const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id;
967d043ef4cSMark Adams           if (glb_v_id < num_vertices) {
968d043ef4cSMark Adams             for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids
969d043ef4cSMark Adams               // for (PetscInt sp = ctx->species_offset[grid], i0 = 0; sp < ctx->species_offset[grid + 1]; sp++, i0++) -- loop over species for Nf > 1 -- TODO
970d043ef4cSMark Adams               PetscErrorCode ierr_t;
971d043ef4cSMark Adams               DM             sw = globSwarmArray[LAND_PACK_IDX(v_id, grid)];
972d043ef4cSMark Adams               ierr_t            = PetscInfo(pack, "makeSwarm %" PetscInt_FMT ".%" PetscInt_FMT ") for block %" PetscInt_FMT "\n", v_id, grid, LAND_PACK_IDX(v_id, grid));
973d043ef4cSMark Adams               ierr_t            = makeSwarm(sw, dim, Np_t[grid][tid], xx_t[grid][tid], yy_t[grid][tid], zz_t[grid][tid]);
974d043ef4cSMark Adams               if (ierr_t) ierr = ierr_t;
975d043ef4cSMark Adams             }
976d043ef4cSMark Adams           } // active
977d043ef4cSMark Adams         } // threads
978d043ef4cSMark Adams         // create particle mass matrices
979d043ef4cSMark Adams         //PetscPragmaOMP(parallel for)
9806497c311SBarry Smith         for (PetscInt tid = 0; tid < numthreads; tid++) {
981d043ef4cSMark Adams           const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id;
982d043ef4cSMark Adams           if (glb_v_id < num_vertices) {
983d043ef4cSMark Adams             for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids
984d043ef4cSMark Adams               PetscErrorCode ierr_t;
985d043ef4cSMark Adams               DM             dm = grid_dm[grid];
986d043ef4cSMark Adams               DM             sw = globSwarmArray[LAND_PACK_IDX(v_id, grid)];
987d043ef4cSMark Adams               ierr_t            = PetscInfo(pack, "createMp %" PetscInt_FMT ".%" PetscInt_FMT ") for block %" PetscInt_FMT "\n", v_id, grid, LAND_PACK_IDX(v_id, grid));
988d043ef4cSMark Adams               ierr_t            = createMp(dm, sw, &globMpArray[LAND_PACK_IDX(v_id, grid)]);
989d043ef4cSMark Adams               if (ierr_t) ierr = ierr_t;
990d043ef4cSMark Adams             }
991d043ef4cSMark Adams           } // active
992d043ef4cSMark Adams         } // threads
993d043ef4cSMark Adams         PetscCheck(!ierr, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Error in OMP loop. ierr = %d", (int)ierr);
994d043ef4cSMark Adams         /* Cleanup */
9956497c311SBarry Smith         for (PetscInt tid = 0; tid < numthreads; tid++) {
996d043ef4cSMark Adams           const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id;
997d043ef4cSMark Adams           if (glb_v_id < num_vertices) {
998d043ef4cSMark Adams             for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids
999d043ef4cSMark Adams               PetscCall(PetscFree4(xx_t[grid][tid], yy_t[grid][tid], wp_t[grid][tid], zz_t[grid][tid]));
1000d043ef4cSMark Adams             }
1001d043ef4cSMark Adams           } // active
1002d043ef4cSMark Adams         } // threads
1003d043ef4cSMark Adams       } // batch
1004d043ef4cSMark Adams       // view
1005d043ef4cSMark Adams       if (v_target >= global_vertex_id_0 && v_target < global_vertex_id_0 + ctx->batch_sz) {
1006d043ef4cSMark Adams         /* Visualize particle field */
1007d043ef4cSMark Adams         DM  sw = globSwarmArray[LAND_PACK_IDX(v_target % ctx->batch_sz, g_target)];
1008d043ef4cSMark Adams         Vec f;
1009d043ef4cSMark Adams         PetscCall(DMSetOutputSequenceNumber(sw, 0, 0.0));
1010d043ef4cSMark Adams         PetscCall(DMViewFromOptions(sw, NULL, "-resampled_weights_sw_view"));
1011d043ef4cSMark Adams         PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &f));
1012d043ef4cSMark Adams         PetscCall(PetscObjectSetName((PetscObject)f, "resampled_weights"));
1013d043ef4cSMark Adams         PetscCall(VecViewFromOptions(f, NULL, "-resampled_weights_vec_view"));
1014d043ef4cSMark Adams         PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &f));
1015d043ef4cSMark Adams         PetscCall(DMSwarmViewXDMF(sw, "resampled.xmf"));
1016d043ef4cSMark Adams       }
1017d043ef4cSMark Adams     } // !uniform
1018d043ef4cSMark Adams     // particles to grid, compute moments and entropy, for target vertex only
10191b431c67SMark Adams     if (v_target >= global_vertex_id_0 && v_target < global_vertex_id_0 + ctx->batch_sz && printCtx->print_entropy) {
1020d043ef4cSMark Adams       PetscReal energy_error_rel;
10213cdbf90fSMark 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));
102257508eceSPierre Jolivet       energy_error_rel = PetscAbsReal(moments_1b[2] - moments_0[2]) / moments_0[2];
1023d043ef4cSMark Adams       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Particle Moments:\t number density      momentum (par)     energy             entropy            negative weights  : # OMP threads %g\n", (double)numthreads));
1024d043ef4cSMark Adams       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\tInitial:         %18.12e %19.12e %18.12e %18.12e %g %%\n", (double)moments_0[0], (double)moments_0[1], (double)moments_0[2], (double)moments_0[3], 100 * (double)(moments_0[4] / moments_0[0])));
1025d043ef4cSMark Adams       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\tCoarse-graining: %18.12e %19.12e %18.12e %18.12e %g %%\n", (double)moments_1a[0], (double)moments_1a[1], (double)moments_1a[2], (double)moments_1a[3], 100 * (double)(moments_1a[4] / moments_0[0])));
1026d043ef4cSMark Adams       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\tLandau:          %18.12e %19.12e %18.12e %18.12e %g %%\n", (double)moments_1b[0], (double)moments_1b[1], (double)moments_1b[2], (double)moments_1b[3], 100 * (double)(moments_1b[4] / moments_0[0])));
10271b431c67SMark 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])));
1028d043ef4cSMark Adams       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "(relative) energy conservation: Coarse-graining = %e, Landau = %e (%g %d)\n", (double)(moments_1a[2] - moments_0[2]) / (double)moments_0[2], (double)energy_error_rel, (double)PetscLog10Real(energy_error_rel), (int)(PetscLog10Real(energy_error_rel) + .5)));
10291b431c67SMark Adams     }
10303cdbf90fSMark Adams     // restore vector
10313cdbf90fSMark Adams     PetscCall(DMCompositeRestoreAccessArray(pack, X, nDMs, NULL, globXArray));
10323cdbf90fSMark Adams     // cleanup
10333cdbf90fSMark Adams     for (PetscInt v_id_0 = 0; v_id_0 < ctx->batch_sz; v_id_0 += numthreads) {
10346497c311SBarry Smith       for (PetscInt tid = 0; tid < numthreads; tid++) {
10353cdbf90fSMark Adams         const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id;
10363cdbf90fSMark Adams         if (glb_v_id < num_vertices) {
10373cdbf90fSMark Adams           for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
10383cdbf90fSMark Adams             PetscCall(DMDestroy(&globSwarmArray[LAND_PACK_IDX(v_id, grid)]));
10393cdbf90fSMark Adams             PetscCall(MatDestroy(&globMpArray[LAND_PACK_IDX(v_id, grid)]));
10403cdbf90fSMark Adams           }
10413cdbf90fSMark Adams         }
10423cdbf90fSMark Adams       }
10433cdbf90fSMark Adams     }
10441b431c67SMark Adams   } // user batch, not used
1045f3237bfeSMark Adams   /* Cleanup */
1046f3237bfeSMark Adams   PetscCall(PetscFree(globXArray));
1047f3237bfeSMark Adams   PetscCall(PetscFree(globSwarmArray));
1048f3237bfeSMark Adams   PetscCall(PetscFree(globMpArray));
10493cdbf90fSMark Adams   PetscCall(PetscFree(printCtx));
1050f3237bfeSMark Adams   // clean up mass matrices
1051f3237bfeSMark Adams   for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids
1052f3237bfeSMark Adams     PetscCall(MatDestroy(&g_Mass[grid]));
10536497c311SBarry Smith     for (PetscInt tid = 0; tid < numthreads; tid++) {
1054f3237bfeSMark Adams       PetscCall(VecDestroy(&t_fhat[grid][tid]));
1055f3237bfeSMark Adams       PetscCall(KSPDestroy(&t_ksp[grid][tid]));
1056f3237bfeSMark Adams     }
1057f3237bfeSMark Adams   }
10583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1059f3237bfeSMark Adams }
1060f3237bfeSMark Adams 
1061d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
1062d71ae5a4SJacob Faibussowitsch {
1063f3237bfeSMark Adams   DM         pack;
1064f3237bfeSMark Adams   Vec        X;
10651b431c67SMark Adams   PetscInt   dim = 2, num_vertices = 1, Np = 10, v_target = 0, g_target = 0;
1066f3237bfeSMark Adams   TS         ts;
1067f3237bfeSMark Adams   Mat        J;
1068f3237bfeSMark Adams   LandauCtx *ctx;
10693cdbf90fSMark Adams   PetscReal  shift                     = 0;
10701b431c67SMark Adams   PetscBool  use_uniform_particle_grid = PETSC_TRUE;
1071f3237bfeSMark Adams 
1072327415f7SBarry Smith   PetscFunctionBeginUser;
1073f3237bfeSMark Adams   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
1074f3237bfeSMark Adams   // process args
1075d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_SELF, "", "Collision Options", "DMPLEX");
1076f3237bfeSMark Adams   PetscCall(PetscOptionsInt("-dim", "Velocity space dimension", "ex30.c", dim, &dim, NULL));
10771b431c67SMark Adams   PetscCall(PetscOptionsInt("-number_spatial_vertices", "Number of user spatial vertices to be batched for Landau", "ex30.c", num_vertices, &num_vertices, NULL));
1078f3237bfeSMark 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));
10791b431c67SMark Adams   PetscCall(PetscOptionsBool("-use_uniform_particle_grid", "Use uniform particle grid", "ex30.c", use_uniform_particle_grid, &use_uniform_particle_grid, NULL));
1080d043ef4cSMark Adams   PetscCall(PetscOptionsInt("-vertex_view_target", "Global vertex for diagnostics", "ex30.c", v_target, &v_target, NULL));
1081d043ef4cSMark Adams   PetscCall(PetscOptionsReal("-e_shift", "Bi-Maxwellian shift", "ex30.c", shift, &shift, NULL));
1082d043ef4cSMark Adams   PetscCheck(v_target < num_vertices, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Batch to view %" PetscInt_FMT " should be < number of vertices %" PetscInt_FMT, v_target, num_vertices);
10831b431c67SMark Adams   PetscCall(PetscOptionsInt("-grid_view_target", "Grid to view with diagnostics", "ex30.c", g_target, &g_target, NULL));
1084d0609cedSBarry Smith   PetscOptionsEnd();
1085f3237bfeSMark Adams   /* Create a mesh */
1086f3237bfeSMark Adams   PetscCall(DMPlexLandauCreateVelocitySpace(PETSC_COMM_SELF, dim, "", &X, &J, &pack));
1087d043ef4cSMark Adams   PetscCall(DMGetApplicationContext(pack, &ctx));
1088f3237bfeSMark Adams   PetscCall(DMSetUp(pack));
1089f3237bfeSMark Adams   PetscCall(DMSetOutputSequenceNumber(pack, 0, 0.0));
10901b431c67SMark 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);
10911b431c67SMark 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);
1092d043ef4cSMark Adams   // PetscCheck(!use_uniform_particle_grid || !ctx->sphere, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Can not use -use_uniform_particle_grid and -dm_landau_sphere");
1093f3237bfeSMark Adams   /* Create timestepping solver context */
1094f3237bfeSMark Adams   PetscCall(TSCreate(PETSC_COMM_SELF, &ts));
1095f3237bfeSMark Adams   PetscCall(TSSetDM(ts, pack));
1096f3237bfeSMark Adams   PetscCall(TSSetIFunction(ts, NULL, DMPlexLandauIFunction, NULL));
1097f3237bfeSMark Adams   PetscCall(TSSetIJacobian(ts, J, J, DMPlexLandauIJacobian, NULL));
1098f3237bfeSMark Adams   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
1099f3237bfeSMark Adams   PetscCall(TSSetFromOptions(ts));
1100f3237bfeSMark Adams   PetscCall(PetscObjectSetName((PetscObject)X, "X"));
11016b664d00Smarkadams4   // do particle advance
11021b431c67SMark Adams   PetscCall(go(ts, X, num_vertices, Np, dim, v_target, g_target, shift, use_uniform_particle_grid));
1103984ed092SMark Adams   PetscCall(MatZeroEntries(J)); // need to zero out so as to not reuse it in Landau's logic
1104f3237bfeSMark Adams   /* clean up */
1105f3237bfeSMark Adams   PetscCall(DMPlexLandauDestroyVelocitySpace(&pack));
1106f3237bfeSMark Adams   PetscCall(TSDestroy(&ts));
1107f3237bfeSMark Adams   PetscCall(VecDestroy(&X));
1108f3237bfeSMark Adams   PetscCall(PetscFinalize());
1109f3237bfeSMark Adams   return 0;
1110f3237bfeSMark Adams }
1111f3237bfeSMark Adams 
1112f3237bfeSMark Adams /*TEST
1113f3237bfeSMark Adams 
1114f3237bfeSMark Adams   build:
11153cdbf90fSMark Adams     requires: !complex
1116f3237bfeSMark Adams 
1117f3237bfeSMark Adams   testset:
1118984ed092SMark Adams     requires: double defined(PETSC_USE_DMLANDAU_2D)
1119f3237bfeSMark Adams     output_file: output/ex30_0.out
1120d043ef4cSMark Adams     args: -dim 2 -petscspace_degree 3 -dm_landau_num_species_grid 1,1,1 -dm_refine 1 -number_particles_per_dimension 20 \
11211b431c67SMark Adams           -dm_landau_batch_size 4 -number_spatial_vertices 6 -vertex_view_target 5 -grid_view_target 1 -dm_landau_batch_view_idx 1 \
1122f3237bfeSMark 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 \
1123d043ef4cSMark 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 -ftop_ksp_error_if_not_converged \
1124d043ef4cSMark Adams           -ksp_type gmres -ksp_error_if_not_converged -dm_landau_verbose 4 -print_entropy \
1125d043ef4cSMark Adams           -ptof_ksp_type cg -ptof_pc_type jacobi -ptof_ksp_rtol 1e-12 -ptof_ksp_error_if_not_converged\
1126d043ef4cSMark Adams           -snes_converged_reason -snes_monitor -snes_rtol 1e-12 -snes_stol 1e-12 \
1127d043ef4cSMark 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
1128f3237bfeSMark Adams     test:
1129f3237bfeSMark Adams       suffix: cpu
1130d043ef4cSMark Adams       args: -dm_landau_device_type cpu -pc_type jacobi
1131f3237bfeSMark Adams     test:
1132f3237bfeSMark Adams       suffix: kokkos
1133bbdffb7fSJunchao Zhang       # failed on Sunspot@ALCF with sycl
1134bbdffb7fSJunchao Zhang       requires: kokkos_kernels !openmp !sycl
11353cdbf90fSMark 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
1136f3237bfeSMark Adams 
1137984ed092SMark Adams   testset:
1138984ed092SMark Adams     requires: double !defined(PETSC_USE_DMLANDAU_2D)
1139984ed092SMark Adams     output_file: output/ex30_3d.out
1140d043ef4cSMark Adams     args: -dim 3 -petscspace_degree 2 -dm_landau_num_species_grid 1,1 -dm_refine 0 -number_particles_per_dimension 10 -dm_plex_hash_location \
11411b431c67SMark Adams           -dm_landau_batch_size 1 -number_spatial_vertices 1 -vertex_view_target 0 -grid_view_target 0 -dm_landau_batch_view_idx 0 \
1142d043ef4cSMark Adams           -dm_landau_n 1.000018,1 -dm_landau_thermal_temps 2,1 -dm_landau_ion_masses 2 -dm_landau_ion_charges 1 \
1143d043ef4cSMark Adams           -ftop_ksp_type cg -ftop_pc_type jacobi -ftop_ksp_rtol 1e-12 -ftop_ksp_error_if_not_converged -ksp_type preonly -pc_type lu -ksp_error_if_not_converged \
1144d043ef4cSMark Adams           -ptof_ksp_type cg -ptof_pc_type jacobi -ptof_ksp_rtol 1e-12 -ptof_ksp_error_if_not_converged \
11453cdbf90fSMark Adams           -snes_converged_reason -snes_monitor -snes_rtol 1e-12 -snes_stol 1e-12 \
1146d043ef4cSMark 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
1147984ed092SMark Adams     test:
1148984ed092SMark Adams       suffix: cpu_3d
1149984ed092SMark Adams       args: -dm_landau_device_type cpu
1150984ed092SMark Adams     test:
1151984ed092SMark Adams       suffix: kokkos_3d
11523cdbf90fSMark Adams       requires: kokkos_kernels !openmp
11533cdbf90fSMark 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
11543cdbf90fSMark Adams 
1155d043ef4cSMark Adams   test:
1156d043ef4cSMark Adams     suffix: conserve
11573cdbf90fSMark Adams     requires: !complex double defined(PETSC_USE_DMLANDAU_2D) !cuda
1158d043ef4cSMark Adams     args: -dm_landau_batch_size 4 -dm_refine 0 -dm_landau_num_species_grid 1 -dm_landau_thermal_temps 1 -petscspace_degree 3 -snes_converged_reason -ts_type beuler -ts_dt .1 \
1159d043ef4cSMark Adams           -ts_max_steps 1 -ksp_type preonly -ksp_error_if_not_converged -snes_rtol 1e-14 -snes_stol 1e-14 -dm_landau_device_type cpu -number_particles_per_dimension 20 \
1160d043ef4cSMark Adams           -ptof_ksp_type cg -ptof_pc_type jacobi -ptof_ksp_rtol 1e-14 -ptof_ksp_error_if_not_converged -pc_type lu -dm_landau_simplex 1 -use_uniform_particle_grid false -dm_landau_sphere -print_entropy -number_particles_per_dimension 50 -ftop_ksp_type cg -ftop_pc_type jacobi -ftop_ksp_rtol 1e-14
1161984ed092SMark Adams 
1162f3237bfeSMark Adams TEST*/
1163