xref: /petsc/src/dm/impls/da/dacreate.c (revision 089b283744364aef00a310a92368c00bc3aa30b8)
1 
2 #include <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   PetscBool      flg;
10   char           typeName[256];
11   DM_DA          *dd = (DM_DA*)da->data;
12   PetscInt       refine = 0;
13   PetscBool      negativeMNP = PETSC_FALSE;
14 
15   PetscFunctionBegin;
16   PetscValidHeaderSpecific(da,DM_CLASSID,1);
17 
18   ierr = PetscOptionsBegin(((PetscObject)da)->comm,((PetscObject)da)->prefix,"DMDA Options","DMDA");CHKERRQ(ierr);
19     /* Handle DMDA grid sizes */
20     if (dd->M < 0) {
21       PetscInt newM = -dd->M;
22       ierr        = PetscOptionsInt("-da_grid_x","Number of grid points in x direction","DMDASetSizes",newM,&newM,PETSC_NULL);CHKERRQ(ierr);
23       dd->M       = newM;
24       negativeMNP = PETSC_TRUE;
25     }
26     if (dd->dim > 1 && dd->N < 0) {
27       PetscInt newN = -dd->N;
28       ierr        = PetscOptionsInt("-da_grid_y","Number of grid points in y direction","DMDASetSizes",newN,&newN,PETSC_NULL);CHKERRQ(ierr);
29       dd->N       = newN;
30       negativeMNP = PETSC_TRUE;
31     }
32     if (dd->dim > 2 && dd->P < 0) {
33       PetscInt newP = -dd->P;
34       ierr        = PetscOptionsInt("-da_grid_z","Number of grid points in z direction","DMDASetSizes",newP,&newP,PETSC_NULL);CHKERRQ(ierr);
35       dd->P       = newP;
36       negativeMNP = PETSC_TRUE;
37     }
38     /* Handle DMDA parallel distibution */
39     ierr = PetscOptionsInt("-da_processors_x","Number of processors in x direction","DMDASetNumProcs",dd->m,&dd->m,PETSC_NULL);CHKERRQ(ierr);
40     if (dd->dim > 1) {ierr = PetscOptionsInt("-da_processors_y","Number of processors in y direction","DMDASetNumProcs",dd->n,&dd->n,PETSC_NULL);CHKERRQ(ierr);}
41     if (dd->dim > 2) {ierr = PetscOptionsInt("-da_processors_z","Number of processors in z direction","DMDASetNumProcs",dd->p,&dd->p,PETSC_NULL);CHKERRQ(ierr);}
42     /* Handle DMDA refinement */
43     ierr = PetscOptionsInt("-da_refine_x","Refinement ratio in x direction","DMDASetRefinementFactor",dd->refine_x,&dd->refine_x,PETSC_NULL);CHKERRQ(ierr);
44     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);}
45     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);}
46 
47     if (!VecRegisterAllCalled) {ierr = VecRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
48     ierr = PetscOptionsList("-da_vec_type","Vector type used for created vectors","DMSetVecType",VecList,da->vectype,typeName,256,&flg);CHKERRQ(ierr);
49     if (flg) {
50       ierr = DMSetVecType(da,typeName);CHKERRQ(ierr);
51     }
52 
53     if (negativeMNP) {
54       ierr = PetscOptionsInt("-da_refine","Uniformly refine DA one or more times","None",refine,&refine,PETSC_NULL);CHKERRQ(ierr);
55       while (refine--) {
56         if (dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){
57           dd->M = dd->refine_x*dd->M;
58         } else {
59           dd->M = 1 + dd->refine_x*(dd->M - 1);
60         }
61         if (dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){
62           dd->N = dd->refine_y*dd->N;
63         } else {
64           dd->N = 1 + dd->refine_y*(dd->N - 1);
65         }
66         if (dd->bz == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){
67           dd->P = dd->refine_z*dd->P;
68         } else {
69           dd->P = 1 + dd->refine_z*(dd->P - 1);
70         }
71         da->levelup++;
72       }
73     }
74 
75 
76     /* process any options handlers added with PetscObjectAddOptionsHandler() */
77     ierr = PetscObjectProcessOptionsHandlers((PetscObject)da);CHKERRQ(ierr);
78   ierr = PetscOptionsEnd();CHKERRQ(ierr);
79 
80   PetscFunctionReturn(0);
81 }
82 
83 extern PetscErrorCode  DMCreateGlobalVector_DA(DM,Vec*);
84 extern PetscErrorCode  DMCreateLocalVector_DA(DM,Vec*);
85 extern PetscErrorCode  DMGlobalToLocalBegin_DA(DM,Vec,InsertMode,Vec);
86 extern PetscErrorCode  DMGlobalToLocalEnd_DA(DM,Vec,InsertMode,Vec);
87 extern PetscErrorCode  DMLocalToGlobalBegin_DA(DM,Vec,InsertMode,Vec);
88 extern PetscErrorCode  DMLocalToGlobalEnd_DA(DM,Vec,InsertMode,Vec);
89 extern PetscErrorCode  DMGetInterpolation_DA(DM,DM,Mat*,Vec*);
90 extern PetscErrorCode  DMGetColoring_DA(DM,ISColoringType,const MatType,ISColoring*);
91 extern PetscErrorCode  DMGetMatrix_DA(DM,const MatType,Mat*);
92 extern PetscErrorCode  DMRefine_DA(DM,MPI_Comm,DM*);
93 extern PetscErrorCode  DMCoarsen_DA(DM,MPI_Comm,DM*);
94 extern PetscErrorCode  DMRefineHierarchy_DA(DM,PetscInt,DM[]);
95 extern PetscErrorCode  DMCoarsenHierarchy_DA(DM,PetscInt,DM[]);
96 extern PetscErrorCode  DMGetInjection_DA(DM,DM,VecScatter*);
97 extern PetscErrorCode  DMGetAggregates_DA(DM,DM,Mat*);
98 extern PetscErrorCode  DMView_DA(DM,PetscViewer);
99 extern PetscErrorCode  DMSetUp_DA(DM);
100 extern PetscErrorCode  DMDestroy_DA(DM);
101 
102 EXTERN_C_BEGIN
103 #undef __FUNCT__
104 #define __FUNCT__ "DMCreate_DA"
105 PetscErrorCode  DMCreate_DA(DM da)
106 {
107   PetscErrorCode ierr;
108   DM_DA          *dd;
109 
110   PetscFunctionBegin;
111   PetscValidPointer(da,1);
112   ierr = PetscNewLog(da,DM_DA,&dd);CHKERRQ(ierr);
113   da->data = dd;
114 
115   dd->dim        = -1;
116   dd->interptype = DMDA_Q1;
117   dd->refine_x   = 2;
118   dd->refine_y   = 2;
119   dd->refine_z   = 2;
120   dd->fieldname  = PETSC_NULL;
121   dd->nlocal     = -1;
122   dd->Nlocal     = -1;
123   dd->M          = -1;
124   dd->N          = -1;
125   dd->P          = -1;
126   dd->m          = -1;
127   dd->n          = -1;
128   dd->p          = -1;
129   dd->w          = -1;
130   dd->s          = -1;
131   dd->xs = -1; dd->xe = -1; dd->ys = -1; dd->ye = -1; dd->zs = -1; dd->ze = -1;
132   dd->Xs = -1; dd->Xe = -1; dd->Ys = -1; dd->Ye = -1; dd->Zs = -1; dd->Ze = -1;
133 
134   dd->gtol         = PETSC_NULL;
135   dd->ltog         = PETSC_NULL;
136   dd->ltol         = PETSC_NULL;
137   dd->ao           = PETSC_NULL;
138   dd->base         = -1;
139   dd->bx         = DMDA_BOUNDARY_NONE;
140   dd->by         = DMDA_BOUNDARY_NONE;
141   dd->bz         = DMDA_BOUNDARY_NONE;
142   dd->stencil_type = DMDA_STENCIL_BOX;
143   dd->interptype   = DMDA_Q1;
144   dd->idx          = PETSC_NULL;
145   dd->Nl           = -1;
146   dd->lx           = PETSC_NULL;
147   dd->ly           = PETSC_NULL;
148   dd->lz           = PETSC_NULL;
149 
150   dd->elementtype  = DMDA_ELEMENT_Q1;
151 
152   ierr = PetscStrallocpy(VECSTANDARD,&da->vectype);CHKERRQ(ierr);
153   da->ops->globaltolocalbegin = DMGlobalToLocalBegin_DA;
154   da->ops->globaltolocalend   = DMGlobalToLocalEnd_DA;
155   da->ops->localtoglobalbegin = DMLocalToGlobalBegin_DA;
156   da->ops->localtoglobalend   = DMLocalToGlobalEnd_DA;
157   da->ops->createglobalvector = DMCreateGlobalVector_DA;
158   da->ops->createlocalvector  = DMCreateLocalVector_DA;
159   da->ops->getinterpolation   = DMGetInterpolation_DA;
160   da->ops->getcoloring        = DMGetColoring_DA;
161   da->ops->getmatrix          = DMGetMatrix_DA;
162   da->ops->refine             = DMRefine_DA;
163   da->ops->coarsen            = DMCoarsen_DA;
164   da->ops->refinehierarchy    = DMRefineHierarchy_DA;
165   da->ops->coarsenhierarchy   = DMCoarsenHierarchy_DA;
166   da->ops->getinjection       = DMGetInjection_DA;
167   da->ops->getaggregates      = DMGetAggregates_DA;
168   da->ops->destroy            = DMDestroy_DA;
169   da->ops->view               = 0;
170   da->ops->setfromoptions     = DMSetFromOptions_DA;
171   da->ops->setup              = DMSetUp_DA;
172   PetscFunctionReturn(0);
173 }
174 EXTERN_C_END
175 
176 #undef __FUNCT__
177 #define __FUNCT__ "DMDACreate"
178 /*@
179   DMDACreate - Creates a DMDA object.
180 
181   Collective on MPI_Comm
182 
183   Input Parameter:
184 . comm - The communicator for the DMDA object
185 
186   Output Parameter:
187 . da  - The DMDA object
188 
189   Level: advanced
190 
191   Developers Note: Since there exists DMDACreate1/2/3d() should this routine even exist?
192 
193 .keywords: DMDA, create
194 .seealso:  DMDASetSizes(), DMDADuplicate(),  DMDACreate1d(), DMDACreate2d(), DMDACreate3d()
195 @*/
196 PetscErrorCode  DMDACreate(MPI_Comm comm, DM *da)
197 {
198   PetscErrorCode ierr;
199 
200   PetscFunctionBegin;
201   PetscValidPointer(da,2);
202   ierr = DMCreate(comm,da);CHKERRQ(ierr);
203   ierr = DMSetType(*da,DMDA);CHKERRQ(ierr);
204   PetscFunctionReturn(0);
205 }
206