xref: /petsc/src/dm/impls/swarm/swarm.c (revision 5f50eb2efa956e5fc10c906e884b1cd686a330cb)
1 
2 #define PETSCDM_DLL
3 #include <petsc/private/dmswarmimpl.h>    /*I   "petscdmswarm.h"   I*/
4 #include "data_bucket.h"
5 
6 //typedef PetscErrorCode (*swarm_project)(DM,DM,Vec) DMSwarmProjectMethod; /* swarm, geometry, result */
7 
8 //typedef enum { PROJECT_DMDA_AQ1=0, PROJECT_DMDA_P0 } DMSwarmDMDAProjectionType;
9 
10 #if 0
11 
12 PetscErrorCode DMSwarmDefineFieldVector(DM dm,const char *fieldnames[])
13 {
14   /* Check all fields are of type PETSC_REAL or PETSC_SCALAR */
15   /* Compute summed block size */
16   /* Set guard */
17   PetscFunctionReturn(0);
18 }
19 
20 PetscErrorCode DMSwarmGetGlobalVectorFromFields(DM dm,Vec *vec)
21 {
22   PetscFunctionReturn(0);
23 }
24 
25 PetscErrorCode DMSwarmRestoreGlobalVectorFromFields(DM dm,Vec *vec)
26 {
27   PetscFunctionReturn(0);
28 }
29 
30 /* requires DMSwarmDefineFieldVector has been called */
31 PetscErrorCode DMCreateGlobalVector_Swarm(void)
32 {
33 
34   PetscFunctionReturn(0);
35 }
36 
37 /* requires DMSwarmDefineFieldVector has been called */
38 PetscErrorCode DMLocalGlobalVector_Swarm(void)
39 {
40 
41   PetscFunctionReturn(0);
42 }
43 
44 /* Defines what the local space will be */
45 PetscErrorCode DMSwarmSetOverlap(void)
46 {
47 
48   PetscFunctionReturn(0);
49 }
50 
51 
52 /* coordinates */
53 /*
54 DMGetCoordinateDM returns self
55 DMGetCoordinates and DMGetCoordinatesLocal return same thing
56 Local view could be used to define overlapping information
57 */
58 
59 #endif
60 
61 #undef __FUNCT__
62 #define __FUNCT__ "DMSwarmInitializeFieldRegister"
63 PETSC_EXTERN PetscErrorCode DMSwarmInitializeFieldRegister(DM dm)
64 {
65   DM_Swarm *swarm = (DM_Swarm*)dm->data;
66   swarm->field_registration_initialized = PETSC_TRUE;
67   PetscFunctionReturn(0);
68 }
69 
70 #undef __FUNCT__
71 #define __FUNCT__ "DMSwarmFinalizeFieldRegister"
72 PETSC_EXTERN PetscErrorCode DMSwarmFinalizeFieldRegister(DM dm)
73 {
74   DM_Swarm *swarm = (DM_Swarm*)dm->data;
75   swarm->field_registration_finalized = PETSC_TRUE;
76   DataBucketFinalize(swarm->db);
77   PetscFunctionReturn(0);
78 }
79 
80 #undef __FUNCT__
81 #define __FUNCT__ "DMSwarmSetLocalSizes"
82 PETSC_EXTERN PetscErrorCode DMSwarmSetLocalSizes(DM dm,PetscInt nlocal,PetscInt buffer)
83 {
84   DM_Swarm *swarm = (DM_Swarm*)dm->data;
85 
86   DataBucketSetSizes(swarm->db,nlocal,buffer);
87 
88   PetscFunctionReturn(0);
89 }
90 
91 #undef __FUNCT__
92 #define __FUNCT__ "DMSwarmRegisterPetscDatatypeField"
93 PETSC_EXTERN PetscErrorCode DMSwarmRegisterPetscDatatypeField(DM dm,const char fieldname[],PetscInt blocksize,PetscDataType type)
94 {
95   DM_Swarm *swarm = (DM_Swarm*)dm->data;
96   size_t size;
97 
98   if (!swarm->field_registration_initialized) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMSwarmInitializeFieldRegister() first");
99   if (swarm->field_registration_finalized) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Cannot register additional fields after calling DMSwarmFinalizeFieldRegister() first");
100 
101   if (type == PETSC_OBJECT) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
102   if (type == PETSC_FUNCTION) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
103   if (type == PETSC_STRING) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
104   if (type == PETSC_STRUCT) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
105   if (type == PETSC_DATATYPE_UNKNOWN) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
106 
107   switch (type) {
108     case PETSC_CHAR:
109       size = sizeof(PetscChar);
110       break;
111     case PETSC_SHORT:
112       size = sizeof(PetscShort);
113       break;
114     case PETSC_INT:
115       size = sizeof(PetscInt);
116       break;
117     case PETSC_LONG:
118       size = sizeof(Petsc64bitInt);
119       break;
120     case PETSC_FLOAT:
121       size = sizeof(PetscFloat);
122       break;
123     case PETSC_DOUBLE:
124       size = sizeof(PetscReal);
125       break;
126 
127     default:
128       SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
129       break;
130   }
131 
132   /* Load a specific data type into data bucket, specifying textual name and its size in bytes */
133 	DataBucketRegisterField(swarm->db,fieldname,size,NULL);
134   swarm->db->field[swarm->db->nfields-1]->petsc_type = type;
135 
136   PetscFunctionReturn(0);
137 }
138 
139 #undef __FUNCT__
140 #define __FUNCT__ "DMSwarmRegisterUserStructField"
141 PETSC_EXTERN PetscErrorCode DMSwarmRegisterUserStructField(DM dm,const char fieldname[],size_t size)
142 {
143   DM_Swarm *swarm = (DM_Swarm*)dm->data;
144 
145 	DataBucketRegisterField(swarm->db,fieldname,size,NULL);
146   swarm->db->field[swarm->db->nfields-1]->petsc_type = PETSC_STRUCT ;
147 
148   PetscFunctionReturn(0);
149 }
150 
151 #undef __FUNCT__
152 #define __FUNCT__ "DMSwarmRegisterUserDatatypeField"
153 PETSC_EXTERN PetscErrorCode DMSwarmRegisterUserDatatypeField(DM dm,const char fieldname[],size_t size)
154 {
155   DM_Swarm *swarm = (DM_Swarm*)dm->data;
156 
157 	DataBucketRegisterField(swarm->db,fieldname,size,NULL);
158   swarm->db->field[swarm->db->nfields-1]->petsc_type = PETSC_DATATYPE_UNKNOWN;
159 
160   PetscFunctionReturn(0);
161 }
162 
163 #undef __FUNCT__
164 #define __FUNCT__ "DMSwarmGetField"
165 PETSC_EXTERN PetscErrorCode DMSwarmGetField(DM dm,const char fieldname[],PetscInt *blocksize,PetscDataType *type,void **data)
166 {
167   DM_Swarm *swarm = (DM_Swarm*)dm->data;
168   DataField gfield;
169 
170   DataBucketGetDataFieldByName(swarm->db,fieldname,&gfield);
171   DataFieldGetAccess(gfield);
172   DataFieldGetEntries(gfield,data);
173 
174   PetscFunctionReturn(0);
175 }
176 
177 #undef __FUNCT__
178 #define __FUNCT__ "DMSwarmRestoreField"
179 PETSC_EXTERN PetscErrorCode DMSwarmRestoreField(DM dm,const char fieldname[],PetscInt *blocksize,PetscDataType *type,void **data)
180 {
181   DM_Swarm *swarm = (DM_Swarm*)dm->data;
182   DataField gfield;
183 
184   DataBucketGetDataFieldByName(swarm->db,fieldname,&gfield);
185   DataFieldRestoreAccess(gfield);
186   if (data) *data = NULL;
187 
188   PetscFunctionReturn(0);
189 }
190 
191 
192 #undef __FUNCT__
193 #define __FUNCT__ "DMDestroy_Swarm"
194 PetscErrorCode DMDestroy_Swarm(DM dm)
195 {
196   DM_Swarm *swarm = (DM_Swarm*)dm->data;
197   PetscErrorCode ierr;
198 
199   PetscFunctionBegin;
200   DataBucketDestroy(&swarm->db);
201   ierr = PetscFree(swarm);CHKERRQ(ierr);
202   PetscFunctionReturn(0);
203 }
204 
205 #undef __FUNCT__
206 #define __FUNCT__ "DMView_Swarm"
207 PetscErrorCode DMView_Swarm(DM dm, PetscViewer viewer)
208 {
209   DM_Swarm *swarm = (DM_Swarm*)dm->data;
210   PetscBool      iascii,ibinary,ishdf5,isvtk;
211   PetscErrorCode ierr;
212 
213   PetscFunctionBegin;
214   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
215   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
216   ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
217   ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY,&ibinary);CHKERRQ(ierr);
218   ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERVTK,   &isvtk);CHKERRQ(ierr);
219   ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5,  &ishdf5);CHKERRQ(ierr);
220   if (iascii) {
221     DataBucketView(PetscObjectComm((PetscObject)dm),swarm->db,NULL,DATABUCKET_VIEW_STDOUT);
222   } else if (ibinary) {
223     SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"NO VTK support");
224   } else if (ishdf5) {
225 #if defined(PETSC_HAVE_HDF5)
226     SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"NO HDF5 support");
227 #else
228     SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"HDF5 not supported. Please reconfigure using --download-hdf5");
229 #endif
230   } else if (isvtk) {
231     SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"NO VTK support");
232   }
233   PetscFunctionReturn(0);
234 }
235 
236 #undef __FUNCT__
237 #define __FUNCT__ "DMCreate_Swarm"
238 PETSC_EXTERN PetscErrorCode DMCreate_Swarm(DM dm)
239 {
240   DM_Swarm      *swarm;
241   PetscErrorCode ierr;
242 
243   PetscFunctionBegin;
244   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
245   ierr     = PetscNewLog(dm,&swarm);CHKERRQ(ierr);
246 
247   DataBucketCreate(&swarm->db);
248 
249   dm->dim  = 0;
250   dm->data = swarm;
251 
252   dm->ops->view                            = DMView_Swarm;
253   dm->ops->load                            = NULL;
254   dm->ops->setfromoptions                  = NULL;
255   dm->ops->clone                           = NULL;
256   dm->ops->setup                           = NULL;
257   dm->ops->createdefaultsection            = NULL;
258   dm->ops->createdefaultconstraints        = NULL;
259   dm->ops->createglobalvector              = NULL;
260   dm->ops->createlocalvector               = NULL;
261   dm->ops->getlocaltoglobalmapping         = NULL;
262   dm->ops->createfieldis                   = NULL;
263   dm->ops->createcoordinatedm              = NULL;
264   dm->ops->getcoloring                     = NULL;
265   dm->ops->creatematrix                    = NULL;
266   dm->ops->createinterpolation             = NULL;
267   dm->ops->getaggregates                   = NULL;
268   dm->ops->getinjection                    = NULL;
269   dm->ops->refine                          = NULL;
270   dm->ops->coarsen                         = NULL;
271   dm->ops->refinehierarchy                 = NULL;
272   dm->ops->coarsenhierarchy                = NULL;
273   dm->ops->globaltolocalbegin              = NULL;
274   dm->ops->globaltolocalend                = NULL;
275   dm->ops->localtoglobalbegin              = NULL;
276   dm->ops->localtoglobalend                = NULL;
277   dm->ops->destroy                         = DMDestroy_Swarm;
278   dm->ops->createsubdm                     = NULL;
279   dm->ops->getdimpoints                    = NULL;
280   dm->ops->locatepoints                    = NULL;
281 
282   PetscFunctionReturn(0);
283 }