xref: /petsc/src/dm/impls/da/daltol.c (revision 467446fbca37e95635331dc8e534c1ae283105d1)
147c6ae99SBarry Smith /*
247c6ae99SBarry Smith   Code for manipulating distributed regular arrays in parallel.
347c6ae99SBarry Smith */
447c6ae99SBarry Smith 
5af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h> /*I   "petscdmda.h"   I*/
647c6ae99SBarry Smith 
747c6ae99SBarry Smith /*
8d78e899eSRichard Tran Mills    DMLocalToLocalCreate_DA - Creates the local to local scatter
947c6ae99SBarry Smith 
1020f4b53cSBarry Smith    Collective
1147c6ae99SBarry Smith 
1247c6ae99SBarry Smith    Input Parameter:
1347c6ae99SBarry Smith .  da - the distributed array
1447c6ae99SBarry Smith 
1547c6ae99SBarry Smith */
16d71ae5a4SJacob Faibussowitsch PetscErrorCode DMLocalToLocalCreate_DA(DM da)
17d71ae5a4SJacob Faibussowitsch {
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 
243ba16761SJacob Faibussowitsch   if (dd->ltol) PetscFunctionReturn(PETSC_SUCCESS);
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
28*467446fbSPierre Jolivet      rather than from the plain array.
2947c6ae99SBarry Smith   */
309566063dSJacob Faibussowitsch   PetscCall(VecScatterCopy(dd->gtol, &dd->ltol));
31c73cfb54SMatthew G. Knepley   if (dim == 1) {
3247c6ae99SBarry Smith     left = dd->xs - dd->Xs;
339566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(dd->xe - dd->xs, &idx));
348865f1eaSKarl Rupp     for (j = 0; j < dd->xe - dd->xs; j++) idx[j] = left + j;
35c73cfb54SMatthew G. Knepley   } else if (dim == 2) {
369371c9d4SSatish Balay     left = dd->xs - dd->Xs;
379371c9d4SSatish Balay     down = dd->ys - dd->Ys;
389371c9d4SSatish Balay     up   = down + dd->ye - dd->ys;
399566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1((dd->xe - dd->xs) * (up - down), &idx));
4047c6ae99SBarry Smith     count = 0;
4147c6ae99SBarry Smith     for (i = down; i < up; i++) {
42ad540459SPierre Jolivet       for (j = 0; j < dd->xe - dd->xs; j++) idx[count++] = left + i * (dd->Xe - dd->Xs) + j;
4347c6ae99SBarry Smith     }
44c73cfb54SMatthew G. Knepley   } else if (dim == 3) {
4547c6ae99SBarry Smith     left   = dd->xs - dd->Xs;
469371c9d4SSatish Balay     bottom = dd->ys - dd->Ys;
479371c9d4SSatish Balay     top    = bottom + dd->ye - dd->ys;
489371c9d4SSatish Balay     down   = dd->zs - dd->Zs;
499371c9d4SSatish Balay     up     = down + dd->ze - dd->zs;
5047c6ae99SBarry Smith     count  = (dd->xe - dd->xs) * (top - bottom) * (up - down);
519566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(count, &idx));
5247c6ae99SBarry Smith     count = 0;
5347c6ae99SBarry Smith     for (i = down; i < up; i++) {
5447c6ae99SBarry Smith       for (j = bottom; j < top; j++) {
55ad540459SPierre 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;
5647c6ae99SBarry Smith       }
5747c6ae99SBarry Smith     }
5863a3b9bcSJacob Faibussowitsch   } else SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_CORRUPT, "DMDA has invalid dimension %" PetscInt_FMT, dim);
5947c6ae99SBarry Smith 
609566063dSJacob Faibussowitsch   PetscCall(VecScatterRemap(dd->ltol, idx, NULL));
619566063dSJacob Faibussowitsch   PetscCall(PetscFree(idx));
623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6347c6ae99SBarry Smith }
6447c6ae99SBarry Smith 
65d71ae5a4SJacob Faibussowitsch PetscErrorCode DMLocalToLocalBegin_DA(DM da, Vec g, InsertMode mode, Vec l)
66d71ae5a4SJacob Faibussowitsch {
6747c6ae99SBarry Smith   DM_DA *dd = (DM_DA *)da->data;
6847c6ae99SBarry Smith 
6947c6ae99SBarry Smith   PetscFunctionBegin;
7047c6ae99SBarry Smith   PetscValidHeaderSpecific(da, DM_CLASSID, 1);
7148a46eb9SPierre Jolivet   if (!dd->ltol) PetscCall(DMLocalToLocalCreate_DA(da));
729566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(dd->ltol, g, l, mode, SCATTER_FORWARD));
733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7447c6ae99SBarry Smith }
7547c6ae99SBarry Smith 
76d71ae5a4SJacob Faibussowitsch PetscErrorCode DMLocalToLocalEnd_DA(DM da, Vec g, InsertMode mode, Vec l)
77d71ae5a4SJacob Faibussowitsch {
7847c6ae99SBarry Smith   DM_DA *dd = (DM_DA *)da->data;
7947c6ae99SBarry Smith 
8047c6ae99SBarry Smith   PetscFunctionBegin;
8147c6ae99SBarry Smith   PetscValidHeaderSpecific(da, DM_CLASSID, 1);
8247c6ae99SBarry Smith   PetscValidHeaderSpecific(g, VEC_CLASSID, 2);
839566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(dd->ltol, g, l, mode, SCATTER_FORWARD));
843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8547c6ae99SBarry Smith }
86