xref: /petsc/src/dm/impls/swarm/tests/ex7.c (revision 2286efddd54511ab18e8e2adb1e023c4bf8f0b92)
13c611800SMark Adams static char help[] = "Example program demonstrating projection between particle and finite element spaces using OpenMP in 2D cylindrical coordinates\n";
23c611800SMark Adams 
346233b44SBarry Smith #include <petscdmplex.h>
446233b44SBarry Smith #include <petscds.h>
546233b44SBarry Smith #include <petscdmswarm.h>
646233b44SBarry Smith #include <petscksp.h>
73c611800SMark Adams #include <petsc/private/petscimpl.h>
83c611800SMark Adams #if defined(PETSC_HAVE_OPENMP) && defined(PETSC_HAVE_THREADSAFETY)
93c611800SMark Adams   #include <omp.h>
103c611800SMark Adams #endif
113c611800SMark Adams 
123c611800SMark Adams typedef struct {
133c611800SMark Adams   Mat MpTrans;
143c611800SMark Adams   Mat Mp;
153c611800SMark Adams   Vec ff;
163c611800SMark Adams   Vec uu;
173c611800SMark Adams } MatShellCtx;
183c611800SMark Adams 
MatMultMtM_SeqAIJ(Mat MtM,Vec xx,Vec yy)19d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultMtM_SeqAIJ(Mat MtM, Vec xx, Vec yy)
20d71ae5a4SJacob Faibussowitsch {
213c611800SMark Adams   MatShellCtx *matshellctx;
223c611800SMark Adams 
233c611800SMark Adams   PetscFunctionBeginUser;
249566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(MtM, &matshellctx));
2528b400f6SJacob Faibussowitsch   PetscCheck(matshellctx, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "No context");
269566063dSJacob Faibussowitsch   PetscCall(MatMult(matshellctx->Mp, xx, matshellctx->ff));
279566063dSJacob Faibussowitsch   PetscCall(MatMult(matshellctx->MpTrans, matshellctx->ff, yy));
283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
293c611800SMark Adams }
303c611800SMark Adams 
MatMultAddMtM_SeqAIJ(Mat MtM,Vec xx,Vec yy,Vec zz)31d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMultAddMtM_SeqAIJ(Mat MtM, Vec xx, Vec yy, Vec zz)
32d71ae5a4SJacob Faibussowitsch {
333c611800SMark Adams   MatShellCtx *matshellctx;
343c611800SMark Adams 
353c611800SMark Adams   PetscFunctionBeginUser;
369566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(MtM, &matshellctx));
3728b400f6SJacob Faibussowitsch   PetscCheck(matshellctx, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "No context");
389566063dSJacob Faibussowitsch   PetscCall(MatMult(matshellctx->Mp, xx, matshellctx->ff));
399566063dSJacob Faibussowitsch   PetscCall(MatMultAdd(matshellctx->MpTrans, matshellctx->ff, yy, zz));
403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
413c611800SMark Adams }
423c611800SMark Adams 
createSwarm(const DM dm,DM * sw)43d71ae5a4SJacob Faibussowitsch PetscErrorCode createSwarm(const DM dm, DM *sw)
44d71ae5a4SJacob Faibussowitsch {
453c611800SMark Adams   PetscInt Nc = 1, dim = 2;
463c611800SMark Adams 
473c611800SMark Adams   PetscFunctionBeginUser;
489566063dSJacob Faibussowitsch   PetscCall(DMCreate(PETSC_COMM_SELF, sw));
499566063dSJacob Faibussowitsch   PetscCall(DMSetType(*sw, DMSWARM));
509566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(*sw, dim));
519566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetType(*sw, DMSWARM_PIC));
529566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetCellDM(*sw, dm));
539566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "w_q", Nc, PETSC_SCALAR));
549566063dSJacob Faibussowitsch   PetscCall(DMSwarmFinalizeFieldRegister(*sw));
559566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*sw));
563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
573c611800SMark Adams }
583c611800SMark Adams 
gridToParticles(const DM dm,DM sw,PetscReal * moments,Vec rhs,Mat M_p)59d71ae5a4SJacob Faibussowitsch PetscErrorCode gridToParticles(const DM dm, DM sw, PetscReal *moments, Vec rhs, Mat M_p)
60d71ae5a4SJacob Faibussowitsch {
613c611800SMark Adams   PetscBool     is_lsqr;
623c611800SMark Adams   KSP           ksp;
633c611800SMark Adams   Mat           PM_p = NULL, MtM, D;
643c611800SMark Adams   Vec           ff;
653c611800SMark Adams   PetscInt      Np, timestep = 0, bs, N, M, nzl;
663c611800SMark Adams   PetscReal     time = 0.0;
673c611800SMark Adams   PetscDataType dtype;
683c611800SMark Adams   MatShellCtx  *matshellctx;
693c611800SMark Adams 
703c611800SMark Adams   PetscFunctionBeginUser;
719566063dSJacob Faibussowitsch   PetscCall(KSPCreate(PETSC_COMM_SELF, &ksp));
729566063dSJacob Faibussowitsch   PetscCall(KSPSetOptionsPrefix(ksp, "ftop_"));
739566063dSJacob Faibussowitsch   PetscCall(KSPSetFromOptions(ksp));
749566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)ksp, KSPLSQR, &is_lsqr));
753c611800SMark Adams   if (!is_lsqr) {
769566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(M_p, &M, &N));
773c611800SMark Adams     if (N > M) {
783c611800SMark Adams       PC pc;
79feadbf81SMark Adams       PetscCall(PetscInfo(ksp, " M (%" PetscInt_FMT ") < N (%" PetscInt_FMT ") -- skip revert to lsqr\n", M, N));
803c611800SMark Adams       is_lsqr = PETSC_TRUE;
819566063dSJacob Faibussowitsch       PetscCall(KSPSetType(ksp, KSPLSQR));
829566063dSJacob Faibussowitsch       PetscCall(KSPGetPC(ksp, &pc));
839566063dSJacob Faibussowitsch       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
843c611800SMark Adams     } else {
859566063dSJacob Faibussowitsch       PetscCall(PetscNew(&matshellctx));
869566063dSJacob Faibussowitsch       PetscCall(MatCreateShell(PetscObjectComm((PetscObject)dm), N, N, PETSC_DECIDE, PETSC_DECIDE, matshellctx, &MtM));
879566063dSJacob Faibussowitsch       PetscCall(MatTranspose(M_p, MAT_INITIAL_MATRIX, &matshellctx->MpTrans));
883c611800SMark Adams       matshellctx->Mp = M_p;
89*57d50842SBarry Smith       PetscCall(MatShellSetOperation(MtM, MATOP_MULT, (PetscErrorCodeFn *)MatMultMtM_SeqAIJ));
90*57d50842SBarry Smith       PetscCall(MatShellSetOperation(MtM, MATOP_MULT_ADD, (PetscErrorCodeFn *)MatMultAddMtM_SeqAIJ));
919566063dSJacob Faibussowitsch       PetscCall(MatCreateVecs(M_p, &matshellctx->uu, &matshellctx->ff));
929566063dSJacob Faibussowitsch       PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, N, N, 1, NULL, &D));
933c611800SMark Adams       for (int i = 0; i < N; i++) {
943c611800SMark Adams         const PetscScalar *vals;
953c611800SMark Adams         const PetscInt    *cols;
963c611800SMark Adams         PetscScalar        dot = 0;
979566063dSJacob Faibussowitsch         PetscCall(MatGetRow(matshellctx->MpTrans, i, &nzl, &cols, &vals));
983c611800SMark Adams         for (int ii = 0; ii < nzl; ii++) dot += PetscSqr(vals[ii]);
9963a3b9bcSJacob Faibussowitsch         PetscCheck(dot != 0.0, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Row %d is empty", i);
1009566063dSJacob Faibussowitsch         PetscCall(MatSetValue(D, i, i, dot, INSERT_VALUES));
1013c611800SMark Adams       }
1029566063dSJacob Faibussowitsch       PetscCall(MatAssemblyBegin(D, MAT_FINAL_ASSEMBLY));
1039566063dSJacob Faibussowitsch       PetscCall(MatAssemblyEnd(D, MAT_FINAL_ASSEMBLY));
1043ba16761SJacob Faibussowitsch       PetscCall(PetscInfo(M_p, "createMtMKSP Have %" PetscInt_FMT " eqs, nzl = %" PetscInt_FMT "\n", N, nzl));
1059566063dSJacob Faibussowitsch       PetscCall(KSPSetOperators(ksp, MtM, D));
1069566063dSJacob Faibussowitsch       PetscCall(MatViewFromOptions(D, NULL, "-ftop2_D_mat_view"));
1079566063dSJacob Faibussowitsch       PetscCall(MatViewFromOptions(M_p, NULL, "-ftop2_Mp_mat_view"));
1089566063dSJacob Faibussowitsch       PetscCall(MatViewFromOptions(matshellctx->MpTrans, NULL, "-ftop2_MpTranspose_mat_view"));
1093c611800SMark Adams     }
1103c611800SMark Adams   }
1113c611800SMark Adams   if (is_lsqr) {
1123c611800SMark Adams     PC        pc;
1133c611800SMark Adams     PetscBool is_bjac;
1149566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(ksp, &pc));
1159566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCBJACOBI, &is_bjac));
1163c611800SMark Adams     if (is_bjac) {
1179566063dSJacob Faibussowitsch       PetscCall(DMSwarmCreateMassMatrixSquare(sw, dm, &PM_p));
1189566063dSJacob Faibussowitsch       PetscCall(KSPSetOperators(ksp, M_p, PM_p));
1193c611800SMark Adams     } else {
1209566063dSJacob Faibussowitsch       PetscCall(KSPSetOperators(ksp, M_p, M_p));
1213c611800SMark Adams     }
1223c611800SMark Adams   }
1239566063dSJacob Faibussowitsch   PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &ff)); // this grabs access !!!!!
1243c611800SMark Adams   if (!is_lsqr) {
1259566063dSJacob Faibussowitsch     PetscCall(KSPSolve(ksp, rhs, matshellctx->uu));
1269566063dSJacob Faibussowitsch     PetscCall(MatMult(M_p, matshellctx->uu, ff));
1279566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&matshellctx->MpTrans));
1289566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&matshellctx->ff));
1299566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&matshellctx->uu));
1309566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&D));
1319566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&MtM));
1329566063dSJacob Faibussowitsch     PetscCall(PetscFree(matshellctx));
1333c611800SMark Adams   } else {
1349566063dSJacob Faibussowitsch     PetscCall(KSPSolveTranspose(ksp, rhs, ff));
1353c611800SMark Adams   }
1369566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&ksp));
1373c611800SMark Adams   /* Visualize particle field */
1389566063dSJacob Faibussowitsch   PetscCall(DMSetOutputSequenceNumber(sw, timestep, time));
1399566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(ff, NULL, "-weights_view"));
1409566063dSJacob Faibussowitsch   PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &ff));
1413c611800SMark Adams 
1423c611800SMark Adams   /* compute energy */
1433c611800SMark Adams   if (moments) {
1443c611800SMark Adams     PetscReal *wq, *coords;
1459566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetLocalSize(sw, &Np));
1469566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void **)&wq));
1479566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
1483c611800SMark Adams     moments[0] = moments[1] = moments[2] = 0;
1493c611800SMark Adams     for (int p = 0; p < Np; p++) {
1503c611800SMark Adams       moments[0] += wq[p];
1513c611800SMark Adams       moments[1] += wq[p] * coords[p * 2 + 0]; // x-momentum
1523c611800SMark Adams       moments[2] += wq[p] * (PetscSqr(coords[p * 2 + 0]) + PetscSqr(coords[p * 2 + 1]));
1533c611800SMark Adams     }
1549566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
1559566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void **)&wq));
1563c611800SMark Adams   }
1579566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&PM_p));
1583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1593c611800SMark Adams }
1603c611800SMark Adams 
particlesToGrid(const DM dm,DM sw,const PetscInt Np,const PetscInt a_tid,const PetscInt dim,const PetscInt target,const PetscReal xx[],const PetscReal yy[],const PetscReal a_wp[],Vec rho,Mat * Mp_out)161d71ae5a4SJacob Faibussowitsch PetscErrorCode particlesToGrid(const DM dm, DM sw, const PetscInt Np, const PetscInt a_tid, const PetscInt dim, const PetscInt target, const PetscReal xx[], const PetscReal yy[], const PetscReal a_wp[], Vec rho, Mat *Mp_out)
162d71ae5a4SJacob Faibussowitsch {
1633c611800SMark Adams   PetscBool     removePoints = PETSC_TRUE;
1643c611800SMark Adams   PetscReal    *wq, *coords;
1653c611800SMark Adams   PetscDataType dtype;
1663c611800SMark Adams   Mat           M_p;
1673c611800SMark Adams   Vec           ff;
1683c611800SMark Adams   PetscInt      bs, p, zero = 0;
1693c611800SMark Adams 
1703c611800SMark Adams   PetscFunctionBeginUser;
1719566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetLocalSizes(sw, Np, zero));
1729566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void **)&wq));
1739566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
1743c611800SMark Adams   for (p = 0; p < Np; p++) {
1753c611800SMark Adams     coords[p * 2 + 0] = xx[p];
1763c611800SMark Adams     coords[p * 2 + 1] = yy[p];
1773c611800SMark Adams     wq[p]             = a_wp[p];
1783c611800SMark Adams   }
1799566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void **)&coords));
1809566063dSJacob Faibussowitsch   PetscCall(DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void **)&wq));
1819566063dSJacob Faibussowitsch   PetscCall(DMSwarmMigrate(sw, removePoints));
1829566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)sw, "Particle Grid"));
1839566063dSJacob Faibussowitsch   if (a_tid == target) PetscCall(DMViewFromOptions(sw, NULL, "-swarm_view"));
1843c611800SMark Adams 
1853c611800SMark Adams   /* Project particles to field */
1863c611800SMark Adams   /* This gives M f = \int_\Omega \phi f, which looks like a rhs for a PDE */
1879566063dSJacob Faibussowitsch   PetscCall(DMCreateMassMatrix(sw, dm, &M_p));
1889566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)rho, "rho"));
1893c611800SMark Adams 
1909566063dSJacob Faibussowitsch   PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &ff)); // this grabs access !!!!!
1919566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)ff, "weights"));
1929566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(M_p, ff, rho));
1939566063dSJacob Faibussowitsch   PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &ff));
1943c611800SMark Adams 
1953c611800SMark Adams   /* Visualize mesh field */
1969566063dSJacob Faibussowitsch   if (a_tid == target) PetscCall(VecViewFromOptions(rho, NULL, "-rho_view"));
1973c611800SMark Adams   // output
1983c611800SMark Adams   *Mp_out = M_p;
1993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2003c611800SMark Adams }
2014d86920dSPierre Jolivet 
maxwellian(PetscInt dim,const PetscReal x[],PetscReal kt_m,PetscReal n,PetscScalar * u)202d71ae5a4SJacob Faibussowitsch static PetscErrorCode maxwellian(PetscInt dim, const PetscReal x[], PetscReal kt_m, PetscReal n, PetscScalar *u)
203d71ae5a4SJacob Faibussowitsch {
2043c611800SMark Adams   PetscInt  i;
2053c611800SMark Adams   PetscReal v2 = 0, theta = 2 * kt_m; /* theta = 2kT/mc^2 */
2063c611800SMark Adams 
2073c611800SMark Adams   PetscFunctionBegin;
2083c611800SMark Adams   /* compute the exponents, v^2 */
2093c611800SMark Adams   for (i = 0; i < dim; ++i) v2 += x[i] * x[i];
2103c611800SMark Adams   /* evaluate the Maxwellian */
2113c611800SMark Adams   u[0] = n * PetscPowReal(PETSC_PI * theta, -1.5) * (PetscExpReal(-v2 / theta)) * 2. * PETSC_PI * x[1]; // radial term for 2D axi-sym.
2123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2133c611800SMark Adams }
21428114503SMark Adams 
2153c611800SMark Adams #define MAX_NUM_THRDS 12
go()216d71ae5a4SJacob Faibussowitsch PetscErrorCode go()
217d71ae5a4SJacob Faibussowitsch {
2183c611800SMark Adams   DM            dm_t[MAX_NUM_THRDS], sw_t[MAX_NUM_THRDS];
2193c611800SMark Adams   PetscFE       fe;
22028114503SMark Adams   PetscInt      dim = 2, Nc = 1, i, faces[3];
2213c611800SMark Adams   PetscInt      Np[2] = {10, 10}, Np2[2], field = 0, target = 0, Np_t[MAX_NUM_THRDS];
22228114503SMark Adams   PetscReal     moments_0[3], moments_1[3], vol = 1;
22328114503SMark Adams   PetscReal     lo[3] = {-5, 0, -5}, hi[3] = {5, 5, 5}, h[3], hp[3], *xx_t[MAX_NUM_THRDS], *yy_t[MAX_NUM_THRDS], *wp_t[MAX_NUM_THRDS];
2243c611800SMark Adams   Vec           rho_t[MAX_NUM_THRDS], rhs_t[MAX_NUM_THRDS];
2253c611800SMark Adams   Mat           M_p_t[MAX_NUM_THRDS];
2263c611800SMark Adams   PetscLogStage stage;
2273c611800SMark Adams   PetscLogEvent swarm_create_ev, solve_ev, solve_loop_ev;
2283c611800SMark Adams #if defined(PETSC_HAVE_OPENMP) && defined(PETSC_HAVE_THREADSAFETY)
2293c611800SMark Adams   PetscInt numthreads = PetscNumOMPThreads;
2303c611800SMark Adams #else
2313c611800SMark Adams   PetscInt numthreads = 1;
2323c611800SMark Adams #endif
2333c611800SMark Adams 
2343c611800SMark Adams   PetscFunctionBeginUser;
2353c611800SMark Adams #if defined(PETSC_HAVE_OPENMP) && defined(PETSC_HAVE_THREADSAFETY)
23663a3b9bcSJacob Faibussowitsch   PetscCheck(numthreads <= MAX_NUM_THRDS, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Too many threads %" PetscInt_FMT " > %d", numthreads, MAX_NUM_THRDS);
23763a3b9bcSJacob Faibussowitsch   PetscCheck(numthreads > 0, PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "No threads %" PetscInt_FMT " > %d", numthreads, MAX_NUM_THRDS);
2383c611800SMark Adams #endif
2393c611800SMark Adams   if (target >= numthreads) target = numthreads - 1;
2409566063dSJacob Faibussowitsch   PetscCall(PetscLogEventRegister("Create Swarm", DM_CLASSID, &swarm_create_ev));
2419566063dSJacob Faibussowitsch   PetscCall(PetscLogEventRegister("Single solve", DM_CLASSID, &solve_ev));
2429566063dSJacob Faibussowitsch   PetscCall(PetscLogEventRegister("Solve loop", DM_CLASSID, &solve_loop_ev));
2439566063dSJacob Faibussowitsch   PetscCall(PetscLogStageRegister("Solve", &stage));
2443c611800SMark Adams   i = dim;
2459566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &i, NULL));
2463c611800SMark Adams   i = dim;
2479566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-np", Np, &i, NULL));
2483c611800SMark Adams   /* Create thread meshes */
2493c611800SMark Adams   for (int tid = 0; tid < numthreads; tid++) {
250917d862eSmarkadams4     PetscBool simplex = PETSC_FALSE;
251917d862eSmarkadams4     PetscCall(PetscOptionsGetBool(NULL, NULL, "-dm_plex_simplex", &simplex, NULL));
2523c611800SMark Adams     // setup mesh dm_t, could use PETSc's Landau create velocity space mesh here to get dm_t[tid]
2539566063dSJacob Faibussowitsch     PetscCall(DMCreate(PETSC_COMM_SELF, &dm_t[tid]));
2549566063dSJacob Faibussowitsch     PetscCall(DMSetType(dm_t[tid], DMPLEX));
2559566063dSJacob Faibussowitsch     PetscCall(DMSetFromOptions(dm_t[tid]));
256917d862eSmarkadams4     PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, Nc, simplex, "", PETSC_DECIDE, &fe));
2579566063dSJacob Faibussowitsch     PetscCall(PetscFESetFromOptions(fe));
2589566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)fe, "fe"));
2599566063dSJacob Faibussowitsch     PetscCall(DMSetField(dm_t[tid], field, NULL, (PetscObject)fe));
2609566063dSJacob Faibussowitsch     PetscCall(DMCreateDS(dm_t[tid]));
2619566063dSJacob Faibussowitsch     PetscCall(PetscFEDestroy(&fe));
2623c611800SMark Adams     // helper vectors
2639566063dSJacob Faibussowitsch     PetscCall(DMCreateGlobalVector(dm_t[tid], &rho_t[tid]));
2649566063dSJacob Faibussowitsch     PetscCall(DMCreateGlobalVector(dm_t[tid], &rhs_t[tid]));
2653c611800SMark Adams     // this mimics application code
2669566063dSJacob Faibussowitsch     PetscCall(DMGetBoundingBox(dm_t[tid], lo, hi));
2673c611800SMark Adams     if (tid == target) {
2689566063dSJacob Faibussowitsch       PetscCall(DMViewFromOptions(dm_t[tid], NULL, "-dm_view"));
2693c611800SMark Adams       for (i = 0, vol = 1; i < dim; i++) {
2703c611800SMark Adams         h[i]  = (hi[i] - lo[i]) / faces[i];
2713c611800SMark Adams         hp[i] = (hi[i] - lo[i]) / Np[i];
2723c611800SMark Adams         vol *= (hi[i] - lo[i]);
27363a3b9bcSJacob Faibussowitsch         PetscCall(PetscInfo(dm_t[tid], " lo = %g hi = %g n = %" PetscInt_FMT " h = %g hp = %g\n", (double)lo[i], (double)hi[i], faces[i], (double)h[i], (double)hp[i]));
2744279555eSSatish Balay         (void)h[i];
2753c611800SMark Adams       }
2763c611800SMark Adams     }
2773c611800SMark Adams   }
2783c611800SMark Adams   // prepare particle data for problems. This mimics application code
2799566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(swarm_create_ev, 0, 0, 0, 0));
2809371c9d4SSatish Balay   Np2[0] = Np[0];
2819371c9d4SSatish Balay   Np2[1] = Np[1];
2823c611800SMark Adams   for (int tid = 0; tid < numthreads; tid++) { // change size of particle list a little
2833c611800SMark Adams     Np_t[tid] = Np2[0] * Np2[1];
2849566063dSJacob Faibussowitsch     PetscCall(PetscMalloc3(Np_t[tid], &xx_t[tid], Np_t[tid], &yy_t[tid], Np_t[tid], &wp_t[tid]));
285ad540459SPierre Jolivet     if (tid == target) moments_0[0] = moments_0[1] = moments_0[2] = 0;
2863c611800SMark Adams     for (int pi = 0, pp = 0; pi < Np2[0]; pi++) {
2873c611800SMark Adams       for (int pj = 0; pj < Np2[1]; pj++, pp++) {
2883c611800SMark Adams         xx_t[tid][pp] = lo[0] + hp[0] / 2. + pi * hp[0];
2893c611800SMark Adams         yy_t[tid][pp] = lo[1] + hp[1] / 2. + pj * hp[1];
2903c611800SMark Adams         {
2913c611800SMark Adams           PetscReal x[] = {xx_t[tid][pp], yy_t[tid][pp]};
2929566063dSJacob Faibussowitsch           PetscCall(maxwellian(2, x, 1.0, vol / (PetscReal)Np_t[tid], &wp_t[tid][pp]));
2933c611800SMark Adams         }
2943c611800SMark Adams         if (tid == target) { //energy_0 += wp_t[tid][pp]*(PetscSqr(xx_t[tid][pp])+PetscSqr(yy_t[tid][pp]));
2953c611800SMark Adams           moments_0[0] += wp_t[tid][pp];
2963c611800SMark Adams           moments_0[1] += wp_t[tid][pp] * xx_t[tid][pp]; // x-momentum
2973c611800SMark Adams           moments_0[2] += wp_t[tid][pp] * (PetscSqr(xx_t[tid][pp]) + PetscSqr(yy_t[tid][pp]));
2983c611800SMark Adams         }
2993c611800SMark Adams       }
3003c611800SMark Adams     }
3019371c9d4SSatish Balay     Np2[0]++;
3029371c9d4SSatish Balay     Np2[1]++;
3033c611800SMark Adams   }
3049566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(swarm_create_ev, 0, 0, 0, 0));
3059566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(solve_ev, 0, 0, 0, 0));
3063c611800SMark Adams   /* Create particle swarm */
3073c611800SMark Adams   PetscPragmaOMP(parallel for)
308f37bacd1SJacob Faibussowitsch   for (int tid = 0; tid < numthreads; tid++) PetscCallAbort(PETSC_COMM_SELF, createSwarm(dm_t[tid], &sw_t[tid]));
3093c611800SMark Adams   PetscPragmaOMP(parallel for)
310f37bacd1SJacob Faibussowitsch   for (int tid = 0; tid < numthreads; tid++) PetscCallAbort(PETSC_COMM_SELF, particlesToGrid(dm_t[tid], sw_t[tid], Np_t[tid], tid, dim, target, xx_t[tid], yy_t[tid], wp_t[tid], rho_t[tid], &M_p_t[tid]));
3113c611800SMark Adams   /* Project field to particles */
3123c611800SMark Adams   /*   This gives f_p = M_p^+ M f */
3133c611800SMark Adams   PetscPragmaOMP(parallel for)
314f37bacd1SJacob Faibussowitsch   for (int tid = 0; tid < numthreads; tid++) PetscCallAbort(PETSC_COMM_SELF, VecCopy(rho_t[tid], rhs_t[tid])); /* Identity: M^1 M rho */
3153c611800SMark Adams   PetscPragmaOMP(parallel for)
316f37bacd1SJacob Faibussowitsch   for (int tid = 0; tid < numthreads; tid++) PetscCallAbort(PETSC_COMM_SELF, gridToParticles(dm_t[tid], sw_t[tid], (tid == target) ? moments_1 : NULL, rhs_t[tid], M_p_t[tid]));
3173c611800SMark Adams   /* Cleanup */
3183c611800SMark Adams   for (int tid = 0; tid < numthreads; tid++) {
3199566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&M_p_t[tid]));
3209566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&sw_t[tid]));
3213c611800SMark Adams   }
3229566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(solve_ev, 0, 0, 0, 0));
3233c611800SMark Adams   //
32463a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF, "Total number density: %20.12e (%20.12e); x-momentum = %g (%g); energy = %g error = %e, %" PetscInt_FMT " particles. Use %" PetscInt_FMT " threads\n", (double)moments_1[0], (double)moments_0[0], (double)moments_1[1], (double)moments_0[1], (double)moments_1[2], (double)((moments_1[2] - moments_0[2]) / moments_0[2]), Np[0] * Np[1], numthreads));
3253c611800SMark Adams   /* Cleanup */
3263c611800SMark Adams   for (int tid = 0; tid < numthreads; tid++) {
3279566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&rho_t[tid]));
3289566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&rhs_t[tid]));
3299566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&dm_t[tid]));
3309566063dSJacob Faibussowitsch     PetscCall(PetscFree3(xx_t[tid], yy_t[tid], wp_t[tid]));
3313c611800SMark Adams   }
3323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3333c611800SMark Adams }
3343c611800SMark Adams 
main(int argc,char ** argv)335d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
336d71ae5a4SJacob Faibussowitsch {
337327415f7SBarry Smith   PetscFunctionBeginUser;
3389884f84cSJunchao Zhang #if defined(PETSC_HAVE_OPENMP) && defined(PETSC_HAVE_THREADSAFETY)
3399884f84cSJunchao Zhang   PETSC_MPI_THREAD_REQUIRED = MPI_THREAD_MULTIPLE; // use thread multiple if multiple threads call petsc
3409884f84cSJunchao Zhang #endif
3419566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
3429566063dSJacob Faibussowitsch   PetscCall(go());
3439566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
344b122ec5aSJacob Faibussowitsch   return 0;
3453c611800SMark Adams }
3463c611800SMark Adams 
3473c611800SMark Adams /*TEST
3483c611800SMark Adams 
3493c611800SMark Adams   build:
3503c611800SMark Adams     requires: !complex
3513c611800SMark Adams 
3523c611800SMark Adams   test:
3533c611800SMark Adams     suffix: 0
3543c611800SMark Adams     requires: double triangle
355feadbf81SMark Adams     args: -dm_plex_simplex 0 -dm_plex_box_faces 8,4 -np 10,10 -dm_plex_box_lower -2.0,0.0 -dm_plex_box_upper 2.0,2.0 -petscspace_degree 2 -ftop_ksp_type lsqr -ftop_pc_type none -dm_view -ftop_ksp_converged_reason -ftop_ksp_rtol 1.e-14
3563c611800SMark Adams     filter: grep -v DM_ | grep -v atomic
357bb928f97SBarry Smith     TODO: broken due to thread safety problems
3583c611800SMark Adams 
3593c611800SMark Adams   test:
3603c611800SMark Adams     suffix: 1
3613c611800SMark Adams     requires: double triangle
362feadbf81SMark Adams     args: -dm_plex_simplex 0 -dm_plex_box_faces 8,4 -np 10,10 -dm_plex_box_lower -2.0,0.0 -dm_plex_box_upper 2.0,2.0 -petscspace_degree 2 -dm_plex_hash_location -ftop_ksp_type lsqr -ftop_pc_type bjacobi -ftop_sub_pc_type lu -ftop_sub_pc_factor_shift_type nonzero -dm_view -ftop_ksp_converged_reason -ftop_ksp_rtol 1.e-14
3633c611800SMark Adams     filter: grep -v DM_ | grep -v atomic
364bb928f97SBarry Smith     TODO: broken due to thread safety problems
3653c611800SMark Adams 
3663c611800SMark Adams   test:
3673c611800SMark Adams     suffix: 2
3683c611800SMark Adams     requires: double triangle
369917d862eSmarkadams4     args: -dm_plex_simplex 1 -dm_plex_box_faces 8,4 -np 15,15 -dm_plex_box_lower -2.0,0.0 -dm_plex_box_upper 2.0,2.0 -petscspace_degree 2 -ftop_ksp_type cg -ftop_pc_type jacobi -dm_view -ftop_ksp_converged_reason -ftop_ksp_rtol 1.e-14
3703c611800SMark Adams     filter: grep -v DM_ | grep -v atomic
371bb928f97SBarry Smith     TODO: broken due to thread safety problems
3723c611800SMark Adams 
3733c611800SMark Adams TEST*/
374