xref: /petsc/src/dm/tutorials/swarm_ex2.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests DMSwarm\n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscdm.h>
5c4762a1bSJed Brown #include <petscdmda.h>
6c4762a1bSJed Brown #include <petscdmswarm.h>
7c4762a1bSJed Brown 
8c4762a1bSJed Brown /*
9c4762a1bSJed Brown  Checks for variable blocksize
10c4762a1bSJed Brown */
11*d71ae5a4SJacob Faibussowitsch PetscErrorCode ex2_1(void)
12*d71ae5a4SJacob Faibussowitsch {
13c4762a1bSJed Brown   DM          dms;
14c4762a1bSJed Brown   Vec         x;
15c4762a1bSJed Brown   PetscMPIInt rank;
16c4762a1bSJed Brown   PetscInt    p, bs, nlocal;
17c4762a1bSJed Brown 
18362febeeSStefano Zampini   PetscFunctionBegin;
199566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
20c4762a1bSJed Brown 
219566063dSJacob Faibussowitsch   PetscCall(DMCreate(PETSC_COMM_WORLD, &dms));
229566063dSJacob Faibussowitsch   PetscCall(DMSetType(dms, DMSWARM));
236a5217c0SMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)dms, "Particles"));
249566063dSJacob Faibussowitsch   PetscCall(DMSwarmInitializeFieldRegister(dms));
259566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(dms, "viscosity", 1, PETSC_REAL));
269566063dSJacob Faibussowitsch   PetscCall(DMSwarmRegisterPetscDatatypeField(dms, "strain", 3, PETSC_REAL));
279566063dSJacob Faibussowitsch   PetscCall(DMSwarmFinalizeFieldRegister(dms));
289566063dSJacob Faibussowitsch   PetscCall(DMSwarmSetLocalSizes(dms, 5 + rank, 4));
299566063dSJacob Faibussowitsch   PetscCall(DMView(dms, PETSC_VIEWER_STDOUT_WORLD));
309566063dSJacob Faibussowitsch   PetscCall(DMSwarmGetLocalSize(dms, &nlocal));
31c4762a1bSJed Brown 
32c4762a1bSJed Brown   {
33c4762a1bSJed Brown     PetscReal *array;
349566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms, "viscosity", &bs, NULL, (void **)&array));
35ad540459SPierre Jolivet     for (p = 0; p < nlocal; p++) array[p] = 11.1 + p * 0.1 + rank * 100.0;
369566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms, "viscosity", &bs, NULL, (void **)&array));
37c4762a1bSJed Brown   }
38c4762a1bSJed Brown 
39c4762a1bSJed Brown   {
40c4762a1bSJed Brown     PetscReal *array;
419566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms, "strain", &bs, NULL, (void **)&array));
42c4762a1bSJed Brown     for (p = 0; p < nlocal; p++) {
43c4762a1bSJed Brown       array[bs * p + 0] = 2.0e-2 + p * 0.001 + rank * 1.0;
44c4762a1bSJed Brown       array[bs * p + 1] = 2.0e-2 + p * 0.002 + rank * 1.0;
45c4762a1bSJed Brown       array[bs * p + 2] = 2.0e-2 + p * 0.003 + rank * 1.0;
46c4762a1bSJed Brown     }
479566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms, "strain", &bs, NULL, (void **)&array));
48c4762a1bSJed Brown   }
49c4762a1bSJed Brown 
509566063dSJacob Faibussowitsch   PetscCall(DMSwarmCreateGlobalVectorFromField(dms, "viscosity", &x));
519566063dSJacob Faibussowitsch   PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD));
529566063dSJacob Faibussowitsch   PetscCall(DMSwarmDestroyGlobalVectorFromField(dms, "viscosity", &x));
53c4762a1bSJed Brown 
549566063dSJacob Faibussowitsch   PetscCall(DMSwarmCreateGlobalVectorFromField(dms, "strain", &x));
559566063dSJacob Faibussowitsch   PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD));
569566063dSJacob Faibussowitsch   PetscCall(DMSwarmDestroyGlobalVectorFromField(dms, "strain", &x));
57c4762a1bSJed Brown 
589566063dSJacob Faibussowitsch   PetscCall(DMSwarmVectorDefineField(dms, "strain"));
599566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(dms, &x));
609566063dSJacob Faibussowitsch   PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD));
619566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
629566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dms));
63c4762a1bSJed Brown   PetscFunctionReturn(0);
64c4762a1bSJed Brown }
65c4762a1bSJed Brown 
66*d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
67*d71ae5a4SJacob Faibussowitsch {
68327415f7SBarry Smith   PetscFunctionBeginUser;
699566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
709566063dSJacob Faibussowitsch   PetscCall(ex2_1());
719566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
72b122ec5aSJacob Faibussowitsch   return 0;
73c4762a1bSJed Brown }
74c4762a1bSJed Brown 
75c4762a1bSJed Brown /*TEST
76c4762a1bSJed Brown 
77c4762a1bSJed Brown    test:
78c4762a1bSJed Brown       requires: !complex double
79c4762a1bSJed Brown       nsize: 3
80c4762a1bSJed Brown       filter: grep -v atomic
81c4762a1bSJed Brown       filter_output: grep -v atomic
82c4762a1bSJed Brown 
83c4762a1bSJed Brown TEST*/
84