1c4762a1bSJed Brown 2c4762a1bSJed Brown const char help[] = "Test DMPlexCoordinatesToReference().\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscdm.h> 5c4762a1bSJed Brown #include <petscdmplex.h> 6c4762a1bSJed Brown 7c4762a1bSJed Brown static PetscErrorCode testIdentity(DM dm, PetscBool dmIsSimplicial, PetscInt cell, PetscRandom randCtx, PetscInt numPoints, PetscReal tol) 8c4762a1bSJed Brown { 9c4762a1bSJed Brown PetscInt i, j, dimC, dimR; 10c4762a1bSJed Brown PetscReal *preimage, *mapped, *inverted; 11c4762a1bSJed Brown PetscErrorCode ierr; 12c4762a1bSJed Brown 13c4762a1bSJed Brown PetscFunctionBegin; 14c4762a1bSJed Brown ierr = DMGetDimension(dm,&dimR);CHKERRQ(ierr); 15c4762a1bSJed Brown ierr = DMGetCoordinateDim(dm,&dimC);CHKERRQ(ierr); 16c4762a1bSJed Brown 17c4762a1bSJed Brown ierr = DMGetWorkArray(dm,dimR * numPoints,MPIU_REAL,&preimage);CHKERRQ(ierr); 18c4762a1bSJed Brown ierr = DMGetWorkArray(dm,dimC * numPoints,MPIU_REAL,&mapped);CHKERRQ(ierr); 19c4762a1bSJed Brown ierr = DMGetWorkArray(dm,dimR * numPoints,MPIU_REAL,&inverted);CHKERRQ(ierr); 20c4762a1bSJed Brown 21c4762a1bSJed Brown for (i = 0; i < dimR * numPoints; i++) { 22c4762a1bSJed Brown ierr = PetscRandomGetValueReal(randCtx, &preimage[i]);CHKERRQ(ierr); 23c4762a1bSJed Brown } 24c4762a1bSJed Brown if (dmIsSimplicial && dimR > 1) { 25c4762a1bSJed Brown for (i = 0; i < numPoints; i++) { 26c4762a1bSJed Brown for (j = 0; j < dimR; j++) { 27c4762a1bSJed Brown PetscReal x = preimage[i * dimR + j]; 28c4762a1bSJed Brown PetscReal y = preimage[i * dimR + ((j + 1) % dimR)]; 29c4762a1bSJed Brown 30c4762a1bSJed Brown preimage[i * dimR + ((j + 1) % dimR)] = -1. + (y + 1.) * 0.5 * (1. - x); 31c4762a1bSJed Brown } 32c4762a1bSJed Brown } 33c4762a1bSJed Brown } 34c4762a1bSJed Brown 35c4762a1bSJed Brown ierr = DMPlexReferenceToCoordinates(dm,cell,numPoints,preimage,mapped);CHKERRQ(ierr); 36c4762a1bSJed Brown ierr = DMPlexCoordinatesToReference(dm,cell,numPoints,mapped,inverted);CHKERRQ(ierr); 37c4762a1bSJed Brown 38c4762a1bSJed Brown for (i = 0; i < numPoints; i++) { 39c4762a1bSJed Brown PetscReal max = 0.; 40c4762a1bSJed Brown 41c4762a1bSJed Brown for (j = 0; j < dimR; j++) { 42c4762a1bSJed Brown max = PetscMax(max,PetscAbsReal(preimage[i * dimR + j] - inverted[i * dimR + j])); 43c4762a1bSJed Brown } 44c4762a1bSJed Brown if (max > tol) { 45c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Bad inversion for cell %D with error %g (tol %g): (",cell,(double)max,(double)tol);CHKERRQ(ierr); 46c4762a1bSJed Brown for (j = 0; j < dimR; j++) { 47c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"%+f",(double) preimage[i * dimR + j]);CHKERRQ(ierr); 48c4762a1bSJed Brown if (j < dimR - 1) { 49c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,",");CHKERRQ(ierr); 50c4762a1bSJed Brown } 51c4762a1bSJed Brown } 52c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,") --> (");CHKERRQ(ierr); 53c4762a1bSJed Brown for (j = 0; j < dimC; j++) { 54c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"%+f",(double) mapped[i * dimC + j]);CHKERRQ(ierr); 55c4762a1bSJed Brown if (j < dimC - 1) { 56c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,",");CHKERRQ(ierr); 57c4762a1bSJed Brown } 58c4762a1bSJed Brown } 59c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,") --> (");CHKERRQ(ierr); 60c4762a1bSJed Brown for (j = 0; j < dimR; j++) { 61c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"%+f",(double) inverted[i * dimR + j]);CHKERRQ(ierr); 62c4762a1bSJed Brown if (j < dimR - 1) { 63c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,",");CHKERRQ(ierr); 64c4762a1bSJed Brown } 65c4762a1bSJed Brown } 66c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,")\n");CHKERRQ(ierr); 67c4762a1bSJed Brown } else { 68c4762a1bSJed Brown char strBuf[BUFSIZ] = {'\0'}; 69c4762a1bSJed Brown size_t offset = 0, count; 70c4762a1bSJed Brown 71c4762a1bSJed Brown ierr = PetscSNPrintfCount(strBuf + offset,BUFSIZ-offset,"Good inversion for cell %D: (",&count,cell);CHKERRQ(ierr); 72c4762a1bSJed Brown offset += count; 73c4762a1bSJed Brown for (j = 0; j < dimR; j++) { 74c4762a1bSJed Brown ierr = PetscSNPrintfCount(strBuf + offset,BUFSIZ-offset,"%+f",&count,(double) preimage[i * dimR + j]);CHKERRQ(ierr); 75c4762a1bSJed Brown offset += count; 76c4762a1bSJed Brown if (j < dimR - 1) { 77c4762a1bSJed Brown ierr = PetscSNPrintfCount(strBuf + offset,BUFSIZ-offset,",",&count);CHKERRQ(ierr); 78c4762a1bSJed Brown offset += count; 79c4762a1bSJed Brown } 80c4762a1bSJed Brown } 81c4762a1bSJed Brown ierr = PetscSNPrintfCount(strBuf + offset,BUFSIZ-offset,") --> (",&count);CHKERRQ(ierr); 82c4762a1bSJed Brown offset += count; 83c4762a1bSJed Brown for (j = 0; j < dimC; j++) { 84c4762a1bSJed Brown ierr = PetscSNPrintfCount(strBuf + offset,BUFSIZ-offset,"%+f",&count,(double) mapped[i * dimC + j]);CHKERRQ(ierr); 85c4762a1bSJed Brown offset += count; 86c4762a1bSJed Brown if (j < dimC - 1) { 87c4762a1bSJed Brown ierr = PetscSNPrintfCount(strBuf + offset,BUFSIZ-offset,",",&count);CHKERRQ(ierr); 88c4762a1bSJed Brown offset += count; 89c4762a1bSJed Brown } 90c4762a1bSJed Brown } 91c4762a1bSJed Brown ierr = PetscSNPrintfCount(strBuf + offset,BUFSIZ-offset,") --> (",&count);CHKERRQ(ierr); 92c4762a1bSJed Brown offset += count; 93c4762a1bSJed Brown for (j = 0; j < dimR; j++) { 94c4762a1bSJed Brown ierr = PetscSNPrintfCount(strBuf + offset,BUFSIZ-offset,"%+f",&count,(double) inverted[i * dimR + j]);CHKERRQ(ierr); 95c4762a1bSJed Brown offset += count; 96c4762a1bSJed Brown if (j < dimR - 1) { 97c4762a1bSJed Brown ierr = PetscSNPrintfCount(strBuf + offset,BUFSIZ-offset,",",&count);CHKERRQ(ierr); 98c4762a1bSJed Brown offset += count; 99c4762a1bSJed Brown } 100c4762a1bSJed Brown } 101c4762a1bSJed Brown ierr = PetscSNPrintfCount(strBuf + offset,BUFSIZ-offset,")\n",&count);CHKERRQ(ierr); 1027d3de750SJacob Faibussowitsch ierr = PetscInfo(dm,"%s",strBuf);CHKERRQ(ierr); 103c4762a1bSJed Brown } 104c4762a1bSJed Brown } 105c4762a1bSJed Brown 106c4762a1bSJed Brown ierr = DMRestoreWorkArray(dm,dimR * numPoints,MPIU_REAL,&inverted);CHKERRQ(ierr); 107c4762a1bSJed Brown ierr = DMRestoreWorkArray(dm,dimC * numPoints,MPIU_REAL,&mapped);CHKERRQ(ierr); 108c4762a1bSJed Brown ierr = DMRestoreWorkArray(dm,dimR * numPoints,MPIU_REAL,&preimage);CHKERRQ(ierr); 109c4762a1bSJed Brown PetscFunctionReturn(0); 110c4762a1bSJed Brown } 111c4762a1bSJed Brown 112c4762a1bSJed Brown static PetscErrorCode identityEmbedding(PetscInt dim, PetscReal time, const PetscReal *x, PetscInt Nf, PetscScalar *u, void *ctx) 113c4762a1bSJed Brown { 114c4762a1bSJed Brown PetscInt i; 115c4762a1bSJed Brown 116c4762a1bSJed Brown for (i = 0; i < dim; i++) {u[i] = x[i];}; 117c4762a1bSJed Brown return 0; 118c4762a1bSJed Brown } 119c4762a1bSJed Brown 120c4762a1bSJed Brown int main(int argc, char **argv) 121c4762a1bSJed Brown { 122c4762a1bSJed Brown PetscRandom randCtx; 123c4762a1bSJed Brown PetscInt dim, dimC, isSimplex, isFE, numTests = 10; 124c4762a1bSJed Brown PetscReal perturb = 0.1, tol = 10. * PETSC_SMALL; 125c4762a1bSJed Brown PetscErrorCode ierr; 126c4762a1bSJed Brown 127c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr; 128c4762a1bSJed Brown ierr = PetscRandomCreate(PETSC_COMM_SELF,&randCtx);CHKERRQ(ierr); 129c4762a1bSJed Brown ierr = PetscRandomSetInterval(randCtx,-1.,1.);CHKERRQ(ierr); 130c4762a1bSJed Brown ierr = PetscRandomSetFromOptions(randCtx);CHKERRQ(ierr); 131c4762a1bSJed Brown ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"ex21",NULL);CHKERRQ(ierr); 132c4762a1bSJed Brown ierr = PetscOptionsReal("-vertex_perturbation","scale of random vertex distortion",NULL,perturb,&perturb,NULL);CHKERRQ(ierr); 133c4762a1bSJed Brown ierr = PetscOptionsBoundedInt("-num_test_points","number of points to test",NULL,numTests,&numTests,NULL,0);CHKERRQ(ierr); 134c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 135c4762a1bSJed Brown for (dim = 1; dim <= 3; dim++) { 136c4762a1bSJed Brown for (dimC = dim; dimC <= PetscMin(3,dim + 1); dimC++) { 137c4762a1bSJed Brown for (isSimplex = 0; isSimplex < 2; isSimplex++) { 138c4762a1bSJed Brown for (isFE = 0; isFE < 2; isFE++) { 139c4762a1bSJed Brown DM dm; 140c4762a1bSJed Brown Vec coords; 141c4762a1bSJed Brown PetscScalar *coordArray; 142c4762a1bSJed Brown PetscReal noise; 143c4762a1bSJed Brown PetscInt i, n, order = 1; 144c4762a1bSJed Brown 14530602db0SMatthew G. Knepley ierr = DMPlexCreateReferenceCell(PETSC_COMM_SELF, DMPolytopeTypeSimpleShape(dim, isSimplex ? PETSC_TRUE : PETSC_FALSE), &dm);CHKERRQ(ierr); 146c4762a1bSJed Brown if (isFE) { 147c4762a1bSJed Brown DM dmCoord; 148c4762a1bSJed Brown PetscSpace sp; 149c4762a1bSJed Brown PetscFE fe; 150c4762a1bSJed Brown Vec localCoords; 151c4762a1bSJed Brown PetscErrorCode (*funcs[1])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *) = {identityEmbedding}; 152c4762a1bSJed Brown void *ctxs[1] = {NULL}; 153c4762a1bSJed Brown 154c4762a1bSJed Brown ierr = PetscFECreateDefault(PetscObjectComm((PetscObject) dm),dim,dim,isSimplex ? PETSC_TRUE : PETSC_FALSE,isSimplex ? NULL : "tensor_" ,PETSC_DEFAULT,&fe);CHKERRQ(ierr); 155c4762a1bSJed Brown ierr = PetscFEGetBasisSpace(fe,&sp);CHKERRQ(ierr); 156c4762a1bSJed Brown ierr = PetscSpaceGetDegree(sp,&order,NULL);CHKERRQ(ierr); 157c4762a1bSJed Brown ierr = DMSetField(dm,0,NULL,(PetscObject)fe);CHKERRQ(ierr); 158c4762a1bSJed Brown ierr = DMCreateDS(dm);CHKERRQ(ierr); 159c4762a1bSJed Brown ierr = DMCreateLocalVector(dm,&localCoords);CHKERRQ(ierr); 160c4762a1bSJed Brown ierr = DMProjectFunctionLocal(dm,0,funcs,ctxs,INSERT_VALUES,localCoords);CHKERRQ(ierr); 161c4762a1bSJed Brown ierr = VecSetDM(localCoords,NULL);CHKERRQ(ierr); /* This is necessary to prevent a reference loop */ 162c4762a1bSJed Brown ierr = DMClone(dm,&dmCoord);CHKERRQ(ierr); 163c4762a1bSJed Brown ierr = DMSetField(dmCoord,0,NULL,(PetscObject)fe);CHKERRQ(ierr); 164c4762a1bSJed Brown ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 165c4762a1bSJed Brown ierr = DMCreateDS(dmCoord);CHKERRQ(ierr); 166c4762a1bSJed Brown ierr = DMSetCoordinateDM(dm,dmCoord);CHKERRQ(ierr); 167c4762a1bSJed Brown ierr = DMDestroy(&dmCoord);CHKERRQ(ierr); 168c4762a1bSJed Brown ierr = DMSetCoordinatesLocal(dm,localCoords);CHKERRQ(ierr); 169c4762a1bSJed Brown ierr = VecDestroy(&localCoords);CHKERRQ(ierr); 170c4762a1bSJed Brown } 1717d3de750SJacob Faibussowitsch ierr = PetscInfo(dm,"Testing %s%D %DD mesh embedded in %DD\n",isSimplex ? "P" : "Q" , order, dim, dimC);CHKERRQ(ierr); 172c4762a1bSJed Brown ierr = DMGetCoordinatesLocal(dm,&coords);CHKERRQ(ierr); 173c4762a1bSJed Brown ierr = VecGetLocalSize(coords,&n);CHKERRQ(ierr); 174c4762a1bSJed Brown if (dimC > dim) { /* reembed in higher dimension */ 175c4762a1bSJed Brown PetscSection sec, newSec; 176c4762a1bSJed Brown PetscInt pStart, pEnd, p, i, newN; 177c4762a1bSJed Brown Vec newVec; 178c4762a1bSJed Brown DM coordDM; 179c4762a1bSJed Brown PetscScalar *newCoordArray; 180c4762a1bSJed Brown 181c4762a1bSJed Brown ierr = DMGetCoordinateSection(dm,&sec);CHKERRQ(ierr); 182c4762a1bSJed Brown ierr = PetscSectionCreate(PetscObjectComm((PetscObject)sec),&newSec);CHKERRQ(ierr); 183c4762a1bSJed Brown ierr = PetscSectionSetNumFields(newSec,1);CHKERRQ(ierr); 184c4762a1bSJed Brown ierr = PetscSectionGetChart(sec,&pStart,&pEnd);CHKERRQ(ierr); 185c4762a1bSJed Brown ierr = PetscSectionSetChart(newSec,pStart,pEnd);CHKERRQ(ierr); 186c4762a1bSJed Brown for (p = pStart; p < pEnd; p++) { 187c4762a1bSJed Brown PetscInt nDof; 188c4762a1bSJed Brown 189c4762a1bSJed Brown ierr = PetscSectionGetDof(sec,p,&nDof);CHKERRQ(ierr); 190*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(nDof % dim,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Coordinate section point %D has %D dofs != 0 mod %D",p,nDof,dim); 191c4762a1bSJed Brown ierr = PetscSectionSetDof(newSec,p,(nDof/dim)*dimC);CHKERRQ(ierr); 192c4762a1bSJed Brown ierr = PetscSectionSetFieldDof(newSec,p,0,(nDof/dim)*dimC);CHKERRQ(ierr); 193c4762a1bSJed Brown } 194c4762a1bSJed Brown ierr = PetscSectionSetUp(newSec);CHKERRQ(ierr); 195c4762a1bSJed Brown ierr = PetscSectionGetStorageSize(newSec,&newN);CHKERRQ(ierr); 196c4762a1bSJed Brown ierr = VecCreateSeq(PETSC_COMM_SELF,newN,&newVec);CHKERRQ(ierr); 197c4762a1bSJed Brown ierr = VecGetArray(newVec,&newCoordArray);CHKERRQ(ierr); 198c4762a1bSJed Brown ierr = VecGetArray(coords,&coordArray);CHKERRQ(ierr); 199c4762a1bSJed Brown for (i = 0; i < n / dim; i++) { 200c4762a1bSJed Brown PetscInt j; 201c4762a1bSJed Brown 202c4762a1bSJed Brown for (j = 0; j < dim; j++) { 203c4762a1bSJed Brown newCoordArray[i * dimC + j] = coordArray[i * dim + j]; 204c4762a1bSJed Brown } 205c4762a1bSJed Brown for (; j < dimC; j++) { 206c4762a1bSJed Brown newCoordArray[i * dimC + j] = 0.; 207c4762a1bSJed Brown } 208c4762a1bSJed Brown } 209c4762a1bSJed Brown ierr = VecRestoreArray(coords,&coordArray);CHKERRQ(ierr); 210c4762a1bSJed Brown ierr = VecRestoreArray(newVec,&newCoordArray);CHKERRQ(ierr); 211c4762a1bSJed Brown ierr = DMSetCoordinateDim(dm,dimC);CHKERRQ(ierr); 212c4762a1bSJed Brown ierr = DMSetCoordinateSection(dm,dimC,newSec);CHKERRQ(ierr); 213c4762a1bSJed Brown ierr = DMSetCoordinatesLocal(dm,newVec);CHKERRQ(ierr); 214c4762a1bSJed Brown ierr = VecDestroy(&newVec);CHKERRQ(ierr); 215c4762a1bSJed Brown ierr = PetscSectionDestroy(&newSec);CHKERRQ(ierr); 216c4762a1bSJed Brown ierr = DMGetCoordinatesLocal(dm,&coords);CHKERRQ(ierr); 217c4762a1bSJed Brown ierr = DMGetCoordinateDM(dm,&coordDM);CHKERRQ(ierr); 218c4762a1bSJed Brown if (isFE) { 219c4762a1bSJed Brown PetscFE fe; 220c4762a1bSJed Brown 221c4762a1bSJed Brown ierr = PetscFECreateDefault(PetscObjectComm((PetscObject) dm),dim,dimC,isSimplex ? PETSC_TRUE : PETSC_FALSE,isSimplex ? NULL : "tensor_",PETSC_DEFAULT,&fe);CHKERRQ(ierr); 222c4762a1bSJed Brown ierr = DMSetField(coordDM,0,NULL,(PetscObject)fe);CHKERRQ(ierr); 223c4762a1bSJed Brown ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 224c4762a1bSJed Brown ierr = DMCreateDS(coordDM);CHKERRQ(ierr); 225c4762a1bSJed Brown } 226c4762a1bSJed Brown ierr = DMSetCoordinateDim(coordDM,dimC);CHKERRQ(ierr); 227c4762a1bSJed Brown ierr = VecGetLocalSize(coords,&n);CHKERRQ(ierr); 228c4762a1bSJed Brown } 229c4762a1bSJed Brown ierr = VecGetArray(coords,&coordArray);CHKERRQ(ierr); 230c4762a1bSJed Brown for (i = 0; i < n; i++) { 231c4762a1bSJed Brown ierr = PetscRandomGetValueReal(randCtx,&noise);CHKERRQ(ierr); 232c4762a1bSJed Brown coordArray[i] += noise * perturb; 233c4762a1bSJed Brown } 234c4762a1bSJed Brown ierr = VecRestoreArray(coords,&coordArray);CHKERRQ(ierr); 235c4762a1bSJed Brown ierr = DMSetCoordinatesLocal(dm, coords);CHKERRQ(ierr); 236c4762a1bSJed Brown ierr = testIdentity(dm, isSimplex ? PETSC_TRUE : PETSC_FALSE, 0, randCtx, numTests, tol);CHKERRQ(ierr); 237c4762a1bSJed Brown ierr = DMDestroy(&dm);CHKERRQ(ierr); 238c4762a1bSJed Brown } 239c4762a1bSJed Brown } 240c4762a1bSJed Brown } 241c4762a1bSJed Brown } 242c4762a1bSJed Brown ierr = PetscRandomDestroy(&randCtx);CHKERRQ(ierr); 243c4762a1bSJed Brown ierr = PetscFinalize(); 244c4762a1bSJed Brown return ierr; 245c4762a1bSJed Brown } 246c4762a1bSJed Brown 247c4762a1bSJed Brown /*TEST 248c4762a1bSJed Brown 249c4762a1bSJed Brown test: 250c4762a1bSJed Brown suffix: 0 251c4762a1bSJed Brown args: -petscspace_degree 2 -tensor_petscspace_degree 2 252c4762a1bSJed Brown 253c4762a1bSJed Brown TEST*/ 254