xref: /petsc/src/dm/impls/da/dacreate.c (revision a790c00499656cdc34e95ec8f15f35cf5f4cdcbb)
1 
2 #include <petsc-private/daimpl.h>    /*I   "petscdmda.h"   I*/
3 
4 #undef __FUNCT__
5 #define __FUNCT__ "DMSetFromOptions_DA"
6 PetscErrorCode  DMSetFromOptions_DA(DM da)
7 {
8   PetscErrorCode ierr;
9   DM_DA          *dd = (DM_DA*)da->data;
10   PetscInt       refine = 0,maxnlevels = 100,*refx,*refy,*refz,n,i;
11   PetscBool      negativeMNP = PETSC_FALSE,bM = PETSC_FALSE,bN = PETSC_FALSE, bP = PETSC_FALSE;
12 
13   PetscFunctionBegin;
14   PetscValidHeaderSpecific(da,DM_CLASSID,1);
15 
16   if (dd->M < 0) {
17     dd->M       = -dd->M;
18     bM          = PETSC_TRUE;
19     negativeMNP = PETSC_TRUE;
20   }
21   if (dd->dim > 1 && dd->N < 0) {
22     dd->N       = -dd->N;
23     bN          = PETSC_TRUE;
24     negativeMNP = PETSC_TRUE;
25   }
26   if (dd->dim > 2 && dd->P < 0) {
27     dd->P       = -dd->P;
28     bP          = PETSC_TRUE;
29     negativeMNP = PETSC_TRUE;
30   }
31 
32   ierr = PetscOptionsHead("DMDA Options");CHKERRQ(ierr);
33     if (bM) {ierr = PetscOptionsInt("-da_grid_x","Number of grid points in x direction","DMDASetSizes",dd->M,&dd->M,PETSC_NULL);CHKERRQ(ierr);}
34     if (bN) {ierr = PetscOptionsInt("-da_grid_y","Number of grid points in y direction","DMDASetSizes",dd->N,&dd->N,PETSC_NULL);CHKERRQ(ierr);}
35     if (bP) {ierr = PetscOptionsInt("-da_grid_z","Number of grid points in z direction","DMDASetSizes",dd->P,&dd->P,PETSC_NULL);CHKERRQ(ierr);}
36     /* Handle DMDA parallel distibution */
37     ierr = PetscOptionsInt("-da_processors_x","Number of processors in x direction","DMDASetNumProcs",dd->m,&dd->m,PETSC_NULL);CHKERRQ(ierr);
38     if (dd->dim > 1) {ierr = PetscOptionsInt("-da_processors_y","Number of processors in y direction","DMDASetNumProcs",dd->n,&dd->n,PETSC_NULL);CHKERRQ(ierr);}
39     if (dd->dim > 2) {ierr = PetscOptionsInt("-da_processors_z","Number of processors in z direction","DMDASetNumProcs",dd->p,&dd->p,PETSC_NULL);CHKERRQ(ierr);}
40     /* Handle DMDA refinement */
41     ierr = PetscOptionsInt("-da_refine_x","Refinement ratio in x direction","DMDASetRefinementFactor",dd->refine_x,&dd->refine_x,PETSC_NULL);CHKERRQ(ierr);
42     if (dd->dim > 1) {ierr = PetscOptionsInt("-da_refine_y","Refinement ratio in y direction","DMDASetRefinementFactor",dd->refine_y,&dd->refine_y,PETSC_NULL);CHKERRQ(ierr);}
43     if (dd->dim > 2) {ierr = PetscOptionsInt("-da_refine_z","Refinement ratio in z direction","DMDASetRefinementFactor",dd->refine_z,&dd->refine_z,PETSC_NULL);CHKERRQ(ierr);}
44     dd->coarsen_x = dd->refine_x; dd->coarsen_y = dd->refine_y; dd->coarsen_z = dd->refine_z;
45 
46     /* Get refinement factors, defaults taken from the coarse DMDA */
47     ierr = PetscMalloc3(maxnlevels,PetscInt,&refx,maxnlevels,PetscInt,&refy,maxnlevels,PetscInt,&refz);CHKERRQ(ierr);
48     ierr = DMDAGetRefinementFactor(da,&refx[0],&refy[0],&refz[0]);CHKERRQ(ierr);
49     for (i=1; i<maxnlevels; i++) {
50       refx[i] = refx[0];
51       refy[i] = refy[0];
52       refz[i] = refz[0];
53     }
54     n = maxnlevels;
55     ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_x",refx,&n,PETSC_NULL);CHKERRQ(ierr);
56     if (da->levelup - da->leveldown >= 0) dd->refine_x = refx[da->levelup - da->leveldown];
57     if (da->levelup - da->leveldown >= 1) dd->coarsen_x = refx[da->levelup - da->leveldown - 1];
58     if (dd->dim > 1) {
59       n = maxnlevels;
60       ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_y",refy,&n,PETSC_NULL);CHKERRQ(ierr);
61       if (da->levelup - da->leveldown >= 0) dd->refine_y = refy[da->levelup - da->leveldown];
62       if (da->levelup - da->leveldown >= 1) dd->coarsen_y = refy[da->levelup - da->leveldown - 1];
63     }
64     if (dd->dim > 2) {
65       n = maxnlevels;
66       ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_z",refz,&n,PETSC_NULL);CHKERRQ(ierr);
67       if (da->levelup - da->leveldown >= 0) dd->refine_z = refz[da->levelup - da->leveldown];
68       if (da->levelup - da->leveldown >= 1) dd->coarsen_z = refz[da->levelup - da->leveldown - 1];
69     }
70 
71     if (negativeMNP) {ierr = PetscOptionsInt("-da_refine","Uniformly refine DA one or more times","None",refine,&refine,PETSC_NULL);CHKERRQ(ierr);}
72   ierr = PetscOptionsTail();CHKERRQ(ierr);
73 
74   while (refine--) {
75     if (da->levelup - da->leveldown >= 0) {
76       dd->refine_x = refx[da->levelup - da->leveldown];
77       dd->refine_y = refy[da->levelup - da->leveldown];
78       dd->refine_z = refz[da->levelup - da->leveldown];
79     }
80     if (dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){
81       dd->M = dd->refine_x*dd->M;
82     } else {
83       dd->M = 1 + dd->refine_x*(dd->M - 1);
84     }
85     if (dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){
86       dd->N = dd->refine_y*dd->N;
87     } else {
88       dd->N = 1 + dd->refine_y*(dd->N - 1);
89     }
90     if (dd->bz == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){
91       dd->P = dd->refine_z*dd->P;
92     } else {
93       dd->P = 1 + dd->refine_z*(dd->P - 1);
94     }
95     da->levelup++;
96   }
97   ierr = PetscFree3(refx,refy,refz);CHKERRQ(ierr);
98   PetscFunctionReturn(0);
99 }
100 
101 extern PetscErrorCode  DMCreateGlobalVector_DA(DM,Vec*);
102 extern PetscErrorCode  DMCreateLocalVector_DA(DM,Vec*);
103 extern PetscErrorCode  DMGlobalToLocalBegin_DA(DM,Vec,InsertMode,Vec);
104 extern PetscErrorCode  DMGlobalToLocalEnd_DA(DM,Vec,InsertMode,Vec);
105 extern PetscErrorCode  DMLocalToGlobalBegin_DA(DM,Vec,InsertMode,Vec);
106 extern PetscErrorCode  DMLocalToGlobalEnd_DA(DM,Vec,InsertMode,Vec);
107 extern PetscErrorCode  DMCreateInterpolation_DA(DM,DM,Mat*,Vec*);
108 extern PetscErrorCode  DMCreateColoring_DA(DM,ISColoringType,const MatType,ISColoring*);
109 extern PetscErrorCode  DMCreateMatrix_DA(DM,const MatType,Mat*);
110 extern PetscErrorCode  DMRefine_DA(DM,MPI_Comm,DM*);
111 extern PetscErrorCode  DMCoarsen_DA(DM,MPI_Comm,DM*);
112 extern PetscErrorCode  DMRefineHierarchy_DA(DM,PetscInt,DM[]);
113 extern PetscErrorCode  DMCoarsenHierarchy_DA(DM,PetscInt,DM[]);
114 extern PetscErrorCode  DMCreateInjection_DA(DM,DM,VecScatter*);
115 extern PetscErrorCode  DMCreateAggregates_DA(DM,DM,Mat*);
116 extern PetscErrorCode  DMView_DA(DM,PetscViewer);
117 extern PetscErrorCode  DMSetUp_DA(DM);
118 extern PetscErrorCode  DMDestroy_DA(DM);
119 
120 #undef __FUNCT__
121 #define __FUNCT__ "DMLoad_DA"
122 PetscErrorCode DMLoad_DA(DM da,PetscViewer viewer)
123 {
124   PetscErrorCode   ierr;
125   PetscInt         dim,m,n,p,dof,swidth;
126   DMDAStencilType  stencil;
127   DMDABoundaryType bx,by,bz;
128   PetscInt         classid = DM_FILE_CLASSID,subclassid = DMDA_FILE_CLASSID;
129   PetscBool        coors;
130   DM               dac;
131   Vec              c;
132 
133   PetscFunctionBegin;
134   ierr = PetscViewerBinaryRead(viewer,&classid,1,PETSC_INT);CHKERRQ(ierr);
135   if (classid != DM_FILE_CLASSID) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONG,"Not DM next in file");
136   ierr = PetscViewerBinaryRead(viewer,&subclassid,1,PETSC_INT);CHKERRQ(ierr);
137   if (subclassid != DMDA_FILE_CLASSID) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONG,"Not DM DA next in file");
138   ierr = PetscViewerBinaryRead(viewer,&dim,1,PETSC_INT);CHKERRQ(ierr);
139   ierr = PetscViewerBinaryRead(viewer,&m,1,PETSC_INT);CHKERRQ(ierr);
140   ierr = PetscViewerBinaryRead(viewer,&n,1,PETSC_INT);CHKERRQ(ierr);
141   ierr = PetscViewerBinaryRead(viewer,&p,1,PETSC_INT);CHKERRQ(ierr);
142   ierr = PetscViewerBinaryRead(viewer,&dof,1,PETSC_INT);CHKERRQ(ierr);
143   ierr = PetscViewerBinaryRead(viewer,&swidth,1,PETSC_INT);CHKERRQ(ierr);
144   ierr = PetscViewerBinaryRead(viewer,&bx,1,PETSC_ENUM);CHKERRQ(ierr);
145   ierr = PetscViewerBinaryRead(viewer,&by,1,PETSC_ENUM);CHKERRQ(ierr);
146   ierr = PetscViewerBinaryRead(viewer,&bz,1,PETSC_ENUM);CHKERRQ(ierr);
147   ierr = PetscViewerBinaryRead(viewer,&stencil,1,PETSC_ENUM);CHKERRQ(ierr);
148 
149   ierr = DMDASetDim(da, dim);CHKERRQ(ierr);
150   ierr = DMDASetSizes(da, m,n,p);CHKERRQ(ierr);
151   ierr = DMDASetBoundaryType(da, bx, by, bz);CHKERRQ(ierr);
152   ierr = DMDASetDof(da, dof);CHKERRQ(ierr);
153   ierr = DMDASetStencilType(da, stencil);CHKERRQ(ierr);
154   ierr = DMDASetStencilWidth(da, swidth);CHKERRQ(ierr);
155   ierr = DMSetUp(da);CHKERRQ(ierr);
156   ierr = PetscViewerBinaryRead(viewer,&coors,1,PETSC_ENUM);CHKERRQ(ierr);
157   if (coors) {
158     ierr = DMDAGetCoordinateDA(da,&dac);CHKERRQ(ierr);
159     ierr = DMCreateGlobalVector(dac,&c);CHKERRQ(ierr);
160     ierr = VecLoad(c,viewer);CHKERRQ(ierr);
161     ierr = DMDASetCoordinates(da,c);CHKERRQ(ierr);
162     ierr = VecDestroy(&c);CHKERRQ(ierr);
163   }
164   PetscFunctionReturn(0);
165 }
166 
167 #undef __FUNCT__
168 #define __FUNCT__ "DMCreateFieldDecomposition_DA"
169 PetscErrorCode DMCreateFieldDecomposition_DA(DM dm, PetscInt *len,char ***namelist, IS **islist, DM** dmlist)
170 {
171   PetscInt       i;
172   PetscErrorCode ierr;
173   DM_DA          *dd = (DM_DA*)dm->data;
174   PetscInt       dof = dd->w;
175 
176   PetscFunctionBegin;
177   if(len) *len = dof;
178   if (islist) {
179     Vec      v;
180     PetscInt rstart,n;
181 
182     ierr = DMGetGlobalVector(dm,&v);CHKERRQ(ierr);
183     ierr = VecGetOwnershipRange(v,&rstart,PETSC_NULL);CHKERRQ(ierr);
184     ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
185     ierr = DMRestoreGlobalVector(dm,&v);CHKERRQ(ierr);
186     ierr = PetscMalloc(dof*sizeof(IS),islist);CHKERRQ(ierr);
187     for (i=0; i<dof; i++) {
188       ierr = ISCreateStride(((PetscObject)dm)->comm,n/dof,rstart+i,dof,&(*islist)[i]);CHKERRQ(ierr);
189     }
190   }
191   if (namelist) {
192     ierr = PetscMalloc(dof*sizeof(const char *), namelist);CHKERRQ(ierr);
193     if (dd->fieldname) {
194       for (i=0; i<dof; i++) {
195         ierr = PetscStrallocpy(dd->fieldname[i],&(*namelist)[i]);CHKERRQ(ierr);
196       }
197     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Currently DMDA must have fieldnames");
198   }
199   if (dmlist) {
200     DM da;
201 
202     ierr = DMDACreate(((PetscObject)dm)->comm, &da);CHKERRQ(ierr);
203     ierr = DMDASetDim(da, dd->dim);CHKERRQ(ierr);
204     ierr = DMDASetSizes(da, dd->M, dd->N, dd->P);CHKERRQ(ierr);
205     ierr = DMDASetNumProcs(da, dd->m, dd->n, dd->p);CHKERRQ(ierr);
206     ierr = DMDASetBoundaryType(da, dd->bx, dd->by, dd->bz);CHKERRQ(ierr);
207     ierr = DMDASetDof(da, 1);CHKERRQ(ierr);
208     ierr = DMDASetStencilType(da, dd->stencil_type);CHKERRQ(ierr);
209     ierr = DMDASetStencilWidth(da, dd->s);CHKERRQ(ierr);
210     ierr = DMSetUp(da);CHKERRQ(ierr);
211     ierr = PetscMalloc(dof*sizeof(DM),dmlist);CHKERRQ(ierr);
212     for (i=0; i<dof-1; i++) {ierr = PetscObjectReference((PetscObject)da);CHKERRQ(ierr);}
213     for (i=0; i<dof; i++) (*dmlist)[i] = da;
214   }
215 
216   PetscFunctionReturn(0);
217 }
218 
219 
220 EXTERN_C_BEGIN
221 #undef __FUNCT__
222 #define __FUNCT__ "DMCreate_DA"
223 PetscErrorCode  DMCreate_DA(DM da)
224 {
225   PetscErrorCode ierr;
226   DM_DA          *dd;
227 
228   PetscFunctionBegin;
229   PetscValidPointer(da,1);
230   ierr = PetscNewLog(da,DM_DA,&dd);CHKERRQ(ierr);
231   da->data = dd;
232 
233   dd->dim        = -1;
234   dd->interptype = DMDA_Q1;
235   dd->refine_x   = 2;
236   dd->refine_y   = 2;
237   dd->refine_z   = 2;
238   dd->coarsen_x  = 2;
239   dd->coarsen_y  = 2;
240   dd->coarsen_z  = 2;
241   dd->fieldname  = PETSC_NULL;
242   dd->nlocal     = -1;
243   dd->Nlocal     = -1;
244   dd->M          = -1;
245   dd->N          = -1;
246   dd->P          = -1;
247   dd->m          = -1;
248   dd->n          = -1;
249   dd->p          = -1;
250   dd->w          = -1;
251   dd->s          = -1;
252   dd->xs = -1; dd->xe = -1; dd->ys = -1; dd->ye = -1; dd->zs = -1; dd->ze = -1;
253   dd->Xs = -1; dd->Xe = -1; dd->Ys = -1; dd->Ye = -1; dd->Zs = -1; dd->Ze = -1;
254 
255   dd->gtol         = PETSC_NULL;
256   dd->ltog         = PETSC_NULL;
257   dd->ltol         = PETSC_NULL;
258   dd->ao           = PETSC_NULL;
259   dd->base         = -1;
260   dd->bx         = DMDA_BOUNDARY_NONE;
261   dd->by         = DMDA_BOUNDARY_NONE;
262   dd->bz         = DMDA_BOUNDARY_NONE;
263   dd->stencil_type = DMDA_STENCIL_BOX;
264   dd->interptype   = DMDA_Q1;
265   dd->idx          = PETSC_NULL;
266   dd->Nl           = -1;
267   dd->lx           = PETSC_NULL;
268   dd->ly           = PETSC_NULL;
269   dd->lz           = PETSC_NULL;
270 
271   dd->elementtype  = DMDA_ELEMENT_Q1;
272 
273   ierr = PetscStrallocpy(VECSTANDARD,&da->vectype);CHKERRQ(ierr);
274   da->ops->globaltolocalbegin  = DMGlobalToLocalBegin_DA;
275   da->ops->globaltolocalend    = DMGlobalToLocalEnd_DA;
276   da->ops->localtoglobalbegin  = DMLocalToGlobalBegin_DA;
277   da->ops->localtoglobalend    = DMLocalToGlobalEnd_DA;
278   da->ops->createglobalvector  = DMCreateGlobalVector_DA;
279   da->ops->createlocalvector   = DMCreateLocalVector_DA;
280   da->ops->createinterpolation = DMCreateInterpolation_DA;
281   da->ops->getcoloring         = DMCreateColoring_DA;
282   da->ops->creatematrix        = DMCreateMatrix_DA;
283   da->ops->refine              = DMRefine_DA;
284   da->ops->coarsen             = DMCoarsen_DA;
285   da->ops->refinehierarchy     = DMRefineHierarchy_DA;
286   da->ops->coarsenhierarchy    = DMCoarsenHierarchy_DA;
287   da->ops->getinjection        = DMCreateInjection_DA;
288   da->ops->getaggregates       = DMCreateAggregates_DA;
289   da->ops->destroy             = DMDestroy_DA;
290   da->ops->view                = 0;
291   da->ops->setfromoptions      = DMSetFromOptions_DA;
292   da->ops->setup               = DMSetUp_DA;
293   da->ops->load                = DMLoad_DA;
294   da->ops->createfielddecomposition = DMCreateFieldDecomposition_DA;
295   PetscFunctionReturn(0);
296 }
297 EXTERN_C_END
298 
299 #undef __FUNCT__
300 #define __FUNCT__ "DMDACreate"
301 /*@
302   DMDACreate - Creates a DMDA object.
303 
304   Collective on MPI_Comm
305 
306   Input Parameter:
307 . comm - The communicator for the DMDA object
308 
309   Output Parameter:
310 . da  - The DMDA object
311 
312   Level: advanced
313 
314   Developers Note: Since there exists DMDACreate1/2/3d() should this routine even exist?
315 
316 .keywords: DMDA, create
317 .seealso:  DMDASetSizes(), DMDADuplicate(),  DMDACreate1d(), DMDACreate2d(), DMDACreate3d()
318 @*/
319 PetscErrorCode  DMDACreate(MPI_Comm comm, DM *da)
320 {
321   PetscErrorCode ierr;
322 
323   PetscFunctionBegin;
324   PetscValidPointer(da,2);
325   ierr = DMCreate(comm,da);CHKERRQ(ierr);
326   ierr = DMSetType(*da,DMDA);CHKERRQ(ierr);
327   PetscFunctionReturn(0);
328 }
329 
330 
331