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