Actual source code: daindex.c
petsc-3.10.5 2019-03-28
2: /*
3: Code for manipulating distributed regular arrays in parallel.
4: */
6: #include <petsc/private/dmdaimpl.h>
8: /*
9: Gets the natural number for each global number on the process.
11: Used by DMDAGetAO() and DMDAGlobalToNatural_Create()
12: */
13: PetscErrorCode DMDAGetNatural_Private(DM da,PetscInt *outNlocal,IS *isnatural)
14: {
16: PetscInt Nlocal,i,j,k,*lidx,lict = 0,dim = da->dim;
17: DM_DA *dd = (DM_DA*)da->data;
20: Nlocal = (dd->xe-dd->xs);
21: if (dim > 1) Nlocal *= (dd->ye-dd->ys);
22: if (dim > 2) Nlocal *= (dd->ze-dd->zs);
24: PetscMalloc1(Nlocal,&lidx);
26: if (dim == 1) {
27: for (i=dd->xs; i<dd->xe; i++) {
28: /* global number in natural ordering */
29: lidx[lict++] = i;
30: }
31: } else if (dim == 2) {
32: for (j=dd->ys; j<dd->ye; j++) {
33: for (i=dd->xs; i<dd->xe; i++) {
34: /* global number in natural ordering */
35: lidx[lict++] = i + j*dd->M*dd->w;
36: }
37: }
38: } else if (dim == 3) {
39: for (k=dd->zs; k<dd->ze; k++) {
40: for (j=dd->ys; j<dd->ye; j++) {
41: for (i=dd->xs; i<dd->xe; i++) {
42: lidx[lict++] = i + j*dd->M*dd->w + k*dd->M*dd->N*dd->w;
43: }
44: }
45: }
46: }
47: *outNlocal = Nlocal;
48: ISCreateGeneral(PetscObjectComm((PetscObject)da),Nlocal,lidx,PETSC_OWN_POINTER,isnatural);
49: return(0);
50: }
52: /*@C
53: DMDASetAOType - Sets the type of application ordering for a distributed array.
55: Collective on DMDA
57: Input Parameter:
58: . da - the distributed array
59: . aotype - type of AO
61: Output Parameters:
63: Level: intermediate
65: Notes:
66: It will generate and error if an AO has already been obtained with a call to DMDAGetAO and the user sets a different AOType
68: .keywords: distributed array, get, global, indices, local-to-global
70: .seealso: DMDACreate2d(), DMDAGetAO(), DMDAGetGhostCorners(), DMDAGetCorners(), DMDALocalToGlocal()
71: DMGlobalToLocalBegin(), DMGlobalToLocalEnd(), DMLocalToLocalBegin(), DMLocalToLocalEnd(), DMDAGetGlobalIndices(), DMDAGetOwnershipRanges(),
72: AO, AOPetscToApplication(), AOApplicationToPetsc()
73: @*/
74: PetscErrorCode DMDASetAOType(DM da,AOType aotype)
75: {
76: DM_DA *dd;
77: PetscBool isdmda;
82: PetscObjectTypeCompare((PetscObject)da,DMDA,&isdmda);
83: if (!isdmda) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Requires a DMDA as input");
84: /* now we can safely dereference */
85: dd = (DM_DA*)da->data;
86: if (dd->ao) { /* check if the already computed AO has the same type as requested */
87: PetscBool match;
88: PetscObjectTypeCompare((PetscObject)dd->ao,aotype,&match);
89: if (!match) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Cannot change AO type");
90: return(0);
91: }
92: PetscFree(dd->aotype);
93: PetscStrallocpy(aotype,(char**)&dd->aotype);
94: return(0);
95: }
97: /*@
98: DMDAGetAO - Gets the application ordering context for a distributed array.
100: Collective on DMDA
102: Input Parameter:
103: . da - the distributed array
105: Output Parameters:
106: . ao - the application ordering context for DMDAs
108: Level: intermediate
110: Notes:
111: In this case, the AO maps to the natural grid ordering that would be used
112: for the DMDA if only 1 processor were employed (ordering most rapidly in the
113: x-direction, then y, then z). Multiple degrees of freedom are numbered
114: for each node (rather than 1 component for the whole grid, then the next
115: component, etc.)
117: Do NOT call AODestroy() on the ao returned by this function.
119: .keywords: distributed array, get, global, indices, local-to-global
121: .seealso: DMDACreate2d(), DMDASetAOType(), DMDAGetGhostCorners(), DMDAGetCorners(), DMDALocalToGlocal()
122: DMGlobalToLocalBegin(), DMGlobalToLocalEnd(), DMLocalToLocalBegin(), DMLocalToLocalEnd(), DMDAGetOwnershipRanges(),
123: AO, AOPetscToApplication(), AOApplicationToPetsc()
124: @*/
125: PetscErrorCode DMDAGetAO(DM da,AO *ao)
126: {
127: DM_DA *dd;
128: PetscBool isdmda;
134: PetscObjectTypeCompare((PetscObject)da,DMDA,&isdmda);
135: if (!isdmda) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Requires a DMDA as input");
136: /* now we can safely dereference */
137: dd = (DM_DA*)da->data;
139: /*
140: Build the natural ordering to PETSc ordering mappings.
141: */
142: if (!dd->ao) {
143: IS ispetsc,isnatural;
145: PetscInt Nlocal;
147: DMDAGetNatural_Private(da,&Nlocal,&isnatural);
148: ISCreateStride(PetscObjectComm((PetscObject)da),Nlocal,dd->base,1,&ispetsc);
149: AOCreate(PetscObjectComm((PetscObject)da),&dd->ao);
150: AOSetIS(dd->ao,isnatural,ispetsc);
151: AOSetType(dd->ao,dd->aotype);
152: PetscLogObjectParent((PetscObject)da,(PetscObject)dd->ao);
153: ISDestroy(&ispetsc);
154: ISDestroy(&isnatural);
155: }
156: *ao = dd->ao;
157: return(0);
158: }