1454e267fSLisandro Dalcin 2af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h> /*I "petscdmda.h" I*/ 3454e267fSLisandro Dalcin 4454e267fSLisandro Dalcin #undef __FUNCT__ 52dde6fd4SLisandro Dalcin #define __FUNCT__ "DMDAGetElements_1D" 62dde6fd4SLisandro Dalcin static PetscErrorCode DMDAGetElements_1D(DM dm,PetscInt *nel,PetscInt *nen,const PetscInt *e[]) 7454e267fSLisandro Dalcin { 8454e267fSLisandro Dalcin PetscErrorCode ierr; 9454e267fSLisandro Dalcin DM_DA *da = (DM_DA*)dm->data; 10454e267fSLisandro Dalcin PetscInt i,xs,xe,Xs,Xe; 11454e267fSLisandro Dalcin PetscInt cnt=0; 125fd66863SKarl Rupp 13454e267fSLisandro Dalcin PetscFunctionBegin; 14454e267fSLisandro Dalcin if (!da->e) { 1598b2d131SBarry Smith if (!da->s) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Cannot get elements for DMDA with zero stencil width"); 16454e267fSLisandro Dalcin ierr = DMDAGetCorners(dm,&xs,0,0,&xe,0,0);CHKERRQ(ierr); 17454e267fSLisandro Dalcin ierr = DMDAGetGhostCorners(dm,&Xs,0,0,&Xe,0,0);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 } 25454e267fSLisandro Dalcin } 26454e267fSLisandro Dalcin *nel = da->ne; 27454e267fSLisandro Dalcin *nen = 2; 28454e267fSLisandro Dalcin *e = da->e; 29454e267fSLisandro Dalcin PetscFunctionReturn(0); 30454e267fSLisandro Dalcin } 31454e267fSLisandro Dalcin 32454e267fSLisandro Dalcin #undef __FUNCT__ 332dde6fd4SLisandro Dalcin #define __FUNCT__ "DMDAGetElements_2D" 342dde6fd4SLisandro Dalcin static PetscErrorCode DMDAGetElements_2D(DM dm,PetscInt *nel,PetscInt *nen,const PetscInt *e[]) 35454e267fSLisandro Dalcin { 36454e267fSLisandro Dalcin PetscErrorCode ierr; 37454e267fSLisandro Dalcin DM_DA *da = (DM_DA*)dm->data; 38454e267fSLisandro Dalcin PetscInt i,xs,xe,Xs,Xe; 39454e267fSLisandro Dalcin PetscInt j,ys,ye,Ys,Ye; 40454e267fSLisandro Dalcin PetscInt cnt=0, cell[4], ns=2, nn=3; 41454e267fSLisandro Dalcin PetscInt c, split[] = {0,1,3, 42454e267fSLisandro Dalcin 2,3,1}; 435fd66863SKarl Rupp 44454e267fSLisandro Dalcin PetscFunctionBegin; 45454e267fSLisandro Dalcin if (!da->e) { 4698b2d131SBarry Smith if (!da->s) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Cannot get elements for DMDA with zero stencil width"); 47454e267fSLisandro Dalcin if (da->elementtype == DMDA_ELEMENT_P1) {ns=2; nn=3;} 48454e267fSLisandro Dalcin if (da->elementtype == DMDA_ELEMENT_Q1) {ns=1; nn=4;} 49454e267fSLisandro Dalcin ierr = DMDAGetCorners(dm,&xs,&ys,0,&xe,&ye,0);CHKERRQ(ierr); 50454e267fSLisandro Dalcin ierr = DMDAGetGhostCorners(dm,&Xs,&Ys,0,&Xe,&Ye,0);CHKERRQ(ierr); 51454e267fSLisandro Dalcin xe += xs; Xe += Xs; if (xs != Xs) xs -= 1; 52454e267fSLisandro Dalcin ye += ys; Ye += Ys; if (ys != Ys) ys -= 1; 53454e267fSLisandro Dalcin da->ne = ns*(xe - xs - 1)*(ye - ys - 1); 54854ce69bSBarry Smith ierr = PetscMalloc1(1 + nn*da->ne,&da->e);CHKERRQ(ierr); 55454e267fSLisandro Dalcin for (j=ys; j<ye-1; j++) { 56454e267fSLisandro Dalcin for (i=xs; i<xe-1; i++) { 57454e267fSLisandro Dalcin cell[0] = (i-Xs ) + (j-Ys )*(Xe-Xs); 58454e267fSLisandro Dalcin cell[1] = (i-Xs+1) + (j-Ys )*(Xe-Xs); 59454e267fSLisandro Dalcin cell[2] = (i-Xs+1) + (j-Ys+1)*(Xe-Xs); 60454e267fSLisandro Dalcin cell[3] = (i-Xs ) + (j-Ys+1)*(Xe-Xs); 61454e267fSLisandro Dalcin if (da->elementtype == DMDA_ELEMENT_P1) { 628865f1eaSKarl Rupp for (c=0; c<ns*nn; c++) da->e[cnt++] = cell[split[c]]; 63454e267fSLisandro Dalcin } 64454e267fSLisandro Dalcin if (da->elementtype == DMDA_ELEMENT_Q1) { 658865f1eaSKarl Rupp for (c=0; c<ns*nn; c++) da->e[cnt++] = cell[c]; 66454e267fSLisandro Dalcin } 67454e267fSLisandro Dalcin } 68454e267fSLisandro Dalcin } 69454e267fSLisandro Dalcin } 70454e267fSLisandro Dalcin *nel = da->ne; 71454e267fSLisandro Dalcin *nen = nn; 72454e267fSLisandro Dalcin *e = da->e; 73454e267fSLisandro Dalcin PetscFunctionReturn(0); 74454e267fSLisandro Dalcin } 75454e267fSLisandro Dalcin 76454e267fSLisandro Dalcin #undef __FUNCT__ 772dde6fd4SLisandro Dalcin #define __FUNCT__ "DMDAGetElements_3D" 782dde6fd4SLisandro Dalcin static PetscErrorCode DMDAGetElements_3D(DM dm,PetscInt *nel,PetscInt *nen,const PetscInt *e[]) 79454e267fSLisandro Dalcin { 80454e267fSLisandro Dalcin PetscErrorCode ierr; 81454e267fSLisandro Dalcin DM_DA *da = (DM_DA*)dm->data; 82454e267fSLisandro Dalcin PetscInt i,xs,xe,Xs,Xe; 83454e267fSLisandro Dalcin PetscInt j,ys,ye,Ys,Ye; 84454e267fSLisandro Dalcin PetscInt k,zs,ze,Zs,Ze; 85454e267fSLisandro Dalcin PetscInt cnt=0, cell[8], ns=6, nn=4; 86454e267fSLisandro Dalcin PetscInt c, split[] = {0,1,3,7, 87454e267fSLisandro Dalcin 0,1,7,4, 88454e267fSLisandro Dalcin 1,2,3,7, 89454e267fSLisandro Dalcin 1,2,7,6, 90454e267fSLisandro Dalcin 1,4,5,7, 91454e267fSLisandro Dalcin 1,5,6,7}; 925fd66863SKarl Rupp 93454e267fSLisandro Dalcin PetscFunctionBegin; 94454e267fSLisandro Dalcin if (!da->e) { 9598b2d131SBarry Smith if (!da->s) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Cannot get elements for DMDA with zero stencil width"); 96454e267fSLisandro Dalcin if (da->elementtype == DMDA_ELEMENT_P1) {ns=6; nn=4;} 97454e267fSLisandro Dalcin if (da->elementtype == DMDA_ELEMENT_Q1) {ns=1; nn=8;} 98454e267fSLisandro Dalcin ierr = DMDAGetCorners(dm,&xs,&ys,&zs,&xe,&ye,&ze);CHKERRQ(ierr); 99454e267fSLisandro Dalcin ierr = DMDAGetGhostCorners(dm,&Xs,&Ys,&Zs,&Xe,&Ye,&Ze);CHKERRQ(ierr); 100454e267fSLisandro Dalcin xe += xs; Xe += Xs; if (xs != Xs) xs -= 1; 101454e267fSLisandro Dalcin ye += ys; Ye += Ys; if (ys != Ys) ys -= 1; 102454e267fSLisandro Dalcin ze += zs; Ze += Zs; if (zs != Zs) zs -= 1; 103454e267fSLisandro Dalcin da->ne = ns*(xe - xs - 1)*(ye - ys - 1)*(ze - zs - 1); 104854ce69bSBarry Smith ierr = PetscMalloc1(1 + nn*da->ne,&da->e);CHKERRQ(ierr); 105454e267fSLisandro Dalcin for (k=zs; k<ze-1; k++) { 106454e267fSLisandro Dalcin for (j=ys; j<ye-1; j++) { 107454e267fSLisandro Dalcin for (i=xs; i<xe-1; i++) { 108454e267fSLisandro Dalcin cell[0] = (i-Xs ) + (j-Ys )*(Xe-Xs) + (k-Zs )*(Xe-Xs)*(Ye-Ys); 109454e267fSLisandro Dalcin cell[1] = (i-Xs+1) + (j-Ys )*(Xe-Xs) + (k-Zs )*(Xe-Xs)*(Ye-Ys); 110454e267fSLisandro Dalcin cell[2] = (i-Xs+1) + (j-Ys+1)*(Xe-Xs) + (k-Zs )*(Xe-Xs)*(Ye-Ys); 111454e267fSLisandro Dalcin cell[3] = (i-Xs ) + (j-Ys+1)*(Xe-Xs) + (k-Zs )*(Xe-Xs)*(Ye-Ys); 112454e267fSLisandro Dalcin cell[4] = (i-Xs ) + (j-Ys )*(Xe-Xs) + (k-Zs+1)*(Xe-Xs)*(Ye-Ys); 113454e267fSLisandro Dalcin cell[5] = (i-Xs+1) + (j-Ys )*(Xe-Xs) + (k-Zs+1)*(Xe-Xs)*(Ye-Ys); 114454e267fSLisandro Dalcin cell[6] = (i-Xs+1) + (j-Ys+1)*(Xe-Xs) + (k-Zs+1)*(Xe-Xs)*(Ye-Ys); 115454e267fSLisandro Dalcin cell[7] = (i-Xs ) + (j-Ys+1)*(Xe-Xs) + (k-Zs+1)*(Xe-Xs)*(Ye-Ys); 116454e267fSLisandro Dalcin if (da->elementtype == DMDA_ELEMENT_P1) { 1178865f1eaSKarl Rupp for (c=0; c<ns*nn; c++) da->e[cnt++] = cell[split[c]]; 118454e267fSLisandro Dalcin } 119454e267fSLisandro Dalcin if (da->elementtype == DMDA_ELEMENT_Q1) { 1208865f1eaSKarl Rupp for (c=0; c<ns*nn; c++) da->e[cnt++] = cell[c]; 121454e267fSLisandro Dalcin } 122454e267fSLisandro Dalcin } 123454e267fSLisandro Dalcin } 124454e267fSLisandro Dalcin } 125454e267fSLisandro Dalcin } 126454e267fSLisandro Dalcin *nel = da->ne; 127454e267fSLisandro Dalcin *nen = nn; 128454e267fSLisandro Dalcin *e = da->e; 129454e267fSLisandro Dalcin PetscFunctionReturn(0); 130454e267fSLisandro Dalcin } 131454e267fSLisandro Dalcin 1322dde6fd4SLisandro Dalcin /*@C 1332dde6fd4SLisandro Dalcin DMDASetElementType - Sets the element type to be returned by DMDAGetElements() 1342dde6fd4SLisandro Dalcin 1352dde6fd4SLisandro Dalcin Not Collective 1362dde6fd4SLisandro Dalcin 1372dde6fd4SLisandro Dalcin Input Parameter: 1382dde6fd4SLisandro Dalcin . da - the DMDA object 1392dde6fd4SLisandro Dalcin 1402dde6fd4SLisandro Dalcin Output Parameters: 1412dde6fd4SLisandro Dalcin . etype - the element type, currently either DMDA_ELEMENT_P1 or ELEMENT_Q1 1422dde6fd4SLisandro Dalcin 1432dde6fd4SLisandro Dalcin Level: intermediate 1442dde6fd4SLisandro Dalcin 1452dde6fd4SLisandro Dalcin .seealso: DMDAElementType, DMDAGetElementType(), DMDAGetElements(), DMDARestoreElements() 1462dde6fd4SLisandro Dalcin @*/ 147454e267fSLisandro Dalcin #undef __FUNCT__ 1482dde6fd4SLisandro Dalcin #define __FUNCT__ "DMDASetElementType" 1492dde6fd4SLisandro Dalcin PetscErrorCode DMDASetElementType(DM da, DMDAElementType etype) 1502dde6fd4SLisandro Dalcin { 1512dde6fd4SLisandro Dalcin DM_DA *dd = (DM_DA*)da->data; 1522dde6fd4SLisandro Dalcin PetscErrorCode ierr; 1532dde6fd4SLisandro Dalcin 1542dde6fd4SLisandro Dalcin PetscFunctionBegin; 1552dde6fd4SLisandro Dalcin PetscValidHeaderSpecific(da,DM_CLASSID,1); 1562dde6fd4SLisandro Dalcin PetscValidLogicalCollectiveEnum(da,etype,2); 1572dde6fd4SLisandro Dalcin if (dd->elementtype != etype) { 1582dde6fd4SLisandro Dalcin ierr = PetscFree(dd->e);CHKERRQ(ierr); 1598865f1eaSKarl Rupp 1602dde6fd4SLisandro Dalcin dd->elementtype = etype; 1612dde6fd4SLisandro Dalcin dd->ne = 0; 1620298fd71SBarry Smith dd->e = NULL; 1632dde6fd4SLisandro Dalcin } 1642dde6fd4SLisandro Dalcin PetscFunctionReturn(0); 1652dde6fd4SLisandro Dalcin } 1662dde6fd4SLisandro Dalcin 1672dde6fd4SLisandro Dalcin /*@C 1682dde6fd4SLisandro Dalcin DMDAGetElementType - Gets the element type to be returned by DMDAGetElements() 1692dde6fd4SLisandro Dalcin 1702dde6fd4SLisandro Dalcin Not Collective 1712dde6fd4SLisandro Dalcin 1722dde6fd4SLisandro Dalcin Input Parameter: 1732dde6fd4SLisandro Dalcin . da - the DMDA object 1742dde6fd4SLisandro Dalcin 1752dde6fd4SLisandro Dalcin Output Parameters: 1762dde6fd4SLisandro Dalcin . etype - the element type, currently either DMDA_ELEMENT_P1 or ELEMENT_Q1 1772dde6fd4SLisandro Dalcin 1782dde6fd4SLisandro Dalcin Level: intermediate 1792dde6fd4SLisandro Dalcin 1802dde6fd4SLisandro Dalcin .seealso: DMDAElementType, DMDASetElementType(), DMDAGetElements(), DMDARestoreElements() 1812dde6fd4SLisandro Dalcin @*/ 1822dde6fd4SLisandro Dalcin #undef __FUNCT__ 1832dde6fd4SLisandro Dalcin #define __FUNCT__ "DMDAGetElementType" 1842dde6fd4SLisandro Dalcin PetscErrorCode DMDAGetElementType(DM da, DMDAElementType *etype) 1852dde6fd4SLisandro Dalcin { 1862dde6fd4SLisandro Dalcin DM_DA *dd = (DM_DA*)da->data; 1875fd66863SKarl Rupp 1882dde6fd4SLisandro Dalcin PetscFunctionBegin; 1892dde6fd4SLisandro Dalcin PetscValidHeaderSpecific(da,DM_CLASSID,1); 1902dde6fd4SLisandro Dalcin PetscValidPointer(etype,2); 1912dde6fd4SLisandro Dalcin *etype = dd->elementtype; 1922dde6fd4SLisandro Dalcin PetscFunctionReturn(0); 1932dde6fd4SLisandro Dalcin } 1942dde6fd4SLisandro Dalcin 1952dde6fd4SLisandro Dalcin /*@C 1962dde6fd4SLisandro Dalcin DMDAGetElements - Gets an array containing the indices (in local coordinates) 1972dde6fd4SLisandro Dalcin of all the local elements 1982dde6fd4SLisandro Dalcin 1992dde6fd4SLisandro Dalcin Not Collective 2002dde6fd4SLisandro Dalcin 2012dde6fd4SLisandro Dalcin Input Parameter: 2022dde6fd4SLisandro Dalcin . dm - the DM object 2032dde6fd4SLisandro Dalcin 2042dde6fd4SLisandro Dalcin Output Parameters: 2052dde6fd4SLisandro Dalcin + nel - number of local elements 2062dde6fd4SLisandro Dalcin . nen - number of element nodes 207c2cd2aa3SJed Brown - e - the local indices of the elements' vertices 2082dde6fd4SLisandro Dalcin 2092dde6fd4SLisandro Dalcin Level: intermediate 2102dde6fd4SLisandro Dalcin 211*009f0576SBarry Smith Notes: Each process uniquely owns a subset of the elements. That is no element is owned by two or more processes. 212*009f0576SBarry Smith 213*009f0576SBarry 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. 214*009f0576SBarry Smith 215c2cd2aa3SJed Brown .seealso: DMDAElementType, DMDASetElementType(), DMDARestoreElements(), VecSetValuesLocal(), MatSetValuesLocal(), DMGlobalToLocalBegin(), DMLocalToGlobalBegin() 2162dde6fd4SLisandro Dalcin @*/ 2172dde6fd4SLisandro Dalcin #undef __FUNCT__ 2182dde6fd4SLisandro Dalcin #define __FUNCT__ "DMDAGetElements" 2192dde6fd4SLisandro Dalcin PetscErrorCode DMDAGetElements(DM dm,PetscInt *nel,PetscInt *nen,const PetscInt *e[]) 220454e267fSLisandro Dalcin { 221c73cfb54SMatthew G. Knepley PetscInt dim; 222454e267fSLisandro Dalcin PetscErrorCode ierr; 2235fd66863SKarl Rupp 224454e267fSLisandro Dalcin PetscFunctionBegin; 225c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 226c73cfb54SMatthew G. Knepley if (dim==-1) { 2270298fd71SBarry Smith *nel = 0; *nen = 0; *e = NULL; 228c73cfb54SMatthew G. Knepley } else if (dim==1) { 2292dde6fd4SLisandro Dalcin ierr = DMDAGetElements_1D(dm,nel,nen,e);CHKERRQ(ierr); 230c73cfb54SMatthew G. Knepley } else if (dim==2) { 2312dde6fd4SLisandro Dalcin ierr = DMDAGetElements_2D(dm,nel,nen,e);CHKERRQ(ierr); 232c73cfb54SMatthew G. Knepley } else if (dim==3) { 2332dde6fd4SLisandro Dalcin ierr = DMDAGetElements_3D(dm,nel,nen,e);CHKERRQ(ierr); 234c73cfb54SMatthew G. Knepley } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim); 235454e267fSLisandro Dalcin PetscFunctionReturn(0); 236454e267fSLisandro Dalcin } 2372dde6fd4SLisandro Dalcin 2382dde6fd4SLisandro Dalcin #undef __FUNCT__ 2392dde6fd4SLisandro Dalcin #define __FUNCT__ "DMDARestoreElements" 2402dde6fd4SLisandro Dalcin /*@C 241*009f0576SBarry Smith DMDARestoreElements - Restores the array obtained with DMDAGetElements() 2422dde6fd4SLisandro Dalcin 2432dde6fd4SLisandro Dalcin Not Collective 2442dde6fd4SLisandro Dalcin 2452dde6fd4SLisandro Dalcin Input Parameter: 2462dde6fd4SLisandro Dalcin + dm - the DM object 2472dde6fd4SLisandro Dalcin . nel - number of local elements 2482dde6fd4SLisandro Dalcin . nen - number of element nodes 249c2cd2aa3SJed Brown - e - the local indices of the elements' vertices 2502dde6fd4SLisandro Dalcin 2512dde6fd4SLisandro Dalcin Level: intermediate 2522dde6fd4SLisandro Dalcin 253*009f0576SBarry Smith Note: You should not access these values after you have called this routine. 254*009f0576SBarry Smith 255*009f0576SBarry Smith This returns signals the DMDA object that you no longer need access to the array information. 256*009f0576SBarry Smith 2572dde6fd4SLisandro Dalcin .seealso: DMDAElementType, DMDASetElementType(), DMDAGetElements() 2582dde6fd4SLisandro Dalcin @*/ 2592dde6fd4SLisandro Dalcin PetscErrorCode DMDARestoreElements(DM dm,PetscInt *nel,PetscInt *nen,const PetscInt *e[]) 2602dde6fd4SLisandro Dalcin { 2612dde6fd4SLisandro Dalcin PetscFunctionBegin; 2622dde6fd4SLisandro Dalcin PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2632dde6fd4SLisandro Dalcin PetscValidIntPointer(nel,2); 2642dde6fd4SLisandro Dalcin PetscValidIntPointer(nen,3); 2652dde6fd4SLisandro Dalcin PetscValidPointer(e,4); 2669cc8bc54SJed Brown *nel = 0; 2679cc8bc54SJed Brown *nen = -1; 2689cc8bc54SJed Brown *e = NULL; 2692dde6fd4SLisandro Dalcin PetscFunctionReturn(0); 2702dde6fd4SLisandro Dalcin } 271