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