1f3237bfeSMark Adams static char help[] = "Grid based Landau collision operator with PIC interface with OpenMP setup. (one species per grid)\n";
2f3237bfeSMark Adams
3f3237bfeSMark Adams /*
4f3237bfeSMark Adams Support 2.5V with axisymmetric coordinates
5f3237bfeSMark Adams - r,z coordinates
6f3237bfeSMark Adams - Domain and species data input by Landau operator
7f3237bfeSMark Adams - "radius" for each grid, normalized with electron thermal velocity
8f3237bfeSMark Adams - Domain: (0,radius) x (-radius,radius), thus first coordinate x[0] is perpendicular velocity and 2pi*x[0] term is added for axisymmetric
9f3237bfeSMark Adams Supports full 3V
10f3237bfeSMark Adams
11f3237bfeSMark Adams */
12f3237bfeSMark Adams
1346233b44SBarry Smith #include <petscdmplex.h>
1446233b44SBarry Smith #include <petscds.h>
1546233b44SBarry Smith #include <petscdmswarm.h>
1646233b44SBarry Smith #include <petscksp.h>
17f3237bfeSMark Adams #include <petsc/private/petscimpl.h>
18f3237bfeSMark Adams #if defined(PETSC_HAVE_OPENMP) && defined(PETSC_HAVE_THREADSAFETY)
19f3237bfeSMark Adams #include <omp.h>
20f3237bfeSMark Adams #endif
21f3237bfeSMark Adams #include <petsclandau.h>
22f3237bfeSMark Adams #include <petscdmcomposite.h>
23f3237bfeSMark Adams
24f3237bfeSMark Adams typedef struct {
25f3237bfeSMark Adams Mat MpTrans;
26f3237bfeSMark Adams Mat Mp;
27f3237bfeSMark Adams Vec ff;
28f3237bfeSMark Adams Vec uu;
29f3237bfeSMark Adams } MatShellCtx;
30f3237bfeSMark Adams
313cdbf90fSMark Adams typedef struct {
323cdbf90fSMark Adams PetscInt v_target;
331b431c67SMark Adams PetscInt g_target;
341b431c67SMark Adams PetscInt global_vertex_id_0;
353cdbf90fSMark Adams DM *globSwarmArray;
363cdbf90fSMark Adams LandauCtx *ctx;
373cdbf90fSMark Adams DM *grid_dm;
383cdbf90fSMark Adams Mat *g_Mass;
393cdbf90fSMark Adams Mat *globMpArray;
403cdbf90fSMark Adams Vec *globXArray;
413cdbf90fSMark Adams PetscBool print;
421b431c67SMark Adams PetscBool print_entropy;
433cdbf90fSMark Adams } PrintCtx;
443cdbf90fSMark Adams
MatMultMtM_SeqAIJ(Mat MtM,Vec xx,Vec yy)45d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultMtM_SeqAIJ(Mat MtM, Vec xx, Vec yy)
46d71ae5a4SJacob Faibussowitsch {
47f3237bfeSMark Adams MatShellCtx *matshellctx;
48f3237bfeSMark Adams
49f3237bfeSMark Adams PetscFunctionBeginUser;
50f3237bfeSMark Adams PetscCall(MatShellGetContext(MtM, &matshellctx));
51f3237bfeSMark Adams PetscCheck(matshellctx, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "No context");
52f3237bfeSMark Adams PetscCall(MatMult(matshellctx->Mp, xx, matshellctx->ff));
53f3237bfeSMark Adams PetscCall(MatMult(matshellctx->MpTrans, matshellctx->ff, yy));
543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
55f3237bfeSMark Adams }
56f3237bfeSMark Adams
MatMultAddMtM_SeqAIJ(Mat MtM,Vec xx,Vec yy,Vec zz)57d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultAddMtM_SeqAIJ(Mat MtM, Vec xx, Vec yy, Vec zz)
58d71ae5a4SJacob Faibussowitsch {
59f3237bfeSMark Adams MatShellCtx *matshellctx;
60f3237bfeSMark Adams
61f3237bfeSMark Adams PetscFunctionBeginUser;
62f3237bfeSMark Adams PetscCall(MatShellGetContext(MtM, &matshellctx));
63f3237bfeSMark Adams PetscCheck(matshellctx, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "No context");
64f3237bfeSMark Adams PetscCall(MatMult(matshellctx->Mp, xx, matshellctx->ff));
65f3237bfeSMark Adams PetscCall(MatMultAdd(matshellctx->MpTrans, matshellctx->ff, yy, zz));
663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
67f3237bfeSMark Adams }
68f3237bfeSMark Adams
createSwarm(const DM dm,PetscInt dim,DM * sw)69d71ae5a4SJacob Faibussowitsch PetscErrorCode createSwarm(const DM dm, PetscInt dim, DM *sw)
70d71ae5a4SJacob Faibussowitsch {
71f3237bfeSMark Adams PetscInt Nc = 1;
72f3237bfeSMark Adams
73f3237bfeSMark Adams PetscFunctionBeginUser;
74f3237bfeSMark Adams PetscCall(DMCreate(PETSC_COMM_SELF, sw));
75f3237bfeSMark Adams PetscCall(DMSetType(*sw, DMSWARM));
76f3237bfeSMark Adams PetscCall(DMSetDimension(*sw, dim));
77f3237bfeSMark Adams PetscCall(DMSwarmSetType(*sw, DMSWARM_PIC));
78f3237bfeSMark Adams PetscCall(DMSwarmSetCellDM(*sw, dm));
793cdbf90fSMark Adams PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "w_q", Nc, PETSC_REAL));
80f3237bfeSMark Adams PetscCall(DMSwarmFinalizeFieldRegister(*sw));
81f3237bfeSMark Adams PetscCall(DMSetFromOptions(*sw));
823cdbf90fSMark Adams PetscCall(PetscObjectSetName((PetscObject)*sw, "Particle Grid"));
833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
84f3237bfeSMark Adams }
85f3237bfeSMark Adams
makeSwarm(DM sw,const PetscInt dim,const PetscInt Np,const PetscReal xx[],const PetscReal yy[],const PetscReal zz[])863cdbf90fSMark Adams static PetscErrorCode makeSwarm(DM sw, const PetscInt dim, const PetscInt Np, const PetscReal xx[], const PetscReal yy[], const PetscReal zz[])
873cdbf90fSMark Adams {
883cdbf90fSMark Adams PetscReal *coords;
893cdbf90fSMark Adams PetscDataType dtype;
903cdbf90fSMark Adams PetscInt bs, p, zero = 0;
913cdbf90fSMark Adams
924d86920dSPierre Jolivet PetscFunctionBeginUser;
933cdbf90fSMark Adams PetscCall(DMSwarmSetLocalSizes(sw, Np, zero));
943cdbf90fSMark Adams PetscCall(DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
953cdbf90fSMark Adams for (p = 0; p < Np; p++) {
963cdbf90fSMark Adams coords[p * dim + 0] = xx[p];
973cdbf90fSMark Adams coords[p * dim + 1] = yy[p];
983cdbf90fSMark Adams if (dim == 3) coords[p * dim + 2] = zz[p];
993cdbf90fSMark Adams }
1003cdbf90fSMark Adams PetscCall(DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
101d52c2f21SMatthew G. Knepley PetscCall(DMSwarmVectorDefineField(sw, "w_q"));
1023cdbf90fSMark Adams PetscFunctionReturn(PETSC_SUCCESS);
1033cdbf90fSMark Adams }
1043cdbf90fSMark Adams
createMp(const DM dm,DM sw,Mat * Mp_out)1053cdbf90fSMark Adams static PetscErrorCode createMp(const DM dm, DM sw, Mat *Mp_out)
1063cdbf90fSMark Adams {
1073cdbf90fSMark Adams PetscBool removePoints = PETSC_TRUE;
1083cdbf90fSMark Adams Mat M_p;
1094d86920dSPierre Jolivet
1103cdbf90fSMark Adams PetscFunctionBeginUser;
1113cdbf90fSMark Adams // migrate after coords are set
1123cdbf90fSMark Adams PetscCall(DMSwarmMigrate(sw, removePoints));
113d043ef4cSMark Adams //
1143cdbf90fSMark Adams PetscCall(PetscObjectSetName((PetscObject)sw, "Particle Grid"));
115d043ef4cSMark Adams
116d043ef4cSMark Adams /* PetscInt N,*count,nmin=10000,nmax=0,ntot=0; */
117d043ef4cSMark Adams /* // count */
118d043ef4cSMark Adams /* PetscCall(DMSwarmCreatePointPerCellCount(sw, &N, &count)); */
119d043ef4cSMark Adams /* for (int i=0, n; i< N ; i++) { */
120d043ef4cSMark Adams /* if ((n=count[i]) > nmax) nmax = n; */
121d043ef4cSMark Adams /* if (n < nmin) nmin = n; */
122d043ef4cSMark Adams /* PetscCall(PetscInfo(dm, " %d) %d particles\n", i, n)); */
123d043ef4cSMark Adams /* ntot += n; */
124d043ef4cSMark Adams /* } */
125d043ef4cSMark Adams /* PetscCall(PetscFree(count)); */
126d043ef4cSMark Adams /* PetscCall(PetscInfo(dm, " %" PetscInt_FMT " max particle / cell, and %" PetscInt_FMT " min, ratio = %g, %" PetscInt_FMT " total\n", nmax, nmin, (double)nmax/(double)nmin,ntot)); */
127d043ef4cSMark Adams
1283cdbf90fSMark Adams /* This gives M f = \int_\Omega \phi f, which looks like a rhs for a PDE */
1293cdbf90fSMark Adams PetscCall(DMCreateMassMatrix(sw, dm, &M_p));
1303cdbf90fSMark Adams PetscCall(DMViewFromOptions(sw, NULL, "-ex30_sw_view"));
1313cdbf90fSMark Adams // output
1323cdbf90fSMark Adams *Mp_out = M_p;
1333cdbf90fSMark Adams PetscFunctionReturn(PETSC_SUCCESS);
1343cdbf90fSMark Adams }
1353cdbf90fSMark Adams
particlesToGrid(const DM dm,DM sw,const PetscInt a_tid,const PetscInt dim,const PetscReal a_wp[],Vec rho,Mat M_p)1361b431c67SMark Adams static PetscErrorCode particlesToGrid(const DM dm, DM sw, const PetscInt a_tid, const PetscInt dim, const PetscReal a_wp[], Vec rho, Mat M_p)
1373cdbf90fSMark Adams {
1383cdbf90fSMark Adams PetscReal *wq;
1393cdbf90fSMark Adams PetscDataType dtype;
1403cdbf90fSMark Adams Vec ff;
1411b431c67SMark Adams PetscInt bs, p, Np;
1423cdbf90fSMark Adams
1433cdbf90fSMark Adams PetscFunctionBeginUser;
1443cdbf90fSMark Adams PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void **)&wq));
1451b431c67SMark Adams PetscCall(DMSwarmGetLocalSize(sw, &Np));
1463cdbf90fSMark Adams for (p = 0; p < Np; p++) wq[p] = a_wp[p];
1473cdbf90fSMark Adams PetscCall(DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void **)&wq));
1483cdbf90fSMark Adams PetscCall(PetscObjectSetName((PetscObject)rho, "rho"));
1493cdbf90fSMark Adams PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &ff));
1503cdbf90fSMark Adams PetscCall(PetscObjectSetName((PetscObject)ff, "weights"));
1513cdbf90fSMark Adams PetscCall(MatMultTranspose(M_p, ff, rho));
1523cdbf90fSMark Adams PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &ff));
1533cdbf90fSMark Adams PetscFunctionReturn(PETSC_SUCCESS);
1543cdbf90fSMark Adams }
1553cdbf90fSMark Adams
1563cdbf90fSMark Adams //
1573cdbf90fSMark Adams // add grid to arg 'sw.w_q'
1583cdbf90fSMark Adams //
gridToParticles(const DM dm,DM sw,const Vec rhs,Vec work_ferhs,Mat M_p,Mat Mass)159d043ef4cSMark Adams PetscErrorCode gridToParticles(const DM dm, DM sw, const Vec rhs, Vec work_ferhs, Mat M_p, Mat Mass)
160d71ae5a4SJacob Faibussowitsch {
161f3237bfeSMark Adams PetscBool is_lsqr;
162f3237bfeSMark Adams KSP ksp;
163d043ef4cSMark Adams Mat PM_p = NULL, MtM, D = NULL;
164f3237bfeSMark Adams Vec ff;
165f3237bfeSMark Adams PetscInt N, M, nzl;
166d043ef4cSMark Adams MatShellCtx *matshellctx = NULL;
1671b431c67SMark Adams PC pc;
168f3237bfeSMark Adams
169f3237bfeSMark Adams PetscFunctionBeginUser;
1703643261fSMark Adams // 1) apply M in, for Moore-Penrose with mass: Mp (Mp' Mp)^-1 M
171d043ef4cSMark Adams PetscCall(MatMult(Mass, rhs, work_ferhs));
1723643261fSMark Adams // 2) pseudo-inverse, first part: (Mp' Mp)^-1
173f3237bfeSMark Adams PetscCall(KSPCreate(PETSC_COMM_SELF, &ksp));
1741b431c67SMark Adams PetscCall(KSPSetType(ksp, KSPCG));
1751b431c67SMark Adams PetscCall(KSPGetPC(ksp, &pc));
1761b431c67SMark Adams PetscCall(PCSetType(pc, PCJACOBI));
177f3237bfeSMark Adams PetscCall(KSPSetOptionsPrefix(ksp, "ftop_"));
178f3237bfeSMark Adams PetscCall(KSPSetFromOptions(ksp));
179f3237bfeSMark Adams PetscCall(PetscObjectTypeCompare((PetscObject)ksp, KSPLSQR, &is_lsqr));
180f3237bfeSMark Adams if (!is_lsqr) {
181f3237bfeSMark Adams PetscCall(MatGetLocalSize(M_p, &M, &N));
182f3237bfeSMark Adams if (N > M) {
1831b431c67SMark Adams PetscCall(PetscInfo(ksp, " M (%" PetscInt_FMT ") < M (%" PetscInt_FMT ") more vertices than particles: revert to lsqr\n", M, N));
184f3237bfeSMark Adams is_lsqr = PETSC_TRUE;
185f3237bfeSMark Adams PetscCall(KSPSetType(ksp, KSPLSQR));
186d043ef4cSMark Adams PetscCall(PCSetType(pc, PCNONE)); // should not happen, but could solve stable (Mp^T Mp), move projection Mp before solve
187f3237bfeSMark Adams } else {
188f3237bfeSMark Adams PetscCall(PetscNew(&matshellctx));
189d043ef4cSMark Adams PetscCall(MatCreateVecs(M_p, &matshellctx->uu, &matshellctx->ff));
190d043ef4cSMark Adams if (0) {
191d043ef4cSMark Adams PetscCall(MatTransposeMatMult(M_p, M_p, MAT_INITIAL_MATRIX, 4, &MtM));
192d043ef4cSMark Adams PetscCall(KSPSetOperators(ksp, MtM, MtM));
193d043ef4cSMark Adams PetscCall(PetscInfo(M_p, "createMtM KSP with explicit Mp'Mp\n"));
194d043ef4cSMark Adams PetscCall(MatViewFromOptions(MtM, NULL, "-ftop2_MtM_mat_view"));
195d043ef4cSMark Adams } else {
196f3237bfeSMark Adams PetscCall(MatCreateShell(PetscObjectComm((PetscObject)dm), N, N, PETSC_DECIDE, PETSC_DECIDE, matshellctx, &MtM));
197f3237bfeSMark Adams PetscCall(MatTranspose(M_p, MAT_INITIAL_MATRIX, &matshellctx->MpTrans));
198f3237bfeSMark Adams matshellctx->Mp = M_p;
19957d50842SBarry Smith PetscCall(MatShellSetOperation(MtM, MATOP_MULT, (PetscErrorCodeFn *)MatMultMtM_SeqAIJ));
20057d50842SBarry Smith PetscCall(MatShellSetOperation(MtM, MATOP_MULT_ADD, (PetscErrorCodeFn *)MatMultAddMtM_SeqAIJ));
201f3237bfeSMark Adams PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, N, N, 1, NULL, &D));
202d043ef4cSMark Adams PetscCall(MatViewFromOptions(matshellctx->MpTrans, NULL, "-ftop2_MpT_mat_view"));
2036497c311SBarry Smith for (PetscInt i = 0; i < N; i++) {
204f3237bfeSMark Adams const PetscScalar *vals;
205f3237bfeSMark Adams const PetscInt *cols;
206f3237bfeSMark Adams PetscScalar dot = 0;
207f3237bfeSMark Adams PetscCall(MatGetRow(matshellctx->MpTrans, i, &nzl, &cols, &vals));
2086497c311SBarry Smith for (PetscInt ii = 0; ii < nzl; ii++) dot += PetscSqr(vals[ii]);
209d043ef4cSMark Adams if (dot < PETSC_MACHINE_EPSILON) {
2106497c311SBarry Smith PetscCall(PetscInfo(ksp, "empty row in pseudo-inverse %d\n", (int)i));
211d043ef4cSMark Adams is_lsqr = PETSC_TRUE; // empty rows
212d043ef4cSMark Adams PetscCall(KSPSetType(ksp, KSPLSQR));
213d043ef4cSMark Adams PetscCall(PCSetType(pc, PCNONE)); // should not happen, but could solve stable (Mp Mp^T), move projection Mp before solve
214d043ef4cSMark Adams // clean up
215d043ef4cSMark Adams PetscCall(MatDestroy(&matshellctx->MpTrans));
216d043ef4cSMark Adams PetscCall(VecDestroy(&matshellctx->ff));
217d043ef4cSMark Adams PetscCall(VecDestroy(&matshellctx->uu));
218d043ef4cSMark Adams PetscCall(MatDestroy(&D));
219d043ef4cSMark Adams PetscCall(MatDestroy(&MtM));
220d043ef4cSMark Adams PetscCall(PetscFree(matshellctx));
221d043ef4cSMark Adams D = NULL;
222d043ef4cSMark Adams break;
223d043ef4cSMark Adams }
224f3237bfeSMark Adams PetscCall(MatSetValue(D, i, i, dot, INSERT_VALUES));
225f3237bfeSMark Adams }
226d043ef4cSMark Adams if (D) {
227f3237bfeSMark Adams PetscCall(MatAssemblyBegin(D, MAT_FINAL_ASSEMBLY));
228f3237bfeSMark Adams PetscCall(MatAssemblyEnd(D, MAT_FINAL_ASSEMBLY));
229f3237bfeSMark Adams PetscCall(PetscInfo(M_p, "createMtMKSP Have %" PetscInt_FMT " eqs, nzl = %" PetscInt_FMT "\n", N, nzl));
230f3237bfeSMark Adams PetscCall(KSPSetOperators(ksp, MtM, D));
231f3237bfeSMark Adams PetscCall(MatViewFromOptions(D, NULL, "-ftop2_D_mat_view"));
232f3237bfeSMark Adams PetscCall(MatViewFromOptions(M_p, NULL, "-ftop2_Mp_mat_view"));
233f3237bfeSMark Adams PetscCall(MatViewFromOptions(matshellctx->MpTrans, NULL, "-ftop2_MpTranspose_mat_view"));
2346b664d00Smarkadams4 PetscCall(MatViewFromOptions(MtM, NULL, "-ftop2_MtM_mat_view"));
235f3237bfeSMark Adams }
236f3237bfeSMark Adams }
237d043ef4cSMark Adams }
238d043ef4cSMark Adams }
239f3237bfeSMark Adams if (is_lsqr) {
240d043ef4cSMark Adams PC pc2;
241f3237bfeSMark Adams PetscBool is_bjac;
242d043ef4cSMark Adams PetscCall(KSPGetPC(ksp, &pc2));
243d043ef4cSMark Adams PetscCall(PetscObjectTypeCompare((PetscObject)pc2, PCBJACOBI, &is_bjac));
244f3237bfeSMark Adams if (is_bjac) {
245f3237bfeSMark Adams PetscCall(DMSwarmCreateMassMatrixSquare(sw, dm, &PM_p));
246f3237bfeSMark Adams PetscCall(KSPSetOperators(ksp, M_p, PM_p));
247f3237bfeSMark Adams } else {
248f3237bfeSMark Adams PetscCall(KSPSetOperators(ksp, M_p, M_p));
249f3237bfeSMark Adams }
250d043ef4cSMark Adams PetscCall(MatViewFromOptions(M_p, NULL, "-ftop2_Mp_mat_view"));
251f3237bfeSMark Adams }
252f3237bfeSMark Adams PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &ff)); // this grabs access
253f3237bfeSMark Adams if (!is_lsqr) {
2544d1837e9SPierre Jolivet PetscCall(KSPSolve(ksp, work_ferhs, matshellctx->uu));
2553643261fSMark Adams // 3) with Moore-Penrose apply Mp: M_p (Mp' Mp)^-1 M
256f3237bfeSMark Adams PetscCall(MatMult(M_p, matshellctx->uu, ff));
257d043ef4cSMark Adams if (D) PetscCall(MatDestroy(&D));
258d043ef4cSMark Adams PetscCall(MatDestroy(&MtM));
259d043ef4cSMark Adams if (matshellctx->MpTrans) PetscCall(MatDestroy(&matshellctx->MpTrans));
260f3237bfeSMark Adams PetscCall(VecDestroy(&matshellctx->ff));
261f3237bfeSMark Adams PetscCall(VecDestroy(&matshellctx->uu));
262f3237bfeSMark Adams PetscCall(PetscFree(matshellctx));
263f3237bfeSMark Adams } else {
2643643261fSMark Adams // finally with LSQR apply M_p^\dagger
2654d1837e9SPierre Jolivet PetscCall(KSPSolveTranspose(ksp, work_ferhs, ff));
266f3237bfeSMark Adams }
267f3237bfeSMark Adams PetscCall(KSPDestroy(&ksp));
268f3237bfeSMark Adams PetscCall(MatDestroy(&PM_p));
269f3237bfeSMark Adams PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &ff));
2703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
271f3237bfeSMark Adams }
272f3237bfeSMark Adams
2733cdbf90fSMark Adams #define EX30_MAX_NUM_THRDS 12
2743cdbf90fSMark Adams #define EX30_MAX_BATCH_SZ 1024
2753cdbf90fSMark Adams //
2763cdbf90fSMark Adams // add grid to arg 'globSwarmArray[].w_q'
2773cdbf90fSMark Adams //
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)2783cdbf90fSMark 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)
279d71ae5a4SJacob Faibussowitsch {
2803cdbf90fSMark Adams PetscErrorCode ierr = (PetscErrorCode)0; // used for inside thread loops
281f3237bfeSMark Adams
282f3237bfeSMark Adams PetscFunctionBeginUser;
2833cdbf90fSMark Adams // map back to particles
2843cdbf90fSMark Adams for (PetscInt v_id_0 = 0; v_id_0 < ctx->batch_sz; v_id_0 += numthreads) {
2853cdbf90fSMark 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));
2863cdbf90fSMark Adams //PetscPragmaOMP(parallel for)
2876497c311SBarry Smith for (PetscInt tid = 0; tid < numthreads; tid++) {
2883cdbf90fSMark Adams const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id + v_id;
2893cdbf90fSMark Adams if (glb_v_id < num_vertices) {
2903cdbf90fSMark Adams for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids
2913cdbf90fSMark Adams PetscErrorCode ierr_t;
2923cdbf90fSMark 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));
2933cdbf90fSMark 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]);
2943cdbf90fSMark Adams if (ierr_t) ierr = ierr_t;
295f3237bfeSMark Adams }
2963cdbf90fSMark Adams }
2973cdbf90fSMark Adams }
298d043ef4cSMark Adams PetscCheck(!ierr, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Error in OMP loop. ierr = %d", (int)ierr);
2993cdbf90fSMark Adams /* Get moments */
3003cdbf90fSMark Adams PetscCall(PetscInfo(grid_dm[0], "Cleanup batches %" PetscInt_FMT " to %" PetscInt_FMT "\n", v_id_0, v_id_0 + numthreads));
3016497c311SBarry Smith for (PetscInt tid = 0; tid < numthreads; tid++) {
3023cdbf90fSMark Adams const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id + v_id;
3033cdbf90fSMark Adams if (glb_v_id == v_target) {
3043cdbf90fSMark Adams for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
3053cdbf90fSMark Adams PetscDataType dtype;
3063cdbf90fSMark Adams PetscReal *wp, *coords;
3073cdbf90fSMark Adams DM sw = globSwarmArray[LAND_PACK_IDX(v_id, grid)];
3083cdbf90fSMark Adams PetscInt npoints, bs = 1;
3093cdbf90fSMark Adams PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void **)&wp)); // take data out here
3103cdbf90fSMark Adams PetscCall(DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
3113cdbf90fSMark Adams PetscCall(DMSwarmGetLocalSize(sw, &npoints));
3126497c311SBarry Smith for (PetscInt p = 0; p < npoints; p++) {
3133cdbf90fSMark 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]];
3146497c311SBarry Smith for (PetscInt i = 0; i < dim; ++i) v2 += PetscSqr(coords[p * dim + i]);
3153cdbf90fSMark Adams moments[0] += w;
3163cdbf90fSMark Adams moments[1] += w * ctx->v_0 * coords[p * dim + 1]; // z-momentum
3171b431c67SMark Adams moments[2] += w * 0.5 * ctx->v_0 * ctx->v_0 * v2;
3183cdbf90fSMark Adams }
3193cdbf90fSMark Adams PetscCall(DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void **)&wp));
320f3237bfeSMark Adams PetscCall(DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
3213cdbf90fSMark Adams }
3223cdbf90fSMark Adams const PetscReal N_inv = 1 / moments[0];
3231b431c67SMark Adams PetscCall(PetscInfo(grid_dm[0], "gridToParticles_private [%" PetscInt_FMT "], n = %g\n", v_id, (double)moments[0]));
3243cdbf90fSMark Adams for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
3253cdbf90fSMark Adams PetscDataType dtype;
3263cdbf90fSMark Adams PetscReal *wp, *coords;
3273cdbf90fSMark Adams DM sw = globSwarmArray[LAND_PACK_IDX(v_id, grid)];
3283cdbf90fSMark Adams PetscInt npoints, bs = 1;
3293cdbf90fSMark Adams PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void **)&wp)); // take data out here
3303cdbf90fSMark Adams PetscCall(DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
3313cdbf90fSMark Adams PetscCall(DMSwarmGetLocalSize(sw, &npoints));
3326497c311SBarry Smith for (PetscInt p = 0; p < npoints; p++) {
3333cdbf90fSMark 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;
3341b431c67SMark Adams if (w > PETSC_REAL_MIN) {
3353cdbf90fSMark Adams moments[3] -= ww * PetscLogReal(ww);
336d043ef4cSMark Adams PetscCheck(ww < 1 - PETSC_MACHINE_EPSILON, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "ww (%g) > 1", (double)ww);
337d043ef4cSMark Adams } else moments[4] -= w; // keep track of density that is lost
3383cdbf90fSMark Adams }
3393cdbf90fSMark Adams PetscCall(DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void **)&wp));
3403cdbf90fSMark Adams PetscCall(DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
3413cdbf90fSMark Adams }
3423cdbf90fSMark Adams }
3433cdbf90fSMark Adams } // thread batch
3443cdbf90fSMark Adams } // batch
3453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
346f3237bfeSMark Adams }
3473cdbf90fSMark Adams
maxwellian(PetscInt dim,const PetscReal x[],PetscReal kt_m,PetscReal n,PetscReal shift,PetscScalar * u)3483cdbf90fSMark Adams static void maxwellian(PetscInt dim, const PetscReal x[], PetscReal kt_m, PetscReal n, PetscReal shift, PetscScalar *u)
349d71ae5a4SJacob Faibussowitsch {
350f3237bfeSMark Adams PetscInt i;
351f3237bfeSMark Adams PetscReal v2 = 0, theta = 2.0 * kt_m; /* theta = 2kT/mc^2 */
352f3237bfeSMark Adams
3533cdbf90fSMark Adams if (shift != 0.) {
3543cdbf90fSMark Adams v2 = 0;
3553cdbf90fSMark Adams for (i = 0; i < dim - 1; ++i) v2 += x[i] * x[i];
3563cdbf90fSMark Adams v2 += (x[dim - 1] - shift) * (x[dim - 1] - shift);
3573cdbf90fSMark Adams /* evaluate the shifted Maxwellian */
3583cdbf90fSMark Adams u[0] += n * PetscPowReal(PETSC_PI * theta, -1.5) * (PetscExpReal(-v2 / theta));
3591b431c67SMark Adams } else {
3601b431c67SMark Adams /* compute the exponents, v^2 */
3611b431c67SMark Adams for (i = 0; i < dim; ++i) v2 += x[i] * x[i];
3621b431c67SMark Adams /* evaluate the Maxwellian */
3631b431c67SMark Adams u[0] += n * PetscPowReal(PETSC_PI * theta, -1.5) * (PetscExpReal(-v2 / theta));
3643cdbf90fSMark Adams }
365f3237bfeSMark Adams }
366f3237bfeSMark Adams
PostStep(TS ts)3673cdbf90fSMark Adams static PetscErrorCode PostStep(TS ts)
3683cdbf90fSMark Adams {
3691b431c67SMark Adams PetscInt n, dim, nDMs, v_id;
3703cdbf90fSMark Adams PetscReal t;
3713cdbf90fSMark Adams LandauCtx *ctx;
3723cdbf90fSMark Adams Vec X;
3733cdbf90fSMark Adams PrintCtx *printCtx;
3743cdbf90fSMark Adams DM pack;
3751b431c67SMark Adams PetscReal moments[5], e_grid[LANDAU_MAX_GRIDS];
3763cdbf90fSMark Adams
3773cdbf90fSMark Adams PetscFunctionBeginUser;
3783cdbf90fSMark Adams PetscCall(TSGetApplicationContext(ts, &printCtx));
3791b431c67SMark Adams if (!printCtx->print && !printCtx->print_entropy) PetscFunctionReturn(PETSC_SUCCESS);
3803cdbf90fSMark Adams ctx = printCtx->ctx;
3811b431c67SMark Adams if (printCtx->v_target < printCtx->global_vertex_id_0 || printCtx->v_target >= printCtx->global_vertex_id_0 + ctx->batch_sz) PetscFunctionReturn(PETSC_SUCCESS);
3826497c311SBarry Smith for (PetscInt i = 0; i < 5; i++) moments[i] = 0;
3836497c311SBarry Smith for (PetscInt i = 0; i < LANDAU_MAX_GRIDS; i++) e_grid[i] = 0;
3841b431c67SMark Adams v_id = printCtx->v_target % ctx->batch_sz;
3853cdbf90fSMark Adams PetscCall(TSGetDM(ts, &pack));
3863cdbf90fSMark Adams PetscCall(DMGetDimension(pack, &dim));
3873cdbf90fSMark Adams PetscCall(DMCompositeGetNumberDM(pack, &nDMs)); // number of vertices * number of grids
3883cdbf90fSMark Adams PetscCall(TSGetSolution(ts, &X));
3893cdbf90fSMark Adams PetscCall(TSGetStepNumber(ts, &n));
3903cdbf90fSMark Adams PetscCall(TSGetTime(ts, &t));
3913cdbf90fSMark Adams PetscCall(DMCompositeGetAccessArray(pack, X, nDMs, NULL, printCtx->globXArray));
392d043ef4cSMark Adams if (printCtx->print_entropy && printCtx->v_target >= 0 && 0) {
3933cdbf90fSMark Adams for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
3943cdbf90fSMark Adams PetscDataType dtype;
3953cdbf90fSMark Adams PetscReal *wp, *coords;
3963cdbf90fSMark Adams DM sw = printCtx->globSwarmArray[LAND_PACK_IDX(v_id, grid)];
3973cdbf90fSMark Adams Vec work, subX = printCtx->globXArray[LAND_PACK_IDX(v_id, grid)];
3981b431c67SMark Adams PetscInt bs, NN;
3993cdbf90fSMark Adams // C-G moments
4003cdbf90fSMark Adams PetscCall(VecDuplicate(subX, &work));
4013cdbf90fSMark Adams PetscCall(gridToParticles(printCtx->grid_dm[grid], sw, subX, work, printCtx->globMpArray[LAND_PACK_IDX(v_id, grid)], printCtx->g_Mass[grid]));
4023cdbf90fSMark Adams PetscCall(VecDestroy(&work));
4033cdbf90fSMark Adams // moments
4043cdbf90fSMark Adams PetscCall(DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
4051b431c67SMark Adams PetscCall(DMSwarmGetLocalSize(sw, &NN));
4061b431c67SMark Adams PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void **)&wp));
4076497c311SBarry Smith for (PetscInt pp = 0; pp < NN; pp++) {
4081b431c67SMark 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]];
4096497c311SBarry Smith for (PetscInt i = 0; i < dim; ++i) v2 += PetscSqr(coords[pp * dim + i]);
4103cdbf90fSMark Adams moments[0] += w;
4113cdbf90fSMark Adams moments[1] += w * ctx->v_0 * coords[pp * dim + 1]; // z-momentum
4121b431c67SMark Adams moments[2] += w * 0.5 * ctx->v_0 * ctx->v_0 * v2;
4131b431c67SMark Adams e_grid[grid] += w * 0.5 * ctx->v_0 * ctx->v_0 * v2;
4143cdbf90fSMark Adams }
4153cdbf90fSMark Adams PetscCall(DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
4163cdbf90fSMark Adams PetscCall(DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void **)&wp));
4173cdbf90fSMark Adams }
4181b431c67SMark Adams // entropy
4191b431c67SMark Adams const PetscReal N_inv = 1 / moments[0];
4201b431c67SMark Adams for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
4211b431c67SMark Adams PetscDataType dtype;
4221b431c67SMark Adams PetscReal *wp, *coords;
4231b431c67SMark Adams DM sw = printCtx->globSwarmArray[LAND_PACK_IDX(v_id, grid)];
4241b431c67SMark Adams PetscInt bs, NN;
4251b431c67SMark Adams PetscCall(DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
4261b431c67SMark Adams PetscCall(DMSwarmGetLocalSize(sw, &NN));
4271b431c67SMark Adams PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void **)&wp));
4286497c311SBarry Smith for (PetscInt pp = 0; pp < NN; pp++) {
4291b431c67SMark Adams PetscReal fact = (dim == 2) ? 2.0 * PETSC_PI * coords[pp * dim + 0] : 1, w = fact * wp[pp] * ctx->n_0 * ctx->masses[ctx->species_offset[grid]], ww = w * N_inv;
4301b431c67SMark Adams if (w > PETSC_REAL_MIN) {
4311b431c67SMark Adams moments[3] -= ww * PetscLogReal(ww);
4321b431c67SMark Adams } else moments[4] -= w;
4331b431c67SMark Adams }
4341b431c67SMark Adams PetscCall(DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
4351b431c67SMark Adams PetscCall(DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void **)&wp));
4361b431c67SMark Adams }
437d043ef4cSMark Adams PetscCall(PetscInfo(X, "%4d) time %e, Landau particle moments: 0: %18.12e 1: %19.12e 2: %18.12e entropy: %e loss %e. energy = %e + %e + %e\n", (int)n, (double)t, (double)moments[0], (double)moments[1], (double)moments[2], (double)moments[3], (double)(moments[4] / moments[0]), (double)e_grid[0], (double)e_grid[1], (double)e_grid[2]));
4381b431c67SMark Adams }
4391b431c67SMark Adams if (printCtx->print && printCtx->g_target >= 0) {
4401b431c67SMark Adams PetscInt grid = printCtx->g_target, id;
4411b431c67SMark Adams static PetscReal last_t = -100000, period = .5;
4421b431c67SMark Adams if (last_t == -100000) last_t = -period + t;
4431b431c67SMark Adams if (t >= last_t + period) {
4441b431c67SMark Adams last_t = t;
4451b431c67SMark Adams PetscCall(DMGetOutputSequenceNumber(ctx->plex[grid], &id, NULL));
4461b431c67SMark Adams PetscCall(DMSetOutputSequenceNumber(ctx->plex[grid], id + 1, t));
4471b431c67SMark Adams PetscCall(VecViewFromOptions(printCtx->globXArray[LAND_PACK_IDX(v_id % ctx->batch_sz, grid)], NULL, "-ex30_vec_view"));
4481b431c67SMark Adams if (ctx->num_grids > grid + 1) {
4491b431c67SMark Adams PetscCall(DMSetOutputSequenceNumber(ctx->plex[grid + 1], id + 1, t));
4501b431c67SMark Adams PetscCall(VecViewFromOptions(printCtx->globXArray[LAND_PACK_IDX(v_id % ctx->batch_sz, grid + 1)], NULL, "-ex30_vec_view2"));
4511b431c67SMark Adams }
4521b431c67SMark Adams PetscCall(PetscInfo(X, "%4d) time %e View\n", (int)n, (double)t));
4531b431c67SMark Adams }
4541b431c67SMark Adams }
4553cdbf90fSMark Adams PetscCall(DMCompositeRestoreAccessArray(pack, X, nDMs, NULL, printCtx->globXArray));
4563cdbf90fSMark Adams PetscFunctionReturn(PETSC_SUCCESS);
4573cdbf90fSMark Adams }
4583cdbf90fSMark Adams
go(TS ts,Vec X,const PetscInt num_vertices,const PetscInt a_Np,const PetscInt dim,const PetscInt v_target,const PetscInt g_target,PetscReal shift,PetscBool use_uniform_particle_grid)4591b431c67SMark Adams PetscErrorCode go(TS ts, Vec X, const PetscInt num_vertices, const PetscInt a_Np, const PetscInt dim, const PetscInt v_target, const PetscInt g_target, PetscReal shift, PetscBool use_uniform_particle_grid)
460d71ae5a4SJacob Faibussowitsch {
461f3237bfeSMark Adams DM pack, *globSwarmArray, grid_dm[LANDAU_MAX_GRIDS];
462f3237bfeSMark Adams Mat *globMpArray, g_Mass[LANDAU_MAX_GRIDS];
4633cdbf90fSMark Adams KSP t_ksp[LANDAU_MAX_GRIDS][EX30_MAX_NUM_THRDS];
4643cdbf90fSMark Adams Vec t_fhat[LANDAU_MAX_GRIDS][EX30_MAX_NUM_THRDS];
4651b431c67SMark Adams PetscInt nDMs;
4663cdbf90fSMark Adams PetscErrorCode ierr = (PetscErrorCode)0; // used for inside thread loops
467f3237bfeSMark Adams #if defined(PETSC_HAVE_OPENMP) && defined(PETSC_HAVE_THREADSAFETY)
468f3237bfeSMark Adams PetscInt numthreads = PetscNumOMPThreads;
469f3237bfeSMark Adams #else
470f3237bfeSMark Adams PetscInt numthreads = 1;
471f3237bfeSMark Adams #endif
472f3237bfeSMark Adams LandauCtx *ctx;
473f3237bfeSMark Adams Vec *globXArray;
4741b431c67SMark Adams PetscReal moments_0[5], moments_1a[5], moments_1b[5], dt_init;
4753cdbf90fSMark Adams PrintCtx *printCtx;
476f3237bfeSMark Adams
477f3237bfeSMark Adams PetscFunctionBeginUser;
4783cdbf90fSMark 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);
4793cdbf90fSMark Adams PetscCheck(numthreads > 0, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Number threads %" PetscInt_FMT " > %d", numthreads, EX30_MAX_NUM_THRDS);
480f3237bfeSMark Adams PetscCall(TSGetDM(ts, &pack));
481f3237bfeSMark Adams PetscCall(DMGetApplicationContext(pack, &ctx));
482f3237bfeSMark 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);
4833cdbf90fSMark Adams PetscCall(DMCompositeGetNumberDM(pack, &nDMs)); // number of vertices * number of grids
4841b431c67SMark Adams PetscCall(PetscInfo(pack, "Have %" PetscInt_FMT " total grids, with %" PetscInt_FMT " Landau local batched and %" PetscInt_FMT " global items (vertices) %d DMs\n", ctx->num_grids, ctx->batch_sz, num_vertices, (int)nDMs));
485f3237bfeSMark Adams PetscCall(PetscMalloc(sizeof(*globXArray) * nDMs, &globXArray));
486f3237bfeSMark Adams PetscCall(PetscMalloc(sizeof(*globMpArray) * nDMs, &globMpArray));
487f3237bfeSMark Adams PetscCall(PetscMalloc(sizeof(*globSwarmArray) * nDMs, &globSwarmArray));
4883cdbf90fSMark Adams // print ctx
4893cdbf90fSMark Adams PetscCall(PetscNew(&printCtx));
4903cdbf90fSMark Adams PetscCall(TSSetApplicationContext(ts, printCtx));
4913cdbf90fSMark Adams printCtx->v_target = v_target;
4921b431c67SMark Adams printCtx->g_target = g_target;
4933cdbf90fSMark Adams printCtx->ctx = ctx;
4943cdbf90fSMark Adams printCtx->globSwarmArray = globSwarmArray;
4953cdbf90fSMark Adams printCtx->grid_dm = grid_dm;
4963cdbf90fSMark Adams printCtx->globMpArray = globMpArray;
4973cdbf90fSMark Adams printCtx->g_Mass = g_Mass;
4983cdbf90fSMark Adams printCtx->globXArray = globXArray;
4991b431c67SMark Adams printCtx->print_entropy = PETSC_FALSE;
50015229ffcSPierre Jolivet PetscOptionsBegin(PETSC_COMM_SELF, "", "Print Options", "DMPLEX");
5011b431c67SMark Adams PetscCall(PetscOptionsBool("-print_entropy", "Print entropy and moments at each time step", "ex30.c", printCtx->print_entropy, &printCtx->print_entropy, NULL));
5021b431c67SMark Adams PetscOptionsEnd();
5033cdbf90fSMark Adams // view
504f3237bfeSMark Adams PetscCall(DMViewFromOptions(ctx->plex[g_target], NULL, "-ex30_dm_view"));
5053a7d0413SPierre Jolivet if (ctx->num_grids > g_target + 1) PetscCall(DMViewFromOptions(ctx->plex[g_target + 1], NULL, "-ex30_dm_view2"));
5063cdbf90fSMark Adams // create mesh mass matrices
5076b664d00Smarkadams4 PetscCall(VecZeroEntries(X));
508f3237bfeSMark Adams PetscCall(DMCompositeGetAccessArray(pack, X, nDMs, NULL, globXArray)); // just to duplicate
509f3237bfeSMark Adams for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids
510f3237bfeSMark Adams Vec subX = globXArray[LAND_PACK_IDX(0, grid)];
511f3237bfeSMark Adams DM dm = ctx->plex[grid];
512f3237bfeSMark Adams PetscSection s;
513f3237bfeSMark Adams grid_dm[grid] = dm;
514f3237bfeSMark Adams PetscCall(DMCreateMassMatrix(dm, dm, &g_Mass[grid]));
515f3237bfeSMark Adams //
516f3237bfeSMark Adams PetscCall(DMGetLocalSection(dm, &s));
517f3237bfeSMark Adams PetscCall(DMPlexCreateClosureIndex(dm, s));
5186497c311SBarry Smith for (PetscInt tid = 0; tid < numthreads; tid++) {
5191b431c67SMark Adams PC pc;
520f3237bfeSMark Adams PetscCall(VecDuplicate(subX, &t_fhat[grid][tid]));
521f3237bfeSMark Adams PetscCall(KSPCreate(PETSC_COMM_SELF, &t_ksp[grid][tid]));
5221b431c67SMark Adams PetscCall(KSPSetType(t_ksp[grid][tid], KSPCG));
5231b431c67SMark Adams PetscCall(KSPGetPC(t_ksp[grid][tid], &pc));
5241b431c67SMark Adams PetscCall(PCSetType(pc, PCJACOBI));
525f3237bfeSMark Adams PetscCall(KSPSetOptionsPrefix(t_ksp[grid][tid], "ptof_"));
526f3237bfeSMark Adams PetscCall(KSPSetOperators(t_ksp[grid][tid], g_Mass[grid], g_Mass[grid]));
527f3237bfeSMark Adams PetscCall(KSPSetFromOptions(t_ksp[grid][tid]));
528f3237bfeSMark Adams }
529f3237bfeSMark Adams }
530f3237bfeSMark Adams PetscCall(DMCompositeRestoreAccessArray(pack, X, nDMs, NULL, globXArray));
531f3237bfeSMark Adams PetscCall(TSGetTimeStep(ts, &dt_init)); // we could have an adaptive time stepper
5323cdbf90fSMark Adams // loop over all vertices in chucks that are batched for TSSolve
5336497c311SBarry Smith for (PetscInt i = 0; i < 5; i++) moments_0[i] = moments_1a[i] = moments_1b[i] = 0;
5343cdbf90fSMark 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
5356b664d00Smarkadams4 PetscCall(TSSetTime(ts, 0));
5366b664d00Smarkadams4 PetscCall(TSSetStepNumber(ts, 0));
5376b664d00Smarkadams4 PetscCall(TSSetTimeStep(ts, dt_init));
538f3237bfeSMark Adams PetscCall(DMCompositeGetAccessArray(pack, X, nDMs, NULL, globXArray));
5391b431c67SMark Adams printCtx->global_vertex_id_0 = global_vertex_id_0;
5403cdbf90fSMark Adams if (v_target >= global_vertex_id_0 && v_target < global_vertex_id_0 + ctx->batch_sz) {
5413cdbf90fSMark Adams PetscCall(PetscObjectSetName((PetscObject)globXArray[LAND_PACK_IDX(v_target % ctx->batch_sz, g_target)], "rho"));
5423cdbf90fSMark Adams printCtx->print = PETSC_TRUE;
5433cdbf90fSMark Adams } else printCtx->print = PETSC_FALSE;
5443cdbf90fSMark Adams // create fake particles in batches with threads
5453cdbf90fSMark Adams for (PetscInt v_id_0 = 0; v_id_0 < ctx->batch_sz; v_id_0 += numthreads) {
546d043ef4cSMark Adams PetscReal *xx_t[LANDAU_MAX_GRIDS][EX30_MAX_NUM_THRDS], *yy_t[LANDAU_MAX_GRIDS][EX30_MAX_NUM_THRDS], *zz_t[LANDAU_MAX_GRIDS][EX30_MAX_NUM_THRDS], *wp_t[LANDAU_MAX_GRIDS][EX30_MAX_NUM_THRDS] /* , radiuses[80000] */;
5473cdbf90fSMark Adams PetscInt Np_t[LANDAU_MAX_GRIDS][EX30_MAX_NUM_THRDS];
548f3237bfeSMark Adams // make particles
5496497c311SBarry Smith for (PetscInt tid = 0; tid < numthreads; tid++) {
5503cdbf90fSMark Adams const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id;
5513cdbf90fSMark Adams if (glb_v_id < num_vertices) { // the ragged edge (in last batch)
5521b431c67SMark Adams PetscInt Npp0 = a_Np + (glb_v_id % (a_Np / 10 + 1)), nTargetP[LANDAU_MAX_GRIDS]; // n of particels in each dim with load imbalance
553f3237bfeSMark Adams for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids
554d043ef4cSMark Adams // for (PetscInt sp = ctx->species_offset[grid], i0 = 0; sp < ctx->species_offset[grid + 1]; sp++, i0++) {
5553cdbf90fSMark 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 */
556f3237bfeSMark 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
557f3237bfeSMark Adams PetscInt Npi = Npp0, Npj = 2 * Npp0, Npk = 1;
5581b431c67SMark Adams PetscRandom rand;
559d043ef4cSMark Adams PetscReal sigma = ctx->thermal_speed[grid] / ctx->thermal_speed[0], p2_shift = grid == 0 ? shift : -shift; // symmetric shift of e vs ions
5601b431c67SMark Adams PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rand));
5611b431c67SMark Adams PetscCall(PetscRandomSetInterval(rand, 0., 1.));
5621b431c67SMark Adams PetscCall(PetscRandomSetFromOptions(rand));
563f3237bfeSMark Adams if (dim == 2) lo[0] = 0; // Landau coordinate (r,z)
564f3237bfeSMark Adams else Npi = Npj = Npk = Npp0;
5653cdbf90fSMark Adams // User: use glb_v_id to index into your data
566d043ef4cSMark Adams const PetscInt NNreal = Npi * Npj * Npk, NN = NNreal + (dim == 2 ? 3 : 6); // make room for bounding box
567f3237bfeSMark Adams Np_t[grid][tid] = NN;
5683cdbf90fSMark Adams if (glb_v_id == v_target) nTargetP[grid] = NN;
569f3237bfeSMark 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]));
570f3237bfeSMark Adams hp[0] = (hi[0] - lo[0]) / Npi;
571f3237bfeSMark Adams hp[1] = (hi[1] - lo[1]) / Npj;
572f3237bfeSMark Adams hp[2] = (hi[2] - lo[2]) / Npk;
573f3237bfeSMark Adams if (dim == 2) hp[2] = 1;
57463a3b9bcSJacob 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
575f3237bfeSMark Adams vole = hp[0] * hp[1] * hp[2] * ctx->n[grid]; // fix for multi-species
5763cdbf90fSMark 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));
5776497c311SBarry Smith for (PetscInt pj = 0, pp = 0; pj < Npj; pj++) {
5786497c311SBarry Smith for (PetscInt pk = 0; pk < Npk; pk++) {
5796497c311SBarry Smith for (PetscInt pi = 0; pi < Npi; pi++, pp++) {
5801b431c67SMark Adams PetscReal p_shift = p2_shift;
5811b431c67SMark Adams wp_t[grid][tid][pp] = 0;
5821b431c67SMark Adams if (use_uniform_particle_grid) {
583f3237bfeSMark Adams xx_t[grid][tid][pp] = lo[0] + hp[0] / 2.0 + pi * hp[0];
584f3237bfeSMark Adams yy_t[grid][tid][pp] = lo[1] + hp[1] / 2.0 + pj * hp[1];
585f3237bfeSMark Adams if (dim == 3) zz_t[grid][tid][pp] = lo[2] + hp[2] / 2.0 + pk * hp[2];
586f3237bfeSMark Adams PetscReal x[] = {xx_t[grid][tid][pp], yy_t[grid][tid][pp], dim == 2 ? 0 : zz_t[grid][tid][pp]};
5871b431c67SMark Adams p_shift *= ctx->thermal_speed[grid] / ctx->v_0;
588d043ef4cSMark Adams if (ctx->sphere && PetscSqrtReal(PetscSqr(xx_t[grid][tid][pp]) + PetscSqr(yy_t[grid][tid][pp])) > 0.92 * hi[0]) {
589d043ef4cSMark Adams wp_t[grid][tid][pp] = 0;
590d043ef4cSMark Adams } else {
5911b431c67SMark Adams maxwellian(dim, x, kT_m, vole, p_shift, &wp_t[grid][tid][pp]);
5921b431c67SMark Adams if (ctx->num_grids == 1 && shift != 0) { // bi-maxwellian, electron plasma
593d043ef4cSMark Adams maxwellian(dim, x, kT_m, vole, -p_shift, &wp_t[grid][tid][pp]); // symmetric shift of electron plasma
594d043ef4cSMark Adams }
5951b431c67SMark Adams }
5961b431c67SMark Adams } else {
5971b431c67SMark Adams PetscReal u1, u2;
5981b431c67SMark Adams do {
5991b431c67SMark Adams do {
6001b431c67SMark Adams PetscCall(PetscRandomGetValueReal(rand, &u1));
6011b431c67SMark Adams } while (u1 == 0);
6021b431c67SMark Adams PetscCall(PetscRandomGetValueReal(rand, &u2));
6031b431c67SMark Adams //compute z0 and z1
604d043ef4cSMark Adams PetscReal mag = sigma * PetscSqrtReal(-2.0 * PetscLogReal(u1)); // is this the same scale grid Maxwellian? t_therm = sigma
6051b431c67SMark Adams xx_t[grid][tid][pp] = mag * PetscCosReal(2.0 * PETSC_PI * u2);
6061b431c67SMark Adams yy_t[grid][tid][pp] = mag * PetscSinReal(2.0 * PETSC_PI * u2);
6071b431c67SMark Adams if (dim == 2 && xx_t[grid][tid][pp] < lo[0]) xx_t[grid][tid][pp] = -xx_t[grid][tid][pp];
6081b431c67SMark Adams if (dim == 3) zz_t[grid][tid][pp] = lo[2] + hp[2] / 2.0 + pk * hp[2];
609d043ef4cSMark Adams if (!ctx->sphere) {
610d043ef4cSMark Adams if (dim == 2 && xx_t[grid][tid][pp] < 0) xx_t[grid][tid][pp] = -xx_t[grid][tid][pp]; // ???
611d043ef4cSMark Adams else if (dim == 3) {
612d043ef4cSMark Adams while (zz_t[grid][tid][pp] >= hi[2] || zz_t[grid][tid][pp] <= lo[2]) zz_t[grid][tid][pp] *= .9;
613d043ef4cSMark Adams }
614d043ef4cSMark Adams while (xx_t[grid][tid][pp] >= hi[0] || xx_t[grid][tid][pp] <= lo[0]) xx_t[grid][tid][pp] *= .9;
615d043ef4cSMark Adams while (yy_t[grid][tid][pp] >= hi[1] || yy_t[grid][tid][pp] <= lo[1]) yy_t[grid][tid][pp] *= .9;
616d043ef4cSMark Adams } else { // 2D
617d043ef4cSMark Adams //if (glb_v_id == v_target && pp < 80000) radiuses[pp] = PetscSqrtReal(PetscSqr(xx_t[grid][tid][pp]) + PetscSqr(yy_t[grid][tid][pp]));
618d043ef4cSMark Adams while (PetscSqrtReal(PetscSqr(xx_t[grid][tid][pp]) + PetscSqr(yy_t[grid][tid][pp])) > 0.92 * hi[0]) { // safety factor for facets of sphere
619d043ef4cSMark Adams xx_t[grid][tid][pp] *= .9;
620d043ef4cSMark Adams yy_t[grid][tid][pp] *= .9;
621d043ef4cSMark Adams }
622d043ef4cSMark Adams }
6231b431c67SMark Adams if (ctx->num_grids == 1 && pp % 2 == 0) p_shift = 0; // one species, split bi-max
6241b431c67SMark Adams p_shift *= ctx->thermal_speed[grid] / ctx->v_0;
6251b431c67SMark Adams if (dim == 3) zz_t[grid][tid][pp] += p_shift;
6261b431c67SMark Adams else yy_t[grid][tid][pp] += p_shift;
627d043ef4cSMark Adams wp_t[grid][tid][pp] += ctx->n[grid] / NNreal * PetscSqrtReal(ctx->masses[ctx->species_offset[grid]] / ctx->masses[0]);
6281b431c67SMark Adams if (p_shift <= 0) break; // add bi-max for electron plasma only
6291b431c67SMark Adams p_shift = -p_shift;
6301b431c67SMark Adams } while (ctx->num_grids == 1); // add bi-max for electron plasma only
6311b431c67SMark Adams }
6321b431c67SMark Adams {
6333cdbf90fSMark Adams if (glb_v_id == v_target) {
6341b431c67SMark Adams PetscReal x[] = {xx_t[grid][tid][pp], yy_t[grid][tid][pp], dim == 2 ? 0 : zz_t[grid][tid][pp]};
6353cdbf90fSMark 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]];
6366497c311SBarry Smith for (PetscInt i = 0; i < dim; ++i) v2 += PetscSqr(x[i]);
6373cdbf90fSMark Adams moments_0[0] += w; // not thread safe
6383cdbf90fSMark Adams moments_0[1] += w * ctx->v_0 * x[1]; // z-momentum
6391b431c67SMark Adams moments_0[2] += w * 0.5 * ctx->v_0 * ctx->v_0 * v2;
640f3237bfeSMark Adams }
641f3237bfeSMark Adams }
642f3237bfeSMark Adams }
643f3237bfeSMark Adams }
644f3237bfeSMark Adams }
645d043ef4cSMark Adams if (dim == 2) { // fix bounding box
6466497c311SBarry Smith PetscInt pp = NNreal;
647d043ef4cSMark Adams wp_t[grid][tid][pp] = 0;
648d043ef4cSMark Adams xx_t[grid][tid][pp] = 1.e-7;
649d043ef4cSMark Adams yy_t[grid][tid][pp++] = hi[1] - 5.e-7;
650d043ef4cSMark Adams wp_t[grid][tid][pp] = 0;
651d043ef4cSMark Adams xx_t[grid][tid][pp] = hi[0] - 5.e-7;
652d043ef4cSMark Adams yy_t[grid][tid][pp++] = 0;
653d043ef4cSMark Adams wp_t[grid][tid][pp] = 0;
654d043ef4cSMark Adams xx_t[grid][tid][pp] = 1.e-7;
655d043ef4cSMark Adams yy_t[grid][tid][pp++] = lo[1] + 5.e-7;
656d043ef4cSMark Adams } else {
6576497c311SBarry Smith const PetscInt p0 = NNreal;
658ac530a7eSPierre Jolivet for (PetscInt pj = 0; pj < 6; pj++) xx_t[grid][tid][p0 + pj] = yy_t[grid][tid][p0 + pj] = zz_t[grid][tid][p0 + pj] = wp_t[grid][tid][p0 + pj] = 0;
659d043ef4cSMark Adams xx_t[grid][tid][p0 + 0] = lo[0];
660d043ef4cSMark Adams xx_t[grid][tid][p0 + 1] = hi[0];
661d043ef4cSMark Adams yy_t[grid][tid][p0 + 2] = lo[1];
662d043ef4cSMark Adams yy_t[grid][tid][p0 + 3] = hi[1];
663d043ef4cSMark Adams zz_t[grid][tid][p0 + 4] = lo[2];
664d043ef4cSMark Adams zz_t[grid][tid][p0 + 5] = hi[2];
665d043ef4cSMark Adams }
6661b431c67SMark Adams PetscCall(PetscRandomDestroy(&rand));
6673cdbf90fSMark Adams }
6681b431c67SMark Adams // entropy init, need global n
6693cdbf90fSMark Adams if (glb_v_id == v_target) {
6701b431c67SMark Adams const PetscReal N_inv = 1 / moments_0[0];
6713cdbf90fSMark Adams PetscCall(PetscInfo(pack, "Target %" PetscInt_FMT " with %" PetscInt_FMT " particels\n", glb_v_id, nTargetP[0]));
6723cdbf90fSMark Adams for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
6731b431c67SMark Adams const PetscInt NN = nTargetP[grid];
6746497c311SBarry Smith for (PetscInt pp = 0; pp < NN; pp++) {
6751b431c67SMark Adams const PetscReal fact = dim == 2 ? 2.0 * PETSC_PI * xx_t[grid][tid][pp] : 1, w = fact * ctx->n_0 * ctx->masses[ctx->species_offset[grid]] * wp_t[grid][tid][pp], ww = w * N_inv;
6761b431c67SMark Adams if (w > PETSC_REAL_MIN) {
6773cdbf90fSMark Adams moments_0[3] -= ww * PetscLogReal(ww);
6783cdbf90fSMark Adams PetscCheck(ww < 1 - PETSC_MACHINE_EPSILON, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "ww (%g) > 1", (double)ww);
6791b431c67SMark Adams } else moments_0[4] -= w;
6803cdbf90fSMark Adams }
681f3237bfeSMark Adams } // grid
6821b431c67SMark Adams } // target
683f3237bfeSMark Adams } // active
6843cdbf90fSMark Adams } // threads
685f3237bfeSMark Adams /* Create particle swarm */
6866497c311SBarry Smith for (PetscInt tid = 0; tid < numthreads; tid++) {
6873cdbf90fSMark Adams const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id;
6883cdbf90fSMark Adams if (glb_v_id < num_vertices) { // the ragged edge of the last batch
689f3237bfeSMark Adams for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids
690f3237bfeSMark Adams PetscSection section;
691f3237bfeSMark Adams PetscInt Nf;
692f3237bfeSMark Adams DM dm = grid_dm[grid];
693d043ef4cSMark Adams PetscCall(DMGetLocalSection(dm, §ion));
694d043ef4cSMark Adams PetscCall(PetscSectionGetNumFields(section, &Nf));
695bd158744SPierre Jolivet PetscCheck(Nf == 1, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Only one species per grid supported -- todo");
696d043ef4cSMark Adams PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
697d043ef4cSMark Adams PetscCall(PetscInfo(pack, "call createSwarm [%" PetscInt_FMT ".%" PetscInt_FMT "] local block index %" PetscInt_FMT "\n", v_id, grid, LAND_PACK_IDX(v_id, grid)));
698d043ef4cSMark Adams PetscCall(createSwarm(dm, dim, &globSwarmArray[LAND_PACK_IDX(v_id, grid)]));
699f3237bfeSMark Adams }
700f3237bfeSMark Adams } // active
7013cdbf90fSMark Adams } // threads
702d043ef4cSMark Adams PetscCheck(ierr != 9999, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Only support one species per grid");
7033cdbf90fSMark Adams // make globMpArray
704f3237bfeSMark Adams PetscPragmaOMP(parallel for)
7056497c311SBarry Smith for (PetscInt tid = 0; tid < numthreads; tid++) {
7063cdbf90fSMark Adams const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id;
7073cdbf90fSMark Adams if (glb_v_id < num_vertices) {
7083cdbf90fSMark Adams for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids
709d043ef4cSMark Adams // for (PetscInt sp = ctx->species_offset[grid], i0 = 0; sp < ctx->species_offset[grid + 1]; sp++, i0++) -- loop over species for Nf > 1 -- TODO
7103cdbf90fSMark Adams PetscErrorCode ierr_t;
7113cdbf90fSMark Adams DM sw = globSwarmArray[LAND_PACK_IDX(v_id, grid)];
7121b431c67SMark Adams ierr_t = PetscInfo(pack, "makeSwarm %" PetscInt_FMT ".%" PetscInt_FMT ") for block %" PetscInt_FMT "\n", v_id, grid, LAND_PACK_IDX(v_id, grid));
7133cdbf90fSMark Adams ierr_t = makeSwarm(sw, dim, Np_t[grid][tid], xx_t[grid][tid], yy_t[grid][tid], zz_t[grid][tid]);
7143cdbf90fSMark Adams if (ierr_t) ierr = ierr_t;
7153cdbf90fSMark Adams }
7163cdbf90fSMark Adams }
7173cdbf90fSMark Adams }
7186497c311SBarry Smith for (PetscInt tid = 0; tid < numthreads; tid++) {
7193cdbf90fSMark Adams const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id;
7203cdbf90fSMark Adams if (glb_v_id < num_vertices) {
721f3237bfeSMark Adams for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids
722f3237bfeSMark Adams DM dm = grid_dm[grid];
7233cdbf90fSMark Adams DM sw = globSwarmArray[LAND_PACK_IDX(v_id, grid)];
724d043ef4cSMark Adams PetscCall(PetscInfo(pack, "createMp %" PetscInt_FMT ".%" PetscInt_FMT ") for block %" PetscInt_FMT "\n", v_id, grid, LAND_PACK_IDX(v_id, grid)));
725d043ef4cSMark Adams PetscCall(createMp(dm, sw, &globMpArray[LAND_PACK_IDX(v_id, grid)]));
726d043ef4cSMark Adams PetscCall(MatViewFromOptions(globMpArray[LAND_PACK_IDX(v_id, grid)], NULL, "-mp_mat_view"));
7273cdbf90fSMark Adams }
7283cdbf90fSMark Adams }
7293cdbf90fSMark Adams }
7303cdbf90fSMark Adams // p --> g: set X
7313cdbf90fSMark Adams // PetscPragmaOMP(parallel for)
7326497c311SBarry Smith for (PetscInt tid = 0; tid < numthreads; tid++) {
7333cdbf90fSMark Adams const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id;
7343cdbf90fSMark Adams if (glb_v_id < num_vertices) {
7353cdbf90fSMark Adams for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids
7363cdbf90fSMark Adams PetscErrorCode ierr_t;
7373cdbf90fSMark Adams DM dm = grid_dm[grid];
7383cdbf90fSMark Adams DM sw = globSwarmArray[LAND_PACK_IDX(v_id, grid)];
7393cdbf90fSMark Adams Vec subX = globXArray[LAND_PACK_IDX(v_id, grid)], work = t_fhat[grid][tid];
7401b431c67SMark Adams ierr_t = PetscInfo(pack, "particlesToGrid %" PetscInt_FMT ".%" PetscInt_FMT ") for block %" PetscInt_FMT "\n", v_id, grid, LAND_PACK_IDX(v_id, grid));
7411b431c67SMark Adams ierr_t = particlesToGrid(dm, sw, tid, dim, wp_t[grid][tid], subX, globMpArray[LAND_PACK_IDX(v_id, grid)]);
742f3237bfeSMark Adams if (ierr_t) ierr = ierr_t;
743f3237bfeSMark Adams // u = M^_1 f_w
744f3237bfeSMark Adams ierr_t = VecCopy(subX, work);
745f3237bfeSMark Adams ierr_t = KSPSolve(t_ksp[grid][tid], work, subX);
746f3237bfeSMark Adams if (ierr_t) ierr = ierr_t;
747f3237bfeSMark Adams }
748f3237bfeSMark Adams }
749f3237bfeSMark Adams }
750d043ef4cSMark Adams PetscCheck(!ierr, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Error in OMP loop. ierr = %d", (int)ierr);
751f3237bfeSMark Adams /* Cleanup */
7526497c311SBarry Smith for (PetscInt tid = 0; tid < numthreads; tid++) {
7533cdbf90fSMark Adams const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id;
7543cdbf90fSMark Adams if (glb_v_id < num_vertices) {
755f3237bfeSMark Adams for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids
756f3237bfeSMark Adams PetscCall(PetscFree4(xx_t[grid][tid], yy_t[grid][tid], wp_t[grid][tid], zz_t[grid][tid]));
757f3237bfeSMark Adams }
758f3237bfeSMark Adams } // active
7593cdbf90fSMark Adams } // threads
7603cdbf90fSMark Adams } // (fake) particle loop
7613cdbf90fSMark Adams // standard view of initial conditions
7623cdbf90fSMark Adams if (v_target >= global_vertex_id_0 && v_target < global_vertex_id_0 + ctx->batch_sz) {
7633cdbf90fSMark Adams PetscCall(DMSetOutputSequenceNumber(ctx->plex[g_target], 0, 0.0));
7643cdbf90fSMark Adams PetscCall(VecViewFromOptions(globXArray[LAND_PACK_IDX(v_target % ctx->batch_sz, g_target)], NULL, "-ex30_vec_view"));
7651b431c67SMark Adams if (ctx->num_grids > g_target + 1) {
7661b431c67SMark Adams PetscCall(DMSetOutputSequenceNumber(ctx->plex[g_target + 1], 0, 0.0));
7671b431c67SMark Adams PetscCall(VecViewFromOptions(globXArray[LAND_PACK_IDX(v_target % ctx->batch_sz, g_target + 1)], NULL, "-ex30_vec_view2"));
7681b431c67SMark Adams }
769d043ef4cSMark Adams PetscCall(MatViewFromOptions(globMpArray[LAND_PACK_IDX(v_target % ctx->batch_sz, g_target)], NULL, "-ex30_mass_mat_view"));
770d043ef4cSMark Adams PetscCall(DMViewFromOptions(globSwarmArray[LAND_PACK_IDX(v_target % ctx->batch_sz, g_target)], NULL, "-ex30_sw_view"));
771d043ef4cSMark Adams PetscCall(DMSwarmViewXDMF(globSwarmArray[LAND_PACK_IDX(v_target % ctx->batch_sz, g_target)], "initial_swarm.xmf")); // writes a file by default!!!
772f3237bfeSMark Adams }
773d043ef4cSMark Adams // coarse graining moments_1a, bring f back from grid before advance
7741b431c67SMark Adams if (v_target >= global_vertex_id_0 && v_target < global_vertex_id_0 + ctx->batch_sz && printCtx->print_entropy) {
7753cdbf90fSMark Adams const PetscInt v_id = v_target % ctx->batch_sz;
776f3237bfeSMark Adams for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
777f3237bfeSMark Adams PetscDataType dtype;
778f3237bfeSMark Adams PetscReal *wp, *coords;
7793cdbf90fSMark Adams DM sw = globSwarmArray[LAND_PACK_IDX(v_id, grid)];
7803cdbf90fSMark Adams Vec work, subX = globXArray[LAND_PACK_IDX(v_id, grid)];
7811b431c67SMark Adams PetscInt bs, NN;
7823cdbf90fSMark Adams // C-G moments
7833cdbf90fSMark Adams PetscCall(VecDuplicate(subX, &work));
7843cdbf90fSMark Adams PetscCall(gridToParticles(grid_dm[grid], sw, subX, work, globMpArray[LAND_PACK_IDX(v_id, grid)], g_Mass[grid]));
7853cdbf90fSMark Adams PetscCall(VecDestroy(&work));
7863cdbf90fSMark Adams // moments
787f3237bfeSMark Adams PetscCall(DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
7881b431c67SMark Adams PetscCall(DMSwarmGetLocalSize(sw, &NN));
7891b431c67SMark Adams PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void **)&wp));
7906497c311SBarry Smith for (PetscInt pp = 0; pp < NN; pp++) {
7911b431c67SMark 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]];
7926497c311SBarry Smith for (PetscInt i = 0; i < dim; ++i) v2 += PetscSqr(coords[pp * dim + i]);
7933cdbf90fSMark Adams moments_1a[0] += w;
7943cdbf90fSMark Adams moments_1a[1] += w * ctx->v_0 * coords[pp * dim + 1]; // z-momentum
7951b431c67SMark Adams moments_1a[2] += w * 0.5 * ctx->v_0 * ctx->v_0 * v2;
796f3237bfeSMark Adams }
797f3237bfeSMark Adams PetscCall(DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
798f3237bfeSMark Adams PetscCall(DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void **)&wp));
799f3237bfeSMark Adams }
8001b431c67SMark Adams // entropy
8011b431c67SMark Adams const PetscReal N_inv = 1 / moments_1a[0];
8021b431c67SMark Adams PetscCall(PetscInfo(pack, "Entropy batch %" PetscInt_FMT " of %" PetscInt_FMT ", n = %g\n", v_target, num_vertices, (double)(1 / N_inv)));
8031b431c67SMark Adams for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
8041b431c67SMark Adams PetscDataType dtype;
8051b431c67SMark Adams PetscReal *wp, *coords;
8061b431c67SMark Adams DM sw = globSwarmArray[LAND_PACK_IDX(v_id, grid)];
8071b431c67SMark Adams PetscInt bs, NN;
8081b431c67SMark Adams PetscCall(DMSwarmGetLocalSize(sw, &NN));
8091b431c67SMark Adams PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void **)&wp));
8101b431c67SMark Adams PetscCall(DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
8116497c311SBarry Smith for (PetscInt pp = 0; pp < NN; pp++) {
8121b431c67SMark Adams PetscReal fact = (dim == 2) ? 2.0 * PETSC_PI * coords[pp * dim + 0] : 1, w = fact * wp[pp] * ctx->n_0 * ctx->masses[ctx->species_offset[grid]], ww = w * N_inv;
8131b431c67SMark Adams if (w > PETSC_REAL_MIN) {
8141b431c67SMark Adams moments_1a[3] -= ww * PetscLogReal(ww);
8151b431c67SMark Adams PetscCheck(ww < 1 - PETSC_MACHINE_EPSILON, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "ww (%g) > 1", (double)ww);
8161b431c67SMark Adams } else moments_1a[4] -= w;
8171b431c67SMark Adams }
8181b431c67SMark Adams PetscCall(DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void **)&wp));
8191b431c67SMark Adams PetscCall(DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
8201b431c67SMark Adams }
821f3237bfeSMark Adams }
8223cdbf90fSMark Adams // restore vector
823f3237bfeSMark Adams PetscCall(DMCompositeRestoreAccessArray(pack, X, nDMs, NULL, globXArray));
824d043ef4cSMark Adams // view initial grid
8253a7d0413SPierre Jolivet if (v_target >= global_vertex_id_0 && v_target < global_vertex_id_0 + ctx->batch_sz) PetscCall(DMPlexLandauPrintNorms(X, 0));
8263cdbf90fSMark Adams // advance
8273cdbf90fSMark Adams PetscCall(TSSetSolution(ts, X));
8281b431c67SMark Adams PetscCall(PetscInfo(pack, "Advance vertex %" PetscInt_FMT " to %" PetscInt_FMT "\n", global_vertex_id_0, global_vertex_id_0 + ctx->batch_sz));
8293cdbf90fSMark Adams PetscCall(TSSetPostStep(ts, PostStep));
8303cdbf90fSMark Adams PetscCall(PostStep(ts));
8313cdbf90fSMark Adams PetscCall(TSSolve(ts, X));
8323cdbf90fSMark Adams // view
8333cdbf90fSMark Adams PetscCall(DMCompositeGetAccessArray(pack, X, nDMs, NULL, globXArray));
8343cdbf90fSMark Adams if (v_target >= global_vertex_id_0 && v_target < global_vertex_id_0 + ctx->batch_sz) {
835d043ef4cSMark Adams /* Visualize original particle field */
8361b431c67SMark Adams DM sw = globSwarmArray[LAND_PACK_IDX(v_target % ctx->batch_sz, g_target)];
8371b431c67SMark Adams Vec f;
8381b431c67SMark Adams PetscCall(DMSetOutputSequenceNumber(sw, 0, 0.0));
8391b431c67SMark Adams PetscCall(DMViewFromOptions(grid_dm[g_target], NULL, "-weights_dm_view"));
8401b431c67SMark Adams PetscCall(DMViewFromOptions(sw, NULL, "-weights_sw_view"));
8411b431c67SMark Adams PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &f));
8421b431c67SMark Adams PetscCall(PetscObjectSetName((PetscObject)f, "weights"));
8431b431c67SMark Adams PetscCall(VecViewFromOptions(f, NULL, "-weights_vec_view"));
8441b431c67SMark Adams PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &f));
845d043ef4cSMark Adams //
846d043ef4cSMark Adams PetscCall(DMPlexLandauPrintNorms(X, 1));
8473cdbf90fSMark Adams }
848d043ef4cSMark Adams if (!use_uniform_particle_grid) { // resample to uniform grid
849d043ef4cSMark Adams for (PetscInt v_id_0 = 0; v_id_0 < ctx->batch_sz; v_id_0 += numthreads) {
850d043ef4cSMark 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];
851d043ef4cSMark Adams PetscInt Np_t[LANDAU_MAX_GRIDS][EX30_MAX_NUM_THRDS];
8526497c311SBarry Smith for (PetscInt tid = 0; tid < numthreads; tid++) {
853d043ef4cSMark Adams const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id;
854d043ef4cSMark Adams if (glb_v_id < num_vertices) {
855d043ef4cSMark Adams // create uniform grid w/o weights & smaller
856d043ef4cSMark Adams PetscInt Npp0 = (a_Np + (glb_v_id % (a_Np / 10 + 1))) / 2, Nv; // 1/2 of uniform particle grid size
857d043ef4cSMark Adams for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
858d043ef4cSMark Adams // for (PetscInt sp = ctx->species_offset[grid], i0 = 0; sp < ctx->species_offset[grid + 1]; sp++, i0++)
859d043ef4cSMark 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];
860d043ef4cSMark Adams PetscInt Npi = Npp0, Npj = 2 * Npp0, Npk = 1, NN;
861d043ef4cSMark Adams // delete old particles and particle mass matrix
862d043ef4cSMark Adams PetscCall(DMDestroy(&globSwarmArray[LAND_PACK_IDX(v_id, grid)]));
863d043ef4cSMark Adams PetscCall(MatDestroy(&globMpArray[LAND_PACK_IDX(v_id, grid)]));
864d043ef4cSMark Adams // create fake particles in batches with threads
865d043ef4cSMark Adams PetscCall(MatGetLocalSize(g_Mass[grid], &Nv, NULL));
866d043ef4cSMark Adams if (dim == 2) lo[0] = 0;
867d043ef4cSMark Adams else Npi = Npj = Npk = Npp0;
868d043ef4cSMark Adams NN = Npi * Npj * Npk + (dim == 2 ? 3 : 6); // make a regular grid of particles Npp x Npp
869d043ef4cSMark Adams while (Npi * Npj * Npk < Nv) { // make stable - no LS
870d043ef4cSMark Adams Npi++;
871d043ef4cSMark Adams Npj++;
872d043ef4cSMark Adams Npk++;
873d043ef4cSMark Adams NN = Npi * Npj * Npk + (dim == 2 ? 3 : 6);
874d043ef4cSMark Adams }
875d043ef4cSMark Adams Np_t[grid][tid] = NN;
876d043ef4cSMark 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]));
877d043ef4cSMark Adams hp[0] = (hi[0] - lo[0]) / Npi;
878d043ef4cSMark Adams hp[1] = (hi[1] - lo[1]) / Npj;
879d043ef4cSMark Adams hp[2] = (hi[2] - lo[2]) / Npk;
880d043ef4cSMark Adams if (dim == 2) hp[2] = 1;
881d043ef4cSMark Adams PetscCall(PetscInfo(pack, "Resampling %d particles, %d vertices\n", (int)NN, (int)Nv)); // temp
8826497c311SBarry Smith for (PetscInt pj = 0, pp = 0; pj < Npj; pj++) {
8836497c311SBarry Smith for (PetscInt pk = 0; pk < Npk; pk++) {
8846497c311SBarry Smith for (PetscInt pi = 0; pi < Npi; pi++, pp++) {
885d043ef4cSMark Adams wp_t[grid][tid][pp] = 0;
886d043ef4cSMark Adams xx_t[grid][tid][pp] = lo[0] + hp[0] / 2.0 + pi * hp[0];
887d043ef4cSMark Adams yy_t[grid][tid][pp] = lo[1] + hp[1] / 2.0 + pj * hp[1];
888d043ef4cSMark Adams if (dim == 3) zz_t[grid][tid][pp] = lo[2] + hp[2] / 2.0 + pk * hp[2];
889d043ef4cSMark Adams }
890d043ef4cSMark Adams }
891d043ef4cSMark Adams }
892d043ef4cSMark Adams if (dim == 2) { // fix bounding box
8936497c311SBarry Smith PetscInt pp = NN - 3;
894d043ef4cSMark Adams wp_t[grid][tid][pp] = 0;
895d043ef4cSMark Adams xx_t[grid][tid][pp] = 1.e-7;
896d043ef4cSMark Adams yy_t[grid][tid][pp++] = hi[1] - 5.e-7;
897d043ef4cSMark Adams wp_t[grid][tid][pp] = 0;
898d043ef4cSMark Adams xx_t[grid][tid][pp] = hi[0] - 5.e-7;
899d043ef4cSMark Adams yy_t[grid][tid][pp++] = 0;
900d043ef4cSMark Adams wp_t[grid][tid][pp] = 0;
901d043ef4cSMark Adams xx_t[grid][tid][pp] = 1.e-7;
902d043ef4cSMark Adams yy_t[grid][tid][pp++] = lo[1] + 5.e-7;
903d043ef4cSMark Adams } else {
9046497c311SBarry Smith const PetscInt p0 = NN - 6;
905ac530a7eSPierre Jolivet for (PetscInt pj = 0; pj < 6; pj++) xx_t[grid][tid][p0 + pj] = yy_t[grid][tid][p0 + pj] = zz_t[grid][tid][p0 + pj] = wp_t[grid][tid][p0 + pj] = 0;
906d043ef4cSMark Adams xx_t[grid][tid][p0 + 0] = lo[0];
907d043ef4cSMark Adams xx_t[grid][tid][p0 + 1] = hi[0];
908d043ef4cSMark Adams yy_t[grid][tid][p0 + 2] = lo[1];
909d043ef4cSMark Adams yy_t[grid][tid][p0 + 3] = hi[1];
910d043ef4cSMark Adams zz_t[grid][tid][p0 + 4] = lo[2];
911d043ef4cSMark Adams zz_t[grid][tid][p0 + 5] = hi[2];
912d043ef4cSMark Adams }
913d043ef4cSMark Adams }
914d043ef4cSMark Adams } // active
915d043ef4cSMark Adams } // threads
916d043ef4cSMark Adams /* Create particle swarm */
9176497c311SBarry Smith for (PetscInt tid = 0; tid < numthreads; tid++) {
918d043ef4cSMark Adams const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id;
919d043ef4cSMark Adams if (glb_v_id < num_vertices) { // the ragged edge of the last batch
920d043ef4cSMark Adams for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids
921d043ef4cSMark Adams // for (PetscInt sp = ctx->species_offset[grid], i0 = 0; sp < ctx->species_offset[grid + 1]; sp++, i0++) -- loop over species for Nf > 1 -- TODO
922d043ef4cSMark Adams PetscErrorCode ierr_t;
923d043ef4cSMark Adams PetscSection section;
924d043ef4cSMark Adams PetscInt Nf;
925d043ef4cSMark Adams DM dm = grid_dm[grid];
926d043ef4cSMark Adams ierr_t = DMGetLocalSection(dm, §ion);
927d043ef4cSMark Adams ierr_t = PetscSectionGetNumFields(section, &Nf);
928d043ef4cSMark Adams if (Nf != 1) ierr_t = (PetscErrorCode)9999;
929d043ef4cSMark Adams else {
930d043ef4cSMark Adams ierr_t = DMViewFromOptions(dm, NULL, "-dm_view");
931d043ef4cSMark Adams ierr_t = PetscInfo(pack, "call createSwarm [%" PetscInt_FMT ".%" PetscInt_FMT "] local block index %" PetscInt_FMT "\n", v_id, grid, LAND_PACK_IDX(v_id, grid));
932d043ef4cSMark Adams ierr_t = createSwarm(dm, dim, &globSwarmArray[LAND_PACK_IDX(v_id, grid)]);
933d043ef4cSMark Adams }
934d043ef4cSMark Adams if (ierr_t) ierr = ierr_t;
935d043ef4cSMark Adams }
936d043ef4cSMark Adams } // active
937d043ef4cSMark Adams } // threads
938d043ef4cSMark Adams PetscCheck(ierr != 9999, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Only support one species per grid");
939d043ef4cSMark Adams PetscCheck(!ierr, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Error in OMP loop. ierr = %d", (int)ierr);
940d043ef4cSMark Adams // make globMpArray
941d043ef4cSMark Adams PetscPragmaOMP(parallel for)
9426497c311SBarry Smith for (PetscInt tid = 0; tid < numthreads; tid++) {
943d043ef4cSMark Adams const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id;
944d043ef4cSMark Adams if (glb_v_id < num_vertices) {
945d043ef4cSMark Adams for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids
946d043ef4cSMark Adams // for (PetscInt sp = ctx->species_offset[grid], i0 = 0; sp < ctx->species_offset[grid + 1]; sp++, i0++) -- loop over species for Nf > 1 -- TODO
947d043ef4cSMark Adams PetscErrorCode ierr_t;
948d043ef4cSMark Adams DM sw = globSwarmArray[LAND_PACK_IDX(v_id, grid)];
949d043ef4cSMark Adams ierr_t = PetscInfo(pack, "makeSwarm %" PetscInt_FMT ".%" PetscInt_FMT ") for block %" PetscInt_FMT "\n", v_id, grid, LAND_PACK_IDX(v_id, grid));
950d043ef4cSMark Adams ierr_t = makeSwarm(sw, dim, Np_t[grid][tid], xx_t[grid][tid], yy_t[grid][tid], zz_t[grid][tid]);
951d043ef4cSMark Adams if (ierr_t) ierr = ierr_t;
952d043ef4cSMark Adams }
953d043ef4cSMark Adams } // active
954d043ef4cSMark Adams } // threads
955d043ef4cSMark Adams // create particle mass matrices
956d043ef4cSMark Adams //PetscPragmaOMP(parallel for)
9576497c311SBarry Smith for (PetscInt tid = 0; tid < numthreads; tid++) {
958d043ef4cSMark Adams const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id;
959d043ef4cSMark Adams if (glb_v_id < num_vertices) {
960d043ef4cSMark Adams for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids
961d043ef4cSMark Adams PetscErrorCode ierr_t;
962d043ef4cSMark Adams DM dm = grid_dm[grid];
963d043ef4cSMark Adams DM sw = globSwarmArray[LAND_PACK_IDX(v_id, grid)];
964d043ef4cSMark Adams ierr_t = PetscInfo(pack, "createMp %" PetscInt_FMT ".%" PetscInt_FMT ") for block %" PetscInt_FMT "\n", v_id, grid, LAND_PACK_IDX(v_id, grid));
965d043ef4cSMark Adams ierr_t = createMp(dm, sw, &globMpArray[LAND_PACK_IDX(v_id, grid)]);
966d043ef4cSMark Adams if (ierr_t) ierr = ierr_t;
967d043ef4cSMark Adams }
968d043ef4cSMark Adams } // active
969d043ef4cSMark Adams } // threads
970d043ef4cSMark Adams PetscCheck(!ierr, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Error in OMP loop. ierr = %d", (int)ierr);
971d043ef4cSMark Adams /* Cleanup */
9726497c311SBarry Smith for (PetscInt tid = 0; tid < numthreads; tid++) {
973d043ef4cSMark Adams const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id;
974d043ef4cSMark Adams if (glb_v_id < num_vertices) {
975d043ef4cSMark Adams for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids
976d043ef4cSMark Adams PetscCall(PetscFree4(xx_t[grid][tid], yy_t[grid][tid], wp_t[grid][tid], zz_t[grid][tid]));
977d043ef4cSMark Adams }
978d043ef4cSMark Adams } // active
979d043ef4cSMark Adams } // threads
980d043ef4cSMark Adams } // batch
981d043ef4cSMark Adams // view
982d043ef4cSMark Adams if (v_target >= global_vertex_id_0 && v_target < global_vertex_id_0 + ctx->batch_sz) {
983d043ef4cSMark Adams /* Visualize particle field */
984d043ef4cSMark Adams DM sw = globSwarmArray[LAND_PACK_IDX(v_target % ctx->batch_sz, g_target)];
985d043ef4cSMark Adams Vec f;
986d043ef4cSMark Adams PetscCall(DMSetOutputSequenceNumber(sw, 0, 0.0));
987d043ef4cSMark Adams PetscCall(DMViewFromOptions(sw, NULL, "-resampled_weights_sw_view"));
988d043ef4cSMark Adams PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &f));
989d043ef4cSMark Adams PetscCall(PetscObjectSetName((PetscObject)f, "resampled_weights"));
990d043ef4cSMark Adams PetscCall(VecViewFromOptions(f, NULL, "-resampled_weights_vec_view"));
991d043ef4cSMark Adams PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &f));
992d043ef4cSMark Adams PetscCall(DMSwarmViewXDMF(sw, "resampled.xmf"));
993d043ef4cSMark Adams }
994d043ef4cSMark Adams } // !uniform
995d043ef4cSMark Adams // particles to grid, compute moments and entropy, for target vertex only
9961b431c67SMark Adams if (v_target >= global_vertex_id_0 && v_target < global_vertex_id_0 + ctx->batch_sz && printCtx->print_entropy) {
997d043ef4cSMark Adams PetscReal energy_error_rel;
9983cdbf90fSMark 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));
99957508eceSPierre Jolivet energy_error_rel = PetscAbsReal(moments_1b[2] - moments_0[2]) / moments_0[2];
1000d043ef4cSMark Adams PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Particle Moments:\t number density momentum (par) energy entropy negative weights : # OMP threads %g\n", (double)numthreads));
1001d043ef4cSMark Adams PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\tInitial: %18.12e %19.12e %18.12e %18.12e %g %%\n", (double)moments_0[0], (double)moments_0[1], (double)moments_0[2], (double)moments_0[3], 100 * (double)(moments_0[4] / moments_0[0])));
1002d043ef4cSMark Adams PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\tCoarse-graining: %18.12e %19.12e %18.12e %18.12e %g %%\n", (double)moments_1a[0], (double)moments_1a[1], (double)moments_1a[2], (double)moments_1a[3], 100 * (double)(moments_1a[4] / moments_0[0])));
1003d043ef4cSMark Adams PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\tLandau: %18.12e %19.12e %18.12e %18.12e %g %%\n", (double)moments_1b[0], (double)moments_1b[1], (double)moments_1b[2], (double)moments_1b[3], 100 * (double)(moments_1b[4] / moments_0[0])));
10041b431c67SMark Adams PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Coarse-graining entropy generation = %e ; Landau entropy generation = %e\n", (double)(moments_1a[3] - moments_0[3]), (double)(moments_1b[3] - moments_0[3])));
1005d043ef4cSMark Adams PetscCall(PetscPrintf(PETSC_COMM_WORLD, "(relative) energy conservation: Coarse-graining = %e, Landau = %e (%g %d)\n", (double)(moments_1a[2] - moments_0[2]) / (double)moments_0[2], (double)energy_error_rel, (double)PetscLog10Real(energy_error_rel), (int)(PetscLog10Real(energy_error_rel) + .5)));
10061b431c67SMark Adams }
10073cdbf90fSMark Adams // restore vector
10083cdbf90fSMark Adams PetscCall(DMCompositeRestoreAccessArray(pack, X, nDMs, NULL, globXArray));
10093cdbf90fSMark Adams // cleanup
10103cdbf90fSMark Adams for (PetscInt v_id_0 = 0; v_id_0 < ctx->batch_sz; v_id_0 += numthreads) {
10116497c311SBarry Smith for (PetscInt tid = 0; tid < numthreads; tid++) {
10123cdbf90fSMark Adams const PetscInt v_id = v_id_0 + tid, glb_v_id = global_vertex_id_0 + v_id;
10133cdbf90fSMark Adams if (glb_v_id < num_vertices) {
10143cdbf90fSMark Adams for (PetscInt grid = 0; grid < ctx->num_grids; grid++) {
10153cdbf90fSMark Adams PetscCall(DMDestroy(&globSwarmArray[LAND_PACK_IDX(v_id, grid)]));
10163cdbf90fSMark Adams PetscCall(MatDestroy(&globMpArray[LAND_PACK_IDX(v_id, grid)]));
10173cdbf90fSMark Adams }
10183cdbf90fSMark Adams }
10193cdbf90fSMark Adams }
10203cdbf90fSMark Adams }
10211b431c67SMark Adams } // user batch, not used
1022f3237bfeSMark Adams /* Cleanup */
1023f3237bfeSMark Adams PetscCall(PetscFree(globXArray));
1024f3237bfeSMark Adams PetscCall(PetscFree(globSwarmArray));
1025f3237bfeSMark Adams PetscCall(PetscFree(globMpArray));
10263cdbf90fSMark Adams PetscCall(PetscFree(printCtx));
1027f3237bfeSMark Adams // clean up mass matrices
1028f3237bfeSMark Adams for (PetscInt grid = 0; grid < ctx->num_grids; grid++) { // add same particels for all grids
1029f3237bfeSMark Adams PetscCall(MatDestroy(&g_Mass[grid]));
10306497c311SBarry Smith for (PetscInt tid = 0; tid < numthreads; tid++) {
1031f3237bfeSMark Adams PetscCall(VecDestroy(&t_fhat[grid][tid]));
1032f3237bfeSMark Adams PetscCall(KSPDestroy(&t_ksp[grid][tid]));
1033f3237bfeSMark Adams }
1034f3237bfeSMark Adams }
10353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1036f3237bfeSMark Adams }
1037f3237bfeSMark Adams
main(int argc,char ** argv)1038d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
1039d71ae5a4SJacob Faibussowitsch {
1040f3237bfeSMark Adams DM pack;
1041f3237bfeSMark Adams Vec X;
10421b431c67SMark Adams PetscInt dim = 2, num_vertices = 1, Np = 10, v_target = 0, g_target = 0;
1043f3237bfeSMark Adams TS ts;
1044f3237bfeSMark Adams Mat J;
1045f3237bfeSMark Adams LandauCtx *ctx;
10463cdbf90fSMark Adams PetscReal shift = 0;
10471b431c67SMark Adams PetscBool use_uniform_particle_grid = PETSC_TRUE;
1048f3237bfeSMark Adams
1049327415f7SBarry Smith PetscFunctionBeginUser;
1050f3237bfeSMark Adams PetscCall(PetscInitialize(&argc, &argv, NULL, help));
1051f3237bfeSMark Adams // process args
1052d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_SELF, "", "Collision Options", "DMPLEX");
1053f3237bfeSMark Adams PetscCall(PetscOptionsInt("-dim", "Velocity space dimension", "ex30.c", dim, &dim, NULL));
10541b431c67SMark Adams PetscCall(PetscOptionsInt("-number_spatial_vertices", "Number of user spatial vertices to be batched for Landau", "ex30.c", num_vertices, &num_vertices, NULL));
1055f3237bfeSMark 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));
10561b431c67SMark Adams PetscCall(PetscOptionsBool("-use_uniform_particle_grid", "Use uniform particle grid", "ex30.c", use_uniform_particle_grid, &use_uniform_particle_grid, NULL));
1057d043ef4cSMark Adams PetscCall(PetscOptionsInt("-vertex_view_target", "Global vertex for diagnostics", "ex30.c", v_target, &v_target, NULL));
1058d043ef4cSMark Adams PetscCall(PetscOptionsReal("-e_shift", "Bi-Maxwellian shift", "ex30.c", shift, &shift, NULL));
1059d043ef4cSMark 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);
10601b431c67SMark Adams PetscCall(PetscOptionsInt("-grid_view_target", "Grid to view with diagnostics", "ex30.c", g_target, &g_target, NULL));
1061d0609cedSBarry Smith PetscOptionsEnd();
1062f3237bfeSMark Adams /* Create a mesh */
1063f3237bfeSMark Adams PetscCall(DMPlexLandauCreateVelocitySpace(PETSC_COMM_SELF, dim, "", &X, &J, &pack));
1064d043ef4cSMark Adams PetscCall(DMGetApplicationContext(pack, &ctx));
1065f3237bfeSMark Adams PetscCall(DMSetUp(pack));
1066f3237bfeSMark Adams PetscCall(DMSetOutputSequenceNumber(pack, 0, 0.0));
10671b431c67SMark Adams PetscCheck(g_target < ctx->num_grids, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Grid to view %" PetscInt_FMT " should be < number of grids %" PetscInt_FMT, g_target, ctx->num_grids);
10681b431c67SMark Adams PetscCheck(ctx->batch_view_idx == v_target % ctx->batch_sz, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Global view index %" PetscInt_FMT " mode batch size %" PetscInt_FMT " != ctx->batch_view_idx %" PetscInt_FMT, v_target, ctx->batch_sz, ctx->batch_view_idx);
1069d043ef4cSMark Adams // PetscCheck(!use_uniform_particle_grid || !ctx->sphere, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Can not use -use_uniform_particle_grid and -dm_landau_sphere");
1070f3237bfeSMark Adams /* Create timestepping solver context */
1071f3237bfeSMark Adams PetscCall(TSCreate(PETSC_COMM_SELF, &ts));
1072f3237bfeSMark Adams PetscCall(TSSetDM(ts, pack));
1073f3237bfeSMark Adams PetscCall(TSSetIFunction(ts, NULL, DMPlexLandauIFunction, NULL));
1074f3237bfeSMark Adams PetscCall(TSSetIJacobian(ts, J, J, DMPlexLandauIJacobian, NULL));
1075f3237bfeSMark Adams PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
1076f3237bfeSMark Adams PetscCall(TSSetFromOptions(ts));
1077f3237bfeSMark Adams PetscCall(PetscObjectSetName((PetscObject)X, "X"));
10786b664d00Smarkadams4 // do particle advance
10791b431c67SMark Adams PetscCall(go(ts, X, num_vertices, Np, dim, v_target, g_target, shift, use_uniform_particle_grid));
1080984ed092SMark Adams PetscCall(MatZeroEntries(J)); // need to zero out so as to not reuse it in Landau's logic
1081f3237bfeSMark Adams /* clean up */
1082f3237bfeSMark Adams PetscCall(DMPlexLandauDestroyVelocitySpace(&pack));
1083f3237bfeSMark Adams PetscCall(TSDestroy(&ts));
1084f3237bfeSMark Adams PetscCall(VecDestroy(&X));
1085f3237bfeSMark Adams PetscCall(PetscFinalize());
1086f3237bfeSMark Adams return 0;
1087f3237bfeSMark Adams }
1088f3237bfeSMark Adams
1089f3237bfeSMark Adams /*TEST
1090f3237bfeSMark Adams
1091f3237bfeSMark Adams build:
10923cdbf90fSMark Adams requires: !complex
1093f3237bfeSMark Adams
1094f3237bfeSMark Adams testset:
1095984ed092SMark Adams requires: double defined(PETSC_USE_DMLANDAU_2D)
1096f3237bfeSMark Adams output_file: output/ex30_0.out
1097d043ef4cSMark Adams args: -dim 2 -petscspace_degree 3 -dm_landau_num_species_grid 1,1,1 -dm_refine 1 -number_particles_per_dimension 20 \
10981b431c67SMark Adams -dm_landau_batch_size 4 -number_spatial_vertices 6 -vertex_view_target 5 -grid_view_target 1 -dm_landau_batch_view_idx 1 \
1099f3237bfeSMark 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 \
1100d043ef4cSMark Adams -ftop_ksp_rtol 1e-10 -ftop_ksp_type lsqr -ftop_pc_type bjacobi -ftop_sub_pc_factor_shift_type nonzero -ftop_sub_pc_type lu -ftop_ksp_error_if_not_converged \
1101d043ef4cSMark Adams -ksp_type gmres -ksp_error_if_not_converged -dm_landau_verbose 4 -print_entropy \
1102d043ef4cSMark Adams -ptof_ksp_type cg -ptof_pc_type jacobi -ptof_ksp_rtol 1e-12 -ptof_ksp_error_if_not_converged\
1103d043ef4cSMark Adams -snes_converged_reason -snes_monitor -snes_rtol 1e-12 -snes_stol 1e-12 \
1104*188af4bfSBarry Smith -ts_time_step 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
1105f3237bfeSMark Adams test:
1106f3237bfeSMark Adams suffix: cpu
1107d043ef4cSMark Adams args: -dm_landau_device_type cpu -pc_type jacobi
1108f3237bfeSMark Adams test:
1109f3237bfeSMark Adams suffix: kokkos
1110bbdffb7fSJunchao Zhang # failed on Sunspot@ALCF with sycl
1111bbdffb7fSJunchao Zhang requires: kokkos_kernels !openmp !sycl
11123cdbf90fSMark 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
1113f3237bfeSMark Adams
1114984ed092SMark Adams testset:
1115984ed092SMark Adams requires: double !defined(PETSC_USE_DMLANDAU_2D)
1116984ed092SMark Adams output_file: output/ex30_3d.out
1117d043ef4cSMark Adams args: -dim 3 -petscspace_degree 2 -dm_landau_num_species_grid 1,1 -dm_refine 0 -number_particles_per_dimension 10 -dm_plex_hash_location \
11181b431c67SMark Adams -dm_landau_batch_size 1 -number_spatial_vertices 1 -vertex_view_target 0 -grid_view_target 0 -dm_landau_batch_view_idx 0 \
1119d043ef4cSMark Adams -dm_landau_n 1.000018,1 -dm_landau_thermal_temps 2,1 -dm_landau_ion_masses 2 -dm_landau_ion_charges 1 \
1120d043ef4cSMark Adams -ftop_ksp_type cg -ftop_pc_type jacobi -ftop_ksp_rtol 1e-12 -ftop_ksp_error_if_not_converged -ksp_type preonly -pc_type lu -ksp_error_if_not_converged \
1121d043ef4cSMark Adams -ptof_ksp_type cg -ptof_pc_type jacobi -ptof_ksp_rtol 1e-12 -ptof_ksp_error_if_not_converged \
11223cdbf90fSMark Adams -snes_converged_reason -snes_monitor -snes_rtol 1e-12 -snes_stol 1e-12 \
1123*188af4bfSBarry Smith -ts_time_step 0.1 -ts_exact_final_time stepover -ts_max_snes_failures -1 -ts_max_steps 1 -ts_monitor -ts_type beuler -print_entropy
1124984ed092SMark Adams test:
1125984ed092SMark Adams suffix: cpu_3d
1126984ed092SMark Adams args: -dm_landau_device_type cpu
1127984ed092SMark Adams test:
1128984ed092SMark Adams suffix: kokkos_3d
11293cdbf90fSMark Adams requires: kokkos_kernels !openmp
11303cdbf90fSMark 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
11313cdbf90fSMark Adams
1132d043ef4cSMark Adams test:
1133d043ef4cSMark Adams suffix: conserve
11343cdbf90fSMark Adams requires: !complex double defined(PETSC_USE_DMLANDAU_2D) !cuda
1135*188af4bfSBarry Smith args: -dm_landau_batch_size 4 -dm_refine 0 -dm_landau_num_species_grid 1 -dm_landau_thermal_temps 1 -petscspace_degree 3 -snes_converged_reason -ts_type beuler -ts_time_step .1 \
1136d043ef4cSMark Adams -ts_max_steps 1 -ksp_type preonly -ksp_error_if_not_converged -snes_rtol 1e-14 -snes_stol 1e-14 -dm_landau_device_type cpu -number_particles_per_dimension 20 \
1137d043ef4cSMark Adams -ptof_ksp_type cg -ptof_pc_type jacobi -ptof_ksp_rtol 1e-14 -ptof_ksp_error_if_not_converged -pc_type lu -dm_landau_simplex 1 -use_uniform_particle_grid false -dm_landau_sphere -print_entropy -number_particles_per_dimension 50 -ftop_ksp_type cg -ftop_pc_type jacobi -ftop_ksp_rtol 1e-14
1138984ed092SMark Adams
1139f3237bfeSMark Adams TEST*/
1140