xref: /petsc/src/dm/impls/da/daltol.c (revision 63a3b9bc7a1f24f247904ccba9383635fe6abade)
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 */
17d78e899eSRichard Tran Mills PetscErrorCode  DMLocalToLocalCreate_DA(DM da)
1847c6ae99SBarry Smith {
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 
2547c6ae99SBarry Smith   if (dd->ltol) PetscFunctionReturn(0);
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));
329566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)da,(PetscObject)dd->ltol));
33c73cfb54SMatthew G. Knepley   if (dim == 1) {
3447c6ae99SBarry Smith     left = dd->xs - dd->Xs;
359566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(dd->xe-dd->xs,&idx));
368865f1eaSKarl Rupp     for (j=0; j<dd->xe-dd->xs; j++) idx[j] = left + j;
37c73cfb54SMatthew G. Knepley   } else if (dim == 2) {
3847c6ae99SBarry Smith     left  = dd->xs - dd->Xs; down  = dd->ys - dd->Ys; 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++) {
4247c6ae99SBarry Smith       for (j=0; j<dd->xe-dd->xs; j++) {
4347c6ae99SBarry Smith         idx[count++] = left + i*(dd->Xe-dd->Xs) + j;
4447c6ae99SBarry Smith       }
4547c6ae99SBarry Smith     }
46c73cfb54SMatthew G. Knepley   } else if (dim == 3) {
4747c6ae99SBarry Smith     left   = dd->xs - dd->Xs;
4847c6ae99SBarry Smith     bottom = dd->ys - dd->Ys; top = bottom + dd->ye-dd->ys;
4947c6ae99SBarry Smith     down   = dd->zs - dd->Zs; 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++) {
5547c6ae99SBarry Smith         for (k=0; k<dd->xe-dd->xs; k++) {
5647c6ae99SBarry Smith           idx[count++] = (left+j*(dd->Xe-dd->Xs))+i*(dd->Xe-dd->Xs)*(dd->Ye-dd->Ys) + k;
5747c6ae99SBarry Smith         }
5847c6ae99SBarry Smith       }
5947c6ae99SBarry Smith     }
60*63a3b9bcSJacob Faibussowitsch   } else SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_CORRUPT,"DMDA has invalid dimension %" PetscInt_FMT,dim);
6147c6ae99SBarry Smith 
629566063dSJacob Faibussowitsch   PetscCall(VecScatterRemap(dd->ltol,idx,NULL));
639566063dSJacob Faibussowitsch   PetscCall(PetscFree(idx));
6447c6ae99SBarry Smith   PetscFunctionReturn(0);
6547c6ae99SBarry Smith }
6647c6ae99SBarry Smith 
67d78e899eSRichard Tran Mills /*
68d78e899eSRichard Tran Mills    DMLocalToLocalBegin_DA - Maps from a local vector (including ghost points
6947c6ae99SBarry Smith    that contain irrelevant values) to another local vector where the ghost
70d78e899eSRichard Tran Mills    points in the second are set correctly. Must be followed by DMLocalToLocalEnd_DA().
7147c6ae99SBarry Smith 
72d083f849SBarry Smith    Neighbor-wise Collective on da
7347c6ae99SBarry Smith 
7447c6ae99SBarry Smith    Input Parameters:
7547c6ae99SBarry Smith +  da - the distributed array context
7647c6ae99SBarry Smith .  g - the original local vector
7747c6ae99SBarry Smith -  mode - one of INSERT_VALUES or ADD_VALUES
7847c6ae99SBarry Smith 
7947c6ae99SBarry Smith    Output Parameter:
8047c6ae99SBarry Smith .  l  - the local vector with correct ghost values
8147c6ae99SBarry Smith 
8247c6ae99SBarry Smith    Notes:
8347c6ae99SBarry Smith    The local vectors used here need not be the same as those
84564755cdSBarry Smith    obtained from DMCreateLocalVector(), BUT they
8547c6ae99SBarry Smith    must have the same parallel data layout; they could, for example, be
86aa219208SBarry Smith    obtained with VecDuplicate() from the DMDA originating vectors.
8747c6ae99SBarry Smith 
88d78e899eSRichard Tran Mills */
89d78e899eSRichard Tran Mills PetscErrorCode  DMLocalToLocalBegin_DA(DM da,Vec g,InsertMode mode,Vec l)
9047c6ae99SBarry Smith {
9147c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
9247c6ae99SBarry Smith 
9347c6ae99SBarry Smith   PetscFunctionBegin;
9447c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
9547c6ae99SBarry Smith   if (!dd->ltol) {
969566063dSJacob Faibussowitsch     PetscCall(DMLocalToLocalCreate_DA(da));
9747c6ae99SBarry Smith   }
989566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(dd->ltol,g,l,mode,SCATTER_FORWARD));
9947c6ae99SBarry Smith   PetscFunctionReturn(0);
10047c6ae99SBarry Smith }
10147c6ae99SBarry Smith 
102d78e899eSRichard Tran Mills /*
103d78e899eSRichard Tran Mills    DMLocalToLocalEnd_DA - Maps from a local vector (including ghost points
10447c6ae99SBarry Smith    that contain irrelevant values) to another local vector where the ghost
105a5b23f4aSJose E. Roman    points in the second are set correctly.  Must be preceded by
106d78e899eSRichard Tran Mills    DMLocalToLocalBegin_DA().
10747c6ae99SBarry Smith 
108d083f849SBarry Smith    Neighbor-wise Collective on da
10947c6ae99SBarry Smith 
11047c6ae99SBarry Smith    Input Parameters:
11147c6ae99SBarry Smith +  da - the distributed array context
11247c6ae99SBarry Smith .  g - the original local vector
11347c6ae99SBarry Smith -  mode - one of INSERT_VALUES or ADD_VALUES
11447c6ae99SBarry Smith 
11547c6ae99SBarry Smith    Output Parameter:
11647c6ae99SBarry Smith .  l  - the local vector with correct ghost values
11747c6ae99SBarry Smith 
11847c6ae99SBarry Smith    Note:
11947c6ae99SBarry Smith    The local vectors used here need not be the same as those
120564755cdSBarry Smith    obtained from DMCreateLocalVector(), BUT they
12147c6ae99SBarry Smith    must have the same parallel data layout; they could, for example, be
122aa219208SBarry Smith    obtained with VecDuplicate() from the DMDA originating vectors.
12347c6ae99SBarry Smith 
124d78e899eSRichard Tran Mills */
125d78e899eSRichard Tran Mills PetscErrorCode  DMLocalToLocalEnd_DA(DM da,Vec g,InsertMode mode,Vec l)
12647c6ae99SBarry Smith {
12747c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
12847c6ae99SBarry Smith 
12947c6ae99SBarry Smith   PetscFunctionBegin;
13047c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
13147c6ae99SBarry Smith   PetscValidHeaderSpecific(g,VEC_CLASSID,2);
1329566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(dd->ltol,g,l,mode,SCATTER_FORWARD));
13347c6ae99SBarry Smith   PetscFunctionReturn(0);
13447c6ae99SBarry Smith }
135