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