xref: /petsc/src/dm/impls/swarm/tutorials/ex1.c (revision 3c61180074fb218d306db6efca283cdc905c977f)
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 {
10*3c611800SMark Adams   DM              dm, sw;
114cd81d59SMatthew G. Knepley   PetscFE         fe;
124cd81d59SMatthew G. Knepley   KSP             ksp;
13*3c611800SMark Adams   PC              pc;
14*3c611800SMark Adams   Mat             M_p, PM_p=NULL;
15*3c611800SMark Adams   Vec             f, rho, rhs;
16*3c611800SMark Adams   PetscInt        dim, Nc = 1, timestep = 0, i, faces[3];
1730602db0SMatthew G. Knepley   PetscInt        Np = 10, p, field = 0, zero = 0, bs;
18*3c611800SMark 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;
24*3c611800SMark 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 */
2830602db0SMatthew G. Knepley   ierr = DMCreate(PETSC_COMM_WORLD, &dm);CHKERRQ(ierr);
2930602db0SMatthew G. Knepley   ierr = DMSetType(dm, DMPLEX);CHKERRQ(ierr);
3030602db0SMatthew G. Knepley   ierr = DMSetFromOptions(dm);CHKERRQ(ierr);
314cd81d59SMatthew G. Knepley   ierr = DMViewFromOptions(dm, NULL, "-dm_view");CHKERRQ(ierr);
3230602db0SMatthew G. Knepley 
3330602db0SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
3430602db0SMatthew G. Knepley   i    = dim;
3530602db0SMatthew G. Knepley   ierr = PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &i, NULL);CHKERRQ(ierr);
36*3c611800SMark Adams   ierr = PetscOptionsGetInt(NULL, NULL, "-np", &Np, NULL);CHKERRQ(ierr);
3730602db0SMatthew G. Knepley   ierr = DMGetBoundingBox(dm, lo, hi);CHKERRQ(ierr);
3830602db0SMatthew G. Knepley   for (i=0;i<dim;i++) {
3930602db0SMatthew G. Knepley     h[i] = (hi[i] - lo[i])/faces[i];
4030602db0SMatthew G. Knepley     ierr = PetscPrintf(PETSC_COMM_SELF," lo = %g hi = %g n = %D h = %g\n",lo[i],hi[i],faces[i],h[i]);CHKERRQ(ierr);
4130602db0SMatthew G. Knepley   }
4230602db0SMatthew G. Knepley 
434cd81d59SMatthew G. Knepley   ierr = PetscFECreateDefault(PETSC_COMM_SELF, dim, Nc, PETSC_FALSE, "", PETSC_DECIDE, &fe);CHKERRQ(ierr);
444cd81d59SMatthew G. Knepley   ierr = PetscFESetFromOptions(fe);CHKERRQ(ierr);
454cd81d59SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject)fe, "fe");CHKERRQ(ierr);
464cd81d59SMatthew G. Knepley   ierr = DMSetField(dm, field, NULL, (PetscObject)fe);CHKERRQ(ierr);
474cd81d59SMatthew G. Knepley   ierr = DMCreateDS(dm);CHKERRQ(ierr);
484cd81d59SMatthew G. Knepley   ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
494cd81d59SMatthew G. Knepley   /* Create particle swarm */
504cd81d59SMatthew G. Knepley   ierr = DMCreate(PETSC_COMM_SELF, &sw);CHKERRQ(ierr);
514cd81d59SMatthew G. Knepley   ierr = DMSetType(sw, DMSWARM);CHKERRQ(ierr);
524cd81d59SMatthew G. Knepley   ierr = DMSetDimension(sw, dim);CHKERRQ(ierr);
534cd81d59SMatthew G. Knepley   ierr = DMSwarmSetType(sw, DMSWARM_PIC);CHKERRQ(ierr);
544cd81d59SMatthew G. Knepley   ierr = DMSwarmSetCellDM(sw, dm);CHKERRQ(ierr);
554cd81d59SMatthew G. Knepley   ierr = DMSwarmRegisterPetscDatatypeField(sw, "w_q", Nc, PETSC_SCALAR);CHKERRQ(ierr);
564cd81d59SMatthew G. Knepley   ierr = DMSwarmFinalizeFieldRegister(sw);CHKERRQ(ierr);
574cd81d59SMatthew G. Knepley   ierr = DMSwarmSetLocalSizes(sw, Np, zero);CHKERRQ(ierr);
584cd81d59SMatthew G. Knepley   ierr = DMSetFromOptions(sw);CHKERRQ(ierr);
594cd81d59SMatthew G. Knepley   ierr = DMSwarmGetField(sw, "w_q", &bs, &dtype, (void**)&wq);CHKERRQ(ierr);
604cd81d59SMatthew G. Knepley   ierr = DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void**)&coords);CHKERRQ(ierr);
61*3c611800SMark 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;
65*3c611800SMark Adams     energy_0 += wq[p]*(PetscSqr(coords[p*2+0])+PetscSqr(coords[p*2+1]));
664cd81d59SMatthew G. Knepley   }
674cd81d59SMatthew G. Knepley   ierr = DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void**)&coords);CHKERRQ(ierr);
684cd81d59SMatthew G. Knepley   ierr = DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void**)&wq);CHKERRQ(ierr);
694cd81d59SMatthew G. Knepley   ierr = DMSwarmMigrate(sw, removePoints);CHKERRQ(ierr);
704cd81d59SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject)sw, "Particle Grid");CHKERRQ(ierr);
714cd81d59SMatthew G. Knepley   ierr = DMViewFromOptions(sw, NULL, "-swarm_view");CHKERRQ(ierr);
72*3c611800SMark 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 */
754cd81d59SMatthew G. Knepley   ierr = DMCreateMassMatrix(sw, dm, &M_p);CHKERRQ(ierr);
764cd81d59SMatthew G. Knepley   ierr = DMCreateGlobalVector(dm, &rho);CHKERRQ(ierr);
774cd81d59SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject)rho, "rho");CHKERRQ(ierr);
78*3c611800SMark Adams 
794cd81d59SMatthew G. Knepley   ierr = DMSwarmCreateGlobalVectorFromField(sw, "w_q", &f);CHKERRQ(ierr);
804cd81d59SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject)f, "weights");CHKERRQ(ierr);
814cd81d59SMatthew G. Knepley   ierr = MatMultTranspose(M_p, f, rho);CHKERRQ(ierr);
824cd81d59SMatthew G. Knepley 
834cd81d59SMatthew G. Knepley   /* Visualize mesh field */
844cd81d59SMatthew G. Knepley   ierr = DMSetOutputSequenceNumber(dm, timestep, time);CHKERRQ(ierr);
854cd81d59SMatthew G. Knepley   ierr = VecViewFromOptions(rho, NULL, "-rho_view");CHKERRQ(ierr);
864cd81d59SMatthew G. Knepley 
874cd81d59SMatthew G. Knepley   /* Project field to particles */
884cd81d59SMatthew G. Knepley   /*   This gives f_p = M_p^+ M f */
894cd81d59SMatthew G. Knepley   ierr = DMCreateGlobalVector(dm, &rhs);CHKERRQ(ierr);
9030602db0SMatthew G. Knepley   ierr = VecCopy(rho, rhs);CHKERRQ(ierr); /* Identity: M^1 M rho */
91*3c611800SMark Adams 
924cd81d59SMatthew G. Knepley   ierr = KSPCreate(PETSC_COMM_WORLD, &ksp);CHKERRQ(ierr);
934cd81d59SMatthew G. Knepley   ierr = KSPSetOptionsPrefix(ksp, "ftop_");CHKERRQ(ierr);
944cd81d59SMatthew G. Knepley   ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr);
95*3c611800SMark Adams   ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
96*3c611800SMark Adams   ierr = PetscObjectTypeCompare((PetscObject)pc,PCBJACOBI,&is_bjac);
97*3c611800SMark Adams   if (is_bjac) {
98*3c611800SMark Adams     ierr = DMSwarmCreateMassMatrixSquare(sw, dm, &PM_p);CHKERRQ(ierr);
99*3c611800SMark Adams     ierr = KSPSetOperators(ksp, M_p, PM_p);CHKERRQ(ierr);
100*3c611800SMark Adams   } else {
1014cd81d59SMatthew G. Knepley     ierr = KSPSetOperators(ksp, M_p, M_p);CHKERRQ(ierr);
102*3c611800SMark Adams   }
1034cd81d59SMatthew G. Knepley   ierr = KSPSolveTranspose(ksp, rhs, f);CHKERRQ(ierr);
1044cd81d59SMatthew G. Knepley   ierr = KSPDestroy(&ksp);CHKERRQ(ierr);
1054cd81d59SMatthew G. Knepley   ierr = VecDestroy(&rhs);CHKERRQ(ierr);
1064cd81d59SMatthew G. Knepley 
1074cd81d59SMatthew G. Knepley   /* Visualize particle field */
1084cd81d59SMatthew G. Knepley   ierr = DMSetOutputSequenceNumber(sw, timestep, time);CHKERRQ(ierr);
1094cd81d59SMatthew G. Knepley   ierr = VecViewFromOptions(f, NULL, "-weights_view");CHKERRQ(ierr);
1104cd81d59SMatthew G. Knepley   ierr = VecNorm(f,NORM_1,&norm);CHKERRQ(ierr);
1114cd81d59SMatthew G. Knepley   ierr = DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &f);CHKERRQ(ierr);
112*3c611800SMark Adams 
113*3c611800SMark Adams   /* compute energy */
114*3c611800SMark Adams   ierr = DMSwarmGetField(sw, "w_q", &bs, &dtype, (void**)&wq);CHKERRQ(ierr);
115*3c611800SMark Adams   ierr = DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void**)&coords);CHKERRQ(ierr);
116*3c611800SMark Adams   for (p=0,energy_1=0;p<Np;p++) {
117*3c611800SMark Adams     energy_1 += wq[p]*(PetscSqr(coords[p*2+0])+PetscSqr(coords[p*2+1]));
118*3c611800SMark Adams   }
119*3c611800SMark Adams   ierr = DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void**)&coords);CHKERRQ(ierr);
120*3c611800SMark Adams   ierr = DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void**)&wq);CHKERRQ(ierr);
121*3c611800SMark Adams   ierr = PetscPrintf(PETSC_COMM_SELF,"Total number = %20.12e. energy = %20.12e error = %20.12e\n", norm, energy_0, (energy_1-energy_0)/energy_0);CHKERRQ(ierr);
122*3c611800SMark Adams   /* Cleanup */
1234cd81d59SMatthew G. Knepley   ierr = MatDestroy(&M_p);CHKERRQ(ierr);
124*3c611800SMark Adams   ierr = MatDestroy(&PM_p);CHKERRQ(ierr);
1254cd81d59SMatthew G. Knepley   ierr = VecDestroy(&rho);CHKERRQ(ierr);
1264cd81d59SMatthew G. Knepley   ierr = DMDestroy(&sw);CHKERRQ(ierr);
1274cd81d59SMatthew G. Knepley   ierr = DMDestroy(&dm);CHKERRQ(ierr);
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
139*3c611800SMark Adams     requires: double triangle
140*3c611800SMark 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
141*3c611800SMark Adams     filter: grep -v DM_ | grep -v atomic
142*3c611800SMark Adams 
143*3c611800SMark Adams   test:
144*3c611800SMark Adams     suffix: bjacobi
145*3c611800SMark Adams     requires: double triangle
146*3c611800SMark 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