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 DM_DA *da = (DM_DA*)dm->data; 7454e267fSLisandro Dalcin PetscInt i,xs,xe,Xs,Xe; 8454e267fSLisandro Dalcin PetscInt cnt=0; 95fd66863SKarl Rupp 10454e267fSLisandro Dalcin PetscFunctionBegin; 11454e267fSLisandro Dalcin if (!da->e) { 12f951a8fcSStefano Zampini PetscInt corners[2]; 13f951a8fcSStefano Zampini 147a8be351SBarry Smith PetscCheck(da->s,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Cannot get elements for DMDA with zero stencil width"); 159566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(dm,&xs,NULL,NULL,&xe,NULL,NULL)); 169566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(dm,&Xs,NULL,NULL,&Xe,NULL,NULL)); 17454e267fSLisandro Dalcin xe += xs; Xe += Xs; if (xs != Xs) xs -= 1; 18454e267fSLisandro Dalcin da->ne = 1*(xe - xs - 1); 199566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(1 + 2*da->ne,&da->e)); 20454e267fSLisandro Dalcin for (i=xs; i<xe-1; i++) { 21454e267fSLisandro Dalcin da->e[cnt++] = (i-Xs); 22454e267fSLisandro Dalcin da->e[cnt++] = (i-Xs+1); 23454e267fSLisandro Dalcin } 24cfbed8a3SStefano Zampini da->nen = 2; 25cfbed8a3SStefano Zampini 26f951a8fcSStefano Zampini corners[0] = (xs -Xs); 27f951a8fcSStefano Zampini corners[1] = (xe-1-Xs); 289566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF,2,corners,PETSC_COPY_VALUES,&da->ecorners)); 29454e267fSLisandro Dalcin } 30454e267fSLisandro Dalcin *nel = da->ne; 31cfbed8a3SStefano Zampini *nen = da->nen; 32454e267fSLisandro Dalcin *e = da->e; 33454e267fSLisandro Dalcin PetscFunctionReturn(0); 34454e267fSLisandro Dalcin } 35454e267fSLisandro Dalcin 362dde6fd4SLisandro Dalcin static PetscErrorCode DMDAGetElements_2D(DM dm,PetscInt *nel,PetscInt *nen,const PetscInt *e[]) 37454e267fSLisandro Dalcin { 38454e267fSLisandro Dalcin DM_DA *da = (DM_DA*)dm->data; 39454e267fSLisandro Dalcin PetscInt i,xs,xe,Xs,Xe; 40454e267fSLisandro Dalcin PetscInt j,ys,ye,Ys,Ye; 41cfbed8a3SStefano Zampini PetscInt cnt=0, cell[4], ns=2; 42454e267fSLisandro Dalcin PetscInt c, split[] = {0,1,3, 43454e267fSLisandro Dalcin 2,3,1}; 445fd66863SKarl Rupp 45454e267fSLisandro Dalcin PetscFunctionBegin; 46454e267fSLisandro Dalcin if (!da->e) { 47cfbed8a3SStefano Zampini PetscInt corners[4],nn = 0; 48f951a8fcSStefano Zampini 497a8be351SBarry Smith PetscCheck(da->s,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Cannot get elements for DMDA with zero stencil width"); 50cfbed8a3SStefano Zampini 51cfbed8a3SStefano Zampini switch (da->elementtype) { 52cfbed8a3SStefano Zampini case DMDA_ELEMENT_Q1: 53cfbed8a3SStefano Zampini da->nen = 4; 54cfbed8a3SStefano Zampini break; 55cfbed8a3SStefano Zampini case DMDA_ELEMENT_P1: 56cfbed8a3SStefano Zampini da->nen = 3; 57cfbed8a3SStefano Zampini break; 58cfbed8a3SStefano Zampini default: 5998921bdaSJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Unknown element type %d",da->elementtype); 60cfbed8a3SStefano Zampini } 61cfbed8a3SStefano Zampini nn = da->nen; 62cfbed8a3SStefano Zampini 63b644c393SDave May if (da->elementtype == DMDA_ELEMENT_P1) {ns=2;} 64b644c393SDave May if (da->elementtype == DMDA_ELEMENT_Q1) {ns=1;} 659566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(dm,&xs,&ys,NULL,&xe,&ye,NULL)); 669566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(dm,&Xs,&Ys,NULL,&Xe,&Ye,NULL)); 67454e267fSLisandro Dalcin xe += xs; Xe += Xs; if (xs != Xs) xs -= 1; 68454e267fSLisandro Dalcin ye += ys; Ye += Ys; if (ys != Ys) ys -= 1; 69454e267fSLisandro Dalcin da->ne = ns*(xe - xs - 1)*(ye - ys - 1); 709566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(1 + nn*da->ne,&da->e)); 71454e267fSLisandro Dalcin for (j=ys; j<ye-1; j++) { 72454e267fSLisandro Dalcin for (i=xs; i<xe-1; i++) { 73454e267fSLisandro Dalcin cell[0] = (i-Xs) + (j-Ys)*(Xe-Xs); 74454e267fSLisandro Dalcin cell[1] = (i-Xs+1) + (j-Ys)*(Xe-Xs); 75454e267fSLisandro Dalcin cell[2] = (i-Xs+1) + (j-Ys+1)*(Xe-Xs); 76454e267fSLisandro Dalcin cell[3] = (i-Xs) + (j-Ys+1)*(Xe-Xs); 77454e267fSLisandro Dalcin if (da->elementtype == DMDA_ELEMENT_P1) { 788865f1eaSKarl Rupp for (c=0; c<ns*nn; c++) da->e[cnt++] = cell[split[c]]; 79454e267fSLisandro Dalcin } 80454e267fSLisandro Dalcin if (da->elementtype == DMDA_ELEMENT_Q1) { 818865f1eaSKarl Rupp for (c=0; c<ns*nn; c++) da->e[cnt++] = cell[c]; 82454e267fSLisandro Dalcin } 83454e267fSLisandro Dalcin } 84454e267fSLisandro Dalcin } 85cfbed8a3SStefano Zampini 86f951a8fcSStefano Zampini corners[0] = (xs -Xs) + (ys -Ys)*(Xe-Xs); 87f951a8fcSStefano Zampini corners[1] = (xe-1-Xs) + (ys -Ys)*(Xe-Xs); 88f951a8fcSStefano Zampini corners[2] = (xs -Xs) + (ye-1-Ys)*(Xe-Xs); 89f951a8fcSStefano Zampini corners[3] = (xe-1-Xs) + (ye-1-Ys)*(Xe-Xs); 909566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF,4,corners,PETSC_COPY_VALUES,&da->ecorners)); 91454e267fSLisandro Dalcin } 92454e267fSLisandro Dalcin *nel = da->ne; 93cfbed8a3SStefano Zampini *nen = da->nen; 94454e267fSLisandro Dalcin *e = da->e; 95454e267fSLisandro Dalcin PetscFunctionReturn(0); 96454e267fSLisandro Dalcin } 97454e267fSLisandro Dalcin 982dde6fd4SLisandro Dalcin static PetscErrorCode DMDAGetElements_3D(DM dm,PetscInt *nel,PetscInt *nen,const PetscInt *e[]) 99454e267fSLisandro Dalcin { 100454e267fSLisandro Dalcin DM_DA *da = (DM_DA*)dm->data; 101454e267fSLisandro Dalcin PetscInt i,xs,xe,Xs,Xe; 102454e267fSLisandro Dalcin PetscInt j,ys,ye,Ys,Ye; 103454e267fSLisandro Dalcin PetscInt k,zs,ze,Zs,Ze; 104cfbed8a3SStefano Zampini PetscInt cnt=0, cell[8], ns=6; 105454e267fSLisandro Dalcin PetscInt c, split[] = {0,1,3,7, 106454e267fSLisandro Dalcin 0,1,7,4, 107454e267fSLisandro Dalcin 1,2,3,7, 108454e267fSLisandro Dalcin 1,2,7,6, 109454e267fSLisandro Dalcin 1,4,5,7, 110454e267fSLisandro Dalcin 1,5,6,7}; 1115fd66863SKarl Rupp 112454e267fSLisandro Dalcin PetscFunctionBegin; 113454e267fSLisandro Dalcin if (!da->e) { 114cfbed8a3SStefano Zampini PetscInt corners[8],nn = 0; 115f951a8fcSStefano Zampini 1167a8be351SBarry Smith PetscCheck(da->s,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Cannot get elements for DMDA with zero stencil width"); 117cfbed8a3SStefano Zampini 118cfbed8a3SStefano Zampini switch (da->elementtype) { 119cfbed8a3SStefano Zampini case DMDA_ELEMENT_Q1: 120cfbed8a3SStefano Zampini da->nen = 8; 121cfbed8a3SStefano Zampini break; 122cfbed8a3SStefano Zampini case DMDA_ELEMENT_P1: 123cfbed8a3SStefano Zampini da->nen = 4; 124cfbed8a3SStefano Zampini break; 125cfbed8a3SStefano Zampini default: 12698921bdaSJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Unknown element type %d",da->elementtype); 127cfbed8a3SStefano Zampini } 128cfbed8a3SStefano Zampini nn = da->nen; 129cfbed8a3SStefano Zampini 130b644c393SDave May if (da->elementtype == DMDA_ELEMENT_P1) {ns=6;} 131b644c393SDave May if (da->elementtype == DMDA_ELEMENT_Q1) {ns=1;} 1329566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(dm,&xs,&ys,&zs,&xe,&ye,&ze)); 1339566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(dm,&Xs,&Ys,&Zs,&Xe,&Ye,&Ze)); 134454e267fSLisandro Dalcin xe += xs; Xe += Xs; if (xs != Xs) xs -= 1; 135454e267fSLisandro Dalcin ye += ys; Ye += Ys; if (ys != Ys) ys -= 1; 136454e267fSLisandro Dalcin ze += zs; Ze += Zs; if (zs != Zs) zs -= 1; 137454e267fSLisandro Dalcin da->ne = ns*(xe - xs - 1)*(ye - ys - 1)*(ze - zs - 1); 1389566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(1 + nn*da->ne,&da->e)); 139454e267fSLisandro Dalcin for (k=zs; k<ze-1; k++) { 140454e267fSLisandro Dalcin for (j=ys; j<ye-1; j++) { 141454e267fSLisandro Dalcin for (i=xs; i<xe-1; i++) { 142454e267fSLisandro Dalcin cell[0] = (i -Xs) + (j -Ys)*(Xe-Xs) + (k -Zs)*(Xe-Xs)*(Ye-Ys); 1432da392ccSBarry Smith cell[1] = (i+1-Xs) + (j -Ys)*(Xe-Xs) + (k -Zs)*(Xe-Xs)*(Ye-Ys); 1442da392ccSBarry Smith cell[2] = (i+1-Xs) + (j+1-Ys)*(Xe-Xs) + (k -Zs)*(Xe-Xs)*(Ye-Ys); 1452da392ccSBarry Smith cell[3] = (i -Xs) + (j+1-Ys)*(Xe-Xs) + (k -Zs)*(Xe-Xs)*(Ye-Ys); 1462da392ccSBarry Smith cell[4] = (i -Xs) + (j -Ys)*(Xe-Xs) + (k+1-Zs)*(Xe-Xs)*(Ye-Ys); 1472da392ccSBarry Smith cell[5] = (i+1-Xs) + (j -Ys)*(Xe-Xs) + (k+1-Zs)*(Xe-Xs)*(Ye-Ys); 1482da392ccSBarry Smith cell[6] = (i+1-Xs) + (j+1-Ys)*(Xe-Xs) + (k+1-Zs)*(Xe-Xs)*(Ye-Ys); 1492da392ccSBarry Smith cell[7] = (i -Xs) + (j+1-Ys)*(Xe-Xs) + (k+1-Zs)*(Xe-Xs)*(Ye-Ys); 150454e267fSLisandro Dalcin if (da->elementtype == DMDA_ELEMENT_P1) { 1518865f1eaSKarl Rupp for (c=0; c<ns*nn; c++) da->e[cnt++] = cell[split[c]]; 152454e267fSLisandro Dalcin } 153454e267fSLisandro Dalcin if (da->elementtype == DMDA_ELEMENT_Q1) { 1548865f1eaSKarl Rupp for (c=0; c<ns*nn; c++) da->e[cnt++] = cell[c]; 155454e267fSLisandro Dalcin } 156454e267fSLisandro Dalcin } 157454e267fSLisandro Dalcin } 158454e267fSLisandro Dalcin } 159cfbed8a3SStefano Zampini 160f951a8fcSStefano Zampini corners[0] = (xs -Xs) + (ys -Ys)*(Xe-Xs) + (zs- Zs)*(Xe-Xs)*(Ye-Ys); 161f951a8fcSStefano Zampini corners[1] = (xe-1-Xs) + (ys -Ys)*(Xe-Xs) + (zs- Zs)*(Xe-Xs)*(Ye-Ys); 162f951a8fcSStefano Zampini corners[2] = (xs -Xs) + (ye-1-Ys)*(Xe-Xs) + (zs- Zs)*(Xe-Xs)*(Ye-Ys); 163f951a8fcSStefano Zampini corners[3] = (xe-1-Xs) + (ye-1-Ys)*(Xe-Xs) + (zs- Zs)*(Xe-Xs)*(Ye-Ys); 164f951a8fcSStefano Zampini corners[4] = (xs -Xs) + (ys -Ys)*(Xe-Xs) + (ze-1-Zs)*(Xe-Xs)*(Ye-Ys); 165f951a8fcSStefano Zampini corners[5] = (xe-1-Xs) + (ys -Ys)*(Xe-Xs) + (ze-1-Zs)*(Xe-Xs)*(Ye-Ys); 166f951a8fcSStefano Zampini corners[6] = (xs -Xs) + (ye-1-Ys)*(Xe-Xs) + (ze-1-Zs)*(Xe-Xs)*(Ye-Ys); 167f951a8fcSStefano Zampini corners[7] = (xe-1-Xs) + (ye-1-Ys)*(Xe-Xs) + (ze-1-Zs)*(Xe-Xs)*(Ye-Ys); 1689566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF,8,corners,PETSC_COPY_VALUES,&da->ecorners)); 169454e267fSLisandro Dalcin } 170454e267fSLisandro Dalcin *nel = da->ne; 171cfbed8a3SStefano Zampini *nen = da->nen; 172454e267fSLisandro Dalcin *e = da->e; 173454e267fSLisandro Dalcin PetscFunctionReturn(0); 174454e267fSLisandro Dalcin } 175454e267fSLisandro Dalcin 176f5f57ec0SBarry Smith /*@ 177d4a6ed37SStefano Zampini DMDAGetElementsCorners - Returns the global (x,y,z) indices of the lower left 178d4a6ed37SStefano Zampini corner of the non-overlapping decomposition identified by DMDAGetElements() 179d4a6ed37SStefano Zampini 180d4a6ed37SStefano Zampini Not Collective 181d4a6ed37SStefano Zampini 182d4a6ed37SStefano Zampini Input Parameter: 183d4a6ed37SStefano Zampini . da - the DM object 184d4a6ed37SStefano Zampini 185d4a6ed37SStefano Zampini Output Parameters: 186d4a6ed37SStefano Zampini + gx - the x index 187d4a6ed37SStefano Zampini . gy - the y index 188d4a6ed37SStefano Zampini - gz - the z index 189d4a6ed37SStefano Zampini 190d4a6ed37SStefano Zampini Level: intermediate 191d4a6ed37SStefano Zampini 192d4a6ed37SStefano Zampini Notes: 193d4a6ed37SStefano Zampini 194d4a6ed37SStefano Zampini .seealso: DMDAElementType, DMDASetElementType(), DMDAGetElements() 195d4a6ed37SStefano Zampini @*/ 196d4a6ed37SStefano Zampini PetscErrorCode DMDAGetElementsCorners(DM da, PetscInt *gx, PetscInt *gy, PetscInt *gz) 197d4a6ed37SStefano Zampini { 198d4a6ed37SStefano Zampini PetscInt xs,Xs; 199d4a6ed37SStefano Zampini PetscInt ys,Ys; 200d4a6ed37SStefano Zampini PetscInt zs,Zs; 201d4a6ed37SStefano Zampini PetscBool isda; 202d4a6ed37SStefano Zampini 203d4a6ed37SStefano Zampini PetscFunctionBegin; 204a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 205d4a6ed37SStefano Zampini if (gx) PetscValidIntPointer(gx,2); 206d4a6ed37SStefano Zampini if (gy) PetscValidIntPointer(gy,3); 207d4a6ed37SStefano Zampini if (gz) PetscValidIntPointer(gz,4); 2089566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)da,DMDA,&isda)); 2097a8be351SBarry Smith PetscCheck(isda,PetscObjectComm((PetscObject)da),PETSC_ERR_USER,"Not for DM type %s",((PetscObject)da)->type_name); 2109566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,&ys,&zs,NULL,NULL,NULL)); 2119566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da,&Xs,&Ys,&Zs,NULL,NULL,NULL)); 212d4a6ed37SStefano Zampini if (xs != Xs) xs -= 1; 213d4a6ed37SStefano Zampini if (ys != Ys) ys -= 1; 214d4a6ed37SStefano Zampini if (zs != Zs) zs -= 1; 215d4a6ed37SStefano Zampini if (gx) *gx = xs; 216d4a6ed37SStefano Zampini if (gy) *gy = ys; 217d4a6ed37SStefano Zampini if (gz) *gz = zs; 218d4a6ed37SStefano Zampini PetscFunctionReturn(0); 219d4a6ed37SStefano Zampini } 220d4a6ed37SStefano Zampini 221d4a6ed37SStefano Zampini /*@ 222d4a6ed37SStefano Zampini DMDAGetElementsSizes - Gets the local number of elements per direction for the non-overlapping decomposition identified by DMDAGetElements() 223d4a6ed37SStefano Zampini 224d4a6ed37SStefano Zampini Not Collective 225d4a6ed37SStefano Zampini 226d4a6ed37SStefano Zampini Input Parameter: 227d4a6ed37SStefano Zampini . da - the DM object 228d4a6ed37SStefano Zampini 229d4a6ed37SStefano Zampini Output Parameters: 230d4a6ed37SStefano Zampini + mx - number of local elements in x-direction 231d4a6ed37SStefano Zampini . my - number of local elements in y-direction 232d4a6ed37SStefano Zampini - mz - number of local elements in z-direction 233d4a6ed37SStefano Zampini 234d4a6ed37SStefano Zampini Level: intermediate 235d4a6ed37SStefano Zampini 23695452b02SPatrick Sanan Notes: 23795452b02SPatrick Sanan It returns the same number of elements, irrespective of the DMDAElementType 238d4a6ed37SStefano Zampini 239d4a6ed37SStefano Zampini .seealso: DMDAElementType, DMDASetElementType(), DMDAGetElements 240d4a6ed37SStefano Zampini @*/ 241d4a6ed37SStefano Zampini PetscErrorCode DMDAGetElementsSizes(DM da, PetscInt *mx, PetscInt *my, PetscInt *mz) 242d4a6ed37SStefano Zampini { 243d4a6ed37SStefano Zampini PetscInt xs,xe,Xs; 244d4a6ed37SStefano Zampini PetscInt ys,ye,Ys; 245d4a6ed37SStefano Zampini PetscInt zs,ze,Zs; 246d4a6ed37SStefano Zampini PetscInt dim; 247d4a6ed37SStefano Zampini PetscBool isda; 248d4a6ed37SStefano Zampini 249d4a6ed37SStefano Zampini PetscFunctionBegin; 250a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 251d4a6ed37SStefano Zampini if (mx) PetscValidIntPointer(mx,2); 252d4a6ed37SStefano Zampini if (my) PetscValidIntPointer(my,3); 253d4a6ed37SStefano Zampini if (mz) PetscValidIntPointer(mz,4); 2549566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)da,DMDA,&isda)); 2557a8be351SBarry Smith PetscCheck(isda,PetscObjectComm((PetscObject)da),PETSC_ERR_USER,"Not for DM type %s",((PetscObject)da)->type_name); 2569566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,&ys,&zs,&xe,&ye,&ze)); 2579566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da,&Xs,&Ys,&Zs,NULL,NULL,NULL)); 258d4a6ed37SStefano Zampini xe += xs; if (xs != Xs) xs -= 1; 259d4a6ed37SStefano Zampini ye += ys; if (ys != Ys) ys -= 1; 260d4a6ed37SStefano Zampini ze += zs; if (zs != Zs) zs -= 1; 261d4a6ed37SStefano Zampini if (mx) *mx = 0; 262d4a6ed37SStefano Zampini if (my) *my = 0; 263d4a6ed37SStefano Zampini if (mz) *mz = 0; 2649566063dSJacob Faibussowitsch PetscCall(DMGetDimension(da,&dim)); 265d4a6ed37SStefano Zampini switch (dim) { 266d4a6ed37SStefano Zampini case 3: 267d4a6ed37SStefano Zampini if (mz) *mz = ze - zs - 1; 268d4a6ed37SStefano Zampini case 2: 269d4a6ed37SStefano Zampini if (my) *my = ye - ys - 1; 270d4a6ed37SStefano Zampini case 1: 271d4a6ed37SStefano Zampini if (mx) *mx = xe - xs - 1; 272d4a6ed37SStefano Zampini break; 273d4a6ed37SStefano Zampini } 274d4a6ed37SStefano Zampini PetscFunctionReturn(0); 275d4a6ed37SStefano Zampini } 276d4a6ed37SStefano Zampini 277d4a6ed37SStefano Zampini /*@ 2782dde6fd4SLisandro Dalcin DMDASetElementType - Sets the element type to be returned by DMDAGetElements() 2792dde6fd4SLisandro Dalcin 2802dde6fd4SLisandro Dalcin Not Collective 2812dde6fd4SLisandro Dalcin 2822dde6fd4SLisandro Dalcin Input Parameter: 2832dde6fd4SLisandro Dalcin . da - the DMDA object 2842dde6fd4SLisandro Dalcin 2852dde6fd4SLisandro Dalcin Output Parameters: 2866420f107SDave May . etype - the element type, currently either DMDA_ELEMENT_P1 or DMDA_ELEMENT_Q1 2872dde6fd4SLisandro Dalcin 2882dde6fd4SLisandro Dalcin Level: intermediate 2892dde6fd4SLisandro Dalcin 2902dde6fd4SLisandro Dalcin .seealso: DMDAElementType, DMDAGetElementType(), DMDAGetElements(), DMDARestoreElements() 2912dde6fd4SLisandro Dalcin @*/ 2922dde6fd4SLisandro Dalcin PetscErrorCode DMDASetElementType(DM da, DMDAElementType etype) 2932dde6fd4SLisandro Dalcin { 2942dde6fd4SLisandro Dalcin DM_DA *dd = (DM_DA*)da->data; 295f951a8fcSStefano Zampini PetscBool isda; 2962dde6fd4SLisandro Dalcin 2972dde6fd4SLisandro Dalcin PetscFunctionBegin; 298a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 2992dde6fd4SLisandro Dalcin PetscValidLogicalCollectiveEnum(da,etype,2); 3009566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)da,DMDA,&isda)); 301f951a8fcSStefano Zampini if (!isda) PetscFunctionReturn(0); 3022dde6fd4SLisandro Dalcin if (dd->elementtype != etype) { 3039566063dSJacob Faibussowitsch PetscCall(PetscFree(dd->e)); 3049566063dSJacob Faibussowitsch PetscCall(ISDestroy(&dd->ecorners)); 3058865f1eaSKarl Rupp 3062dde6fd4SLisandro Dalcin dd->elementtype = etype; 3072dde6fd4SLisandro Dalcin dd->ne = 0; 308cfbed8a3SStefano Zampini dd->nen = 0; 3090298fd71SBarry Smith dd->e = NULL; 3102dde6fd4SLisandro Dalcin } 3112dde6fd4SLisandro Dalcin PetscFunctionReturn(0); 3122dde6fd4SLisandro Dalcin } 3132dde6fd4SLisandro Dalcin 314f5f57ec0SBarry Smith /*@ 3152dde6fd4SLisandro Dalcin DMDAGetElementType - Gets the element type to be returned by DMDAGetElements() 3162dde6fd4SLisandro Dalcin 3172dde6fd4SLisandro Dalcin Not Collective 3182dde6fd4SLisandro Dalcin 3192dde6fd4SLisandro Dalcin Input Parameter: 3202dde6fd4SLisandro Dalcin . da - the DMDA object 3212dde6fd4SLisandro Dalcin 3222dde6fd4SLisandro Dalcin Output Parameters: 3236420f107SDave May . etype - the element type, currently either DMDA_ELEMENT_P1 or DMDA_ELEMENT_Q1 3242dde6fd4SLisandro Dalcin 3252dde6fd4SLisandro Dalcin Level: intermediate 3262dde6fd4SLisandro Dalcin 3272dde6fd4SLisandro Dalcin .seealso: DMDAElementType, DMDASetElementType(), DMDAGetElements(), DMDARestoreElements() 3282dde6fd4SLisandro Dalcin @*/ 3292dde6fd4SLisandro Dalcin PetscErrorCode DMDAGetElementType(DM da, DMDAElementType *etype) 3302dde6fd4SLisandro Dalcin { 3312dde6fd4SLisandro Dalcin DM_DA *dd = (DM_DA*)da->data; 332f951a8fcSStefano Zampini PetscBool isda; 3335fd66863SKarl Rupp 3342dde6fd4SLisandro Dalcin PetscFunctionBegin; 335a9a02de4SBarry Smith PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA); 3362dde6fd4SLisandro Dalcin PetscValidPointer(etype,2); 3379566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)da,DMDA,&isda)); 3387a8be351SBarry Smith PetscCheck(isda,PetscObjectComm((PetscObject)da),PETSC_ERR_USER,"Not for DM type %s",((PetscObject)da)->type_name); 3392dde6fd4SLisandro Dalcin *etype = dd->elementtype; 3402dde6fd4SLisandro Dalcin PetscFunctionReturn(0); 3412dde6fd4SLisandro Dalcin } 3422dde6fd4SLisandro Dalcin 3432dde6fd4SLisandro Dalcin /*@C 3442dde6fd4SLisandro Dalcin DMDAGetElements - Gets an array containing the indices (in local coordinates) 3452dde6fd4SLisandro Dalcin of all the local elements 3462dde6fd4SLisandro Dalcin 3472dde6fd4SLisandro Dalcin Not Collective 3482dde6fd4SLisandro Dalcin 3492dde6fd4SLisandro Dalcin Input Parameter: 3502dde6fd4SLisandro Dalcin . dm - the DM object 3512dde6fd4SLisandro Dalcin 3522dde6fd4SLisandro Dalcin Output Parameters: 3532dde6fd4SLisandro Dalcin + nel - number of local elements 3542dde6fd4SLisandro Dalcin . nen - number of element nodes 355c2cd2aa3SJed Brown - e - the local indices of the elements' vertices 3562dde6fd4SLisandro Dalcin 3572dde6fd4SLisandro Dalcin Level: intermediate 3582dde6fd4SLisandro Dalcin 35993386343SJed Brown Notes: 36093386343SJed Brown Call DMDARestoreElements() once you have finished accessing the elements. 36193386343SJed Brown 36245950774SSatish Balay Each process uniquely owns a subset of the elements. That is no element is owned by two or more processes. 363009f0576SBarry Smith 364009f0576SBarry 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. 365009f0576SBarry Smith 366f5f57ec0SBarry Smith Not supported in Fortran 367f5f57ec0SBarry Smith 36893386343SJed Brown .seealso: DMDAElementType, DMDASetElementType(), VecSetValuesLocal(), MatSetValuesLocal(), DMGlobalToLocalBegin(), DMLocalToGlobalBegin() 3692dde6fd4SLisandro Dalcin @*/ 3702dde6fd4SLisandro Dalcin PetscErrorCode DMDAGetElements(DM dm,PetscInt *nel,PetscInt *nen,const PetscInt *e[]) 371454e267fSLisandro Dalcin { 372c73cfb54SMatthew G. Knepley PetscInt dim; 37399c57625SBarry Smith DM_DA *dd = (DM_DA*)dm->data; 374f951a8fcSStefano Zampini PetscBool isda; 3755fd66863SKarl Rupp 376454e267fSLisandro Dalcin PetscFunctionBegin; 377a9a02de4SBarry Smith PetscValidHeaderSpecificType(dm,DM_CLASSID,1,DMDA); 378f951a8fcSStefano Zampini PetscValidIntPointer(nel,2); 379f951a8fcSStefano Zampini PetscValidIntPointer(nen,3); 380f951a8fcSStefano Zampini PetscValidPointer(e,4); 3819566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm,DMDA,&isda)); 3827a8be351SBarry Smith PetscCheck(isda,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for DM type %s",((PetscObject)dm)->type_name); 383*08401ef6SPierre Jolivet PetscCheck(dd->stencil_type != DMDA_STENCIL_STAR,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"DMDAGetElements() requires you use a stencil type of DMDA_STENCIL_BOX"); 3849566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 385cfbed8a3SStefano Zampini if (dd->e) { 386cfbed8a3SStefano Zampini *nel = dd->ne; 387cfbed8a3SStefano Zampini *nen = dd->nen; 388cfbed8a3SStefano Zampini *e = dd->e; 389cfbed8a3SStefano Zampini PetscFunctionReturn(0); 390cfbed8a3SStefano Zampini } 391c73cfb54SMatthew G. Knepley if (dim==-1) { 3920298fd71SBarry Smith *nel = 0; *nen = 0; *e = NULL; 393c73cfb54SMatthew G. Knepley } else if (dim==1) { 3949566063dSJacob Faibussowitsch PetscCall(DMDAGetElements_1D(dm,nel,nen,e)); 395c73cfb54SMatthew G. Knepley } else if (dim==2) { 3969566063dSJacob Faibussowitsch PetscCall(DMDAGetElements_2D(dm,nel,nen,e)); 397c73cfb54SMatthew G. Knepley } else if (dim==3) { 3989566063dSJacob Faibussowitsch PetscCall(DMDAGetElements_3D(dm,nel,nen,e)); 39998921bdaSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D",dim); 400454e267fSLisandro Dalcin PetscFunctionReturn(0); 401454e267fSLisandro Dalcin } 4022dde6fd4SLisandro Dalcin 403f951a8fcSStefano Zampini /*@ 404d4a6ed37SStefano Zampini DMDAGetSubdomainCornersIS - Gets an index set containing the corner indices (in local coordinates) 405f951a8fcSStefano Zampini of the non-overlapping decomposition identified by DMDAGetElements 406f951a8fcSStefano Zampini 407f951a8fcSStefano Zampini Not Collective 408f951a8fcSStefano Zampini 409f951a8fcSStefano Zampini Input Parameter: 410f951a8fcSStefano Zampini . dm - the DM object 411f951a8fcSStefano Zampini 412f951a8fcSStefano Zampini Output Parameters: 413f951a8fcSStefano Zampini . is - the index set 414f951a8fcSStefano Zampini 415f951a8fcSStefano Zampini Level: intermediate 416f951a8fcSStefano Zampini 41795452b02SPatrick Sanan Notes: 41895452b02SPatrick Sanan Call DMDARestoreSubdomainCornersIS() once you have finished accessing the index set. 419f951a8fcSStefano Zampini 420f951a8fcSStefano Zampini .seealso: DMDAElementType, DMDASetElementType(), DMDAGetElements(), DMDARestoreElementsCornersIS() 421f951a8fcSStefano Zampini @*/ 422d4a6ed37SStefano Zampini PetscErrorCode DMDAGetSubdomainCornersIS(DM dm,IS *is) 423f951a8fcSStefano Zampini { 424f951a8fcSStefano Zampini DM_DA *dd = (DM_DA*)dm->data; 425f951a8fcSStefano Zampini PetscBool isda; 426f951a8fcSStefano Zampini 427f951a8fcSStefano Zampini PetscFunctionBegin; 428a9a02de4SBarry Smith PetscValidHeaderSpecificType(dm,DM_CLASSID,1,DMDA); 429f951a8fcSStefano Zampini PetscValidPointer(is,2); 4309566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm,DMDA,&isda)); 4317a8be351SBarry Smith PetscCheck(isda,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Not for DM type %s",((PetscObject)dm)->type_name); 432*08401ef6SPierre Jolivet PetscCheck(dd->stencil_type != DMDA_STENCIL_STAR,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"DMDAGetElement() requires you use a stencil type of DMDA_STENCIL_BOX"); 433f951a8fcSStefano Zampini if (!dd->ecorners) { /* compute elements if not yet done */ 434f951a8fcSStefano Zampini const PetscInt *e; 435f951a8fcSStefano Zampini PetscInt nel,nen; 436f951a8fcSStefano Zampini 4379566063dSJacob Faibussowitsch PetscCall(DMDAGetElements(dm,&nel,&nen,&e)); 4389566063dSJacob Faibussowitsch PetscCall(DMDARestoreElements(dm,&nel,&nen,&e)); 439f951a8fcSStefano Zampini } 440f951a8fcSStefano Zampini *is = dd->ecorners; 441f951a8fcSStefano Zampini PetscFunctionReturn(0); 442f951a8fcSStefano Zampini } 443f951a8fcSStefano Zampini 4442dde6fd4SLisandro Dalcin /*@C 445009f0576SBarry Smith DMDARestoreElements - Restores the array obtained with DMDAGetElements() 4462dde6fd4SLisandro Dalcin 4472dde6fd4SLisandro Dalcin Not Collective 4482dde6fd4SLisandro Dalcin 449d8d19677SJose E. Roman Input Parameters: 4502dde6fd4SLisandro Dalcin + dm - the DM object 4512dde6fd4SLisandro Dalcin . nel - number of local elements 4522dde6fd4SLisandro Dalcin . nen - number of element nodes 453c2cd2aa3SJed Brown - e - the local indices of the elements' vertices 4542dde6fd4SLisandro Dalcin 4552dde6fd4SLisandro Dalcin Level: intermediate 4562dde6fd4SLisandro Dalcin 457009f0576SBarry Smith Note: You should not access these values after you have called this routine. 458009f0576SBarry Smith 459c67aec63SBarry Smith This restore signals the DMDA object that you no longer need access to the array information. 460009f0576SBarry Smith 461f5f57ec0SBarry Smith Not supported in Fortran 462f5f57ec0SBarry Smith 4632dde6fd4SLisandro Dalcin .seealso: DMDAElementType, DMDASetElementType(), DMDAGetElements() 4642dde6fd4SLisandro Dalcin @*/ 4652dde6fd4SLisandro Dalcin PetscErrorCode DMDARestoreElements(DM dm,PetscInt *nel,PetscInt *nen,const PetscInt *e[]) 4662dde6fd4SLisandro Dalcin { 4672dde6fd4SLisandro Dalcin PetscFunctionBegin; 468a9a02de4SBarry Smith PetscValidHeaderSpecificType(dm,DM_CLASSID,1,DMDA); 4692dde6fd4SLisandro Dalcin PetscValidIntPointer(nel,2); 4702dde6fd4SLisandro Dalcin PetscValidIntPointer(nen,3); 4712dde6fd4SLisandro Dalcin PetscValidPointer(e,4); 4729cc8bc54SJed Brown *nel = 0; 4739cc8bc54SJed Brown *nen = -1; 4749cc8bc54SJed Brown *e = NULL; 4752dde6fd4SLisandro Dalcin PetscFunctionReturn(0); 4762dde6fd4SLisandro Dalcin } 477f951a8fcSStefano Zampini 478f951a8fcSStefano Zampini /*@ 479d4a6ed37SStefano Zampini DMDARestoreSubdomainCornersIS - Restores the IS obtained with DMDAGetSubdomainCornersIS() 480f951a8fcSStefano Zampini 481f951a8fcSStefano Zampini Not Collective 482f951a8fcSStefano Zampini 483d8d19677SJose E. Roman Input Parameters: 484f951a8fcSStefano Zampini + dm - the DM object 485f951a8fcSStefano Zampini - is - the index set 486f951a8fcSStefano Zampini 487f951a8fcSStefano Zampini Level: intermediate 488f951a8fcSStefano Zampini 489f951a8fcSStefano Zampini Note: 490f951a8fcSStefano Zampini 491d4a6ed37SStefano Zampini .seealso: DMDAElementType, DMDASetElementType(), DMDAGetSubdomainCornersIS() 492f951a8fcSStefano Zampini @*/ 493d4a6ed37SStefano Zampini PetscErrorCode DMDARestoreSubdomainCornersIS(DM dm,IS *is) 494f951a8fcSStefano Zampini { 495f951a8fcSStefano Zampini PetscFunctionBegin; 496a9a02de4SBarry Smith PetscValidHeaderSpecificType(dm,DM_CLASSID,1,DMDA); 497f951a8fcSStefano Zampini PetscValidHeaderSpecific(*is,IS_CLASSID,2); 498f951a8fcSStefano Zampini *is = NULL; 499f951a8fcSStefano Zampini PetscFunctionReturn(0); 500f951a8fcSStefano Zampini } 501