xref: /petsc/src/dm/impls/swarm/tests/ex7.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
13c611800SMark Adams static char help[] = "Example program demonstrating projection between particle and finite element spaces using OpenMP in 2D cylindrical coordinates\n";
23c611800SMark Adams 
33c611800SMark Adams #include "petscdmplex.h"
43c611800SMark Adams #include "petscds.h"
53c611800SMark Adams #include "petscdmswarm.h"
63c611800SMark Adams #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 
193c611800SMark Adams PetscErrorCode MatMultMtM_SeqAIJ(Mat MtM,Vec xx,Vec yy)
203c611800SMark Adams {
213c611800SMark Adams   MatShellCtx    *matshellctx;
223c611800SMark Adams 
233c611800SMark Adams   PetscFunctionBeginUser;
24*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShellGetContext(MtM,&matshellctx));
252c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!matshellctx,PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "No context");
26*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(matshellctx->Mp, xx, matshellctx->ff));
27*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(matshellctx->MpTrans, matshellctx->ff, yy));
283c611800SMark Adams   PetscFunctionReturn(0);
293c611800SMark Adams }
303c611800SMark Adams 
313c611800SMark Adams PetscErrorCode MatMultAddMtM_SeqAIJ(Mat MtM,Vec xx, Vec yy, Vec zz)
323c611800SMark Adams {
333c611800SMark Adams   MatShellCtx    *matshellctx;
343c611800SMark Adams 
353c611800SMark Adams   PetscFunctionBeginUser;
36*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShellGetContext(MtM,&matshellctx));
372c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!matshellctx,PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "No context");
38*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(matshellctx->Mp, xx, matshellctx->ff));
39*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultAdd(matshellctx->MpTrans, matshellctx->ff, yy, zz));
403c611800SMark Adams   PetscFunctionReturn(0);
413c611800SMark Adams }
423c611800SMark Adams 
433c611800SMark Adams PetscErrorCode createSwarm(const DM dm, DM *sw)
443c611800SMark Adams {
453c611800SMark Adams   PetscInt Nc = 1, dim = 2;
463c611800SMark Adams 
473c611800SMark Adams   PetscFunctionBeginUser;
48*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreate(PETSC_COMM_SELF, sw));
49*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetType(*sw, DMSWARM));
50*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetDimension(*sw, dim));
51*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmSetType(*sw, DMSWARM_PIC));
52*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmSetCellDM(*sw, dm));
53*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmRegisterPetscDatatypeField(*sw, "w_q", Nc, PETSC_SCALAR));
54*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmFinalizeFieldRegister(*sw));
55*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(*sw));
563c611800SMark Adams   PetscFunctionReturn(0);
573c611800SMark Adams }
583c611800SMark Adams 
593c611800SMark Adams PetscErrorCode gridToParticles(const DM dm, DM sw, PetscReal *moments, Vec rhs, Mat M_p)
603c611800SMark Adams {
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;
71*5f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPCreate(PETSC_COMM_SELF, &ksp));
72*5f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPSetOptionsPrefix(ksp, "ftop_"));
73*5f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPSetFromOptions(ksp));
74*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectTypeCompare((PetscObject)ksp,KSPLSQR,&is_lsqr));
753c611800SMark Adams   if (!is_lsqr) {
76*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatGetLocalSize(M_p, &M, &N));
773c611800SMark Adams     if (N>M) {
783c611800SMark Adams       PC        pc;
79*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscInfo(ksp, " M (%D) < M (%D) -- skip revert to lsqr\n",M,N));
803c611800SMark Adams       is_lsqr = PETSC_TRUE;
81*5f80ce2aSJacob Faibussowitsch       CHKERRQ(KSPSetType(ksp,KSPLSQR));
82*5f80ce2aSJacob Faibussowitsch       CHKERRQ(KSPGetPC(ksp,&pc));
83*5f80ce2aSJacob Faibussowitsch       CHKERRQ(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 {
85*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscNew(&matshellctx));
86*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatCreateShell(PetscObjectComm((PetscObject)dm),N,N,PETSC_DECIDE,PETSC_DECIDE,matshellctx,&MtM));
87*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatTranspose(M_p,MAT_INITIAL_MATRIX,&matshellctx->MpTrans));
883c611800SMark Adams       matshellctx->Mp = M_p;
89*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatShellSetOperation(MtM, MATOP_MULT, (void (*)(void))MatMultMtM_SeqAIJ));
90*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatShellSetOperation(MtM, MATOP_MULT_ADD, (void (*)(void))MatMultAddMtM_SeqAIJ));
91*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatCreateVecs(M_p,&matshellctx->uu,&matshellctx->ff));
92*5f80ce2aSJacob Faibussowitsch       CHKERRQ(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;
97*5f80ce2aSJacob Faibussowitsch         CHKERRQ(MatGetRow(matshellctx->MpTrans,i,&nzl,&cols,&vals));
983c611800SMark Adams         for (int ii=0 ; ii<nzl ; ii++) dot += PetscSqr(vals[ii]);
992c71b3e2SJacob Faibussowitsch         PetscCheckFalse(dot==0.0,PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Row %D is empty", i);
100*5f80ce2aSJacob Faibussowitsch         CHKERRQ(MatSetValue(D,i,i,dot,INSERT_VALUES));
1013c611800SMark Adams       }
102*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatAssemblyBegin(D, MAT_FINAL_ASSEMBLY));
103*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatAssemblyEnd(D, MAT_FINAL_ASSEMBLY));
1047d3de750SJacob Faibussowitsch       PetscInfo(M_p,"createMtMKSP Have %D eqs, nzl = %D\n",N,nzl);
105*5f80ce2aSJacob Faibussowitsch       CHKERRQ(KSPSetOperators(ksp, MtM, D));
106*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatViewFromOptions(D,NULL,"-ftop2_D_mat_view"));
107*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatViewFromOptions(M_p,NULL,"-ftop2_Mp_mat_view"));
108*5f80ce2aSJacob Faibussowitsch       CHKERRQ(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;
114*5f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPGetPC(ksp,&pc));
115*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectTypeCompare((PetscObject)pc,PCBJACOBI,&is_bjac));
1163c611800SMark Adams     if (is_bjac) {
117*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMSwarmCreateMassMatrixSquare(sw, dm, &PM_p));
118*5f80ce2aSJacob Faibussowitsch       CHKERRQ(KSPSetOperators(ksp, M_p, PM_p));
1193c611800SMark Adams     } else {
120*5f80ce2aSJacob Faibussowitsch       CHKERRQ(KSPSetOperators(ksp, M_p, M_p));
1213c611800SMark Adams     }
1223c611800SMark Adams   }
123*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &ff)); // this grabs access !!!!!
1243c611800SMark Adams   if (!is_lsqr) {
125*5f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPSolve(ksp, rhs, matshellctx->uu));
126*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMult(M_p, matshellctx->uu, ff));
127*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&matshellctx->MpTrans));
128*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&matshellctx->ff));
129*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&matshellctx->uu));
130*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&D));
131*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&MtM));
132*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(matshellctx));
1333c611800SMark Adams   } else {
134*5f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPSolveTranspose(ksp, rhs, ff));
1353c611800SMark Adams   }
136*5f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPDestroy(&ksp));
1373c611800SMark Adams   /* Visualize particle field */
138*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetOutputSequenceNumber(sw, timestep, time));
139*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecViewFromOptions(ff, NULL, "-weights_view"));
140*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &ff));
1413c611800SMark Adams 
1423c611800SMark Adams   /* compute energy */
1433c611800SMark Adams   if (moments) {
1443c611800SMark Adams     PetscReal *wq, *coords;
145*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmGetLocalSize(sw,&Np));
146*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void**)&wq));
147*5f80ce2aSJacob Faibussowitsch     CHKERRQ(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     }
154*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void**)&coords));
155*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void**)&wq));
1563c611800SMark Adams   }
157*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&PM_p));
1583c611800SMark Adams   PetscFunctionReturn(0);
1593c611800SMark Adams }
1603c611800SMark Adams 
1613c611800SMark Adams PetscErrorCode particlesToGrid(const DM dm, DM sw, const PetscInt Np, const PetscInt a_tid, const PetscInt dim, const PetscInt target,
1623c611800SMark Adams                                const PetscReal xx[], const PetscReal yy[], const PetscReal a_wp[], Vec rho, Mat *Mp_out)
1633c611800SMark Adams {
1643c611800SMark Adams 
1653c611800SMark Adams   PetscBool      removePoints = PETSC_TRUE;
1663c611800SMark Adams   PetscReal      *wq, *coords;
1673c611800SMark Adams   PetscDataType  dtype;
1683c611800SMark Adams   Mat            M_p;
1693c611800SMark Adams   Vec            ff;
1703c611800SMark Adams   PetscInt       bs,p,zero=0;
1713c611800SMark Adams 
1723c611800SMark Adams   PetscFunctionBeginUser;
173*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmSetLocalSizes(sw, Np, zero));
174*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void**)&wq));
175*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void**)&coords));
1763c611800SMark Adams   for (p=0;p<Np;p++) {
1773c611800SMark Adams     coords[p*2+0]  = xx[p];
1783c611800SMark Adams     coords[p*2+1]  = yy[p];
1793c611800SMark Adams     wq[p]          = a_wp[p];
1803c611800SMark Adams   }
181*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void**)&coords));
182*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void**)&wq));
183*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmMigrate(sw, removePoints));
184*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject)sw, "Particle Grid"));
185*5f80ce2aSJacob Faibussowitsch   if (a_tid==target) CHKERRQ(DMViewFromOptions(sw, NULL, "-swarm_view"));
1863c611800SMark Adams 
1873c611800SMark Adams   /* Project particles to field */
1883c611800SMark Adams   /* This gives M f = \int_\Omega \phi f, which looks like a rhs for a PDE */
189*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateMassMatrix(sw, dm, &M_p));
190*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject)rho, "rho"));
1913c611800SMark Adams 
192*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &ff)); // this grabs access !!!!!
193*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject)ff, "weights"));
194*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultTranspose(M_p, ff, rho));
195*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &ff));
1963c611800SMark Adams 
1973c611800SMark Adams   /* Visualize mesh field */
198*5f80ce2aSJacob Faibussowitsch   if (a_tid==target) CHKERRQ(VecViewFromOptions(rho, NULL, "-rho_view"));
1993c611800SMark Adams   // output
2003c611800SMark Adams   *Mp_out = M_p;
2013c611800SMark Adams 
2023c611800SMark Adams   PetscFunctionReturn(0);
2033c611800SMark Adams }
2043c611800SMark Adams static PetscErrorCode maxwellian(PetscInt dim, const PetscReal x[], PetscReal kt_m, PetscReal n, PetscScalar *u)
2053c611800SMark Adams {
2063c611800SMark Adams   PetscInt      i;
2073c611800SMark Adams   PetscReal     v2 = 0, theta = 2*kt_m; /* theta = 2kT/mc^2 */
2083c611800SMark Adams 
2093c611800SMark Adams   PetscFunctionBegin;
2103c611800SMark Adams   /* compute the exponents, v^2 */
2113c611800SMark Adams   for (i = 0; i < dim; ++i) v2 += x[i]*x[i];
2123c611800SMark Adams   /* evaluate the Maxwellian */
2133c611800SMark Adams   u[0] = n*PetscPowReal(PETSC_PI*theta,-1.5)*(PetscExpReal(-v2/theta)) * 2.*PETSC_PI*x[1]; // radial term for 2D axi-sym.
2143c611800SMark Adams   PetscFunctionReturn(0);
2153c611800SMark Adams }
2163c611800SMark Adams #define NUM_SOLVE_LOOPS 100
2173c611800SMark Adams #define MAX_NUM_THRDS 12
2183c611800SMark Adams PetscErrorCode go()
2193c611800SMark Adams {
2203c611800SMark Adams   DM              dm_t[MAX_NUM_THRDS], sw_t[MAX_NUM_THRDS];
2213c611800SMark Adams   PetscFE         fe;
2223c611800SMark Adams   PetscInt        dim = 2, Nc = 1, timestep = 0, i, faces[3];
2233c611800SMark Adams   PetscInt        Np[2] = {10,10}, Np2[2], field = 0, target = 0, Np_t[MAX_NUM_THRDS];
2243c611800SMark Adams   PetscReal       time = 0.0, moments_0[3], moments_1[3], vol;
2253c611800SMark 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], solve_time = 0;
2263c611800SMark Adams   Vec             rho_t[MAX_NUM_THRDS], rhs_t[MAX_NUM_THRDS];
2273c611800SMark Adams   Mat             M_p_t[MAX_NUM_THRDS];
2283c611800SMark Adams #if defined PETSC_USE_LOG
2293c611800SMark Adams   PetscLogStage   stage;
2303c611800SMark Adams   PetscLogEvent   swarm_create_ev, solve_ev, solve_loop_ev;
2313c611800SMark Adams #endif
2323c611800SMark Adams #if defined(PETSC_HAVE_OPENMP) && defined(PETSC_HAVE_THREADSAFETY)
2333c611800SMark Adams   PetscInt        numthreads = PetscNumOMPThreads;
2343c611800SMark Adams   double          starttime, endtime;
2353c611800SMark Adams #else
2363c611800SMark Adams   PetscInt        numthreads = 1;
2373c611800SMark Adams #endif
2383c611800SMark Adams 
2393c611800SMark Adams   PetscFunctionBeginUser;
2403c611800SMark Adams #if defined(PETSC_HAVE_OPENMP) && defined(PETSC_HAVE_THREADSAFETY)
2412c71b3e2SJacob Faibussowitsch   PetscCheckFalse(numthreads>MAX_NUM_THRDS,PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Too many threads %D > %D", numthreads, MAX_NUM_THRDS);
2422c71b3e2SJacob Faibussowitsch   PetscCheckFalse(numthreads<=0,PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "No threads %D > %D ", numthreads,  MAX_NUM_THRDS);
2433c611800SMark Adams #endif
2443c611800SMark Adams   if (target >= numthreads) target = numthreads-1;
245*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventRegister("Create Swarm", DM_CLASSID, &swarm_create_ev));
246*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventRegister("Single solve", DM_CLASSID, &solve_ev));
247*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventRegister("Solve loop", DM_CLASSID, &solve_loop_ev));
248*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogStageRegister("Solve", &stage));
2493c611800SMark Adams   i    = dim;
250*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &i, NULL));
2513c611800SMark Adams   i    = dim;
252*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetIntArray(NULL, NULL, "-np", Np,  &i, NULL));
2533c611800SMark Adams   /* Create thread meshes */
2543c611800SMark Adams   for (int tid=0; tid<numthreads; tid++) {
2553c611800SMark Adams     // setup mesh dm_t, could use PETSc's Landau create velocity space mesh here to get dm_t[tid]
256*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMCreate(PETSC_COMM_SELF, &dm_t[tid]));
257*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSetType(dm_t[tid], DMPLEX));
258*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSetFromOptions(dm_t[tid]));
259*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFECreateDefault(PETSC_COMM_SELF, dim, Nc, PETSC_FALSE, "", PETSC_DECIDE, &fe));
260*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFESetFromOptions(fe));
261*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectSetName((PetscObject)fe, "fe"));
262*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSetField(dm_t[tid], field, NULL, (PetscObject)fe));
263*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMCreateDS(dm_t[tid]));
264*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFEDestroy(&fe));
2653c611800SMark Adams     // helper vectors
266*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSetOutputSequenceNumber(dm_t[tid], timestep, time)); // not used
267*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMCreateGlobalVector(dm_t[tid], &rho_t[tid]));
268*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMCreateGlobalVector(dm_t[tid], &rhs_t[tid]));
2693c611800SMark Adams     // this mimics application code
270*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetBoundingBox(dm_t[tid], lo, hi));
2713c611800SMark Adams     if (tid==target) {
272*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMViewFromOptions(dm_t[tid], NULL, "-dm_view"));
2733c611800SMark Adams       for (i=0,vol=1;i<dim;i++) {
2743c611800SMark Adams         h[i] = (hi[i] - lo[i])/faces[i];
2753c611800SMark Adams         hp[i] = (hi[i] - lo[i])/Np[i];
2763c611800SMark Adams         vol *= (hi[i] - lo[i]);
277*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscInfo(dm_t[tid]," lo = %g hi = %g n = %D h = %g hp = %g\n",lo[i],hi[i],faces[i],h[i],hp[i]));
2783c611800SMark Adams       }
2793c611800SMark Adams     }
2803c611800SMark Adams   }
2813c611800SMark Adams   // prepare particle data for problems. This mimics application code
282*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventBegin(swarm_create_ev,0,0,0,0));
2833c611800SMark Adams   Np2[0] = Np[0]; Np2[1] = Np[1];
2843c611800SMark Adams   for (int tid=0; tid<numthreads; tid++) { // change size of particle list a little
2853c611800SMark Adams     Np_t[tid] = Np2[0]*Np2[1];
286*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc3(Np_t[tid],&xx_t[tid],Np_t[tid],&yy_t[tid],Np_t[tid],&wp_t[tid]));
2873c611800SMark Adams     if (tid==target) {moments_0[0] = moments_0[1] = moments_0[2] = 0;}
2883c611800SMark Adams     for (int pi=0, pp=0;pi<Np2[0];pi++) {
2893c611800SMark Adams       for (int pj=0;pj<Np2[1];pj++,pp++) {
2903c611800SMark Adams         xx_t[tid][pp] = lo[0] + hp[0]/2. + pi*hp[0];
2913c611800SMark Adams         yy_t[tid][pp] = lo[1] + hp[1]/2. + pj*hp[1];
2923c611800SMark Adams         {
2933c611800SMark Adams           PetscReal x[] = {xx_t[tid][pp],yy_t[tid][pp]};
294*5f80ce2aSJacob Faibussowitsch           CHKERRQ(maxwellian(2, x, 1.0, vol/(PetscReal)Np_t[tid], &wp_t[tid][pp]));
2953c611800SMark Adams         }
2963c611800SMark Adams         if (tid==target) { //energy_0 += wp_t[tid][pp]*(PetscSqr(xx_t[tid][pp])+PetscSqr(yy_t[tid][pp]));
2973c611800SMark Adams           moments_0[0] += wp_t[tid][pp];
2983c611800SMark Adams           moments_0[1] += wp_t[tid][pp] * xx_t[tid][pp]; // x-momentum
2993c611800SMark Adams           moments_0[2] += wp_t[tid][pp] * (PetscSqr(xx_t[tid][pp]) + PetscSqr(yy_t[tid][pp]));
3003c611800SMark Adams         }
3013c611800SMark Adams       }
3023c611800SMark Adams     }
3033c611800SMark Adams     Np2[0]++; Np2[1]++;
3043c611800SMark Adams   }
305*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventEnd(swarm_create_ev,0,0,0,0));
306*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventBegin(solve_ev,0,0,0,0));
3073c611800SMark Adams   /* Create particle swarm */
3083c611800SMark Adams   PetscPragmaOMP(parallel for)
3093c611800SMark Adams   for (int tid=0; tid<numthreads; tid++) {
310*5f80ce2aSJacob Faibussowitsch     CHKERRABORT(PETSC_COMM_SELF,createSwarm(dm_t[tid], &sw_t[tid]));
3113c611800SMark Adams   }
3123c611800SMark Adams   PetscPragmaOMP(parallel for)
3133c611800SMark Adams   for (int tid=0; tid<numthreads; tid++) {
314*5f80ce2aSJacob Faibussowitsch     CHKERRABORT(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]));
3153c611800SMark Adams   }
3163c611800SMark Adams   /* Project field to particles */
3173c611800SMark Adams   /*   This gives f_p = M_p^+ M f */
3183c611800SMark Adams   PetscPragmaOMP(parallel for)
3193c611800SMark Adams   for (int tid=0; tid<numthreads; tid++) {
320*5f80ce2aSJacob Faibussowitsch     CHKERRABORT(PETSC_COMM_SELF,VecCopy(rho_t[tid], rhs_t[tid])); /* Identity: M^1 M rho */
3213c611800SMark Adams   }
3223c611800SMark Adams   PetscPragmaOMP(parallel for)
3233c611800SMark Adams   for (int tid=0; tid<numthreads; tid++) {
324*5f80ce2aSJacob Faibussowitsch     CHKERRABORT(PETSC_COMM_SELF,gridToParticles(dm_t[tid], sw_t[tid], (tid==target) ?  moments_1 : NULL, rhs_t[tid], M_p_t[tid]));
3253c611800SMark Adams   }
3263c611800SMark Adams   /* Cleanup */
3273c611800SMark Adams   for (int tid=0; tid<numthreads; tid++) {
328*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&M_p_t[tid]));
329*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDestroy(&sw_t[tid]));
3303c611800SMark Adams   }
331*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventEnd(solve_ev,0,0,0,0));
3323c611800SMark Adams   /* for timing */
333*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsClearValue(NULL,"-ftop_ksp_converged_reason"));
334*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsClearValue(NULL,"-ftop_ksp_monitor"));
335*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsClearValue(NULL,"-ftop_ksp_view"));
336*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsClearValue(NULL,"-info"));
3373c611800SMark Adams   // repeat solve many times to get warmed up data
338*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogStagePush(stage));
3393c611800SMark Adams #if defined(PETSC_HAVE_OPENMP) && defined(PETSC_HAVE_THREADSAFETY)
3403c611800SMark Adams   starttime = MPI_Wtime();
3413c611800SMark Adams #endif
342*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventBegin(solve_loop_ev,0,0,0,0));
3433c611800SMark Adams   for (int d=0; d<NUM_SOLVE_LOOPS; d++) {
3443c611800SMark Adams   /* Create particle swarm */
3453c611800SMark Adams     PetscPragmaOMP(parallel for)
3463c611800SMark Adams     for (int tid=0; tid<numthreads; tid++) {
347*5f80ce2aSJacob Faibussowitsch       CHKERRABORT(PETSC_COMM_SELF,createSwarm(dm_t[tid], &sw_t[tid]));
3483c611800SMark Adams     }
3493c611800SMark Adams     PetscPragmaOMP(parallel for)
3503c611800SMark Adams     for (int tid=0; tid<numthreads; tid++) {
351*5f80ce2aSJacob Faibussowitsch       CHKERRABORT(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]));
3523c611800SMark Adams     }
3533c611800SMark Adams     PetscPragmaOMP(parallel for)
3543c611800SMark Adams     for (int tid=0; tid<numthreads; tid++) {
355*5f80ce2aSJacob Faibussowitsch       CHKERRABORT(PETSC_COMM_SELF,VecCopy(rho_t[tid], rhs_t[tid])); /* Identity: M^1 M rho */
3563c611800SMark Adams     }
3573c611800SMark Adams     PetscPragmaOMP(parallel for)
3583c611800SMark Adams     for (int tid=0; tid<numthreads; tid++) {
359*5f80ce2aSJacob Faibussowitsch       CHKERRABORT(PETSC_COMM_SELF,gridToParticles(dm_t[tid], sw_t[tid], NULL, rhs_t[tid], M_p_t[tid]));
3603c611800SMark Adams     }
3613c611800SMark Adams     /* Cleanup */
3623c611800SMark Adams     for (int tid=0; tid<numthreads; tid++) {
363*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatDestroy(&M_p_t[tid]));
364*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMDestroy(&sw_t[tid]));
3653c611800SMark Adams     }
3663c611800SMark Adams   }
367*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventEnd(solve_loop_ev,0,0,0,0));
3683c611800SMark Adams #if defined(PETSC_HAVE_OPENMP) && defined(PETSC_HAVE_THREADSAFETY)
3693c611800SMark Adams   endtime = MPI_Wtime();
3703c611800SMark Adams   solve_time += (endtime - starttime);
3713c611800SMark Adams #endif
372*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogStagePop());
3733c611800SMark Adams   //
374*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Total number density: %20.12e (%20.12e); x-momentum = %g (%g); energy = %g error = %e, %D particles. Use %D threads, Solve time: %g\n", moments_1[0], moments_0[0], moments_1[1], moments_0[1], moments_1[2], (moments_1[2]-moments_0[2])/moments_0[2],Np[0]*Np[1],numthreads,solve_time));
3753c611800SMark Adams   /* Cleanup */
3763c611800SMark Adams   for (int tid=0; tid<numthreads; tid++) {
377*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&rho_t[tid]));
378*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&rhs_t[tid]));
379*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDestroy(&dm_t[tid]));
380*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree3(xx_t[tid],yy_t[tid],wp_t[tid]));
3813c611800SMark Adams   }
3823c611800SMark Adams   PetscFunctionReturn(0);
3833c611800SMark Adams }
3843c611800SMark Adams 
3853c611800SMark Adams int main(int argc, char **argv)
3863c611800SMark Adams {
3873c611800SMark Adams   PetscErrorCode  ierr;
3883c611800SMark Adams   ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
389*5f80ce2aSJacob Faibussowitsch   CHKERRQ(go());
3903c611800SMark Adams   ierr = PetscFinalize();
3913c611800SMark Adams   return ierr;
3923c611800SMark Adams }
3933c611800SMark Adams 
3943c611800SMark Adams /*TEST
3953c611800SMark Adams 
3963c611800SMark Adams   build:
3973c611800SMark Adams     requires: !complex
3983c611800SMark Adams 
3993c611800SMark Adams   test:
4003c611800SMark Adams     suffix: 0
4013c611800SMark Adams     requires: double triangle
4023c611800SMark Adams     args: -dm_plex_simplex 0 -dm_plex_box_faces 4,2 -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
4033c611800SMark Adams     filter: grep -v DM_ | grep -v atomic
4043c611800SMark Adams 
4053c611800SMark Adams   test:
4063c611800SMark Adams     suffix: 1
4073c611800SMark Adams     requires: double triangle
4083c611800SMark Adams     args: -dm_plex_simplex 0 -dm_plex_box_faces 4,2 -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
4093c611800SMark Adams     filter: grep -v DM_ | grep -v atomic
4103c611800SMark Adams 
4113c611800SMark Adams   test:
4123c611800SMark Adams     suffix: 2
4133c611800SMark Adams     requires: double triangle
4143c611800SMark Adams     args: -dm_plex_simplex 0 -dm_plex_box_faces 4,2 -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 cg -ftop_pc_type jacobi -dm_view -ftop_ksp_converged_reason -ftop_ksp_rtol 1.e-14
4153c611800SMark Adams     filter: grep -v DM_ | grep -v atomic
4163c611800SMark Adams 
4173c611800SMark Adams TEST*/
418