14cd81d59SMatthew G. Knepley static char help[] = "Example program demonstrating projection between particle and finite element spaces\n\n"; 24cd81d59SMatthew G. Knepley 34cd81d59SMatthew G. Knepley #include "petscdmplex.h" 44cd81d59SMatthew G. Knepley #include "petscds.h" 54cd81d59SMatthew G. Knepley #include "petscdmswarm.h" 64cd81d59SMatthew G. Knepley #include "petscksp.h" 74cd81d59SMatthew G. Knepley 84cd81d59SMatthew G. Knepley int main(int argc, char **argv) 94cd81d59SMatthew G. Knepley { 103c611800SMark Adams DM dm, sw; 114cd81d59SMatthew G. Knepley PetscFE fe; 124cd81d59SMatthew G. Knepley KSP ksp; 133c611800SMark Adams PC pc; 143c611800SMark Adams Mat M_p, PM_p=NULL; 153c611800SMark Adams Vec f, rho, rhs; 163c611800SMark Adams PetscInt dim, Nc = 1, timestep = 0, i, faces[3]; 1730602db0SMatthew G. Knepley PetscInt Np = 10, p, field = 0, zero = 0, bs; 183c611800SMark Adams PetscReal time = 0.0, norm, energy_0, energy_1; 1930602db0SMatthew G. Knepley PetscReal lo[3], hi[3], h[3]; 2030602db0SMatthew G. Knepley PetscBool removePoints = PETSC_TRUE; 2130602db0SMatthew G. Knepley PetscReal *wq, *coords; 224cd81d59SMatthew G. Knepley PetscDataType dtype; 2330602db0SMatthew G. Knepley PetscErrorCode ierr; 243c611800SMark Adams PetscBool is_bjac; 254cd81d59SMatthew G. Knepley 264cd81d59SMatthew G. Knepley ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr; 274cd81d59SMatthew G. Knepley /* Create a mesh */ 28*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreate(PETSC_COMM_WORLD, &dm)); 29*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetType(dm, DMPLEX)); 30*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(dm)); 31*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(dm, NULL, "-dm_view")); 3230602db0SMatthew G. Knepley 33*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(dm, &dim)); 3430602db0SMatthew G. Knepley i = dim; 35*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &i, NULL)); 36*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL, NULL, "-np", &Np, NULL)); 37*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetBoundingBox(dm, lo, hi)); 3830602db0SMatthew G. Knepley for (i=0;i<dim;i++) { 3930602db0SMatthew G. Knepley h[i] = (hi[i] - lo[i])/faces[i]; 40*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF," lo = %g hi = %g n = %D h = %g\n",lo[i],hi[i],faces[i],h[i])); 4130602db0SMatthew G. Knepley } 4230602db0SMatthew G. Knepley 43*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFECreateDefault(PETSC_COMM_SELF, dim, Nc, PETSC_FALSE, "", PETSC_DECIDE, &fe)); 44*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFESetFromOptions(fe)); 45*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)fe, "fe")); 46*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetField(dm, field, NULL, (PetscObject)fe)); 47*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateDS(dm)); 48*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFEDestroy(&fe)); 494cd81d59SMatthew G. Knepley /* Create particle swarm */ 50*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreate(PETSC_COMM_SELF, &sw)); 51*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetType(sw, DMSWARM)); 52*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetDimension(sw, dim)); 53*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmSetType(sw, DMSWARM_PIC)); 54*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmSetCellDM(sw, dm)); 55*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmRegisterPetscDatatypeField(sw, "w_q", Nc, PETSC_SCALAR)); 56*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmFinalizeFieldRegister(sw)); 57*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmSetLocalSizes(sw, Np, zero)); 58*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(sw)); 59*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void**)&wq)); 60*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void**)&coords)); 613c611800SMark Adams for (p=0,energy_0=0;p<Np;p++) { 624cd81d59SMatthew G. Knepley coords[p*2+0] = -PetscCosReal((PetscReal)(p+1)/(PetscReal)(Np+1) * PETSC_PI); 634cd81d59SMatthew G. Knepley coords[p*2+1] = PetscSinReal((PetscReal)(p+1)/(PetscReal)(Np+1) * PETSC_PI); 644cd81d59SMatthew G. Knepley wq[p] = 1.0; 653c611800SMark Adams energy_0 += wq[p]*(PetscSqr(coords[p*2+0])+PetscSqr(coords[p*2+1])); 664cd81d59SMatthew G. Knepley } 67*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void**)&coords)); 68*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void**)&wq)); 69*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmMigrate(sw, removePoints)); 70*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)sw, "Particle Grid")); 71*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(sw, NULL, "-swarm_view")); 723c611800SMark Adams 734cd81d59SMatthew G. Knepley /* Project particles to field */ 744cd81d59SMatthew G. Knepley /* This gives M f = \int_\Omega \phi f, which looks like a rhs for a PDE */ 75*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateMassMatrix(sw, dm, &M_p)); 76*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(dm, &rho)); 77*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)rho, "rho")); 783c611800SMark Adams 79*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &f)); 80*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)f, "weights")); 81*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTranspose(M_p, f, rho)); 824cd81d59SMatthew G. Knepley 834cd81d59SMatthew G. Knepley /* Visualize mesh field */ 84*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetOutputSequenceNumber(dm, timestep, time)); 85*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecViewFromOptions(rho, NULL, "-rho_view")); 864cd81d59SMatthew G. Knepley 874cd81d59SMatthew G. Knepley /* Project field to particles */ 884cd81d59SMatthew G. Knepley /* This gives f_p = M_p^+ M f */ 89*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(dm, &rhs)); 90*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(rho, rhs)); /* Identity: M^1 M rho */ 913c611800SMark Adams 92*5f80ce2aSJacob Faibussowitsch CHKERRQ(KSPCreate(PETSC_COMM_WORLD, &ksp)); 93*5f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetOptionsPrefix(ksp, "ftop_")); 94*5f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetFromOptions(ksp)); 95*5f80ce2aSJacob Faibussowitsch CHKERRQ(KSPGetPC(ksp,&pc)); 963c611800SMark Adams ierr = PetscObjectTypeCompare((PetscObject)pc,PCBJACOBI,&is_bjac); 973c611800SMark Adams if (is_bjac) { 98*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmCreateMassMatrixSquare(sw, dm, &PM_p)); 99*5f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetOperators(ksp, M_p, PM_p)); 1003c611800SMark Adams } else { 101*5f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetOperators(ksp, M_p, M_p)); 1023c611800SMark Adams } 103*5f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSolveTranspose(ksp, rhs, f)); 104*5f80ce2aSJacob Faibussowitsch CHKERRQ(KSPDestroy(&ksp)); 105*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&rhs)); 1064cd81d59SMatthew G. Knepley 1074cd81d59SMatthew G. Knepley /* Visualize particle field */ 108*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetOutputSequenceNumber(sw, timestep, time)); 109*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecViewFromOptions(f, NULL, "-weights_view")); 110*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(f,NORM_1,&norm)); 111*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &f)); 1123c611800SMark Adams 1133c611800SMark Adams /* compute energy */ 114*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void**)&wq)); 115*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void**)&coords)); 1163c611800SMark Adams for (p=0,energy_1=0;p<Np;p++) { 1173c611800SMark Adams energy_1 += wq[p]*(PetscSqr(coords[p*2+0])+PetscSqr(coords[p*2+1])); 1183c611800SMark Adams } 119*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void**)&coords)); 120*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void**)&wq)); 121*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Total number = %20.12e. energy = %20.12e error = %20.12e\n", norm, energy_0, (energy_1-energy_0)/energy_0)); 1223c611800SMark Adams /* Cleanup */ 123*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&M_p)); 124*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&PM_p)); 125*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&rho)); 126*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&sw)); 127*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dm)); 1284cd81d59SMatthew G. Knepley ierr = PetscFinalize(); 1294cd81d59SMatthew G. Knepley return ierr; 1304cd81d59SMatthew G. Knepley } 1314cd81d59SMatthew G. Knepley 1324cd81d59SMatthew G. Knepley /*TEST 1334cd81d59SMatthew G. Knepley 1344cd81d59SMatthew G. Knepley build: 1354cd81d59SMatthew G. Knepley requires: !complex 1364cd81d59SMatthew G. Knepley 1374cd81d59SMatthew G. Knepley test: 1384cd81d59SMatthew G. Knepley suffix: 0 1393c611800SMark Adams requires: double triangle 1403c611800SMark Adams args: -dm_plex_simplex 0 -dm_plex_box_faces 4,2 -np 50 -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 -swarm_view -ftop_ksp_rtol 1.e-14 1413c611800SMark Adams filter: grep -v DM_ | grep -v atomic 1423c611800SMark Adams 1433c611800SMark Adams test: 1443c611800SMark Adams suffix: bjacobi 1453c611800SMark Adams requires: double triangle 1463c611800SMark Adams args: -dm_plex_simplex 0 -dm_plex_box_faces 4,2 -np 50 -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 -swarm_view -ftop_ksp_rtol 1.e-14 1474cd81d59SMatthew G. Knepley filter: grep -v DM_ | grep -v atomic 1484cd81d59SMatthew G. Knepley 1494cd81d59SMatthew G. Knepley TEST*/ 150