147c6ae99SBarry Smith 247c6ae99SBarry Smith /* 347c6ae99SBarry Smith Code for manipulating distributed regular arrays in parallel. 447c6ae99SBarry Smith */ 547c6ae99SBarry Smith 6af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h> /*I "petscdmda.h" I*/ 747c6ae99SBarry Smith 847c6ae99SBarry Smith /* 9d78e899eSRichard Tran Mills DMLocalToLocalCreate_DA - Creates the local to local scatter 1047c6ae99SBarry Smith 11d083f849SBarry Smith Collective on da 1247c6ae99SBarry Smith 1347c6ae99SBarry Smith Input Parameter: 1447c6ae99SBarry Smith . da - the distributed array 1547c6ae99SBarry Smith 1647c6ae99SBarry Smith */ 179371c9d4SSatish Balay PetscErrorCode DMLocalToLocalCreate_DA(DM da) { 18c73cfb54SMatthew G. Knepley PetscInt *idx, left, j, count, up, down, i, bottom, top, k, dim = da->dim; 1947c6ae99SBarry Smith DM_DA *dd = (DM_DA *)da->data; 2047c6ae99SBarry Smith 2147c6ae99SBarry Smith PetscFunctionBegin; 2247c6ae99SBarry Smith PetscValidHeaderSpecific(da, DM_CLASSID, 1); 2347c6ae99SBarry Smith 2447c6ae99SBarry Smith if (dd->ltol) PetscFunctionReturn(0); 2547c6ae99SBarry Smith /* 2647c6ae99SBarry Smith We simply remap the values in the from part of 2747c6ae99SBarry Smith global to local to read from an array with the ghost values 2847c6ae99SBarry Smith rather then from the plain array. 2947c6ae99SBarry Smith */ 309566063dSJacob Faibussowitsch PetscCall(VecScatterCopy(dd->gtol, &dd->ltol)); 319566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)da, (PetscObject)dd->ltol)); 32c73cfb54SMatthew G. Knepley if (dim == 1) { 3347c6ae99SBarry Smith left = dd->xs - dd->Xs; 349566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dd->xe - dd->xs, &idx)); 358865f1eaSKarl Rupp for (j = 0; j < dd->xe - dd->xs; j++) idx[j] = left + j; 36c73cfb54SMatthew G. Knepley } else if (dim == 2) { 379371c9d4SSatish Balay left = dd->xs - dd->Xs; 389371c9d4SSatish Balay down = dd->ys - dd->Ys; 399371c9d4SSatish Balay up = down + dd->ye - dd->ys; 409566063dSJacob Faibussowitsch PetscCall(PetscMalloc1((dd->xe - dd->xs) * (up - down), &idx)); 4147c6ae99SBarry Smith count = 0; 4247c6ae99SBarry Smith for (i = down; i < up; i++) { 439371c9d4SSatish Balay for (j = 0; j < dd->xe - dd->xs; j++) { idx[count++] = left + i * (dd->Xe - dd->Xs) + j; } 4447c6ae99SBarry Smith } 45c73cfb54SMatthew G. Knepley } else if (dim == 3) { 4647c6ae99SBarry Smith left = dd->xs - dd->Xs; 479371c9d4SSatish Balay bottom = dd->ys - dd->Ys; 489371c9d4SSatish Balay top = bottom + dd->ye - dd->ys; 499371c9d4SSatish Balay down = dd->zs - dd->Zs; 509371c9d4SSatish Balay up = down + dd->ze - dd->zs; 5147c6ae99SBarry Smith count = (dd->xe - dd->xs) * (top - bottom) * (up - down); 529566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(count, &idx)); 5347c6ae99SBarry Smith count = 0; 5447c6ae99SBarry Smith for (i = down; i < up; i++) { 5547c6ae99SBarry Smith for (j = bottom; j < top; j++) { 569371c9d4SSatish Balay for (k = 0; k < dd->xe - dd->xs; k++) { idx[count++] = (left + j * (dd->Xe - dd->Xs)) + i * (dd->Xe - dd->Xs) * (dd->Ye - dd->Ys) + k; } 5747c6ae99SBarry Smith } 5847c6ae99SBarry Smith } 5963a3b9bcSJacob Faibussowitsch } else SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_CORRUPT, "DMDA has invalid dimension %" PetscInt_FMT, dim); 6047c6ae99SBarry Smith 619566063dSJacob Faibussowitsch PetscCall(VecScatterRemap(dd->ltol, idx, NULL)); 629566063dSJacob Faibussowitsch PetscCall(PetscFree(idx)); 6347c6ae99SBarry Smith PetscFunctionReturn(0); 6447c6ae99SBarry Smith } 6547c6ae99SBarry Smith 66d78e899eSRichard Tran Mills /* 67d78e899eSRichard Tran Mills DMLocalToLocalBegin_DA - Maps from a local vector (including ghost points 6847c6ae99SBarry Smith that contain irrelevant values) to another local vector where the ghost 69d78e899eSRichard Tran Mills points in the second are set correctly. Must be followed by DMLocalToLocalEnd_DA(). 7047c6ae99SBarry Smith 71d083f849SBarry Smith Neighbor-wise Collective on da 7247c6ae99SBarry Smith 7347c6ae99SBarry Smith Input Parameters: 7447c6ae99SBarry Smith + da - the distributed array context 7547c6ae99SBarry Smith . g - the original local vector 7647c6ae99SBarry Smith - mode - one of INSERT_VALUES or ADD_VALUES 7747c6ae99SBarry Smith 7847c6ae99SBarry Smith Output Parameter: 7947c6ae99SBarry Smith . l - the local vector with correct ghost values 8047c6ae99SBarry Smith 8147c6ae99SBarry Smith Notes: 8247c6ae99SBarry Smith The local vectors used here need not be the same as those 83564755cdSBarry Smith obtained from DMCreateLocalVector(), BUT they 8447c6ae99SBarry Smith must have the same parallel data layout; they could, for example, be 85aa219208SBarry Smith obtained with VecDuplicate() from the DMDA originating vectors. 8647c6ae99SBarry Smith 87d78e899eSRichard Tran Mills */ 889371c9d4SSatish Balay PetscErrorCode DMLocalToLocalBegin_DA(DM da, Vec g, InsertMode mode, Vec l) { 8947c6ae99SBarry Smith DM_DA *dd = (DM_DA *)da->data; 9047c6ae99SBarry Smith 9147c6ae99SBarry Smith PetscFunctionBegin; 9247c6ae99SBarry Smith PetscValidHeaderSpecific(da, DM_CLASSID, 1); 93*48a46eb9SPierre Jolivet if (!dd->ltol) PetscCall(DMLocalToLocalCreate_DA(da)); 949566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(dd->ltol, g, l, mode, SCATTER_FORWARD)); 9547c6ae99SBarry Smith PetscFunctionReturn(0); 9647c6ae99SBarry Smith } 9747c6ae99SBarry Smith 98d78e899eSRichard Tran Mills /* 99d78e899eSRichard Tran Mills DMLocalToLocalEnd_DA - Maps from a local vector (including ghost points 10047c6ae99SBarry Smith that contain irrelevant values) to another local vector where the ghost 101a5b23f4aSJose E. Roman points in the second are set correctly. Must be preceded by 102d78e899eSRichard Tran Mills DMLocalToLocalBegin_DA(). 10347c6ae99SBarry Smith 104d083f849SBarry Smith Neighbor-wise Collective on da 10547c6ae99SBarry Smith 10647c6ae99SBarry Smith Input Parameters: 10747c6ae99SBarry Smith + da - the distributed array context 10847c6ae99SBarry Smith . g - the original local vector 10947c6ae99SBarry Smith - mode - one of INSERT_VALUES or ADD_VALUES 11047c6ae99SBarry Smith 11147c6ae99SBarry Smith Output Parameter: 11247c6ae99SBarry Smith . l - the local vector with correct ghost values 11347c6ae99SBarry Smith 11447c6ae99SBarry Smith Note: 11547c6ae99SBarry Smith The local vectors used here need not be the same as those 116564755cdSBarry Smith obtained from DMCreateLocalVector(), BUT they 11747c6ae99SBarry Smith must have the same parallel data layout; they could, for example, be 118aa219208SBarry Smith obtained with VecDuplicate() from the DMDA originating vectors. 11947c6ae99SBarry Smith 120d78e899eSRichard Tran Mills */ 1219371c9d4SSatish Balay PetscErrorCode DMLocalToLocalEnd_DA(DM da, Vec g, InsertMode mode, Vec l) { 12247c6ae99SBarry Smith DM_DA *dd = (DM_DA *)da->data; 12347c6ae99SBarry Smith 12447c6ae99SBarry Smith PetscFunctionBegin; 12547c6ae99SBarry Smith PetscValidHeaderSpecific(da, DM_CLASSID, 1); 12647c6ae99SBarry Smith PetscValidHeaderSpecific(g, VEC_CLASSID, 2); 1279566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(dd->ltol, g, l, mode, SCATTER_FORWARD)); 12847c6ae99SBarry Smith PetscFunctionReturn(0); 12947c6ae99SBarry Smith } 130