1454e267fSLisandro Dalcin 2af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h> /*I "petscdmda.h" I*/ 3454e267fSLisandro Dalcin 42dde6fd4SLisandro Dalcin static PetscErrorCode DMDAGetElements_1D(DM dm,PetscInt *nel,PetscInt *nen,const PetscInt *e[]) 5454e267fSLisandro Dalcin { 6454e267fSLisandro Dalcin PetscErrorCode ierr; 7454e267fSLisandro Dalcin DM_DA *da = (DM_DA*)dm->data; 8454e267fSLisandro Dalcin PetscInt i,xs,xe,Xs,Xe; 9454e267fSLisandro Dalcin PetscInt cnt=0; 105fd66863SKarl Rupp 11454e267fSLisandro Dalcin PetscFunctionBegin; 12454e267fSLisandro Dalcin if (!da->e) { 13f951a8fcSStefano Zampini PetscInt corners[2]; 14f951a8fcSStefano Zampini 15*7a8be351SBarry Smith PetscCheck(da->s,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Cannot get elements for DMDA with zero stencil width"); 16ea78f98cSLisandro Dalcin ierr = DMDAGetCorners(dm,&xs,NULL,NULL,&xe,NULL,NULL);CHKERRQ(ierr); 17ea78f98cSLisandro Dalcin ierr = DMDAGetGhostCorners(dm,&Xs,NULL,NULL,&Xe,NULL,NULL);CHKERRQ(ierr); 18454e267fSLisandro Dalcin xe += xs; Xe += Xs; if (xs != Xs) xs -= 1; 19454e267fSLisandro Dalcin da->ne = 1*(xe - xs - 1); 20854ce69bSBarry Smith ierr = PetscMalloc1(1 + 2*da->ne,&da->e);CHKERRQ(ierr); 21454e267fSLisandro Dalcin for (i=xs; i<xe-1; i++) { 22454e267fSLisandro Dalcin da->e[cnt++] = (i-Xs); 23454e267fSLisandro Dalcin da->e[cnt++] = (i-Xs+1); 24454e267fSLisandro Dalcin } 25cfbed8a3SStefano Zampini da->nen = 2; 26cfbed8a3SStefano Zampini 27f951a8fcSStefano Zampini corners[0] = (xs -Xs); 28f951a8fcSStefano Zampini corners[1] = (xe-1-Xs); 29f951a8fcSStefano Zampini ierr = ISCreateGeneral(PETSC_COMM_SELF,2,corners,PETSC_COPY_VALUES,&da->ecorners);CHKERRQ(ierr); 30454e267fSLisandro Dalcin } 31454e267fSLisandro Dalcin *nel = da->ne; 32cfbed8a3SStefano Zampini *nen = da->nen; 33454e267fSLisandro Dalcin *e = da->e; 34454e267fSLisandro Dalcin PetscFunctionReturn(0); 35454e267fSLisandro Dalcin } 36454e267fSLisandro Dalcin 372dde6fd4SLisandro Dalcin static PetscErrorCode DMDAGetElements_2D(DM dm,PetscInt *nel,PetscInt *nen,const PetscInt *e[]) 38454e267fSLisandro Dalcin { 39454e267fSLisandro Dalcin PetscErrorCode ierr; 40454e267fSLisandro Dalcin DM_DA *da = (DM_DA*)dm->data; 41454e267fSLisandro Dalcin PetscInt i,xs,xe,Xs,Xe; 42454e267fSLisandro Dalcin PetscInt j,ys,ye,Ys,Ye; 43cfbed8a3SStefano Zampini PetscInt cnt=0, cell[4], ns=2; 44454e267fSLisandro Dalcin PetscInt c, split[] = {0,1,3, 45454e267fSLisandro Dalcin 2,3,1}; 465fd66863SKarl Rupp 47454e267fSLisandro Dalcin PetscFunctionBegin; 48454e267fSLisandro Dalcin if (!da->e) { 49cfbed8a3SStefano Zampini PetscInt corners[4],nn = 0; 50f951a8fcSStefano Zampini 51*7a8be351SBarry Smith PetscCheck(da->s,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Cannot get elements for DMDA with zero stencil width"); 52cfbed8a3SStefano Zampini 53cfbed8a3SStefano Zampini switch (da->elementtype) { 54cfbed8a3SStefano Zampini case DMDA_ELEMENT_Q1: 55cfbed8a3SStefano Zampini da->nen = 4; 56cfbed8a3SStefano Zampini break; 57cfbed8a3SStefano Zampini case DMDA_ELEMENT_P1: 58cfbed8a3SStefano Zampini da->nen = 3; 59cfbed8a3SStefano Zampini break; 60cfbed8a3SStefano Zampini default: 6198921bdaSJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Unknown element type %d",da->elementtype); 62cfbed8a3SStefano Zampini } 63cfbed8a3SStefano Zampini nn = da->nen; 64cfbed8a3SStefano Zampini 65b644c393SDave May if (da->elementtype == DMDA_ELEMENT_P1) {ns=2;} 66b644c393SDave May if (da->elementtype == DMDA_ELEMENT_Q1) {ns=1;} 67ea78f98cSLisandro Dalcin ierr = DMDAGetCorners(dm,&xs,&ys,NULL,&xe,&ye,NULL);CHKERRQ(ierr); 68ea78f98cSLisandro Dalcin ierr = DMDAGetGhostCorners(dm,&Xs,&Ys,NULL,&Xe,&Ye,NULL);CHKERRQ(ierr); 69454e267fSLisandro Dalcin xe += xs; Xe += Xs; if (xs != Xs) xs -= 1; 70454e267fSLisandro Dalcin ye += ys; Ye += Ys; if (ys != Ys) ys -= 1; 71454e267fSLisandro Dalcin da->ne = ns*(xe - xs - 1)*(ye - ys - 1); 72854ce69bSBarry Smith ierr = PetscMalloc1(1 + nn*da->ne,&da->e);CHKERRQ(ierr); 73454e267fSLisandro Dalcin for (j=ys; j<ye-1; j++) { 74454e267fSLisandro Dalcin for (i=xs; i<xe-1; i++) { 75454e267fSLisandro Dalcin cell[0] = (i-Xs) + (j-Ys)*(Xe-Xs); 76454e267fSLisandro Dalcin cell[1] = (i-Xs+1) + (j-Ys)*(Xe-Xs); 77454e267fSLisandro Dalcin cell[2] = (i-Xs+1) + (j-Ys+1)*(Xe-Xs); 78454e267fSLisandro Dalcin cell[3] = (i-Xs) + (j-Ys+1)*(Xe-Xs); 79454e267fSLisandro Dalcin if (da->elementtype == DMDA_ELEMENT_P1) { 808865f1eaSKarl Rupp for (c=0; c<ns*nn; c++) da->e[cnt++] = cell[split[c]]; 81454e267fSLisandro Dalcin } 82454e267fSLisandro Dalcin if (da->elementtype == DMDA_ELEMENT_Q1) { 838865f1eaSKarl Rupp for (c=0; c<ns*nn; c++) da->e[cnt++] = cell[c]; 84454e267fSLisandro Dalcin } 85454e267fSLisandro Dalcin } 86454e267fSLisandro Dalcin } 87cfbed8a3SStefano Zampini 88f951a8fcSStefano Zampini corners[0] = (xs -Xs) + (ys -Ys)*(Xe-Xs); 89f951a8fcSStefano Zampini corners[1] = (xe-1-Xs) + (ys -Ys)*(Xe-Xs); 90f951a8fcSStefano Zampini corners[2] = (xs -Xs) + (ye-1-Ys)*(Xe-Xs); 91f951a8fcSStefano Zampini corners[3] = (xe-1-Xs) + (ye-1-Ys)*(Xe-Xs); 92f951a8fcSStefano Zampini ierr = ISCreateGeneral(PETSC_COMM_SELF,4,corners,PETSC_COPY_VALUES,&da->ecorners);CHKERRQ(ierr); 93454e267fSLisandro Dalcin } 94454e267fSLisandro Dalcin *nel = da->ne; 95cfbed8a3SStefano Zampini *nen = da->nen; 96454e267fSLisandro Dalcin *e = da->e; 97454e267fSLisandro Dalcin PetscFunctionReturn(0); 98454e267fSLisandro Dalcin } 99454e267fSLisandro Dalcin 1002dde6fd4SLisandro Dalcin static PetscErrorCode DMDAGetElements_3D(DM dm,PetscInt *nel,PetscInt *nen,const PetscInt *e[]) 101454e267fSLisandro Dalcin { 102454e267fSLisandro Dalcin PetscErrorCode ierr; 103454e267fSLisandro Dalcin DM_DA *da = (DM_DA*)dm->data; 104454e267fSLisandro Dalcin PetscInt i,xs,xe,Xs,Xe; 105454e267fSLisandro Dalcin PetscInt j,ys,ye,Ys,Ye; 106454e267fSLisandro Dalcin PetscInt k,zs,ze,Zs,Ze; 107cfbed8a3SStefano Zampini PetscInt cnt=0, cell[8], ns=6; 108454e267fSLisandro Dalcin PetscInt c, split[] = {0,1,3,7, 109454e267fSLisandro Dalcin 0,1,7,4, 110454e267fSLisandro Dalcin 1,2,3,7, 111454e267fSLisandro Dalcin 1,2,7,6, 112454e267fSLisandro Dalcin 1,4,5,7, 113454e267fSLisandro Dalcin 1,5,6,7}; 1145fd66863SKarl Rupp 115454e267fSLisandro Dalcin PetscFunctionBegin; 116454e267fSLisandro Dalcin if (!da->e) { 117cfbed8a3SStefano Zampini PetscInt corners[8],nn = 0; 118f951a8fcSStefano Zampini 119*7a8be351SBarry Smith PetscCheck(da->s,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Cannot get elements for DMDA with zero stencil width"); 120cfbed8a3SStefano Zampini 121cfbed8a3SStefano Zampini switch (da->elementtype) { 122cfbed8a3SStefano Zampini case DMDA_ELEMENT_Q1: 123cfbed8a3SStefano Zampini da->nen = 8; 124cfbed8a3SStefano Zampini break; 125cfbed8a3SStefano Zampini case DMDA_ELEMENT_P1: 126cfbed8a3SStefano Zampini da->nen = 4; 127cfbed8a3SStefano Zampini break; 128cfbed8a3SStefano Zampini default: 12998921bdaSJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Unknown element type %d",da->elementtype); 130cfbed8a3SStefano Zampini } 131cfbed8a3SStefano Zampini nn = da->nen; 132cfbed8a3SStefano Zampini 133b644c393SDave May if (da->elementtype == DMDA_ELEMENT_P1) {ns=6;} 134b644c393SDave May if (da->elementtype == DMDA_ELEMENT_Q1) {ns=1;} 135454e267fSLisandro Dalcin ierr = DMDAGetCorners(dm,&xs,&ys,&zs,&xe,&ye,&ze);CHKERRQ(ierr); 136454e267fSLisandro Dalcin ierr = DMDAGetGhostCorners(dm,&Xs,&Ys,&Zs,&Xe,&Ye,&Ze);CHKERRQ(ierr); 137454e267fSLisandro Dalcin xe += xs; Xe += Xs; if (xs != Xs) xs -= 1; 138454e267fSLisandro Dalcin ye += ys; Ye += Ys; if (ys != Ys) ys -= 1; 139454e267fSLisandro Dalcin ze += zs; Ze += Zs; if (zs != Zs) zs -= 1; 140454e267fSLisandro Dalcin da->ne = ns*(xe - xs - 1)*(ye - ys - 1)*(ze - zs - 1); 141854ce69bSBarry Smith ierr = PetscMalloc1(1 + nn*da->ne,&da->e);CHKERRQ(ierr); 142454e267fSLisandro Dalcin for (k=zs; k<ze-1; k++) { 143454e267fSLisandro Dalcin for (j=ys; j<ye-1; j++) { 144454e267fSLisandro Dalcin for (i=xs; i<xe-1; i++) { 145454e267fSLisandro Dalcin cell[0] = (i -Xs) + (j -Ys)*(Xe-Xs) + (k -Zs)*(Xe-Xs)*(Ye-Ys); 1462da392ccSBarry Smith cell[1] = (i+1-Xs) + (j -Ys)*(Xe-Xs) + (k -Zs)*(Xe-Xs)*(Ye-Ys); 1472da392ccSBarry Smith cell[2] = (i+1-Xs) + (j+1-Ys)*(Xe-Xs) + (k -Zs)*(Xe-Xs)*(Ye-Ys); 1482da392ccSBarry Smith cell[3] = (i -Xs) + (j+1-Ys)*(Xe-Xs) + (k -Zs)*(Xe-Xs)*(Ye-Ys); 1492da392ccSBarry Smith cell[4] = (i -Xs) + (j -Ys)*(Xe-Xs) + (k+1-Zs)*(Xe-Xs)*(Ye-Ys); 1502da392ccSBarry Smith cell[5] = (i+1-Xs) + (j -Ys)*(Xe-Xs) + (k+1-Zs)*(Xe-Xs)*(Ye-Ys); 1512da392ccSBarry Smith cell[6] = (i+1-Xs) + (j+1-Ys)*(Xe-Xs) + (k+1-Zs)*(Xe-Xs)*(Ye-Ys); 1522da392ccSBarry Smith cell[7] = (i -Xs) + (j+1-Ys)*(Xe-Xs) + (k+1-Zs)*(Xe-Xs)*(Ye-Ys); 153454e267fSLisandro Dalcin if (da->elementtype == DMDA_ELEMENT_P1) { 1548865f1eaSKarl Rupp for (c=0; c<ns*nn; c++) da->e[cnt++] = cell[split[c]]; 155454e267fSLisandro Dalcin } 156454e267fSLisandro Dalcin if (da->elementtype == DMDA_ELEMENT_Q1) { 1578865f1eaSKarl Rupp for (c=0; c<ns*nn; c++) da->e[cnt++] = cell[c]; 158454e267fSLisandro Dalcin } 159454e267fSLisandro Dalcin } 160454e267fSLisandro Dalcin } 161454e267fSLisandro Dalcin } 162cfbed8a3SStefano Zampini 163f951a8fcSStefano Zampini corners[0] = (xs -Xs) + (ys -Ys)*(Xe-Xs) + (zs- Zs)*(Xe-Xs)*(Ye-Ys); 164f951a8fcSStefano Zampini corners[1] = (xe-1-Xs) + (ys -Ys)*(Xe-Xs) + (zs- Zs)*(Xe-Xs)*(Ye-Ys); 165f951a8fcSStefano Zampini corners[2] = (xs -Xs) + (ye-1-Ys)*(Xe-Xs) + (zs- Zs)*(Xe-Xs)*(Ye-Ys); 166f951a8fcSStefano Zampini corners[3] = (xe-1-Xs) + (ye-1-Ys)*(Xe-Xs) + (zs- Zs)*(Xe-Xs)*(Ye-Ys); 167f951a8fcSStefano Zampini corners[4] = (xs -Xs) + (ys -Ys)*(Xe-Xs) + (ze-1-Zs)*(Xe-Xs)*(Ye-Ys); 168f951a8fcSStefano Zampini corners[5] = (xe-1-Xs) + (ys -Ys)*(Xe-Xs) + (ze-1-Zs)*(Xe-Xs)*(Ye-Ys); 169f951a8fcSStefano Zampini corners[6] = (xs -Xs) + (ye-1-Ys)*(Xe-Xs) + (ze-1-Zs)*(Xe-Xs)*(Ye-Ys); 170f951a8fcSStefano Zampini corners[7] = (xe-1-Xs) + (ye-1-Ys)*(Xe-Xs) + (ze-1-Zs)*(Xe-Xs)*(Ye-Ys); 171f951a8fcSStefano Zampini ierr = ISCreateGeneral(PETSC_COMM_SELF,8,corners,PETSC_COPY_VALUES,&da->ecorners);CHKERRQ(ierr); 172454e267fSLisandro Dalcin } 173454e267fSLisandro Dalcin *nel = da->ne; 174cfbed8a3SStefano Zampini *nen = da->nen; 175454e267fSLisandro Dalcin *e = da->e; 176454e267fSLisandro Dalcin PetscFunctionReturn(0); 177454e267fSLisandro Dalcin } 178454e267fSLisandro Dalcin 179f5f57ec0SBarry Smith /*@ 180d4a6ed37SStefano Zampini DMDAGetElementsCorners - Returns the global (x,y,z) indices of the lower left 181d4a6ed37SStefano Zampini corner of the non-overlapping decomposition identified by DMDAGetElements() 182d4a6ed37SStefano Zampini 183d4a6ed37SStefano Zampini Not Collective 184d4a6ed37SStefano Zampini 185d4a6ed37SStefano Zampini Input Parameter: 186d4a6ed37SStefano Zampini . da - the DM object 187d4a6ed37SStefano Zampini 188d4a6ed37SStefano Zampini Output Parameters: 189d4a6ed37SStefano Zampini + gx - the x index 190d4a6ed37SStefano Zampini . gy - the y index 191d4a6ed37SStefano Zampini - gz - the z index 192d4a6ed37SStefano Zampini 193d4a6ed37SStefano Zampini Level: intermediate 194d4a6ed37SStefano Zampini 195d4a6ed37SStefano Zampini Notes: 196d4a6ed37SStefano Zampini 197d4a6ed37SStefano Zampini .seealso: DMDAElementType, DMDASetElementType(), DMDAGetElements() 198d4a6ed37SStefano Zampini @*/ 199d4a6ed37SStefano Zampini PetscErrorCode DMDAGetElementsCorners(DM da, PetscInt *gx, PetscInt *gy, PetscInt *gz) 200d4a6ed37SStefano Zampini { 201d4a6ed37SStefano Zampini PetscInt xs,Xs; 202d4a6ed37SStefano Zampini PetscInt ys,Ys; 203d4a6ed37SStefano Zampini PetscInt zs,Zs; 204d4a6ed37SStefano Zampini PetscBool isda; 205d4a6ed37SStefano Zampini PetscErrorCode ierr; 206d4a6ed37SStefano Zampini 207d4a6ed37SStefano Zampini PetscFunctionBegin; 208a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 209d4a6ed37SStefano Zampini if (gx) PetscValidIntPointer(gx,2); 210d4a6ed37SStefano Zampini if (gy) PetscValidIntPointer(gy,3); 211d4a6ed37SStefano Zampini if (gz) PetscValidIntPointer(gz,4); 212d4a6ed37SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)da,DMDA,&isda);CHKERRQ(ierr); 213*7a8be351SBarry Smith PetscCheck(isda,PetscObjectComm((PetscObject)da),PETSC_ERR_USER,"Not for DM type %s",((PetscObject)da)->type_name); 214d4a6ed37SStefano Zampini ierr = DMDAGetCorners(da,&xs,&ys,&zs,NULL,NULL,NULL);CHKERRQ(ierr); 215d4a6ed37SStefano Zampini ierr = DMDAGetGhostCorners(da,&Xs,&Ys,&Zs,NULL,NULL,NULL);CHKERRQ(ierr); 216d4a6ed37SStefano Zampini if (xs != Xs) xs -= 1; 217d4a6ed37SStefano Zampini if (ys != Ys) ys -= 1; 218d4a6ed37SStefano Zampini if (zs != Zs) zs -= 1; 219d4a6ed37SStefano Zampini if (gx) *gx = xs; 220d4a6ed37SStefano Zampini if (gy) *gy = ys; 221d4a6ed37SStefano Zampini if (gz) *gz = zs; 222d4a6ed37SStefano Zampini PetscFunctionReturn(0); 223d4a6ed37SStefano Zampini } 224d4a6ed37SStefano Zampini 225d4a6ed37SStefano Zampini /*@ 226d4a6ed37SStefano Zampini DMDAGetElementsSizes - Gets the local number of elements per direction for the non-overlapping decomposition identified by DMDAGetElements() 227d4a6ed37SStefano Zampini 228d4a6ed37SStefano Zampini Not Collective 229d4a6ed37SStefano Zampini 230d4a6ed37SStefano Zampini Input Parameter: 231d4a6ed37SStefano Zampini . da - the DM object 232d4a6ed37SStefano Zampini 233d4a6ed37SStefano Zampini Output Parameters: 234d4a6ed37SStefano Zampini + mx - number of local elements in x-direction 235d4a6ed37SStefano Zampini . my - number of local elements in y-direction 236d4a6ed37SStefano Zampini - mz - number of local elements in z-direction 237d4a6ed37SStefano Zampini 238d4a6ed37SStefano Zampini Level: intermediate 239d4a6ed37SStefano Zampini 24095452b02SPatrick Sanan Notes: 24195452b02SPatrick Sanan It returns the same number of elements, irrespective of the DMDAElementType 242d4a6ed37SStefano Zampini 243d4a6ed37SStefano Zampini .seealso: DMDAElementType, DMDASetElementType(), DMDAGetElements 244d4a6ed37SStefano Zampini @*/ 245d4a6ed37SStefano Zampini PetscErrorCode DMDAGetElementsSizes(DM da, PetscInt *mx, PetscInt *my, PetscInt *mz) 246d4a6ed37SStefano Zampini { 247d4a6ed37SStefano Zampini PetscInt xs,xe,Xs; 248d4a6ed37SStefano Zampini PetscInt ys,ye,Ys; 249d4a6ed37SStefano Zampini PetscInt zs,ze,Zs; 250d4a6ed37SStefano Zampini PetscInt dim; 251d4a6ed37SStefano Zampini PetscBool isda; 252d4a6ed37SStefano Zampini PetscErrorCode ierr; 253d4a6ed37SStefano Zampini 254d4a6ed37SStefano Zampini PetscFunctionBegin; 255a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 256d4a6ed37SStefano Zampini if (mx) PetscValidIntPointer(mx,2); 257d4a6ed37SStefano Zampini if (my) PetscValidIntPointer(my,3); 258d4a6ed37SStefano Zampini if (mz) PetscValidIntPointer(mz,4); 259d4a6ed37SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)da,DMDA,&isda);CHKERRQ(ierr); 260*7a8be351SBarry Smith PetscCheck(isda,PetscObjectComm((PetscObject)da),PETSC_ERR_USER,"Not for DM type %s",((PetscObject)da)->type_name); 261d4a6ed37SStefano Zampini ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xe,&ye,&ze);CHKERRQ(ierr); 262d4a6ed37SStefano Zampini ierr = DMDAGetGhostCorners(da,&Xs,&Ys,&Zs,NULL,NULL,NULL);CHKERRQ(ierr); 263d4a6ed37SStefano Zampini xe += xs; if (xs != Xs) xs -= 1; 264d4a6ed37SStefano Zampini ye += ys; if (ys != Ys) ys -= 1; 265d4a6ed37SStefano Zampini ze += zs; if (zs != Zs) zs -= 1; 266d4a6ed37SStefano Zampini if (mx) *mx = 0; 267d4a6ed37SStefano Zampini if (my) *my = 0; 268d4a6ed37SStefano Zampini if (mz) *mz = 0; 269d4a6ed37SStefano Zampini ierr = DMGetDimension(da,&dim);CHKERRQ(ierr); 270d4a6ed37SStefano Zampini switch (dim) { 271d4a6ed37SStefano Zampini case 3: 272d4a6ed37SStefano Zampini if (mz) *mz = ze - zs - 1; 273d4a6ed37SStefano Zampini case 2: 274d4a6ed37SStefano Zampini if (my) *my = ye - ys - 1; 275d4a6ed37SStefano Zampini case 1: 276d4a6ed37SStefano Zampini if (mx) *mx = xe - xs - 1; 277d4a6ed37SStefano Zampini break; 278d4a6ed37SStefano Zampini } 279d4a6ed37SStefano Zampini PetscFunctionReturn(0); 280d4a6ed37SStefano Zampini } 281d4a6ed37SStefano Zampini 282d4a6ed37SStefano Zampini /*@ 2832dde6fd4SLisandro Dalcin DMDASetElementType - Sets the element type to be returned by DMDAGetElements() 2842dde6fd4SLisandro Dalcin 2852dde6fd4SLisandro Dalcin Not Collective 2862dde6fd4SLisandro Dalcin 2872dde6fd4SLisandro Dalcin Input Parameter: 2882dde6fd4SLisandro Dalcin . da - the DMDA object 2892dde6fd4SLisandro Dalcin 2902dde6fd4SLisandro Dalcin Output Parameters: 2916420f107SDave May . etype - the element type, currently either DMDA_ELEMENT_P1 or DMDA_ELEMENT_Q1 2922dde6fd4SLisandro Dalcin 2932dde6fd4SLisandro Dalcin Level: intermediate 2942dde6fd4SLisandro Dalcin 2952dde6fd4SLisandro Dalcin .seealso: DMDAElementType, DMDAGetElementType(), DMDAGetElements(), DMDARestoreElements() 2962dde6fd4SLisandro Dalcin @*/ 2972dde6fd4SLisandro Dalcin PetscErrorCode DMDASetElementType(DM da, DMDAElementType etype) 2982dde6fd4SLisandro Dalcin { 2992dde6fd4SLisandro Dalcin DM_DA *dd = (DM_DA*)da->data; 3002dde6fd4SLisandro Dalcin PetscErrorCode ierr; 301f951a8fcSStefano Zampini PetscBool isda; 3022dde6fd4SLisandro Dalcin 3032dde6fd4SLisandro Dalcin PetscFunctionBegin; 304a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 3052dde6fd4SLisandro Dalcin PetscValidLogicalCollectiveEnum(da,etype,2); 306f951a8fcSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)da,DMDA,&isda);CHKERRQ(ierr); 307f951a8fcSStefano Zampini if (!isda) PetscFunctionReturn(0); 3082dde6fd4SLisandro Dalcin if (dd->elementtype != etype) { 3092dde6fd4SLisandro Dalcin ierr = PetscFree(dd->e);CHKERRQ(ierr); 310f951a8fcSStefano Zampini ierr = ISDestroy(&dd->ecorners);CHKERRQ(ierr); 3118865f1eaSKarl Rupp 3122dde6fd4SLisandro Dalcin dd->elementtype = etype; 3132dde6fd4SLisandro Dalcin dd->ne = 0; 314cfbed8a3SStefano Zampini dd->nen = 0; 3150298fd71SBarry Smith dd->e = NULL; 3162dde6fd4SLisandro Dalcin } 3172dde6fd4SLisandro Dalcin PetscFunctionReturn(0); 3182dde6fd4SLisandro Dalcin } 3192dde6fd4SLisandro Dalcin 320f5f57ec0SBarry Smith /*@ 3212dde6fd4SLisandro Dalcin DMDAGetElementType - Gets the element type to be returned by DMDAGetElements() 3222dde6fd4SLisandro Dalcin 3232dde6fd4SLisandro Dalcin Not Collective 3242dde6fd4SLisandro Dalcin 3252dde6fd4SLisandro Dalcin Input Parameter: 3262dde6fd4SLisandro Dalcin . da - the DMDA object 3272dde6fd4SLisandro Dalcin 3282dde6fd4SLisandro Dalcin Output Parameters: 3296420f107SDave May . etype - the element type, currently either DMDA_ELEMENT_P1 or DMDA_ELEMENT_Q1 3302dde6fd4SLisandro Dalcin 3312dde6fd4SLisandro Dalcin Level: intermediate 3322dde6fd4SLisandro Dalcin 3332dde6fd4SLisandro Dalcin .seealso: DMDAElementType, DMDASetElementType(), DMDAGetElements(), DMDARestoreElements() 3342dde6fd4SLisandro Dalcin @*/ 3352dde6fd4SLisandro Dalcin PetscErrorCode DMDAGetElementType(DM da, DMDAElementType *etype) 3362dde6fd4SLisandro Dalcin { 3372dde6fd4SLisandro Dalcin DM_DA *dd = (DM_DA*)da->data; 338f951a8fcSStefano Zampini PetscErrorCode ierr; 339f951a8fcSStefano Zampini PetscBool isda; 3405fd66863SKarl Rupp 3412dde6fd4SLisandro Dalcin PetscFunctionBegin; 342a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 3432dde6fd4SLisandro Dalcin PetscValidPointer(etype,2); 344f951a8fcSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)da,DMDA,&isda);CHKERRQ(ierr); 345*7a8be351SBarry Smith PetscCheck(isda,PetscObjectComm((PetscObject)da),PETSC_ERR_USER,"Not for DM type %s",((PetscObject)da)->type_name); 3462dde6fd4SLisandro Dalcin *etype = dd->elementtype; 3472dde6fd4SLisandro Dalcin PetscFunctionReturn(0); 3482dde6fd4SLisandro Dalcin } 3492dde6fd4SLisandro Dalcin 3502dde6fd4SLisandro Dalcin /*@C 3512dde6fd4SLisandro Dalcin DMDAGetElements - Gets an array containing the indices (in local coordinates) 3522dde6fd4SLisandro Dalcin of all the local elements 3532dde6fd4SLisandro Dalcin 3542dde6fd4SLisandro Dalcin Not Collective 3552dde6fd4SLisandro Dalcin 3562dde6fd4SLisandro Dalcin Input Parameter: 3572dde6fd4SLisandro Dalcin . dm - the DM object 3582dde6fd4SLisandro Dalcin 3592dde6fd4SLisandro Dalcin Output Parameters: 3602dde6fd4SLisandro Dalcin + nel - number of local elements 3612dde6fd4SLisandro Dalcin . nen - number of element nodes 362c2cd2aa3SJed Brown - e - the local indices of the elements' vertices 3632dde6fd4SLisandro Dalcin 3642dde6fd4SLisandro Dalcin Level: intermediate 3652dde6fd4SLisandro Dalcin 36693386343SJed Brown Notes: 36793386343SJed Brown Call DMDARestoreElements() once you have finished accessing the elements. 36893386343SJed Brown 36945950774SSatish Balay Each process uniquely owns a subset of the elements. That is no element is owned by two or more processes. 370009f0576SBarry Smith 371009f0576SBarry Smith If on each process you integrate over its owned elements and use ADD_VALUES in Vec/MatSetValuesLocal() then you'll obtain the correct result. 372009f0576SBarry Smith 373f5f57ec0SBarry Smith Not supported in Fortran 374f5f57ec0SBarry Smith 37593386343SJed Brown .seealso: DMDAElementType, DMDASetElementType(), VecSetValuesLocal(), MatSetValuesLocal(), DMGlobalToLocalBegin(), DMLocalToGlobalBegin() 3762dde6fd4SLisandro Dalcin @*/ 3772dde6fd4SLisandro Dalcin PetscErrorCode DMDAGetElements(DM dm,PetscInt *nel,PetscInt *nen,const PetscInt *e[]) 378454e267fSLisandro Dalcin { 379c73cfb54SMatthew G. Knepley PetscInt dim; 380454e267fSLisandro Dalcin PetscErrorCode ierr; 38199c57625SBarry Smith DM_DA *dd = (DM_DA*)dm->data; 382f951a8fcSStefano Zampini PetscBool isda; 3835fd66863SKarl Rupp 384454e267fSLisandro Dalcin PetscFunctionBegin; 385a9a02de4SBarry Smith PetscValidHeaderSpecificType(dm,DM_CLASSID,1,DMDA); 386f951a8fcSStefano Zampini PetscValidIntPointer(nel,2); 387f951a8fcSStefano Zampini PetscValidIntPointer(nen,3); 388f951a8fcSStefano Zampini PetscValidPointer(e,4); 389f951a8fcSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)dm,DMDA,&isda);CHKERRQ(ierr); 390*7a8be351SBarry Smith PetscCheck(isda,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for DM type %s",((PetscObject)dm)->type_name); 3912c71b3e2SJacob Faibussowitsch PetscCheckFalse(dd->stencil_type == DMDA_STENCIL_STAR,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"DMDAGetElements() requires you use a stencil type of DMDA_STENCIL_BOX"); 392c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 393cfbed8a3SStefano Zampini if (dd->e) { 394cfbed8a3SStefano Zampini *nel = dd->ne; 395cfbed8a3SStefano Zampini *nen = dd->nen; 396cfbed8a3SStefano Zampini *e = dd->e; 397cfbed8a3SStefano Zampini PetscFunctionReturn(0); 398cfbed8a3SStefano Zampini } 399c73cfb54SMatthew G. Knepley if (dim==-1) { 4000298fd71SBarry Smith *nel = 0; *nen = 0; *e = NULL; 401c73cfb54SMatthew G. Knepley } else if (dim==1) { 4022dde6fd4SLisandro Dalcin ierr = DMDAGetElements_1D(dm,nel,nen,e);CHKERRQ(ierr); 403c73cfb54SMatthew G. Knepley } else if (dim==2) { 4042dde6fd4SLisandro Dalcin ierr = DMDAGetElements_2D(dm,nel,nen,e);CHKERRQ(ierr); 405c73cfb54SMatthew G. Knepley } else if (dim==3) { 4062dde6fd4SLisandro Dalcin ierr = DMDAGetElements_3D(dm,nel,nen,e);CHKERRQ(ierr); 40798921bdaSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D",dim); 408454e267fSLisandro Dalcin PetscFunctionReturn(0); 409454e267fSLisandro Dalcin } 4102dde6fd4SLisandro Dalcin 411f951a8fcSStefano Zampini /*@ 412d4a6ed37SStefano Zampini DMDAGetSubdomainCornersIS - Gets an index set containing the corner indices (in local coordinates) 413f951a8fcSStefano Zampini of the non-overlapping decomposition identified by DMDAGetElements 414f951a8fcSStefano Zampini 415f951a8fcSStefano Zampini Not Collective 416f951a8fcSStefano Zampini 417f951a8fcSStefano Zampini Input Parameter: 418f951a8fcSStefano Zampini . dm - the DM object 419f951a8fcSStefano Zampini 420f951a8fcSStefano Zampini Output Parameters: 421f951a8fcSStefano Zampini . is - the index set 422f951a8fcSStefano Zampini 423f951a8fcSStefano Zampini Level: intermediate 424f951a8fcSStefano Zampini 42595452b02SPatrick Sanan Notes: 42695452b02SPatrick Sanan Call DMDARestoreSubdomainCornersIS() once you have finished accessing the index set. 427f951a8fcSStefano Zampini 428f951a8fcSStefano Zampini .seealso: DMDAElementType, DMDASetElementType(), DMDAGetElements(), DMDARestoreElementsCornersIS() 429f951a8fcSStefano Zampini @*/ 430d4a6ed37SStefano Zampini PetscErrorCode DMDAGetSubdomainCornersIS(DM dm,IS *is) 431f951a8fcSStefano Zampini { 432f951a8fcSStefano Zampini PetscErrorCode ierr; 433f951a8fcSStefano Zampini DM_DA *dd = (DM_DA*)dm->data; 434f951a8fcSStefano Zampini PetscBool isda; 435f951a8fcSStefano Zampini 436f951a8fcSStefano Zampini PetscFunctionBegin; 437a9a02de4SBarry Smith PetscValidHeaderSpecificType(dm,DM_CLASSID,1,DMDA); 438f951a8fcSStefano Zampini PetscValidPointer(is,2); 439f951a8fcSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)dm,DMDA,&isda);CHKERRQ(ierr); 440*7a8be351SBarry Smith PetscCheck(isda,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for DM type %s",((PetscObject)dm)->type_name); 4412c71b3e2SJacob Faibussowitsch PetscCheckFalse(dd->stencil_type == DMDA_STENCIL_STAR,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"DMDAGetElement() requires you use a stencil type of DMDA_STENCIL_BOX"); 442f951a8fcSStefano Zampini if (!dd->ecorners) { /* compute elements if not yet done */ 443f951a8fcSStefano Zampini const PetscInt *e; 444f951a8fcSStefano Zampini PetscInt nel,nen; 445f951a8fcSStefano Zampini 446f951a8fcSStefano Zampini ierr = DMDAGetElements(dm,&nel,&nen,&e);CHKERRQ(ierr); 447f951a8fcSStefano Zampini ierr = DMDARestoreElements(dm,&nel,&nen,&e);CHKERRQ(ierr); 448f951a8fcSStefano Zampini } 449f951a8fcSStefano Zampini *is = dd->ecorners; 450f951a8fcSStefano Zampini PetscFunctionReturn(0); 451f951a8fcSStefano Zampini } 452f951a8fcSStefano Zampini 4532dde6fd4SLisandro Dalcin /*@C 454009f0576SBarry Smith DMDARestoreElements - Restores the array obtained with DMDAGetElements() 4552dde6fd4SLisandro Dalcin 4562dde6fd4SLisandro Dalcin Not Collective 4572dde6fd4SLisandro Dalcin 458d8d19677SJose E. Roman Input Parameters: 4592dde6fd4SLisandro Dalcin + dm - the DM object 4602dde6fd4SLisandro Dalcin . nel - number of local elements 4612dde6fd4SLisandro Dalcin . nen - number of element nodes 462c2cd2aa3SJed Brown - e - the local indices of the elements' vertices 4632dde6fd4SLisandro Dalcin 4642dde6fd4SLisandro Dalcin Level: intermediate 4652dde6fd4SLisandro Dalcin 466009f0576SBarry Smith Note: You should not access these values after you have called this routine. 467009f0576SBarry Smith 468c67aec63SBarry Smith This restore signals the DMDA object that you no longer need access to the array information. 469009f0576SBarry Smith 470f5f57ec0SBarry Smith Not supported in Fortran 471f5f57ec0SBarry Smith 4722dde6fd4SLisandro Dalcin .seealso: DMDAElementType, DMDASetElementType(), DMDAGetElements() 4732dde6fd4SLisandro Dalcin @*/ 4742dde6fd4SLisandro Dalcin PetscErrorCode DMDARestoreElements(DM dm,PetscInt *nel,PetscInt *nen,const PetscInt *e[]) 4752dde6fd4SLisandro Dalcin { 4762dde6fd4SLisandro Dalcin PetscFunctionBegin; 477a9a02de4SBarry Smith PetscValidHeaderSpecificType(dm,DM_CLASSID,1,DMDA); 4782dde6fd4SLisandro Dalcin PetscValidIntPointer(nel,2); 4792dde6fd4SLisandro Dalcin PetscValidIntPointer(nen,3); 4802dde6fd4SLisandro Dalcin PetscValidPointer(e,4); 4819cc8bc54SJed Brown *nel = 0; 4829cc8bc54SJed Brown *nen = -1; 4839cc8bc54SJed Brown *e = NULL; 4842dde6fd4SLisandro Dalcin PetscFunctionReturn(0); 4852dde6fd4SLisandro Dalcin } 486f951a8fcSStefano Zampini 487f951a8fcSStefano Zampini /*@ 488d4a6ed37SStefano Zampini DMDARestoreSubdomainCornersIS - Restores the IS obtained with DMDAGetSubdomainCornersIS() 489f951a8fcSStefano Zampini 490f951a8fcSStefano Zampini Not Collective 491f951a8fcSStefano Zampini 492d8d19677SJose E. Roman Input Parameters: 493f951a8fcSStefano Zampini + dm - the DM object 494f951a8fcSStefano Zampini - is - the index set 495f951a8fcSStefano Zampini 496f951a8fcSStefano Zampini Level: intermediate 497f951a8fcSStefano Zampini 498f951a8fcSStefano Zampini Note: 499f951a8fcSStefano Zampini 500d4a6ed37SStefano Zampini .seealso: DMDAElementType, DMDASetElementType(), DMDAGetSubdomainCornersIS() 501f951a8fcSStefano Zampini @*/ 502d4a6ed37SStefano Zampini PetscErrorCode DMDARestoreSubdomainCornersIS(DM dm,IS *is) 503f951a8fcSStefano Zampini { 504f951a8fcSStefano Zampini PetscFunctionBegin; 505a9a02de4SBarry Smith PetscValidHeaderSpecificType(dm,DM_CLASSID,1,DMDA); 506f951a8fcSStefano Zampini PetscValidHeaderSpecific(*is,IS_CLASSID,2); 507f951a8fcSStefano Zampini *is = NULL; 508f951a8fcSStefano Zampini PetscFunctionReturn(0); 509f951a8fcSStefano Zampini } 510