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; 233c611800SMark Adams PetscBool is_bjac; 244cd81d59SMatthew G. Knepley 259566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL,help)); 264cd81d59SMatthew G. Knepley /* Create a mesh */ 279566063dSJacob Faibussowitsch PetscCall(DMCreate(PETSC_COMM_WORLD, &dm)); 289566063dSJacob Faibussowitsch PetscCall(DMSetType(dm, DMPLEX)); 299566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dm)); 309566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 3130602db0SMatthew G. Knepley 329566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 3330602db0SMatthew G. Knepley i = dim; 349566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &i, NULL)); 359566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-np", &Np, NULL)); 369566063dSJacob Faibussowitsch PetscCall(DMGetBoundingBox(dm, lo, hi)); 3730602db0SMatthew G. Knepley for (i=0;i<dim;i++) { 3830602db0SMatthew G. Knepley h[i] = (hi[i] - lo[i])/faces[i]; 39*63a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF," lo = %g hi = %g n = %" PetscInt_FMT " h = %g\n",(double)lo[i],(double)hi[i],faces[i],(double)h[i])); 4030602db0SMatthew G. Knepley } 4130602db0SMatthew G. Knepley 429566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, Nc, PETSC_FALSE, "", PETSC_DECIDE, &fe)); 439566063dSJacob Faibussowitsch PetscCall(PetscFESetFromOptions(fe)); 449566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe, "fe")); 459566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, field, NULL, (PetscObject)fe)); 469566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 479566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 484cd81d59SMatthew G. Knepley /* Create particle swarm */ 499566063dSJacob Faibussowitsch PetscCall(DMCreate(PETSC_COMM_SELF, &sw)); 509566063dSJacob Faibussowitsch PetscCall(DMSetType(sw, DMSWARM)); 519566063dSJacob Faibussowitsch PetscCall(DMSetDimension(sw, dim)); 529566063dSJacob Faibussowitsch PetscCall(DMSwarmSetType(sw, DMSWARM_PIC)); 539566063dSJacob Faibussowitsch PetscCall(DMSwarmSetCellDM(sw, dm)); 549566063dSJacob Faibussowitsch PetscCall(DMSwarmRegisterPetscDatatypeField(sw, "w_q", Nc, PETSC_SCALAR)); 559566063dSJacob Faibussowitsch PetscCall(DMSwarmFinalizeFieldRegister(sw)); 569566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(sw, Np, zero)); 579566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(sw)); 589566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void**)&wq)); 599566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void**)&coords)); 603c611800SMark Adams for (p=0,energy_0=0;p<Np;p++) { 614cd81d59SMatthew G. Knepley coords[p*2+0] = -PetscCosReal((PetscReal)(p+1)/(PetscReal)(Np+1) * PETSC_PI); 624cd81d59SMatthew G. Knepley coords[p*2+1] = PetscSinReal((PetscReal)(p+1)/(PetscReal)(Np+1) * PETSC_PI); 634cd81d59SMatthew G. Knepley wq[p] = 1.0; 643c611800SMark Adams energy_0 += wq[p]*(PetscSqr(coords[p*2+0])+PetscSqr(coords[p*2+1])); 654cd81d59SMatthew G. Knepley } 669566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void**)&coords)); 679566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void**)&wq)); 689566063dSJacob Faibussowitsch PetscCall(DMSwarmMigrate(sw, removePoints)); 699566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)sw, "Particle Grid")); 709566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(sw, NULL, "-swarm_view")); 713c611800SMark Adams 724cd81d59SMatthew G. Knepley /* Project particles to field */ 734cd81d59SMatthew G. Knepley /* This gives M f = \int_\Omega \phi f, which looks like a rhs for a PDE */ 749566063dSJacob Faibussowitsch PetscCall(DMCreateMassMatrix(sw, dm, &M_p)); 759566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &rho)); 769566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)rho, "rho")); 773c611800SMark Adams 789566063dSJacob Faibussowitsch PetscCall(DMSwarmCreateGlobalVectorFromField(sw, "w_q", &f)); 799566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)f, "weights")); 809566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(M_p, f, rho)); 814cd81d59SMatthew G. Knepley 824cd81d59SMatthew G. Knepley /* Visualize mesh field */ 839566063dSJacob Faibussowitsch PetscCall(DMSetOutputSequenceNumber(dm, timestep, time)); 849566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(rho, NULL, "-rho_view")); 854cd81d59SMatthew G. Knepley 864cd81d59SMatthew G. Knepley /* Project field to particles */ 874cd81d59SMatthew G. Knepley /* This gives f_p = M_p^+ M f */ 889566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &rhs)); 899566063dSJacob Faibussowitsch PetscCall(VecCopy(rho, rhs)); /* Identity: M^1 M rho */ 903c611800SMark Adams 919566063dSJacob Faibussowitsch PetscCall(KSPCreate(PETSC_COMM_WORLD, &ksp)); 929566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(ksp, "ftop_")); 939566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(ksp)); 949566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp,&pc)); 959566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc,PCBJACOBI,&is_bjac)); 963c611800SMark Adams if (is_bjac) { 979566063dSJacob Faibussowitsch PetscCall(DMSwarmCreateMassMatrixSquare(sw, dm, &PM_p)); 989566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(ksp, M_p, PM_p)); 993c611800SMark Adams } else { 1009566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(ksp, M_p, M_p)); 1013c611800SMark Adams } 1029566063dSJacob Faibussowitsch PetscCall(KSPSolveTranspose(ksp, rhs, f)); 1039566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&ksp)); 1049566063dSJacob Faibussowitsch PetscCall(VecDestroy(&rhs)); 1054cd81d59SMatthew G. Knepley 1064cd81d59SMatthew G. Knepley /* Visualize particle field */ 1079566063dSJacob Faibussowitsch PetscCall(DMSetOutputSequenceNumber(sw, timestep, time)); 1089566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(f, NULL, "-weights_view")); 1099566063dSJacob Faibussowitsch PetscCall(VecNorm(f,NORM_1,&norm)); 1109566063dSJacob Faibussowitsch PetscCall(DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &f)); 1113c611800SMark Adams 1123c611800SMark Adams /* compute energy */ 1139566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(sw, "w_q", &bs, &dtype, (void**)&wq)); 1149566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void**)&coords)); 1153c611800SMark Adams for (p=0,energy_1=0;p<Np;p++) { 1163c611800SMark Adams energy_1 += wq[p]*(PetscSqr(coords[p*2+0])+PetscSqr(coords[p*2+1])); 1173c611800SMark Adams } 1189566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void**)&coords)); 1199566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void**)&wq)); 120*63a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Total number = %20.12e. energy = %20.12e error = %20.12e\n", (double)norm, (double)energy_0, (double)((energy_1-energy_0)/energy_0))); 1213c611800SMark Adams /* Cleanup */ 1229566063dSJacob Faibussowitsch PetscCall(MatDestroy(&M_p)); 1239566063dSJacob Faibussowitsch PetscCall(MatDestroy(&PM_p)); 1249566063dSJacob Faibussowitsch PetscCall(VecDestroy(&rho)); 1259566063dSJacob Faibussowitsch PetscCall(DMDestroy(&sw)); 1269566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 1279566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 128b122ec5aSJacob Faibussowitsch return 0; 1294cd81d59SMatthew G. Knepley } 1304cd81d59SMatthew G. Knepley 1314cd81d59SMatthew G. Knepley /*TEST 1324cd81d59SMatthew G. Knepley 1334cd81d59SMatthew G. Knepley build: 1344cd81d59SMatthew G. Knepley requires: !complex 1354cd81d59SMatthew G. Knepley 1364cd81d59SMatthew G. Knepley test: 1374cd81d59SMatthew G. Knepley suffix: 0 1383c611800SMark Adams requires: double triangle 1393c611800SMark 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 1403c611800SMark Adams filter: grep -v DM_ | grep -v atomic 1413c611800SMark Adams 1423c611800SMark Adams test: 1433c611800SMark Adams suffix: bjacobi 1443c611800SMark Adams requires: double triangle 1453c611800SMark 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 1464cd81d59SMatthew G. Knepley filter: grep -v DM_ | grep -v atomic 1474cd81d59SMatthew G. Knepley 1484cd81d59SMatthew G. Knepley TEST*/ 149