xref: /petsc/src/dm/impls/swarm/swarm.c (revision cb1d1399e4986c428df25dfe271dbdc012ec953d)
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 
657795646SDave May //typedef PetscErrorCode (*swarm_project)(DM,DM,Vec) DMSwarmProjectMethod; /* swarm, geometry, result */
757795646SDave May 
857795646SDave May //typedef enum { PROJECT_DMDA_AQ1=0, PROJECT_DMDA_P0 } DMSwarmDMDAProjectionType;
957795646SDave May 
1057795646SDave May #if 0
1157795646SDave May 
1257795646SDave May /* Defines what the local space will be */
13b62e03f8SDave May PetscErrorCode DMSwarmSetOverlap(void)
1457795646SDave May {
1557795646SDave May 
16b62e03f8SDave May   PetscFunctionReturn(0);
1757795646SDave May }
1857795646SDave May 
1957795646SDave May 
2057795646SDave May /* coordinates */
2157795646SDave May /*
2257795646SDave May DMGetCoordinateDM returns self
2357795646SDave May DMGetCoordinates and DMGetCoordinatesLocal return same thing
2457795646SDave May Local view could be used to define overlapping information
2557795646SDave May */
2657795646SDave May 
2757795646SDave May #endif
2857795646SDave May 
2957795646SDave May #undef __FUNCT__
30b5bcf523SDave May #define __FUNCT__ "DMSwarmVectorDefineField"
31b5bcf523SDave May PETSC_EXTERN PetscErrorCode DMSwarmVectorDefineField(DM dm,const char fieldname[])
32b5bcf523SDave May {
33b5bcf523SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
34b5bcf523SDave May   PetscErrorCode ierr;
35b5bcf523SDave May   PetscInt bs,n;
36b5bcf523SDave May   PetscScalar *array;
37b5bcf523SDave May   PetscDataType type;
38b5bcf523SDave May 
396845f8f5SDave May   ierr = DataBucketGetSizes(swarm->db,&n,NULL,NULL);CHKERRQ(ierr);
40b5bcf523SDave May   ierr = DMSwarmGetField(dm,fieldname,&bs,&type,(void**)&array);CHKERRQ(ierr);
41b5bcf523SDave May 
42b5bcf523SDave May   /* Check all fields are of type PETSC_REAL or PETSC_SCALAR */
43b5bcf523SDave May   if (type != PETSC_REAL) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only valid for PETSC_REAL");
44b5bcf523SDave May 
45b5bcf523SDave May   PetscSNPrintf(swarm->vec_field_name,PETSC_MAX_PATH_LEN-1,"%s",fieldname);
46b5bcf523SDave May   swarm->vec_field_set = PETSC_TRUE;
47b5bcf523SDave May   swarm->vec_field_bs = 1;//bs;
48b5bcf523SDave May   swarm->vec_field_nlocal = n;
49b5bcf523SDave May 
50b5bcf523SDave May   PetscFunctionReturn(0);
51b5bcf523SDave May }
52b5bcf523SDave May 
53b5bcf523SDave May /* requires DMSwarmDefineFieldVector has been called */
54b5bcf523SDave May #undef __FUNCT__
55b5bcf523SDave May #define __FUNCT__ "DMCreateGlobalVector_Swarm"
56b5bcf523SDave May PetscErrorCode DMCreateGlobalVector_Swarm(DM dm,Vec *vec)
57b5bcf523SDave May {
58b5bcf523SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
59b5bcf523SDave May   PetscErrorCode ierr;
60b5bcf523SDave May   Vec x;
61b5bcf523SDave May   char name[PETSC_MAX_PATH_LEN];
62b5bcf523SDave May 
63b5bcf523SDave May   if (!swarm->vec_field_set) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMSwarmVectorDefineField first");
64b5bcf523SDave May   PetscSNPrintf(name,PETSC_MAX_PATH_LEN-1,"DMSwarmField_%s",swarm->vec_field_name);
65b5bcf523SDave May   ierr = VecCreate(PetscObjectComm((PetscObject)dm),&x);CHKERRQ(ierr);
66b5bcf523SDave May   ierr = PetscObjectSetName((PetscObject)x,name);CHKERRQ(ierr);
67b5bcf523SDave May   ierr = VecSetSizes(x,swarm->vec_field_nlocal,PETSC_DETERMINE);CHKERRQ(ierr);
68b5bcf523SDave May   ierr = VecSetBlockSize(x,swarm->vec_field_bs);CHKERRQ(ierr);
69b5bcf523SDave May   ierr = VecSetFromOptions(x);CHKERRQ(ierr);
70b5bcf523SDave May   *vec = x;
71b5bcf523SDave May 
72b5bcf523SDave May   PetscFunctionReturn(0);
73b5bcf523SDave May }
74b5bcf523SDave May 
75b5bcf523SDave May /* requires DMSwarmDefineFieldVector has been called */
76b5bcf523SDave May #undef __FUNCT__
77b5bcf523SDave May #define __FUNCT__ "DMCreateLocalVector_Swarm"
78b5bcf523SDave May PetscErrorCode DMCreateLocalVector_Swarm(DM dm,Vec *vec)
79b5bcf523SDave May {
80b5bcf523SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
81b5bcf523SDave May   PetscErrorCode ierr;
82b5bcf523SDave May   Vec x;
83b5bcf523SDave May 
84b5bcf523SDave May   char name[PETSC_MAX_PATH_LEN];
85b5bcf523SDave May 
86b5bcf523SDave May   if (!swarm->vec_field_set) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMSwarmVectorDefineField first");
87b5bcf523SDave May   PetscSNPrintf(name,PETSC_MAX_PATH_LEN-1,"DMSwarmField_%s",swarm->vec_field_name);
88b5bcf523SDave May   ierr = VecCreate(PETSC_COMM_SELF,&x);CHKERRQ(ierr);
89b5bcf523SDave May   ierr = PetscObjectSetName((PetscObject)x,name);CHKERRQ(ierr);
90b5bcf523SDave May   ierr = VecSetSizes(x,swarm->vec_field_nlocal,swarm->vec_field_nlocal);CHKERRQ(ierr);
91b5bcf523SDave May   ierr = VecSetBlockSize(x,swarm->vec_field_bs);CHKERRQ(ierr);
92b5bcf523SDave May   ierr = VecSetFromOptions(x);CHKERRQ(ierr);
93b5bcf523SDave May   *vec = x;
94b5bcf523SDave May 
95b5bcf523SDave May   PetscFunctionReturn(0);
96b5bcf523SDave May }
97b5bcf523SDave May 
98b5bcf523SDave May #undef __FUNCT__
99b5bcf523SDave May #define __FUNCT__ "DMSwarmCreateGlobalVectorFromField"
100b5bcf523SDave May PETSC_EXTERN PetscErrorCode DMSwarmCreateGlobalVectorFromField(DM dm,const char fieldname[],Vec *vec)
101b5bcf523SDave May {
102b5bcf523SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
103b5bcf523SDave May   PetscErrorCode ierr;
104b5bcf523SDave May   PetscInt bs,n;
105b5bcf523SDave May   PetscScalar *array;
106b5bcf523SDave May   Vec x;
107b5bcf523SDave May   PetscDataType type;
108b5bcf523SDave May   char name[PETSC_MAX_PATH_LEN];
109b5bcf523SDave May 
1106845f8f5SDave May   ierr = DataBucketGetSizes(swarm->db,&n,NULL,NULL);CHKERRQ(ierr);
111b5bcf523SDave May   ierr = DMSwarmGetField(dm,fieldname,&bs,&type,(void**)&array);CHKERRQ(ierr);
112b5bcf523SDave May 
113b5bcf523SDave May   /* Check all fields are of type PETSC_REAL or PETSC_SCALAR */
114b5bcf523SDave May   if (type != PETSC_REAL) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Only valid for PETSC_REAL");
115b5bcf523SDave May 
116b5bcf523SDave May   ierr = VecCreateSeqWithArray(PetscObjectComm((PetscObject)dm),1,n,array,&x);CHKERRQ(ierr);
117b5bcf523SDave May 
118b5bcf523SDave May   PetscSNPrintf(name,PETSC_MAX_PATH_LEN-1,"DMSwarm_VecFieldInPlace_%s",fieldname);
119b5bcf523SDave May   ierr = PetscObjectComposeFunction((PetscObject)x,name,DMSwarmDestroyGlobalVectorFromField);CHKERRQ(ierr);
120b5bcf523SDave May 
121b5bcf523SDave May   /* Set guard */
122b5bcf523SDave May 
123b5bcf523SDave May   *vec = x;
124b5bcf523SDave May   PetscFunctionReturn(0);
125b5bcf523SDave May }
126b5bcf523SDave May 
127b5bcf523SDave May #undef __FUNCT__
128b5bcf523SDave May #define __FUNCT__ "DMSwarmDestroyGlobalVectorFromField"
129b5bcf523SDave May PETSC_EXTERN PetscErrorCode DMSwarmDestroyGlobalVectorFromField(DM dm,const char fieldname[],Vec *vec)
130b5bcf523SDave May {
131b5bcf523SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
132b5bcf523SDave May   PetscErrorCode ierr;
133b5bcf523SDave May   DataField gfield;
134b5bcf523SDave May   char name[PETSC_MAX_PATH_LEN];
135b5bcf523SDave May   void (*fptr)(void);
136b5bcf523SDave May 
137b5bcf523SDave May   /* get data field */
1382eac95f8SDave May   ierr = DataBucketGetDataFieldByName(swarm->db,fieldname,&gfield);CHKERRQ(ierr);
139b5bcf523SDave May 
140b5bcf523SDave May   /* check vector is an inplace array */
141b5bcf523SDave May   PetscSNPrintf(name,PETSC_MAX_PATH_LEN-1,"DMSwarm_VecFieldInPlace_%s",fieldname);
142b5bcf523SDave May   ierr = PetscObjectQueryFunction((PetscObject)(*vec),name,&fptr);CHKERRQ(ierr);
143b5bcf523SDave May   if (!fptr) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Vector being destroyed was not created from DMSwarm field(%s)",fieldname);
144b5bcf523SDave May 
145b5bcf523SDave May   /* restore data field */
1466845f8f5SDave May   ierr = DataFieldRestoreAccess(gfield);CHKERRQ(ierr);
147b5bcf523SDave May 
148b5bcf523SDave May   ierr = VecDestroy(vec);CHKERRQ(ierr);
149b5bcf523SDave May 
150b5bcf523SDave May   PetscFunctionReturn(0);
151b5bcf523SDave May }
152b5bcf523SDave May 
153b5bcf523SDave May /*
154b5bcf523SDave May PETSC_EXTERN PetscErrorCode DMSwarmCreateGlobalVectorFromFields(DM dm,const PetscInt nf,const char *fieldnames[],Vec *vec)
155b5bcf523SDave May {
156b5bcf523SDave May   PetscFunctionReturn(0);
157b5bcf523SDave May }
158b5bcf523SDave May 
159b5bcf523SDave May PETSC_EXTERN PetscErrorCode DMSwarmRestoreGlobalVectorFromFields(DM dm,Vec *vec)
160b5bcf523SDave May {
161b5bcf523SDave May   PetscFunctionReturn(0);
162b5bcf523SDave May }
163b5bcf523SDave May */
164b5bcf523SDave May 
165b5bcf523SDave May #undef __FUNCT__
1665f50eb2eSDave May #define __FUNCT__ "DMSwarmInitializeFieldRegister"
1675f50eb2eSDave May PETSC_EXTERN PetscErrorCode DMSwarmInitializeFieldRegister(DM dm)
1685f50eb2eSDave May {
1695f50eb2eSDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
1705f50eb2eSDave May   swarm->field_registration_initialized = PETSC_TRUE;
1715f50eb2eSDave May   PetscFunctionReturn(0);
1725f50eb2eSDave May }
1735f50eb2eSDave May 
1745f50eb2eSDave May #undef __FUNCT__
1755f50eb2eSDave May #define __FUNCT__ "DMSwarmFinalizeFieldRegister"
1765f50eb2eSDave May PETSC_EXTERN PetscErrorCode DMSwarmFinalizeFieldRegister(DM dm)
1775f50eb2eSDave May {
1785f50eb2eSDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
1796845f8f5SDave May   PetscErrorCode ierr;
1806845f8f5SDave May 
1815f50eb2eSDave May   swarm->field_registration_finalized = PETSC_TRUE;
1826845f8f5SDave May   ierr = DataBucketFinalize(swarm->db);CHKERRQ(ierr);
1835f50eb2eSDave May   PetscFunctionReturn(0);
1845f50eb2eSDave May }
1855f50eb2eSDave May 
1865f50eb2eSDave May #undef __FUNCT__
1875f50eb2eSDave May #define __FUNCT__ "DMSwarmSetLocalSizes"
1885f50eb2eSDave May PETSC_EXTERN PetscErrorCode DMSwarmSetLocalSizes(DM dm,PetscInt nlocal,PetscInt buffer)
1895f50eb2eSDave May {
1905f50eb2eSDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
1916845f8f5SDave May   PetscErrorCode ierr;
1925f50eb2eSDave May 
1936845f8f5SDave May   ierr = DataBucketSetSizes(swarm->db,nlocal,buffer);CHKERRQ(ierr);
1945f50eb2eSDave May 
1955f50eb2eSDave May   PetscFunctionReturn(0);
1965f50eb2eSDave May }
1975f50eb2eSDave May 
1985f50eb2eSDave May #undef __FUNCT__
199b62e03f8SDave May #define __FUNCT__ "DMSwarmRegisterPetscDatatypeField"
2005f50eb2eSDave May PETSC_EXTERN PetscErrorCode DMSwarmRegisterPetscDatatypeField(DM dm,const char fieldname[],PetscInt blocksize,PetscDataType type)
201b62e03f8SDave May {
2022eac95f8SDave May   PetscErrorCode ierr;
203b62e03f8SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
204b62e03f8SDave May   size_t size;
205b62e03f8SDave May 
2065f50eb2eSDave May   if (!swarm->field_registration_initialized) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMSwarmInitializeFieldRegister() first");
2075f50eb2eSDave May   if (swarm->field_registration_finalized) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Cannot register additional fields after calling DMSwarmFinalizeFieldRegister() first");
2085f50eb2eSDave May 
2095f50eb2eSDave May   if (type == PETSC_OBJECT) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
2105f50eb2eSDave May   if (type == PETSC_FUNCTION) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
2115f50eb2eSDave May   if (type == PETSC_STRING) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
2125f50eb2eSDave May   if (type == PETSC_STRUCT) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
2135f50eb2eSDave May   if (type == PETSC_DATATYPE_UNKNOWN) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
214b62e03f8SDave May 
215b62e03f8SDave May   switch (type) {
216b62e03f8SDave May     case PETSC_CHAR:
217b62e03f8SDave May       size = sizeof(PetscChar);
218b62e03f8SDave May       break;
219b62e03f8SDave May     case PETSC_SHORT:
220b62e03f8SDave May       size = sizeof(PetscShort);
221b62e03f8SDave May       break;
222b62e03f8SDave May     case PETSC_INT:
223b62e03f8SDave May       size = sizeof(PetscInt);
224b62e03f8SDave May       break;
225b62e03f8SDave May     case PETSC_LONG:
226b62e03f8SDave May       size = sizeof(Petsc64bitInt);
227b62e03f8SDave May       break;
228b62e03f8SDave May     case PETSC_FLOAT:
229b62e03f8SDave May       size = sizeof(PetscFloat);
230b62e03f8SDave May       break;
231b62e03f8SDave May     case PETSC_DOUBLE:
232b62e03f8SDave May       size = sizeof(PetscReal);
233b62e03f8SDave May       break;
234b62e03f8SDave May 
235b62e03f8SDave May     default:
2365f50eb2eSDave May       SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
237b62e03f8SDave May       break;
238b62e03f8SDave May   }
239b62e03f8SDave May 
240b62e03f8SDave May   /* Load a specific data type into data bucket, specifying textual name and its size in bytes */
2412eac95f8SDave May 	ierr = DataBucketRegisterField(swarm->db,"DMSwarmRegisterPetscDatatypeField",fieldname,size,NULL);CHKERRQ(ierr);
242b62e03f8SDave May   swarm->db->field[swarm->db->nfields-1]->petsc_type = type;
243b62e03f8SDave May 
244b62e03f8SDave May   PetscFunctionReturn(0);
245b62e03f8SDave May }
246b62e03f8SDave May 
247b62e03f8SDave May #undef __FUNCT__
248b62e03f8SDave May #define __FUNCT__ "DMSwarmRegisterUserStructField"
2495f50eb2eSDave May PETSC_EXTERN PetscErrorCode DMSwarmRegisterUserStructField(DM dm,const char fieldname[],size_t size)
250b62e03f8SDave May {
2512eac95f8SDave May   PetscErrorCode ierr;
252b62e03f8SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
253b62e03f8SDave May 
2542eac95f8SDave May 	ierr = DataBucketRegisterField(swarm->db,"DMSwarmRegisterUserStructField",fieldname,size,NULL);CHKERRQ(ierr);
255b62e03f8SDave May   swarm->db->field[swarm->db->nfields-1]->petsc_type = PETSC_STRUCT ;
256b62e03f8SDave May 
257b62e03f8SDave May   PetscFunctionReturn(0);
258b62e03f8SDave May }
259b62e03f8SDave May 
260b62e03f8SDave May #undef __FUNCT__
261b62e03f8SDave May #define __FUNCT__ "DMSwarmRegisterUserDatatypeField"
2625f50eb2eSDave May PETSC_EXTERN PetscErrorCode DMSwarmRegisterUserDatatypeField(DM dm,const char fieldname[],size_t size)
263b62e03f8SDave May {
264b62e03f8SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
2656845f8f5SDave May   PetscErrorCode ierr;
266b62e03f8SDave May 
2676845f8f5SDave May 	ierr = DataBucketRegisterField(swarm->db,"DMSwarmRegisterUserDatatypeField",fieldname,size,NULL);CHKERRQ(ierr);
268b62e03f8SDave May   swarm->db->field[swarm->db->nfields-1]->petsc_type = PETSC_DATATYPE_UNKNOWN;
269b62e03f8SDave May 
270b62e03f8SDave May   PetscFunctionReturn(0);
271b62e03f8SDave May }
272b62e03f8SDave May 
273b62e03f8SDave May #undef __FUNCT__
274b62e03f8SDave May #define __FUNCT__ "DMSwarmGetField"
2755f50eb2eSDave May PETSC_EXTERN PetscErrorCode DMSwarmGetField(DM dm,const char fieldname[],PetscInt *blocksize,PetscDataType *type,void **data)
276b62e03f8SDave May {
277b62e03f8SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
278b62e03f8SDave May   DataField gfield;
2792eac95f8SDave May   PetscErrorCode ierr;
280b62e03f8SDave May 
2812eac95f8SDave May   ierr = DataBucketGetDataFieldByName(swarm->db,fieldname,&gfield);CHKERRQ(ierr);
2826845f8f5SDave May   ierr = DataFieldGetAccess(gfield);CHKERRQ(ierr);
2836845f8f5SDave May   ierr = DataFieldGetEntries(gfield,data);CHKERRQ(ierr);
284b5bcf523SDave May   if (blocksize) {}
285b5bcf523SDave May   if (type) { *type = gfield->petsc_type; }
286b62e03f8SDave May 
287b62e03f8SDave May   PetscFunctionReturn(0);
288b62e03f8SDave May }
289b62e03f8SDave May 
290b62e03f8SDave May #undef __FUNCT__
291b62e03f8SDave May #define __FUNCT__ "DMSwarmRestoreField"
2925f50eb2eSDave May PETSC_EXTERN PetscErrorCode DMSwarmRestoreField(DM dm,const char fieldname[],PetscInt *blocksize,PetscDataType *type,void **data)
293b62e03f8SDave May {
294b62e03f8SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
295b62e03f8SDave May   DataField gfield;
2962eac95f8SDave May   PetscErrorCode ierr;
297b62e03f8SDave May 
2982eac95f8SDave May   ierr = DataBucketGetDataFieldByName(swarm->db,fieldname,&gfield);CHKERRQ(ierr);
2996845f8f5SDave May   ierr = DataFieldRestoreAccess(gfield);CHKERRQ(ierr);
300b62e03f8SDave May   if (data) *data = NULL;
301b62e03f8SDave May 
302b62e03f8SDave May   PetscFunctionReturn(0);
303b62e03f8SDave May }
304b62e03f8SDave May 
305*cb1d1399SDave May #undef __FUNCT__
306*cb1d1399SDave May #define __FUNCT__ "DMSwarmAddPoint"
307*cb1d1399SDave May PETSC_EXTERN PetscErrorCode DMSwarmAddPoint(DM dm)
308*cb1d1399SDave May {
309*cb1d1399SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
310*cb1d1399SDave May   PetscErrorCode ierr;
311*cb1d1399SDave May 
312*cb1d1399SDave May   ierr = DataBucketAddPoint(swarm->db);CHKERRQ(ierr);
313*cb1d1399SDave May   PetscFunctionReturn(0);
314*cb1d1399SDave May }
315*cb1d1399SDave May 
316*cb1d1399SDave May #undef __FUNCT__
317*cb1d1399SDave May #define __FUNCT__ "DMSwarmAddNPoints"
318*cb1d1399SDave May PETSC_EXTERN PetscErrorCode DMSwarmAddNPoints(DM dm,PetscInt npoints)
319*cb1d1399SDave May {
320*cb1d1399SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
321*cb1d1399SDave May   PetscErrorCode ierr;
322*cb1d1399SDave May   PetscInt nlocal;
323*cb1d1399SDave May 
324*cb1d1399SDave May   ierr = DataBucketGetSizes(swarm->db,&nlocal,NULL,NULL);CHKERRQ(ierr);
325*cb1d1399SDave May   nlocal = nlocal + npoints;
326*cb1d1399SDave May   ierr = DataBucketSetSizes(swarm->db,nlocal,-1);CHKERRQ(ierr);
327*cb1d1399SDave May   PetscFunctionReturn(0);
328*cb1d1399SDave May }
329*cb1d1399SDave May 
330*cb1d1399SDave May #undef __FUNCT__
331*cb1d1399SDave May #define __FUNCT__ "DMSwarmRemovePoint"
332*cb1d1399SDave May PETSC_EXTERN PetscErrorCode DMSwarmRemovePoint(DM dm)
333*cb1d1399SDave May {
334*cb1d1399SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
335*cb1d1399SDave May   PetscErrorCode ierr;
336*cb1d1399SDave May 
337*cb1d1399SDave May   ierr = DataBucketRemovePoint(swarm->db);CHKERRQ(ierr);
338*cb1d1399SDave May   PetscFunctionReturn(0);
339*cb1d1399SDave May }
340*cb1d1399SDave May 
341*cb1d1399SDave May #undef __FUNCT__
342*cb1d1399SDave May #define __FUNCT__ "DMSwarmRemovePointAtIndex"
343*cb1d1399SDave May PETSC_EXTERN PetscErrorCode DMSwarmRemovePointAtIndex(DM dm,PetscInt idx)
344*cb1d1399SDave May {
345*cb1d1399SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
346*cb1d1399SDave May   PetscErrorCode ierr;
347*cb1d1399SDave May 
348*cb1d1399SDave May   ierr = DataBucketRemovePointAtIndex(swarm->db,idx);CHKERRQ(ierr);
349*cb1d1399SDave May   PetscFunctionReturn(0);
350*cb1d1399SDave May }
351b62e03f8SDave May 
352b62e03f8SDave May #undef __FUNCT__
35357795646SDave May #define __FUNCT__ "DMDestroy_Swarm"
35457795646SDave May PetscErrorCode DMDestroy_Swarm(DM dm)
35557795646SDave May {
35657795646SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
35757795646SDave May   PetscErrorCode ierr;
35857795646SDave May 
35957795646SDave May   PetscFunctionBegin;
3606845f8f5SDave May   ierr = DataBucketDestroy(&swarm->db);CHKERRQ(ierr);
36157795646SDave May   ierr = PetscFree(swarm);CHKERRQ(ierr);
36257795646SDave May   PetscFunctionReturn(0);
36357795646SDave May }
36457795646SDave May 
36557795646SDave May #undef __FUNCT__
3665f50eb2eSDave May #define __FUNCT__ "DMView_Swarm"
3675f50eb2eSDave May PetscErrorCode DMView_Swarm(DM dm, PetscViewer viewer)
3685f50eb2eSDave May {
3695f50eb2eSDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
3705f50eb2eSDave May   PetscBool      iascii,ibinary,ishdf5,isvtk;
3715f50eb2eSDave May   PetscErrorCode ierr;
3725f50eb2eSDave May 
3735f50eb2eSDave May   PetscFunctionBegin;
3745f50eb2eSDave May   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3755f50eb2eSDave May   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
3765f50eb2eSDave May   ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
3775f50eb2eSDave May   ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY,&ibinary);CHKERRQ(ierr);
3785f50eb2eSDave May   ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERVTK,   &isvtk);CHKERRQ(ierr);
3795f50eb2eSDave May   ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5,  &ishdf5);CHKERRQ(ierr);
3805f50eb2eSDave May   if (iascii) {
3816845f8f5SDave May     ierr = DataBucketView(PetscObjectComm((PetscObject)dm),swarm->db,NULL,DATABUCKET_VIEW_STDOUT);CHKERRQ(ierr);
3825f50eb2eSDave May   } else if (ibinary) {
3835f50eb2eSDave May     SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"NO VTK support");
3845f50eb2eSDave May   } else if (ishdf5) {
3855f50eb2eSDave May #if defined(PETSC_HAVE_HDF5)
3865f50eb2eSDave May     SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"NO HDF5 support");
3875f50eb2eSDave May #else
3885f50eb2eSDave May     SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"HDF5 not supported. Please reconfigure using --download-hdf5");
3895f50eb2eSDave May #endif
3905f50eb2eSDave May   } else if (isvtk) {
3915f50eb2eSDave May     SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"NO VTK support");
3925f50eb2eSDave May   }
3935f50eb2eSDave May   PetscFunctionReturn(0);
3945f50eb2eSDave May }
3955f50eb2eSDave May 
3965f50eb2eSDave May #undef __FUNCT__
39757795646SDave May #define __FUNCT__ "DMCreate_Swarm"
39857795646SDave May PETSC_EXTERN PetscErrorCode DMCreate_Swarm(DM dm)
39957795646SDave May {
40057795646SDave May   DM_Swarm      *swarm;
40157795646SDave May   PetscErrorCode ierr;
40257795646SDave May 
40357795646SDave May   PetscFunctionBegin;
40457795646SDave May   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
40557795646SDave May   ierr     = PetscNewLog(dm,&swarm);CHKERRQ(ierr);
40657795646SDave May 
4076845f8f5SDave May   ierr = DataBucketCreate(&swarm->db);CHKERRQ(ierr);
408b5bcf523SDave May   swarm->vec_field_set = PETSC_FALSE;
409b62e03f8SDave May 
41057795646SDave May   dm->dim  = 0;
41157795646SDave May   dm->data = swarm;
41257795646SDave May 
4135f50eb2eSDave May   dm->ops->view                            = DMView_Swarm;
41457795646SDave May   dm->ops->load                            = NULL;
41557795646SDave May   dm->ops->setfromoptions                  = NULL;
41657795646SDave May   dm->ops->clone                           = NULL;
41757795646SDave May   dm->ops->setup                           = NULL;
41857795646SDave May   dm->ops->createdefaultsection            = NULL;
41957795646SDave May   dm->ops->createdefaultconstraints        = NULL;
420b5bcf523SDave May   dm->ops->createglobalvector              = DMCreateGlobalVector_Swarm;
421b5bcf523SDave May   dm->ops->createlocalvector               = DMCreateLocalVector_Swarm;
42257795646SDave May   dm->ops->getlocaltoglobalmapping         = NULL;
42357795646SDave May   dm->ops->createfieldis                   = NULL;
42457795646SDave May   dm->ops->createcoordinatedm              = NULL;
42557795646SDave May   dm->ops->getcoloring                     = NULL;
42657795646SDave May   dm->ops->creatematrix                    = NULL;
42757795646SDave May   dm->ops->createinterpolation             = NULL;
42857795646SDave May   dm->ops->getaggregates                   = NULL;
42957795646SDave May   dm->ops->getinjection                    = NULL;
43057795646SDave May   dm->ops->refine                          = NULL;
43157795646SDave May   dm->ops->coarsen                         = NULL;
43257795646SDave May   dm->ops->refinehierarchy                 = NULL;
43357795646SDave May   dm->ops->coarsenhierarchy                = NULL;
43457795646SDave May   dm->ops->globaltolocalbegin              = NULL;
43557795646SDave May   dm->ops->globaltolocalend                = NULL;
43657795646SDave May   dm->ops->localtoglobalbegin              = NULL;
43757795646SDave May   dm->ops->localtoglobalend                = NULL;
43857795646SDave May   dm->ops->destroy                         = DMDestroy_Swarm;
43957795646SDave May   dm->ops->createsubdm                     = NULL;
44057795646SDave May   dm->ops->getdimpoints                    = NULL;
44157795646SDave May   dm->ops->locatepoints                    = NULL;
44257795646SDave May 
44357795646SDave May   PetscFunctionReturn(0);
44457795646SDave May }