xref: /petsc/src/dm/impls/swarm/swarm.c (revision 5f50eb2efa956e5fc10c906e884b1cd686a330cb)
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 PetscErrorCode DMSwarmDefineFieldVector(DM dm,const char *fieldnames[])
1357795646SDave May {
1457795646SDave May   /* Check all fields are of type PETSC_REAL or PETSC_SCALAR */
1557795646SDave May   /* Compute summed block size */
1657795646SDave May   /* Set guard */
1757795646SDave May   PetscFunctionReturn(0);
1857795646SDave May }
1957795646SDave May 
2057795646SDave May PetscErrorCode DMSwarmGetGlobalVectorFromFields(DM dm,Vec *vec)
2157795646SDave May {
2257795646SDave May   PetscFunctionReturn(0);
2357795646SDave May }
24b62e03f8SDave May 
2557795646SDave May PetscErrorCode DMSwarmRestoreGlobalVectorFromFields(DM dm,Vec *vec)
2657795646SDave May {
2757795646SDave May   PetscFunctionReturn(0);
2857795646SDave May }
2957795646SDave May 
3057795646SDave May /* requires DMSwarmDefineFieldVector has been called */
31b62e03f8SDave May PetscErrorCode DMCreateGlobalVector_Swarm(void)
3257795646SDave May {
3357795646SDave May 
3457795646SDave May   PetscFunctionReturn(0);
3557795646SDave May }
3657795646SDave May 
3757795646SDave May /* requires DMSwarmDefineFieldVector has been called */
38b62e03f8SDave May PetscErrorCode DMLocalGlobalVector_Swarm(void)
3957795646SDave May {
4057795646SDave May 
41b62e03f8SDave May   PetscFunctionReturn(0);
4257795646SDave May }
4357795646SDave May 
4457795646SDave May /* Defines what the local space will be */
45b62e03f8SDave May PetscErrorCode DMSwarmSetOverlap(void)
4657795646SDave May {
4757795646SDave May 
48b62e03f8SDave May   PetscFunctionReturn(0);
4957795646SDave May }
5057795646SDave May 
5157795646SDave May 
5257795646SDave May /* coordinates */
5357795646SDave May /*
5457795646SDave May DMGetCoordinateDM returns self
5557795646SDave May DMGetCoordinates and DMGetCoordinatesLocal return same thing
5657795646SDave May Local view could be used to define overlapping information
5757795646SDave May */
5857795646SDave May 
5957795646SDave May #endif
6057795646SDave May 
6157795646SDave May #undef __FUNCT__
62*5f50eb2eSDave May #define __FUNCT__ "DMSwarmInitializeFieldRegister"
63*5f50eb2eSDave May PETSC_EXTERN PetscErrorCode DMSwarmInitializeFieldRegister(DM dm)
64*5f50eb2eSDave May {
65*5f50eb2eSDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
66*5f50eb2eSDave May   swarm->field_registration_initialized = PETSC_TRUE;
67*5f50eb2eSDave May   PetscFunctionReturn(0);
68*5f50eb2eSDave May }
69*5f50eb2eSDave May 
70*5f50eb2eSDave May #undef __FUNCT__
71*5f50eb2eSDave May #define __FUNCT__ "DMSwarmFinalizeFieldRegister"
72*5f50eb2eSDave May PETSC_EXTERN PetscErrorCode DMSwarmFinalizeFieldRegister(DM dm)
73*5f50eb2eSDave May {
74*5f50eb2eSDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
75*5f50eb2eSDave May   swarm->field_registration_finalized = PETSC_TRUE;
76*5f50eb2eSDave May   DataBucketFinalize(swarm->db);
77*5f50eb2eSDave May   PetscFunctionReturn(0);
78*5f50eb2eSDave May }
79*5f50eb2eSDave May 
80*5f50eb2eSDave May #undef __FUNCT__
81*5f50eb2eSDave May #define __FUNCT__ "DMSwarmSetLocalSizes"
82*5f50eb2eSDave May PETSC_EXTERN PetscErrorCode DMSwarmSetLocalSizes(DM dm,PetscInt nlocal,PetscInt buffer)
83*5f50eb2eSDave May {
84*5f50eb2eSDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
85*5f50eb2eSDave May 
86*5f50eb2eSDave May   DataBucketSetSizes(swarm->db,nlocal,buffer);
87*5f50eb2eSDave May 
88*5f50eb2eSDave May   PetscFunctionReturn(0);
89*5f50eb2eSDave May }
90*5f50eb2eSDave May 
91*5f50eb2eSDave May #undef __FUNCT__
92b62e03f8SDave May #define __FUNCT__ "DMSwarmRegisterPetscDatatypeField"
93*5f50eb2eSDave May PETSC_EXTERN PetscErrorCode DMSwarmRegisterPetscDatatypeField(DM dm,const char fieldname[],PetscInt blocksize,PetscDataType type)
94b62e03f8SDave May {
95b62e03f8SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
96b62e03f8SDave May   size_t size;
97b62e03f8SDave May 
98*5f50eb2eSDave May   if (!swarm->field_registration_initialized) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Must call DMSwarmInitializeFieldRegister() first");
99*5f50eb2eSDave May   if (swarm->field_registration_finalized) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Cannot register additional fields after calling DMSwarmFinalizeFieldRegister() first");
100*5f50eb2eSDave May 
101*5f50eb2eSDave May   if (type == PETSC_OBJECT) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
102*5f50eb2eSDave May   if (type == PETSC_FUNCTION) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
103*5f50eb2eSDave May   if (type == PETSC_STRING) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
104*5f50eb2eSDave May   if (type == PETSC_STRUCT) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
105*5f50eb2eSDave May   if (type == PETSC_DATATYPE_UNKNOWN) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
106b62e03f8SDave May 
107b62e03f8SDave May   switch (type) {
108b62e03f8SDave May     case PETSC_CHAR:
109b62e03f8SDave May       size = sizeof(PetscChar);
110b62e03f8SDave May       break;
111b62e03f8SDave May     case PETSC_SHORT:
112b62e03f8SDave May       size = sizeof(PetscShort);
113b62e03f8SDave May       break;
114b62e03f8SDave May     case PETSC_INT:
115b62e03f8SDave May       size = sizeof(PetscInt);
116b62e03f8SDave May       break;
117b62e03f8SDave May     case PETSC_LONG:
118b62e03f8SDave May       size = sizeof(Petsc64bitInt);
119b62e03f8SDave May       break;
120b62e03f8SDave May     case PETSC_FLOAT:
121b62e03f8SDave May       size = sizeof(PetscFloat);
122b62e03f8SDave May       break;
123b62e03f8SDave May     case PETSC_DOUBLE:
124b62e03f8SDave May       size = sizeof(PetscReal);
125b62e03f8SDave May       break;
126b62e03f8SDave May 
127b62e03f8SDave May     default:
128*5f50eb2eSDave May       SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
129b62e03f8SDave May       break;
130b62e03f8SDave May   }
131b62e03f8SDave May 
132b62e03f8SDave May   /* Load a specific data type into data bucket, specifying textual name and its size in bytes */
133b62e03f8SDave May 	DataBucketRegisterField(swarm->db,fieldname,size,NULL);
134b62e03f8SDave May   swarm->db->field[swarm->db->nfields-1]->petsc_type = type;
135b62e03f8SDave May 
136b62e03f8SDave May   PetscFunctionReturn(0);
137b62e03f8SDave May }
138b62e03f8SDave May 
139b62e03f8SDave May #undef __FUNCT__
140b62e03f8SDave May #define __FUNCT__ "DMSwarmRegisterUserStructField"
141*5f50eb2eSDave May PETSC_EXTERN PetscErrorCode DMSwarmRegisterUserStructField(DM dm,const char fieldname[],size_t size)
142b62e03f8SDave May {
143b62e03f8SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
144b62e03f8SDave May 
145b62e03f8SDave May 	DataBucketRegisterField(swarm->db,fieldname,size,NULL);
146b62e03f8SDave May   swarm->db->field[swarm->db->nfields-1]->petsc_type = PETSC_STRUCT ;
147b62e03f8SDave May 
148b62e03f8SDave May   PetscFunctionReturn(0);
149b62e03f8SDave May }
150b62e03f8SDave May 
151b62e03f8SDave May #undef __FUNCT__
152b62e03f8SDave May #define __FUNCT__ "DMSwarmRegisterUserDatatypeField"
153*5f50eb2eSDave May PETSC_EXTERN PetscErrorCode DMSwarmRegisterUserDatatypeField(DM dm,const char fieldname[],size_t size)
154b62e03f8SDave May {
155b62e03f8SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
156b62e03f8SDave May 
157b62e03f8SDave May 	DataBucketRegisterField(swarm->db,fieldname,size,NULL);
158b62e03f8SDave May   swarm->db->field[swarm->db->nfields-1]->petsc_type = PETSC_DATATYPE_UNKNOWN;
159b62e03f8SDave May 
160b62e03f8SDave May   PetscFunctionReturn(0);
161b62e03f8SDave May }
162b62e03f8SDave May 
163b62e03f8SDave May #undef __FUNCT__
164b62e03f8SDave May #define __FUNCT__ "DMSwarmGetField"
165*5f50eb2eSDave May PETSC_EXTERN PetscErrorCode DMSwarmGetField(DM dm,const char fieldname[],PetscInt *blocksize,PetscDataType *type,void **data)
166b62e03f8SDave May {
167b62e03f8SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
168b62e03f8SDave May   DataField gfield;
169b62e03f8SDave May 
170b62e03f8SDave May   DataBucketGetDataFieldByName(swarm->db,fieldname,&gfield);
171b62e03f8SDave May   DataFieldGetAccess(gfield);
172b62e03f8SDave May   DataFieldGetEntries(gfield,data);
173b62e03f8SDave May 
174b62e03f8SDave May   PetscFunctionReturn(0);
175b62e03f8SDave May }
176b62e03f8SDave May 
177b62e03f8SDave May #undef __FUNCT__
178b62e03f8SDave May #define __FUNCT__ "DMSwarmRestoreField"
179*5f50eb2eSDave May PETSC_EXTERN PetscErrorCode DMSwarmRestoreField(DM dm,const char fieldname[],PetscInt *blocksize,PetscDataType *type,void **data)
180b62e03f8SDave May {
181b62e03f8SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
182b62e03f8SDave May   DataField gfield;
183b62e03f8SDave May 
184b62e03f8SDave May   DataBucketGetDataFieldByName(swarm->db,fieldname,&gfield);
185b62e03f8SDave May   DataFieldRestoreAccess(gfield);
186b62e03f8SDave May   if (data) *data = NULL;
187b62e03f8SDave May 
188b62e03f8SDave May   PetscFunctionReturn(0);
189b62e03f8SDave May }
190b62e03f8SDave May 
191b62e03f8SDave May 
192b62e03f8SDave May #undef __FUNCT__
19357795646SDave May #define __FUNCT__ "DMDestroy_Swarm"
19457795646SDave May PetscErrorCode DMDestroy_Swarm(DM dm)
19557795646SDave May {
19657795646SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
19757795646SDave May   PetscErrorCode ierr;
19857795646SDave May 
19957795646SDave May   PetscFunctionBegin;
200b62e03f8SDave May   DataBucketDestroy(&swarm->db);
20157795646SDave May   ierr = PetscFree(swarm);CHKERRQ(ierr);
20257795646SDave May   PetscFunctionReturn(0);
20357795646SDave May }
20457795646SDave May 
20557795646SDave May #undef __FUNCT__
206*5f50eb2eSDave May #define __FUNCT__ "DMView_Swarm"
207*5f50eb2eSDave May PetscErrorCode DMView_Swarm(DM dm, PetscViewer viewer)
208*5f50eb2eSDave May {
209*5f50eb2eSDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
210*5f50eb2eSDave May   PetscBool      iascii,ibinary,ishdf5,isvtk;
211*5f50eb2eSDave May   PetscErrorCode ierr;
212*5f50eb2eSDave May 
213*5f50eb2eSDave May   PetscFunctionBegin;
214*5f50eb2eSDave May   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
215*5f50eb2eSDave May   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
216*5f50eb2eSDave May   ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
217*5f50eb2eSDave May   ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY,&ibinary);CHKERRQ(ierr);
218*5f50eb2eSDave May   ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERVTK,   &isvtk);CHKERRQ(ierr);
219*5f50eb2eSDave May   ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5,  &ishdf5);CHKERRQ(ierr);
220*5f50eb2eSDave May   if (iascii) {
221*5f50eb2eSDave May     DataBucketView(PetscObjectComm((PetscObject)dm),swarm->db,NULL,DATABUCKET_VIEW_STDOUT);
222*5f50eb2eSDave May   } else if (ibinary) {
223*5f50eb2eSDave May     SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"NO VTK support");
224*5f50eb2eSDave May   } else if (ishdf5) {
225*5f50eb2eSDave May #if defined(PETSC_HAVE_HDF5)
226*5f50eb2eSDave May     SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"NO HDF5 support");
227*5f50eb2eSDave May #else
228*5f50eb2eSDave May     SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"HDF5 not supported. Please reconfigure using --download-hdf5");
229*5f50eb2eSDave May #endif
230*5f50eb2eSDave May   } else if (isvtk) {
231*5f50eb2eSDave May     SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"NO VTK support");
232*5f50eb2eSDave May   }
233*5f50eb2eSDave May   PetscFunctionReturn(0);
234*5f50eb2eSDave May }
235*5f50eb2eSDave May 
236*5f50eb2eSDave May #undef __FUNCT__
23757795646SDave May #define __FUNCT__ "DMCreate_Swarm"
23857795646SDave May PETSC_EXTERN PetscErrorCode DMCreate_Swarm(DM dm)
23957795646SDave May {
24057795646SDave May   DM_Swarm      *swarm;
24157795646SDave May   PetscErrorCode ierr;
24257795646SDave May 
24357795646SDave May   PetscFunctionBegin;
24457795646SDave May   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
24557795646SDave May   ierr     = PetscNewLog(dm,&swarm);CHKERRQ(ierr);
24657795646SDave May 
247b62e03f8SDave May   DataBucketCreate(&swarm->db);
248b62e03f8SDave May 
24957795646SDave May   dm->dim  = 0;
25057795646SDave May   dm->data = swarm;
25157795646SDave May 
252*5f50eb2eSDave May   dm->ops->view                            = DMView_Swarm;
25357795646SDave May   dm->ops->load                            = NULL;
25457795646SDave May   dm->ops->setfromoptions                  = NULL;
25557795646SDave May   dm->ops->clone                           = NULL;
25657795646SDave May   dm->ops->setup                           = NULL;
25757795646SDave May   dm->ops->createdefaultsection            = NULL;
25857795646SDave May   dm->ops->createdefaultconstraints        = NULL;
25957795646SDave May   dm->ops->createglobalvector              = NULL;
26057795646SDave May   dm->ops->createlocalvector               = NULL;
26157795646SDave May   dm->ops->getlocaltoglobalmapping         = NULL;
26257795646SDave May   dm->ops->createfieldis                   = NULL;
26357795646SDave May   dm->ops->createcoordinatedm              = NULL;
26457795646SDave May   dm->ops->getcoloring                     = NULL;
26557795646SDave May   dm->ops->creatematrix                    = NULL;
26657795646SDave May   dm->ops->createinterpolation             = NULL;
26757795646SDave May   dm->ops->getaggregates                   = NULL;
26857795646SDave May   dm->ops->getinjection                    = NULL;
26957795646SDave May   dm->ops->refine                          = NULL;
27057795646SDave May   dm->ops->coarsen                         = NULL;
27157795646SDave May   dm->ops->refinehierarchy                 = NULL;
27257795646SDave May   dm->ops->coarsenhierarchy                = NULL;
27357795646SDave May   dm->ops->globaltolocalbegin              = NULL;
27457795646SDave May   dm->ops->globaltolocalend                = NULL;
27557795646SDave May   dm->ops->localtoglobalbegin              = NULL;
27657795646SDave May   dm->ops->localtoglobalend                = NULL;
27757795646SDave May   dm->ops->destroy                         = DMDestroy_Swarm;
27857795646SDave May   dm->ops->createsubdm                     = NULL;
27957795646SDave May   dm->ops->getdimpoints                    = NULL;
28057795646SDave May   dm->ops->locatepoints                    = NULL;
28157795646SDave May 
28257795646SDave May   PetscFunctionReturn(0);
28357795646SDave May }