xref: /petsc/src/ts/tests/ex30.c (revision a259d2ee2575540d930166e4f236b09b11c05c1d)
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 
13f3237bfeSMark Adams #include "petscdmplex.h"
14f3237bfeSMark Adams #include "petscds.h"
15f3237bfeSMark Adams #include "petscdmswarm.h"
16f3237bfeSMark Adams #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 
31d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultMtM_SeqAIJ(Mat MtM, Vec xx, Vec yy)
32d71ae5a4SJacob Faibussowitsch {
33f3237bfeSMark Adams   MatShellCtx *matshellctx;
34f3237bfeSMark Adams 
35f3237bfeSMark Adams   PetscFunctionBeginUser;
36f3237bfeSMark Adams   PetscCall(MatShellGetContext(MtM, &matshellctx));
37f3237bfeSMark Adams   PetscCheck(matshellctx, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "No context");
38f3237bfeSMark Adams   PetscCall(MatMult(matshellctx->Mp, xx, matshellctx->ff));
39f3237bfeSMark Adams   PetscCall(MatMult(matshellctx->MpTrans, matshellctx->ff, yy));
403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
41f3237bfeSMark Adams }
42f3237bfeSMark Adams 
43d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultAddMtM_SeqAIJ(Mat MtM, Vec xx, Vec yy, Vec zz)
44d71ae5a4SJacob Faibussowitsch {
45f3237bfeSMark Adams   MatShellCtx *matshellctx;
46f3237bfeSMark Adams 
47f3237bfeSMark Adams   PetscFunctionBeginUser;
48f3237bfeSMark Adams   PetscCall(MatShellGetContext(MtM, &matshellctx));
49f3237bfeSMark Adams   PetscCheck(matshellctx, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "No context");
50f3237bfeSMark Adams   PetscCall(MatMult(matshellctx->Mp, xx, matshellctx->ff));
51f3237bfeSMark Adams   PetscCall(MatMultAdd(matshellctx->MpTrans, matshellctx->ff, yy, zz));
523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
53f3237bfeSMark Adams }
54f3237bfeSMark Adams 
55d71ae5a4SJacob Faibussowitsch PetscErrorCode createSwarm(const DM dm, PetscInt dim, DM *sw)
56d71ae5a4SJacob Faibussowitsch {
57f3237bfeSMark Adams   PetscInt Nc = 1;
58f3237bfeSMark Adams 
59f3237bfeSMark Adams   PetscFunctionBeginUser;
60f3237bfeSMark Adams   PetscCall(DMCreate(PETSC_COMM_SELF, sw));
61f3237bfeSMark Adams   PetscCall(DMSetType(*sw, DMSWARM));
62f3237bfeSMark Adams   PetscCall(DMSetDimension(*sw, dim));
63f3237bfeSMark Adams   PetscCall(DMSwarmSetType(*sw, DMSWARM_PIC));
64f3237bfeSMark Adams   PetscCall(DMSwarmSetCellDM(*sw, dm));
65f3237bfeSMark Adams   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "w_q", Nc, PETSC_SCALAR));
66f3237bfeSMark Adams   PetscCall(DMSwarmFinalizeFieldRegister(*sw));
67f3237bfeSMark Adams   PetscCall(DMSetFromOptions(*sw));
683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
69f3237bfeSMark Adams }
70f3237bfeSMark Adams 
71d71ae5a4SJacob Faibussowitsch PetscErrorCode gridToParticles(const DM dm, DM sw, Vec rhs, Vec work, Mat M_p, Mat Mass)
72d71ae5a4SJacob Faibussowitsch {
73f3237bfeSMark Adams   PetscBool    is_lsqr;
74f3237bfeSMark Adams   KSP          ksp;
75f3237bfeSMark Adams   Mat          PM_p = NULL, MtM, D;
76f3237bfeSMark Adams   Vec          ff;
77f3237bfeSMark Adams   PetscInt     N, M, nzl;
78f3237bfeSMark Adams   MatShellCtx *matshellctx;
79f3237bfeSMark Adams 
80f3237bfeSMark Adams   PetscFunctionBeginUser;
81f3237bfeSMark Adams   PetscCall(MatMult(Mass, rhs, work));
82f3237bfeSMark Adams   PetscCall(VecCopy(work, rhs));
83f3237bfeSMark Adams   // pseudo-inverse
84f3237bfeSMark Adams   PetscCall(KSPCreate(PETSC_COMM_SELF, &ksp));
85f3237bfeSMark Adams   PetscCall(KSPSetOptionsPrefix(ksp, "ftop_"));
86f3237bfeSMark Adams   PetscCall(KSPSetFromOptions(ksp));
87f3237bfeSMark Adams   PetscCall(PetscObjectTypeCompare((PetscObject)ksp, KSPLSQR, &is_lsqr));
88f3237bfeSMark Adams   if (!is_lsqr) {
89f3237bfeSMark Adams     PetscCall(MatGetLocalSize(M_p, &M, &N));
90f3237bfeSMark Adams     if (N > M) {
91f3237bfeSMark Adams       PC pc;
92f3237bfeSMark Adams       PetscCall(PetscInfo(ksp, " M (%" PetscInt_FMT ") < M (%" PetscInt_FMT ") -- skip revert to lsqr\n", M, N));
93f3237bfeSMark Adams       is_lsqr = PETSC_TRUE;
94f3237bfeSMark Adams       PetscCall(KSPSetType(ksp, KSPLSQR));
95f3237bfeSMark Adams       PetscCall(KSPGetPC(ksp, &pc));
96f3237bfeSMark Adams       PetscCall(PCSetType(pc, PCNONE)); // could put in better solver -ftop_pc_type bjacobi -ftop_sub_pc_type lu -ftop_sub_pc_factor_shift_type nonzero
97f3237bfeSMark Adams     } else {
98f3237bfeSMark Adams       PetscCall(PetscNew(&matshellctx));
99f3237bfeSMark Adams       PetscCall(MatCreateShell(PetscObjectComm((PetscObject)dm), N, N, PETSC_DECIDE, PETSC_DECIDE, matshellctx, &MtM));
100f3237bfeSMark Adams       PetscCall(MatTranspose(M_p, MAT_INITIAL_MATRIX, &matshellctx->MpTrans));
101f3237bfeSMark Adams       matshellctx->Mp = M_p;
102f3237bfeSMark Adams       PetscCall(MatShellSetOperation(MtM, MATOP_MULT, (void (*)(void))MatMultMtM_SeqAIJ));
103f3237bfeSMark Adams       PetscCall(MatShellSetOperation(MtM, MATOP_MULT_ADD, (void (*)(void))MatMultAddMtM_SeqAIJ));
104f3237bfeSMark Adams       PetscCall(MatCreateVecs(M_p, &matshellctx->uu, &matshellctx->ff));
105f3237bfeSMark Adams       PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, N, N, 1, NULL, &D));
106f3237bfeSMark Adams       PetscCall(MatViewFromOptions(matshellctx->MpTrans, NULL, "-ftop2_Mp_mat_view"));
107f3237bfeSMark Adams       for (int i = 0; i < N; i++) {
108f3237bfeSMark Adams         const PetscScalar *vals;
109f3237bfeSMark Adams         const PetscInt    *cols;
110f3237bfeSMark Adams         PetscScalar        dot = 0;
111f3237bfeSMark Adams         PetscCall(MatGetRow(matshellctx->MpTrans, i, &nzl, &cols, &vals));
112f3237bfeSMark Adams         for (int ii = 0; ii < nzl; ii++) dot += PetscSqr(vals[ii]);
11363a3b9bcSJacob Faibussowitsch         PetscCheck(dot != 0.0, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Row %d is empty", i);
114f3237bfeSMark Adams         PetscCall(MatSetValue(D, i, i, dot, INSERT_VALUES));
115f3237bfeSMark Adams       }
116f3237bfeSMark Adams       PetscCall(MatAssemblyBegin(D, MAT_FINAL_ASSEMBLY));
117f3237bfeSMark Adams       PetscCall(MatAssemblyEnd(D, MAT_FINAL_ASSEMBLY));
118f3237bfeSMark Adams       PetscCall(PetscInfo(M_p, "createMtMKSP Have %" PetscInt_FMT " eqs, nzl = %" PetscInt_FMT "\n", N, nzl));
119f3237bfeSMark Adams       PetscCall(KSPSetOperators(ksp, MtM, D));
120f3237bfeSMark Adams       PetscCall(MatViewFromOptions(D, NULL, "-ftop2_D_mat_view"));
121f3237bfeSMark Adams       PetscCall(MatViewFromOptions(M_p, NULL, "-ftop2_Mp_mat_view"));
122f3237bfeSMark Adams       PetscCall(MatViewFromOptions(matshellctx->MpTrans, NULL, "-ftop2_MpTranspose_mat_view"));
1236b664d00Smarkadams4       PetscCall(MatViewFromOptions(MtM, NULL, "-ftop2_MtM_mat_view"));
124f3237bfeSMark Adams     }
125f3237bfeSMark Adams   }
126f3237bfeSMark Adams   if (is_lsqr) {
127f3237bfeSMark Adams     PC        pc;
128f3237bfeSMark Adams     PetscBool is_bjac;
129f3237bfeSMark Adams     PetscCall(KSPGetPC(ksp, &pc));
130f3237bfeSMark Adams     PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCBJACOBI, &is_bjac));
131f3237bfeSMark Adams     if (is_bjac) {
132f3237bfeSMark Adams       PetscCall(DMSwarmCreateMassMatrixSquare(sw, dm, &PM_p));
133f3237bfeSMark Adams       PetscCall(KSPSetOperators(ksp, M_p, PM_p));
134f3237bfeSMark Adams     } else {
135f3237bfeSMark Adams       PetscCall(KSPSetOperators(ksp, M_p, M_p));
136f3237bfeSMark Adams     }
137f3237bfeSMark Adams   }
138f3237bfeSMark Adams   PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &ff)); // this grabs access
139f3237bfeSMark Adams   if (!is_lsqr) {
140f3237bfeSMark Adams     PetscCall(KSPSolve(ksp, rhs, matshellctx->uu));
141f3237bfeSMark Adams     PetscCall(MatMult(M_p, matshellctx->uu, ff));
142f3237bfeSMark Adams     PetscCall(MatDestroy(&matshellctx->MpTrans));
143f3237bfeSMark Adams     PetscCall(VecDestroy(&matshellctx->ff));
144f3237bfeSMark Adams     PetscCall(VecDestroy(&matshellctx->uu));
145f3237bfeSMark Adams     PetscCall(MatDestroy(&D));
146f3237bfeSMark Adams     PetscCall(MatDestroy(&MtM));
147f3237bfeSMark Adams     PetscCall(PetscFree(matshellctx));
148f3237bfeSMark Adams   } else {
149f3237bfeSMark Adams     PetscCall(KSPSolveTranspose(ksp, rhs, ff));
150f3237bfeSMark Adams   }
151f3237bfeSMark Adams   PetscCall(KSPDestroy(&ksp));
152f3237bfeSMark Adams   /* Visualize particle field */
153f3237bfeSMark Adams   PetscCall(VecViewFromOptions(ff, NULL, "-weights_view"));
154f3237bfeSMark Adams   PetscCall(MatDestroy(&PM_p));
155f3237bfeSMark Adams   PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &ff));
156f3237bfeSMark Adams 
1573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
158f3237bfeSMark Adams }
159f3237bfeSMark Adams 
160d71ae5a4SJacob Faibussowitsch PetscErrorCode particlesToGrid(const DM dm, DM sw, const PetscInt Np, const PetscInt a_tid, const PetscInt dim, const PetscReal xx[], const PetscReal yy[], const PetscReal zz[], const PetscReal a_wp[], Vec rho, Mat *Mp_out)
161d71ae5a4SJacob Faibussowitsch {
162f3237bfeSMark Adams   PetscBool     removePoints = PETSC_TRUE;
163f3237bfeSMark Adams   PetscReal    *wq, *coords;
164f3237bfeSMark Adams   PetscDataType dtype;
165f3237bfeSMark Adams   Mat           M_p;
166f3237bfeSMark Adams   Vec           ff;
167f3237bfeSMark Adams   PetscInt      bs, p, zero = 0;
168f3237bfeSMark Adams 
169f3237bfeSMark Adams   PetscFunctionBeginUser;
170f3237bfeSMark Adams   PetscCall(DMSwarmSetLocalSizes(sw, Np, zero));
171f3237bfeSMark Adams   PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void **)&wq));
172f3237bfeSMark Adams   PetscCall(DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
173f3237bfeSMark Adams   for (p = 0; p < Np; p++) {
174f3237bfeSMark Adams     coords[p * dim + 0] = xx[p];
175f3237bfeSMark Adams     coords[p * dim + 1] = yy[p];
176f3237bfeSMark Adams     wq[p]               = a_wp[p];
177f3237bfeSMark Adams     if (dim == 3) coords[p * dim + 2] = zz[p];
178f3237bfeSMark Adams   }
179f3237bfeSMark Adams   PetscCall(DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
180f3237bfeSMark Adams   PetscCall(DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void **)&wq));
181f3237bfeSMark Adams   PetscCall(DMSwarmMigrate(sw, removePoints));
182f3237bfeSMark Adams   PetscCall(PetscObjectSetName((PetscObject)sw, "Particle Grid"));
183f3237bfeSMark Adams 
184f3237bfeSMark Adams   /* This gives M f = \int_\Omega \phi f, which looks like a rhs for a PDE */
185f3237bfeSMark Adams   PetscCall(DMCreateMassMatrix(sw, dm, &M_p));
186f3237bfeSMark Adams 
187f3237bfeSMark Adams   PetscCall(PetscObjectSetName((PetscObject)rho, "rho"));
188f3237bfeSMark Adams   PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &ff));
189f3237bfeSMark Adams   PetscCall(PetscObjectSetName((PetscObject)ff, "weights"));
190f3237bfeSMark Adams   PetscCall(MatMultTranspose(M_p, ff, rho));
191f3237bfeSMark Adams   PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &ff));
192f3237bfeSMark Adams 
193f3237bfeSMark Adams   // output
194f3237bfeSMark Adams   *Mp_out = M_p;
195f3237bfeSMark Adams 
1963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
197f3237bfeSMark Adams }
198d71ae5a4SJacob Faibussowitsch static void maxwellian(PetscInt dim, const PetscReal x[], PetscReal kt_m, PetscReal n, PetscScalar *u)
199d71ae5a4SJacob Faibussowitsch {
200f3237bfeSMark Adams   PetscInt  i;
201f3237bfeSMark Adams   PetscReal v2 = 0, theta = 2.0 * kt_m; /* theta = 2kT/mc^2 */
202f3237bfeSMark Adams 
203f3237bfeSMark Adams   /* compute the exponents, v^2 */
204f3237bfeSMark Adams   for (i = 0; i < dim; ++i) v2 += x[i] * x[i];
205f3237bfeSMark Adams   /* evaluate the Maxwellian */
206f3237bfeSMark Adams   u[0] = n * PetscPowReal(PETSC_PI * theta, -1.5) * (PetscExpReal(-v2 / theta));
207f3237bfeSMark Adams }
208f3237bfeSMark Adams 
209f3237bfeSMark Adams #define MAX_NUM_THRDS 12
210d71ae5a4SJacob Faibussowitsch PetscErrorCode go(TS ts, Vec X, const PetscInt NUserV, const PetscInt a_Np, const PetscInt dim, const PetscInt b_target, const PetscInt g_target)
211d71ae5a4SJacob Faibussowitsch {
212f3237bfeSMark Adams   DM             pack, *globSwarmArray, grid_dm[LANDAU_MAX_GRIDS];
213f3237bfeSMark Adams   Mat           *globMpArray, g_Mass[LANDAU_MAX_GRIDS];
214f3237bfeSMark Adams   KSP            t_ksp[LANDAU_MAX_GRIDS][MAX_NUM_THRDS];
215f3237bfeSMark Adams   Vec            t_fhat[LANDAU_MAX_GRIDS][MAX_NUM_THRDS];
216f3237bfeSMark Adams   PetscInt       nDMs, glb_b_id, nTargetP = 0;
2176b664d00Smarkadams4   PetscErrorCode ierr = 0; // used for inside thread loops
218f3237bfeSMark Adams #if defined(PETSC_HAVE_OPENMP) && defined(PETSC_HAVE_THREADSAFETY)
219f3237bfeSMark Adams   PetscInt numthreads = PetscNumOMPThreads;
220f3237bfeSMark Adams #else
221f3237bfeSMark Adams   PetscInt numthreads = 1;
222f3237bfeSMark Adams #endif
223f3237bfeSMark Adams   LandauCtx *ctx;
224f3237bfeSMark Adams   Vec       *globXArray;
225f3237bfeSMark Adams   PetscReal  moments_0[3], moments_1[3], dt_init;
226f3237bfeSMark Adams 
227f3237bfeSMark Adams   PetscFunctionBeginUser;
22863a3b9bcSJacob Faibussowitsch   PetscCheck(numthreads <= MAX_NUM_THRDS, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Too many threads %" PetscInt_FMT " > %d", numthreads, MAX_NUM_THRDS);
22963a3b9bcSJacob Faibussowitsch   PetscCheck(numthreads > 0, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Number threads %" PetscInt_FMT " > %d", numthreads, MAX_NUM_THRDS);
230f3237bfeSMark Adams   PetscCall(TSGetDM(ts, &pack));
231f3237bfeSMark Adams   PetscCall(DMGetApplicationContext(pack, &ctx));
232f3237bfeSMark 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);
233f3237bfeSMark Adams   PetscCall(DMCompositeGetNumberDM(pack, &nDMs));
234f3237bfeSMark Adams   PetscCall(PetscInfo(pack, "Have %" PetscInt_FMT " total grids, with %" PetscInt_FMT " Landau local batched and %" PetscInt_FMT " global items (vertices)\n", ctx->num_grids, ctx->batch_sz, NUserV));
235f3237bfeSMark Adams   PetscCall(PetscMalloc(sizeof(*globXArray) * nDMs, &globXArray));
236f3237bfeSMark Adams   PetscCall(PetscMalloc(sizeof(*globMpArray) * nDMs, &globMpArray));
237f3237bfeSMark Adams   PetscCall(PetscMalloc(sizeof(*globSwarmArray) * nDMs, &globSwarmArray));
238f3237bfeSMark Adams   PetscCall(DMViewFromOptions(ctx->plex[g_target], NULL, "-ex30_dm_view"));
239f3237bfeSMark Adams   // create mass matrices
2406b664d00Smarkadams4   PetscCall(VecZeroEntries(X));
241f3237bfeSMark Adams   PetscCall(DMCompositeGetAccessArray(pack, X, nDMs, NULL, globXArray)); // just to duplicate
242f3237bfeSMark Adams   for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {               // add same particels for all grids
243f3237bfeSMark Adams     Vec          subX = globXArray[LAND_PACK_IDX(0, grid)];
244f3237bfeSMark Adams     DM           dm   = ctx->plex[grid];
245f3237bfeSMark Adams     PetscSection s;
246f3237bfeSMark Adams     grid_dm[grid] = dm;
247f3237bfeSMark Adams     PetscCall(DMCreateMassMatrix(dm, dm, &g_Mass[grid]));
248f3237bfeSMark Adams     //
249f3237bfeSMark Adams     PetscCall(DMGetLocalSection(dm, &s));
250f3237bfeSMark Adams     PetscCall(DMPlexCreateClosureIndex(dm, s));
251f3237bfeSMark Adams     for (int tid = 0; tid < numthreads; tid++) {
252f3237bfeSMark Adams       PetscCall(VecDuplicate(subX, &t_fhat[grid][tid]));
253f3237bfeSMark Adams       PetscCall(KSPCreate(PETSC_COMM_SELF, &t_ksp[grid][tid]));
254f3237bfeSMark Adams       PetscCall(KSPSetOptionsPrefix(t_ksp[grid][tid], "ptof_"));
255f3237bfeSMark Adams       PetscCall(KSPSetOperators(t_ksp[grid][tid], g_Mass[grid], g_Mass[grid]));
256f3237bfeSMark Adams       PetscCall(KSPSetFromOptions(t_ksp[grid][tid]));
257f3237bfeSMark Adams     }
258f3237bfeSMark Adams   }
259f3237bfeSMark Adams   PetscCall(DMCompositeRestoreAccessArray(pack, X, nDMs, NULL, globXArray));
260f3237bfeSMark Adams   // create particle raw data. could use OMP with a thread safe malloc, but this is just the fake user
261f3237bfeSMark Adams   for (int i = 0; i < 3; i++) moments_0[i] = moments_1[i] = 0;
262f3237bfeSMark Adams   PetscCall(TSGetTimeStep(ts, &dt_init)); // we could have an adaptive time stepper
263f3237bfeSMark Adams   for (PetscInt global_batch_id = 0; global_batch_id < NUserV; global_batch_id += ctx->batch_sz) {
2646b664d00Smarkadams4     PetscCall(TSSetTime(ts, 0));
2656b664d00Smarkadams4     PetscCall(TSSetStepNumber(ts, 0));
2666b664d00Smarkadams4     PetscCall(TSSetTimeStep(ts, dt_init));
267f3237bfeSMark Adams     PetscCall(DMCompositeGetAccessArray(pack, X, nDMs, NULL, globXArray));
26848a46eb9SPierre Jolivet     if (b_target >= global_batch_id && b_target < global_batch_id + ctx->batch_sz) PetscCall(PetscObjectSetName((PetscObject)globXArray[LAND_PACK_IDX(b_target % ctx->batch_sz, g_target)], "rho"));
269f3237bfeSMark Adams     // create fake particles
270f3237bfeSMark Adams     for (PetscInt b_id_0 = 0; b_id_0 < ctx->batch_sz; b_id_0 += numthreads) {
271f3237bfeSMark Adams       PetscReal *xx_t[LANDAU_MAX_GRIDS][MAX_NUM_THRDS], *yy_t[LANDAU_MAX_GRIDS][MAX_NUM_THRDS], *zz_t[LANDAU_MAX_GRIDS][MAX_NUM_THRDS], *wp_t[LANDAU_MAX_GRIDS][MAX_NUM_THRDS];
272f3237bfeSMark Adams       PetscInt   Np_t[LANDAU_MAX_GRIDS][MAX_NUM_THRDS];
273f3237bfeSMark Adams       // make particles
274f3237bfeSMark Adams       for (int tid = 0; tid < numthreads; tid++) {
275f3237bfeSMark Adams         const PetscInt b_id = b_id_0 + tid;
276f3237bfeSMark Adams         if ((glb_b_id = global_batch_id + b_id) < NUserV) {        // the ragged edge of the last batch
277f3237bfeSMark Adams           PetscInt Npp0 = a_Np + (glb_b_id % a_Np), NN;            // fake user: number of particels in each dimension with add some load imbalance and diff (<2x)
278f3237bfeSMark Adams           for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids
2799371c9d4SSatish Balay             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 -- TODO */
2809371c9d4SSatish Balay             ;
281f3237bfeSMark 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
282f3237bfeSMark Adams             PetscInt  Npi = Npp0, Npj = 2 * Npp0, Npk = 1;
283f3237bfeSMark Adams             if (dim == 2) lo[0] = 0; // Landau coordinate (r,z)
284f3237bfeSMark Adams             else Npi = Npj = Npk = Npp0;
285f3237bfeSMark Adams             // User: use glb_b_id to index into your data
286f3237bfeSMark Adams             NN = Npi * Npj * Npk; // make a regular grid of particles Npp x Npp
287f3237bfeSMark Adams             if (glb_b_id == b_target) {
288f3237bfeSMark Adams               nTargetP = NN;
289f3237bfeSMark Adams               PetscCall(PetscInfo(pack, "Target %" PetscInt_FMT " with %" PetscInt_FMT " particels\n", glb_b_id, NN));
290f3237bfeSMark Adams             }
291f3237bfeSMark Adams             Np_t[grid][tid] = NN;
292f3237bfeSMark 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]));
293f3237bfeSMark Adams             hp[0] = (hi[0] - lo[0]) / Npi;
294f3237bfeSMark Adams             hp[1] = (hi[1] - lo[1]) / Npj;
295f3237bfeSMark Adams             hp[2] = (hi[2] - lo[2]) / Npk;
296f3237bfeSMark Adams             if (dim == 2) hp[2] = 1;
29763a3b9bcSJacob 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
298f3237bfeSMark Adams             vole = hp[0] * hp[1] * hp[2] * ctx->n[grid];                                                                                                                           // fix for multi-species
299f3237bfeSMark Adams             PetscCall(PetscInfo(pack, "Vertex %" PetscInt_FMT ", grid %" PetscInt_FMT " with %" PetscInt_FMT " particles (diagnostic target = %" PetscInt_FMT ")\n", glb_b_id, grid, NN, b_target));
300f3237bfeSMark Adams             for (int pj = 0, pp = 0; pj < Npj; pj++) {
301f3237bfeSMark Adams               for (int pk = 0; pk < Npk; pk++) {
302f3237bfeSMark Adams                 for (int pi = 0; pi < Npi; pi++, pp++) {
303f3237bfeSMark Adams                   xx_t[grid][tid][pp] = lo[0] + hp[0] / 2.0 + pi * hp[0];
304f3237bfeSMark Adams                   yy_t[grid][tid][pp] = lo[1] + hp[1] / 2.0 + pj * hp[1];
305f3237bfeSMark Adams                   if (dim == 3) zz_t[grid][tid][pp] = lo[2] + hp[2] / 2.0 + pk * hp[2];
306f3237bfeSMark Adams                   {
307f3237bfeSMark Adams                     PetscReal x[] = {xx_t[grid][tid][pp], yy_t[grid][tid][pp], dim == 2 ? 0 : zz_t[grid][tid][pp]};
308f3237bfeSMark Adams                     maxwellian(dim, x, kT_m, vole, &wp_t[grid][tid][pp]);
309f3237bfeSMark Adams                     //PetscCall(PetscInfo(pack,"%" PetscInt_FMT ") x = %14.7e, %14.7e, %14.7e, n = %14.7e, w = %14.7e\n", pp, x[0], x[1], dim==2 ? 0 : x[2], ctx->n[grid], wp_t[grid][tid][pp])); // temp
310f3237bfeSMark Adams                     if (glb_b_id == b_target) {
311f3237bfeSMark Adams                       PetscReal v2 = 0, fact = dim == 2 ? 2.0 * PETSC_PI * x[0] : 1;
312f3237bfeSMark Adams                       for (int i = 0; i < dim; ++i) v2 += PetscSqr(x[i]);
313f3237bfeSMark Adams                       moments_0[0] += fact * wp_t[grid][tid][pp] * ctx->n_0 * ctx->masses[ctx->species_offset[grid]];
314f3237bfeSMark Adams                       moments_0[1] += fact * wp_t[grid][tid][pp] * ctx->n_0 * ctx->v_0 * ctx->masses[ctx->species_offset[grid]] * x[1]; // z-momentum
315f3237bfeSMark Adams                       moments_0[2] += fact * wp_t[grid][tid][pp] * ctx->n_0 * ctx->v_0 * ctx->v_0 * ctx->masses[ctx->species_offset[grid]] * v2;
316f3237bfeSMark Adams                     }
317f3237bfeSMark Adams                   }
318f3237bfeSMark Adams                 }
319f3237bfeSMark Adams               }
320f3237bfeSMark Adams             }
321f3237bfeSMark Adams           } // grid
322f3237bfeSMark Adams         }   // active
323f3237bfeSMark Adams       }     // fake threads
324f3237bfeSMark Adams       /* Create particle swarm */
325f3237bfeSMark Adams       PetscPragmaOMP(parallel for)
326f37bacd1SJacob Faibussowitsch       for (int tid = 0; tid < numthreads; tid++) {
327f3237bfeSMark Adams         const PetscInt b_id = b_id_0 + tid;
328f3237bfeSMark Adams         if ((glb_b_id = global_batch_id + b_id) < NUserV) { // the ragged edge of the last batch
329f3237bfeSMark Adams           //PetscCall(PetscInfo(pack,"Create swarms for 'glob' index %" PetscInt_FMT " create swarm\n",glb_b_id));
330f3237bfeSMark Adams           for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids
331f3237bfeSMark Adams             PetscErrorCode ierr_t;
332f3237bfeSMark Adams             PetscSection   section;
333f3237bfeSMark Adams             PetscInt       Nf;
334f3237bfeSMark Adams             DM             dm = grid_dm[grid];
335f3237bfeSMark Adams             ierr_t            = DMGetLocalSection(dm, &section);
336f3237bfeSMark Adams             ierr_t            = PetscSectionGetNumFields(section, &Nf);
337f3237bfeSMark Adams             if (Nf != 1) ierr_t = 9999;
338f3237bfeSMark Adams             else {
339f3237bfeSMark Adams               ierr_t = DMViewFromOptions(dm, NULL, "-dm_view");
340f3237bfeSMark Adams               ierr_t = PetscInfo(pack, "call createSwarm [%" PetscInt_FMT ".%" PetscInt_FMT "] local batch index %" PetscInt_FMT "\n", b_id, grid, LAND_PACK_IDX(b_id, grid));
341f3237bfeSMark Adams               ierr_t = createSwarm(dm, dim, &globSwarmArray[LAND_PACK_IDX(b_id, grid)]);
342f3237bfeSMark Adams             }
343f3237bfeSMark Adams             if (ierr_t) ierr = ierr_t;
344f3237bfeSMark Adams           }
345f3237bfeSMark Adams         } // active
346f3237bfeSMark Adams       }
347cad9d221SBarry Smith       PetscCheck(ierr != 9999, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Only support one species per grid");
3486b664d00Smarkadams4       PetscCheck(!ierr, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Error in OMP loop. ierr = %d", (int)ierr);
349f3237bfeSMark Adams       // p --> g: make globMpArray & set X
350f3237bfeSMark Adams       PetscPragmaOMP(parallel for)
351f37bacd1SJacob Faibussowitsch       for (int tid = 0; tid < numthreads; tid++) {
352f3237bfeSMark Adams         const PetscInt b_id = b_id_0 + tid;
353f3237bfeSMark Adams         if ((glb_b_id = global_batch_id + b_id) < NUserV) {
354f3237bfeSMark Adams           for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids
355f3237bfeSMark Adams             PetscErrorCode ierr_t;
356f3237bfeSMark Adams             DM             dm   = grid_dm[grid];
357f3237bfeSMark Adams             DM             sw   = globSwarmArray[LAND_PACK_IDX(b_id, grid)];
358f3237bfeSMark Adams             Vec            subX = globXArray[LAND_PACK_IDX(b_id, grid)], work = t_fhat[grid][tid];
3597864358aSSatish Balay             ierr_t = PetscInfo(pack, "particlesToGrid %" PetscInt_FMT ".%" PetscInt_FMT ") particlesToGrid for local batch %" PetscInt_FMT "\n", global_batch_id, grid, LAND_PACK_IDX(b_id, grid));
360f3237bfeSMark Adams             ierr_t = particlesToGrid(dm, sw, Np_t[grid][tid], tid, dim, xx_t[grid][tid], yy_t[grid][tid], zz_t[grid][tid], wp_t[grid][tid], subX, &globMpArray[LAND_PACK_IDX(b_id, grid)]);
361f3237bfeSMark Adams             if (ierr_t) ierr = ierr_t;
362f3237bfeSMark Adams             // u = M^_1 f_w
363f3237bfeSMark Adams             ierr_t = VecCopy(subX, work);
364f3237bfeSMark Adams             ierr_t = KSPSolve(t_ksp[grid][tid], work, subX);
365f3237bfeSMark Adams             if (ierr_t) ierr = ierr_t;
366f3237bfeSMark Adams           }
367f3237bfeSMark Adams         }
368f3237bfeSMark Adams       }
3696b664d00Smarkadams4       PetscCheck(!ierr, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Error in OMP loop. ierr = %d", (int)ierr);
370f3237bfeSMark Adams       /* Cleanup */
371f3237bfeSMark Adams       for (int tid = 0; tid < numthreads; tid++) {
372f3237bfeSMark Adams         const PetscInt b_id = b_id_0 + tid;
373f3237bfeSMark Adams         if ((glb_b_id = global_batch_id + b_id) < NUserV) {
374f3237bfeSMark Adams           PetscCall(PetscInfo(pack, "Free for global batch %" PetscInt_FMT " of %" PetscInt_FMT "\n", glb_b_id + 1, NUserV));
375f3237bfeSMark Adams           for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids
376f3237bfeSMark Adams             PetscCall(PetscFree4(xx_t[grid][tid], yy_t[grid][tid], wp_t[grid][tid], zz_t[grid][tid]));
377f3237bfeSMark Adams           }
378f3237bfeSMark Adams         } // active
379f3237bfeSMark Adams       }
3806b664d00Smarkadams4     } // batches
38148a46eb9SPierre Jolivet     if (b_target >= global_batch_id && b_target < global_batch_id + ctx->batch_sz) PetscCall(VecViewFromOptions(globXArray[LAND_PACK_IDX(b_target % ctx->batch_sz, g_target)], NULL, "-ex30_vec_view"));
382f3237bfeSMark Adams     PetscCall(DMCompositeRestoreAccessArray(pack, X, nDMs, NULL, globXArray));
383f3237bfeSMark Adams     PetscCall(DMPlexLandauPrintNorms(X, 0));
384f3237bfeSMark Adams     // advance
385f3237bfeSMark Adams     PetscCall(TSSetSolution(ts, X));
386f3237bfeSMark Adams     PetscCall(PetscInfo(pack, "Advance vertex %" PetscInt_FMT " to %" PetscInt_FMT " (with padding)\n", global_batch_id, global_batch_id + ctx->batch_sz));
387f3237bfeSMark Adams     PetscCall(TSSolve(ts, X));
388f3237bfeSMark Adams     PetscCall(DMPlexLandauPrintNorms(X, 1));
389f3237bfeSMark Adams     PetscCall(DMCompositeGetAccessArray(pack, X, nDMs, NULL, globXArray));
390f3237bfeSMark Adams     // map back to particles
391f3237bfeSMark Adams     for (PetscInt b_id_0 = 0; b_id_0 < ctx->batch_sz; b_id_0 += numthreads) {
392f3237bfeSMark Adams       PetscCall(PetscInfo(pack, "g2p: global batch %" PetscInt_FMT " of %" PetscInt_FMT ", Landau batch %" PetscInt_FMT " of %" PetscInt_FMT ": map back to particles\n", global_batch_id + 1, NUserV, b_id_0 + 1, ctx->batch_sz));
393f3237bfeSMark Adams       PetscPragmaOMP(parallel for)
394f37bacd1SJacob Faibussowitsch       for (int tid = 0; tid < numthreads; tid++) {
395f3237bfeSMark Adams         const PetscInt b_id = b_id_0 + tid;
396f3237bfeSMark Adams         if ((glb_b_id = global_batch_id + b_id) < NUserV) {
397f3237bfeSMark Adams           for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids
398f3237bfeSMark Adams             PetscErrorCode ierr_t;
3997864358aSSatish Balay             ierr_t = PetscInfo(pack, "gridToParticles: global batch %" PetscInt_FMT ", local batch b=%" PetscInt_FMT ", grid g=%" PetscInt_FMT ", index(b,g) %" PetscInt_FMT "\n", global_batch_id, b_id, grid, LAND_PACK_IDX(b_id, grid));
400f3237bfeSMark Adams             ierr_t = gridToParticles(grid_dm[grid], globSwarmArray[LAND_PACK_IDX(b_id, grid)], globXArray[LAND_PACK_IDX(b_id, grid)], t_fhat[grid][tid], globMpArray[LAND_PACK_IDX(b_id, grid)], g_Mass[grid]);
401f3237bfeSMark Adams             if (ierr_t) ierr = ierr_t;
402f3237bfeSMark Adams           }
403f3237bfeSMark Adams         }
404f3237bfeSMark Adams       }
4056b664d00Smarkadams4       PetscCheck(!ierr, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Error in OMP loop. ierr = %d", (int)ierr);
406f3237bfeSMark Adams       /* Cleanup, and get data */
407f3237bfeSMark Adams       PetscCall(PetscInfo(pack, "Cleanup batches %" PetscInt_FMT " to %" PetscInt_FMT "\n", b_id_0, b_id_0 + numthreads));
408f3237bfeSMark Adams       for (int tid = 0; tid < numthreads; tid++) {
409f3237bfeSMark Adams         const PetscInt b_id = b_id_0 + tid;
410f3237bfeSMark Adams         if ((glb_b_id = global_batch_id + b_id) < NUserV) {
411f3237bfeSMark Adams           for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
412f3237bfeSMark Adams             PetscDataType dtype;
413f3237bfeSMark Adams             PetscReal    *wp, *coords;
414f3237bfeSMark Adams             DM            sw = globSwarmArray[LAND_PACK_IDX(b_id, grid)];
415f3237bfeSMark Adams             PetscInt      npoints, bs = 1;
416f3237bfeSMark Adams             PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void **)&wp)); // take data out here
417f3237bfeSMark Adams             if (glb_b_id == b_target) {
418f3237bfeSMark Adams               PetscCall(DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
419f3237bfeSMark Adams               PetscCall(DMSwarmGetLocalSize(sw, &npoints));
420f3237bfeSMark Adams               for (int p = 0; p < npoints; p++) {
421f3237bfeSMark Adams                 PetscReal v2 = 0, fact = dim == 2 ? 2.0 * PETSC_PI * coords[p * dim + 0] : 1;
422f3237bfeSMark Adams                 for (int i = 0; i < dim; ++i) v2 += PetscSqr(coords[p * dim + i]);
423f3237bfeSMark Adams                 moments_1[0] += fact * wp[p] * ctx->n_0 * ctx->masses[ctx->species_offset[grid]];
424f3237bfeSMark Adams                 moments_1[1] += fact * wp[p] * ctx->n_0 * ctx->v_0 * ctx->masses[ctx->species_offset[grid]] * coords[p * dim + 1]; // z-momentum
425f3237bfeSMark Adams                 moments_1[2] += fact * wp[p] * ctx->n_0 * ctx->v_0 * ctx->v_0 * ctx->masses[ctx->species_offset[grid]] * v2;
426f3237bfeSMark Adams               }
427f3237bfeSMark Adams               PetscCall(DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
428f3237bfeSMark Adams             }
429f3237bfeSMark Adams             PetscCall(DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void **)&wp));
430f3237bfeSMark Adams             PetscCall(DMDestroy(&globSwarmArray[LAND_PACK_IDX(b_id, grid)]));
431f3237bfeSMark Adams             PetscCall(MatDestroy(&globMpArray[LAND_PACK_IDX(b_id, grid)]));
432f3237bfeSMark Adams           }
433f3237bfeSMark Adams         }
434f3237bfeSMark Adams       }
435f3237bfeSMark Adams     } // thread batch
436f3237bfeSMark Adams     PetscCall(DMCompositeRestoreAccessArray(pack, X, nDMs, NULL, globXArray));
437f3237bfeSMark Adams   } // user batch
438f3237bfeSMark Adams   /* Cleanup */
439f3237bfeSMark Adams   PetscCall(PetscFree(globXArray));
440f3237bfeSMark Adams   PetscCall(PetscFree(globSwarmArray));
441f3237bfeSMark Adams   PetscCall(PetscFree(globMpArray));
442f3237bfeSMark Adams   // clean up mass matrices
443f3237bfeSMark Adams   for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids
444f3237bfeSMark Adams     PetscCall(MatDestroy(&g_Mass[grid]));
445f3237bfeSMark Adams     for (int tid = 0; tid < numthreads; tid++) {
446f3237bfeSMark Adams       PetscCall(VecDestroy(&t_fhat[grid][tid]));
447f3237bfeSMark Adams       PetscCall(KSPDestroy(&t_ksp[grid][tid]));
448f3237bfeSMark Adams     }
449f3237bfeSMark Adams   }
450*a259d2eeSmarkadams4   {
451*a259d2eeSmarkadams4     PetscReal relerr = PetscAbsReal((moments_1[2] - moments_0[2]) / moments_0[2]), logEerr = PetscLog10Real(relerr);
452*a259d2eeSmarkadams4     PetscCall(PetscInfo(X, "Total number density: %20.12e (%20.12e); x-momentum = %20.12e (%20.12e); energy = %20.12e (%20.12e) rerr = %e (log10(rerr) = %e, %" PetscInt_FMT "), %" PetscInt_FMT " particles. Use %" PetscInt_FMT " threads\n", (double)moments_1[0], (double)moments_0[0], (double)moments_1[1], (double)moments_0[1], (double)moments_1[2], (double)moments_0[2], (double)relerr, (double)logEerr, -((-(PetscInt)logEerr + 1) / 2) * 2, nTargetP, numthreads));
453*a259d2eeSmarkadams4     PetscFunctionReturn(PETSC_SUCCESS);
454*a259d2eeSmarkadams4   }
4553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
456f3237bfeSMark Adams }
457f3237bfeSMark Adams 
458d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
459d71ae5a4SJacob Faibussowitsch {
460f3237bfeSMark Adams   DM         pack;
461f3237bfeSMark Adams   Vec        X;
462f3237bfeSMark Adams   PetscInt   dim = 2, nvert = 1, Np = 10, btarget = 0, gtarget = 0;
463f3237bfeSMark Adams   TS         ts;
464f3237bfeSMark Adams   Mat        J;
465f3237bfeSMark Adams   LandauCtx *ctx;
466f3237bfeSMark Adams 
467327415f7SBarry Smith   PetscFunctionBeginUser;
468f3237bfeSMark Adams   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
469f3237bfeSMark Adams   // process args
470d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_SELF, "", "Collision Options", "DMPLEX");
471f3237bfeSMark Adams   PetscCall(PetscOptionsInt("-number_spatial_vertices", "Number of user spatial vertices to be batched for Landau", "ex30.c", nvert, &nvert, NULL));
472f3237bfeSMark Adams   PetscCall(PetscOptionsInt("-dim", "Velocity space dimension", "ex30.c", dim, &dim, NULL));
473f3237bfeSMark 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));
474f3237bfeSMark Adams   PetscCall(PetscOptionsInt("-view_vertex_target", "Batch to view with diagnostics", "ex30.c", btarget, &btarget, NULL));
475f3237bfeSMark Adams   PetscCheck(btarget < nvert, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Batch to view %" PetscInt_FMT " should be < number of vertices %" PetscInt_FMT, btarget, nvert);
476f3237bfeSMark Adams   PetscCall(PetscOptionsInt("-view_grid_target", "Grid to view with diagnostics", "ex30.c", gtarget, &gtarget, NULL));
477d0609cedSBarry Smith   PetscOptionsEnd();
478f3237bfeSMark Adams   /* Create a mesh */
479f3237bfeSMark Adams   PetscCall(DMPlexLandauCreateVelocitySpace(PETSC_COMM_SELF, dim, "", &X, &J, &pack));
480f3237bfeSMark Adams   PetscCall(DMSetUp(pack));
481f3237bfeSMark Adams   PetscCall(DMSetOutputSequenceNumber(pack, 0, 0.0));
482f3237bfeSMark Adams   PetscCall(DMGetApplicationContext(pack, &ctx));
483f3237bfeSMark Adams   PetscCheck(gtarget < ctx->num_grids, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Grid to view %" PetscInt_FMT " should be < number of grids %" PetscInt_FMT, gtarget, ctx->num_grids);
484f3237bfeSMark Adams   PetscCheck(nvert >= ctx->batch_sz, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Number of vertices %" PetscInt_FMT " should be <= batch size %" PetscInt_FMT, nvert, ctx->batch_sz);
485f3237bfeSMark Adams   /* Create timestepping solver context */
486f3237bfeSMark Adams   PetscCall(TSCreate(PETSC_COMM_SELF, &ts));
487f3237bfeSMark Adams   PetscCall(TSSetDM(ts, pack));
488f3237bfeSMark Adams   PetscCall(TSSetIFunction(ts, NULL, DMPlexLandauIFunction, NULL));
489f3237bfeSMark Adams   PetscCall(TSSetIJacobian(ts, J, J, DMPlexLandauIJacobian, NULL));
490f3237bfeSMark Adams   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
491f3237bfeSMark Adams   PetscCall(TSSetFromOptions(ts));
492f3237bfeSMark Adams   PetscCall(PetscObjectSetName((PetscObject)X, "X"));
4936b664d00Smarkadams4   // do particle advance
494f3237bfeSMark Adams   PetscCall(go(ts, X, nvert, Np, dim, btarget, gtarget));
495984ed092SMark Adams   PetscCall(MatZeroEntries(J)); // need to zero out so as to not reuse it in Landau's logic
496f3237bfeSMark Adams   /* clean up */
497f3237bfeSMark Adams   PetscCall(DMPlexLandauDestroyVelocitySpace(&pack));
498f3237bfeSMark Adams   PetscCall(TSDestroy(&ts));
499f3237bfeSMark Adams   PetscCall(VecDestroy(&X));
500f3237bfeSMark Adams   PetscCall(PetscFinalize());
501f3237bfeSMark Adams   return 0;
502f3237bfeSMark Adams }
503f3237bfeSMark Adams 
504f3237bfeSMark Adams /*TEST
505f3237bfeSMark Adams 
506f3237bfeSMark Adams   build:
507f3237bfeSMark Adams     requires: !complex p4est
508f3237bfeSMark Adams 
509f3237bfeSMark Adams   testset:
510984ed092SMark Adams     requires: double defined(PETSC_USE_DMLANDAU_2D)
511f3237bfeSMark Adams     output_file: output/ex30_0.out
512f3237bfeSMark Adams     args: -dim 2 -petscspace_degree 3 -dm_landau_type p4est -dm_landau_num_species_grid 1,1,1 -dm_landau_amr_levels_max 0,0,0 \
513f3237bfeSMark Adams           -dm_landau_amr_post_refine 1 -number_particles_per_dimension 10 -dm_plex_hash_location \
514f3237bfeSMark Adams           -dm_landau_batch_size 2 -number_spatial_vertices 3 -dm_landau_batch_view_idx 1 -view_vertex_target 2 -view_grid_target 1 \
515f3237bfeSMark 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 \
516f3237bfeSMark Adams           -ftop_ksp_converged_reason -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 \
5176b664d00Smarkadams4           -ksp_type preonly -pc_type lu -dm_landau_verbose 4 \
518f3237bfeSMark Adams           -ptof_ksp_type cg -ptof_pc_type jacobi -ptof_ksp_converged_reason -ptof_ksp_rtol 1e-12\
519f3237bfeSMark Adams           -snes_converged_reason -snes_monitor -snes_rtol 1e-14 -snes_stol 1e-14\
520f3237bfeSMark 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 -info :vec
521f3237bfeSMark Adams 
522f3237bfeSMark Adams     test:
523f3237bfeSMark Adams       suffix: cpu
524f3237bfeSMark Adams       args: -dm_landau_device_type cpu
525f3237bfeSMark Adams     test:
526f3237bfeSMark Adams       suffix: kokkos
527f3237bfeSMark Adams       requires: kokkos_kernels
528f3237bfeSMark Adams       args: -dm_landau_device_type kokkos -dm_mat_type aijkokkos -dm_vec_type kokkos
529f3237bfeSMark Adams     test:
530f3237bfeSMark Adams       suffix: cuda
531f3237bfeSMark Adams       requires: cuda
532f3237bfeSMark Adams       args: -dm_landau_device_type cuda -dm_mat_type aijcusparse -dm_vec_type cuda
533f3237bfeSMark Adams 
534984ed092SMark Adams   testset:
535984ed092SMark Adams     requires: double !defined(PETSC_USE_DMLANDAU_2D)
536984ed092SMark Adams     output_file: output/ex30_3d.out
537984ed092SMark Adams     args: -dim 3 -petscspace_degree 2 -dm_landau_type p8est -dm_landau_num_species_grid 1,1,1 -dm_landau_amr_levels_max 0,0,0 \
538984ed092SMark Adams           -dm_landau_amr_post_refine 0 -number_particles_per_dimension 5 -dm_plex_hash_location \
539984ed092SMark Adams           -dm_landau_batch_size 1 -number_spatial_vertices 1 -dm_landau_batch_view_idx 0 -view_vertex_target 0 -view_grid_target 0 \
540984ed092SMark 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 \
541984ed092SMark Adams           -ftop_ksp_converged_reason -ftop_ksp_rtol 1e-12 -ftop_ksp_type cg -ftop_pc_type jacobi \
542984ed092SMark Adams           -ksp_type preonly -pc_type lu \
5436b664d00Smarkadams4           -ptof_ksp_type cg -ptof_pc_type jacobi -ptof_ksp_converged_reason -ptof_ksp_rtol 1e-6\
5446b664d00Smarkadams4           -snes_converged_reason -snes_monitor -snes_rtol 1e-9 -snes_stol 1e-9\
545984ed092SMark Adams           -ts_dt 0.1 -ts_exact_final_time stepover -ts_max_snes_failures -1 -ts_max_steps 1 -ts_monitor -ts_type beuler -info :vec
546984ed092SMark Adams 
547984ed092SMark Adams     test:
548984ed092SMark Adams       suffix: cpu_3d
549984ed092SMark Adams       args: -dm_landau_device_type cpu
550984ed092SMark Adams     test:
551984ed092SMark Adams       suffix: kokkos_3d
552984ed092SMark Adams       requires: kokkos_kernels
553984ed092SMark Adams       args: -dm_landau_device_type kokkos -dm_mat_type aijkokkos -dm_vec_type kokkos
554984ed092SMark Adams     test:
555984ed092SMark Adams       suffix: cuda_3d
556984ed092SMark Adams       requires: cuda
557984ed092SMark Adams       args: -dm_landau_device_type cuda -dm_mat_type aijcusparse -dm_vec_type cuda
558984ed092SMark Adams 
559f3237bfeSMark Adams TEST*/
560