1 2 #include <petsc-private/dmdaimpl.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,flg; 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,NULL);CHKERRQ(ierr);} 34 if (bN) {ierr = PetscOptionsInt("-da_grid_y","Number of grid points in y direction","DMDASetSizes",dd->N,&dd->N,NULL);CHKERRQ(ierr);} 35 if (bP) {ierr = PetscOptionsInt("-da_grid_z","Number of grid points in z direction","DMDASetSizes",dd->P,&dd->P,NULL);CHKERRQ(ierr);} 36 37 ierr = PetscOptionsInt("-da_overlap","Decomposition overlap in all directions","DMDASetOverlap",dd->xol,&dd->xol,&flg);CHKERRQ(ierr); 38 if (flg) {ierr = DMDASetOverlap(da,dd->xol,dd->xol,dd->xol);CHKERRQ(ierr);} 39 ierr = PetscOptionsInt("-da_overlap_x","Decomposition overlap in x direction","DMDASetOverlap",dd->xol,&dd->xol,NULL);CHKERRQ(ierr); 40 if (dd->dim > 1) {ierr = PetscOptionsInt("-da_overlap_y","Decomposition overlap in y direction","DMDASetOverlap",dd->yol,&dd->yol,NULL);CHKERRQ(ierr);} 41 if (dd->dim > 2) {ierr = PetscOptionsInt("-da_overlap_z","Decomposition overlap in z direction","DMDASetOverlap",dd->zol,&dd->zol,NULL);CHKERRQ(ierr);} 42 43 ierr = PetscOptionsInt("-da_local_subdomains","","DMDASetNumLocalSubdomains",dd->Nsub,&dd->Nsub,&flg);CHKERRQ(ierr); 44 if (flg) {ierr = DMDASetNumLocalSubDomains(da,dd->Nsub);CHKERRQ(ierr);} 45 46 /* Handle DMDA parallel distibution */ 47 ierr = PetscOptionsInt("-da_processors_x","Number of processors in x direction","DMDASetNumProcs",dd->m,&dd->m,NULL);CHKERRQ(ierr); 48 if (dd->dim > 1) {ierr = PetscOptionsInt("-da_processors_y","Number of processors in y direction","DMDASetNumProcs",dd->n,&dd->n,NULL);CHKERRQ(ierr);} 49 if (dd->dim > 2) {ierr = PetscOptionsInt("-da_processors_z","Number of processors in z direction","DMDASetNumProcs",dd->p,&dd->p,NULL);CHKERRQ(ierr);} 50 /* Handle DMDA refinement */ 51 ierr = PetscOptionsInt("-da_refine_x","Refinement ratio in x direction","DMDASetRefinementFactor",dd->refine_x,&dd->refine_x,NULL);CHKERRQ(ierr); 52 if (dd->dim > 1) {ierr = PetscOptionsInt("-da_refine_y","Refinement ratio in y direction","DMDASetRefinementFactor",dd->refine_y,&dd->refine_y,NULL);CHKERRQ(ierr);} 53 if (dd->dim > 2) {ierr = PetscOptionsInt("-da_refine_z","Refinement ratio in z direction","DMDASetRefinementFactor",dd->refine_z,&dd->refine_z,NULL);CHKERRQ(ierr);} 54 dd->coarsen_x = dd->refine_x; dd->coarsen_y = dd->refine_y; dd->coarsen_z = dd->refine_z; 55 56 /* Get refinement factors, defaults taken from the coarse DMDA */ 57 ierr = PetscMalloc3(maxnlevels,PetscInt,&refx,maxnlevels,PetscInt,&refy,maxnlevels,PetscInt,&refz);CHKERRQ(ierr); 58 ierr = DMDAGetRefinementFactor(da,&refx[0],&refy[0],&refz[0]);CHKERRQ(ierr); 59 for (i=1; i<maxnlevels; i++) { 60 refx[i] = refx[0]; 61 refy[i] = refy[0]; 62 refz[i] = refz[0]; 63 } 64 n = maxnlevels; 65 ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_x",refx,&n,NULL);CHKERRQ(ierr); 66 if (da->levelup - da->leveldown >= 0) dd->refine_x = refx[da->levelup - da->leveldown]; 67 if (da->levelup - da->leveldown >= 1) dd->coarsen_x = refx[da->levelup - da->leveldown - 1]; 68 if (dd->dim > 1) { 69 n = maxnlevels; 70 ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_y",refy,&n,NULL);CHKERRQ(ierr); 71 if (da->levelup - da->leveldown >= 0) dd->refine_y = refy[da->levelup - da->leveldown]; 72 if (da->levelup - da->leveldown >= 1) dd->coarsen_y = refy[da->levelup - da->leveldown - 1]; 73 } 74 if (dd->dim > 2) { 75 n = maxnlevels; 76 ierr = PetscOptionsGetIntArray(((PetscObject)da)->prefix,"-da_refine_hierarchy_z",refz,&n,NULL);CHKERRQ(ierr); 77 if (da->levelup - da->leveldown >= 0) dd->refine_z = refz[da->levelup - da->leveldown]; 78 if (da->levelup - da->leveldown >= 1) dd->coarsen_z = refz[da->levelup - da->leveldown - 1]; 79 } 80 81 if (negativeMNP) {ierr = PetscOptionsInt("-da_refine","Uniformly refine DA one or more times","None",refine,&refine,NULL);CHKERRQ(ierr);} 82 ierr = PetscOptionsTail();CHKERRQ(ierr); 83 84 while (refine--) { 85 if (dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 86 dd->M = dd->refine_x*dd->M; 87 } else { 88 dd->M = 1 + dd->refine_x*(dd->M - 1); 89 } 90 if (dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 91 dd->N = dd->refine_y*dd->N; 92 } else { 93 dd->N = 1 + dd->refine_y*(dd->N - 1); 94 } 95 if (dd->bz == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0) { 96 dd->P = dd->refine_z*dd->P; 97 } else { 98 dd->P = 1 + dd->refine_z*(dd->P - 1); 99 } 100 da->levelup++; 101 if (da->levelup - da->leveldown >= 0) { 102 dd->refine_x = refx[da->levelup - da->leveldown]; 103 dd->refine_y = refy[da->levelup - da->leveldown]; 104 dd->refine_z = refz[da->levelup - da->leveldown]; 105 } 106 if (da->levelup - da->leveldown >= 1) { 107 dd->coarsen_x = refx[da->levelup - da->leveldown - 1]; 108 dd->coarsen_y = refy[da->levelup - da->leveldown - 1]; 109 dd->coarsen_z = refz[da->levelup - da->leveldown - 1]; 110 } 111 } 112 ierr = PetscFree3(refx,refy,refz);CHKERRQ(ierr); 113 PetscFunctionReturn(0); 114 } 115 116 extern PetscErrorCode DMCreateGlobalVector_DA(DM,Vec*); 117 extern PetscErrorCode DMCreateLocalVector_DA(DM,Vec*); 118 extern PetscErrorCode DMGlobalToLocalBegin_DA(DM,Vec,InsertMode,Vec); 119 extern PetscErrorCode DMGlobalToLocalEnd_DA(DM,Vec,InsertMode,Vec); 120 extern PetscErrorCode DMLocalToGlobalBegin_DA(DM,Vec,InsertMode,Vec); 121 extern PetscErrorCode DMLocalToGlobalEnd_DA(DM,Vec,InsertMode,Vec); 122 extern PetscErrorCode DMLocalToLocalBegin_DA(DM,Vec,InsertMode,Vec); 123 extern PetscErrorCode DMLocalToLocalEnd_DA(DM,Vec,InsertMode,Vec); 124 extern PetscErrorCode DMCreateInterpolation_DA(DM,DM,Mat*,Vec*); 125 extern PetscErrorCode DMCreateColoring_DA(DM,ISColoringType,MatType,ISColoring*); 126 extern PetscErrorCode DMCreateMatrix_DA(DM,MatType,Mat*); 127 extern PetscErrorCode DMCreateCoordinateDM_DA(DM,DM*); 128 extern PetscErrorCode DMRefine_DA(DM,MPI_Comm,DM*); 129 extern PetscErrorCode DMCoarsen_DA(DM,MPI_Comm,DM*); 130 extern PetscErrorCode DMRefineHierarchy_DA(DM,PetscInt,DM[]); 131 extern PetscErrorCode DMCoarsenHierarchy_DA(DM,PetscInt,DM[]); 132 extern PetscErrorCode DMCreateInjection_DA(DM,DM,VecScatter*); 133 extern PetscErrorCode DMCreateAggregates_DA(DM,DM,Mat*); 134 extern PetscErrorCode DMView_DA(DM,PetscViewer); 135 extern PetscErrorCode DMSetUp_DA(DM); 136 extern PetscErrorCode DMDestroy_DA(DM); 137 extern PetscErrorCode DMCreateDomainDecomposition_DA(DM,PetscInt*,char***,IS**,IS**,DM**); 138 extern PetscErrorCode DMCreateDomainDecompositionScatters_DA(DM,PetscInt,DM*,VecScatter**,VecScatter**,VecScatter**); 139 140 #undef __FUNCT__ 141 #define __FUNCT__ "DMLoad_DA" 142 PetscErrorCode DMLoad_DA(DM da,PetscViewer viewer) 143 { 144 PetscErrorCode ierr; 145 PetscInt dim,m,n,p,dof,swidth; 146 DMDAStencilType stencil; 147 DMDABoundaryType bx,by,bz; 148 PetscBool coors; 149 DM dac; 150 Vec c; 151 152 PetscFunctionBegin; 153 ierr = PetscViewerBinaryRead(viewer,&dim,1,PETSC_INT);CHKERRQ(ierr); 154 ierr = PetscViewerBinaryRead(viewer,&m,1,PETSC_INT);CHKERRQ(ierr); 155 ierr = PetscViewerBinaryRead(viewer,&n,1,PETSC_INT);CHKERRQ(ierr); 156 ierr = PetscViewerBinaryRead(viewer,&p,1,PETSC_INT);CHKERRQ(ierr); 157 ierr = PetscViewerBinaryRead(viewer,&dof,1,PETSC_INT);CHKERRQ(ierr); 158 ierr = PetscViewerBinaryRead(viewer,&swidth,1,PETSC_INT);CHKERRQ(ierr); 159 ierr = PetscViewerBinaryRead(viewer,&bx,1,PETSC_ENUM);CHKERRQ(ierr); 160 ierr = PetscViewerBinaryRead(viewer,&by,1,PETSC_ENUM);CHKERRQ(ierr); 161 ierr = PetscViewerBinaryRead(viewer,&bz,1,PETSC_ENUM);CHKERRQ(ierr); 162 ierr = PetscViewerBinaryRead(viewer,&stencil,1,PETSC_ENUM);CHKERRQ(ierr); 163 164 ierr = DMDASetDim(da, dim);CHKERRQ(ierr); 165 ierr = DMDASetSizes(da, m,n,p);CHKERRQ(ierr); 166 ierr = DMDASetBoundaryType(da, bx, by, bz);CHKERRQ(ierr); 167 ierr = DMDASetDof(da, dof);CHKERRQ(ierr); 168 ierr = DMDASetStencilType(da, stencil);CHKERRQ(ierr); 169 ierr = DMDASetStencilWidth(da, swidth);CHKERRQ(ierr); 170 ierr = DMSetUp(da);CHKERRQ(ierr); 171 ierr = PetscViewerBinaryRead(viewer,&coors,1,PETSC_ENUM);CHKERRQ(ierr); 172 if (coors) { 173 ierr = DMGetCoordinateDM(da,&dac);CHKERRQ(ierr); 174 ierr = DMCreateGlobalVector(dac,&c);CHKERRQ(ierr); 175 ierr = VecLoad(c,viewer);CHKERRQ(ierr); 176 ierr = DMSetCoordinates(da,c);CHKERRQ(ierr); 177 ierr = VecDestroy(&c);CHKERRQ(ierr); 178 } 179 PetscFunctionReturn(0); 180 } 181 182 #undef __FUNCT__ 183 #define __FUNCT__ "DMCreateSubDM_DA" 184 PetscErrorCode DMCreateSubDM_DA(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm) 185 { 186 DM_DA *da = (DM_DA*) dm->data; 187 PetscSection section; 188 PetscErrorCode ierr; 189 190 PetscFunctionBegin; 191 if (subdm) { 192 ierr = DMClone(dm, subdm);CHKERRQ(ierr); 193 ierr = DMDASetDof(*subdm, numFields);CHKERRQ(ierr); 194 } 195 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 196 if (section) { 197 ierr = DMCreateSubDM_Section_Private(dm, numFields, fields, is, subdm);CHKERRQ(ierr); 198 } else { 199 if (is) { 200 PetscInt *indices, cnt = 0, dof = da->w, i, j; 201 202 ierr = PetscMalloc(da->Nlocal*numFields/dof * sizeof(PetscInt), &indices);CHKERRQ(ierr); 203 for (i = da->base/dof; i < (da->base+da->Nlocal)/dof; ++i) { 204 for (j = 0; j < numFields; ++j) { 205 indices[cnt++] = dof*i + fields[j]; 206 } 207 } 208 if (cnt != da->Nlocal*numFields/dof) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Count %d does not equal expected value %d", cnt, da->Nlocal*numFields/dof); 209 ierr = ISCreateGeneral(PetscObjectComm((PetscObject) dm), cnt, indices, PETSC_OWN_POINTER, is);CHKERRQ(ierr); 210 } 211 } 212 PetscFunctionReturn(0); 213 } 214 215 #undef __FUNCT__ 216 #define __FUNCT__ "DMCreateFieldDecomposition_DA" 217 PetscErrorCode DMCreateFieldDecomposition_DA(DM dm, PetscInt *len,char ***namelist, IS **islist, DM **dmlist) 218 { 219 PetscInt i; 220 PetscErrorCode ierr; 221 DM_DA *dd = (DM_DA*)dm->data; 222 PetscInt dof = dd->w; 223 224 PetscFunctionBegin; 225 if (len) *len = dof; 226 if (islist) { 227 Vec v; 228 PetscInt rstart,n; 229 230 ierr = DMGetGlobalVector(dm,&v);CHKERRQ(ierr); 231 ierr = VecGetOwnershipRange(v,&rstart,NULL);CHKERRQ(ierr); 232 ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr); 233 ierr = DMRestoreGlobalVector(dm,&v);CHKERRQ(ierr); 234 ierr = PetscMalloc(dof*sizeof(IS),islist);CHKERRQ(ierr); 235 for (i=0; i<dof; i++) { 236 ierr = ISCreateStride(PetscObjectComm((PetscObject)dm),n/dof,rstart+i,dof,&(*islist)[i]);CHKERRQ(ierr); 237 } 238 } 239 if (namelist) { 240 ierr = PetscMalloc(dof*sizeof(const char*), namelist);CHKERRQ(ierr); 241 if (dd->fieldname) { 242 for (i=0; i<dof; i++) { 243 ierr = PetscStrallocpy(dd->fieldname[i],&(*namelist)[i]);CHKERRQ(ierr); 244 } 245 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Currently DMDA must have fieldnames"); 246 } 247 if (dmlist) { 248 DM da; 249 250 ierr = DMDACreate(PetscObjectComm((PetscObject)dm), &da);CHKERRQ(ierr); 251 ierr = DMDASetDim(da, dd->dim);CHKERRQ(ierr); 252 ierr = DMDASetSizes(da, dd->M, dd->N, dd->P);CHKERRQ(ierr); 253 ierr = DMDASetNumProcs(da, dd->m, dd->n, dd->p);CHKERRQ(ierr); 254 ierr = DMDASetBoundaryType(da, dd->bx, dd->by, dd->bz);CHKERRQ(ierr); 255 ierr = DMDASetDof(da, 1);CHKERRQ(ierr); 256 ierr = DMDASetStencilType(da, dd->stencil_type);CHKERRQ(ierr); 257 ierr = DMDASetStencilWidth(da, dd->s);CHKERRQ(ierr); 258 ierr = DMSetUp(da);CHKERRQ(ierr); 259 ierr = PetscMalloc(dof*sizeof(DM),dmlist);CHKERRQ(ierr); 260 for (i=0; i<dof-1; i++) {ierr = PetscObjectReference((PetscObject)da);CHKERRQ(ierr);} 261 for (i=0; i<dof; i++) (*dmlist)[i] = da; 262 } 263 PetscFunctionReturn(0); 264 } 265 266 #undef __FUNCT__ 267 #define __FUNCT__ "DMClone_DA" 268 PetscErrorCode DMClone_DA(DM dm, DM *newdm) 269 { 270 DM_DA *da = (DM_DA *) dm->data; 271 PetscErrorCode ierr; 272 273 PetscFunctionBegin; 274 ierr = DMSetType(*newdm, DMDA);CHKERRQ(ierr); 275 ierr = DMDASetDim(*newdm, da->dim);CHKERRQ(ierr); 276 ierr = DMDASetSizes(*newdm, da->M, da->N, da->P);CHKERRQ(ierr); 277 ierr = DMDASetNumProcs(*newdm, da->m, da->n, da->p);CHKERRQ(ierr); 278 ierr = DMDASetBoundaryType(*newdm, da->bx, da->by, da->bz);CHKERRQ(ierr); 279 ierr = DMDASetDof(*newdm, da->w);CHKERRQ(ierr); 280 ierr = DMDASetStencilType(*newdm, da->stencil_type);CHKERRQ(ierr); 281 ierr = DMDASetStencilWidth(*newdm, da->s);CHKERRQ(ierr); 282 ierr = DMDASetOwnershipRanges(*newdm, da->lx, da->ly, da->lz);CHKERRQ(ierr); 283 PetscFunctionReturn(0); 284 } 285 286 /*MC 287 DMDA = "da" - A DM object that is used to manage data for a structured grid in 1, 2, or 3 dimensions. 288 In the global representation of the vector each process stores a non-overlapping rectangular (or slab in 3d) portion of the grid points. 289 In the local representation these rectangular regions (slabs) are extended in all directions by a stencil width. 290 291 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 292 vertex centered. 293 294 295 Level: intermediate 296 297 .seealso: DMType, DMCOMPOSITE, DMDACreate(), DMCreate(), DMSetType() 298 M*/ 299 300 301 #undef __FUNCT__ 302 #define __FUNCT__ "DMCreate_DA" 303 PETSC_EXTERN PetscErrorCode DMCreate_DA(DM da) 304 { 305 PetscErrorCode ierr; 306 DM_DA *dd; 307 308 PetscFunctionBegin; 309 PetscValidPointer(da,1); 310 ierr = PetscNewLog(da,DM_DA,&dd);CHKERRQ(ierr); 311 da->data = dd; 312 313 dd->dim = -1; 314 dd->interptype = DMDA_Q1; 315 dd->refine_x = 2; 316 dd->refine_y = 2; 317 dd->refine_z = 2; 318 dd->coarsen_x = 2; 319 dd->coarsen_y = 2; 320 dd->coarsen_z = 2; 321 dd->fieldname = NULL; 322 dd->nlocal = -1; 323 dd->Nlocal = -1; 324 dd->M = -1; 325 dd->N = -1; 326 dd->P = -1; 327 dd->m = -1; 328 dd->n = -1; 329 dd->p = -1; 330 dd->w = -1; 331 dd->s = -1; 332 333 dd->xs = -1; dd->xe = -1; dd->ys = -1; dd->ye = -1; dd->zs = -1; dd->ze = -1; 334 dd->Xs = -1; dd->Xe = -1; dd->Ys = -1; dd->Ye = -1; dd->Zs = -1; dd->Ze = -1; 335 336 dd->Nsub = 1; 337 dd->xol = 0; 338 dd->yol = 0; 339 dd->zol = 0; 340 dd->xo = 0; 341 dd->yo = 0; 342 dd->zo = 0; 343 dd->Mo = -1; 344 dd->No = -1; 345 dd->Po = -1; 346 347 dd->gtol = NULL; 348 dd->ltog = NULL; 349 dd->ltol = NULL; 350 dd->ao = NULL; 351 dd->base = -1; 352 dd->bx = DMDA_BOUNDARY_NONE; 353 dd->by = DMDA_BOUNDARY_NONE; 354 dd->bz = DMDA_BOUNDARY_NONE; 355 dd->stencil_type = DMDA_STENCIL_BOX; 356 dd->interptype = DMDA_Q1; 357 dd->idx = NULL; 358 dd->Nl = -1; 359 dd->lx = NULL; 360 dd->ly = NULL; 361 dd->lz = NULL; 362 363 dd->elementtype = DMDA_ELEMENT_Q1; 364 365 ierr = PetscStrallocpy(VECSTANDARD,(char**)&da->vectype);CHKERRQ(ierr); 366 367 da->ops->globaltolocalbegin = DMGlobalToLocalBegin_DA; 368 da->ops->globaltolocalend = DMGlobalToLocalEnd_DA; 369 da->ops->localtoglobalbegin = DMLocalToGlobalBegin_DA; 370 da->ops->localtoglobalend = DMLocalToGlobalEnd_DA; 371 da->ops->localtolocalbegin = DMLocalToLocalBegin_DA; 372 da->ops->localtolocalend = DMLocalToLocalEnd_DA; 373 da->ops->createglobalvector = DMCreateGlobalVector_DA; 374 da->ops->createlocalvector = DMCreateLocalVector_DA; 375 da->ops->createinterpolation = DMCreateInterpolation_DA; 376 da->ops->getcoloring = DMCreateColoring_DA; 377 da->ops->creatematrix = DMCreateMatrix_DA; 378 da->ops->refine = DMRefine_DA; 379 da->ops->coarsen = DMCoarsen_DA; 380 da->ops->refinehierarchy = DMRefineHierarchy_DA; 381 da->ops->coarsenhierarchy = DMCoarsenHierarchy_DA; 382 da->ops->getinjection = DMCreateInjection_DA; 383 da->ops->getaggregates = DMCreateAggregates_DA; 384 da->ops->destroy = DMDestroy_DA; 385 da->ops->view = 0; 386 da->ops->setfromoptions = DMSetFromOptions_DA; 387 da->ops->setup = DMSetUp_DA; 388 da->ops->clone = DMClone_DA; 389 da->ops->load = DMLoad_DA; 390 da->ops->createcoordinatedm = DMCreateCoordinateDM_DA; 391 da->ops->createsubdm = DMCreateSubDM_DA; 392 da->ops->createfielddecomposition = DMCreateFieldDecomposition_DA; 393 da->ops->createdomaindecomposition = DMCreateDomainDecomposition_DA; 394 da->ops->createddscatters = DMCreateDomainDecompositionScatters_DA; 395 PetscFunctionReturn(0); 396 } 397 398 #undef __FUNCT__ 399 #define __FUNCT__ "DMDACreate" 400 /*@ 401 DMDACreate - Creates a DMDA object. 402 403 Collective on MPI_Comm 404 405 Input Parameter: 406 . comm - The communicator for the DMDA object 407 408 Output Parameter: 409 . da - The DMDA object 410 411 Level: advanced 412 413 Developers Note: Since there exists DMDACreate1/2/3d() should this routine even exist? 414 415 .keywords: DMDA, create 416 .seealso: DMDASetSizes(), DMDADuplicate(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d() 417 @*/ 418 PetscErrorCode DMDACreate(MPI_Comm comm, DM *da) 419 { 420 PetscErrorCode ierr; 421 422 PetscFunctionBegin; 423 PetscValidPointer(da,2); 424 ierr = DMCreate(comm,da);CHKERRQ(ierr); 425 ierr = DMSetType(*da,DMDA);CHKERRQ(ierr); 426 PetscFunctionReturn(0); 427 } 428 429 430