static char help[] = "Example program demonstrating projection between particle and finite element spaces using OpenMP in 2D cylindrical coordinates\n"; #include "petscdmplex.h" #include "petscds.h" #include "petscdmswarm.h" #include "petscksp.h" #include #if defined(PETSC_HAVE_OPENMP) && defined(PETSC_HAVE_THREADSAFETY) #include #endif typedef struct { Mat MpTrans; Mat Mp; Vec ff; Vec uu; } MatShellCtx; PetscErrorCode MatMultMtM_SeqAIJ(Mat MtM,Vec xx,Vec yy) { MatShellCtx *matshellctx; PetscFunctionBeginUser; PetscCall(MatShellGetContext(MtM,&matshellctx)); PetscCheck(matshellctx,PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "No context"); PetscCall(MatMult(matshellctx->Mp, xx, matshellctx->ff)); PetscCall(MatMult(matshellctx->MpTrans, matshellctx->ff, yy)); PetscFunctionReturn(0); } PetscErrorCode MatMultAddMtM_SeqAIJ(Mat MtM,Vec xx, Vec yy, Vec zz) { MatShellCtx *matshellctx; PetscFunctionBeginUser; PetscCall(MatShellGetContext(MtM,&matshellctx)); PetscCheck(matshellctx,PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "No context"); PetscCall(MatMult(matshellctx->Mp, xx, matshellctx->ff)); PetscCall(MatMultAdd(matshellctx->MpTrans, matshellctx->ff, yy, zz)); PetscFunctionReturn(0); } PetscErrorCode createSwarm(const DM dm, DM *sw) { PetscInt Nc = 1, dim = 2; PetscFunctionBeginUser; PetscCall(DMCreate(PETSC_COMM_SELF, sw)); PetscCall(DMSetType(*sw, DMSWARM)); PetscCall(DMSetDimension(*sw, dim)); PetscCall(DMSwarmSetType(*sw, DMSWARM_PIC)); PetscCall(DMSwarmSetCellDM(*sw, dm)); PetscCall(DMSwarmRegisterPetscDatatypeField(*sw, "w_q", Nc, PETSC_SCALAR)); PetscCall(DMSwarmFinalizeFieldRegister(*sw)); PetscCall(DMSetFromOptions(*sw)); PetscFunctionReturn(0); } PetscErrorCode gridToParticles(const DM dm, DM sw, PetscReal *moments, Vec rhs, Mat M_p) { PetscBool is_lsqr; KSP ksp; Mat PM_p=NULL,MtM,D; Vec ff; PetscInt Np, timestep = 0, bs, N, M, nzl; PetscReal time = 0.0; PetscDataType dtype; MatShellCtx *matshellctx; PetscFunctionBeginUser; PetscCall(KSPCreate(PETSC_COMM_SELF, &ksp)); PetscCall(KSPSetOptionsPrefix(ksp, "ftop_")); PetscCall(KSPSetFromOptions(ksp)); PetscCall(PetscObjectTypeCompare((PetscObject)ksp,KSPLSQR,&is_lsqr)); if (!is_lsqr) { PetscCall(MatGetLocalSize(M_p, &M, &N)); if (N>M) { PC pc; PetscCall(PetscInfo(ksp, " M (%D) < M (%D) -- skip revert to lsqr\n",M,N)); is_lsqr = PETSC_TRUE; PetscCall(KSPSetType(ksp,KSPLSQR)); PetscCall(KSPGetPC(ksp,&pc)); 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 } else { PetscCall(PetscNew(&matshellctx)); PetscCall(MatCreateShell(PetscObjectComm((PetscObject)dm),N,N,PETSC_DECIDE,PETSC_DECIDE,matshellctx,&MtM)); PetscCall(MatTranspose(M_p,MAT_INITIAL_MATRIX,&matshellctx->MpTrans)); matshellctx->Mp = M_p; PetscCall(MatShellSetOperation(MtM, MATOP_MULT, (void (*)(void))MatMultMtM_SeqAIJ)); PetscCall(MatShellSetOperation(MtM, MATOP_MULT_ADD, (void (*)(void))MatMultAddMtM_SeqAIJ)); PetscCall(MatCreateVecs(M_p,&matshellctx->uu,&matshellctx->ff)); PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF,N,N,1,NULL,&D)); for (int i=0 ; iMpTrans,i,&nzl,&cols,&vals)); for (int ii=0 ; iiMpTrans,NULL,"-ftop2_MpTranspose_mat_view")); } } if (is_lsqr) { PC pc; PetscBool is_bjac; PetscCall(KSPGetPC(ksp,&pc)); PetscCall(PetscObjectTypeCompare((PetscObject)pc,PCBJACOBI,&is_bjac)); if (is_bjac) { PetscCall(DMSwarmCreateMassMatrixSquare(sw, dm, &PM_p)); PetscCall(KSPSetOperators(ksp, M_p, PM_p)); } else { PetscCall(KSPSetOperators(ksp, M_p, M_p)); } } PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &ff)); // this grabs access !!!!! if (!is_lsqr) { PetscCall(KSPSolve(ksp, rhs, matshellctx->uu)); PetscCall(MatMult(M_p, matshellctx->uu, ff)); PetscCall(MatDestroy(&matshellctx->MpTrans)); PetscCall(VecDestroy(&matshellctx->ff)); PetscCall(VecDestroy(&matshellctx->uu)); PetscCall(MatDestroy(&D)); PetscCall(MatDestroy(&MtM)); PetscCall(PetscFree(matshellctx)); } else { PetscCall(KSPSolveTranspose(ksp, rhs, ff)); } PetscCall(KSPDestroy(&ksp)); /* Visualize particle field */ PetscCall(DMSetOutputSequenceNumber(sw, timestep, time)); PetscCall(VecViewFromOptions(ff, NULL, "-weights_view")); PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &ff)); /* compute energy */ if (moments) { PetscReal *wq, *coords; PetscCall(DMSwarmGetLocalSize(sw,&Np)); PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void**)&wq)); PetscCall(DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void**)&coords)); moments[0] = moments[1] = moments[2] = 0; for (int p=0;pMAX_NUM_THRDS,PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Too many threads %D > %D", numthreads, MAX_NUM_THRDS); PetscCheckFalse(numthreads<=0,PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "No threads %D > %D ", numthreads, MAX_NUM_THRDS); #endif if (target >= numthreads) target = numthreads-1; PetscCall(PetscLogEventRegister("Create Swarm", DM_CLASSID, &swarm_create_ev)); PetscCall(PetscLogEventRegister("Single solve", DM_CLASSID, &solve_ev)); PetscCall(PetscLogEventRegister("Solve loop", DM_CLASSID, &solve_loop_ev)); PetscCall(PetscLogStageRegister("Solve", &stage)); i = dim; PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &i, NULL)); i = dim; PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-np", Np, &i, NULL)); /* Create thread meshes */ for (int tid=0; tid