xref: /petsc/src/dm/impls/da/daltol.c (revision 20f4b53cbb5e9bd9ef12b76a8697d60d197cda17)
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 
11*20f4b53cSBarry Smith    Collective
1247c6ae99SBarry Smith 
1347c6ae99SBarry Smith    Input Parameter:
1447c6ae99SBarry Smith .  da - the distributed array
1547c6ae99SBarry Smith 
1647c6ae99SBarry Smith */
17d71ae5a4SJacob Faibussowitsch PetscErrorCode DMLocalToLocalCreate_DA(DM da)
18d71ae5a4SJacob Faibussowitsch {
19c73cfb54SMatthew G. Knepley   PetscInt *idx, left, j, count, up, down, i, bottom, top, k, dim = da->dim;
2047c6ae99SBarry Smith   DM_DA    *dd = (DM_DA *)da->data;
2147c6ae99SBarry Smith 
2247c6ae99SBarry Smith   PetscFunctionBegin;
2347c6ae99SBarry Smith   PetscValidHeaderSpecific(da, DM_CLASSID, 1);
2447c6ae99SBarry Smith 
253ba16761SJacob Faibussowitsch   if (dd->ltol) PetscFunctionReturn(PETSC_SUCCESS);
2647c6ae99SBarry Smith   /*
2747c6ae99SBarry Smith      We simply remap the values in the from part of
2847c6ae99SBarry Smith      global to local to read from an array with the ghost values
2947c6ae99SBarry Smith      rather then from the plain array.
3047c6ae99SBarry Smith   */
319566063dSJacob Faibussowitsch   PetscCall(VecScatterCopy(dd->gtol, &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++) {
43ad540459SPierre Jolivet       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++) {
56ad540459SPierre Jolivet         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));
633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6447c6ae99SBarry Smith }
6547c6ae99SBarry Smith 
66d71ae5a4SJacob Faibussowitsch PetscErrorCode DMLocalToLocalBegin_DA(DM da, Vec g, InsertMode mode, Vec l)
67d71ae5a4SJacob Faibussowitsch {
6847c6ae99SBarry Smith   DM_DA *dd = (DM_DA *)da->data;
6947c6ae99SBarry Smith 
7047c6ae99SBarry Smith   PetscFunctionBegin;
7147c6ae99SBarry Smith   PetscValidHeaderSpecific(da, DM_CLASSID, 1);
7248a46eb9SPierre Jolivet   if (!dd->ltol) PetscCall(DMLocalToLocalCreate_DA(da));
739566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(dd->ltol, g, l, mode, SCATTER_FORWARD));
743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7547c6ae99SBarry Smith }
7647c6ae99SBarry Smith 
77d71ae5a4SJacob Faibussowitsch PetscErrorCode DMLocalToLocalEnd_DA(DM da, Vec g, InsertMode mode, Vec l)
78d71ae5a4SJacob Faibussowitsch {
7947c6ae99SBarry Smith   DM_DA *dd = (DM_DA *)da->data;
8047c6ae99SBarry Smith 
8147c6ae99SBarry Smith   PetscFunctionBegin;
8247c6ae99SBarry Smith   PetscValidHeaderSpecific(da, DM_CLASSID, 1);
8347c6ae99SBarry Smith   PetscValidHeaderSpecific(g, VEC_CLASSID, 2);
849566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(dd->ltol, g, l, mode, SCATTER_FORWARD));
853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8647c6ae99SBarry Smith }
87