xref: /petsc/src/ts/tests/ex30.c (revision 3cdbf90f99184cddde97da155a0864051da5c72d)
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 
31*3cdbf90fSMark Adams typedef struct {
32*3cdbf90fSMark Adams   PetscInt   v_target;
33*3cdbf90fSMark Adams   DM        *globSwarmArray;
34*3cdbf90fSMark Adams   LandauCtx *ctx;
35*3cdbf90fSMark Adams   PetscInt  *nTargetP;
36*3cdbf90fSMark Adams   PetscReal  N_inv;
37*3cdbf90fSMark Adams   DM        *grid_dm;
38*3cdbf90fSMark Adams   Mat       *g_Mass;
39*3cdbf90fSMark Adams   Mat       *globMpArray;
40*3cdbf90fSMark Adams   Vec       *globXArray;
41*3cdbf90fSMark Adams   PetscBool  print;
42*3cdbf90fSMark Adams } PrintCtx;
43*3cdbf90fSMark Adams 
44d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultMtM_SeqAIJ(Mat MtM, Vec xx, Vec yy)
45d71ae5a4SJacob Faibussowitsch {
46f3237bfeSMark Adams   MatShellCtx *matshellctx;
47f3237bfeSMark Adams 
48f3237bfeSMark Adams   PetscFunctionBeginUser;
49f3237bfeSMark Adams   PetscCall(MatShellGetContext(MtM, &matshellctx));
50f3237bfeSMark Adams   PetscCheck(matshellctx, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "No context");
51f3237bfeSMark Adams   PetscCall(MatMult(matshellctx->Mp, xx, matshellctx->ff));
52f3237bfeSMark Adams   PetscCall(MatMult(matshellctx->MpTrans, matshellctx->ff, yy));
533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
54f3237bfeSMark Adams }
55f3237bfeSMark Adams 
56d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultAddMtM_SeqAIJ(Mat MtM, Vec xx, Vec yy, Vec zz)
57d71ae5a4SJacob Faibussowitsch {
58f3237bfeSMark Adams   MatShellCtx *matshellctx;
59f3237bfeSMark Adams 
60f3237bfeSMark Adams   PetscFunctionBeginUser;
61f3237bfeSMark Adams   PetscCall(MatShellGetContext(MtM, &matshellctx));
62f3237bfeSMark Adams   PetscCheck(matshellctx, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "No context");
63f3237bfeSMark Adams   PetscCall(MatMult(matshellctx->Mp, xx, matshellctx->ff));
64f3237bfeSMark Adams   PetscCall(MatMultAdd(matshellctx->MpTrans, matshellctx->ff, yy, zz));
653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
66f3237bfeSMark Adams }
67f3237bfeSMark Adams 
68d71ae5a4SJacob Faibussowitsch PetscErrorCode createSwarm(const DM dm, PetscInt dim, DM *sw)
69d71ae5a4SJacob Faibussowitsch {
70f3237bfeSMark Adams   PetscInt Nc = 1;
71f3237bfeSMark Adams 
72f3237bfeSMark Adams   PetscFunctionBeginUser;
73f3237bfeSMark Adams   PetscCall(DMCreate(PETSC_COMM_SELF, sw));
74f3237bfeSMark Adams   PetscCall(DMSetType(*sw, DMSWARM));
75f3237bfeSMark Adams   PetscCall(DMSetDimension(*sw, dim));
76f3237bfeSMark Adams   PetscCall(DMSwarmSetType(*sw, DMSWARM_PIC));
77f3237bfeSMark Adams   PetscCall(DMSwarmSetCellDM(*sw, dm));
78*3cdbf90fSMark Adams   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "w_q", Nc, PETSC_REAL));
79f3237bfeSMark Adams   PetscCall(DMSwarmFinalizeFieldRegister(*sw));
80f3237bfeSMark Adams   PetscCall(DMSetFromOptions(*sw));
81*3cdbf90fSMark Adams   PetscCall(PetscObjectSetName((PetscObject)*sw, "Particle Grid"));
823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
83f3237bfeSMark Adams }
84f3237bfeSMark Adams 
85*3cdbf90fSMark Adams static PetscErrorCode makeSwarm(DM sw, const PetscInt dim, const PetscInt Np, const PetscReal xx[], const PetscReal yy[], const PetscReal zz[])
86*3cdbf90fSMark Adams {
87*3cdbf90fSMark Adams   PetscReal    *coords;
88*3cdbf90fSMark Adams   PetscDataType dtype;
89*3cdbf90fSMark Adams   PetscInt      bs, p, zero = 0;
90*3cdbf90fSMark Adams   PetscFunctionBeginUser;
91*3cdbf90fSMark Adams 
92*3cdbf90fSMark Adams   PetscCall(DMSwarmSetLocalSizes(sw, Np, zero));
93*3cdbf90fSMark Adams   PetscCall(DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
94*3cdbf90fSMark Adams   for (p = 0; p < Np; p++) {
95*3cdbf90fSMark Adams     coords[p * dim + 0] = xx[p];
96*3cdbf90fSMark Adams     coords[p * dim + 1] = yy[p];
97*3cdbf90fSMark Adams     if (dim == 3) coords[p * dim + 2] = zz[p];
98*3cdbf90fSMark Adams   }
99*3cdbf90fSMark Adams   PetscCall(DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
100*3cdbf90fSMark Adams   PetscFunctionReturn(PETSC_SUCCESS);
101*3cdbf90fSMark Adams }
102*3cdbf90fSMark Adams 
103*3cdbf90fSMark Adams static PetscErrorCode createMp(const DM dm, DM sw, Mat *Mp_out)
104*3cdbf90fSMark Adams {
105*3cdbf90fSMark Adams   PetscBool removePoints = PETSC_TRUE;
106*3cdbf90fSMark Adams   Mat       M_p;
107*3cdbf90fSMark Adams   PetscFunctionBeginUser;
108*3cdbf90fSMark Adams   // migrate after coords are set
109*3cdbf90fSMark Adams   PetscCall(DMSwarmMigrate(sw, removePoints));
110*3cdbf90fSMark Adams   PetscCall(PetscObjectSetName((PetscObject)sw, "Particle Grid"));
111*3cdbf90fSMark Adams   /* This gives M f = \int_\Omega \phi f, which looks like a rhs for a PDE */
112*3cdbf90fSMark Adams   PetscCall(DMCreateMassMatrix(sw, dm, &M_p));
113*3cdbf90fSMark Adams   PetscCall(DMViewFromOptions(sw, NULL, "-ex30_sw_view"));
114*3cdbf90fSMark Adams   // output
115*3cdbf90fSMark Adams   *Mp_out = M_p;
116*3cdbf90fSMark Adams   PetscFunctionReturn(PETSC_SUCCESS);
117*3cdbf90fSMark Adams }
118*3cdbf90fSMark Adams 
119*3cdbf90fSMark Adams static PetscErrorCode particlesToGrid(const DM dm, DM sw, const PetscInt Np, const PetscInt a_tid, const PetscInt dim, const PetscReal a_wp[], Vec rho, Mat M_p)
120*3cdbf90fSMark Adams {
121*3cdbf90fSMark Adams   PetscReal    *wq;
122*3cdbf90fSMark Adams   PetscDataType dtype;
123*3cdbf90fSMark Adams   Vec           ff;
124*3cdbf90fSMark Adams   PetscInt      bs, p;
125*3cdbf90fSMark Adams 
126*3cdbf90fSMark Adams   PetscFunctionBeginUser;
127*3cdbf90fSMark Adams   PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void **)&wq));
128*3cdbf90fSMark Adams   for (p = 0; p < Np; p++) wq[p] = a_wp[p];
129*3cdbf90fSMark Adams   PetscCall(DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void **)&wq));
130*3cdbf90fSMark Adams   PetscCall(PetscObjectSetName((PetscObject)rho, "rho"));
131*3cdbf90fSMark Adams   PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &ff));
132*3cdbf90fSMark Adams   PetscCall(PetscObjectSetName((PetscObject)ff, "weights"));
133*3cdbf90fSMark Adams   PetscCall(MatMultTranspose(M_p, ff, rho));
134*3cdbf90fSMark Adams   PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &ff));
135*3cdbf90fSMark Adams   PetscFunctionReturn(PETSC_SUCCESS);
136*3cdbf90fSMark Adams }
137*3cdbf90fSMark Adams 
138*3cdbf90fSMark Adams //
139*3cdbf90fSMark Adams // add grid to arg 'sw.w_q'
140*3cdbf90fSMark Adams //
141*3cdbf90fSMark Adams PetscErrorCode gridToParticles(const DM dm, DM sw, const Vec rhs, Vec work, Mat M_p, Mat Mass)
142d71ae5a4SJacob Faibussowitsch {
143f3237bfeSMark Adams   PetscBool    is_lsqr;
144f3237bfeSMark Adams   KSP          ksp;
145f3237bfeSMark Adams   Mat          PM_p = NULL, MtM, D;
146f3237bfeSMark Adams   Vec          ff;
147f3237bfeSMark Adams   PetscInt     N, M, nzl;
148f3237bfeSMark Adams   MatShellCtx *matshellctx;
149f3237bfeSMark Adams 
150f3237bfeSMark Adams   PetscFunctionBeginUser;
151f3237bfeSMark Adams   PetscCall(MatMult(Mass, rhs, work));
152f3237bfeSMark Adams   // pseudo-inverse
153f3237bfeSMark Adams   PetscCall(KSPCreate(PETSC_COMM_SELF, &ksp));
154f3237bfeSMark Adams   PetscCall(KSPSetOptionsPrefix(ksp, "ftop_"));
155f3237bfeSMark Adams   PetscCall(KSPSetFromOptions(ksp));
156f3237bfeSMark Adams   PetscCall(PetscObjectTypeCompare((PetscObject)ksp, KSPLSQR, &is_lsqr));
157f3237bfeSMark Adams   if (!is_lsqr) {
158f3237bfeSMark Adams     PetscCall(MatGetLocalSize(M_p, &M, &N));
159f3237bfeSMark Adams     if (N > M) {
160f3237bfeSMark Adams       PC pc;
161f3237bfeSMark Adams       PetscCall(PetscInfo(ksp, " M (%" PetscInt_FMT ") < M (%" PetscInt_FMT ") -- skip revert to lsqr\n", M, N));
162f3237bfeSMark Adams       is_lsqr = PETSC_TRUE;
163f3237bfeSMark Adams       PetscCall(KSPSetType(ksp, KSPLSQR));
164f3237bfeSMark Adams       PetscCall(KSPGetPC(ksp, &pc));
165f3237bfeSMark 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
166f3237bfeSMark Adams     } else {
167f3237bfeSMark Adams       PetscCall(PetscNew(&matshellctx));
168f3237bfeSMark Adams       PetscCall(MatCreateShell(PetscObjectComm((PetscObject)dm), N, N, PETSC_DECIDE, PETSC_DECIDE, matshellctx, &MtM));
169f3237bfeSMark Adams       PetscCall(MatTranspose(M_p, MAT_INITIAL_MATRIX, &matshellctx->MpTrans));
170f3237bfeSMark Adams       matshellctx->Mp = M_p;
171f3237bfeSMark Adams       PetscCall(MatShellSetOperation(MtM, MATOP_MULT, (void (*)(void))MatMultMtM_SeqAIJ));
172f3237bfeSMark Adams       PetscCall(MatShellSetOperation(MtM, MATOP_MULT_ADD, (void (*)(void))MatMultAddMtM_SeqAIJ));
173f3237bfeSMark Adams       PetscCall(MatCreateVecs(M_p, &matshellctx->uu, &matshellctx->ff));
174f3237bfeSMark Adams       PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, N, N, 1, NULL, &D));
175f3237bfeSMark Adams       PetscCall(MatViewFromOptions(matshellctx->MpTrans, NULL, "-ftop2_Mp_mat_view"));
176f3237bfeSMark Adams       for (int i = 0; i < N; i++) {
177f3237bfeSMark Adams         const PetscScalar *vals;
178f3237bfeSMark Adams         const PetscInt    *cols;
179f3237bfeSMark Adams         PetscScalar        dot = 0;
180f3237bfeSMark Adams         PetscCall(MatGetRow(matshellctx->MpTrans, i, &nzl, &cols, &vals));
181f3237bfeSMark Adams         for (int ii = 0; ii < nzl; ii++) dot += PetscSqr(vals[ii]);
18263a3b9bcSJacob Faibussowitsch         PetscCheck(dot != 0.0, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Row %d is empty", i);
183f3237bfeSMark Adams         PetscCall(MatSetValue(D, i, i, dot, INSERT_VALUES));
184f3237bfeSMark Adams       }
185f3237bfeSMark Adams       PetscCall(MatAssemblyBegin(D, MAT_FINAL_ASSEMBLY));
186f3237bfeSMark Adams       PetscCall(MatAssemblyEnd(D, MAT_FINAL_ASSEMBLY));
187f3237bfeSMark Adams       PetscCall(PetscInfo(M_p, "createMtMKSP Have %" PetscInt_FMT " eqs, nzl = %" PetscInt_FMT "\n", N, nzl));
188f3237bfeSMark Adams       PetscCall(KSPSetOperators(ksp, MtM, D));
189f3237bfeSMark Adams       PetscCall(MatViewFromOptions(D, NULL, "-ftop2_D_mat_view"));
190f3237bfeSMark Adams       PetscCall(MatViewFromOptions(M_p, NULL, "-ftop2_Mp_mat_view"));
191f3237bfeSMark Adams       PetscCall(MatViewFromOptions(matshellctx->MpTrans, NULL, "-ftop2_MpTranspose_mat_view"));
1926b664d00Smarkadams4       PetscCall(MatViewFromOptions(MtM, NULL, "-ftop2_MtM_mat_view"));
193f3237bfeSMark Adams     }
194f3237bfeSMark Adams   }
195f3237bfeSMark Adams   if (is_lsqr) {
196f3237bfeSMark Adams     PC        pc;
197f3237bfeSMark Adams     PetscBool is_bjac;
198f3237bfeSMark Adams     PetscCall(KSPGetPC(ksp, &pc));
199f3237bfeSMark Adams     PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCBJACOBI, &is_bjac));
200f3237bfeSMark Adams     if (is_bjac) {
201f3237bfeSMark Adams       PetscCall(DMSwarmCreateMassMatrixSquare(sw, dm, &PM_p));
202f3237bfeSMark Adams       PetscCall(KSPSetOperators(ksp, M_p, PM_p));
203f3237bfeSMark Adams     } else {
204f3237bfeSMark Adams       PetscCall(KSPSetOperators(ksp, M_p, M_p));
205f3237bfeSMark Adams     }
206f3237bfeSMark Adams   }
207f3237bfeSMark Adams   PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &ff)); // this grabs access
208f3237bfeSMark Adams   if (!is_lsqr) {
209*3cdbf90fSMark Adams     PetscCall(KSPSolve(ksp, work, matshellctx->uu));
210f3237bfeSMark Adams     PetscCall(MatMult(M_p, matshellctx->uu, ff));
211f3237bfeSMark Adams     PetscCall(MatDestroy(&matshellctx->MpTrans));
212f3237bfeSMark Adams     PetscCall(VecDestroy(&matshellctx->ff));
213f3237bfeSMark Adams     PetscCall(VecDestroy(&matshellctx->uu));
214f3237bfeSMark Adams     PetscCall(MatDestroy(&D));
215f3237bfeSMark Adams     PetscCall(MatDestroy(&MtM));
216f3237bfeSMark Adams     PetscCall(PetscFree(matshellctx));
217f3237bfeSMark Adams   } else {
218*3cdbf90fSMark Adams     PetscCall(KSPSolveTranspose(ksp, work, ff));
219f3237bfeSMark Adams   }
220f3237bfeSMark Adams   PetscCall(KSPDestroy(&ksp));
221f3237bfeSMark Adams   /* Visualize particle field */
222f3237bfeSMark Adams   PetscCall(VecViewFromOptions(ff, NULL, "-weights_view"));
223f3237bfeSMark Adams   PetscCall(MatDestroy(&PM_p));
224f3237bfeSMark Adams   PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &ff));
225f3237bfeSMark Adams 
2263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
227f3237bfeSMark Adams }
228f3237bfeSMark Adams 
229*3cdbf90fSMark Adams #define EX30_MAX_NUM_THRDS 12
230*3cdbf90fSMark Adams #define EX30_MAX_BATCH_SZ  1024
231*3cdbf90fSMark Adams //
232*3cdbf90fSMark Adams // add grid to arg 'globSwarmArray[].w_q'
233*3cdbf90fSMark Adams //
234*3cdbf90fSMark 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)
235d71ae5a4SJacob Faibussowitsch {
236*3cdbf90fSMark Adams   PetscErrorCode ierr = (PetscErrorCode)0; // used for inside thread loops
237f3237bfeSMark Adams 
238f3237bfeSMark Adams   PetscFunctionBeginUser;
239*3cdbf90fSMark Adams   // map back to particles
240*3cdbf90fSMark Adams   for (PetscInt v_id_0 = 0; v_id_0 < ctx->batch_sz; v_id_0 += numthreads) {
241*3cdbf90fSMark 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));
242*3cdbf90fSMark Adams     //PetscPragmaOMP(parallel for)
243*3cdbf90fSMark Adams     for (int tid = 0; tid < numthreads; tid++) {
244*3cdbf90fSMark Adams       const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id + v_id;
245*3cdbf90fSMark Adams       if (glb_v_id < num_vertices) {
246*3cdbf90fSMark Adams         for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids
247*3cdbf90fSMark Adams           PetscErrorCode ierr_t;
248*3cdbf90fSMark 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));
249*3cdbf90fSMark 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]);
250*3cdbf90fSMark Adams           if (ierr_t) ierr = ierr_t;
251f3237bfeSMark Adams         }
252*3cdbf90fSMark Adams       }
253*3cdbf90fSMark Adams     }
254*3cdbf90fSMark Adams     PetscCheck(!ierr, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Error in OMP loop. ierr = %d", (int)ierr);
255*3cdbf90fSMark Adams     /* Get moments */
256*3cdbf90fSMark Adams     PetscCall(PetscInfo(grid_dm[0], "Cleanup batches %" PetscInt_FMT " to %" PetscInt_FMT "\n", v_id_0, v_id_0 + numthreads));
257*3cdbf90fSMark Adams     for (int tid = 0; tid < numthreads; tid++) {
258*3cdbf90fSMark Adams       const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id + v_id;
259*3cdbf90fSMark Adams       if (glb_v_id == v_target) {
260*3cdbf90fSMark Adams         for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
261*3cdbf90fSMark Adams           PetscDataType dtype;
262*3cdbf90fSMark Adams           PetscReal    *wp, *coords;
263*3cdbf90fSMark Adams           DM            sw = globSwarmArray[LAND_PACK_IDX(v_id, grid)];
264*3cdbf90fSMark Adams           PetscInt      npoints, bs = 1;
265*3cdbf90fSMark Adams           PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void **)&wp)); // take data out here
266*3cdbf90fSMark Adams           PetscCall(DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
267*3cdbf90fSMark Adams           PetscCall(DMSwarmGetLocalSize(sw, &npoints));
268*3cdbf90fSMark Adams           for (int p = 0; p < npoints; p++) {
269*3cdbf90fSMark 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]];
270*3cdbf90fSMark Adams             for (int i = 0; i < dim; ++i) v2 += PetscSqr(coords[p * dim + i]);
271*3cdbf90fSMark Adams             moments[0] += w;
272*3cdbf90fSMark Adams             moments[1] += w * ctx->v_0 * coords[p * dim + 1]; // z-momentum
273*3cdbf90fSMark Adams             moments[2] += w * ctx->v_0 * ctx->v_0 * v2;
274*3cdbf90fSMark Adams           }
275*3cdbf90fSMark Adams           PetscCall(DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void **)&wp));
276f3237bfeSMark Adams           PetscCall(DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
277*3cdbf90fSMark Adams         }
278*3cdbf90fSMark Adams         const PetscReal N_inv = 1 / moments[0];
279*3cdbf90fSMark Adams         for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
280*3cdbf90fSMark Adams           PetscDataType dtype;
281*3cdbf90fSMark Adams           PetscReal    *wp, *coords;
282*3cdbf90fSMark Adams           DM            sw = globSwarmArray[LAND_PACK_IDX(v_id, grid)];
283*3cdbf90fSMark Adams           PetscInt      npoints, bs = 1;
284*3cdbf90fSMark Adams           PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void **)&wp)); // take data out here
285*3cdbf90fSMark Adams           PetscCall(DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
286*3cdbf90fSMark Adams           PetscCall(DMSwarmGetLocalSize(sw, &npoints));
287*3cdbf90fSMark Adams           for (int p = 0; p < npoints; p++) {
288*3cdbf90fSMark 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;
289*3cdbf90fSMark Adams             if (ww > PETSC_REAL_MIN) {
290*3cdbf90fSMark Adams               moments[3] -= ww * PetscLogReal(ww);
291*3cdbf90fSMark Adams               PetscCheck(ww < 1 - PETSC_MACHINE_EPSILON, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "ww (%g) > 1", (double)ww);
292*3cdbf90fSMark Adams             }
293*3cdbf90fSMark Adams           }
294*3cdbf90fSMark Adams           PetscCall(DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void **)&wp));
295*3cdbf90fSMark Adams           PetscCall(DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
296*3cdbf90fSMark Adams         }
297*3cdbf90fSMark Adams       }
298*3cdbf90fSMark Adams     } // thread batch
299*3cdbf90fSMark Adams   }   // batch
3003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
301f3237bfeSMark Adams }
302*3cdbf90fSMark Adams 
303*3cdbf90fSMark Adams static void maxwellian(PetscInt dim, const PetscReal x[], PetscReal kt_m, PetscReal n, PetscReal shift, PetscScalar *u)
304d71ae5a4SJacob Faibussowitsch {
305f3237bfeSMark Adams   PetscInt  i;
306f3237bfeSMark Adams   PetscReal v2 = 0, theta = 2.0 * kt_m; /* theta = 2kT/mc^2 */
307f3237bfeSMark Adams 
308f3237bfeSMark Adams   /* compute the exponents, v^2 */
309f3237bfeSMark Adams   for (i = 0; i < dim; ++i) v2 += x[i] * x[i];
310f3237bfeSMark Adams   /* evaluate the Maxwellian */
311f3237bfeSMark Adams   u[0] = n * PetscPowReal(PETSC_PI * theta, -1.5) * (PetscExpReal(-v2 / theta));
312*3cdbf90fSMark Adams   if (shift != 0.) {
313*3cdbf90fSMark Adams     v2 = 0;
314*3cdbf90fSMark Adams     for (i = 0; i < dim - 1; ++i) v2 += x[i] * x[i];
315*3cdbf90fSMark Adams     v2 += (x[dim - 1] - shift) * (x[dim - 1] - shift);
316*3cdbf90fSMark Adams     /* evaluate the shifted Maxwellian */
317*3cdbf90fSMark Adams     u[0] += n * PetscPowReal(PETSC_PI * theta, -1.5) * (PetscExpReal(-v2 / theta));
318*3cdbf90fSMark Adams   }
319f3237bfeSMark Adams }
320f3237bfeSMark Adams 
321*3cdbf90fSMark Adams static PetscErrorCode PostStep(TS ts)
322*3cdbf90fSMark Adams {
323*3cdbf90fSMark Adams   PetscInt   n, dim, nDMs;
324*3cdbf90fSMark Adams   PetscReal  t;
325*3cdbf90fSMark Adams   LandauCtx *ctx;
326*3cdbf90fSMark Adams   Vec        X;
327*3cdbf90fSMark Adams   PrintCtx  *printCtx;
328*3cdbf90fSMark Adams   PetscReal  moments[4];
329*3cdbf90fSMark Adams   DM         pack;
330*3cdbf90fSMark Adams 
331*3cdbf90fSMark Adams   PetscFunctionBeginUser;
332*3cdbf90fSMark Adams   PetscCall(TSGetApplicationContext(ts, &printCtx));
333*3cdbf90fSMark Adams   if (!printCtx->print) PetscFunctionReturn(PETSC_SUCCESS);
334*3cdbf90fSMark Adams 
335*3cdbf90fSMark Adams   for (int i = 0; i < 4; i++) moments[i] = 0;
336*3cdbf90fSMark Adams   ctx = printCtx->ctx;
337*3cdbf90fSMark Adams   PetscCall(TSGetDM(ts, &pack));
338*3cdbf90fSMark Adams   PetscCall(DMGetDimension(pack, &dim));
339*3cdbf90fSMark Adams   PetscCall(DMCompositeGetNumberDM(pack, &nDMs)); // number of vertices * number of grids
340*3cdbf90fSMark Adams   PetscCall(TSGetSolution(ts, &X));
341*3cdbf90fSMark Adams   const PetscInt v_id = printCtx->v_target % ctx->batch_sz;
342*3cdbf90fSMark Adams   PetscCall(TSGetSolution(ts, &X));
343*3cdbf90fSMark Adams   PetscCall(TSGetStepNumber(ts, &n));
344*3cdbf90fSMark Adams   PetscCall(TSGetTime(ts, &t));
345*3cdbf90fSMark Adams   PetscCall(DMCompositeGetAccessArray(pack, X, nDMs, NULL, printCtx->globXArray));
346*3cdbf90fSMark Adams   for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
347*3cdbf90fSMark Adams     PetscDataType dtype;
348*3cdbf90fSMark Adams     PetscReal    *wp, *coords;
349*3cdbf90fSMark Adams     DM            sw = printCtx->globSwarmArray[LAND_PACK_IDX(v_id, grid)];
350*3cdbf90fSMark Adams     Vec           work, subX = printCtx->globXArray[LAND_PACK_IDX(v_id, grid)];
351*3cdbf90fSMark Adams     PetscInt      bs, NN     = printCtx->nTargetP[grid];
352*3cdbf90fSMark Adams     // C-G moments
353*3cdbf90fSMark Adams     PetscCall(VecDuplicate(subX, &work));
354*3cdbf90fSMark Adams     PetscCall(gridToParticles(printCtx->grid_dm[grid], sw, subX, work, printCtx->globMpArray[LAND_PACK_IDX(v_id, grid)], printCtx->g_Mass[grid]));
355*3cdbf90fSMark Adams     PetscCall(VecDestroy(&work));
356*3cdbf90fSMark Adams     // moments
357*3cdbf90fSMark Adams     PetscCall(DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
358*3cdbf90fSMark Adams     PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void **)&wp)); // could get NN from sw - todo
359*3cdbf90fSMark Adams     for (int pp = 0; pp < NN; pp++) {
360*3cdbf90fSMark 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]], ww = w * printCtx->N_inv;
361*3cdbf90fSMark Adams       for (int i = 0; i < dim; ++i) v2 += PetscSqr(coords[pp * dim + i]);
362*3cdbf90fSMark Adams       moments[0] += w;
363*3cdbf90fSMark Adams       moments[1] += w * ctx->v_0 * coords[pp * dim + 1]; // z-momentum
364*3cdbf90fSMark Adams       moments[2] += w * ctx->v_0 * ctx->v_0 * v2;
365*3cdbf90fSMark Adams       if (ww > PETSC_REAL_MIN) {
366*3cdbf90fSMark Adams         moments[3] -= ww * PetscLogReal(ww);
367*3cdbf90fSMark Adams         PetscCheck(ww < 1 - PETSC_MACHINE_EPSILON, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "ww (%g) > 1", (double)ww);
368*3cdbf90fSMark Adams       }
369*3cdbf90fSMark Adams     }
370*3cdbf90fSMark Adams     PetscCall(DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
371*3cdbf90fSMark Adams     PetscCall(DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void **)&wp));
372*3cdbf90fSMark Adams   }
373*3cdbf90fSMark Adams   PetscCall(DMCompositeRestoreAccessArray(pack, X, nDMs, NULL, printCtx->globXArray));
374*3cdbf90fSMark Adams   PetscCall(PetscInfo(X, "%4d) time %e, Landau moments: %18.12e %19.12e %18.12e %e\n", (int)n, (double)t, (double)moments[0], (double)moments[1], (double)moments[2], (double)moments[3]));
375*3cdbf90fSMark Adams   PetscFunctionReturn(PETSC_SUCCESS);
376*3cdbf90fSMark Adams }
377*3cdbf90fSMark Adams 
378*3cdbf90fSMark 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)
379d71ae5a4SJacob Faibussowitsch {
380f3237bfeSMark Adams   DM             pack, *globSwarmArray, grid_dm[LANDAU_MAX_GRIDS];
381f3237bfeSMark Adams   Mat           *globMpArray, g_Mass[LANDAU_MAX_GRIDS];
382*3cdbf90fSMark Adams   KSP            t_ksp[LANDAU_MAX_GRIDS][EX30_MAX_NUM_THRDS];
383*3cdbf90fSMark Adams   Vec            t_fhat[LANDAU_MAX_GRIDS][EX30_MAX_NUM_THRDS];
384*3cdbf90fSMark Adams   PetscInt       nDMs, nTargetP[LANDAU_MAX_GRIDS];
385*3cdbf90fSMark Adams   PetscErrorCode ierr = (PetscErrorCode)0; // used for inside thread loops
386f3237bfeSMark Adams #if defined(PETSC_HAVE_OPENMP) && defined(PETSC_HAVE_THREADSAFETY)
387f3237bfeSMark Adams   PetscInt numthreads = PetscNumOMPThreads;
388f3237bfeSMark Adams #else
389f3237bfeSMark Adams   PetscInt numthreads = 1;
390f3237bfeSMark Adams #endif
391f3237bfeSMark Adams   LandauCtx *ctx;
392f3237bfeSMark Adams   Vec       *globXArray;
393*3cdbf90fSMark Adams   PetscReal  moments_0[4], moments_1a[4], moments_1b[4], dt_init;
394*3cdbf90fSMark Adams   PrintCtx  *printCtx;
395f3237bfeSMark Adams 
396f3237bfeSMark Adams   PetscFunctionBeginUser;
397*3cdbf90fSMark 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);
398*3cdbf90fSMark Adams   PetscCheck(numthreads > 0, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Number threads %" PetscInt_FMT " > %d", numthreads, EX30_MAX_NUM_THRDS);
399f3237bfeSMark Adams   PetscCall(TSGetDM(ts, &pack));
400f3237bfeSMark Adams   PetscCall(DMGetApplicationContext(pack, &ctx));
401f3237bfeSMark 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);
402*3cdbf90fSMark Adams   PetscCall(DMCompositeGetNumberDM(pack, &nDMs)); // number of vertices * number of grids
403*3cdbf90fSMark 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, num_vertices));
404f3237bfeSMark Adams   PetscCall(PetscMalloc(sizeof(*globXArray) * nDMs, &globXArray));
405f3237bfeSMark Adams   PetscCall(PetscMalloc(sizeof(*globMpArray) * nDMs, &globMpArray));
406f3237bfeSMark Adams   PetscCall(PetscMalloc(sizeof(*globSwarmArray) * nDMs, &globSwarmArray));
407*3cdbf90fSMark Adams   // print ctx
408*3cdbf90fSMark Adams   PetscCall(PetscNew(&printCtx));
409*3cdbf90fSMark Adams   PetscCall(TSSetApplicationContext(ts, printCtx));
410*3cdbf90fSMark Adams   printCtx->v_target       = v_target;
411*3cdbf90fSMark Adams   printCtx->ctx            = ctx;
412*3cdbf90fSMark Adams   printCtx->nTargetP       = nTargetP;
413*3cdbf90fSMark Adams   printCtx->globSwarmArray = globSwarmArray;
414*3cdbf90fSMark Adams   printCtx->grid_dm        = grid_dm;
415*3cdbf90fSMark Adams   printCtx->globMpArray    = globMpArray;
416*3cdbf90fSMark Adams   printCtx->g_Mass         = g_Mass;
417*3cdbf90fSMark Adams   printCtx->globXArray     = globXArray;
418*3cdbf90fSMark Adams   // view
419f3237bfeSMark Adams   PetscCall(DMViewFromOptions(ctx->plex[g_target], NULL, "-ex30_dm_view"));
420*3cdbf90fSMark Adams   // create mesh mass matrices
4216b664d00Smarkadams4   PetscCall(VecZeroEntries(X));
422f3237bfeSMark Adams   PetscCall(DMCompositeGetAccessArray(pack, X, nDMs, NULL, globXArray)); // just to duplicate
423f3237bfeSMark Adams   for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {               // add same particels for all grids
424f3237bfeSMark Adams     Vec          subX = globXArray[LAND_PACK_IDX(0, grid)];
425f3237bfeSMark Adams     DM           dm   = ctx->plex[grid];
426f3237bfeSMark Adams     PetscSection s;
427f3237bfeSMark Adams     grid_dm[grid] = dm;
428f3237bfeSMark Adams     PetscCall(DMCreateMassMatrix(dm, dm, &g_Mass[grid]));
429f3237bfeSMark Adams     //
430f3237bfeSMark Adams     PetscCall(DMGetLocalSection(dm, &s));
431f3237bfeSMark Adams     PetscCall(DMPlexCreateClosureIndex(dm, s));
432f3237bfeSMark Adams     for (int tid = 0; tid < numthreads; tid++) {
433f3237bfeSMark Adams       PetscCall(VecDuplicate(subX, &t_fhat[grid][tid]));
434f3237bfeSMark Adams       PetscCall(KSPCreate(PETSC_COMM_SELF, &t_ksp[grid][tid]));
435f3237bfeSMark Adams       PetscCall(KSPSetOptionsPrefix(t_ksp[grid][tid], "ptof_"));
436f3237bfeSMark Adams       PetscCall(KSPSetOperators(t_ksp[grid][tid], g_Mass[grid], g_Mass[grid]));
437f3237bfeSMark Adams       PetscCall(KSPSetFromOptions(t_ksp[grid][tid]));
438f3237bfeSMark Adams     }
439f3237bfeSMark Adams   }
440f3237bfeSMark Adams   PetscCall(DMCompositeRestoreAccessArray(pack, X, nDMs, NULL, globXArray));
441*3cdbf90fSMark Adams   for (int i = 0; i < 4; i++) moments_0[i] = moments_1a[i] = moments_1b[i] = 0;
442f3237bfeSMark Adams   PetscCall(TSGetTimeStep(ts, &dt_init)); // we could have an adaptive time stepper
443*3cdbf90fSMark Adams   // loop over all vertices in chucks that are batched for TSSolve
444*3cdbf90fSMark 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
4456b664d00Smarkadams4     PetscCall(TSSetTime(ts, 0));
4466b664d00Smarkadams4     PetscCall(TSSetStepNumber(ts, 0));
4476b664d00Smarkadams4     PetscCall(TSSetTimeStep(ts, dt_init));
448f3237bfeSMark Adams     PetscCall(DMCompositeGetAccessArray(pack, X, nDMs, NULL, globXArray));
449*3cdbf90fSMark Adams     if (v_target >= global_vertex_id_0 && v_target < global_vertex_id_0 + ctx->batch_sz) {
450*3cdbf90fSMark Adams       PetscCall(PetscObjectSetName((PetscObject)globXArray[LAND_PACK_IDX(v_target % ctx->batch_sz, g_target)], "rho"));
451*3cdbf90fSMark Adams       printCtx->print = PETSC_TRUE;
452*3cdbf90fSMark Adams     } else printCtx->print = PETSC_FALSE;
453*3cdbf90fSMark Adams     // create fake particles in batches with threads
454*3cdbf90fSMark Adams     for (PetscInt v_id_0 = 0; v_id_0 < ctx->batch_sz; v_id_0 += numthreads) {
455*3cdbf90fSMark 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];
456*3cdbf90fSMark Adams       PetscInt   Np_t[LANDAU_MAX_GRIDS][EX30_MAX_NUM_THRDS];
457f3237bfeSMark Adams       // make particles
458f3237bfeSMark Adams       for (int tid = 0; tid < numthreads; tid++) {
459*3cdbf90fSMark Adams         const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id;
460*3cdbf90fSMark Adams         if (glb_v_id < num_vertices) {                                                                                                                                            // the ragged edge (in last batch)
461*3cdbf90fSMark Adams           PetscInt Npp0 = a_Np + (glb_v_id % (a_Np / 10 + 1)), NN;                                                                                                                // number of particels in each dimension with add some load imbalance
462f3237bfeSMark Adams           for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {                                                                                                                // add same particels for all grids
463*3cdbf90fSMark 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 */
464f3237bfeSMark 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
465f3237bfeSMark Adams             PetscInt        Npi = Npp0, Npj = 2 * Npp0, Npk = 1;
466f3237bfeSMark Adams             if (dim == 2) lo[0] = 0; // Landau coordinate (r,z)
467f3237bfeSMark Adams             else Npi = Npj = Npk = Npp0;
468*3cdbf90fSMark Adams             // User: use glb_v_id to index into your data
469f3237bfeSMark Adams             NN              = Npi * Npj * Npk; // make a regular grid of particles Npp x Npp
470f3237bfeSMark Adams             Np_t[grid][tid] = NN;
471*3cdbf90fSMark Adams             if (glb_v_id == v_target) nTargetP[grid] = NN;
472f3237bfeSMark 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]));
473f3237bfeSMark Adams             hp[0] = (hi[0] - lo[0]) / Npi;
474f3237bfeSMark Adams             hp[1] = (hi[1] - lo[1]) / Npj;
475f3237bfeSMark Adams             hp[2] = (hi[2] - lo[2]) / Npk;
476f3237bfeSMark Adams             if (dim == 2) hp[2] = 1;
47763a3b9bcSJacob 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
478f3237bfeSMark Adams             vole = hp[0] * hp[1] * hp[2] * ctx->n[grid];                                                                                                                           // fix for multi-species
479*3cdbf90fSMark 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));
480f3237bfeSMark Adams             for (int pj = 0, pp = 0; pj < Npj; pj++) {
481f3237bfeSMark Adams               for (int pk = 0; pk < Npk; pk++) {
482f3237bfeSMark Adams                 for (int pi = 0; pi < Npi; pi++, pp++) {
483f3237bfeSMark Adams                   xx_t[grid][tid][pp] = lo[0] + hp[0] / 2.0 + pi * hp[0];
484f3237bfeSMark Adams                   yy_t[grid][tid][pp] = lo[1] + hp[1] / 2.0 + pj * hp[1];
485f3237bfeSMark Adams                   if (dim == 3) zz_t[grid][tid][pp] = lo[2] + hp[2] / 2.0 + pk * hp[2];
486f3237bfeSMark Adams                   {
487f3237bfeSMark Adams                     PetscReal x[] = {xx_t[grid][tid][pp], yy_t[grid][tid][pp], dim == 2 ? 0 : zz_t[grid][tid][pp]};
488*3cdbf90fSMark Adams                     maxwellian(dim, x, kT_m, vole, shift, &wp_t[grid][tid][pp]);
489f3237bfeSMark 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
490*3cdbf90fSMark Adams                     if (glb_v_id == v_target) {
491*3cdbf90fSMark 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]];
492f3237bfeSMark Adams                       for (int i = 0; i < dim; ++i) v2 += PetscSqr(x[i]);
493*3cdbf90fSMark Adams                       moments_0[0] += w;                   // not thread safe
494*3cdbf90fSMark Adams                       moments_0[1] += w * ctx->v_0 * x[1]; // z-momentum
495*3cdbf90fSMark Adams                       moments_0[2] += w * ctx->v_0 * ctx->v_0 * v2;
496f3237bfeSMark Adams                     }
497f3237bfeSMark Adams                   }
498f3237bfeSMark Adams                 }
499f3237bfeSMark Adams               }
500f3237bfeSMark Adams             }
501*3cdbf90fSMark Adams           }
502*3cdbf90fSMark Adams           // entropy
503*3cdbf90fSMark Adams           if (glb_v_id == v_target) {
504*3cdbf90fSMark Adams             printCtx->N_inv = 1 / moments_0[0];
505*3cdbf90fSMark Adams             PetscCall(PetscInfo(pack, "Target %" PetscInt_FMT " with %" PetscInt_FMT " particels\n", glb_v_id, nTargetP[0]));
506*3cdbf90fSMark Adams             for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
507*3cdbf90fSMark Adams               NN = nTargetP[grid];
508*3cdbf90fSMark Adams               for (int pp = 0; pp < NN; pp++) {
509*3cdbf90fSMark 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 * printCtx->N_inv;
510*3cdbf90fSMark Adams                 if (ww > PETSC_REAL_MIN) {
511*3cdbf90fSMark Adams                   moments_0[3] -= ww * PetscLogReal(ww);
512*3cdbf90fSMark Adams                   PetscCheck(ww < 1 - PETSC_MACHINE_EPSILON, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "ww (%g) > 1", (double)ww);
513*3cdbf90fSMark Adams                 }
514*3cdbf90fSMark Adams               }
515*3cdbf90fSMark Adams             } // diagnostics
516f3237bfeSMark Adams           }   // grid
517f3237bfeSMark Adams         }     // active
518*3cdbf90fSMark Adams       }       // threads
519f3237bfeSMark Adams       /* Create particle swarm */
520f37bacd1SJacob Faibussowitsch       for (int tid = 0; tid < numthreads; tid++) {
521*3cdbf90fSMark Adams         const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id;
522*3cdbf90fSMark Adams         if (glb_v_id < num_vertices) {                             // the ragged edge of the last batch
523f3237bfeSMark Adams           for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids
524f3237bfeSMark Adams             PetscErrorCode ierr_t;
525f3237bfeSMark Adams             PetscSection   section;
526f3237bfeSMark Adams             PetscInt       Nf;
527f3237bfeSMark Adams             DM             dm = grid_dm[grid];
528f3237bfeSMark Adams             ierr_t            = DMGetLocalSection(dm, &section);
529f3237bfeSMark Adams             ierr_t            = PetscSectionGetNumFields(section, &Nf);
530*3cdbf90fSMark Adams             if (Nf != 1) ierr_t = (PetscErrorCode)9999;
531f3237bfeSMark Adams             else {
532f3237bfeSMark Adams               ierr_t = DMViewFromOptions(dm, NULL, "-dm_view");
533*3cdbf90fSMark Adams               ierr_t = PetscInfo(pack, "call createSwarm [%" PetscInt_FMT ".%" PetscInt_FMT "] local batch index %" PetscInt_FMT "\n", v_id, grid, LAND_PACK_IDX(v_id, grid));
534*3cdbf90fSMark Adams               ierr_t = createSwarm(dm, dim, &globSwarmArray[LAND_PACK_IDX(v_id, grid)]);
535f3237bfeSMark Adams             }
536f3237bfeSMark Adams             if (ierr_t) ierr = ierr_t;
537f3237bfeSMark Adams           }
538f3237bfeSMark Adams         } // active
539*3cdbf90fSMark Adams       }   // threads
540cad9d221SBarry Smith       PetscCheck(ierr != 9999, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Only support one species per grid");
5416b664d00Smarkadams4       PetscCheck(!ierr, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Error in OMP loop. ierr = %d", (int)ierr);
542*3cdbf90fSMark Adams       // make globMpArray
543f3237bfeSMark Adams       PetscPragmaOMP(parallel for)
544f37bacd1SJacob Faibussowitsch       for (int tid = 0; tid < numthreads; tid++) {
545*3cdbf90fSMark Adams         const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id;
546*3cdbf90fSMark Adams         if (glb_v_id < num_vertices) {
547*3cdbf90fSMark Adams           for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids
548*3cdbf90fSMark Adams             PetscErrorCode ierr_t;
549*3cdbf90fSMark Adams             DM             sw = globSwarmArray[LAND_PACK_IDX(v_id, grid)];
550*3cdbf90fSMark Adams             ierr_t            = PetscInfo(pack, "makeSwarm %" PetscInt_FMT ".%" PetscInt_FMT ") for batch %" PetscInt_FMT "\n", global_vertex_id_0, grid, LAND_PACK_IDX(v_id, grid));
551*3cdbf90fSMark Adams             ierr_t            = makeSwarm(sw, dim, Np_t[grid][tid], xx_t[grid][tid], yy_t[grid][tid], zz_t[grid][tid]);
552*3cdbf90fSMark Adams             if (ierr_t) ierr = ierr_t;
553*3cdbf90fSMark Adams           }
554*3cdbf90fSMark Adams         }
555*3cdbf90fSMark Adams       }
556*3cdbf90fSMark Adams       //PetscPragmaOMP(parallel for)
557*3cdbf90fSMark Adams       for (int tid = 0; tid < numthreads; tid++) {
558*3cdbf90fSMark Adams         const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id;
559*3cdbf90fSMark Adams         if (glb_v_id < num_vertices) {
560f3237bfeSMark Adams           for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids
561f3237bfeSMark Adams             PetscErrorCode ierr_t;
562f3237bfeSMark Adams             DM             dm = grid_dm[grid];
563*3cdbf90fSMark Adams             DM             sw = globSwarmArray[LAND_PACK_IDX(v_id, grid)];
564*3cdbf90fSMark Adams             ierr_t            = PetscInfo(pack, "createMp %" PetscInt_FMT ".%" PetscInt_FMT ") for batch %" PetscInt_FMT "\n", global_vertex_id_0, grid, LAND_PACK_IDX(v_id, grid));
565*3cdbf90fSMark Adams             ierr_t            = createMp(dm, sw, &globMpArray[LAND_PACK_IDX(v_id, grid)]);
566*3cdbf90fSMark Adams             if (ierr_t) ierr = ierr_t;
567*3cdbf90fSMark Adams           }
568*3cdbf90fSMark Adams         }
569*3cdbf90fSMark Adams       }
570*3cdbf90fSMark Adams       PetscCheck(!ierr, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Error in OMP loop. ierr = %d", (int)ierr);
571*3cdbf90fSMark Adams       // p --> g: set X
572*3cdbf90fSMark Adams       // PetscPragmaOMP(parallel for)
573*3cdbf90fSMark Adams       for (int tid = 0; tid < numthreads; tid++) {
574*3cdbf90fSMark Adams         const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id;
575*3cdbf90fSMark Adams         if (glb_v_id < num_vertices) {
576*3cdbf90fSMark Adams           for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids
577*3cdbf90fSMark Adams             PetscErrorCode ierr_t;
578*3cdbf90fSMark Adams             DM             dm   = grid_dm[grid];
579*3cdbf90fSMark Adams             DM             sw   = globSwarmArray[LAND_PACK_IDX(v_id, grid)];
580*3cdbf90fSMark Adams             Vec            subX = globXArray[LAND_PACK_IDX(v_id, grid)], work = t_fhat[grid][tid];
581*3cdbf90fSMark Adams             ierr_t = PetscInfo(pack, "particlesToGrid %" PetscInt_FMT ".%" PetscInt_FMT ") for local batch %" PetscInt_FMT "\n", global_vertex_id_0, grid, LAND_PACK_IDX(v_id, grid));
582*3cdbf90fSMark Adams             ierr_t = particlesToGrid(dm, sw, Np_t[grid][tid], tid, dim, wp_t[grid][tid], subX, globMpArray[LAND_PACK_IDX(v_id, grid)]);
583f3237bfeSMark Adams             if (ierr_t) ierr = ierr_t;
584f3237bfeSMark Adams             // u = M^_1 f_w
585f3237bfeSMark Adams             ierr_t = VecCopy(subX, work);
586f3237bfeSMark Adams             ierr_t = KSPSolve(t_ksp[grid][tid], work, subX);
587f3237bfeSMark Adams             if (ierr_t) ierr = ierr_t;
588f3237bfeSMark Adams           }
589f3237bfeSMark Adams         }
590f3237bfeSMark Adams       }
5916b664d00Smarkadams4       PetscCheck(!ierr, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Error in OMP loop. ierr = %d", (int)ierr);
592f3237bfeSMark Adams       /* Cleanup */
593f3237bfeSMark Adams       for (int tid = 0; tid < numthreads; tid++) {
594*3cdbf90fSMark Adams         const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id;
595*3cdbf90fSMark Adams         if (glb_v_id < num_vertices) {
596*3cdbf90fSMark Adams           PetscCall(PetscInfo(pack, "Free for global batch %" PetscInt_FMT " of %" PetscInt_FMT "\n", glb_v_id + 1, num_vertices));
597f3237bfeSMark Adams           for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids
598f3237bfeSMark Adams             PetscCall(PetscFree4(xx_t[grid][tid], yy_t[grid][tid], wp_t[grid][tid], zz_t[grid][tid]));
599f3237bfeSMark Adams           }
600f3237bfeSMark Adams         } // active
601*3cdbf90fSMark Adams       }   // threads
602*3cdbf90fSMark Adams     }     // (fake) particle loop
603*3cdbf90fSMark Adams     // standard view of initial conditions
604*3cdbf90fSMark Adams     if (v_target >= global_vertex_id_0 && v_target < global_vertex_id_0 + ctx->batch_sz) {
605*3cdbf90fSMark Adams       PetscCall(DMSetOutputSequenceNumber(ctx->plex[g_target], 0, 0.0));
606*3cdbf90fSMark Adams       PetscCall(VecViewFromOptions(globXArray[LAND_PACK_IDX(v_target % ctx->batch_sz, g_target)], NULL, "-ex30_vec_view"));
607f3237bfeSMark Adams     }
608*3cdbf90fSMark Adams     // coarse graining moments
609*3cdbf90fSMark Adams     if (v_target >= global_vertex_id_0 && v_target < global_vertex_id_0 + ctx->batch_sz) {
610*3cdbf90fSMark Adams       const PetscInt v_id = v_target % ctx->batch_sz;
611f3237bfeSMark Adams       for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
612f3237bfeSMark Adams         PetscDataType dtype;
613f3237bfeSMark Adams         PetscReal    *wp, *coords;
614*3cdbf90fSMark Adams         DM            sw = globSwarmArray[LAND_PACK_IDX(v_id, grid)];
615*3cdbf90fSMark Adams         Vec           work, subX = globXArray[LAND_PACK_IDX(v_id, grid)];
616*3cdbf90fSMark Adams         PetscInt      bs, NN     = nTargetP[grid];
617*3cdbf90fSMark Adams         // C-G moments
618*3cdbf90fSMark Adams         PetscCall(VecDuplicate(subX, &work));
619*3cdbf90fSMark Adams         PetscCall(gridToParticles(grid_dm[grid], sw, subX, work, globMpArray[LAND_PACK_IDX(v_id, grid)], g_Mass[grid]));
620*3cdbf90fSMark Adams         PetscCall(VecDestroy(&work));
621*3cdbf90fSMark Adams         // moments
622f3237bfeSMark Adams         PetscCall(DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
623*3cdbf90fSMark Adams         PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void **)&wp)); // could get NN from sw - todo
624*3cdbf90fSMark Adams         for (int pp = 0; pp < NN; pp++) {
625*3cdbf90fSMark 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]], ww = w * printCtx->N_inv;
626*3cdbf90fSMark Adams           for (int i = 0; i < dim; ++i) v2 += PetscSqr(coords[pp * dim + i]);
627*3cdbf90fSMark Adams           moments_1a[0] += w;
628*3cdbf90fSMark Adams           moments_1a[1] += w * ctx->v_0 * coords[pp * dim + 1]; // z-momentum
629*3cdbf90fSMark Adams           moments_1a[2] += w * ctx->v_0 * ctx->v_0 * v2;
630*3cdbf90fSMark Adams           if (ww > PETSC_REAL_MIN) {
631*3cdbf90fSMark Adams             moments_1a[3] -= ww * PetscLogReal(ww);
632*3cdbf90fSMark Adams             PetscCheck(ww < 1 - PETSC_MACHINE_EPSILON, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "ww (%g) > 1", (double)ww);
633*3cdbf90fSMark Adams           }
634f3237bfeSMark Adams         }
635f3237bfeSMark Adams         PetscCall(DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
636f3237bfeSMark Adams         PetscCall(DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void **)&wp));
637f3237bfeSMark Adams       }
638f3237bfeSMark Adams     }
639*3cdbf90fSMark Adams     // restore vector
640f3237bfeSMark Adams     PetscCall(DMCompositeRestoreAccessArray(pack, X, nDMs, NULL, globXArray));
641*3cdbf90fSMark Adams     // view
642*3cdbf90fSMark Adams     PetscCall(DMPlexLandauPrintNorms(X, 0));
643*3cdbf90fSMark Adams     // advance
644*3cdbf90fSMark Adams     PetscCall(TSSetSolution(ts, X));
645*3cdbf90fSMark Adams     PetscCall(PetscInfo(pack, "Advance vertex %" PetscInt_FMT " to %" PetscInt_FMT " (with padding)\n", global_vertex_id_0, global_vertex_id_0 + ctx->batch_sz));
646*3cdbf90fSMark Adams     PetscCall(TSSetPostStep(ts, PostStep));
647*3cdbf90fSMark Adams     PetscCall(PostStep(ts));
648*3cdbf90fSMark Adams     PetscCall(TSSolve(ts, X));
649*3cdbf90fSMark Adams     // view
650*3cdbf90fSMark Adams     PetscCall(DMPlexLandauPrintNorms(X, 1));
651*3cdbf90fSMark Adams     PetscCall(DMCompositeGetAccessArray(pack, X, nDMs, NULL, globXArray));
652*3cdbf90fSMark Adams     if (v_target >= global_vertex_id_0 && v_target < global_vertex_id_0 + ctx->batch_sz) {
653*3cdbf90fSMark Adams       PetscCall(DMSetOutputSequenceNumber(ctx->plex[g_target], 1, dt_init));
654*3cdbf90fSMark Adams       PetscCall(VecViewFromOptions(globXArray[LAND_PACK_IDX(v_target % ctx->batch_sz, g_target)], NULL, "-ex30_vec_view"));
655*3cdbf90fSMark Adams     }
656*3cdbf90fSMark Adams     // particles to grid, compute moments and entropy
657*3cdbf90fSMark 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));
658*3cdbf90fSMark Adams     // restore vector
659*3cdbf90fSMark Adams     PetscCall(DMCompositeRestoreAccessArray(pack, X, nDMs, NULL, globXArray));
660*3cdbf90fSMark Adams     // cleanup
661*3cdbf90fSMark Adams     for (PetscInt v_id_0 = 0; v_id_0 < ctx->batch_sz; v_id_0 += numthreads) {
662*3cdbf90fSMark Adams       for (int tid = 0; tid < numthreads; tid++) {
663*3cdbf90fSMark Adams         const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id;
664*3cdbf90fSMark Adams         if (glb_v_id < num_vertices) {
665*3cdbf90fSMark Adams           for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
666*3cdbf90fSMark Adams             PetscCall(DMDestroy(&globSwarmArray[LAND_PACK_IDX(v_id, grid)]));
667*3cdbf90fSMark Adams             PetscCall(MatDestroy(&globMpArray[LAND_PACK_IDX(v_id, grid)]));
668*3cdbf90fSMark Adams           }
669*3cdbf90fSMark Adams         }
670*3cdbf90fSMark Adams       }
671*3cdbf90fSMark Adams     }
672f3237bfeSMark Adams   } // user batch
673f3237bfeSMark Adams   /* Cleanup */
674f3237bfeSMark Adams   PetscCall(PetscFree(globXArray));
675f3237bfeSMark Adams   PetscCall(PetscFree(globSwarmArray));
676f3237bfeSMark Adams   PetscCall(PetscFree(globMpArray));
677*3cdbf90fSMark Adams   PetscCall(PetscFree(printCtx));
678f3237bfeSMark Adams   // clean up mass matrices
679f3237bfeSMark Adams   for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids
680f3237bfeSMark Adams     PetscCall(MatDestroy(&g_Mass[grid]));
681f3237bfeSMark Adams     for (int tid = 0; tid < numthreads; tid++) {
682f3237bfeSMark Adams       PetscCall(VecDestroy(&t_fhat[grid][tid]));
683f3237bfeSMark Adams       PetscCall(KSPDestroy(&t_ksp[grid][tid]));
684f3237bfeSMark Adams     }
685f3237bfeSMark Adams   }
686*3cdbf90fSMark Adams   PetscCall(PetscInfo(X, "Moments:\t         number density      x-momentum          energy             entropy : # OMP threads %g\n", (double)numthreads));
687*3cdbf90fSMark Adams   PetscCall(PetscInfo(X, "\tInitial:         %18.12e %19.12e %18.12e %e\n", (double)moments_0[0], (double)moments_0[1], (double)moments_0[2], (double)moments_0[3]));
688*3cdbf90fSMark Adams   PetscCall(PetscInfo(X, "\tCoarse-graining: %18.12e %19.12e %18.12e %e\n", (double)moments_1a[0], (double)moments_1a[1], (double)moments_1a[2], (double)moments_1a[3]));
689*3cdbf90fSMark Adams   PetscCall(PetscInfo(X, "\tLandau:          %18.12e %19.12e %18.12e %e\n", (double)moments_1b[0], (double)moments_1b[1], (double)moments_1b[2], (double)moments_1b[3]));
690*3cdbf90fSMark Adams   PetscCall(PetscInfo(X, "Coarse-graining entropy generation = %e ; Landau entropy generation = %e\n", (double)(moments_1a[3] - moments_0[3]), (double)(moments_1b[3] - moments_0[3])));
691*3cdbf90fSMark Adams   PetscCall(PetscInfo(X, "(relative) energy conservation: Coarse-graining = %e ; Landau = %e\n", (double)(moments_1a[2] - moments_0[2]) / (double)moments_0[2], (double)(moments_1b[2] - moments_0[2]) / (double)moments_0[2]));
692*3cdbf90fSMark Adams 
6933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
694f3237bfeSMark Adams }
695f3237bfeSMark Adams 
696d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
697d71ae5a4SJacob Faibussowitsch {
698f3237bfeSMark Adams   DM         pack;
699f3237bfeSMark Adams   Vec        X;
700*3cdbf90fSMark Adams   PetscInt   dim = 2, num_vertices = 1, Np = 10, v_target = 0, gtarget = 0;
701f3237bfeSMark Adams   TS         ts;
702f3237bfeSMark Adams   Mat        J;
703f3237bfeSMark Adams   LandauCtx *ctx;
704*3cdbf90fSMark Adams   PetscReal  shift = 0;
705f3237bfeSMark Adams 
706327415f7SBarry Smith   PetscFunctionBeginUser;
707f3237bfeSMark Adams   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
708f3237bfeSMark Adams   // process args
709d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_SELF, "", "Collision Options", "DMPLEX");
710*3cdbf90fSMark Adams   PetscCall(PetscOptionsInt("-number_spatial_vertices", "Number of user spatial vertices to be batched for Landau", "ex30.c", num_vertices, &num_vertices, NULL));
711f3237bfeSMark Adams   PetscCall(PetscOptionsInt("-dim", "Velocity space dimension", "ex30.c", dim, &dim, NULL));
712f3237bfeSMark 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));
713*3cdbf90fSMark Adams   PetscCall(PetscOptionsInt("-vertex_view_target", "Vertex to view with diagnostics", "ex30.c", v_target, &v_target, NULL));
714*3cdbf90fSMark Adams   PetscCall(PetscOptionsReal("-e_shift", "Bim-Maxwellian shift", "ex30.c", shift, &shift, NULL));
715*3cdbf90fSMark 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);
716*3cdbf90fSMark Adams   PetscCall(PetscOptionsInt("-grid_view_target", "Grid to view with diagnostics", "ex30.c", gtarget, &gtarget, NULL));
717d0609cedSBarry Smith   PetscOptionsEnd();
718f3237bfeSMark Adams   /* Create a mesh */
719f3237bfeSMark Adams   PetscCall(DMPlexLandauCreateVelocitySpace(PETSC_COMM_SELF, dim, "", &X, &J, &pack));
720f3237bfeSMark Adams   PetscCall(DMSetUp(pack));
721f3237bfeSMark Adams   PetscCall(DMSetOutputSequenceNumber(pack, 0, 0.0));
722f3237bfeSMark Adams   PetscCall(DMGetApplicationContext(pack, &ctx));
723f3237bfeSMark 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);
724*3cdbf90fSMark Adams   PetscCheck(num_vertices >= ctx->batch_sz, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Number of vertices %" PetscInt_FMT " should be <= batch size %" PetscInt_FMT, num_vertices, ctx->batch_sz);
725f3237bfeSMark Adams   /* Create timestepping solver context */
726f3237bfeSMark Adams   PetscCall(TSCreate(PETSC_COMM_SELF, &ts));
727f3237bfeSMark Adams   PetscCall(TSSetDM(ts, pack));
728f3237bfeSMark Adams   PetscCall(TSSetIFunction(ts, NULL, DMPlexLandauIFunction, NULL));
729f3237bfeSMark Adams   PetscCall(TSSetIJacobian(ts, J, J, DMPlexLandauIJacobian, NULL));
730f3237bfeSMark Adams   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
731f3237bfeSMark Adams   PetscCall(TSSetFromOptions(ts));
732f3237bfeSMark Adams   PetscCall(PetscObjectSetName((PetscObject)X, "X"));
7336b664d00Smarkadams4   // do particle advance
734*3cdbf90fSMark Adams   PetscCall(go(ts, X, num_vertices, Np, dim, v_target, gtarget, shift));
735984ed092SMark Adams   PetscCall(MatZeroEntries(J)); // need to zero out so as to not reuse it in Landau's logic
736f3237bfeSMark Adams   /* clean up */
737f3237bfeSMark Adams   PetscCall(DMPlexLandauDestroyVelocitySpace(&pack));
738f3237bfeSMark Adams   PetscCall(TSDestroy(&ts));
739f3237bfeSMark Adams   PetscCall(VecDestroy(&X));
740f3237bfeSMark Adams   PetscCall(PetscFinalize());
741f3237bfeSMark Adams   return 0;
742f3237bfeSMark Adams }
743f3237bfeSMark Adams 
744f3237bfeSMark Adams /*TEST
745f3237bfeSMark Adams 
746f3237bfeSMark Adams   build:
747*3cdbf90fSMark Adams     requires: !complex
748f3237bfeSMark Adams 
749f3237bfeSMark Adams   testset:
750984ed092SMark Adams     requires: double defined(PETSC_USE_DMLANDAU_2D)
751f3237bfeSMark Adams     output_file: output/ex30_0.out
752*3cdbf90fSMark Adams     args: -dim 2 -petscspace_degree 3 -dm_landau_num_species_grid 1,1,1 -dm_refine 1 -number_particles_per_dimension 10 -dm_plex_hash_location \
753*3cdbf90fSMark Adams           -dm_landau_batch_size 4 -number_spatial_vertices 5 -dm_landau_batch_view_idx 1 -vertex_view_target 2 -grid_view_target 1 \
754f3237bfeSMark 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 \
755*3cdbf90fSMark 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 \
7566b664d00Smarkadams4           -ksp_type preonly -pc_type lu -dm_landau_verbose 4 \
757*3cdbf90fSMark Adams           -ptof_ksp_type cg -ptof_pc_type jacobi -ptof_ksp_rtol 1e-12\
758f3237bfeSMark Adams           -snes_converged_reason -snes_monitor -snes_rtol 1e-14 -snes_stol 1e-14\
759*3cdbf90fSMark 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
760f3237bfeSMark Adams 
761f3237bfeSMark Adams     test:
762f3237bfeSMark Adams       suffix: cpu
763f3237bfeSMark Adams       args: -dm_landau_device_type cpu
764f3237bfeSMark Adams     test:
765f3237bfeSMark Adams       suffix: kokkos
766*3cdbf90fSMark Adams       requires: kokkos_kernels !openmp
767*3cdbf90fSMark 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
768f3237bfeSMark Adams 
769984ed092SMark Adams   testset:
770984ed092SMark Adams     requires: double !defined(PETSC_USE_DMLANDAU_2D)
771984ed092SMark Adams     output_file: output/ex30_3d.out
772*3cdbf90fSMark Adams     args: -dim 3 -petscspace_degree 2 -dm_landau_num_species_grid 1,1,1 -dm_refine 0 -number_particles_per_dimension 5 -dm_plex_hash_location \
773*3cdbf90fSMark Adams           -dm_landau_batch_size 1 -number_spatial_vertices 1 -dm_landau_batch_view_idx 0 -vertex_view_target 0 -grid_view_target 0 \
774984ed092SMark 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 \
775*3cdbf90fSMark Adams           -ftop_ksp_rtol 1e-12 -ftop_ksp_type cg -ftop_pc_type jacobi \
776984ed092SMark Adams           -ksp_type preonly -pc_type lu \
777*3cdbf90fSMark Adams           -ptof_ksp_type cg -ptof_pc_type jacobi -ptof_ksp_rtol 1e-12\
778*3cdbf90fSMark Adams           -snes_converged_reason -snes_monitor -snes_rtol 1e-12 -snes_stol 1e-12\
779*3cdbf90fSMark Adams           -ts_dt 0.1 -ts_exact_final_time stepover -ts_max_snes_failures -1 -ts_max_steps 1 -ts_monitor -ts_type beuler
780984ed092SMark Adams 
781984ed092SMark Adams     test:
782984ed092SMark Adams       suffix: cpu_3d
783984ed092SMark Adams       args: -dm_landau_device_type cpu
784984ed092SMark Adams     test:
785984ed092SMark Adams       suffix: kokkos_3d
786*3cdbf90fSMark Adams       requires: kokkos_kernels !openmp
787*3cdbf90fSMark 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
788*3cdbf90fSMark Adams 
789*3cdbf90fSMark Adams   testset:
790*3cdbf90fSMark Adams     requires: !complex double defined(PETSC_USE_DMLANDAU_2D) !cuda
791*3cdbf90fSMark Adams     args: -dm_landau_domain_radius 6 -dm_refine 2 -dm_landau_num_species_grid 1 -dm_landau_thermal_temps 1 -petscspace_degree 3 -snes_converged_reason -ts_type beuler -ts_dt 1 -ts_max_steps 1 -ksp_type preonly -pc_type lu -snes_rtol 1e-12 -snes_stol 1e-12 -dm_landau_device_type cpu -number_particles_per_dimension 30 -e_shift 3 -ftop_ksp_rtol 1e-12 -ptof_ksp_rtol 1e-12 -dm_landau_batch_size 4 -number_spatial_vertices 4 -grid_view_target 0 -vertex_view_target 1
792*3cdbf90fSMark Adams     test:
793*3cdbf90fSMark Adams       suffix: simple
794*3cdbf90fSMark Adams       args: -ex30_dm_view
795*3cdbf90fSMark Adams     test:
796*3cdbf90fSMark Adams       requires: hdf5
797*3cdbf90fSMark Adams       suffix: simple_hdf5
798*3cdbf90fSMark Adams       args: -ex30_dm_view hdf5:sol_e.h5 -ex30_vec_view hdf5:sol_e.h5::append
799984ed092SMark Adams 
800f3237bfeSMark Adams TEST*/
801