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