147c6ae99SBarry Smith 247c6ae99SBarry Smith /* 347c6ae99SBarry Smith Code for manipulating distributed regular arrays in parallel. 447c6ae99SBarry Smith */ 547c6ae99SBarry Smith 64035e84dSBarry Smith #include <petsc-private/dmdaimpl.h> /*I "petscdmda.h" I*/ 747c6ae99SBarry Smith 847c6ae99SBarry Smith #undef __FUNCT__ 9d78e899eSRichard Tran Mills #define __FUNCT__ "DMLocalToLocalCreate_DA" 1047c6ae99SBarry Smith /* 11d78e899eSRichard Tran Mills DMLocalToLocalCreate_DA - Creates the local to local scatter 1247c6ae99SBarry Smith 13aa219208SBarry Smith Collective on DMDA 1447c6ae99SBarry Smith 1547c6ae99SBarry Smith Input Parameter: 1647c6ae99SBarry Smith . da - the distributed array 1747c6ae99SBarry Smith 1847c6ae99SBarry Smith */ 19d78e899eSRichard Tran Mills PetscErrorCode DMLocalToLocalCreate_DA(DM da) 2047c6ae99SBarry Smith { 2147c6ae99SBarry Smith PetscErrorCode ierr; 22*c73cfb54SMatthew G. Knepley PetscInt *idx,left,j,count,up,down,i,bottom,top,k,dim=da->dim; 2347c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 2447c6ae99SBarry Smith 2547c6ae99SBarry Smith PetscFunctionBegin; 2647c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 2747c6ae99SBarry Smith 2847c6ae99SBarry Smith if (dd->ltol) PetscFunctionReturn(0); 2947c6ae99SBarry Smith /* 3047c6ae99SBarry Smith We simply remap the values in the from part of 3147c6ae99SBarry Smith global to local to read from an array with the ghost values 3247c6ae99SBarry Smith rather then from the plain array. 3347c6ae99SBarry Smith */ 3447c6ae99SBarry Smith ierr = VecScatterCopy(dd->gtol,&dd->ltol);CHKERRQ(ierr); 353bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)da,(PetscObject)dd->ltol);CHKERRQ(ierr); 36*c73cfb54SMatthew G. Knepley if (dim == 1) { 3747c6ae99SBarry Smith left = dd->xs - dd->Xs; 38785e854fSJed Brown ierr = PetscMalloc1((dd->xe-dd->xs),&idx);CHKERRQ(ierr); 398865f1eaSKarl Rupp for (j=0; j<dd->xe-dd->xs; j++) idx[j] = left + j; 40*c73cfb54SMatthew G. Knepley } else if (dim == 2) { 4147c6ae99SBarry Smith left = dd->xs - dd->Xs; down = dd->ys - dd->Ys; up = down + dd->ye-dd->ys; 42785e854fSJed Brown ierr = PetscMalloc1((dd->xe-dd->xs)*(up - down),&idx);CHKERRQ(ierr); 4347c6ae99SBarry Smith count = 0; 4447c6ae99SBarry Smith for (i=down; i<up; i++) { 4547c6ae99SBarry Smith for (j=0; j<dd->xe-dd->xs; j++) { 4647c6ae99SBarry Smith idx[count++] = left + i*(dd->Xe-dd->Xs) + j; 4747c6ae99SBarry Smith } 4847c6ae99SBarry Smith } 49*c73cfb54SMatthew G. Knepley } else if (dim == 3) { 5047c6ae99SBarry Smith left = dd->xs - dd->Xs; 5147c6ae99SBarry Smith bottom = dd->ys - dd->Ys; top = bottom + dd->ye-dd->ys; 5247c6ae99SBarry Smith down = dd->zs - dd->Zs; up = down + dd->ze-dd->zs; 5347c6ae99SBarry Smith count = (dd->xe-dd->xs)*(top-bottom)*(up-down); 54785e854fSJed Brown ierr = PetscMalloc1(count,&idx);CHKERRQ(ierr); 5547c6ae99SBarry Smith count = 0; 5647c6ae99SBarry Smith for (i=down; i<up; i++) { 5747c6ae99SBarry Smith for (j=bottom; j<top; j++) { 5847c6ae99SBarry Smith for (k=0; k<dd->xe-dd->xs; k++) { 5947c6ae99SBarry Smith idx[count++] = (left+j*(dd->Xe-dd->Xs))+i*(dd->Xe-dd->Xs)*(dd->Ye-dd->Ys) + k; 6047c6ae99SBarry Smith } 6147c6ae99SBarry Smith } 6247c6ae99SBarry Smith } 63*c73cfb54SMatthew G. Knepley } else SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_CORRUPT,"DMDA has invalid dimension %D",dim); 6447c6ae99SBarry Smith 650298fd71SBarry Smith ierr = VecScatterRemap(dd->ltol,idx,NULL);CHKERRQ(ierr); 6647c6ae99SBarry Smith ierr = PetscFree(idx);CHKERRQ(ierr); 6747c6ae99SBarry Smith PetscFunctionReturn(0); 6847c6ae99SBarry Smith } 6947c6ae99SBarry Smith 7047c6ae99SBarry Smith #undef __FUNCT__ 71d78e899eSRichard Tran Mills #define __FUNCT__ "DMLocalToLocalBegin_DA" 72d78e899eSRichard Tran Mills /* 73d78e899eSRichard Tran Mills DMLocalToLocalBegin_DA - Maps from a local vector (including ghost points 7447c6ae99SBarry Smith that contain irrelevant values) to another local vector where the ghost 75d78e899eSRichard Tran Mills points in the second are set correctly. Must be followed by DMLocalToLocalEnd_DA(). 7647c6ae99SBarry Smith 77aa219208SBarry Smith Neighbor-wise Collective on DMDA and Vec 7847c6ae99SBarry Smith 7947c6ae99SBarry Smith Input Parameters: 8047c6ae99SBarry Smith + da - the distributed array context 8147c6ae99SBarry Smith . g - the original local vector 8247c6ae99SBarry Smith - mode - one of INSERT_VALUES or ADD_VALUES 8347c6ae99SBarry Smith 8447c6ae99SBarry Smith Output Parameter: 8547c6ae99SBarry Smith . l - the local vector with correct ghost values 8647c6ae99SBarry Smith 8747c6ae99SBarry Smith Notes: 8847c6ae99SBarry Smith The local vectors used here need not be the same as those 89564755cdSBarry Smith obtained from DMCreateLocalVector(), BUT they 9047c6ae99SBarry Smith must have the same parallel data layout; they could, for example, be 91aa219208SBarry Smith obtained with VecDuplicate() from the DMDA originating vectors. 9247c6ae99SBarry Smith 93d78e899eSRichard Tran Mills */ 94d78e899eSRichard Tran Mills PetscErrorCode DMLocalToLocalBegin_DA(DM da,Vec g,InsertMode mode,Vec l) 9547c6ae99SBarry Smith { 9647c6ae99SBarry Smith PetscErrorCode ierr; 9747c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 9847c6ae99SBarry Smith 9947c6ae99SBarry Smith PetscFunctionBegin; 10047c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 10147c6ae99SBarry Smith if (!dd->ltol) { 102d78e899eSRichard Tran Mills ierr = DMLocalToLocalCreate_DA(da);CHKERRQ(ierr); 10347c6ae99SBarry Smith } 10447c6ae99SBarry Smith ierr = VecScatterBegin(dd->ltol,g,l,mode,SCATTER_FORWARD);CHKERRQ(ierr); 10547c6ae99SBarry Smith PetscFunctionReturn(0); 10647c6ae99SBarry Smith } 10747c6ae99SBarry Smith 10847c6ae99SBarry Smith #undef __FUNCT__ 109d78e899eSRichard Tran Mills #define __FUNCT__ "DMLocalToLocalEnd_DA" 110d78e899eSRichard Tran Mills /* 111d78e899eSRichard Tran Mills DMLocalToLocalEnd_DA - Maps from a local vector (including ghost points 11247c6ae99SBarry Smith that contain irrelevant values) to another local vector where the ghost 11347c6ae99SBarry Smith points in the second are set correctly. Must be preceeded by 114d78e899eSRichard Tran Mills DMLocalToLocalBegin_DA(). 11547c6ae99SBarry Smith 116aa219208SBarry Smith Neighbor-wise Collective on DMDA and Vec 11747c6ae99SBarry Smith 11847c6ae99SBarry Smith Input Parameters: 11947c6ae99SBarry Smith + da - the distributed array context 12047c6ae99SBarry Smith . g - the original local vector 12147c6ae99SBarry Smith - mode - one of INSERT_VALUES or ADD_VALUES 12247c6ae99SBarry Smith 12347c6ae99SBarry Smith Output Parameter: 12447c6ae99SBarry Smith . l - the local vector with correct ghost values 12547c6ae99SBarry Smith 12647c6ae99SBarry Smith Note: 12747c6ae99SBarry Smith The local vectors used here need not be the same as those 128564755cdSBarry Smith obtained from DMCreateLocalVector(), BUT they 12947c6ae99SBarry Smith must have the same parallel data layout; they could, for example, be 130aa219208SBarry Smith obtained with VecDuplicate() from the DMDA originating vectors. 13147c6ae99SBarry Smith 132d78e899eSRichard Tran Mills */ 133d78e899eSRichard Tran Mills PetscErrorCode DMLocalToLocalEnd_DA(DM da,Vec g,InsertMode mode,Vec l) 13447c6ae99SBarry Smith { 13547c6ae99SBarry Smith PetscErrorCode ierr; 13647c6ae99SBarry Smith DM_DA *dd = (DM_DA*)da->data; 13747c6ae99SBarry Smith 13847c6ae99SBarry Smith PetscFunctionBegin; 13947c6ae99SBarry Smith PetscValidHeaderSpecific(da,DM_CLASSID,1); 14047c6ae99SBarry Smith PetscValidHeaderSpecific(g,VEC_CLASSID,2); 14147c6ae99SBarry Smith PetscValidHeaderSpecific(g,VEC_CLASSID,4); 14247c6ae99SBarry Smith ierr = VecScatterEnd(dd->ltol,g,l,mode,SCATTER_FORWARD);CHKERRQ(ierr); 14347c6ae99SBarry Smith PetscFunctionReturn(0); 14447c6ae99SBarry Smith } 14547c6ae99SBarry Smith 146