xref: /petsc/src/dm/impls/da/daltol.c (revision fe998a80077c9ee0917a39496df43fc256e1b478)
1 
2 /*
3   Code for manipulating distributed regular arrays in parallel.
4 */
5 
6 #include <petsc-private/dmdaimpl.h>    /*I   "petscdmda.h"   I*/
7 
8 #undef __FUNCT__
9 #define __FUNCT__ "DMLocalToLocalCreate_DA"
10 /*
11    DMLocalToLocalCreate_DA - Creates the local to local scatter
12 
13    Collective on DMDA
14 
15    Input Parameter:
16 .  da - the distributed array
17 
18 */
19 PetscErrorCode  DMLocalToLocalCreate_DA(DM da)
20 {
21   PetscErrorCode ierr;
22   PetscInt       *idx,left,j,count,up,down,i,bottom,top,k,dim=da->dim;
23   DM_DA          *dd = (DM_DA*)da->data;
24 
25   PetscFunctionBegin;
26   PetscValidHeaderSpecific(da,DM_CLASSID,1);
27 
28   if (dd->ltol) PetscFunctionReturn(0);
29   /*
30      We simply remap the values in the from part of
31      global to local to read from an array with the ghost values
32      rather then from the plain array.
33   */
34   ierr = VecScatterCopy(dd->gtol,&dd->ltol);CHKERRQ(ierr);
35   ierr = PetscLogObjectParent((PetscObject)da,(PetscObject)dd->ltol);CHKERRQ(ierr);
36   if (dim == 1) {
37     left = dd->xs - dd->Xs;
38     ierr = PetscMalloc1(dd->xe-dd->xs,&idx);CHKERRQ(ierr);
39     for (j=0; j<dd->xe-dd->xs; j++) idx[j] = left + j;
40   } else if (dim == 2) {
41     left  = dd->xs - dd->Xs; down  = dd->ys - dd->Ys; up    = down + dd->ye-dd->ys;
42     ierr  = PetscMalloc1((dd->xe-dd->xs)*(up - down),&idx);CHKERRQ(ierr);
43     count = 0;
44     for (i=down; i<up; i++) {
45       for (j=0; j<dd->xe-dd->xs; j++) {
46         idx[count++] = left + i*(dd->Xe-dd->Xs) + j;
47       }
48     }
49   } else if (dim == 3) {
50     left   = dd->xs - dd->Xs;
51     bottom = dd->ys - dd->Ys; top = bottom + dd->ye-dd->ys;
52     down   = dd->zs - dd->Zs; up  = down + dd->ze-dd->zs;
53     count  = (dd->xe-dd->xs)*(top-bottom)*(up-down);
54     ierr   = PetscMalloc1(count,&idx);CHKERRQ(ierr);
55     count  = 0;
56     for (i=down; i<up; i++) {
57       for (j=bottom; j<top; j++) {
58         for (k=0; k<dd->xe-dd->xs; k++) {
59           idx[count++] = (left+j*(dd->Xe-dd->Xs))+i*(dd->Xe-dd->Xs)*(dd->Ye-dd->Ys) + k;
60         }
61       }
62     }
63   } else SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_CORRUPT,"DMDA has invalid dimension %D",dim);
64 
65   ierr = VecScatterRemap(dd->ltol,idx,NULL);CHKERRQ(ierr);
66   ierr = PetscFree(idx);CHKERRQ(ierr);
67   PetscFunctionReturn(0);
68 }
69 
70 #undef __FUNCT__
71 #define __FUNCT__ "DMLocalToLocalBegin_DA"
72 /*
73    DMLocalToLocalBegin_DA - Maps from a local vector (including ghost points
74    that contain irrelevant values) to another local vector where the ghost
75    points in the second are set correctly. Must be followed by DMLocalToLocalEnd_DA().
76 
77    Neighbor-wise Collective on DMDA and Vec
78 
79    Input Parameters:
80 +  da - the distributed array context
81 .  g - the original local vector
82 -  mode - one of INSERT_VALUES or ADD_VALUES
83 
84    Output Parameter:
85 .  l  - the local vector with correct ghost values
86 
87    Notes:
88    The local vectors used here need not be the same as those
89    obtained from DMCreateLocalVector(), BUT they
90    must have the same parallel data layout; they could, for example, be
91    obtained with VecDuplicate() from the DMDA originating vectors.
92 
93 */
94 PetscErrorCode  DMLocalToLocalBegin_DA(DM da,Vec g,InsertMode mode,Vec l)
95 {
96   PetscErrorCode ierr;
97   DM_DA          *dd = (DM_DA*)da->data;
98 
99   PetscFunctionBegin;
100   PetscValidHeaderSpecific(da,DM_CLASSID,1);
101   if (!dd->ltol) {
102     ierr = DMLocalToLocalCreate_DA(da);CHKERRQ(ierr);
103   }
104   ierr = VecScatterBegin(dd->ltol,g,l,mode,SCATTER_FORWARD);CHKERRQ(ierr);
105   PetscFunctionReturn(0);
106 }
107 
108 #undef __FUNCT__
109 #define __FUNCT__ "DMLocalToLocalEnd_DA"
110 /*
111    DMLocalToLocalEnd_DA - Maps from a local vector (including ghost points
112    that contain irrelevant values) to another local vector where the ghost
113    points in the second are set correctly.  Must be preceeded by
114    DMLocalToLocalBegin_DA().
115 
116    Neighbor-wise Collective on DMDA and Vec
117 
118    Input Parameters:
119 +  da - the distributed array context
120 .  g - the original local vector
121 -  mode - one of INSERT_VALUES or ADD_VALUES
122 
123    Output Parameter:
124 .  l  - the local vector with correct ghost values
125 
126    Note:
127    The local vectors used here need not be the same as those
128    obtained from DMCreateLocalVector(), BUT they
129    must have the same parallel data layout; they could, for example, be
130    obtained with VecDuplicate() from the DMDA originating vectors.
131 
132 */
133 PetscErrorCode  DMLocalToLocalEnd_DA(DM da,Vec g,InsertMode mode,Vec l)
134 {
135   PetscErrorCode ierr;
136   DM_DA          *dd = (DM_DA*)da->data;
137 
138   PetscFunctionBegin;
139   PetscValidHeaderSpecific(da,DM_CLASSID,1);
140   PetscValidHeaderSpecific(g,VEC_CLASSID,2);
141   PetscValidHeaderSpecific(g,VEC_CLASSID,4);
142   ierr = VecScatterEnd(dd->ltol,g,l,mode,SCATTER_FORWARD);CHKERRQ(ierr);
143   PetscFunctionReturn(0);
144 }
145 
146