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; PetscErrorCode ierr; PetscFunctionBeginUser; ierr = MatShellGetContext(MtM,&matshellctx);CHKERRQ(ierr); if (!matshellctx) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "No context"); ierr = MatMult(matshellctx->Mp, xx, matshellctx->ff);CHKERRQ(ierr); ierr = MatMult(matshellctx->MpTrans, matshellctx->ff, yy);CHKERRQ(ierr); PetscFunctionReturn(0); } PetscErrorCode MatMultAddMtM_SeqAIJ(Mat MtM,Vec xx, Vec yy, Vec zz) { MatShellCtx *matshellctx; PetscErrorCode ierr; PetscFunctionBeginUser; ierr = MatShellGetContext(MtM,&matshellctx);CHKERRQ(ierr); if (!matshellctx) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "No context"); ierr = MatMult(matshellctx->Mp, xx, matshellctx->ff);CHKERRQ(ierr); ierr = MatMultAdd(matshellctx->MpTrans, matshellctx->ff, yy, zz);CHKERRQ(ierr); PetscFunctionReturn(0); } PetscErrorCode createSwarm(const DM dm, DM *sw) { PetscErrorCode ierr; PetscInt Nc = 1, dim = 2; PetscFunctionBeginUser; ierr = DMCreate(PETSC_COMM_SELF, sw);CHKERRQ(ierr); ierr = DMSetType(*sw, DMSWARM);CHKERRQ(ierr); ierr = DMSetDimension(*sw, dim);CHKERRQ(ierr); ierr = DMSwarmSetType(*sw, DMSWARM_PIC);CHKERRQ(ierr); ierr = DMSwarmSetCellDM(*sw, dm);CHKERRQ(ierr); ierr = DMSwarmRegisterPetscDatatypeField(*sw, "w_q", Nc, PETSC_SCALAR);CHKERRQ(ierr); ierr = DMSwarmFinalizeFieldRegister(*sw);CHKERRQ(ierr); ierr = DMSetFromOptions(*sw);CHKERRQ(ierr); 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; PetscErrorCode ierr; PetscInt Np, timestep = 0, bs, N, M, nzl; PetscReal time = 0.0; PetscDataType dtype; MatShellCtx *matshellctx; PetscFunctionBeginUser; ierr = KSPCreate(PETSC_COMM_SELF, &ksp);CHKERRQ(ierr); ierr = KSPSetOptionsPrefix(ksp, "ftop_");CHKERRQ(ierr); ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr); ierr = PetscObjectTypeCompare((PetscObject)ksp,KSPLSQR,&is_lsqr); if (!is_lsqr) { ierr = MatGetLocalSize(M_p, &M, &N);CHKERRQ(ierr); if (N>M) { PC pc; ierr = PetscInfo2(ksp, " M (%D) < M (%D) -- skip revert to lsqr\n",M,N);CHKERRQ(ierr); is_lsqr = PETSC_TRUE; ierr = KSPSetType(ksp,KSPLSQR);CHKERRQ(ierr); ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); // could put in better solver -ftop_pc_type bjacobi -ftop_sub_pc_type lu -ftop_sub_pc_factor_shift_type nonzero } else { ierr = PetscNew(&matshellctx);CHKERRQ(ierr); ierr = MatCreateShell(PetscObjectComm((PetscObject)dm),N,N,PETSC_DECIDE,PETSC_DECIDE,matshellctx,&MtM);CHKERRQ(ierr); ierr = MatTranspose(M_p,MAT_INITIAL_MATRIX,&matshellctx->MpTrans);CHKERRQ(ierr); matshellctx->Mp = M_p; ierr = MatShellSetOperation(MtM, MATOP_MULT, (void (*)(void))MatMultMtM_SeqAIJ);CHKERRQ(ierr); ierr = MatShellSetOperation(MtM, MATOP_MULT_ADD, (void (*)(void))MatMultAddMtM_SeqAIJ);CHKERRQ(ierr); ierr = MatCreateVecs(M_p,&matshellctx->uu,&matshellctx->ff);CHKERRQ(ierr); ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,N,N,1,NULL,&D);CHKERRQ(ierr); for (int i=0 ; iMpTrans,i,&nzl,&cols,&vals);CHKERRQ(ierr); for (int ii=0 ; iiMpTrans,NULL,"-ftop2_MpTranspose_mat_view");CHKERRQ(ierr); } } if (is_lsqr) { PC pc; PetscBool is_bjac; ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); ierr = PetscObjectTypeCompare((PetscObject)pc,PCBJACOBI,&is_bjac); if (is_bjac) { ierr = DMSwarmCreateMassMatrixSquare(sw, dm, &PM_p);CHKERRQ(ierr); ierr = KSPSetOperators(ksp, M_p, PM_p);CHKERRQ(ierr); } else { ierr = KSPSetOperators(ksp, M_p, M_p);CHKERRQ(ierr); } } ierr = DMSwarmCreateGlobalVectorFromField(sw, "w_q", &ff);CHKERRQ(ierr); // this grabs access !!!!! if (!is_lsqr) { ierr = KSPSolve(ksp, rhs, matshellctx->uu);CHKERRQ(ierr); ierr = MatMult(M_p, matshellctx->uu, ff);CHKERRQ(ierr); ierr = MatDestroy(&matshellctx->MpTrans);CHKERRQ(ierr); ierr = VecDestroy(&matshellctx->ff);CHKERRQ(ierr); ierr = VecDestroy(&matshellctx->uu);CHKERRQ(ierr); ierr = MatDestroy(&D);CHKERRQ(ierr); ierr = MatDestroy(&MtM);CHKERRQ(ierr); ierr = PetscFree(matshellctx);CHKERRQ(ierr); } else { ierr = KSPSolveTranspose(ksp, rhs, ff);CHKERRQ(ierr); } ierr = KSPDestroy(&ksp);CHKERRQ(ierr); /* Visualize particle field */ ierr = DMSetOutputSequenceNumber(sw, timestep, time);CHKERRQ(ierr); ierr = VecViewFromOptions(ff, NULL, "-weights_view");CHKERRQ(ierr); ierr = DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &ff);CHKERRQ(ierr); /* compute energy */ if (moments) { PetscReal *wq, *coords; ierr = DMSwarmGetLocalSize(sw,&Np);CHKERRQ(ierr); ierr = DMSwarmGetField(sw, "w_q", &bs, &dtype, (void**)&wq);CHKERRQ(ierr); ierr = DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void**)&coords);CHKERRQ(ierr); moments[0] = moments[1] = moments[2] = 0; for (int p=0;pMAX_NUM_THRDS) SETERRQ2(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Too many threads %D > %D", numthreads, MAX_NUM_THRDS); if (numthreads<=0) SETERRQ2(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "No threads %D > %D ", numthreads, MAX_NUM_THRDS); #endif if (target >= numthreads) target = numthreads-1; ierr = PetscLogEventRegister("Create Swarm", DM_CLASSID, &swarm_create_ev);CHKERRQ(ierr); ierr = PetscLogEventRegister("Single solve", DM_CLASSID, &solve_ev);CHKERRQ(ierr); ierr = PetscLogEventRegister("Solve loop", DM_CLASSID, &solve_loop_ev);CHKERRQ(ierr); ierr = PetscLogStageRegister("Solve", &stage);CHKERRQ(ierr); i = dim; ierr = PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &i, NULL);CHKERRQ(ierr); i = dim; ierr = PetscOptionsGetIntArray(NULL, NULL, "-np", Np, &i, NULL);CHKERRQ(ierr); /* Create thread meshes */ for (int tid=0; tid