xref: /petsc/src/dm/impls/da/daltol.c (revision 9a42bb27a39f0cdf3306a1e22d33cd9809484eaa)
147c6ae99SBarry Smith #define PETSCDM_DLL
247c6ae99SBarry Smith 
347c6ae99SBarry Smith /*
447c6ae99SBarry Smith   Code for manipulating distributed regular arrays in parallel.
547c6ae99SBarry Smith */
647c6ae99SBarry Smith 
747c6ae99SBarry Smith #include "private/daimpl.h"    /*I   "petscda.h"   I*/
847c6ae99SBarry Smith 
947c6ae99SBarry Smith #undef __FUNCT__
10*9a42bb27SBarry Smith #define __FUNCT__ "DMDALocalToLocalCreate"
1147c6ae99SBarry Smith /*
12*9a42bb27SBarry Smith    DMDALocalToLocalCreate - Creates the local to local scatter
1347c6ae99SBarry Smith 
1447c6ae99SBarry Smith    Collective on DA
1547c6ae99SBarry Smith 
1647c6ae99SBarry Smith    Input Parameter:
1747c6ae99SBarry Smith .  da - the distributed array
1847c6ae99SBarry Smith 
1947c6ae99SBarry Smith */
20*9a42bb27SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDALocalToLocalCreate(DM da)
2147c6ae99SBarry Smith {
2247c6ae99SBarry Smith   PetscErrorCode ierr;
2347c6ae99SBarry Smith   PetscInt       *idx,left,j,count,up,down,i,bottom,top,k;
2447c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
2547c6ae99SBarry Smith 
2647c6ae99SBarry Smith   PetscFunctionBegin;
2747c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
2847c6ae99SBarry Smith 
2947c6ae99SBarry Smith   if (dd->ltol) PetscFunctionReturn(0);
3047c6ae99SBarry Smith   /*
3147c6ae99SBarry Smith      We simply remap the values in the from part of
3247c6ae99SBarry Smith      global to local to read from an array with the ghost values
3347c6ae99SBarry Smith      rather then from the plain array.
3447c6ae99SBarry Smith   */
3547c6ae99SBarry Smith   ierr = VecScatterCopy(dd->gtol,&dd->ltol);CHKERRQ(ierr);
3647c6ae99SBarry Smith   ierr = PetscLogObjectParent(da,dd->ltol);CHKERRQ(ierr);
3747c6ae99SBarry Smith   if (dd->dim == 1) {
3847c6ae99SBarry Smith     left = dd->xs - dd->Xs;
3947c6ae99SBarry Smith     ierr = PetscMalloc((dd->xe-dd->xs)*sizeof(PetscInt),&idx);CHKERRQ(ierr);
4047c6ae99SBarry Smith     for (j=0; j<dd->xe-dd->xs; j++) {
4147c6ae99SBarry Smith       idx[j] = left + j;
4247c6ae99SBarry Smith     }
4347c6ae99SBarry Smith   } else if (dd->dim == 2) {
4447c6ae99SBarry Smith     left  = dd->xs - dd->Xs; down  = dd->ys - dd->Ys; up    = down + dd->ye-dd->ys;
4547c6ae99SBarry Smith     ierr = PetscMalloc((dd->xe-dd->xs)*(up - down)*sizeof(PetscInt),&idx);CHKERRQ(ierr);
4647c6ae99SBarry Smith     count = 0;
4747c6ae99SBarry Smith     for (i=down; i<up; i++) {
4847c6ae99SBarry Smith       for (j=0; j<dd->xe-dd->xs; j++) {
4947c6ae99SBarry Smith 	idx[count++] = left + i*(dd->Xe-dd->Xs) + j;
5047c6ae99SBarry Smith       }
5147c6ae99SBarry Smith     }
5247c6ae99SBarry Smith   } else if (dd->dim == 3) {
5347c6ae99SBarry Smith     left   = dd->xs - dd->Xs;
5447c6ae99SBarry Smith     bottom = dd->ys - dd->Ys; top = bottom + dd->ye-dd->ys ;
5547c6ae99SBarry Smith     down   = dd->zs - dd->Zs; up  = down + dd->ze-dd->zs;
5647c6ae99SBarry Smith     count  = (dd->xe-dd->xs)*(top-bottom)*(up-down);
5747c6ae99SBarry Smith     ierr = PetscMalloc(count*sizeof(PetscInt),&idx);CHKERRQ(ierr);
5847c6ae99SBarry Smith     count  = 0;
5947c6ae99SBarry Smith     for (i=down; i<up; i++) {
6047c6ae99SBarry Smith       for (j=bottom; j<top; j++) {
6147c6ae99SBarry Smith 	for (k=0; k<dd->xe-dd->xs; k++) {
6247c6ae99SBarry Smith 	  idx[count++] = (left+j*(dd->Xe-dd->Xs))+i*(dd->Xe-dd->Xs)*(dd->Ye-dd->Ys) + k;
6347c6ae99SBarry Smith 	}
6447c6ae99SBarry Smith       }
6547c6ae99SBarry Smith     }
6647c6ae99SBarry Smith   } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_CORRUPT,"DA has invalid dimension %D",dd->dim);
6747c6ae99SBarry Smith 
6847c6ae99SBarry Smith   ierr = VecScatterRemap(dd->ltol,idx,PETSC_NULL);CHKERRQ(ierr);
6947c6ae99SBarry Smith   ierr = PetscFree(idx);CHKERRQ(ierr);
7047c6ae99SBarry Smith   PetscFunctionReturn(0);
7147c6ae99SBarry Smith }
7247c6ae99SBarry Smith 
7347c6ae99SBarry Smith #undef __FUNCT__
7447c6ae99SBarry Smith #define __FUNCT__ "DALocalToLocalBegin"
7547c6ae99SBarry Smith /*@
7647c6ae99SBarry Smith    DALocalToLocalBegin - Maps from a local vector (including ghost points
7747c6ae99SBarry Smith    that contain irrelevant values) to another local vector where the ghost
7847c6ae99SBarry Smith    points in the second are set correctly. Must be followed by DALocalToLocalEnd().
7947c6ae99SBarry Smith 
8047c6ae99SBarry Smith    Neighbor-wise Collective on DA and Vec
8147c6ae99SBarry Smith 
8247c6ae99SBarry Smith    Input Parameters:
8347c6ae99SBarry Smith +  da - the distributed array context
8447c6ae99SBarry Smith .  g - the original local vector
8547c6ae99SBarry Smith -  mode - one of INSERT_VALUES or ADD_VALUES
8647c6ae99SBarry Smith 
8747c6ae99SBarry Smith    Output Parameter:
8847c6ae99SBarry Smith .  l  - the local vector with correct ghost values
8947c6ae99SBarry Smith 
9047c6ae99SBarry Smith    Level: intermediate
9147c6ae99SBarry Smith 
9247c6ae99SBarry Smith    Notes:
9347c6ae99SBarry Smith    The local vectors used here need not be the same as those
9447c6ae99SBarry Smith    obtained from DACreateLocalVector(), BUT they
9547c6ae99SBarry Smith    must have the same parallel data layout; they could, for example, be
9647c6ae99SBarry Smith    obtained with VecDuplicate() from the DA originating vectors.
9747c6ae99SBarry Smith 
9847c6ae99SBarry Smith .keywords: distributed array, local-to-local, begin
9947c6ae99SBarry Smith 
100*9a42bb27SBarry Smith .seealso: DALocalToLocalEnd(), DMLocalToGlobalBegin()
10147c6ae99SBarry Smith @*/
102*9a42bb27SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DALocalToLocalBegin(DM da,Vec g,InsertMode mode,Vec l)
10347c6ae99SBarry Smith {
10447c6ae99SBarry Smith   PetscErrorCode ierr;
10547c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
10647c6ae99SBarry Smith 
10747c6ae99SBarry Smith   PetscFunctionBegin;
10847c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
10947c6ae99SBarry Smith   if (!dd->ltol) {
110*9a42bb27SBarry Smith     ierr = DMDALocalToLocalCreate(da);CHKERRQ(ierr);
11147c6ae99SBarry Smith   }
11247c6ae99SBarry Smith   ierr = VecScatterBegin(dd->ltol,g,l,mode,SCATTER_FORWARD);CHKERRQ(ierr);
11347c6ae99SBarry Smith   PetscFunctionReturn(0);
11447c6ae99SBarry Smith }
11547c6ae99SBarry Smith 
11647c6ae99SBarry Smith #undef __FUNCT__
11747c6ae99SBarry Smith #define __FUNCT__ "DALocalToLocalEnd"
11847c6ae99SBarry Smith /*@
11947c6ae99SBarry Smith    DALocalToLocalEnd - Maps from a local vector (including ghost points
12047c6ae99SBarry Smith    that contain irrelevant values) to another local vector where the ghost
12147c6ae99SBarry Smith    points in the second are set correctly.  Must be preceeded by
12247c6ae99SBarry Smith    DALocalToLocalBegin().
12347c6ae99SBarry Smith 
12447c6ae99SBarry Smith    Neighbor-wise Collective on DA and Vec
12547c6ae99SBarry Smith 
12647c6ae99SBarry Smith    Input Parameters:
12747c6ae99SBarry Smith +  da - the distributed array context
12847c6ae99SBarry Smith .  g - the original local vector
12947c6ae99SBarry Smith -  mode - one of INSERT_VALUES or ADD_VALUES
13047c6ae99SBarry Smith 
13147c6ae99SBarry Smith    Output Parameter:
13247c6ae99SBarry Smith .  l  - the local vector with correct ghost values
13347c6ae99SBarry Smith 
13447c6ae99SBarry Smith    Level: intermediate
13547c6ae99SBarry Smith 
13647c6ae99SBarry Smith    Note:
13747c6ae99SBarry Smith    The local vectors used here need not be the same as those
13847c6ae99SBarry Smith    obtained from DACreateLocalVector(), BUT they
13947c6ae99SBarry Smith    must have the same parallel data layout; they could, for example, be
14047c6ae99SBarry Smith    obtained with VecDuplicate() from the DA originating vectors.
14147c6ae99SBarry Smith 
14247c6ae99SBarry Smith .keywords: distributed array, local-to-local, end
14347c6ae99SBarry Smith 
144*9a42bb27SBarry Smith .seealso: DALocalToLocalBegin(), DMLocalToGlobalBegin()
14547c6ae99SBarry Smith @*/
146*9a42bb27SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DALocalToLocalEnd(DM da,Vec g,InsertMode mode,Vec l)
14747c6ae99SBarry Smith {
14847c6ae99SBarry Smith   PetscErrorCode ierr;
14947c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
15047c6ae99SBarry Smith 
15147c6ae99SBarry Smith   PetscFunctionBegin;
15247c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
15347c6ae99SBarry Smith   PetscValidHeaderSpecific(g,VEC_CLASSID,2);
15447c6ae99SBarry Smith   PetscValidHeaderSpecific(g,VEC_CLASSID,4);
15547c6ae99SBarry Smith   ierr = VecScatterEnd(dd->ltol,g,l,mode,SCATTER_FORWARD);CHKERRQ(ierr);
15647c6ae99SBarry Smith   PetscFunctionReturn(0);
15747c6ae99SBarry Smith }
15847c6ae99SBarry Smith 
159