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; 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)); 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; 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)); 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; 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)); 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; 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; 799566063dSJacob Faibussowitsch PetscCall(PetscInfo(ksp, " M (%D) < M (%D) -- 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; 899566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(MtM, MATOP_MULT, (void (*)(void))MatMultMtM_SeqAIJ)); 909566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(MtM, MATOP_MULT_ADD, (void (*)(void))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]); 9908401ef6SPierre Jolivet 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)); 1047d3de750SJacob Faibussowitsch PetscInfo(M_p,"createMtMKSP Have %D eqs, nzl = %D\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)); 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; 1739566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(sw, Np, zero)); 1749566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void**)&wq)); 1759566063dSJacob Faibussowitsch PetscCall(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 } 1819566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void**)&coords)); 1829566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void**)&wq)); 1839566063dSJacob Faibussowitsch PetscCall(DMSwarmMigrate(sw, removePoints)); 1849566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)sw, "Particle Grid")); 1859566063dSJacob Faibussowitsch if (a_tid==target) PetscCall(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 */ 1899566063dSJacob Faibussowitsch PetscCall(DMCreateMassMatrix(sw, dm, &M_p)); 1909566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)rho, "rho")); 1913c611800SMark Adams 1929566063dSJacob Faibussowitsch PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &ff)); // this grabs access !!!!! 1939566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)ff, "weights")); 1949566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(M_p, ff, rho)); 1959566063dSJacob Faibussowitsch PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &ff)); 1963c611800SMark Adams 1973c611800SMark Adams /* Visualize mesh field */ 1989566063dSJacob Faibussowitsch if (a_tid==target) PetscCall(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 } 216*28114503SMark Adams 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; 222*28114503SMark Adams PetscInt dim = 2, Nc = 1, i, faces[3]; 2233c611800SMark Adams PetscInt Np[2] = {10,10}, Np2[2], field = 0, target = 0, Np_t[MAX_NUM_THRDS]; 224*28114503SMark Adams PetscReal moments_0[3], moments_1[3], vol = 1; 225*28114503SMark 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]; 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 #else 2353c611800SMark Adams PetscInt numthreads = 1; 2363c611800SMark Adams #endif 2373c611800SMark Adams 2383c611800SMark Adams PetscFunctionBeginUser; 2393c611800SMark Adams #if defined(PETSC_HAVE_OPENMP) && defined(PETSC_HAVE_THREADSAFETY) 24008401ef6SPierre Jolivet PetscCheck(numthreads<=MAX_NUM_THRDS,PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Too many threads %D > %D", numthreads, MAX_NUM_THRDS); 24108401ef6SPierre Jolivet PetscCheck(numthreads>0,PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "No threads %D > %D ", numthreads, MAX_NUM_THRDS); 2423c611800SMark Adams #endif 2433c611800SMark Adams if (target >= numthreads) target = numthreads-1; 2449566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister("Create Swarm", DM_CLASSID, &swarm_create_ev)); 2459566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister("Single solve", DM_CLASSID, &solve_ev)); 2469566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister("Solve loop", DM_CLASSID, &solve_loop_ev)); 2479566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("Solve", &stage)); 2483c611800SMark Adams i = dim; 2499566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &i, NULL)); 2503c611800SMark Adams i = dim; 2519566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-np", Np, &i, NULL)); 2523c611800SMark Adams /* Create thread meshes */ 2533c611800SMark Adams for (int tid=0; tid<numthreads; tid++) { 2543c611800SMark Adams // setup mesh dm_t, could use PETSc's Landau create velocity space mesh here to get dm_t[tid] 2559566063dSJacob Faibussowitsch PetscCall(DMCreate(PETSC_COMM_SELF, &dm_t[tid])); 2569566063dSJacob Faibussowitsch PetscCall(DMSetType(dm_t[tid], DMPLEX)); 2579566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dm_t[tid])); 2589566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, Nc, PETSC_FALSE, "", PETSC_DECIDE, &fe)); 2599566063dSJacob Faibussowitsch PetscCall(PetscFESetFromOptions(fe)); 2609566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe, "fe")); 2619566063dSJacob Faibussowitsch PetscCall(DMSetField(dm_t[tid], field, NULL, (PetscObject)fe)); 2629566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm_t[tid])); 2639566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 2643c611800SMark Adams // helper vectors 2659566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm_t[tid], &rho_t[tid])); 2669566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm_t[tid], &rhs_t[tid])); 2673c611800SMark Adams // this mimics application code 2689566063dSJacob Faibussowitsch PetscCall(DMGetBoundingBox(dm_t[tid], lo, hi)); 2693c611800SMark Adams if (tid==target) { 2709566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm_t[tid], NULL, "-dm_view")); 2713c611800SMark Adams for (i=0,vol=1;i<dim;i++) { 2723c611800SMark Adams h[i] = (hi[i] - lo[i])/faces[i]; 2733c611800SMark Adams hp[i] = (hi[i] - lo[i])/Np[i]; 2743c611800SMark Adams vol *= (hi[i] - lo[i]); 2759566063dSJacob Faibussowitsch PetscCall(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])); 2763c611800SMark Adams } 2773c611800SMark Adams } 2783c611800SMark Adams } 2793c611800SMark Adams // prepare particle data for problems. This mimics application code 2809566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(swarm_create_ev,0,0,0,0)); 2813c611800SMark Adams Np2[0] = Np[0]; 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])); 2853c611800SMark Adams 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 } 3013c611800SMark Adams Np2[0]++; Np2[1]++; 3023c611800SMark Adams } 3039566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(swarm_create_ev,0,0,0,0)); 3049566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(solve_ev,0,0,0,0)); 3053c611800SMark Adams /* Create particle swarm */ 3063c611800SMark Adams PetscPragmaOMP(parallel for) 3073c611800SMark Adams for (int tid=0; tid<numthreads; tid++) { 3089566063dSJacob Faibussowitsch PetscCallAbort(PETSC_COMM_SELF,createSwarm(dm_t[tid], &sw_t[tid])); 3093c611800SMark Adams } 3103c611800SMark Adams PetscPragmaOMP(parallel for) 3113c611800SMark Adams for (int tid=0; tid<numthreads; tid++) { 3129566063dSJacob Faibussowitsch 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])); 3133c611800SMark Adams } 3143c611800SMark Adams /* Project field to particles */ 3153c611800SMark Adams /* This gives f_p = M_p^+ M f */ 3163c611800SMark Adams PetscPragmaOMP(parallel for) 3173c611800SMark Adams for (int tid=0; tid<numthreads; tid++) { 3189566063dSJacob Faibussowitsch PetscCallAbort(PETSC_COMM_SELF,VecCopy(rho_t[tid], rhs_t[tid])); /* Identity: M^1 M rho */ 3193c611800SMark Adams } 3203c611800SMark Adams PetscPragmaOMP(parallel for) 3213c611800SMark Adams for (int tid=0; tid<numthreads; tid++) { 3229566063dSJacob Faibussowitsch PetscCallAbort(PETSC_COMM_SELF,gridToParticles(dm_t[tid], sw_t[tid], (tid==target) ? moments_1 : NULL, rhs_t[tid], M_p_t[tid])); 3233c611800SMark Adams } 3243c611800SMark Adams /* Cleanup */ 3253c611800SMark Adams for (int tid=0; tid<numthreads; tid++) { 3269566063dSJacob Faibussowitsch PetscCall(MatDestroy(&M_p_t[tid])); 3279566063dSJacob Faibussowitsch PetscCall(DMDestroy(&sw_t[tid])); 3283c611800SMark Adams } 3299566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(solve_ev,0,0,0,0)); 3303c611800SMark Adams // 331*28114503SMark Adams PetscCall(PetscPrintf(PETSC_COMM_SELF,"Total number density: %20.12e (%20.12e); x-momentum = %g (%g); energy = %g error = %e, %D particles. Use %D threads\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)); 3323c611800SMark Adams /* Cleanup */ 3333c611800SMark Adams for (int tid=0; tid<numthreads; tid++) { 3349566063dSJacob Faibussowitsch PetscCall(VecDestroy(&rho_t[tid])); 3359566063dSJacob Faibussowitsch PetscCall(VecDestroy(&rhs_t[tid])); 3369566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm_t[tid])); 3379566063dSJacob Faibussowitsch PetscCall(PetscFree3(xx_t[tid],yy_t[tid],wp_t[tid])); 3383c611800SMark Adams } 3393c611800SMark Adams PetscFunctionReturn(0); 3403c611800SMark Adams } 3413c611800SMark Adams 3423c611800SMark Adams int main(int argc, char **argv) 3433c611800SMark Adams { 3449566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL,help)); 3459566063dSJacob Faibussowitsch PetscCall(go()); 3469566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 347b122ec5aSJacob Faibussowitsch return 0; 3483c611800SMark Adams } 3493c611800SMark Adams 3503c611800SMark Adams /*TEST 3513c611800SMark Adams 3523c611800SMark Adams build: 3533c611800SMark Adams requires: !complex 3543c611800SMark Adams 3553c611800SMark Adams test: 3563c611800SMark Adams suffix: 0 3573c611800SMark Adams requires: double triangle 358*28114503SMark Adams args: -dm_plex_simplex 0 -dm_plex_box_faces 8,4 -np 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 3593c611800SMark Adams filter: grep -v DM_ | grep -v atomic 3603c611800SMark Adams 3613c611800SMark Adams test: 3623c611800SMark Adams suffix: 1 3633c611800SMark Adams requires: double triangle 364*28114503SMark Adams args: -dm_plex_simplex 0 -dm_plex_box_faces 8,4 -np 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 3653c611800SMark Adams filter: grep -v DM_ | grep -v atomic 3663c611800SMark Adams 3673c611800SMark Adams test: 3683c611800SMark Adams suffix: 2 3693c611800SMark Adams requires: double triangle 370*28114503SMark Adams args: -dm_plex_simplex 0 -dm_plex_box_faces 8,4 -np 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 cg -ftop_pc_type jacobi -dm_view -ftop_ksp_converged_reason -ftop_ksp_rtol 1.e-14 3713c611800SMark Adams filter: grep -v DM_ | grep -v atomic 3723c611800SMark Adams 3733c611800SMark Adams TEST*/ 374