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__ "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++) idx[j] = left + j; 40 } else if (dd->dim == 2) { 41 left = dd->xs - dd->Xs; down = dd->ys - dd->Ys; up = down + dd->ye-dd->ys; 42 ierr = PetscMalloc((dd->xe-dd->xs)*(up - down)*sizeof(PetscInt),&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 (dd->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 = PetscMalloc(count*sizeof(PetscInt),&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",dd->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__ "DMDALocalToLocalBegin" 72 /*@ 73 DMDALocalToLocalBegin - 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 DMDALocalToLocalEnd(). 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 Level: intermediate 88 89 Notes: 90 The local vectors used here need not be the same as those 91 obtained from DMCreateLocalVector(), BUT they 92 must have the same parallel data layout; they could, for example, be 93 obtained with VecDuplicate() from the DMDA originating vectors. 94 95 .keywords: distributed array, local-to-local, begin 96 97 .seealso: DMDALocalToLocalEnd(), DMLocalToGlobalBegin() 98 @*/ 99 PetscErrorCode DMDALocalToLocalBegin(DM da,Vec g,InsertMode mode,Vec l) 100 { 101 PetscErrorCode ierr; 102 DM_DA *dd = (DM_DA*)da->data; 103 104 PetscFunctionBegin; 105 PetscValidHeaderSpecific(da,DM_CLASSID,1); 106 if (!dd->ltol) { 107 ierr = DMDALocalToLocalCreate(da);CHKERRQ(ierr); 108 } 109 ierr = VecScatterBegin(dd->ltol,g,l,mode,SCATTER_FORWARD);CHKERRQ(ierr); 110 PetscFunctionReturn(0); 111 } 112 113 #undef __FUNCT__ 114 #define __FUNCT__ "DMDALocalToLocalEnd" 115 /*@ 116 DMDALocalToLocalEnd - Maps from a local vector (including ghost points 117 that contain irrelevant values) to another local vector where the ghost 118 points in the second are set correctly. Must be preceeded by 119 DMDALocalToLocalBegin(). 120 121 Neighbor-wise Collective on DMDA and Vec 122 123 Input Parameters: 124 + da - the distributed array context 125 . g - the original local vector 126 - mode - one of INSERT_VALUES or ADD_VALUES 127 128 Output Parameter: 129 . l - the local vector with correct ghost values 130 131 Level: intermediate 132 133 Note: 134 The local vectors used here need not be the same as those 135 obtained from DMCreateLocalVector(), BUT they 136 must have the same parallel data layout; they could, for example, be 137 obtained with VecDuplicate() from the DMDA originating vectors. 138 139 .keywords: distributed array, local-to-local, end 140 141 .seealso: DMDALocalToLocalBegin(), DMLocalToGlobalBegin() 142 @*/ 143 PetscErrorCode DMDALocalToLocalEnd(DM da,Vec g,InsertMode mode,Vec l) 144 { 145 PetscErrorCode ierr; 146 DM_DA *dd = (DM_DA*)da->data; 147 148 PetscFunctionBegin; 149 PetscValidHeaderSpecific(da,DM_CLASSID,1); 150 PetscValidHeaderSpecific(g,VEC_CLASSID,2); 151 PetscValidHeaderSpecific(g,VEC_CLASSID,4); 152 ierr = VecScatterEnd(dd->ltol,g,l,mode,SCATTER_FORWARD);CHKERRQ(ierr); 153 PetscFunctionReturn(0); 154 } 155 156