xref: /petsc/src/dm/impls/da/daltol.c (revision d083f849a86f1f43e18d534ee43954e2786cb29a)
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*d083f849SBarry 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 {
1947c6ae99SBarry Smith   PetscErrorCode ierr;
20c73cfb54SMatthew G. Knepley   PetscInt       *idx,left,j,count,up,down,i,bottom,top,k,dim=da->dim;
2147c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
2247c6ae99SBarry Smith 
2347c6ae99SBarry Smith   PetscFunctionBegin;
2447c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
2547c6ae99SBarry Smith 
2647c6ae99SBarry Smith   if (dd->ltol) PetscFunctionReturn(0);
2747c6ae99SBarry Smith   /*
2847c6ae99SBarry Smith      We simply remap the values in the from part of
2947c6ae99SBarry Smith      global to local to read from an array with the ghost values
3047c6ae99SBarry Smith      rather then from the plain array.
3147c6ae99SBarry Smith   */
3247c6ae99SBarry Smith   ierr = VecScatterCopy(dd->gtol,&dd->ltol);CHKERRQ(ierr);
333bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)da,(PetscObject)dd->ltol);CHKERRQ(ierr);
34c73cfb54SMatthew G. Knepley   if (dim == 1) {
3547c6ae99SBarry Smith     left = dd->xs - dd->Xs;
36854ce69bSBarry Smith     ierr = PetscMalloc1(dd->xe-dd->xs,&idx);CHKERRQ(ierr);
378865f1eaSKarl Rupp     for (j=0; j<dd->xe-dd->xs; j++) idx[j] = left + j;
38c73cfb54SMatthew G. Knepley   } else if (dim == 2) {
3947c6ae99SBarry Smith     left  = dd->xs - dd->Xs; down  = dd->ys - dd->Ys; up    = down + dd->ye-dd->ys;
40785e854fSJed Brown     ierr  = PetscMalloc1((dd->xe-dd->xs)*(up - down),&idx);CHKERRQ(ierr);
4147c6ae99SBarry Smith     count = 0;
4247c6ae99SBarry Smith     for (i=down; i<up; i++) {
4347c6ae99SBarry Smith       for (j=0; j<dd->xe-dd->xs; j++) {
4447c6ae99SBarry Smith         idx[count++] = left + i*(dd->Xe-dd->Xs) + j;
4547c6ae99SBarry Smith       }
4647c6ae99SBarry Smith     }
47c73cfb54SMatthew G. Knepley   } else if (dim == 3) {
4847c6ae99SBarry Smith     left   = dd->xs - dd->Xs;
4947c6ae99SBarry Smith     bottom = dd->ys - dd->Ys; top = bottom + dd->ye-dd->ys;
5047c6ae99SBarry Smith     down   = dd->zs - dd->Zs; up  = down + dd->ze-dd->zs;
5147c6ae99SBarry Smith     count  = (dd->xe-dd->xs)*(top-bottom)*(up-down);
52785e854fSJed Brown     ierr   = PetscMalloc1(count,&idx);CHKERRQ(ierr);
5347c6ae99SBarry Smith     count  = 0;
5447c6ae99SBarry Smith     for (i=down; i<up; i++) {
5547c6ae99SBarry Smith       for (j=bottom; j<top; j++) {
5647c6ae99SBarry Smith         for (k=0; k<dd->xe-dd->xs; k++) {
5747c6ae99SBarry Smith           idx[count++] = (left+j*(dd->Xe-dd->Xs))+i*(dd->Xe-dd->Xs)*(dd->Ye-dd->Ys) + k;
5847c6ae99SBarry Smith         }
5947c6ae99SBarry Smith       }
6047c6ae99SBarry Smith     }
61c73cfb54SMatthew G. Knepley   } else SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_CORRUPT,"DMDA has invalid dimension %D",dim);
6247c6ae99SBarry Smith 
630298fd71SBarry Smith   ierr = VecScatterRemap(dd->ltol,idx,NULL);CHKERRQ(ierr);
6447c6ae99SBarry Smith   ierr = PetscFree(idx);CHKERRQ(ierr);
6547c6ae99SBarry Smith   PetscFunctionReturn(0);
6647c6ae99SBarry Smith }
6747c6ae99SBarry Smith 
68d78e899eSRichard Tran Mills /*
69d78e899eSRichard Tran Mills    DMLocalToLocalBegin_DA - Maps from a local vector (including ghost points
7047c6ae99SBarry Smith    that contain irrelevant values) to another local vector where the ghost
71d78e899eSRichard Tran Mills    points in the second are set correctly. Must be followed by DMLocalToLocalEnd_DA().
7247c6ae99SBarry Smith 
73*d083f849SBarry Smith    Neighbor-wise Collective on da
7447c6ae99SBarry Smith 
7547c6ae99SBarry Smith    Input Parameters:
7647c6ae99SBarry Smith +  da - the distributed array context
7747c6ae99SBarry Smith .  g - the original local vector
7847c6ae99SBarry Smith -  mode - one of INSERT_VALUES or ADD_VALUES
7947c6ae99SBarry Smith 
8047c6ae99SBarry Smith    Output Parameter:
8147c6ae99SBarry Smith .  l  - the local vector with correct ghost values
8247c6ae99SBarry Smith 
8347c6ae99SBarry Smith    Notes:
8447c6ae99SBarry Smith    The local vectors used here need not be the same as those
85564755cdSBarry Smith    obtained from DMCreateLocalVector(), BUT they
8647c6ae99SBarry Smith    must have the same parallel data layout; they could, for example, be
87aa219208SBarry Smith    obtained with VecDuplicate() from the DMDA originating vectors.
8847c6ae99SBarry Smith 
89d78e899eSRichard Tran Mills */
90d78e899eSRichard Tran Mills PetscErrorCode  DMLocalToLocalBegin_DA(DM da,Vec g,InsertMode mode,Vec l)
9147c6ae99SBarry Smith {
9247c6ae99SBarry Smith   PetscErrorCode ierr;
9347c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
9447c6ae99SBarry Smith 
9547c6ae99SBarry Smith   PetscFunctionBegin;
9647c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
9747c6ae99SBarry Smith   if (!dd->ltol) {
98d78e899eSRichard Tran Mills     ierr = DMLocalToLocalCreate_DA(da);CHKERRQ(ierr);
9947c6ae99SBarry Smith   }
10047c6ae99SBarry Smith   ierr = VecScatterBegin(dd->ltol,g,l,mode,SCATTER_FORWARD);CHKERRQ(ierr);
10147c6ae99SBarry Smith   PetscFunctionReturn(0);
10247c6ae99SBarry Smith }
10347c6ae99SBarry Smith 
104d78e899eSRichard Tran Mills /*
105d78e899eSRichard Tran Mills    DMLocalToLocalEnd_DA - Maps from a local vector (including ghost points
10647c6ae99SBarry Smith    that contain irrelevant values) to another local vector where the ghost
10747c6ae99SBarry Smith    points in the second are set correctly.  Must be preceeded by
108d78e899eSRichard Tran Mills    DMLocalToLocalBegin_DA().
10947c6ae99SBarry Smith 
110*d083f849SBarry Smith    Neighbor-wise Collective on da
11147c6ae99SBarry Smith 
11247c6ae99SBarry Smith    Input Parameters:
11347c6ae99SBarry Smith +  da - the distributed array context
11447c6ae99SBarry Smith .  g - the original local vector
11547c6ae99SBarry Smith -  mode - one of INSERT_VALUES or ADD_VALUES
11647c6ae99SBarry Smith 
11747c6ae99SBarry Smith    Output Parameter:
11847c6ae99SBarry Smith .  l  - the local vector with correct ghost values
11947c6ae99SBarry Smith 
12047c6ae99SBarry Smith    Note:
12147c6ae99SBarry Smith    The local vectors used here need not be the same as those
122564755cdSBarry Smith    obtained from DMCreateLocalVector(), BUT they
12347c6ae99SBarry Smith    must have the same parallel data layout; they could, for example, be
124aa219208SBarry Smith    obtained with VecDuplicate() from the DMDA originating vectors.
12547c6ae99SBarry Smith 
126d78e899eSRichard Tran Mills */
127d78e899eSRichard Tran Mills PetscErrorCode  DMLocalToLocalEnd_DA(DM da,Vec g,InsertMode mode,Vec l)
12847c6ae99SBarry Smith {
12947c6ae99SBarry Smith   PetscErrorCode ierr;
13047c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
13147c6ae99SBarry Smith 
13247c6ae99SBarry Smith   PetscFunctionBegin;
13347c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
13447c6ae99SBarry Smith   PetscValidHeaderSpecific(g,VEC_CLASSID,2);
13547c6ae99SBarry Smith   PetscValidHeaderSpecific(g,VEC_CLASSID,4);
13647c6ae99SBarry Smith   ierr = VecScatterEnd(dd->ltol,g,l,mode,SCATTER_FORWARD);CHKERRQ(ierr);
13747c6ae99SBarry Smith   PetscFunctionReturn(0);
13847c6ae99SBarry Smith }
13947c6ae99SBarry Smith 
140