xref: /petsc/src/dm/impls/swarm/swarm.c (revision b62e03f8c0ccd15bcac32967d2e4773ecbab4f95)
157795646SDave May 
257795646SDave May #define PETSCDM_DLL
357795646SDave May #include <petsc/private/dmswarmimpl.h>    /*I   "petscdmswarm.h"   I*/
4*b62e03f8SDave 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 }
24*b62e03f8SDave 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 */
31*b62e03f8SDave May PetscErrorCode DMCreateGlobalVector_Swarm(void)
3257795646SDave May {
3357795646SDave May 
3457795646SDave May   PetscFunctionReturn(0);
3557795646SDave May }
3657795646SDave May 
3757795646SDave May /* requires DMSwarmDefineFieldVector has been called */
38*b62e03f8SDave May PetscErrorCode DMLocalGlobalVector_Swarm(void)
3957795646SDave May {
4057795646SDave May 
41*b62e03f8SDave May   PetscFunctionReturn(0);
4257795646SDave May }
4357795646SDave May 
4457795646SDave May /* Defines what the local space will be */
45*b62e03f8SDave May PetscErrorCode DMSwarmSetOverlap(void)
4657795646SDave May {
4757795646SDave May 
48*b62e03f8SDave 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*b62e03f8SDave May #define __FUNCT__ "DMSwarmRegisterPetscDatatypeField"
63*b62e03f8SDave May PetscErrorCode DMSwarmRegisterPetscDatatypeField(DM dm,const char fieldname[],PetscInt blocksize,PetscDataType type)
64*b62e03f8SDave May {
65*b62e03f8SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
66*b62e03f8SDave May   size_t size;
67*b62e03f8SDave May 
68*b62e03f8SDave May   if (type == PETSC_OBJECT) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
69*b62e03f8SDave May   if (type == PETSC_FUNCTION) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
70*b62e03f8SDave May   if (type == PETSC_STRING) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
71*b62e03f8SDave May   if (type == PETSC_STRUCT) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
72*b62e03f8SDave May   if (type == PETSC_DATATYPE_UNKNOWN) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
73*b62e03f8SDave May 
74*b62e03f8SDave May   switch (type) {
75*b62e03f8SDave May     case PETSC_CHAR:
76*b62e03f8SDave May       size = sizeof(PetscChar);
77*b62e03f8SDave May       break;
78*b62e03f8SDave May     case PETSC_SHORT:
79*b62e03f8SDave May       size = sizeof(PetscShort);
80*b62e03f8SDave May       break;
81*b62e03f8SDave May     case PETSC_INT:
82*b62e03f8SDave May       size = sizeof(PetscInt);
83*b62e03f8SDave May       break;
84*b62e03f8SDave May     case PETSC_LONG:
85*b62e03f8SDave May       size = sizeof(Petsc64bitInt);
86*b62e03f8SDave May       break;
87*b62e03f8SDave May     case PETSC_FLOAT:
88*b62e03f8SDave May       size = sizeof(PetscFloat);
89*b62e03f8SDave May       break;
90*b62e03f8SDave May     case PETSC_DOUBLE:
91*b62e03f8SDave May       size = sizeof(PetscReal);
92*b62e03f8SDave May       break;
93*b62e03f8SDave May 
94*b62e03f8SDave May     default:
95*b62e03f8SDave May       SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
96*b62e03f8SDave May       break;
97*b62e03f8SDave May   }
98*b62e03f8SDave May 
99*b62e03f8SDave May   /* Load a specific data type into data bucket, specifying textual name and its size in bytes */
100*b62e03f8SDave May 	DataBucketRegisterField(swarm->db,fieldname,size,NULL);
101*b62e03f8SDave May   swarm->db->field[swarm->db->nfields-1]->petsc_type = type;
102*b62e03f8SDave May 
103*b62e03f8SDave May   PetscFunctionReturn(0);
104*b62e03f8SDave May }
105*b62e03f8SDave May 
106*b62e03f8SDave May #undef __FUNCT__
107*b62e03f8SDave May #define __FUNCT__ "DMSwarmRegisterUserStructField"
108*b62e03f8SDave May PetscErrorCode DMSwarmRegisterUserStructField(DM dm,const char fieldname[],size_t size)
109*b62e03f8SDave May {
110*b62e03f8SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
111*b62e03f8SDave May 
112*b62e03f8SDave May 	DataBucketRegisterField(swarm->db,fieldname,size,NULL);
113*b62e03f8SDave May   swarm->db->field[swarm->db->nfields-1]->petsc_type = PETSC_STRUCT ;
114*b62e03f8SDave May 
115*b62e03f8SDave May   PetscFunctionReturn(0);
116*b62e03f8SDave May }
117*b62e03f8SDave May 
118*b62e03f8SDave May #undef __FUNCT__
119*b62e03f8SDave May #define __FUNCT__ "DMSwarmRegisterUserDatatypeField"
120*b62e03f8SDave May PetscErrorCode DMSwarmRegisterUserDatatypeField(DM dm,const char fieldname[],size_t size)
121*b62e03f8SDave May {
122*b62e03f8SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
123*b62e03f8SDave May 
124*b62e03f8SDave May 	DataBucketRegisterField(swarm->db,fieldname,size,NULL);
125*b62e03f8SDave May   swarm->db->field[swarm->db->nfields-1]->petsc_type = PETSC_DATATYPE_UNKNOWN;
126*b62e03f8SDave May 
127*b62e03f8SDave May   PetscFunctionReturn(0);
128*b62e03f8SDave May }
129*b62e03f8SDave May 
130*b62e03f8SDave May #undef __FUNCT__
131*b62e03f8SDave May #define __FUNCT__ "DMSwarmGetField"
132*b62e03f8SDave May PetscErrorCode DMSwarmGetField(DM dm,const char fieldname[],PetscInt *blocksize,PetscDataType *type,void **data)
133*b62e03f8SDave May {
134*b62e03f8SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
135*b62e03f8SDave May   DataField gfield;
136*b62e03f8SDave May 
137*b62e03f8SDave May   DataBucketGetDataFieldByName(swarm->db,fieldname,&gfield);
138*b62e03f8SDave May   DataFieldGetAccess(gfield);
139*b62e03f8SDave May   DataFieldGetEntries(gfield,data);
140*b62e03f8SDave May 
141*b62e03f8SDave May   PetscFunctionReturn(0);
142*b62e03f8SDave May }
143*b62e03f8SDave May 
144*b62e03f8SDave May #undef __FUNCT__
145*b62e03f8SDave May #define __FUNCT__ "DMSwarmRestoreField"
146*b62e03f8SDave May PetscErrorCode DMSwarmRestoreField(DM dm,const char fieldname[],PetscInt *blocksize,PetscDataType *type,void **data)
147*b62e03f8SDave May {
148*b62e03f8SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
149*b62e03f8SDave May   DataField gfield;
150*b62e03f8SDave May 
151*b62e03f8SDave May   DataBucketGetDataFieldByName(swarm->db,fieldname,&gfield);
152*b62e03f8SDave May   DataFieldRestoreAccess(gfield);
153*b62e03f8SDave May   if (data) *data = NULL;
154*b62e03f8SDave May 
155*b62e03f8SDave May   PetscFunctionReturn(0);
156*b62e03f8SDave May }
157*b62e03f8SDave May 
158*b62e03f8SDave May 
159*b62e03f8SDave May #undef __FUNCT__
16057795646SDave May #define __FUNCT__ "DMDestroy_Swarm"
16157795646SDave May PetscErrorCode DMDestroy_Swarm(DM dm)
16257795646SDave May {
16357795646SDave May   DM_Swarm *swarm = (DM_Swarm*)dm->data;
16457795646SDave May   PetscErrorCode ierr;
16557795646SDave May 
16657795646SDave May   PetscFunctionBegin;
167*b62e03f8SDave May   DataBucketDestroy(&swarm->db);
16857795646SDave May   ierr = PetscFree(swarm);CHKERRQ(ierr);
16957795646SDave May   PetscFunctionReturn(0);
17057795646SDave May }
17157795646SDave May 
17257795646SDave May #undef __FUNCT__
17357795646SDave May #define __FUNCT__ "DMCreate_Swarm"
17457795646SDave May PETSC_EXTERN PetscErrorCode DMCreate_Swarm(DM dm)
17557795646SDave May {
17657795646SDave May   DM_Swarm      *swarm;
17757795646SDave May   PetscErrorCode ierr;
17857795646SDave May 
17957795646SDave May   PetscFunctionBegin;
18057795646SDave May   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
18157795646SDave May   ierr     = PetscNewLog(dm,&swarm);CHKERRQ(ierr);
18257795646SDave May 
183*b62e03f8SDave May   DataBucketCreate(&swarm->db);
184*b62e03f8SDave May 
18557795646SDave May   dm->dim  = 0;
18657795646SDave May   dm->data = swarm;
18757795646SDave May 
18857795646SDave May   dm->ops->view                            = NULL;
18957795646SDave May   dm->ops->load                            = NULL;
19057795646SDave May   dm->ops->setfromoptions                  = NULL;
19157795646SDave May   dm->ops->clone                           = NULL;
19257795646SDave May   dm->ops->setup                           = NULL;
19357795646SDave May   dm->ops->createdefaultsection            = NULL;
19457795646SDave May   dm->ops->createdefaultconstraints        = NULL;
19557795646SDave May   dm->ops->createglobalvector              = NULL;
19657795646SDave May   dm->ops->createlocalvector               = NULL;
19757795646SDave May   dm->ops->getlocaltoglobalmapping         = NULL;
19857795646SDave May   dm->ops->createfieldis                   = NULL;
19957795646SDave May   dm->ops->createcoordinatedm              = NULL;
20057795646SDave May   dm->ops->getcoloring                     = NULL;
20157795646SDave May   dm->ops->creatematrix                    = NULL;
20257795646SDave May   dm->ops->createinterpolation             = NULL;
20357795646SDave May   dm->ops->getaggregates                   = NULL;
20457795646SDave May   dm->ops->getinjection                    = NULL;
20557795646SDave May   dm->ops->refine                          = NULL;
20657795646SDave May   dm->ops->coarsen                         = NULL;
20757795646SDave May   dm->ops->refinehierarchy                 = NULL;
20857795646SDave May   dm->ops->coarsenhierarchy                = NULL;
20957795646SDave May   dm->ops->globaltolocalbegin              = NULL;
21057795646SDave May   dm->ops->globaltolocalend                = NULL;
21157795646SDave May   dm->ops->localtoglobalbegin              = NULL;
21257795646SDave May   dm->ops->localtoglobalend                = NULL;
21357795646SDave May   dm->ops->destroy                         = DMDestroy_Swarm;
21457795646SDave May   dm->ops->createsubdm                     = NULL;
21557795646SDave May   dm->ops->getdimpoints                    = NULL;
21657795646SDave May   dm->ops->locatepoints                    = NULL;
21757795646SDave May 
21857795646SDave May   PetscFunctionReturn(0);
21957795646SDave May }