1*47c6ae99SBarry Smith #define PETSCDM_DLL 2*47c6ae99SBarry Smith 3*47c6ae99SBarry Smith /* 4*47c6ae99SBarry Smith Code for manipulating distributed regular arrays in parallel. 5*47c6ae99SBarry Smith */ 6*47c6ae99SBarry Smith 7*47c6ae99SBarry Smith #include "private/daimpl.h" /*I "petscda.h" I*/ 8*47c6ae99SBarry Smith 9*47c6ae99SBarry Smith #undef __FUNCT__ 10*47c6ae99SBarry Smith #define __FUNCT__ "DAGetGlobalIndices" 11*47c6ae99SBarry Smith /*@C 12*47c6ae99SBarry Smith DAGetGlobalIndices - Returns the global node number of all local nodes, 13*47c6ae99SBarry Smith including ghost nodes. 14*47c6ae99SBarry Smith 15*47c6ae99SBarry Smith Not Collective 16*47c6ae99SBarry Smith 17*47c6ae99SBarry Smith Input Parameter: 18*47c6ae99SBarry Smith . da - the distributed array 19*47c6ae99SBarry Smith 20*47c6ae99SBarry Smith Output Parameters: 21*47c6ae99SBarry Smith + n - the number of local elements, including ghost nodes (or PETSC_NULL) 22*47c6ae99SBarry Smith - idx - the global indices 23*47c6ae99SBarry Smith 24*47c6ae99SBarry Smith Level: intermediate 25*47c6ae99SBarry Smith 26*47c6ae99SBarry Smith Note: 27*47c6ae99SBarry Smith For DA_STENCIL_STAR stencils the inactive corner ghost nodes are also included 28*47c6ae99SBarry Smith in the list of local indices (even though those nodes are not updated 29*47c6ae99SBarry Smith during calls to DAXXXToXXX(). 30*47c6ae99SBarry Smith 31*47c6ae99SBarry Smith Essentially the same data is returned in the form of a local-to-global mapping 32*47c6ae99SBarry Smith with the routine DAGetISLocalToGlobalMapping(); 33*47c6ae99SBarry Smith 34*47c6ae99SBarry Smith Fortran Note: 35*47c6ae99SBarry Smith This routine is used differently from Fortran 36*47c6ae99SBarry Smith .vb 37*47c6ae99SBarry Smith DA da 38*47c6ae99SBarry Smith integer n,da_array(1) 39*47c6ae99SBarry Smith PetscOffset i_da 40*47c6ae99SBarry Smith integer ierr 41*47c6ae99SBarry Smith call DAGetGlobalIndices(da,n,da_array,i_da,ierr) 42*47c6ae99SBarry Smith 43*47c6ae99SBarry Smith C Access first local entry in list 44*47c6ae99SBarry Smith value = da_array(i_da + 1) 45*47c6ae99SBarry Smith .ve 46*47c6ae99SBarry Smith 47*47c6ae99SBarry Smith See the <A href="../../docs/manual.pdf#nameddest=ch_fortran">Fortran chapter</A> of the users manual for details. 48*47c6ae99SBarry Smith 49*47c6ae99SBarry Smith .keywords: distributed array, get, global, indices, local-to-global 50*47c6ae99SBarry Smith 51*47c6ae99SBarry Smith .seealso: DACreate2d(), DAGetGhostCorners(), DAGetCorners(), DALocalToGlobal() 52*47c6ae99SBarry Smith DAGlobalToLocalBegin(), DAGlobalToLocalEnd(), DALocalToLocalBegin(), DAGetAO(), DAGetGlobalIndicesF90() 53*47c6ae99SBarry Smith DAGetISLocalToGlobalMapping(), DACreate3d(), DACreate1d(), DALocalToLocalEnd(), DAGetOwnershipRanges() 54*47c6ae99SBarry Smith @*/ 55*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DAGetGlobalIndices(DA da,PetscInt *n,PetscInt **idx) 56*47c6ae99SBarry Smith { 57*47c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 58*47c6ae99SBarry Smith 59*47c6ae99SBarry Smith PetscFunctionBegin; 60*47c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 61*47c6ae99SBarry Smith if (n) *n = dd->Nl; 62*47c6ae99SBarry Smith if (idx) *idx = dd->idx; 63*47c6ae99SBarry Smith PetscFunctionReturn(0); 64*47c6ae99SBarry Smith } 65*47c6ae99SBarry Smith 66*47c6ae99SBarry Smith #undef __FUNCT__ 67*47c6ae99SBarry Smith #define __FUNCT__ "DAGetNatural_Private" 68*47c6ae99SBarry Smith /* 69*47c6ae99SBarry Smith Gets the natural number for each global number on the process. 70*47c6ae99SBarry Smith 71*47c6ae99SBarry Smith Used by DAGetAO() and DAGlobalToNatural_Create() 72*47c6ae99SBarry Smith */ 73*47c6ae99SBarry Smith PetscErrorCode DAGetNatural_Private(DA da,PetscInt *outNlocal,IS *isnatural) 74*47c6ae99SBarry Smith { 75*47c6ae99SBarry Smith PetscErrorCode ierr; 76*47c6ae99SBarry Smith PetscInt Nlocal,i,j,k,*lidx,lict = 0; 77*47c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 78*47c6ae99SBarry Smith 79*47c6ae99SBarry Smith PetscFunctionBegin; 80*47c6ae99SBarry Smith Nlocal = (dd->xe-dd->xs); 81*47c6ae99SBarry Smith if (dd->dim > 1) { 82*47c6ae99SBarry Smith Nlocal *= (dd->ye-dd->ys); 83*47c6ae99SBarry Smith } 84*47c6ae99SBarry Smith if (dd->dim > 2) { 85*47c6ae99SBarry Smith Nlocal *= (dd->ze-dd->zs); 86*47c6ae99SBarry Smith } 87*47c6ae99SBarry Smith 88*47c6ae99SBarry Smith ierr = PetscMalloc(Nlocal*sizeof(PetscInt),&lidx);CHKERRQ(ierr); 89*47c6ae99SBarry Smith 90*47c6ae99SBarry Smith if (dd->dim == 1) { 91*47c6ae99SBarry Smith for (i=dd->xs; i<dd->xe; i++) { 92*47c6ae99SBarry Smith /* global number in natural ordering */ 93*47c6ae99SBarry Smith lidx[lict++] = i; 94*47c6ae99SBarry Smith } 95*47c6ae99SBarry Smith } else if (dd->dim == 2) { 96*47c6ae99SBarry Smith for (j=dd->ys; j<dd->ye; j++) { 97*47c6ae99SBarry Smith for (i=dd->xs; i<dd->xe; i++) { 98*47c6ae99SBarry Smith /* global number in natural ordering */ 99*47c6ae99SBarry Smith lidx[lict++] = i + j*dd->M*dd->w; 100*47c6ae99SBarry Smith } 101*47c6ae99SBarry Smith } 102*47c6ae99SBarry Smith } else if (dd->dim == 3) { 103*47c6ae99SBarry Smith for (k=dd->zs; k<dd->ze; k++) { 104*47c6ae99SBarry Smith for (j=dd->ys; j<dd->ye; j++) { 105*47c6ae99SBarry Smith for (i=dd->xs; i<dd->xe; i++) { 106*47c6ae99SBarry Smith lidx[lict++] = i + j*dd->M*dd->w + k*dd->M*dd->N*dd->w; 107*47c6ae99SBarry Smith } 108*47c6ae99SBarry Smith } 109*47c6ae99SBarry Smith } 110*47c6ae99SBarry Smith } 111*47c6ae99SBarry Smith *outNlocal = Nlocal; 112*47c6ae99SBarry Smith ierr = ISCreateGeneral(((PetscObject)da)->comm,Nlocal,lidx,PETSC_OWN_POINTER,isnatural);CHKERRQ(ierr); 113*47c6ae99SBarry Smith PetscFunctionReturn(0); 114*47c6ae99SBarry Smith } 115*47c6ae99SBarry Smith 116*47c6ae99SBarry Smith #undef __FUNCT__ 117*47c6ae99SBarry Smith #define __FUNCT__ "DAGetAO" 118*47c6ae99SBarry Smith /*@ 119*47c6ae99SBarry Smith DAGetAO - Gets the application ordering context for a distributed array. 120*47c6ae99SBarry Smith 121*47c6ae99SBarry Smith Collective on DA 122*47c6ae99SBarry Smith 123*47c6ae99SBarry Smith Input Parameter: 124*47c6ae99SBarry Smith . da - the distributed array 125*47c6ae99SBarry Smith 126*47c6ae99SBarry Smith Output Parameters: 127*47c6ae99SBarry Smith . ao - the application ordering context for DAs 128*47c6ae99SBarry Smith 129*47c6ae99SBarry Smith Level: intermediate 130*47c6ae99SBarry Smith 131*47c6ae99SBarry Smith Notes: 132*47c6ae99SBarry Smith In this case, the AO maps to the natural grid ordering that would be used 133*47c6ae99SBarry Smith for the DA if only 1 processor were employed (ordering most rapidly in the 134*47c6ae99SBarry Smith x-direction, then y, then z). Multiple degrees of freedom are numbered 135*47c6ae99SBarry Smith for each node (rather than 1 component for the whole grid, then the next 136*47c6ae99SBarry Smith component, etc.) 137*47c6ae99SBarry Smith 138*47c6ae99SBarry Smith .keywords: distributed array, get, global, indices, local-to-global 139*47c6ae99SBarry Smith 140*47c6ae99SBarry Smith .seealso: DACreate2d(), DAGetGhostCorners(), DAGetCorners(), DALocalToGlocal() 141*47c6ae99SBarry Smith DAGlobalToLocalBegin(), DAGlobalToLocalEnd(), DALocalToLocalBegin(), DALocalToLocalEnd(), DAGetGlobalIndices(), DAGetOwnershipRanges(), 142*47c6ae99SBarry Smith AO, AOPetscToApplication(), AOApplicationToPetsc() 143*47c6ae99SBarry Smith @*/ 144*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DAGetAO(DA da,AO *ao) 145*47c6ae99SBarry Smith { 146*47c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 147*47c6ae99SBarry Smith 148*47c6ae99SBarry Smith PetscFunctionBegin; 149*47c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 150*47c6ae99SBarry Smith PetscValidPointer(ao,2); 151*47c6ae99SBarry Smith 152*47c6ae99SBarry Smith /* 153*47c6ae99SBarry Smith Build the natural ordering to PETSc ordering mappings. 154*47c6ae99SBarry Smith */ 155*47c6ae99SBarry Smith if (!dd->ao) { 156*47c6ae99SBarry Smith IS ispetsc,isnatural; 157*47c6ae99SBarry Smith PetscErrorCode ierr; 158*47c6ae99SBarry Smith PetscInt Nlocal; 159*47c6ae99SBarry Smith 160*47c6ae99SBarry Smith ierr = DAGetNatural_Private(da,&Nlocal,&isnatural);CHKERRQ(ierr); 161*47c6ae99SBarry Smith ierr = ISCreateStride(((PetscObject)da)->comm,Nlocal,dd->base,1,&ispetsc);CHKERRQ(ierr); 162*47c6ae99SBarry Smith ierr = AOCreateBasicIS(isnatural,ispetsc,&dd->ao);CHKERRQ(ierr); 163*47c6ae99SBarry Smith ierr = PetscLogObjectParent(da,dd->ao);CHKERRQ(ierr); 164*47c6ae99SBarry Smith ierr = ISDestroy(ispetsc);CHKERRQ(ierr); 165*47c6ae99SBarry Smith ierr = ISDestroy(isnatural);CHKERRQ(ierr); 166*47c6ae99SBarry Smith } 167*47c6ae99SBarry Smith *ao = dd->ao; 168*47c6ae99SBarry Smith PetscFunctionReturn(0); 169*47c6ae99SBarry Smith } 170*47c6ae99SBarry Smith 171*47c6ae99SBarry Smith /*MC 172*47c6ae99SBarry Smith DAGetGlobalIndicesF90 - Returns a Fortran90 pointer to the list of 173*47c6ae99SBarry Smith global indices (global node number of all local nodes, including 174*47c6ae99SBarry Smith ghost nodes). 175*47c6ae99SBarry Smith 176*47c6ae99SBarry Smith Synopsis: 177*47c6ae99SBarry Smith DAGetGlobalIndicesF90(DA da,integer n,{integer, pointer :: idx(:)},integer ierr) 178*47c6ae99SBarry Smith 179*47c6ae99SBarry Smith Not Collective 180*47c6ae99SBarry Smith 181*47c6ae99SBarry Smith Input Parameter: 182*47c6ae99SBarry Smith . da - the distributed array 183*47c6ae99SBarry Smith 184*47c6ae99SBarry Smith Output Parameters: 185*47c6ae99SBarry Smith + n - the number of local elements, including ghost nodes (or PETSC_NULL) 186*47c6ae99SBarry Smith . idx - the Fortran90 pointer to the global indices 187*47c6ae99SBarry Smith - ierr - error code 188*47c6ae99SBarry Smith 189*47c6ae99SBarry Smith Level: intermediate 190*47c6ae99SBarry Smith 191*47c6ae99SBarry Smith Notes: 192*47c6ae99SBarry Smith Not yet supported for all F90 compilers 193*47c6ae99SBarry Smith 194*47c6ae99SBarry Smith .keywords: distributed array, get, global, indices, local-to-global, f90 195*47c6ae99SBarry Smith 196*47c6ae99SBarry Smith .seealso: DAGetGlobalIndices() 197*47c6ae99SBarry Smith M*/ 198