1 static char help[] = "Example program demonstrating projection between particle and finite element spaces\n\n"; 2 3 #include "petscdmplex.h" 4 #include "petscds.h" 5 #include "petscdmswarm.h" 6 #include "petscksp.h" 7 8 int main(int argc, char **argv) 9 { 10 DM dm, sw; 11 PetscFE fe; 12 KSP ksp; 13 PC pc; 14 Mat M_p, PM_p=NULL; 15 Vec f, rho, rhs; 16 PetscInt dim, Nc = 1, timestep = 0, i, faces[3]; 17 PetscInt Np = 10, p, field = 0, zero = 0, bs; 18 PetscReal time = 0.0, norm, energy_0, energy_1; 19 PetscReal lo[3], hi[3], h[3]; 20 PetscBool removePoints = PETSC_TRUE; 21 PetscReal *wq, *coords; 22 PetscDataType dtype; 23 PetscErrorCode ierr; 24 PetscBool is_bjac; 25 26 ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr; 27 /* Create a mesh */ 28 CHKERRQ(DMCreate(PETSC_COMM_WORLD, &dm)); 29 CHKERRQ(DMSetType(dm, DMPLEX)); 30 CHKERRQ(DMSetFromOptions(dm)); 31 CHKERRQ(DMViewFromOptions(dm, NULL, "-dm_view")); 32 33 CHKERRQ(DMGetDimension(dm, &dim)); 34 i = dim; 35 CHKERRQ(PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &i, NULL)); 36 CHKERRQ(PetscOptionsGetInt(NULL, NULL, "-np", &Np, NULL)); 37 CHKERRQ(DMGetBoundingBox(dm, lo, hi)); 38 for (i=0;i<dim;i++) { 39 h[i] = (hi[i] - lo[i])/faces[i]; 40 CHKERRQ(PetscPrintf(PETSC_COMM_SELF," lo = %g hi = %g n = %D h = %g\n",lo[i],hi[i],faces[i],h[i])); 41 } 42 43 CHKERRQ(PetscFECreateDefault(PETSC_COMM_SELF, dim, Nc, PETSC_FALSE, "", PETSC_DECIDE, &fe)); 44 CHKERRQ(PetscFESetFromOptions(fe)); 45 CHKERRQ(PetscObjectSetName((PetscObject)fe, "fe")); 46 CHKERRQ(DMSetField(dm, field, NULL, (PetscObject)fe)); 47 CHKERRQ(DMCreateDS(dm)); 48 CHKERRQ(PetscFEDestroy(&fe)); 49 /* Create particle swarm */ 50 CHKERRQ(DMCreate(PETSC_COMM_SELF, &sw)); 51 CHKERRQ(DMSetType(sw, DMSWARM)); 52 CHKERRQ(DMSetDimension(sw, dim)); 53 CHKERRQ(DMSwarmSetType(sw, DMSWARM_PIC)); 54 CHKERRQ(DMSwarmSetCellDM(sw, dm)); 55 CHKERRQ(DMSwarmRegisterPetscDatatypeField(sw, "w_q", Nc, PETSC_SCALAR)); 56 CHKERRQ(DMSwarmFinalizeFieldRegister(sw)); 57 CHKERRQ(DMSwarmSetLocalSizes(sw, Np, zero)); 58 CHKERRQ(DMSetFromOptions(sw)); 59 CHKERRQ(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void**)&wq)); 60 CHKERRQ(DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void**)&coords)); 61 for (p=0,energy_0=0;p<Np;p++) { 62 coords[p*2+0] = -PetscCosReal((PetscReal)(p+1)/(PetscReal)(Np+1) * PETSC_PI); 63 coords[p*2+1] = PetscSinReal((PetscReal)(p+1)/(PetscReal)(Np+1) * PETSC_PI); 64 wq[p] = 1.0; 65 energy_0 += wq[p]*(PetscSqr(coords[p*2+0])+PetscSqr(coords[p*2+1])); 66 } 67 CHKERRQ(DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void**)&coords)); 68 CHKERRQ(DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void**)&wq)); 69 CHKERRQ(DMSwarmMigrate(sw, removePoints)); 70 CHKERRQ(PetscObjectSetName((PetscObject)sw, "Particle Grid")); 71 CHKERRQ(DMViewFromOptions(sw, NULL, "-swarm_view")); 72 73 /* Project particles to field */ 74 /* This gives M f = \int_\Omega \phi f, which looks like a rhs for a PDE */ 75 CHKERRQ(DMCreateMassMatrix(sw, dm, &M_p)); 76 CHKERRQ(DMCreateGlobalVector(dm, &rho)); 77 CHKERRQ(PetscObjectSetName((PetscObject)rho, "rho")); 78 79 CHKERRQ(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &f)); 80 CHKERRQ(PetscObjectSetName((PetscObject)f, "weights")); 81 CHKERRQ(MatMultTranspose(M_p, f, rho)); 82 83 /* Visualize mesh field */ 84 CHKERRQ(DMSetOutputSequenceNumber(dm, timestep, time)); 85 CHKERRQ(VecViewFromOptions(rho, NULL, "-rho_view")); 86 87 /* Project field to particles */ 88 /* This gives f_p = M_p^+ M f */ 89 CHKERRQ(DMCreateGlobalVector(dm, &rhs)); 90 CHKERRQ(VecCopy(rho, rhs)); /* Identity: M^1 M rho */ 91 92 CHKERRQ(KSPCreate(PETSC_COMM_WORLD, &ksp)); 93 CHKERRQ(KSPSetOptionsPrefix(ksp, "ftop_")); 94 CHKERRQ(KSPSetFromOptions(ksp)); 95 CHKERRQ(KSPGetPC(ksp,&pc)); 96 ierr = PetscObjectTypeCompare((PetscObject)pc,PCBJACOBI,&is_bjac); 97 if (is_bjac) { 98 CHKERRQ(DMSwarmCreateMassMatrixSquare(sw, dm, &PM_p)); 99 CHKERRQ(KSPSetOperators(ksp, M_p, PM_p)); 100 } else { 101 CHKERRQ(KSPSetOperators(ksp, M_p, M_p)); 102 } 103 CHKERRQ(KSPSolveTranspose(ksp, rhs, f)); 104 CHKERRQ(KSPDestroy(&ksp)); 105 CHKERRQ(VecDestroy(&rhs)); 106 107 /* Visualize particle field */ 108 CHKERRQ(DMSetOutputSequenceNumber(sw, timestep, time)); 109 CHKERRQ(VecViewFromOptions(f, NULL, "-weights_view")); 110 CHKERRQ(VecNorm(f,NORM_1,&norm)); 111 CHKERRQ(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &f)); 112 113 /* compute energy */ 114 CHKERRQ(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void**)&wq)); 115 CHKERRQ(DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void**)&coords)); 116 for (p=0,energy_1=0;p<Np;p++) { 117 energy_1 += wq[p]*(PetscSqr(coords[p*2+0])+PetscSqr(coords[p*2+1])); 118 } 119 CHKERRQ(DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void**)&coords)); 120 CHKERRQ(DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void**)&wq)); 121 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)); 122 /* Cleanup */ 123 CHKERRQ(MatDestroy(&M_p)); 124 CHKERRQ(MatDestroy(&PM_p)); 125 CHKERRQ(VecDestroy(&rho)); 126 CHKERRQ(DMDestroy(&sw)); 127 CHKERRQ(DMDestroy(&dm)); 128 ierr = PetscFinalize(); 129 return ierr; 130 } 131 132 /*TEST 133 134 build: 135 requires: !complex 136 137 test: 138 suffix: 0 139 requires: double triangle 140 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 141 filter: grep -v DM_ | grep -v atomic 142 143 test: 144 suffix: bjacobi 145 requires: double triangle 146 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 147 filter: grep -v DM_ | grep -v atomic 148 149 TEST*/ 150