xref: /petsc/src/dm/tutorials/swarm_ex2.c (revision 6a5217c03994f2d95bb2e6dbd8bed42381aeb015)
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 */
11c4762a1bSJed Brown PetscErrorCode ex2_1(void)
12c4762a1bSJed Brown {
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));
23*6a5217c0SMatthew 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));
35c4762a1bSJed Brown     for (p=0; p<nlocal; p++) {
36c4762a1bSJed Brown       array[p] = 11.1 + p*0.1 + rank*100.0;
37c4762a1bSJed Brown     }
389566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms,"viscosity",&bs,NULL,(void**)&array));
39c4762a1bSJed Brown   }
40c4762a1bSJed Brown 
41c4762a1bSJed Brown   {
42c4762a1bSJed Brown     PetscReal *array;
439566063dSJacob Faibussowitsch     PetscCall(DMSwarmGetField(dms,"strain",&bs,NULL,(void**)&array));
44c4762a1bSJed Brown     for (p=0; p<nlocal; p++) {
45c4762a1bSJed Brown       array[bs*p+0] = 2.0e-2 + p*0.001 + rank*1.0;
46c4762a1bSJed Brown       array[bs*p+1] = 2.0e-2 + p*0.002 + rank*1.0;
47c4762a1bSJed Brown       array[bs*p+2] = 2.0e-2 + p*0.003 + rank*1.0;
48c4762a1bSJed Brown     }
499566063dSJacob Faibussowitsch     PetscCall(DMSwarmRestoreField(dms,"strain",&bs,NULL,(void**)&array));
50c4762a1bSJed Brown   }
51c4762a1bSJed Brown 
529566063dSJacob Faibussowitsch   PetscCall(DMSwarmCreateGlobalVectorFromField(dms,"viscosity",&x));
539566063dSJacob Faibussowitsch   PetscCall(VecView(x,PETSC_VIEWER_STDOUT_WORLD));
549566063dSJacob Faibussowitsch   PetscCall(DMSwarmDestroyGlobalVectorFromField(dms,"viscosity",&x));
55c4762a1bSJed Brown 
569566063dSJacob Faibussowitsch   PetscCall(DMSwarmCreateGlobalVectorFromField(dms,"strain",&x));
579566063dSJacob Faibussowitsch   PetscCall(VecView(x,PETSC_VIEWER_STDOUT_WORLD));
589566063dSJacob Faibussowitsch   PetscCall(DMSwarmDestroyGlobalVectorFromField(dms,"strain",&x));
59c4762a1bSJed Brown 
609566063dSJacob Faibussowitsch   PetscCall(DMSwarmVectorDefineField(dms,"strain"));
619566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(dms,&x));
629566063dSJacob Faibussowitsch   PetscCall(VecView(x,PETSC_VIEWER_STDOUT_WORLD));
639566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
649566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dms));
65c4762a1bSJed Brown   PetscFunctionReturn(0);
66c4762a1bSJed Brown }
67c4762a1bSJed Brown 
68c4762a1bSJed Brown int main(int argc,char **argv)
69c4762a1bSJed Brown {
70c4762a1bSJed Brown 
719566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
729566063dSJacob Faibussowitsch   PetscCall(ex2_1());
739566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
74b122ec5aSJacob Faibussowitsch   return 0;
75c4762a1bSJed Brown }
76c4762a1bSJed Brown 
77c4762a1bSJed Brown /*TEST
78c4762a1bSJed Brown 
79c4762a1bSJed Brown    test:
80c4762a1bSJed Brown       requires: !complex double
81c4762a1bSJed Brown       nsize: 3
82c4762a1bSJed Brown       filter: grep -v atomic
83c4762a1bSJed Brown       filter_output: grep -v atomic
84c4762a1bSJed Brown 
85c4762a1bSJed Brown TEST*/
86