xref: /petsc/src/dm/impls/swarm/swarm.c (revision f0cdbbba46eb2dee2dd444b383259dc1f084f526)
157795646SDave May 
257795646SDave May #define PETSCDM_DLL
357795646SDave May #include <petsc/private/dmswarmimpl.h>    /*I   "petscdmswarm.h"   I*/
4b62e03f8SDave May #include "data_bucket.h"
557795646SDave May 
6*f0cdbbbaSDave May const char* DMSwarmTypeNames[] = { "basic", "pic", 0 };
7*f0cdbbbaSDave May const char* DMSwarmMigrateTypeNames[] = { "basic", "dmcellnscatter", "dmcellexact", "user", 0 };
8*f0cdbbbaSDave May const char* DMSwarmCollectTypeNames[] = { "basic", "boundingbox", "general", "user", 0 };
9*f0cdbbbaSDave May 
10*f0cdbbbaSDave May const char DMSwarmField_pid[] = "DMSwarm_pid";
11*f0cdbbbaSDave May const char DMSwarmField_rank[] = "DMSwarm_rank";
12*f0cdbbbaSDave May const char DMSwarmPICField_coor[] = "DMSwarmPIC_coor";
13*f0cdbbbaSDave May 
14*f0cdbbbaSDave May 
15dcf43ee8SDave May PetscErrorCode DMSwarmMigrate_Push_Basic(DM dm,PetscBool remove_sent_points);
16dcf43ee8SDave May 
17dcf43ee8SDave May 
1857795646SDave May //typedef PetscErrorCode (*swarm_project)(DM,DM,Vec) DMSwarmProjectMethod; /* swarm, geometry, result */
1957795646SDave May 
2057795646SDave May //typedef enum { PROJECT_DMDA_AQ1=0, PROJECT_DMDA_P0 } DMSwarmDMDAProjectionType;
2157795646SDave May 
2257795646SDave May #if 0
2357795646SDave May 
2457795646SDave May /* Defines what the local space will be */
25b62e03f8SDave May PetscErrorCode DMSwarmSetOverlap(void)
2657795646SDave May {
2757795646SDave May 
28b62e03f8SDave May   PetscFunctionReturn(0);
2957795646SDave May }
3057795646SDave May 
3157795646SDave May 
3257795646SDave May /* coordinates */
3357795646SDave May /*
3457795646SDave May DMGetCoordinateDM returns self
3557795646SDave May DMGetCoordinates and DMGetCoordinatesLocal return same thing
3657795646SDave May Local view could be used to define overlapping information
3757795646SDave May */
3857795646SDave May 
3957795646SDave May #endif
4057795646SDave May 
4157795646SDave May #undef __FUNCT__
42b5bcf523SDave May #define __FUNCT__ "DMSwarmVectorDefineField"
43b5bcf523SDave May PETSC_EXTERN PetscErrorCode DMSwarmVectorDefineField(DM dm,const char fieldname[])
44b5bcf523SDave May {
45b5bcf523SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
46b5bcf523SDave May   PetscErrorCode ierr;
47b5bcf523SDave May   PetscInt bs,n;
48b5bcf523SDave May   PetscScalar *array;
49b5bcf523SDave May   PetscDataType type;
50b5bcf523SDave May 
513454631fSDave May   if (!swarm->issetup) { ierr = DMSetUp(dm);CHKERRQ(ierr); }
526845f8f5SDave May   ierr = DataBucketGetSizes(swarm->db,&n,NULL,NULL);CHKERRQ(ierr);
53b5bcf523SDave May   ierr = DMSwarmGetField(dm,fieldname,&bs,&type,(void**)&array);CHKERRQ(ierr);
54b5bcf523SDave May 
55b5bcf523SDave May   /* Check all fields are of type PETSC_REAL or PETSC_SCALAR */
56b5bcf523SDave May   if (type != PETSC_REAL) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only valid for PETSC_REAL");
57b5bcf523SDave May 
58b5bcf523SDave May   PetscSNPrintf(swarm->vec_field_name,PETSC_MAX_PATH_LEN-1,"%s",fieldname);
59b5bcf523SDave May   swarm->vec_field_set = PETSC_TRUE;
601b1ea282SDave May   swarm->vec_field_bs = bs;
61b5bcf523SDave May   swarm->vec_field_nlocal = n;
62dcf43ee8SDave May   ierr = DMSwarmRestoreField(dm,fieldname,&bs,&type,(void**)&array);CHKERRQ(ierr);
63b5bcf523SDave May 
64b5bcf523SDave May   PetscFunctionReturn(0);
65b5bcf523SDave May }
66b5bcf523SDave May 
67b5bcf523SDave May /* requires DMSwarmDefineFieldVector has been called */
68b5bcf523SDave May #undef __FUNCT__
69b5bcf523SDave May #define __FUNCT__ "DMCreateGlobalVector_Swarm"
70b5bcf523SDave May PetscErrorCode DMCreateGlobalVector_Swarm(DM dm,Vec *vec)
71b5bcf523SDave May {
72b5bcf523SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
73b5bcf523SDave May   PetscErrorCode ierr;
74b5bcf523SDave May   Vec x;
75b5bcf523SDave May   char name[PETSC_MAX_PATH_LEN];
76b5bcf523SDave May 
773454631fSDave May   if (!swarm->issetup) { ierr = DMSetUp(dm);CHKERRQ(ierr); }
78b5bcf523SDave May   if (!swarm->vec_field_set) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMSwarmVectorDefineField first");
79b5bcf523SDave May   PetscSNPrintf(name,PETSC_MAX_PATH_LEN-1,"DMSwarmField_%s",swarm->vec_field_name);
80b5bcf523SDave May   ierr = VecCreate(PetscObjectComm((PetscObject)dm),&x);CHKERRQ(ierr);
81b5bcf523SDave May   ierr = PetscObjectSetName((PetscObject)x,name);CHKERRQ(ierr);
821b1ea282SDave May   ierr = VecSetSizes(x,swarm->db->L*swarm->vec_field_bs,PETSC_DETERMINE);CHKERRQ(ierr);
83b5bcf523SDave May   ierr = VecSetBlockSize(x,swarm->vec_field_bs);CHKERRQ(ierr);
84b5bcf523SDave May   ierr = VecSetFromOptions(x);CHKERRQ(ierr);
85b5bcf523SDave May   *vec = x;
86b5bcf523SDave May 
87b5bcf523SDave May   PetscFunctionReturn(0);
88b5bcf523SDave May }
89b5bcf523SDave May 
90b5bcf523SDave May /* requires DMSwarmDefineFieldVector has been called */
91b5bcf523SDave May #undef __FUNCT__
92b5bcf523SDave May #define __FUNCT__ "DMCreateLocalVector_Swarm"
93b5bcf523SDave May PetscErrorCode DMCreateLocalVector_Swarm(DM dm,Vec *vec)
94b5bcf523SDave May {
95b5bcf523SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
96b5bcf523SDave May   PetscErrorCode ierr;
97b5bcf523SDave May   Vec x;
98b5bcf523SDave May   char name[PETSC_MAX_PATH_LEN];
99b5bcf523SDave May 
1003454631fSDave May   if (!swarm->issetup) { ierr = DMSetUp(dm);CHKERRQ(ierr); }
101b5bcf523SDave May   if (!swarm->vec_field_set) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMSwarmVectorDefineField first");
102b5bcf523SDave May   PetscSNPrintf(name,PETSC_MAX_PATH_LEN-1,"DMSwarmField_%s",swarm->vec_field_name);
103b5bcf523SDave May   ierr = VecCreate(PETSC_COMM_SELF,&x);CHKERRQ(ierr);
104b5bcf523SDave May   ierr = PetscObjectSetName((PetscObject)x,name);CHKERRQ(ierr);
1051b1ea282SDave May   ierr = VecSetSizes(x,swarm->db->L*swarm->vec_field_bs,swarm->db->L);CHKERRQ(ierr);
106b5bcf523SDave May   ierr = VecSetBlockSize(x,swarm->vec_field_bs);CHKERRQ(ierr);
107b5bcf523SDave May   ierr = VecSetFromOptions(x);CHKERRQ(ierr);
108b5bcf523SDave May   *vec = x;
109b5bcf523SDave May 
110b5bcf523SDave May   PetscFunctionReturn(0);
111b5bcf523SDave May }
112b5bcf523SDave May 
113b5bcf523SDave May #undef __FUNCT__
114b5bcf523SDave May #define __FUNCT__ "DMSwarmCreateGlobalVectorFromField"
115b5bcf523SDave May PETSC_EXTERN PetscErrorCode DMSwarmCreateGlobalVectorFromField(DM dm,const char fieldname[],Vec *vec)
116b5bcf523SDave May {
117b5bcf523SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
118b5bcf523SDave May   PetscErrorCode ierr;
119b5bcf523SDave May   PetscInt bs,n;
120b5bcf523SDave May   PetscScalar *array;
121b5bcf523SDave May   Vec x;
122b5bcf523SDave May   PetscDataType type;
123b5bcf523SDave May   char name[PETSC_MAX_PATH_LEN];
1243454631fSDave May   PetscMPIInt commsize;
125b5bcf523SDave May 
1263454631fSDave May   if (!swarm->issetup) { ierr = DMSetUp(dm);CHKERRQ(ierr); }
1276845f8f5SDave May   ierr = DataBucketGetSizes(swarm->db,&n,NULL,NULL);CHKERRQ(ierr);
128b5bcf523SDave May   ierr = DMSwarmGetField(dm,fieldname,&bs,&type,(void**)&array);CHKERRQ(ierr);
129b5bcf523SDave May 
130b5bcf523SDave May   /* Check all fields are of type PETSC_REAL or PETSC_SCALAR */
131b5bcf523SDave May   if (type != PETSC_REAL) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only valid for PETSC_REAL");
132b5bcf523SDave May 
1333454631fSDave May   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm),&commsize);CHKERRQ(ierr);
1343454631fSDave May   if (commsize == 1) {
1351b1ea282SDave May     ierr = VecCreateSeqWithArray(PetscObjectComm((PetscObject)dm),bs,n*bs,array,&x);CHKERRQ(ierr);
1363454631fSDave May   } else {
1371b1ea282SDave May     ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)dm),bs,n*bs,PETSC_DETERMINE,array,&x);CHKERRQ(ierr);
1383454631fSDave May   }
1396864fe11SDave May   PetscSNPrintf(name,PETSC_MAX_PATH_LEN-1,"DMSwarmSharedField_%s",fieldname);
140dcf43ee8SDave May   ierr = PetscObjectSetName((PetscObject)x,name);CHKERRQ(ierr);
141b5bcf523SDave May 
142b5bcf523SDave May   /* Set guard */
143dcf43ee8SDave May   PetscSNPrintf(name,PETSC_MAX_PATH_LEN-1,"DMSwarm_VecFieldInPlace_%s",fieldname);
144dcf43ee8SDave May   ierr = PetscObjectComposeFunction((PetscObject)x,name,DMSwarmDestroyGlobalVectorFromField);CHKERRQ(ierr);
145b5bcf523SDave May 
146b5bcf523SDave May   *vec = x;
147b5bcf523SDave May   PetscFunctionReturn(0);
148b5bcf523SDave May }
149b5bcf523SDave May 
150b5bcf523SDave May #undef __FUNCT__
151b5bcf523SDave May #define __FUNCT__ "DMSwarmDestroyGlobalVectorFromField"
152b5bcf523SDave May PETSC_EXTERN PetscErrorCode DMSwarmDestroyGlobalVectorFromField(DM dm,const char fieldname[],Vec *vec)
153b5bcf523SDave May {
154b5bcf523SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
155b5bcf523SDave May   PetscErrorCode ierr;
156b5bcf523SDave May   DataField gfield;
157b5bcf523SDave May   char name[PETSC_MAX_PATH_LEN];
158b5bcf523SDave May   void (*fptr)(void);
159b5bcf523SDave May 
160b5bcf523SDave May   /* get data field */
1612eac95f8SDave May   ierr = DataBucketGetDataFieldByName(swarm->db,fieldname,&gfield);CHKERRQ(ierr);
162b5bcf523SDave May 
163b5bcf523SDave May   /* check vector is an inplace array */
164b5bcf523SDave May   PetscSNPrintf(name,PETSC_MAX_PATH_LEN-1,"DMSwarm_VecFieldInPlace_%s",fieldname);
165b5bcf523SDave May   ierr = PetscObjectQueryFunction((PetscObject)(*vec),name,&fptr);CHKERRQ(ierr);
166b5bcf523SDave May   if (!fptr) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Vector being destroyed was not created from DMSwarm field(%s)",fieldname);
167b5bcf523SDave May 
168b5bcf523SDave May   /* restore data field */
1696845f8f5SDave May   ierr = DataFieldRestoreAccess(gfield);CHKERRQ(ierr);
170b5bcf523SDave May 
171b5bcf523SDave May   ierr = VecDestroy(vec);CHKERRQ(ierr);
172b5bcf523SDave May 
173b5bcf523SDave May   PetscFunctionReturn(0);
174b5bcf523SDave May }
175b5bcf523SDave May 
176b5bcf523SDave May /*
177b5bcf523SDave May PETSC_EXTERN PetscErrorCode DMSwarmCreateGlobalVectorFromFields(DM dm,const PetscInt nf,const char *fieldnames[],Vec *vec)
178b5bcf523SDave May {
179b5bcf523SDave May   PetscFunctionReturn(0);
180b5bcf523SDave May }
181b5bcf523SDave May 
182b5bcf523SDave May PETSC_EXTERN PetscErrorCode DMSwarmRestoreGlobalVectorFromFields(DM dm,Vec *vec)
183b5bcf523SDave May {
184b5bcf523SDave May   PetscFunctionReturn(0);
185b5bcf523SDave May }
186b5bcf523SDave May */
187b5bcf523SDave May 
188b5bcf523SDave May #undef __FUNCT__
1895f50eb2eSDave May #define __FUNCT__ "DMSwarmInitializeFieldRegister"
1905f50eb2eSDave May PETSC_EXTERN PetscErrorCode DMSwarmInitializeFieldRegister(DM dm)
1915f50eb2eSDave May {
1925f50eb2eSDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
1933454631fSDave May   PetscErrorCode ierr;
1943454631fSDave May 
1955f50eb2eSDave May   swarm->field_registration_initialized = PETSC_TRUE;
1963454631fSDave May 
197*f0cdbbbaSDave May   ierr = DMSwarmRegisterPetscDatatypeField(dm,DMSwarmField_pid,1,PETSC_LONG);CHKERRQ(ierr); /* unique identifer */
198*f0cdbbbaSDave May   ierr = DMSwarmRegisterPetscDatatypeField(dm,DMSwarmField_rank,1,PETSC_INT);CHKERRQ(ierr); /* used for communication */
1993454631fSDave May 
2005f50eb2eSDave May   PetscFunctionReturn(0);
2015f50eb2eSDave May }
2025f50eb2eSDave May 
2035f50eb2eSDave May #undef __FUNCT__
2045f50eb2eSDave May #define __FUNCT__ "DMSwarmFinalizeFieldRegister"
2055f50eb2eSDave May PETSC_EXTERN PetscErrorCode DMSwarmFinalizeFieldRegister(DM dm)
2065f50eb2eSDave May {
2075f50eb2eSDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
2086845f8f5SDave May   PetscErrorCode ierr;
2096845f8f5SDave May 
210*f0cdbbbaSDave May   if (!swarm->field_registration_finalized) {
2116845f8f5SDave May     ierr = DataBucketFinalize(swarm->db);CHKERRQ(ierr);
212*f0cdbbbaSDave May   }
213*f0cdbbbaSDave May   swarm->field_registration_finalized = PETSC_TRUE;
2145f50eb2eSDave May   PetscFunctionReturn(0);
2155f50eb2eSDave May }
2165f50eb2eSDave May 
2175f50eb2eSDave May #undef __FUNCT__
2185f50eb2eSDave May #define __FUNCT__ "DMSwarmSetLocalSizes"
2195f50eb2eSDave May PETSC_EXTERN PetscErrorCode DMSwarmSetLocalSizes(DM dm,PetscInt nlocal,PetscInt buffer)
2205f50eb2eSDave May {
2215f50eb2eSDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
2226845f8f5SDave May   PetscErrorCode ierr;
2235f50eb2eSDave May 
2246845f8f5SDave May   ierr = DataBucketSetSizes(swarm->db,nlocal,buffer);CHKERRQ(ierr);
2255f50eb2eSDave May 
2265f50eb2eSDave May   PetscFunctionReturn(0);
2275f50eb2eSDave May }
2285f50eb2eSDave May 
2295f50eb2eSDave May #undef __FUNCT__
230b16650c8SDave May #define __FUNCT__ "DMSwarmSetCellDM"
231b16650c8SDave May PETSC_EXTERN PetscErrorCode DMSwarmSetCellDM(DM dm,DM dmcell)
232b16650c8SDave May {
233b16650c8SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
234b16650c8SDave May   swarm->dmcell = dmcell;
235b16650c8SDave May   PetscFunctionReturn(0);
236b16650c8SDave May }
237b16650c8SDave May 
238b16650c8SDave May #undef __FUNCT__
239fe39f135SDave May #define __FUNCT__ "DMSwarmGetCellDM"
240fe39f135SDave May PETSC_EXTERN PetscErrorCode DMSwarmGetCellDM(DM dm,DM *dmcell)
241fe39f135SDave May {
242fe39f135SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
243fe39f135SDave May   *dmcell = swarm->dmcell;
244fe39f135SDave May   PetscFunctionReturn(0);
245fe39f135SDave May }
246fe39f135SDave May 
247fe39f135SDave May #undef __FUNCT__
248dcf43ee8SDave May #define __FUNCT__ "DMSwarmGetLocalSize"
249dcf43ee8SDave May PETSC_EXTERN PetscErrorCode DMSwarmGetLocalSize(DM dm,PetscInt *nlocal)
250dcf43ee8SDave May {
251dcf43ee8SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
252dcf43ee8SDave May   PetscErrorCode ierr;
253dcf43ee8SDave May 
254dcf43ee8SDave May   if (nlocal) {
255dcf43ee8SDave May     ierr = DataBucketGetSizes(swarm->db,nlocal,NULL,NULL);CHKERRQ(ierr);
256dcf43ee8SDave May   }
257dcf43ee8SDave May 
258dcf43ee8SDave May   PetscFunctionReturn(0);
259dcf43ee8SDave May }
260dcf43ee8SDave May 
261dcf43ee8SDave May #undef __FUNCT__
262dcf43ee8SDave May #define __FUNCT__ "DMSwarmGetSize"
263dcf43ee8SDave May PETSC_EXTERN PetscErrorCode DMSwarmGetSize(DM dm,PetscInt *n)
264dcf43ee8SDave May {
265dcf43ee8SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
266dcf43ee8SDave May   PetscErrorCode ierr;
267dcf43ee8SDave May   PetscInt nlocal,ng;
268dcf43ee8SDave May 
269dcf43ee8SDave May   ierr = DataBucketGetSizes(swarm->db,&nlocal,NULL,NULL);CHKERRQ(ierr);
270dcf43ee8SDave May   ierr = MPI_Allreduce(&nlocal,&ng,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
271dcf43ee8SDave May   if (n) { *n = ng; }
272dcf43ee8SDave May   PetscFunctionReturn(0);
273dcf43ee8SDave May }
274dcf43ee8SDave May 
275dcf43ee8SDave May #undef __FUNCT__
276b62e03f8SDave May #define __FUNCT__ "DMSwarmRegisterPetscDatatypeField"
2775f50eb2eSDave May PETSC_EXTERN PetscErrorCode DMSwarmRegisterPetscDatatypeField(DM dm,const char fieldname[],PetscInt blocksize,PetscDataType type)
278b62e03f8SDave May {
2792eac95f8SDave May   PetscErrorCode ierr;
280b62e03f8SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
281b62e03f8SDave May   size_t size;
282b62e03f8SDave May 
2835f50eb2eSDave May   if (!swarm->field_registration_initialized) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMSwarmInitializeFieldRegister() first");
2845f50eb2eSDave May   if (swarm->field_registration_finalized) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Cannot register additional fields after calling DMSwarmFinalizeFieldRegister() first");
2855f50eb2eSDave May 
2865f50eb2eSDave May   if (type == PETSC_OBJECT) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
2875f50eb2eSDave May   if (type == PETSC_FUNCTION) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
2885f50eb2eSDave May   if (type == PETSC_STRING) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
2895f50eb2eSDave May   if (type == PETSC_STRUCT) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
2905f50eb2eSDave May   if (type == PETSC_DATATYPE_UNKNOWN) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
291b62e03f8SDave May 
292b62e03f8SDave May   switch (type) {
293b62e03f8SDave May     case PETSC_CHAR:
294b62e03f8SDave May       size = sizeof(PetscChar);
295b62e03f8SDave May       break;
296b62e03f8SDave May     case PETSC_SHORT:
297b62e03f8SDave May       size = sizeof(PetscShort);
298b62e03f8SDave May       break;
299b62e03f8SDave May     case PETSC_INT:
300b62e03f8SDave May       size = sizeof(PetscInt);
301b62e03f8SDave May       break;
302b62e03f8SDave May     case PETSC_LONG:
303b62e03f8SDave May       size = sizeof(Petsc64bitInt);
304b62e03f8SDave May       break;
305b62e03f8SDave May     case PETSC_FLOAT:
306b62e03f8SDave May       size = sizeof(PetscFloat);
307b62e03f8SDave May       break;
308b62e03f8SDave May     case PETSC_DOUBLE:
309b62e03f8SDave May       size = sizeof(PetscReal);
310b62e03f8SDave May       break;
311b62e03f8SDave May 
312b62e03f8SDave May     default:
3135f50eb2eSDave May       SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
314b62e03f8SDave May       break;
315b62e03f8SDave May   }
316b62e03f8SDave May 
317b62e03f8SDave May   /* Load a specific data type into data bucket, specifying textual name and its size in bytes */
31852c3ed93SDave May 	ierr = DataBucketRegisterField(swarm->db,"DMSwarmRegisterPetscDatatypeField",fieldname,blocksize*size,NULL);CHKERRQ(ierr);
31952c3ed93SDave May   {
32052c3ed93SDave May     DataField gfield;
32152c3ed93SDave May 
32252c3ed93SDave May     ierr = DataBucketGetDataFieldByName(swarm->db,fieldname,&gfield);CHKERRQ(ierr);
32352c3ed93SDave May     ierr = DataFieldSetBlockSize(gfield,blocksize);CHKERRQ(ierr);
32452c3ed93SDave May   }
325b62e03f8SDave May   swarm->db->field[swarm->db->nfields-1]->petsc_type = type;
326b62e03f8SDave May 
327b62e03f8SDave May   PetscFunctionReturn(0);
328b62e03f8SDave May }
329b62e03f8SDave May 
330b62e03f8SDave May #undef __FUNCT__
331b62e03f8SDave May #define __FUNCT__ "DMSwarmRegisterUserStructField"
3325f50eb2eSDave May PETSC_EXTERN PetscErrorCode DMSwarmRegisterUserStructField(DM dm,const char fieldname[],size_t size)
333b62e03f8SDave May {
3342eac95f8SDave May   PetscErrorCode ierr;
335b62e03f8SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
336b62e03f8SDave May 
3372eac95f8SDave May 	ierr = DataBucketRegisterField(swarm->db,"DMSwarmRegisterUserStructField",fieldname,size,NULL);CHKERRQ(ierr);
338b62e03f8SDave May   swarm->db->field[swarm->db->nfields-1]->petsc_type = PETSC_STRUCT ;
339b62e03f8SDave May 
340b62e03f8SDave May   PetscFunctionReturn(0);
341b62e03f8SDave May }
342b62e03f8SDave May 
343b62e03f8SDave May #undef __FUNCT__
344b62e03f8SDave May #define __FUNCT__ "DMSwarmRegisterUserDatatypeField"
345320740a0SDave May PETSC_EXTERN PetscErrorCode DMSwarmRegisterUserDatatypeField(DM dm,const char fieldname[],size_t size,PetscInt blocksize)
346b62e03f8SDave May {
347b62e03f8SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
3486845f8f5SDave May   PetscErrorCode ierr;
349b62e03f8SDave May 
350320740a0SDave May 	ierr = DataBucketRegisterField(swarm->db,"DMSwarmRegisterUserDatatypeField",fieldname,blocksize*size,NULL);CHKERRQ(ierr);
351320740a0SDave May   {
352320740a0SDave May     DataField gfield;
353320740a0SDave May 
354320740a0SDave May     ierr = DataBucketGetDataFieldByName(swarm->db,fieldname,&gfield);CHKERRQ(ierr);
355320740a0SDave May     ierr = DataFieldSetBlockSize(gfield,blocksize);CHKERRQ(ierr);
356320740a0SDave May   }
357b62e03f8SDave May   swarm->db->field[swarm->db->nfields-1]->petsc_type = PETSC_DATATYPE_UNKNOWN;
358b62e03f8SDave May 
359b62e03f8SDave May   PetscFunctionReturn(0);
360b62e03f8SDave May }
361b62e03f8SDave May 
362b62e03f8SDave May #undef __FUNCT__
363b62e03f8SDave May #define __FUNCT__ "DMSwarmGetField"
3645f50eb2eSDave May PETSC_EXTERN PetscErrorCode DMSwarmGetField(DM dm,const char fieldname[],PetscInt *blocksize,PetscDataType *type,void **data)
365b62e03f8SDave May {
366b62e03f8SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
367b62e03f8SDave May   DataField gfield;
3682eac95f8SDave May   PetscErrorCode ierr;
369b62e03f8SDave May 
3703454631fSDave May   if (!swarm->issetup) { ierr = DMSetUp(dm);CHKERRQ(ierr); }
3713454631fSDave May 
3722eac95f8SDave May   ierr = DataBucketGetDataFieldByName(swarm->db,fieldname,&gfield);CHKERRQ(ierr);
3736845f8f5SDave May   ierr = DataFieldGetAccess(gfield);CHKERRQ(ierr);
3746845f8f5SDave May   ierr = DataFieldGetEntries(gfield,data);CHKERRQ(ierr);
3751b1ea282SDave May   if (blocksize) {*blocksize = gfield->bs; }
376b5bcf523SDave May   if (type) { *type = gfield->petsc_type; }
377b62e03f8SDave May 
378b62e03f8SDave May   PetscFunctionReturn(0);
379b62e03f8SDave May }
380b62e03f8SDave May 
381b62e03f8SDave May #undef __FUNCT__
382b62e03f8SDave May #define __FUNCT__ "DMSwarmRestoreField"
3835f50eb2eSDave May PETSC_EXTERN PetscErrorCode DMSwarmRestoreField(DM dm,const char fieldname[],PetscInt *blocksize,PetscDataType *type,void **data)
384b62e03f8SDave May {
385b62e03f8SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
386b62e03f8SDave May   DataField gfield;
3872eac95f8SDave May   PetscErrorCode ierr;
388b62e03f8SDave May 
3892eac95f8SDave May   ierr = DataBucketGetDataFieldByName(swarm->db,fieldname,&gfield);CHKERRQ(ierr);
3906845f8f5SDave May   ierr = DataFieldRestoreAccess(gfield);CHKERRQ(ierr);
391b62e03f8SDave May   if (data) *data = NULL;
392b62e03f8SDave May 
393b62e03f8SDave May   PetscFunctionReturn(0);
394b62e03f8SDave May }
395b62e03f8SDave May 
396cb1d1399SDave May #undef __FUNCT__
397cb1d1399SDave May #define __FUNCT__ "DMSwarmAddPoint"
398cb1d1399SDave May PETSC_EXTERN PetscErrorCode DMSwarmAddPoint(DM dm)
399cb1d1399SDave May {
400cb1d1399SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
401cb1d1399SDave May   PetscErrorCode ierr;
402cb1d1399SDave May 
4033454631fSDave May   if (!swarm->issetup) { ierr = DMSetUp(dm);CHKERRQ(ierr); }
404cb1d1399SDave May   ierr = DataBucketAddPoint(swarm->db);CHKERRQ(ierr);
405cb1d1399SDave May   PetscFunctionReturn(0);
406cb1d1399SDave May }
407cb1d1399SDave May 
408cb1d1399SDave May #undef __FUNCT__
409cb1d1399SDave May #define __FUNCT__ "DMSwarmAddNPoints"
410cb1d1399SDave May PETSC_EXTERN PetscErrorCode DMSwarmAddNPoints(DM dm,PetscInt npoints)
411cb1d1399SDave May {
412cb1d1399SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
413cb1d1399SDave May   PetscErrorCode ierr;
414cb1d1399SDave May   PetscInt nlocal;
415cb1d1399SDave May 
416cb1d1399SDave May   ierr = DataBucketGetSizes(swarm->db,&nlocal,NULL,NULL);CHKERRQ(ierr);
417cb1d1399SDave May   nlocal = nlocal + npoints;
418cb1d1399SDave May   ierr = DataBucketSetSizes(swarm->db,nlocal,-1);CHKERRQ(ierr);
419cb1d1399SDave May   PetscFunctionReturn(0);
420cb1d1399SDave May }
421cb1d1399SDave May 
422cb1d1399SDave May #undef __FUNCT__
423cb1d1399SDave May #define __FUNCT__ "DMSwarmRemovePoint"
424cb1d1399SDave May PETSC_EXTERN PetscErrorCode DMSwarmRemovePoint(DM dm)
425cb1d1399SDave May {
426cb1d1399SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
427cb1d1399SDave May   PetscErrorCode ierr;
428cb1d1399SDave May 
429cb1d1399SDave May   ierr = DataBucketRemovePoint(swarm->db);CHKERRQ(ierr);
430cb1d1399SDave May   PetscFunctionReturn(0);
431cb1d1399SDave May }
432cb1d1399SDave May 
433cb1d1399SDave May #undef __FUNCT__
434cb1d1399SDave May #define __FUNCT__ "DMSwarmRemovePointAtIndex"
435cb1d1399SDave May PETSC_EXTERN PetscErrorCode DMSwarmRemovePointAtIndex(DM dm,PetscInt idx)
436cb1d1399SDave May {
437cb1d1399SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
438cb1d1399SDave May   PetscErrorCode ierr;
439cb1d1399SDave May 
440cb1d1399SDave May   ierr = DataBucketRemovePointAtIndex(swarm->db,idx);CHKERRQ(ierr);
441cb1d1399SDave May   PetscFunctionReturn(0);
442cb1d1399SDave May }
443b62e03f8SDave May 
4443454631fSDave May #undef __FUNCT__
445095059a4SDave May #define __FUNCT__ "DMSwarmMigrate_Basic"
446095059a4SDave May PetscErrorCode DMSwarmMigrate_Basic(DM dm,PetscBool remove_sent_points)
4473454631fSDave May {
448dcf43ee8SDave May   PetscErrorCode ierr;
449dcf43ee8SDave May   ierr = DMSwarmMigrate_Push_Basic(dm,remove_sent_points);CHKERRQ(ierr);
4503454631fSDave May   PetscFunctionReturn(0);
4513454631fSDave May }
4523454631fSDave May 
453*f0cdbbbaSDave May PetscErrorCode DMSwarmMigrate_CellDMScatter(DM dm,PetscBool remove_sent_points);
454*f0cdbbbaSDave May PetscErrorCode DMSwarmMigrate_CellDMExact(DM dm,PetscBool remove_sent_points);
455*f0cdbbbaSDave May 
4563454631fSDave May #undef __FUNCT__
457095059a4SDave May #define __FUNCT__ "DMSwarmMigrate"
458095059a4SDave May PETSC_EXTERN PetscErrorCode DMSwarmMigrate(DM dm,PetscBool remove_sent_points)
4593454631fSDave May {
460*f0cdbbbaSDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
4613454631fSDave May   PetscErrorCode ierr;
462*f0cdbbbaSDave May 
463*f0cdbbbaSDave May   switch (swarm->migrate_type) {
464*f0cdbbbaSDave May 
465*f0cdbbbaSDave May     case DMSWARM_MIGRATE_BASIC:
466095059a4SDave May       ierr = DMSwarmMigrate_Basic(dm,remove_sent_points);CHKERRQ(ierr);
467*f0cdbbbaSDave May       break;
468*f0cdbbbaSDave May 
469*f0cdbbbaSDave May     case DMSWARM_MIGRATE_DMCELLNSCATTER:
470*f0cdbbbaSDave May       ierr = DMSwarmMigrate_CellDMScatter(dm,remove_sent_points);CHKERRQ(ierr);
471*f0cdbbbaSDave May       break;
472*f0cdbbbaSDave May 
473*f0cdbbbaSDave May     case DMSWARM_MIGRATE_DMCELLEXACT:
474*f0cdbbbaSDave May       SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DMSWARM_MIGRATE_DMCELLEXACT not implemented");
475*f0cdbbbaSDave May       //ierr = DMSwarmMigrate_CellDMExact(dm,remove_sent_points);CHKERRQ(ierr);
476*f0cdbbbaSDave May       break;
477*f0cdbbbaSDave May 
478*f0cdbbbaSDave May     case DMSWARM_MIGRATE_USER:
479*f0cdbbbaSDave May       SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DMSWARM_MIGRATE_USER not implemented");
480*f0cdbbbaSDave May       //ierr = swarm->migrate(dm,remove_sent_points);CHKERRQ(ierr);
481*f0cdbbbaSDave May       break;
482*f0cdbbbaSDave May 
483*f0cdbbbaSDave May     default:
484*f0cdbbbaSDave May       SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DMSWARM_MIGRATE type unknown");
485*f0cdbbbaSDave May       break;
486*f0cdbbbaSDave May   }
4873454631fSDave May   PetscFunctionReturn(0);
4883454631fSDave May }
4893454631fSDave May 
490*f0cdbbbaSDave May PetscErrorCode DMSwarmMigrate_GlobalToLocal_Basic(DM dm,PetscInt *globalsize);
491*f0cdbbbaSDave May 
4922712d1f2SDave May #undef __FUNCT__
493fe39f135SDave May #define __FUNCT__ "DMSwarmCollectViewCreate"
494fe39f135SDave May PETSC_EXTERN PetscErrorCode DMSwarmCollectViewCreate(DM dm)
4952712d1f2SDave May {
4962712d1f2SDave May   PetscErrorCode ierr;
4972712d1f2SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
4982712d1f2SDave May   PetscInt ng;
4992712d1f2SDave May 
500480eef7bSDave May   if (swarm->collect_view_active) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"CollectView currently active");
5012712d1f2SDave May 
502480eef7bSDave May   ierr = DMSwarmGetLocalSize(dm,&ng);CHKERRQ(ierr);
503480eef7bSDave May   switch (swarm->collect_type) {
504*f0cdbbbaSDave May 
505480eef7bSDave May     case DMSWARM_COLLECT_BASIC:
5062712d1f2SDave May       ierr = DMSwarmMigrate_GlobalToLocal_Basic(dm,&ng);CHKERRQ(ierr);
507480eef7bSDave May       break;
508*f0cdbbbaSDave May 
509480eef7bSDave May     case DMSWARM_COLLECT_DMDABOUNDINGBOX:
510*f0cdbbbaSDave May       SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DMSWARM_COLLECT_DMDABOUNDINGBOX not implemented");
511fe39f135SDave May       //ierr = DMSwarmCollect_DMDABoundingBox(dm,&ng);CHKERRQ(ierr);
512480eef7bSDave May       break;
513*f0cdbbbaSDave May 
514480eef7bSDave May     case DMSWARM_COLLECT_GENERAL:
515*f0cdbbbaSDave May       SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DMSWARM_COLLECT_GENERAL not implemented");
516fe39f135SDave May       //ierr = DMSwarmCollect_General(dm,..,,..,&ng);CHKERRQ(ierr);
517480eef7bSDave May       break;
518480eef7bSDave May 
519480eef7bSDave May     default:
520*f0cdbbbaSDave May       SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"DMSWARM_COLLECT type unknown");
521480eef7bSDave May       break;
522480eef7bSDave May   }
523480eef7bSDave May 
524480eef7bSDave May   swarm->collect_view_active = PETSC_TRUE;
525480eef7bSDave May   swarm->collect_view_reset_nlocal = ng;
5262712d1f2SDave May 
5272712d1f2SDave May   PetscFunctionReturn(0);
5282712d1f2SDave May }
5292712d1f2SDave May 
5302712d1f2SDave May #undef __FUNCT__
531fe39f135SDave May #define __FUNCT__ "DMSwarmCollectViewDestroy"
532fe39f135SDave May PETSC_EXTERN PetscErrorCode DMSwarmCollectViewDestroy(DM dm)
5332712d1f2SDave May {
5342712d1f2SDave May   PetscErrorCode ierr;
5352712d1f2SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
5362712d1f2SDave May 
537480eef7bSDave May   if (!swarm->collect_view_active) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"CollectView is currently not active");
538480eef7bSDave May   ierr = DMSwarmSetLocalSizes(dm,swarm->collect_view_reset_nlocal,-1);CHKERRQ(ierr);
539480eef7bSDave May   swarm->collect_view_active = PETSC_FALSE;
5402712d1f2SDave May 
5412712d1f2SDave May   PetscFunctionReturn(0);
5422712d1f2SDave May }
5433454631fSDave May 
5443454631fSDave May #undef __FUNCT__
545*f0cdbbbaSDave May #define __FUNCT__ "DMSwarmSetUpPIC"
546*f0cdbbbaSDave May PetscErrorCode DMSwarmSetUpPIC(DM dm)
547*f0cdbbbaSDave May {
548*f0cdbbbaSDave May   PetscInt dim;
549*f0cdbbbaSDave May   PetscErrorCode ierr;
550*f0cdbbbaSDave May 
551*f0cdbbbaSDave May   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
552*f0cdbbbaSDave May   if (dim < 1) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Dimension must be 1,2,3 - found %D",dim);
553*f0cdbbbaSDave May   if (dim > 3) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Dimension must be 1,2,3 - found %D",dim);
554*f0cdbbbaSDave May   ierr = DMSwarmRegisterPetscDatatypeField(dm,DMSwarmPICField_coor,dim,PETSC_DOUBLE);CHKERRQ(ierr);
555*f0cdbbbaSDave May   PetscFunctionReturn(0);
556*f0cdbbbaSDave May }
557*f0cdbbbaSDave May 
558*f0cdbbbaSDave May #undef __FUNCT__
559*f0cdbbbaSDave May #define __FUNCT__ "DMSwarmSetType"
560*f0cdbbbaSDave May PETSC_EXTERN PetscErrorCode DMSwarmSetType(DM dm,DMSwarmType stype)
561*f0cdbbbaSDave May {
562*f0cdbbbaSDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
563*f0cdbbbaSDave May   PetscErrorCode ierr;
564*f0cdbbbaSDave May 
565*f0cdbbbaSDave May   swarm->swarm_type = stype;
566*f0cdbbbaSDave May   if (swarm->swarm_type == DMSWARM_PIC) {
567*f0cdbbbaSDave May     ierr = DMSwarmSetUpPIC(dm);CHKERRQ(ierr);
568*f0cdbbbaSDave May   }
569*f0cdbbbaSDave May   PetscFunctionReturn(0);
570*f0cdbbbaSDave May }
571*f0cdbbbaSDave May 
572*f0cdbbbaSDave May #undef __FUNCT__
5733454631fSDave May #define __FUNCT__ "DMSetup_Swarm"
5743454631fSDave May PetscErrorCode DMSetup_Swarm(DM dm)
5753454631fSDave May {
5763454631fSDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
5773454631fSDave May   PetscErrorCode ierr;
5783454631fSDave May   PetscMPIInt rank;
5793454631fSDave May   PetscInt p,npoints,*rankval;
5803454631fSDave May 
5813454631fSDave May   if (swarm->issetup) PetscFunctionReturn(0);
5823454631fSDave May 
5833454631fSDave May   swarm->issetup = PETSC_TRUE;
5843454631fSDave May 
585*f0cdbbbaSDave May   if (swarm->swarm_type == DMSWARM_PIC) {
586*f0cdbbbaSDave May     /* check dmcell exists */
587*f0cdbbbaSDave May     if (!swarm->dmcell) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"DMSWARM_PIC requires you call DMSwarmSetCellDM");
588*f0cdbbbaSDave May 
589*f0cdbbbaSDave May     if (swarm->dmcell->ops->locatepointssubdomain) {
590*f0cdbbbaSDave May       /* check methods exists for exact ownership identificiation */
591*f0cdbbbaSDave May       PetscPrintf(PetscObjectComm((PetscObject)dm),"  DMSWARM_PIC: Using method CellDM->ops->LocatePointsSubdomain\n");
592*f0cdbbbaSDave May       swarm->migrate_type = DMSWARM_MIGRATE_DMCELLEXACT;
593*f0cdbbbaSDave May     } else {
594*f0cdbbbaSDave May       /* check methods exist for point location AND rank neighbor identification */
595*f0cdbbbaSDave May       if (swarm->dmcell->ops->locatepoints) {
596*f0cdbbbaSDave May         PetscPrintf(PetscObjectComm((PetscObject)dm),"  DMSWARM_PIC: Using method CellDM->LocatePoints\n");
597*f0cdbbbaSDave May       } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"DMSWARM_PIC requires the method CellDM->ops->locatepoints be defined");
598*f0cdbbbaSDave May 
599*f0cdbbbaSDave May       if (swarm->dmcell->ops->getneighbors) {
600*f0cdbbbaSDave May         PetscPrintf(PetscObjectComm((PetscObject)dm),"  DMSWARM_PIC: Using method CellDM->GetNeigbors\n");
601*f0cdbbbaSDave May       } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"DMSWARM_PIC requires the method CellDM->ops->getneighbors be defined");
602*f0cdbbbaSDave May 
603*f0cdbbbaSDave May       swarm->migrate_type = DMSWARM_MIGRATE_DMCELLNSCATTER;
604*f0cdbbbaSDave May     }
605*f0cdbbbaSDave May   }
606*f0cdbbbaSDave May 
607*f0cdbbbaSDave May   ierr = DMSwarmFinalizeFieldRegister(dm);CHKERRQ(ierr);
608*f0cdbbbaSDave May 
6093454631fSDave May   /* check some fields were registered */
6103454631fSDave May   if (swarm->db->nfields <= 2) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"At least one field user must be registered via DMSwarmRegisterXXX()");
6113454631fSDave May 
6123454631fSDave May   /* check local sizes were set */
6133454631fSDave May   if (swarm->db->L == -1) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Local sizes must be set via DMSwarmSetLocalSizes()");
6143454631fSDave May 
6153454631fSDave May   /* initialize values in pid and rank placeholders */
6163454631fSDave May   /* TODO: [pid - use MPI_Scan] */
6173454631fSDave May 
6183454631fSDave May   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
6193454631fSDave May   ierr = DataBucketGetSizes(swarm->db,&npoints,NULL,NULL);CHKERRQ(ierr);
620*f0cdbbbaSDave May   ierr = DMSwarmGetField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
6213454631fSDave May   for (p=0; p<npoints; p++) {
6223454631fSDave May     rankval[p] = (PetscInt)rank;
6233454631fSDave May   }
624*f0cdbbbaSDave May   ierr = DMSwarmRestoreField(dm,DMSwarmField_rank,NULL,NULL,(void**)&rankval);CHKERRQ(ierr);
6253454631fSDave May 
6263454631fSDave May   PetscFunctionReturn(0);
6273454631fSDave May }
6283454631fSDave May 
629b62e03f8SDave May #undef __FUNCT__
63057795646SDave May #define __FUNCT__ "DMDestroy_Swarm"
63157795646SDave May PetscErrorCode DMDestroy_Swarm(DM dm)
63257795646SDave May {
63357795646SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
63457795646SDave May   PetscErrorCode ierr;
63557795646SDave May 
63657795646SDave May   PetscFunctionBegin;
6376845f8f5SDave May   ierr = DataBucketDestroy(&swarm->db);CHKERRQ(ierr);
63857795646SDave May   ierr = PetscFree(swarm);CHKERRQ(ierr);
63957795646SDave May   PetscFunctionReturn(0);
64057795646SDave May }
64157795646SDave May 
64257795646SDave May #undef __FUNCT__
6435f50eb2eSDave May #define __FUNCT__ "DMView_Swarm"
6445f50eb2eSDave May PetscErrorCode DMView_Swarm(DM dm, PetscViewer viewer)
6455f50eb2eSDave May {
6465f50eb2eSDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
6475f50eb2eSDave May   PetscBool      iascii,ibinary,ishdf5,isvtk;
6485f50eb2eSDave May   PetscErrorCode ierr;
6495f50eb2eSDave May 
6505f50eb2eSDave May   PetscFunctionBegin;
6515f50eb2eSDave May   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
6525f50eb2eSDave May   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
6535f50eb2eSDave May   ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
6545f50eb2eSDave May   ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY,&ibinary);CHKERRQ(ierr);
6555f50eb2eSDave May   ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERVTK,   &isvtk);CHKERRQ(ierr);
6565f50eb2eSDave May   ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5,  &ishdf5);CHKERRQ(ierr);
6575f50eb2eSDave May   if (iascii) {
6586845f8f5SDave May     ierr = DataBucketView(PetscObjectComm((PetscObject)dm),swarm->db,NULL,DATABUCKET_VIEW_STDOUT);CHKERRQ(ierr);
6595f50eb2eSDave May   } else if (ibinary) {
6605f50eb2eSDave May     SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"NO VTK support");
6615f50eb2eSDave May   } else if (ishdf5) {
6625f50eb2eSDave May #if defined(PETSC_HAVE_HDF5)
6635f50eb2eSDave May     SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"NO HDF5 support");
6645f50eb2eSDave May #else
6655f50eb2eSDave May     SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"HDF5 not supported. Please reconfigure using --download-hdf5");
6665f50eb2eSDave May #endif
6675f50eb2eSDave May   } else if (isvtk) {
6685f50eb2eSDave May     SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"NO VTK support");
6695f50eb2eSDave May   }
6705f50eb2eSDave May   PetscFunctionReturn(0);
6715f50eb2eSDave May }
6725f50eb2eSDave May 
6735f50eb2eSDave May #undef __FUNCT__
67457795646SDave May #define __FUNCT__ "DMCreate_Swarm"
67557795646SDave May PETSC_EXTERN PetscErrorCode DMCreate_Swarm(DM dm)
67657795646SDave May {
67757795646SDave May   DM_Swarm      *swarm;
67857795646SDave May   PetscErrorCode ierr;
67957795646SDave May 
68057795646SDave May   PetscFunctionBegin;
68157795646SDave May   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
68257795646SDave May   ierr     = PetscNewLog(dm,&swarm);CHKERRQ(ierr);
683*f0cdbbbaSDave May   dm->data = swarm;
68457795646SDave May 
6856845f8f5SDave May   ierr = DataBucketCreate(&swarm->db);CHKERRQ(ierr);
686*f0cdbbbaSDave May   ierr = DMSwarmInitializeFieldRegister(dm);CHKERRQ(ierr);
687*f0cdbbbaSDave May 
688b5bcf523SDave May   swarm->vec_field_set = PETSC_FALSE;
6893454631fSDave May   swarm->issetup = PETSC_FALSE;
690480eef7bSDave May   swarm->swarm_type = DMSWARM_BASIC;
691480eef7bSDave May   swarm->migrate_type = DMSWARM_MIGRATE_BASIC;
692480eef7bSDave May   swarm->collect_type = DMSWARM_COLLECT_BASIC;
69340c453e9SDave May   swarm->migrate_error_on_missing_point = PETSC_FALSE;
694b62e03f8SDave May 
695*f0cdbbbaSDave May   swarm->dmcell = NULL;
696*f0cdbbbaSDave May   swarm->collect_view_active = PETSC_FALSE;
697*f0cdbbbaSDave May   swarm->collect_view_reset_nlocal = -1;
69857795646SDave May 
699*f0cdbbbaSDave May   dm->dim  = 0;
7005f50eb2eSDave May   dm->ops->view                            = DMView_Swarm;
70157795646SDave May   dm->ops->load                            = NULL;
70257795646SDave May   dm->ops->setfromoptions                  = NULL;
70357795646SDave May   dm->ops->clone                           = NULL;
7043454631fSDave May   dm->ops->setup                           = DMSetup_Swarm;
70557795646SDave May   dm->ops->createdefaultsection            = NULL;
70657795646SDave May   dm->ops->createdefaultconstraints        = NULL;
707b5bcf523SDave May   dm->ops->createglobalvector              = DMCreateGlobalVector_Swarm;
708b5bcf523SDave May   dm->ops->createlocalvector               = DMCreateLocalVector_Swarm;
70957795646SDave May   dm->ops->getlocaltoglobalmapping         = NULL;
71057795646SDave May   dm->ops->createfieldis                   = NULL;
71157795646SDave May   dm->ops->createcoordinatedm              = NULL;
71257795646SDave May   dm->ops->getcoloring                     = NULL;
71357795646SDave May   dm->ops->creatematrix                    = NULL;
71457795646SDave May   dm->ops->createinterpolation             = NULL;
71557795646SDave May   dm->ops->getaggregates                   = NULL;
71657795646SDave May   dm->ops->getinjection                    = NULL;
71757795646SDave May   dm->ops->refine                          = NULL;
71857795646SDave May   dm->ops->coarsen                         = NULL;
71957795646SDave May   dm->ops->refinehierarchy                 = NULL;
72057795646SDave May   dm->ops->coarsenhierarchy                = NULL;
72157795646SDave May   dm->ops->globaltolocalbegin              = NULL;
72257795646SDave May   dm->ops->globaltolocalend                = NULL;
72357795646SDave May   dm->ops->localtoglobalbegin              = NULL;
72457795646SDave May   dm->ops->localtoglobalend                = NULL;
72557795646SDave May   dm->ops->destroy                         = DMDestroy_Swarm;
72657795646SDave May   dm->ops->createsubdm                     = NULL;
72757795646SDave May   dm->ops->getdimpoints                    = NULL;
72857795646SDave May   dm->ops->locatepoints                    = NULL;
72957795646SDave May 
73057795646SDave May   PetscFunctionReturn(0);
73157795646SDave May }