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__ "DALocalToLocalCreate" 11*47c6ae99SBarry Smith /* 12*47c6ae99SBarry Smith DALocalToLocalCreate - Creates the local to local scatter 13*47c6ae99SBarry Smith 14*47c6ae99SBarry Smith Collective on DA 15*47c6ae99SBarry Smith 16*47c6ae99SBarry Smith Input Parameter: 17*47c6ae99SBarry Smith . da - the distributed array 18*47c6ae99SBarry Smith 19*47c6ae99SBarry Smith */ 20*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DALocalToLocalCreate(DA da) 21*47c6ae99SBarry Smith { 22*47c6ae99SBarry Smith PetscErrorCode ierr; 23*47c6ae99SBarry Smith PetscInt *idx,left,j,count,up,down,i,bottom,top,k; 24*47c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 25*47c6ae99SBarry Smith 26*47c6ae99SBarry Smith PetscFunctionBegin; 27*47c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 28*47c6ae99SBarry Smith 29*47c6ae99SBarry Smith if (dd->ltol) PetscFunctionReturn(0); 30*47c6ae99SBarry Smith /* 31*47c6ae99SBarry Smith We simply remap the values in the from part of 32*47c6ae99SBarry Smith global to local to read from an array with the ghost values 33*47c6ae99SBarry Smith rather then from the plain array. 34*47c6ae99SBarry Smith */ 35*47c6ae99SBarry Smith ierr = VecScatterCopy(dd->gtol,&dd->ltol);CHKERRQ(ierr); 36*47c6ae99SBarry Smith ierr = PetscLogObjectParent(da,dd->ltol);CHKERRQ(ierr); 37*47c6ae99SBarry Smith if (dd->dim == 1) { 38*47c6ae99SBarry Smith left = dd->xs - dd->Xs; 39*47c6ae99SBarry Smith ierr = PetscMalloc((dd->xe-dd->xs)*sizeof(PetscInt),&idx);CHKERRQ(ierr); 40*47c6ae99SBarry Smith for (j=0; j<dd->xe-dd->xs; j++) { 41*47c6ae99SBarry Smith idx[j] = left + j; 42*47c6ae99SBarry Smith } 43*47c6ae99SBarry Smith } else if (dd->dim == 2) { 44*47c6ae99SBarry Smith left = dd->xs - dd->Xs; down = dd->ys - dd->Ys; up = down + dd->ye-dd->ys; 45*47c6ae99SBarry Smith ierr = PetscMalloc((dd->xe-dd->xs)*(up - down)*sizeof(PetscInt),&idx);CHKERRQ(ierr); 46*47c6ae99SBarry Smith count = 0; 47*47c6ae99SBarry Smith for (i=down; i<up; i++) { 48*47c6ae99SBarry Smith for (j=0; j<dd->xe-dd->xs; j++) { 49*47c6ae99SBarry Smith idx[count++] = left + i*(dd->Xe-dd->Xs) + j; 50*47c6ae99SBarry Smith } 51*47c6ae99SBarry Smith } 52*47c6ae99SBarry Smith } else if (dd->dim == 3) { 53*47c6ae99SBarry Smith left = dd->xs - dd->Xs; 54*47c6ae99SBarry Smith bottom = dd->ys - dd->Ys; top = bottom + dd->ye-dd->ys ; 55*47c6ae99SBarry Smith down = dd->zs - dd->Zs; up = down + dd->ze-dd->zs; 56*47c6ae99SBarry Smith count = (dd->xe-dd->xs)*(top-bottom)*(up-down); 57*47c6ae99SBarry Smith ierr = PetscMalloc(count*sizeof(PetscInt),&idx);CHKERRQ(ierr); 58*47c6ae99SBarry Smith count = 0; 59*47c6ae99SBarry Smith for (i=down; i<up; i++) { 60*47c6ae99SBarry Smith for (j=bottom; j<top; j++) { 61*47c6ae99SBarry Smith for (k=0; k<dd->xe-dd->xs; k++) { 62*47c6ae99SBarry Smith idx[count++] = (left+j*(dd->Xe-dd->Xs))+i*(dd->Xe-dd->Xs)*(dd->Ye-dd->Ys) + k; 63*47c6ae99SBarry Smith } 64*47c6ae99SBarry Smith } 65*47c6ae99SBarry Smith } 66*47c6ae99SBarry Smith } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_CORRUPT,"DA has invalid dimension %D",dd->dim); 67*47c6ae99SBarry Smith 68*47c6ae99SBarry Smith ierr = VecScatterRemap(dd->ltol,idx,PETSC_NULL);CHKERRQ(ierr); 69*47c6ae99SBarry Smith ierr = PetscFree(idx);CHKERRQ(ierr); 70*47c6ae99SBarry Smith PetscFunctionReturn(0); 71*47c6ae99SBarry Smith } 72*47c6ae99SBarry Smith 73*47c6ae99SBarry Smith #undef __FUNCT__ 74*47c6ae99SBarry Smith #define __FUNCT__ "DALocalToLocalBegin" 75*47c6ae99SBarry Smith /*@ 76*47c6ae99SBarry Smith DALocalToLocalBegin - Maps from a local vector (including ghost points 77*47c6ae99SBarry Smith that contain irrelevant values) to another local vector where the ghost 78*47c6ae99SBarry Smith points in the second are set correctly. Must be followed by DALocalToLocalEnd(). 79*47c6ae99SBarry Smith 80*47c6ae99SBarry Smith Neighbor-wise Collective on DA and Vec 81*47c6ae99SBarry Smith 82*47c6ae99SBarry Smith Input Parameters: 83*47c6ae99SBarry Smith + da - the distributed array context 84*47c6ae99SBarry Smith . g - the original local vector 85*47c6ae99SBarry Smith - mode - one of INSERT_VALUES or ADD_VALUES 86*47c6ae99SBarry Smith 87*47c6ae99SBarry Smith Output Parameter: 88*47c6ae99SBarry Smith . l - the local vector with correct ghost values 89*47c6ae99SBarry Smith 90*47c6ae99SBarry Smith Level: intermediate 91*47c6ae99SBarry Smith 92*47c6ae99SBarry Smith Notes: 93*47c6ae99SBarry Smith The local vectors used here need not be the same as those 94*47c6ae99SBarry Smith obtained from DACreateLocalVector(), BUT they 95*47c6ae99SBarry Smith must have the same parallel data layout; they could, for example, be 96*47c6ae99SBarry Smith obtained with VecDuplicate() from the DA originating vectors. 97*47c6ae99SBarry Smith 98*47c6ae99SBarry Smith .keywords: distributed array, local-to-local, begin 99*47c6ae99SBarry Smith 100*47c6ae99SBarry Smith .seealso: DALocalToLocalEnd(), DALocalToGlobal() 101*47c6ae99SBarry Smith @*/ 102*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DALocalToLocalBegin(DA da,Vec g,InsertMode mode,Vec l) 103*47c6ae99SBarry Smith { 104*47c6ae99SBarry Smith PetscErrorCode ierr; 105*47c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 106*47c6ae99SBarry Smith 107*47c6ae99SBarry Smith PetscFunctionBegin; 108*47c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 109*47c6ae99SBarry Smith if (!dd->ltol) { 110*47c6ae99SBarry Smith ierr = DALocalToLocalCreate(da);CHKERRQ(ierr); 111*47c6ae99SBarry Smith } 112*47c6ae99SBarry Smith ierr = VecScatterBegin(dd->ltol,g,l,mode,SCATTER_FORWARD);CHKERRQ(ierr); 113*47c6ae99SBarry Smith PetscFunctionReturn(0); 114*47c6ae99SBarry Smith } 115*47c6ae99SBarry Smith 116*47c6ae99SBarry Smith #undef __FUNCT__ 117*47c6ae99SBarry Smith #define __FUNCT__ "DALocalToLocalEnd" 118*47c6ae99SBarry Smith /*@ 119*47c6ae99SBarry Smith DALocalToLocalEnd - Maps from a local vector (including ghost points 120*47c6ae99SBarry Smith that contain irrelevant values) to another local vector where the ghost 121*47c6ae99SBarry Smith points in the second are set correctly. Must be preceeded by 122*47c6ae99SBarry Smith DALocalToLocalBegin(). 123*47c6ae99SBarry Smith 124*47c6ae99SBarry Smith Neighbor-wise Collective on DA and Vec 125*47c6ae99SBarry Smith 126*47c6ae99SBarry Smith Input Parameters: 127*47c6ae99SBarry Smith + da - the distributed array context 128*47c6ae99SBarry Smith . g - the original local vector 129*47c6ae99SBarry Smith - mode - one of INSERT_VALUES or ADD_VALUES 130*47c6ae99SBarry Smith 131*47c6ae99SBarry Smith Output Parameter: 132*47c6ae99SBarry Smith . l - the local vector with correct ghost values 133*47c6ae99SBarry Smith 134*47c6ae99SBarry Smith Level: intermediate 135*47c6ae99SBarry Smith 136*47c6ae99SBarry Smith Note: 137*47c6ae99SBarry Smith The local vectors used here need not be the same as those 138*47c6ae99SBarry Smith obtained from DACreateLocalVector(), BUT they 139*47c6ae99SBarry Smith must have the same parallel data layout; they could, for example, be 140*47c6ae99SBarry Smith obtained with VecDuplicate() from the DA originating vectors. 141*47c6ae99SBarry Smith 142*47c6ae99SBarry Smith .keywords: distributed array, local-to-local, end 143*47c6ae99SBarry Smith 144*47c6ae99SBarry Smith .seealso: DALocalToLocalBegin(), DALocalToGlobal() 145*47c6ae99SBarry Smith @*/ 146*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DALocalToLocalEnd(DA da,Vec g,InsertMode mode,Vec l) 147*47c6ae99SBarry Smith { 148*47c6ae99SBarry Smith PetscErrorCode ierr; 149*47c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 150*47c6ae99SBarry Smith 151*47c6ae99SBarry Smith PetscFunctionBegin; 152*47c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 153*47c6ae99SBarry Smith PetscValidHeaderSpecific(g,VEC_CLASSID,2); 154*47c6ae99SBarry Smith PetscValidHeaderSpecific(g,VEC_CLASSID,4); 155*47c6ae99SBarry Smith ierr = VecScatterEnd(dd->ltol,g,l,mode,SCATTER_FORWARD);CHKERRQ(ierr); 156*47c6ae99SBarry Smith PetscFunctionReturn(0); 157*47c6ae99SBarry Smith } 158*47c6ae99SBarry Smith 159