xref: /petsc/src/dm/tutorials/swarm_ex2.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1 
2 static char help[] = "Tests DMSwarm\n\n";
3 
4 #include <petscdm.h>
5 #include <petscdmda.h>
6 #include <petscdmswarm.h>
7 
8 /*
9  Checks for variable blocksize
10 */
11 PetscErrorCode ex2_1(void)
12 {
13   DM             dms;
14   PetscErrorCode ierr;
15   Vec            x;
16   PetscMPIInt    rank;
17   PetscInt       p,bs,nlocal;
18 
19   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
20 
21   ierr = DMCreate(PETSC_COMM_WORLD,&dms);CHKERRQ(ierr);
22   ierr = DMSetType(dms,DMSWARM);CHKERRQ(ierr);
23   ierr = DMSwarmInitializeFieldRegister(dms);CHKERRQ(ierr);
24   ierr = DMSwarmRegisterPetscDatatypeField(dms,"viscosity",1,PETSC_REAL);CHKERRQ(ierr);
25   ierr = DMSwarmRegisterPetscDatatypeField(dms,"strain",3,PETSC_REAL);CHKERRQ(ierr);
26   ierr = DMSwarmFinalizeFieldRegister(dms);CHKERRQ(ierr);
27   ierr = DMSwarmSetLocalSizes(dms,5+rank,4);CHKERRQ(ierr);
28   ierr = DMView(dms,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
29   ierr = DMSwarmGetLocalSize(dms,&nlocal);CHKERRQ(ierr);
30 
31   {
32     PetscReal *array;
33     ierr = DMSwarmGetField(dms,"viscosity",&bs,NULL,(void**)&array);CHKERRQ(ierr);
34     for (p=0; p<nlocal; p++) {
35       array[p] = 11.1 + p*0.1 + rank*100.0;
36     }
37     ierr = DMSwarmRestoreField(dms,"viscosity",&bs,NULL,(void**)&array);CHKERRQ(ierr);
38   }
39 
40   {
41     PetscReal *array;
42     ierr = DMSwarmGetField(dms,"strain",&bs,NULL,(void**)&array);CHKERRQ(ierr);
43     for (p=0; p<nlocal; p++) {
44       array[bs*p+0] = 2.0e-2 + p*0.001 + rank*1.0;
45       array[bs*p+1] = 2.0e-2 + p*0.002 + rank*1.0;
46       array[bs*p+2] = 2.0e-2 + p*0.003 + rank*1.0;
47     }
48     ierr = DMSwarmRestoreField(dms,"strain",&bs,NULL,(void**)&array);CHKERRQ(ierr);
49   }
50 
51   ierr = DMSwarmCreateGlobalVectorFromField(dms,"viscosity",&x);CHKERRQ(ierr);
52   ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
53   ierr = DMSwarmDestroyGlobalVectorFromField(dms,"viscosity",&x);CHKERRQ(ierr);
54 
55   ierr = DMSwarmCreateGlobalVectorFromField(dms,"strain",&x);CHKERRQ(ierr);
56   ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
57   ierr = DMSwarmDestroyGlobalVectorFromField(dms,"strain",&x);CHKERRQ(ierr);
58 
59   ierr = DMSwarmVectorDefineField(dms,"strain");CHKERRQ(ierr);
60   ierr = DMCreateGlobalVector(dms,&x);CHKERRQ(ierr);
61   ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
62   ierr = VecDestroy(&x);CHKERRQ(ierr);
63   ierr = DMDestroy(&dms);CHKERRQ(ierr);
64   PetscFunctionReturn(0);
65 }
66 
67 
68 int main(int argc,char **argv)
69 {
70   PetscErrorCode ierr;
71 
72   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
73   ierr = ex2_1();CHKERRQ(ierr);
74   ierr = PetscFinalize();
75   return ierr;
76 }
77 
78 /*TEST
79 
80    test:
81       requires: !complex double
82       nsize: 3
83       filter: grep -v atomic
84       filter_output: grep -v atomic
85 
86 TEST*/
87