1 2 /* 3 Code for manipulating distributed regular arrays in parallel. 4 */ 5 6 #include <petsc/private/dmdaimpl.h> /*I "petscdmda.h" I*/ 7 #include <petscdmfield.h> 8 9 PetscErrorCode DMCreateCoordinateDM_DA(DM dm, DM *cdm) 10 { 11 PetscFunctionBegin; 12 PetscCall(DMDACreateCompatibleDMDA(dm,dm->dim,cdm)); 13 PetscFunctionReturn(0); 14 } 15 16 PetscErrorCode DMCreateCoordinateField_DA(DM dm, DMField *field) 17 { 18 PetscReal gmin[3], gmax[3]; 19 PetscScalar corners[24]; 20 PetscInt dim; 21 PetscInt i, j; 22 DM cdm; 23 24 PetscFunctionBegin; 25 PetscCall(DMGetDimension(dm,&dim)); 26 /* TODO: this is wrong if coordinates are not rectilinear */ 27 PetscCall(DMGetBoundingBox(dm,gmin,gmax)); 28 for (i = 0; i < (1 << dim); i++) { 29 for (j = 0; j < dim; j++) { 30 corners[i*dim + j] = (i & (1 << j)) ? gmax[j] : gmin[j]; 31 } 32 } 33 PetscCall(DMClone(dm,&cdm)); 34 PetscCall(DMFieldCreateDA(cdm,dim,corners,field)); 35 PetscCall(DMDestroy(&cdm)); 36 PetscFunctionReturn(0); 37 } 38 39 /*@C 40 DMDASetFieldName - Sets the names of individual field components in multicomponent 41 vectors associated with a DMDA. 42 43 Logically collective; name must contain a common value 44 45 Input Parameters: 46 + da - the distributed array 47 . nf - field number for the DMDA (0, 1, ... dof-1), where dof indicates the 48 number of degrees of freedom per node within the DMDA 49 - names - the name of the field (component) 50 51 Notes: 52 It must be called after having called DMSetUp(). 53 54 Level: intermediate 55 56 .seealso: `DMDAGetFieldName()`, `DMDASetCoordinateName()`, `DMDAGetCoordinateName()`, `DMDASetFieldNames()`, `DMSetUp()` 57 @*/ 58 PetscErrorCode DMDASetFieldName(DM da,PetscInt nf,const char name[]) 59 { 60 DM_DA *dd = (DM_DA*)da->data; 61 62 PetscFunctionBegin; 63 PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 64 PetscCheck(nf >= 0 && nf < dd->w,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Invalid field number: %" PetscInt_FMT,nf); 65 PetscCheck(dd->fieldname,PetscObjectComm((PetscObject)da),PETSC_ERR_ORDER,"You should call DMSetUp() first"); 66 PetscCall(PetscFree(dd->fieldname[nf])); 67 PetscCall(PetscStrallocpy(name,&dd->fieldname[nf])); 68 PetscFunctionReturn(0); 69 } 70 71 /*@C 72 DMDAGetFieldNames - Gets the name of each component in the vector associated with the DMDA 73 74 Not collective; names will contain a common value 75 76 Input Parameter: 77 . dm - the DMDA object 78 79 Output Parameter: 80 . names - the names of the components, final string is NULL, will have the same number of entries as the dof used in creating the DMDA 81 82 Level: intermediate 83 84 Not supported from Fortran, use DMDAGetFieldName() 85 86 .seealso: `DMDAGetFieldName()`, `DMDASetCoordinateName()`, `DMDAGetCoordinateName()`, `DMDASetFieldName()`, `DMDASetFieldNames()` 87 @*/ 88 PetscErrorCode DMDAGetFieldNames(DM da,const char * const **names) 89 { 90 DM_DA *dd = (DM_DA*)da->data; 91 92 PetscFunctionBegin; 93 *names = (const char * const *) dd->fieldname; 94 PetscFunctionReturn(0); 95 } 96 97 /*@C 98 DMDASetFieldNames - Sets the name of each component in the vector associated with the DMDA 99 100 Logically collective; names must contain a common value 101 102 Input Parameters: 103 + dm - the DMDA object 104 - names - the names of the components, final string must be NULL, must have the same number of entries as the dof used in creating the DMDA 105 106 Notes: 107 It must be called after having called DMSetUp(). 108 109 Level: intermediate 110 111 Not supported from Fortran, use DMDASetFieldName() 112 113 .seealso: `DMDAGetFieldName()`, `DMDASetCoordinateName()`, `DMDAGetCoordinateName()`, `DMDASetFieldName()`, `DMSetUp()` 114 @*/ 115 PetscErrorCode DMDASetFieldNames(DM da,const char * const *names) 116 { 117 DM_DA *dd = (DM_DA*)da->data; 118 char **fieldname; 119 PetscInt nf = 0; 120 121 PetscFunctionBegin; 122 PetscCheck(dd->fieldname,PetscObjectComm((PetscObject)da),PETSC_ERR_ORDER,"You should call DMSetUp() first"); 123 while (names[nf++]) {}; 124 PetscCheck(nf == dd->w+1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Invalid number of fields %" PetscInt_FMT,nf-1); 125 PetscCall(PetscStrArrayallocpy(names,&fieldname)); 126 PetscCall(PetscStrArrayDestroy(&dd->fieldname)); 127 dd->fieldname = fieldname; 128 PetscFunctionReturn(0); 129 } 130 131 /*@C 132 DMDAGetFieldName - Gets the names of individual field components in multicomponent 133 vectors associated with a DMDA. 134 135 Not collective; name will contain a common value 136 137 Input Parameters: 138 + da - the distributed array 139 - nf - field number for the DMDA (0, 1, ... dof-1), where dof indicates the 140 number of degrees of freedom per node within the DMDA 141 142 Output Parameter: 143 . names - the name of the field (component) 144 145 Notes: 146 It must be called after having called DMSetUp(). 147 148 Level: intermediate 149 150 .seealso: `DMDASetFieldName()`, `DMDASetCoordinateName()`, `DMDAGetCoordinateName()`, `DMSetUp()` 151 @*/ 152 PetscErrorCode DMDAGetFieldName(DM da,PetscInt nf,const char **name) 153 { 154 DM_DA *dd = (DM_DA*)da->data; 155 156 PetscFunctionBegin; 157 PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 158 PetscValidPointer(name,3); 159 PetscCheck(nf >= 0 && nf < dd->w,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Invalid field number: %" PetscInt_FMT,nf); 160 PetscCheck(dd->fieldname,PetscObjectComm((PetscObject)da),PETSC_ERR_ORDER,"You should call DMSetUp() first"); 161 *name = dd->fieldname[nf]; 162 PetscFunctionReturn(0); 163 } 164 165 /*@C 166 DMDASetCoordinateName - Sets the name of the coordinate directions associated with a DMDA, for example "x" or "y" 167 168 Logically collective; name must contain a common value 169 170 Input Parameters: 171 + dm - the DM 172 . nf - coordinate number for the DMDA (0, 1, ... dim-1), 173 - name - the name of the coordinate 174 175 Notes: 176 It must be called after having called DMSetUp(). 177 178 Level: intermediate 179 180 Not supported from Fortran 181 182 .seealso: `DMDAGetCoordinateName()`, `DMDASetFieldName()`, `DMDAGetFieldName()`, `DMSetUp()` 183 @*/ 184 PetscErrorCode DMDASetCoordinateName(DM dm,PetscInt nf,const char name[]) 185 { 186 DM_DA *dd = (DM_DA*)dm->data; 187 188 PetscFunctionBegin; 189 PetscValidHeaderSpecificType(dm,DM_CLASSID,1,DMDA); 190 PetscCheck(nf >= 0 && nf < dm->dim,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Invalid coordinate number: %" PetscInt_FMT,nf); 191 PetscCheck(dd->coordinatename,PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"You should call DMSetUp() first"); 192 PetscCall(PetscFree(dd->coordinatename[nf])); 193 PetscCall(PetscStrallocpy(name,&dd->coordinatename[nf])); 194 PetscFunctionReturn(0); 195 } 196 197 /*@C 198 DMDAGetCoordinateName - Gets the name of a coodinate direction associated with a DMDA. 199 200 Not collective; name will contain a common value 201 202 Input Parameters: 203 + dm - the DM 204 - nf - number for the DMDA (0, 1, ... dim-1) 205 206 Output Parameter: 207 . names - the name of the coordinate direction 208 209 Notes: 210 It must be called after having called DMSetUp(). 211 212 Level: intermediate 213 214 Not supported from Fortran 215 216 .seealso: `DMDASetCoordinateName()`, `DMDASetFieldName()`, `DMDAGetFieldName()`, `DMSetUp()` 217 @*/ 218 PetscErrorCode DMDAGetCoordinateName(DM dm,PetscInt nf,const char **name) 219 { 220 DM_DA *dd = (DM_DA*)dm->data; 221 222 PetscFunctionBegin; 223 PetscValidHeaderSpecificType(dm,DM_CLASSID,1,DMDA); 224 PetscValidPointer(name,3); 225 PetscCheck(nf >= 0 && nf < dm->dim,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Invalid coordinate number: %" PetscInt_FMT,nf); 226 PetscCheck(dd->coordinatename,PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"You should call DMSetUp() first"); 227 *name = dd->coordinatename[nf]; 228 PetscFunctionReturn(0); 229 } 230 231 /*@C 232 DMDAGetCorners - Returns the global (x,y,z) indices of the lower left 233 corner and size of the local region, excluding ghost points. 234 235 Not collective 236 237 Input Parameter: 238 . da - the distributed array 239 240 Output Parameters: 241 + x - the corner index for the first dimension 242 . y - the corner index for the second dimension (only used in 2D and 3D problems) 243 . z - the corner index for the third dimension (only used in 3D problems) 244 . m - the width in the first dimension 245 . n - the width in the second dimension (only used in 2D and 3D problems) 246 - p - the width in the third dimension (only used in 3D problems) 247 248 Note: 249 The corner information is independent of the number of degrees of 250 freedom per node set with the DMDACreateXX() routine. Thus the x, y, z, and 251 m, n, p can be thought of as coordinates on a logical grid, where each 252 grid point has (potentially) several degrees of freedom. 253 Any of y, z, n, and p can be passed in as NULL if not needed. 254 255 Level: beginner 256 257 .seealso: `DMDAGetGhostCorners()`, `DMDAGetOwnershipRanges()`, `DMStagGetCorners()` 258 @*/ 259 PetscErrorCode DMDAGetCorners(DM da,PetscInt *x,PetscInt *y,PetscInt *z,PetscInt *m,PetscInt *n,PetscInt *p) 260 { 261 PetscInt w; 262 DM_DA *dd = (DM_DA*)da->data; 263 264 PetscFunctionBegin; 265 PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 266 /* since the xs, xe ... have all been multiplied by the number of degrees 267 of freedom per cell, w = dd->w, we divide that out before returning.*/ 268 w = dd->w; 269 if (x) *x = dd->xs/w + dd->xo; 270 /* the y and z have NOT been multiplied by w */ 271 if (y) *y = dd->ys + dd->yo; 272 if (z) *z = dd->zs + dd->zo; 273 if (m) *m = (dd->xe - dd->xs)/w; 274 if (n) *n = (dd->ye - dd->ys); 275 if (p) *p = (dd->ze - dd->zs); 276 PetscFunctionReturn(0); 277 } 278 279 PetscErrorCode DMGetLocalBoundingIndices_DMDA(DM dm, PetscReal lmin[], PetscReal lmax[]) 280 { 281 DMDALocalInfo info; 282 283 PetscFunctionBegin; 284 PetscCall(DMDAGetLocalInfo(dm, &info)); 285 lmin[0] = info.xs; 286 lmin[1] = info.ys; 287 lmin[2] = info.zs; 288 lmax[0] = info.xs + info.xm-1; 289 lmax[1] = info.ys + info.ym-1; 290 lmax[2] = info.zs + info.zm-1; 291 PetscFunctionReturn(0); 292 } 293 294 /*@ 295 DMDAGetReducedDMDA - Deprecated; use DMDACreateCompatibleDMDA() 296 297 Level: deprecated 298 @*/ 299 PetscErrorCode DMDAGetReducedDMDA(DM da,PetscInt nfields,DM *nda) 300 { 301 PetscFunctionBegin; 302 PetscCall(DMDACreateCompatibleDMDA(da,nfields,nda)); 303 PetscFunctionReturn(0); 304 } 305 306 /*@ 307 DMDACreateCompatibleDMDA - Creates a DMDA with the same layout but with fewer or more fields 308 309 Collective 310 311 Input Parameters: 312 + da - the distributed array 313 - nfields - number of fields in new DMDA 314 315 Output Parameter: 316 . nda - the new DMDA 317 318 Level: intermediate 319 320 .seealso: `DMDAGetGhostCorners()`, `DMSetCoordinates()`, `DMDASetUniformCoordinates()`, `DMGetCoordinates()`, `DMDAGetGhostedCoordinates()`, `DMStagCreateCompatibleDMStag()` 321 @*/ 322 PetscErrorCode DMDACreateCompatibleDMDA(DM da,PetscInt nfields,DM *nda) 323 { 324 DM_DA *dd = (DM_DA*)da->data; 325 PetscInt s,m,n,p,M,N,P,dim,Mo,No,Po; 326 const PetscInt *lx,*ly,*lz; 327 DMBoundaryType bx,by,bz; 328 DMDAStencilType stencil_type; 329 Vec coords; 330 PetscInt ox,oy,oz; 331 PetscInt cl,rl; 332 333 PetscFunctionBegin; 334 dim = da->dim; 335 M = dd->M; 336 N = dd->N; 337 P = dd->P; 338 m = dd->m; 339 n = dd->n; 340 p = dd->p; 341 s = dd->s; 342 bx = dd->bx; 343 by = dd->by; 344 bz = dd->bz; 345 346 stencil_type = dd->stencil_type; 347 348 PetscCall(DMDAGetOwnershipRanges(da,&lx,&ly,&lz)); 349 if (dim == 1) { 350 PetscCall(DMDACreate1d(PetscObjectComm((PetscObject)da),bx,M,nfields,s,dd->lx,nda)); 351 } else if (dim == 2) { 352 PetscCall(DMDACreate2d(PetscObjectComm((PetscObject)da),bx,by,stencil_type,M,N,m,n,nfields,s,lx,ly,nda)); 353 } else if (dim == 3) { 354 PetscCall(DMDACreate3d(PetscObjectComm((PetscObject)da),bx,by,bz,stencil_type,M,N,P,m,n,p,nfields,s,lx,ly,lz,nda)); 355 } 356 PetscCall(DMSetUp(*nda)); 357 PetscCall(DMGetCoordinates(da, &coords)); 358 PetscCall(DMSetCoordinates(*nda, coords)); 359 360 /* allow for getting a reduced DA corresponding to a domain decomposition */ 361 PetscCall(DMDAGetOffset(da,&ox,&oy,&oz,&Mo,&No,&Po)); 362 PetscCall(DMDASetOffset(*nda,ox,oy,oz,Mo,No,Po)); 363 364 /* allow for getting a reduced DA corresponding to a coarsened DA */ 365 PetscCall(DMGetCoarsenLevel(da,&cl)); 366 PetscCall(DMGetRefineLevel(da,&rl)); 367 368 (*nda)->levelup = rl; 369 (*nda)->leveldown = cl; 370 PetscFunctionReturn(0); 371 } 372 373 /*@C 374 DMDAGetCoordinateArray - Gets an array containing the coordinates of the DMDA 375 376 Not collective 377 378 Input Parameter: 379 . dm - the DM 380 381 Output Parameter: 382 . xc - the coordinates 383 384 Level: intermediate 385 386 Not supported from Fortran 387 388 .seealso: `DMDASetCoordinateName()`, `DMDASetFieldName()`, `DMDAGetFieldName()`, `DMDARestoreCoordinateArray()` 389 @*/ 390 PetscErrorCode DMDAGetCoordinateArray(DM dm,void *xc) 391 { 392 DM cdm; 393 Vec x; 394 395 PetscFunctionBegin; 396 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 397 PetscCall(DMGetCoordinates(dm,&x)); 398 PetscCall(DMGetCoordinateDM(dm,&cdm)); 399 PetscCall(DMDAVecGetArray(cdm,x,xc)); 400 PetscFunctionReturn(0); 401 } 402 403 /*@C 404 DMDARestoreCoordinateArray - Sets an array containing the coordinates of the DMDA 405 406 Not collective 407 408 Input Parameters: 409 + dm - the DM 410 - xc - the coordinates 411 412 Level: intermediate 413 414 Not supported from Fortran 415 416 .seealso: `DMDASetCoordinateName()`, `DMDASetFieldName()`, `DMDAGetFieldName()`, `DMDAGetCoordinateArray()` 417 @*/ 418 PetscErrorCode DMDARestoreCoordinateArray(DM dm,void *xc) 419 { 420 DM cdm; 421 Vec x; 422 423 PetscFunctionBegin; 424 PetscValidHeaderSpecific(dm,DM_CLASSID,1); 425 PetscCall(DMGetCoordinates(dm,&x)); 426 PetscCall(DMGetCoordinateDM(dm,&cdm)); 427 PetscCall(DMDAVecRestoreArray(cdm,x,xc)); 428 PetscFunctionReturn(0); 429 } 430