xref: /petsc/src/dm/tutorials/swarm_ex2.c (revision 362febeeeb69b91ebadcb4b2dc0a22cb6dfc4097)
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   PetscErrorCode ierr;
15c4762a1bSJed Brown   Vec            x;
16c4762a1bSJed Brown   PetscMPIInt    rank;
17c4762a1bSJed Brown   PetscInt       p,bs,nlocal;
18c4762a1bSJed Brown 
19*362febeeSStefano Zampini   PetscFunctionBegin;
20ffc4695bSBarry Smith   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr);
21c4762a1bSJed Brown 
22c4762a1bSJed Brown   ierr = DMCreate(PETSC_COMM_WORLD,&dms);CHKERRQ(ierr);
23c4762a1bSJed Brown   ierr = DMSetType(dms,DMSWARM);CHKERRQ(ierr);
24c4762a1bSJed Brown   ierr = DMSwarmInitializeFieldRegister(dms);CHKERRQ(ierr);
25c4762a1bSJed Brown   ierr = DMSwarmRegisterPetscDatatypeField(dms,"viscosity",1,PETSC_REAL);CHKERRQ(ierr);
26c4762a1bSJed Brown   ierr = DMSwarmRegisterPetscDatatypeField(dms,"strain",3,PETSC_REAL);CHKERRQ(ierr);
27c4762a1bSJed Brown   ierr = DMSwarmFinalizeFieldRegister(dms);CHKERRQ(ierr);
28c4762a1bSJed Brown   ierr = DMSwarmSetLocalSizes(dms,5+rank,4);CHKERRQ(ierr);
29c4762a1bSJed Brown   ierr = DMView(dms,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
30c4762a1bSJed Brown   ierr = DMSwarmGetLocalSize(dms,&nlocal);CHKERRQ(ierr);
31c4762a1bSJed Brown 
32c4762a1bSJed Brown   {
33c4762a1bSJed Brown     PetscReal *array;
34c4762a1bSJed Brown     ierr = DMSwarmGetField(dms,"viscosity",&bs,NULL,(void**)&array);CHKERRQ(ierr);
35c4762a1bSJed Brown     for (p=0; p<nlocal; p++) {
36c4762a1bSJed Brown       array[p] = 11.1 + p*0.1 + rank*100.0;
37c4762a1bSJed Brown     }
38c4762a1bSJed Brown     ierr = DMSwarmRestoreField(dms,"viscosity",&bs,NULL,(void**)&array);CHKERRQ(ierr);
39c4762a1bSJed Brown   }
40c4762a1bSJed Brown 
41c4762a1bSJed Brown   {
42c4762a1bSJed Brown     PetscReal *array;
43c4762a1bSJed Brown     ierr = DMSwarmGetField(dms,"strain",&bs,NULL,(void**)&array);CHKERRQ(ierr);
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     }
49c4762a1bSJed Brown     ierr = DMSwarmRestoreField(dms,"strain",&bs,NULL,(void**)&array);CHKERRQ(ierr);
50c4762a1bSJed Brown   }
51c4762a1bSJed Brown 
52c4762a1bSJed Brown   ierr = DMSwarmCreateGlobalVectorFromField(dms,"viscosity",&x);CHKERRQ(ierr);
53c4762a1bSJed Brown   ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
54c4762a1bSJed Brown   ierr = DMSwarmDestroyGlobalVectorFromField(dms,"viscosity",&x);CHKERRQ(ierr);
55c4762a1bSJed Brown 
56c4762a1bSJed Brown   ierr = DMSwarmCreateGlobalVectorFromField(dms,"strain",&x);CHKERRQ(ierr);
57c4762a1bSJed Brown   ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
58c4762a1bSJed Brown   ierr = DMSwarmDestroyGlobalVectorFromField(dms,"strain",&x);CHKERRQ(ierr);
59c4762a1bSJed Brown 
60c4762a1bSJed Brown   ierr = DMSwarmVectorDefineField(dms,"strain");CHKERRQ(ierr);
61c4762a1bSJed Brown   ierr = DMCreateGlobalVector(dms,&x);CHKERRQ(ierr);
62c4762a1bSJed Brown   ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
63c4762a1bSJed Brown   ierr = VecDestroy(&x);CHKERRQ(ierr);
64c4762a1bSJed Brown   ierr = DMDestroy(&dms);CHKERRQ(ierr);
65c4762a1bSJed Brown   PetscFunctionReturn(0);
66c4762a1bSJed Brown }
67c4762a1bSJed Brown 
68c4762a1bSJed Brown int main(int argc,char **argv)
69c4762a1bSJed Brown {
70c4762a1bSJed Brown   PetscErrorCode ierr;
71c4762a1bSJed Brown 
72c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
73c4762a1bSJed Brown   ierr = ex2_1();CHKERRQ(ierr);
74c4762a1bSJed Brown   ierr = PetscFinalize();
75c4762a1bSJed Brown   return ierr;
76c4762a1bSJed Brown }
77c4762a1bSJed Brown 
78c4762a1bSJed Brown /*TEST
79c4762a1bSJed Brown 
80c4762a1bSJed Brown    test:
81c4762a1bSJed Brown       requires: !complex double
82c4762a1bSJed Brown       nsize: 3
83c4762a1bSJed Brown       filter: grep -v atomic
84c4762a1bSJed Brown       filter_output: grep -v atomic
85c4762a1bSJed Brown 
86c4762a1bSJed Brown TEST*/
87