xref: /petsc/src/dm/impls/swarm/swarm.c (revision b62e03f8c0ccd15bcac32967d2e4773ecbab4f95)
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__ "DMSwarmRegisterPetscDatatypeField"
63 PetscErrorCode DMSwarmRegisterPetscDatatypeField(DM dm,const char fieldname[],PetscInt blocksize,PetscDataType type)
64 {
65   DM_Swarm *swarm = (DM_Swarm*)dm->data;
66   size_t size;
67 
68   if (type == PETSC_OBJECT) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
69   if (type == PETSC_FUNCTION) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
70   if (type == PETSC_STRING) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
71   if (type == PETSC_STRUCT) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
72   if (type == PETSC_DATATYPE_UNKNOWN) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
73 
74   switch (type) {
75     case PETSC_CHAR:
76       size = sizeof(PetscChar);
77       break;
78     case PETSC_SHORT:
79       size = sizeof(PetscShort);
80       break;
81     case PETSC_INT:
82       size = sizeof(PetscInt);
83       break;
84     case PETSC_LONG:
85       size = sizeof(Petsc64bitInt);
86       break;
87     case PETSC_FLOAT:
88       size = sizeof(PetscFloat);
89       break;
90     case PETSC_DOUBLE:
91       size = sizeof(PetscReal);
92       break;
93 
94     default:
95       SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Valid for {char,short,int,long,float,double}");
96       break;
97   }
98 
99   /* Load a specific data type into data bucket, specifying textual name and its size in bytes */
100 	DataBucketRegisterField(swarm->db,fieldname,size,NULL);
101   swarm->db->field[swarm->db->nfields-1]->petsc_type = type;
102 
103   PetscFunctionReturn(0);
104 }
105 
106 #undef __FUNCT__
107 #define __FUNCT__ "DMSwarmRegisterUserStructField"
108 PetscErrorCode DMSwarmRegisterUserStructField(DM dm,const char fieldname[],size_t size)
109 {
110   DM_Swarm *swarm = (DM_Swarm*)dm->data;
111 
112 	DataBucketRegisterField(swarm->db,fieldname,size,NULL);
113   swarm->db->field[swarm->db->nfields-1]->petsc_type = PETSC_STRUCT ;
114 
115   PetscFunctionReturn(0);
116 }
117 
118 #undef __FUNCT__
119 #define __FUNCT__ "DMSwarmRegisterUserDatatypeField"
120 PetscErrorCode DMSwarmRegisterUserDatatypeField(DM dm,const char fieldname[],size_t size)
121 {
122   DM_Swarm *swarm = (DM_Swarm*)dm->data;
123 
124 	DataBucketRegisterField(swarm->db,fieldname,size,NULL);
125   swarm->db->field[swarm->db->nfields-1]->petsc_type = PETSC_DATATYPE_UNKNOWN;
126 
127   PetscFunctionReturn(0);
128 }
129 
130 #undef __FUNCT__
131 #define __FUNCT__ "DMSwarmGetField"
132 PetscErrorCode DMSwarmGetField(DM dm,const char fieldname[],PetscInt *blocksize,PetscDataType *type,void **data)
133 {
134   DM_Swarm *swarm = (DM_Swarm*)dm->data;
135   DataField gfield;
136 
137   DataBucketGetDataFieldByName(swarm->db,fieldname,&gfield);
138   DataFieldGetAccess(gfield);
139   DataFieldGetEntries(gfield,data);
140 
141   PetscFunctionReturn(0);
142 }
143 
144 #undef __FUNCT__
145 #define __FUNCT__ "DMSwarmRestoreField"
146 PetscErrorCode DMSwarmRestoreField(DM dm,const char fieldname[],PetscInt *blocksize,PetscDataType *type,void **data)
147 {
148   DM_Swarm *swarm = (DM_Swarm*)dm->data;
149   DataField gfield;
150 
151   DataBucketGetDataFieldByName(swarm->db,fieldname,&gfield);
152   DataFieldRestoreAccess(gfield);
153   if (data) *data = NULL;
154 
155   PetscFunctionReturn(0);
156 }
157 
158 
159 #undef __FUNCT__
160 #define __FUNCT__ "DMDestroy_Swarm"
161 PetscErrorCode DMDestroy_Swarm(DM dm)
162 {
163   DM_Swarm *swarm = (DM_Swarm*)dm->data;
164   PetscErrorCode ierr;
165 
166   PetscFunctionBegin;
167   DataBucketDestroy(&swarm->db);
168   ierr = PetscFree(swarm);CHKERRQ(ierr);
169   PetscFunctionReturn(0);
170 }
171 
172 #undef __FUNCT__
173 #define __FUNCT__ "DMCreate_Swarm"
174 PETSC_EXTERN PetscErrorCode DMCreate_Swarm(DM dm)
175 {
176   DM_Swarm      *swarm;
177   PetscErrorCode ierr;
178 
179   PetscFunctionBegin;
180   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
181   ierr     = PetscNewLog(dm,&swarm);CHKERRQ(ierr);
182 
183   DataBucketCreate(&swarm->db);
184 
185   dm->dim  = 0;
186   dm->data = swarm;
187 
188   dm->ops->view                            = NULL;
189   dm->ops->load                            = NULL;
190   dm->ops->setfromoptions                  = NULL;
191   dm->ops->clone                           = NULL;
192   dm->ops->setup                           = NULL;
193   dm->ops->createdefaultsection            = NULL;
194   dm->ops->createdefaultconstraints        = NULL;
195   dm->ops->createglobalvector              = NULL;
196   dm->ops->createlocalvector               = NULL;
197   dm->ops->getlocaltoglobalmapping         = NULL;
198   dm->ops->createfieldis                   = NULL;
199   dm->ops->createcoordinatedm              = NULL;
200   dm->ops->getcoloring                     = NULL;
201   dm->ops->creatematrix                    = NULL;
202   dm->ops->createinterpolation             = NULL;
203   dm->ops->getaggregates                   = NULL;
204   dm->ops->getinjection                    = NULL;
205   dm->ops->refine                          = NULL;
206   dm->ops->coarsen                         = NULL;
207   dm->ops->refinehierarchy                 = NULL;
208   dm->ops->coarsenhierarchy                = NULL;
209   dm->ops->globaltolocalbegin              = NULL;
210   dm->ops->globaltolocalend                = NULL;
211   dm->ops->localtoglobalbegin              = NULL;
212   dm->ops->localtoglobalend                = NULL;
213   dm->ops->destroy                         = DMDestroy_Swarm;
214   dm->ops->createsubdm                     = NULL;
215   dm->ops->getdimpoints                    = NULL;
216   dm->ops->locatepoints                    = NULL;
217 
218   PetscFunctionReturn(0);
219 }