xref: /petsc/src/dm/impls/da/dacreate.c (revision 8865f1ea4ea560cd84ab8db62e98b7095cdff96f)
17d0a6c19SBarry Smith 
2b45d2f2cSJed Brown #include <petsc-private/daimpl.h>    /*I   "petscdmda.h"   I*/
347c6ae99SBarry Smith 
447c6ae99SBarry Smith #undef __FUNCT__
59a42bb27SBarry Smith #define __FUNCT__ "DMSetFromOptions_DA"
67087cfbeSBarry Smith PetscErrorCode  DMSetFromOptions_DA(DM da)
747c6ae99SBarry Smith {
847c6ae99SBarry Smith   PetscErrorCode ierr;
947c6ae99SBarry Smith   DM_DA          *dd         = (DM_DA*)da->data;
10397b6216SJed Brown   PetscInt       refine      = 0,maxnlevels = 100,*refx,*refy,*refz,n,i;
11235683edSBarry Smith   PetscBool      negativeMNP = PETSC_FALSE,bM = PETSC_FALSE,bN = PETSC_FALSE, bP = PETSC_FALSE;
1247c6ae99SBarry Smith 
1347c6ae99SBarry Smith   PetscFunctionBegin;
1447c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
1547c6ae99SBarry Smith 
1647c6ae99SBarry Smith   if (dd->M < 0) {
17235683edSBarry Smith     dd->M       = -dd->M;
18235683edSBarry Smith     bM          = PETSC_TRUE;
19e0f5d30fSBarry Smith     negativeMNP = PETSC_TRUE;
2047c6ae99SBarry Smith   }
2147c6ae99SBarry Smith   if (dd->dim > 1 && dd->N < 0) {
22235683edSBarry Smith     dd->N       = -dd->N;
23235683edSBarry Smith     bN          = PETSC_TRUE;
24e0f5d30fSBarry Smith     negativeMNP = PETSC_TRUE;
2547c6ae99SBarry Smith   }
2647c6ae99SBarry Smith   if (dd->dim > 2 && dd->P < 0) {
27235683edSBarry Smith     dd->P       = -dd->P;
28235683edSBarry Smith     bP          = PETSC_TRUE;
29e0f5d30fSBarry Smith     negativeMNP = PETSC_TRUE;
3047c6ae99SBarry Smith   }
31235683edSBarry Smith 
32f9ba7244SBarry Smith   ierr = PetscOptionsHead("DMDA Options");CHKERRQ(ierr);
33235683edSBarry Smith   if (bM) {ierr = PetscOptionsInt("-da_grid_x","Number of grid points in x direction","DMDASetSizes",dd->M,&dd->M,PETSC_NULL);CHKERRQ(ierr);}
34235683edSBarry Smith   if (bN) {ierr = PetscOptionsInt("-da_grid_y","Number of grid points in y direction","DMDASetSizes",dd->N,&dd->N,PETSC_NULL);CHKERRQ(ierr);}
35235683edSBarry Smith   if (bP) {ierr = PetscOptionsInt("-da_grid_z","Number of grid points in z direction","DMDASetSizes",dd->P,&dd->P,PETSC_NULL);CHKERRQ(ierr);}
3688661749SPeter Brune   ierr = PetscOptionsInt("-da_overlap","Overlap between local grids","DMDASetOverlap",dd->overlap,&dd->overlap,PETSC_NULL);CHKERRQ(ierr);
37aa219208SBarry Smith   /* Handle DMDA parallel distibution */
38aa219208SBarry Smith   ierr = PetscOptionsInt("-da_processors_x","Number of processors in x direction","DMDASetNumProcs",dd->m,&dd->m,PETSC_NULL);CHKERRQ(ierr);
39aa219208SBarry Smith   if (dd->dim > 1) {ierr = PetscOptionsInt("-da_processors_y","Number of processors in y direction","DMDASetNumProcs",dd->n,&dd->n,PETSC_NULL);CHKERRQ(ierr);}
40aa219208SBarry Smith   if (dd->dim > 2) {ierr = PetscOptionsInt("-da_processors_z","Number of processors in z direction","DMDASetNumProcs",dd->p,&dd->p,PETSC_NULL);CHKERRQ(ierr);}
41aa219208SBarry Smith   /* Handle DMDA refinement */
42aa219208SBarry Smith   ierr = PetscOptionsInt("-da_refine_x","Refinement ratio in x direction","DMDASetRefinementFactor",dd->refine_x,&dd->refine_x,PETSC_NULL);CHKERRQ(ierr);
43aa219208SBarry Smith   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);}
44aa219208SBarry Smith   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);}
4538fb6c91SJed Brown   dd->coarsen_x = dd->refine_x; dd->coarsen_y = dd->refine_y; dd->coarsen_z = dd->refine_z;
4647c6ae99SBarry Smith 
47397b6216SJed Brown   /* Get refinement factors, defaults taken from the coarse DMDA */
48397b6216SJed Brown   ierr = PetscMalloc3(maxnlevels,PetscInt,&refx,maxnlevels,PetscInt,&refy,maxnlevels,PetscInt,&refz);CHKERRQ(ierr);
49397b6216SJed Brown   ierr = DMDAGetRefinementFactor(da,&refx[0],&refy[0],&refz[0]);CHKERRQ(ierr);
50397b6216SJed Brown   for (i=1; i<maxnlevels; i++) {
51397b6216SJed Brown     refx[i] = refx[0];
52397b6216SJed Brown     refy[i] = refy[0];
53397b6216SJed Brown     refz[i] = refz[0];
54397b6216SJed Brown   }
55397b6216SJed Brown   n    = maxnlevels;
56397b6216SJed Brown   ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_x",refx,&n,PETSC_NULL);CHKERRQ(ierr);
57397b6216SJed Brown   if (da->levelup - da->leveldown >= 0) dd->refine_x = refx[da->levelup - da->leveldown];
58397b6216SJed Brown   if (da->levelup - da->leveldown >= 1) dd->coarsen_x = refx[da->levelup - da->leveldown - 1];
59397b6216SJed Brown   if (dd->dim > 1) {
60397b6216SJed Brown     n    = maxnlevels;
61397b6216SJed Brown     ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_y",refy,&n,PETSC_NULL);CHKERRQ(ierr);
62397b6216SJed Brown     if (da->levelup - da->leveldown >= 0) dd->refine_y = refy[da->levelup - da->leveldown];
63397b6216SJed Brown     if (da->levelup - da->leveldown >= 1) dd->coarsen_y = refy[da->levelup - da->leveldown - 1];
64397b6216SJed Brown   }
65397b6216SJed Brown   if (dd->dim > 2) {
66397b6216SJed Brown     n    = maxnlevels;
67397b6216SJed Brown     ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_z",refz,&n,PETSC_NULL);CHKERRQ(ierr);
68397b6216SJed Brown     if (da->levelup - da->leveldown >= 0) dd->refine_z = refz[da->levelup - da->leveldown];
69397b6216SJed Brown     if (da->levelup - da->leveldown >= 1) dd->coarsen_z = refz[da->levelup - da->leveldown - 1];
70397b6216SJed Brown   }
71397b6216SJed Brown 
72235683edSBarry Smith   if (negativeMNP) {ierr = PetscOptionsInt("-da_refine","Uniformly refine DA one or more times","None",refine,&refine,PETSC_NULL);CHKERRQ(ierr);}
73f9ba7244SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
74235683edSBarry Smith 
75e0f5d30fSBarry Smith   while (refine--) {
76e0f5d30fSBarry Smith     if (dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
77e0f5d30fSBarry Smith       dd->M = dd->refine_x*dd->M;
78e0f5d30fSBarry Smith     } else {
79e0f5d30fSBarry Smith       dd->M = 1 + dd->refine_x*(dd->M - 1);
80e0f5d30fSBarry Smith     }
81e0f5d30fSBarry Smith     if (dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
82e0f5d30fSBarry Smith       dd->N = dd->refine_y*dd->N;
83e0f5d30fSBarry Smith     } else {
84e0f5d30fSBarry Smith       dd->N = 1 + dd->refine_y*(dd->N - 1);
85e0f5d30fSBarry Smith     }
86e0f5d30fSBarry Smith     if (dd->bz == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) {
87e0f5d30fSBarry Smith       dd->P = dd->refine_z*dd->P;
88e0f5d30fSBarry Smith     } else {
89e0f5d30fSBarry Smith       dd->P = 1 + dd->refine_z*(dd->P - 1);
90e0f5d30fSBarry Smith     }
91e0f5d30fSBarry Smith     da->levelup++;
9285fc4b34SJed Brown     if (da->levelup - da->leveldown >= 0) {
9385fc4b34SJed Brown       dd->refine_x = refx[da->levelup - da->leveldown];
9485fc4b34SJed Brown       dd->refine_y = refy[da->levelup - da->leveldown];
9585fc4b34SJed Brown       dd->refine_z = refz[da->levelup - da->leveldown];
9685fc4b34SJed Brown     }
9785fc4b34SJed Brown     if (da->levelup - da->leveldown >= 1) {
9885fc4b34SJed Brown       dd->coarsen_x = refx[da->levelup - da->leveldown - 1];
9985fc4b34SJed Brown       dd->coarsen_y = refy[da->levelup - da->leveldown - 1];
10085fc4b34SJed Brown       dd->coarsen_z = refz[da->levelup - da->leveldown - 1];
10185fc4b34SJed Brown     }
102e0f5d30fSBarry Smith   }
103397b6216SJed Brown   ierr = PetscFree3(refx,refy,refz);CHKERRQ(ierr);
10447c6ae99SBarry Smith   PetscFunctionReturn(0);
10547c6ae99SBarry Smith }
10647c6ae99SBarry Smith 
1077087cfbeSBarry Smith extern PetscErrorCode  DMCreateGlobalVector_DA(DM,Vec*);
1087087cfbeSBarry Smith extern PetscErrorCode  DMCreateLocalVector_DA(DM,Vec*);
1097087cfbeSBarry Smith extern PetscErrorCode  DMGlobalToLocalBegin_DA(DM,Vec,InsertMode,Vec);
1107087cfbeSBarry Smith extern PetscErrorCode  DMGlobalToLocalEnd_DA(DM,Vec,InsertMode,Vec);
1117087cfbeSBarry Smith extern PetscErrorCode  DMLocalToGlobalBegin_DA(DM,Vec,InsertMode,Vec);
1127087cfbeSBarry Smith extern PetscErrorCode  DMLocalToGlobalEnd_DA(DM,Vec,InsertMode,Vec);
113e727c939SJed Brown extern PetscErrorCode  DMCreateInterpolation_DA(DM,DM,Mat*,Vec*);
11419fd82e9SBarry Smith extern PetscErrorCode  DMCreateColoring_DA(DM,ISColoringType,MatType,ISColoring*);
11519fd82e9SBarry Smith extern PetscErrorCode  DMCreateMatrix_DA(DM,MatType,Mat*);
1166636e97aSMatthew G Knepley extern PetscErrorCode  DMCreateCoordinateDM_DA(DM,DM*);
1177087cfbeSBarry Smith extern PetscErrorCode  DMRefine_DA(DM,MPI_Comm,DM*);
1187087cfbeSBarry Smith extern PetscErrorCode  DMCoarsen_DA(DM,MPI_Comm,DM*);
1197087cfbeSBarry Smith extern PetscErrorCode  DMRefineHierarchy_DA(DM,PetscInt,DM[]);
1207087cfbeSBarry Smith extern PetscErrorCode  DMCoarsenHierarchy_DA(DM,PetscInt,DM[]);
121e727c939SJed Brown extern PetscErrorCode  DMCreateInjection_DA(DM,DM,VecScatter*);
122e727c939SJed Brown extern PetscErrorCode  DMCreateAggregates_DA(DM,DM,Mat*);
1237087cfbeSBarry Smith extern PetscErrorCode  DMView_DA(DM,PetscViewer);
1247087cfbeSBarry Smith extern PetscErrorCode  DMSetUp_DA(DM);
1257087cfbeSBarry Smith extern PetscErrorCode  DMDestroy_DA(DM);
126e30e807fSPeter Brune extern PetscErrorCode  DMCreateDomainDecomposition_DA(DM,PetscInt*,char***,IS**,IS**,DM**);
1270a696f66SPeter Brune extern PetscErrorCode  DMCreateDomainDecompositionDM_DA(DM,const char*,DM*);
128e30e807fSPeter Brune extern PetscErrorCode  DMCreateDomainDecompositionScatters_DA(DM,PetscInt,DM*,VecScatter**,VecScatter**,VecScatter**);
1299a42bb27SBarry Smith 
130b859378eSBarry Smith #undef __FUNCT__
131b859378eSBarry Smith #define __FUNCT__ "DMLoad_DA"
132b859378eSBarry Smith PetscErrorCode DMLoad_DA(DM da,PetscViewer viewer)
133b859378eSBarry Smith {
134b859378eSBarry Smith   PetscErrorCode   ierr;
135bc2bf880SBarry Smith   PetscInt         dim,m,n,p,dof,swidth;
136b859378eSBarry Smith   DMDAStencilType  stencil;
137b859378eSBarry Smith   DMDABoundaryType bx,by,bz;
138bc2bf880SBarry Smith   PetscBool        coors;
139bc2bf880SBarry Smith   DM               dac;
140bc2bf880SBarry Smith   Vec              c;
141b859378eSBarry Smith 
142b859378eSBarry Smith   PetscFunctionBegin;
143b859378eSBarry Smith   ierr = PetscViewerBinaryRead(viewer,&dim,1,PETSC_INT);CHKERRQ(ierr);
144b859378eSBarry Smith   ierr = PetscViewerBinaryRead(viewer,&m,1,PETSC_INT);CHKERRQ(ierr);
145b859378eSBarry Smith   ierr = PetscViewerBinaryRead(viewer,&n,1,PETSC_INT);CHKERRQ(ierr);
146b859378eSBarry Smith   ierr = PetscViewerBinaryRead(viewer,&p,1,PETSC_INT);CHKERRQ(ierr);
147b859378eSBarry Smith   ierr = PetscViewerBinaryRead(viewer,&dof,1,PETSC_INT);CHKERRQ(ierr);
148b859378eSBarry Smith   ierr = PetscViewerBinaryRead(viewer,&swidth,1,PETSC_INT);CHKERRQ(ierr);
149b859378eSBarry Smith   ierr = PetscViewerBinaryRead(viewer,&bx,1,PETSC_ENUM);CHKERRQ(ierr);
150b859378eSBarry Smith   ierr = PetscViewerBinaryRead(viewer,&by,1,PETSC_ENUM);CHKERRQ(ierr);
151b859378eSBarry Smith   ierr = PetscViewerBinaryRead(viewer,&bz,1,PETSC_ENUM);CHKERRQ(ierr);
152b859378eSBarry Smith   ierr = PetscViewerBinaryRead(viewer,&stencil,1,PETSC_ENUM);CHKERRQ(ierr);
153b859378eSBarry Smith 
154bc2bf880SBarry Smith   ierr = DMDASetDim(da, dim);CHKERRQ(ierr);
155bc2bf880SBarry Smith   ierr = DMDASetSizes(da, m,n,p);CHKERRQ(ierr);
156b859378eSBarry Smith   ierr = DMDASetBoundaryType(da, bx, by, bz);CHKERRQ(ierr);
157b859378eSBarry Smith   ierr = DMDASetDof(da, dof);CHKERRQ(ierr);
158b859378eSBarry Smith   ierr = DMDASetStencilType(da, stencil);CHKERRQ(ierr);
159b859378eSBarry Smith   ierr = DMDASetStencilWidth(da, swidth);CHKERRQ(ierr);
160b859378eSBarry Smith   ierr = DMSetUp(da);CHKERRQ(ierr);
161bc2bf880SBarry Smith   ierr = PetscViewerBinaryRead(viewer,&coors,1,PETSC_ENUM);CHKERRQ(ierr);
162bc2bf880SBarry Smith   if (coors) {
1636636e97aSMatthew G Knepley     ierr = DMGetCoordinateDM(da,&dac);CHKERRQ(ierr);
164bc2bf880SBarry Smith     ierr = DMCreateGlobalVector(dac,&c);CHKERRQ(ierr);
165bc2bf880SBarry Smith     ierr = VecLoad(c,viewer);CHKERRQ(ierr);
1666636e97aSMatthew G Knepley     ierr = DMSetCoordinates(da,c);CHKERRQ(ierr);
167bc2bf880SBarry Smith     ierr = VecDestroy(&c);CHKERRQ(ierr);
168bc2bf880SBarry Smith   }
169b859378eSBarry Smith   PetscFunctionReturn(0);
170b859378eSBarry Smith }
171b859378eSBarry Smith 
1721c3fb106SBarry Smith #undef __FUNCT__
173d724dfffSBarry Smith #define __FUNCT__ "DMCreateSubDM_DA"
174d724dfffSBarry Smith PetscErrorCode DMCreateSubDM_DA(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm)
175d724dfffSBarry Smith {
176d724dfffSBarry Smith   PetscErrorCode ierr;
177d724dfffSBarry Smith   DM_DA          *da = (DM_DA*)dm->data;
178d724dfffSBarry Smith 
179d724dfffSBarry Smith   PetscFunctionBegin;
180d724dfffSBarry Smith   if (da->dim != 2) SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Support only implemented for 2d");
181d724dfffSBarry Smith   if (subdm) {
182d724dfffSBarry Smith     ierr = DMDACreate2d(((PetscObject)dm)->comm,da->bx,da->by,da->stencil_type,da->M,da->N,da->m,da->n,numFields,da->s,da->lx,da->ly,subdm);CHKERRQ(ierr);
183d724dfffSBarry Smith   }
184d724dfffSBarry Smith   if (is) {
185d724dfffSBarry Smith     PetscInt *indices,cnt = 0, dof = da->w,i,j;
186d724dfffSBarry Smith     ierr = PetscMalloc(sizeof(PetscInt)*da->Nlocal*numFields/dof,&indices);CHKERRQ(ierr);
187d724dfffSBarry Smith     for (i=da->base/dof; i<(da->base+da->Nlocal)/dof; i++) {
188d724dfffSBarry Smith       for (j=0; j<numFields; j++) {
189d724dfffSBarry Smith         indices[cnt++] = dof*i + fields[j];
190d724dfffSBarry Smith       }
191d724dfffSBarry Smith     }
192d724dfffSBarry Smith     if (cnt != da->Nlocal*numFields/dof) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count does not equal expected value");
193d724dfffSBarry Smith     ierr = ISCreateGeneral(((PetscObject)dm)->comm,da->Nlocal*numFields/dof,indices,PETSC_OWN_POINTER,is);CHKERRQ(ierr);
194d724dfffSBarry Smith   }
195d724dfffSBarry Smith   PetscFunctionReturn(0);
196d724dfffSBarry Smith }
197d724dfffSBarry Smith 
198d724dfffSBarry Smith #undef __FUNCT__
19916621825SDmitry Karpeev #define __FUNCT__ "DMCreateFieldDecomposition_DA"
20016621825SDmitry Karpeev PetscErrorCode DMCreateFieldDecomposition_DA(DM dm, PetscInt *len,char ***namelist, IS **islist, DM **dmlist)
2011c3fb106SBarry Smith {
2021c3fb106SBarry Smith   PetscInt       i;
2031c3fb106SBarry Smith   PetscErrorCode ierr;
2041c3fb106SBarry Smith   DM_DA          *dd = (DM_DA*)dm->data;
2051c3fb106SBarry Smith   PetscInt       dof = dd->w;
2061c3fb106SBarry Smith 
2071c3fb106SBarry Smith   PetscFunctionBegin;
208731c8d9eSDmitry Karpeev   if (len) *len = dof;
2091c3fb106SBarry Smith   if (islist) {
2101c3fb106SBarry Smith     Vec      v;
2111c3fb106SBarry Smith     PetscInt rstart,n;
2121c3fb106SBarry Smith 
2131c3fb106SBarry Smith     ierr = DMGetGlobalVector(dm,&v);CHKERRQ(ierr);
2141c3fb106SBarry Smith     ierr = VecGetOwnershipRange(v,&rstart,PETSC_NULL);CHKERRQ(ierr);
2151c3fb106SBarry Smith     ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
2161c3fb106SBarry Smith     ierr = DMRestoreGlobalVector(dm,&v);CHKERRQ(ierr);
2171c3fb106SBarry Smith     ierr = PetscMalloc(dof*sizeof(IS),islist);CHKERRQ(ierr);
2181c3fb106SBarry Smith     for (i=0; i<dof; i++) {
2191c3fb106SBarry Smith       ierr = ISCreateStride(((PetscObject)dm)->comm,n/dof,rstart+i,dof,&(*islist)[i]);CHKERRQ(ierr);
2201c3fb106SBarry Smith     }
2211c3fb106SBarry Smith   }
2221c3fb106SBarry Smith   if (namelist) {
2231c3fb106SBarry Smith     ierr = PetscMalloc(dof*sizeof(const char*), namelist);CHKERRQ(ierr);
2241c3fb106SBarry Smith     if (dd->fieldname) {
2251c3fb106SBarry Smith       for (i=0; i<dof; i++) {
2261c3fb106SBarry Smith         ierr = PetscStrallocpy(dd->fieldname[i],&(*namelist)[i]);CHKERRQ(ierr);
2271c3fb106SBarry Smith       }
2281c3fb106SBarry Smith     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Currently DMDA must have fieldnames");
2291c3fb106SBarry Smith   }
2301c3fb106SBarry Smith   if (dmlist) {
2311c3fb106SBarry Smith     DM da;
2321c3fb106SBarry Smith 
2331c3fb106SBarry Smith     ierr = DMDACreate(((PetscObject)dm)->comm, &da);CHKERRQ(ierr);
2341c3fb106SBarry Smith     ierr = DMDASetDim(da, dd->dim);CHKERRQ(ierr);
2351c3fb106SBarry Smith     ierr = DMDASetSizes(da, dd->M, dd->N, dd->P);CHKERRQ(ierr);
2361c3fb106SBarry Smith     ierr = DMDASetNumProcs(da, dd->m, dd->n, dd->p);CHKERRQ(ierr);
2371c3fb106SBarry Smith     ierr = DMDASetBoundaryType(da, dd->bx, dd->by, dd->bz);CHKERRQ(ierr);
2381c3fb106SBarry Smith     ierr = DMDASetDof(da, 1);CHKERRQ(ierr);
2391c3fb106SBarry Smith     ierr = DMDASetStencilType(da, dd->stencil_type);CHKERRQ(ierr);
2401c3fb106SBarry Smith     ierr = DMDASetStencilWidth(da, dd->s);CHKERRQ(ierr);
2411c3fb106SBarry Smith     ierr = DMSetUp(da);CHKERRQ(ierr);
2421c3fb106SBarry Smith     ierr = PetscMalloc(dof*sizeof(DM),dmlist);CHKERRQ(ierr);
243c893992eSBarry Smith     for (i=0; i<dof-1; i++) {ierr = PetscObjectReference((PetscObject)da);CHKERRQ(ierr);}
2441c3fb106SBarry Smith     for (i=0; i<dof; i++) (*dmlist)[i] = da;
2451c3fb106SBarry Smith   }
2461c3fb106SBarry Smith   PetscFunctionReturn(0);
2471c3fb106SBarry Smith }
2481c3fb106SBarry Smith 
2493efe6655SBarry Smith /*MC
2503efe6655SBarry Smith    DMDA = "da" - A DM object that is used to manage data for a structured grid in 1, 2, or 3 dimensions.
2513efe6655SBarry Smith          In the global representation of the vector each process stores a non-overlapping rectangular (or slab in 3d) portion of the grid points.
2523efe6655SBarry Smith          In the local representation these rectangular regions (slabs) are extended in all directions by a stencil width.
2533efe6655SBarry Smith 
2543efe6655SBarry Smith          The vectors can be thought of as either cell centered or vertex centered on the mesh. But some variables cannot be cell centered and others
2553efe6655SBarry Smith          vertex centered.
2563efe6655SBarry Smith 
2573efe6655SBarry Smith 
2583efe6655SBarry Smith   Level: intermediate
2593efe6655SBarry Smith 
2603efe6655SBarry Smith .seealso: DMType, DMCOMPOSITE, DMDACreate(), DMCreate(), DMSetType()
2613efe6655SBarry Smith M*/
2623efe6655SBarry Smith 
2631c3fb106SBarry Smith 
264a4121054SBarry Smith EXTERN_C_BEGIN
26547c6ae99SBarry Smith #undef __FUNCT__
266a4121054SBarry Smith #define __FUNCT__ "DMCreate_DA"
2677087cfbeSBarry Smith PetscErrorCode  DMCreate_DA(DM da)
26847c6ae99SBarry Smith {
26947c6ae99SBarry Smith   PetscErrorCode ierr;
27047c6ae99SBarry Smith   DM_DA          *dd;
27147c6ae99SBarry Smith 
27247c6ae99SBarry Smith   PetscFunctionBegin;
273a4121054SBarry Smith   PetscValidPointer(da,1);
274a4121054SBarry Smith   ierr     = PetscNewLog(da,DM_DA,&dd);CHKERRQ(ierr);
275a4121054SBarry Smith   da->data = dd;
27647c6ae99SBarry Smith 
27747c6ae99SBarry Smith   dd->dim        = -1;
278aa219208SBarry Smith   dd->interptype = DMDA_Q1;
27947c6ae99SBarry Smith   dd->refine_x   = 2;
28047c6ae99SBarry Smith   dd->refine_y   = 2;
28147c6ae99SBarry Smith   dd->refine_z   = 2;
28281c108dcSJed Brown   dd->coarsen_x  = 2;
28381c108dcSJed Brown   dd->coarsen_y  = 2;
28481c108dcSJed Brown   dd->coarsen_z  = 2;
28547c6ae99SBarry Smith   dd->fieldname  = PETSC_NULL;
28647c6ae99SBarry Smith   dd->nlocal     = -1;
28747c6ae99SBarry Smith   dd->Nlocal     = -1;
28847c6ae99SBarry Smith   dd->M          = -1;
28947c6ae99SBarry Smith   dd->N          = -1;
29047c6ae99SBarry Smith   dd->P          = -1;
29147c6ae99SBarry Smith   dd->m          = -1;
29247c6ae99SBarry Smith   dd->n          = -1;
29347c6ae99SBarry Smith   dd->p          = -1;
29447c6ae99SBarry Smith   dd->w          = -1;
29547c6ae99SBarry Smith   dd->s          = -1;
296*8865f1eaSKarl Rupp 
29747c6ae99SBarry Smith   dd->xs = -1; dd->xe = -1; dd->ys = -1; dd->ye = -1; dd->zs = -1; dd->ze = -1;
29847c6ae99SBarry Smith   dd->Xs = -1; dd->Xe = -1; dd->Ys = -1; dd->Ye = -1; dd->Zs = -1; dd->Ze = -1;
29947c6ae99SBarry Smith 
3000a696f66SPeter Brune   dd->decompositiondm = PETSC_FALSE;
30188661749SPeter Brune   dd->overlap         = 0;
302d886c4f4SPeter Brune   dd->xo              = 0;
303d886c4f4SPeter Brune   dd->yo              = 0;
304d886c4f4SPeter Brune   dd->zo              = 0;
305e30e807fSPeter Brune   dd->Mo              = -1;
306e30e807fSPeter Brune   dd->No              = -1;
307e30e807fSPeter Brune   dd->Po              = -1;
30888661749SPeter Brune 
30947c6ae99SBarry Smith   dd->gtol         = PETSC_NULL;
31047c6ae99SBarry Smith   dd->ltog         = PETSC_NULL;
31147c6ae99SBarry Smith   dd->ltol         = PETSC_NULL;
31247c6ae99SBarry Smith   dd->ao           = PETSC_NULL;
31347c6ae99SBarry Smith   dd->base         = -1;
3141321219cSEthan Coon   dd->bx           = DMDA_BOUNDARY_NONE;
3151321219cSEthan Coon   dd->by           = DMDA_BOUNDARY_NONE;
3161321219cSEthan Coon   dd->bz           = DMDA_BOUNDARY_NONE;
317aa219208SBarry Smith   dd->stencil_type = DMDA_STENCIL_BOX;
318aa219208SBarry Smith   dd->interptype   = DMDA_Q1;
31947c6ae99SBarry Smith   dd->idx          = PETSC_NULL;
32047c6ae99SBarry Smith   dd->Nl           = -1;
32147c6ae99SBarry Smith   dd->lx           = PETSC_NULL;
32247c6ae99SBarry Smith   dd->ly           = PETSC_NULL;
32347c6ae99SBarry Smith   dd->lz           = PETSC_NULL;
32447c6ae99SBarry Smith 
325454e267fSLisandro Dalcin   dd->elementtype = DMDA_ELEMENT_Q1;
326454e267fSLisandro Dalcin 
32719fd82e9SBarry Smith   ierr = PetscStrallocpy(VECSTANDARD,(char**)&da->vectype);CHKERRQ(ierr);
328*8865f1eaSKarl Rupp 
329a4121054SBarry Smith   da->ops->globaltolocalbegin          = DMGlobalToLocalBegin_DA;
330a4121054SBarry Smith   da->ops->globaltolocalend            = DMGlobalToLocalEnd_DA;
331a4121054SBarry Smith   da->ops->localtoglobalbegin          = DMLocalToGlobalBegin_DA;
332a4121054SBarry Smith   da->ops->localtoglobalend            = DMLocalToGlobalEnd_DA;
333a4121054SBarry Smith   da->ops->createglobalvector          = DMCreateGlobalVector_DA;
334a4121054SBarry Smith   da->ops->createlocalvector           = DMCreateLocalVector_DA;
33525296bd5SBarry Smith   da->ops->createinterpolation         = DMCreateInterpolation_DA;
336e727c939SJed Brown   da->ops->getcoloring                 = DMCreateColoring_DA;
33725296bd5SBarry Smith   da->ops->creatematrix                = DMCreateMatrix_DA;
338a4121054SBarry Smith   da->ops->refine                      = DMRefine_DA;
339a4121054SBarry Smith   da->ops->coarsen                     = DMCoarsen_DA;
340a4121054SBarry Smith   da->ops->refinehierarchy             = DMRefineHierarchy_DA;
341a4121054SBarry Smith   da->ops->coarsenhierarchy            = DMCoarsenHierarchy_DA;
342e727c939SJed Brown   da->ops->getinjection                = DMCreateInjection_DA;
343e727c939SJed Brown   da->ops->getaggregates               = DMCreateAggregates_DA;
344a4121054SBarry Smith   da->ops->destroy                     = DMDestroy_DA;
345a4121054SBarry Smith   da->ops->view                        = 0;
346a4121054SBarry Smith   da->ops->setfromoptions              = DMSetFromOptions_DA;
347a4121054SBarry Smith   da->ops->setup                       = DMSetUp_DA;
348b859378eSBarry Smith   da->ops->load                        = DMLoad_DA;
3496636e97aSMatthew G Knepley   da->ops->createcoordinatedm          = DMCreateCoordinateDM_DA;
350d724dfffSBarry Smith   da->ops->createsubdm                 = DMCreateSubDM_DA;
35116621825SDmitry Karpeev   da->ops->createfielddecomposition    = DMCreateFieldDecomposition_DA;
352e30e807fSPeter Brune   da->ops->createdomaindecomposition   = DMCreateDomainDecomposition_DA;
3530a696f66SPeter Brune   da->ops->createdomaindecompositiondm = DMCreateDomainDecompositionDM_DA;
354e30e807fSPeter Brune   da->ops->createddscatters            = DMCreateDomainDecompositionScatters_DA;
355a4121054SBarry Smith   PetscFunctionReturn(0);
356a4121054SBarry Smith }
357a4121054SBarry Smith EXTERN_C_END
35847c6ae99SBarry Smith 
359a4121054SBarry Smith #undef __FUNCT__
360a4121054SBarry Smith #define __FUNCT__ "DMDACreate"
361a4121054SBarry Smith /*@
362a4121054SBarry Smith   DMDACreate - Creates a DMDA object.
363a4121054SBarry Smith 
364a4121054SBarry Smith   Collective on MPI_Comm
365a4121054SBarry Smith 
366a4121054SBarry Smith   Input Parameter:
367a4121054SBarry Smith . comm - The communicator for the DMDA object
368a4121054SBarry Smith 
369a4121054SBarry Smith   Output Parameter:
370a4121054SBarry Smith . da  - The DMDA object
371a4121054SBarry Smith 
372e0f5d30fSBarry Smith   Level: advanced
373e0f5d30fSBarry Smith 
374e0f5d30fSBarry Smith   Developers Note: Since there exists DMDACreate1/2/3d() should this routine even exist?
375a4121054SBarry Smith 
376a4121054SBarry Smith .keywords: DMDA, create
377e0f5d30fSBarry Smith .seealso:  DMDASetSizes(), DMDADuplicate(),  DMDACreate1d(), DMDACreate2d(), DMDACreate3d()
378a4121054SBarry Smith @*/
3797087cfbeSBarry Smith PetscErrorCode  DMDACreate(MPI_Comm comm, DM *da)
380a4121054SBarry Smith {
381a4121054SBarry Smith   PetscErrorCode ierr;
382a4121054SBarry Smith 
383a4121054SBarry Smith   PetscFunctionBegin;
384a4121054SBarry Smith   PetscValidPointer(da,2);
385a4121054SBarry Smith   ierr = DMCreate(comm,da);CHKERRQ(ierr);
386a4121054SBarry Smith   ierr = DMSetType(*da,DMDA);CHKERRQ(ierr);
38747c6ae99SBarry Smith   PetscFunctionReturn(0);
38847c6ae99SBarry Smith }
389b859378eSBarry Smith 
390b859378eSBarry Smith 
391