1c4762a1bSJed Brown const char help[] = "Test DMPlexCoordinatesToReference().\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscdm.h> 4c4762a1bSJed Brown #include <petscdmplex.h> 5c4762a1bSJed Brown 6d71ae5a4SJacob Faibussowitsch static PetscErrorCode testIdentity(DM dm, PetscBool dmIsSimplicial, PetscInt cell, PetscRandom randCtx, PetscInt numPoints, PetscReal tol) 7d71ae5a4SJacob Faibussowitsch { 8c4762a1bSJed Brown PetscInt i, j, dimC, dimR; 9c4762a1bSJed Brown PetscReal *preimage, *mapped, *inverted; 10c4762a1bSJed Brown 11c4762a1bSJed Brown PetscFunctionBegin; 129566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dimR)); 139566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &dimC)); 14c4762a1bSJed Brown 159566063dSJacob Faibussowitsch PetscCall(DMGetWorkArray(dm, dimR * numPoints, MPIU_REAL, &preimage)); 169566063dSJacob Faibussowitsch PetscCall(DMGetWorkArray(dm, dimC * numPoints, MPIU_REAL, &mapped)); 179566063dSJacob Faibussowitsch PetscCall(DMGetWorkArray(dm, dimR * numPoints, MPIU_REAL, &inverted)); 18c4762a1bSJed Brown 1948a46eb9SPierre Jolivet for (i = 0; i < dimR * numPoints; i++) PetscCall(PetscRandomGetValueReal(randCtx, &preimage[i])); 20c4762a1bSJed Brown if (dmIsSimplicial && dimR > 1) { 21c4762a1bSJed Brown for (i = 0; i < numPoints; i++) { 22c4762a1bSJed Brown for (j = 0; j < dimR; j++) { 23c4762a1bSJed Brown PetscReal x = preimage[i * dimR + j]; 24c4762a1bSJed Brown PetscReal y = preimage[i * dimR + ((j + 1) % dimR)]; 25c4762a1bSJed Brown 26c4762a1bSJed Brown preimage[i * dimR + ((j + 1) % dimR)] = -1. + (y + 1.) * 0.5 * (1. - x); 27c4762a1bSJed Brown } 28c4762a1bSJed Brown } 29c4762a1bSJed Brown } 30c4762a1bSJed Brown 319566063dSJacob Faibussowitsch PetscCall(DMPlexReferenceToCoordinates(dm, cell, numPoints, preimage, mapped)); 329566063dSJacob Faibussowitsch PetscCall(DMPlexCoordinatesToReference(dm, cell, numPoints, mapped, inverted)); 33c4762a1bSJed Brown 34c4762a1bSJed Brown for (i = 0; i < numPoints; i++) { 35c4762a1bSJed Brown PetscReal max = 0.; 36c4762a1bSJed Brown 37ad540459SPierre Jolivet for (j = 0; j < dimR; j++) max = PetscMax(max, PetscAbsReal(preimage[i * dimR + j] - inverted[i * dimR + j])); 38c4762a1bSJed Brown if (max > tol) { 3963a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "Bad inversion for cell %" PetscInt_FMT " with error %g (tol %g): (", cell, (double)max, (double)tol)); 40c4762a1bSJed Brown for (j = 0; j < dimR; j++) { 419566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "%+f", (double)preimage[i * dimR + j])); 4248a46eb9SPierre Jolivet if (j < dimR - 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, ",")); 43c4762a1bSJed Brown } 449566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, ") --> (")); 45c4762a1bSJed Brown for (j = 0; j < dimC; j++) { 469566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "%+f", (double)mapped[i * dimC + j])); 4748a46eb9SPierre Jolivet if (j < dimC - 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, ",")); 48c4762a1bSJed Brown } 499566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, ") --> (")); 50c4762a1bSJed Brown for (j = 0; j < dimR; j++) { 519566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "%+f", (double)inverted[i * dimR + j])); 5248a46eb9SPierre Jolivet if (j < dimR - 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, ",")); 53c4762a1bSJed Brown } 549566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, ")\n")); 55c4762a1bSJed Brown } else { 56c4762a1bSJed Brown char strBuf[BUFSIZ] = {'\0'}; 57c4762a1bSJed Brown size_t offset = 0, count; 58c4762a1bSJed Brown 5963a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintfCount(strBuf + offset, BUFSIZ - offset, "Good inversion for cell %" PetscInt_FMT ": (", &count, cell)); 60c4762a1bSJed Brown offset += count; 61c4762a1bSJed Brown for (j = 0; j < dimR; j++) { 629566063dSJacob Faibussowitsch PetscCall(PetscSNPrintfCount(strBuf + offset, BUFSIZ - offset, "%+f", &count, (double)preimage[i * dimR + j])); 63c4762a1bSJed Brown offset += count; 64c4762a1bSJed Brown if (j < dimR - 1) { 659566063dSJacob Faibussowitsch PetscCall(PetscSNPrintfCount(strBuf + offset, BUFSIZ - offset, ",", &count)); 66c4762a1bSJed Brown offset += count; 67c4762a1bSJed Brown } 68c4762a1bSJed Brown } 699566063dSJacob Faibussowitsch PetscCall(PetscSNPrintfCount(strBuf + offset, BUFSIZ - offset, ") --> (", &count)); 70c4762a1bSJed Brown offset += count; 71c4762a1bSJed Brown for (j = 0; j < dimC; j++) { 729566063dSJacob Faibussowitsch PetscCall(PetscSNPrintfCount(strBuf + offset, BUFSIZ - offset, "%+f", &count, (double)mapped[i * dimC + j])); 73c4762a1bSJed Brown offset += count; 74c4762a1bSJed Brown if (j < dimC - 1) { 759566063dSJacob Faibussowitsch PetscCall(PetscSNPrintfCount(strBuf + offset, BUFSIZ - offset, ",", &count)); 76c4762a1bSJed Brown offset += count; 77c4762a1bSJed Brown } 78c4762a1bSJed Brown } 799566063dSJacob Faibussowitsch PetscCall(PetscSNPrintfCount(strBuf + offset, BUFSIZ - offset, ") --> (", &count)); 80c4762a1bSJed Brown offset += count; 81c4762a1bSJed Brown for (j = 0; j < dimR; j++) { 829566063dSJacob Faibussowitsch PetscCall(PetscSNPrintfCount(strBuf + offset, BUFSIZ - offset, "%+f", &count, (double)inverted[i * dimR + j])); 83c4762a1bSJed Brown offset += count; 84c4762a1bSJed Brown if (j < dimR - 1) { 859566063dSJacob Faibussowitsch PetscCall(PetscSNPrintfCount(strBuf + offset, BUFSIZ - offset, ",", &count)); 86c4762a1bSJed Brown offset += count; 87c4762a1bSJed Brown } 88c4762a1bSJed Brown } 899d3446b2SPierre Jolivet PetscCall(PetscSNPrintfCount(strBuf + offset, BUFSIZ - offset, ")", &count)); 909d3446b2SPierre Jolivet PetscCall(PetscInfo(dm, "%s\n", strBuf)); 91c4762a1bSJed Brown } 92c4762a1bSJed Brown } 93c4762a1bSJed Brown 949566063dSJacob Faibussowitsch PetscCall(DMRestoreWorkArray(dm, dimR * numPoints, MPIU_REAL, &inverted)); 959566063dSJacob Faibussowitsch PetscCall(DMRestoreWorkArray(dm, dimC * numPoints, MPIU_REAL, &mapped)); 969566063dSJacob Faibussowitsch PetscCall(DMRestoreWorkArray(dm, dimR * numPoints, MPIU_REAL, &preimage)); 973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 98c4762a1bSJed Brown } 99c4762a1bSJed Brown 100*2a8381b2SBarry Smith static PetscErrorCode identityEmbedding(PetscInt dim, PetscReal time, const PetscReal *x, PetscInt Nf, PetscScalar *u, PetscCtx ctx) 101d71ae5a4SJacob Faibussowitsch { 102c4762a1bSJed Brown PetscInt i; 103c4762a1bSJed Brown 104ad540459SPierre Jolivet for (i = 0; i < dim; i++) u[i] = x[i]; 1053ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 106c4762a1bSJed Brown } 107c4762a1bSJed Brown 108d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 109d71ae5a4SJacob Faibussowitsch { 110c4762a1bSJed Brown PetscRandom randCtx; 111c4762a1bSJed Brown PetscInt dim, dimC, isSimplex, isFE, numTests = 10; 112c4762a1bSJed Brown PetscReal perturb = 0.1, tol = 10. * PETSC_SMALL; 113c4762a1bSJed Brown 114327415f7SBarry Smith PetscFunctionBeginUser; 1159566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 1169566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &randCtx)); 1179566063dSJacob Faibussowitsch PetscCall(PetscRandomSetInterval(randCtx, -1., 1.)); 1189566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(randCtx)); 119d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "ex21", NULL); 1209566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-vertex_perturbation", "scale of random vertex distortion", NULL, perturb, &perturb, NULL)); 1219566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-num_test_points", "number of points to test", NULL, numTests, &numTests, NULL, 0)); 122d0609cedSBarry Smith PetscOptionsEnd(); 123c4762a1bSJed Brown for (dim = 1; dim <= 3; dim++) { 124c4762a1bSJed Brown for (dimC = dim; dimC <= PetscMin(3, dim + 1); dimC++) { 125c4762a1bSJed Brown for (isSimplex = 0; isSimplex < 2; isSimplex++) { 126c4762a1bSJed Brown for (isFE = 0; isFE < 2; isFE++) { 127c4762a1bSJed Brown DM dm; 128c4762a1bSJed Brown Vec coords; 129c4762a1bSJed Brown PetscScalar *coordArray; 130c4762a1bSJed Brown PetscReal noise; 131c4762a1bSJed Brown PetscInt i, n, order = 1; 132c4762a1bSJed Brown 1339566063dSJacob Faibussowitsch PetscCall(DMPlexCreateReferenceCell(PETSC_COMM_SELF, DMPolytopeTypeSimpleShape(dim, isSimplex ? PETSC_TRUE : PETSC_FALSE), &dm)); 134c4762a1bSJed Brown if (isFE) { 135c4762a1bSJed Brown DM dmCoord; 136c4762a1bSJed Brown PetscSpace sp; 137c4762a1bSJed Brown PetscFE fe; 138c4762a1bSJed Brown Vec localCoords; 139c4762a1bSJed Brown PetscErrorCode (*funcs[1])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *) = {identityEmbedding}; 140*2a8381b2SBarry Smith PetscCtx ctxs[1] = {NULL}; 141c4762a1bSJed Brown 1429566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)dm), dim, dim, isSimplex ? PETSC_TRUE : PETSC_FALSE, isSimplex ? NULL : "tensor_", PETSC_DEFAULT, &fe)); 1439566063dSJacob Faibussowitsch PetscCall(PetscFEGetBasisSpace(fe, &sp)); 1449566063dSJacob Faibussowitsch PetscCall(PetscSpaceGetDegree(sp, &order, NULL)); 1459566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe)); 1469566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 1479566063dSJacob Faibussowitsch PetscCall(DMCreateLocalVector(dm, &localCoords)); 1489566063dSJacob Faibussowitsch PetscCall(DMProjectFunctionLocal(dm, 0, funcs, ctxs, INSERT_VALUES, localCoords)); 1499566063dSJacob Faibussowitsch PetscCall(VecSetDM(localCoords, NULL)); /* This is necessary to prevent a reference loop */ 1509566063dSJacob Faibussowitsch PetscCall(DMClone(dm, &dmCoord)); 1519566063dSJacob Faibussowitsch PetscCall(DMSetField(dmCoord, 0, NULL, (PetscObject)fe)); 1529566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 1539566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dmCoord)); 1549566063dSJacob Faibussowitsch PetscCall(DMSetCoordinateDM(dm, dmCoord)); 1559566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmCoord)); 1569566063dSJacob Faibussowitsch PetscCall(DMSetCoordinatesLocal(dm, localCoords)); 1579566063dSJacob Faibussowitsch PetscCall(VecDestroy(&localCoords)); 158c4762a1bSJed Brown } 15963a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(dm, "Testing %s%" PetscInt_FMT " %" PetscInt_FMT "D mesh embedded in %" PetscInt_FMT "D\n", isSimplex ? "P" : "Q", order, dim, dimC)); 1609566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coords)); 1619566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(coords, &n)); 162c4762a1bSJed Brown if (dimC > dim) { /* reembed in higher dimension */ 163c4762a1bSJed Brown PetscSection sec, newSec; 164c4762a1bSJed Brown PetscInt pStart, pEnd, p, i, newN; 165c4762a1bSJed Brown Vec newVec; 166c4762a1bSJed Brown DM coordDM; 167c4762a1bSJed Brown PetscScalar *newCoordArray; 168c4762a1bSJed Brown 1699566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(dm, &sec)); 1709566063dSJacob Faibussowitsch PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)sec), &newSec)); 1719566063dSJacob Faibussowitsch PetscCall(PetscSectionSetNumFields(newSec, 1)); 1729566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(sec, &pStart, &pEnd)); 1739566063dSJacob Faibussowitsch PetscCall(PetscSectionSetChart(newSec, pStart, pEnd)); 174c4762a1bSJed Brown for (p = pStart; p < pEnd; p++) { 175c4762a1bSJed Brown PetscInt nDof; 176c4762a1bSJed Brown 1779566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(sec, p, &nDof)); 17863a3b9bcSJacob Faibussowitsch PetscCheck(nDof % dim == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Coordinate section point %" PetscInt_FMT " has %" PetscInt_FMT " dofs != 0 mod %" PetscInt_FMT, p, nDof, dim); 1799566063dSJacob Faibussowitsch PetscCall(PetscSectionSetDof(newSec, p, (nDof / dim) * dimC)); 1809566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldDof(newSec, p, 0, (nDof / dim) * dimC)); 181c4762a1bSJed Brown } 1829566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(newSec)); 1839566063dSJacob Faibussowitsch PetscCall(PetscSectionGetStorageSize(newSec, &newN)); 1849566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, newN, &newVec)); 1859566063dSJacob Faibussowitsch PetscCall(VecGetArray(newVec, &newCoordArray)); 1869566063dSJacob Faibussowitsch PetscCall(VecGetArray(coords, &coordArray)); 187c4762a1bSJed Brown for (i = 0; i < n / dim; i++) { 188c4762a1bSJed Brown PetscInt j; 189c4762a1bSJed Brown 190ad540459SPierre Jolivet for (j = 0; j < dim; j++) newCoordArray[i * dimC + j] = coordArray[i * dim + j]; 191ad540459SPierre Jolivet for (; j < dimC; j++) newCoordArray[i * dimC + j] = 0.; 192c4762a1bSJed Brown } 1939566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(coords, &coordArray)); 1949566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(newVec, &newCoordArray)); 1959566063dSJacob Faibussowitsch PetscCall(DMSetCoordinateDim(dm, dimC)); 1969566063dSJacob Faibussowitsch PetscCall(DMSetCoordinateSection(dm, dimC, newSec)); 1979566063dSJacob Faibussowitsch PetscCall(DMSetCoordinatesLocal(dm, newVec)); 1989566063dSJacob Faibussowitsch PetscCall(VecDestroy(&newVec)); 1999566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(&newSec)); 2009566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coords)); 2019566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dm, &coordDM)); 202c4762a1bSJed Brown if (isFE) { 203c4762a1bSJed Brown PetscFE fe; 204c4762a1bSJed Brown 2059566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)dm), dim, dimC, isSimplex ? PETSC_TRUE : PETSC_FALSE, isSimplex ? NULL : "tensor_", PETSC_DEFAULT, &fe)); 2069566063dSJacob Faibussowitsch PetscCall(DMSetField(coordDM, 0, NULL, (PetscObject)fe)); 2079566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 2089566063dSJacob Faibussowitsch PetscCall(DMCreateDS(coordDM)); 209c4762a1bSJed Brown } 2109566063dSJacob Faibussowitsch PetscCall(DMSetCoordinateDim(coordDM, dimC)); 2119566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(coords, &n)); 212c4762a1bSJed Brown } 2139566063dSJacob Faibussowitsch PetscCall(VecGetArray(coords, &coordArray)); 214c4762a1bSJed Brown for (i = 0; i < n; i++) { 2159566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValueReal(randCtx, &noise)); 216c4762a1bSJed Brown coordArray[i] += noise * perturb; 217c4762a1bSJed Brown } 2189566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(coords, &coordArray)); 2199566063dSJacob Faibussowitsch PetscCall(DMSetCoordinatesLocal(dm, coords)); 2209566063dSJacob Faibussowitsch PetscCall(testIdentity(dm, isSimplex ? PETSC_TRUE : PETSC_FALSE, 0, randCtx, numTests, tol)); 2219566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 222c4762a1bSJed Brown } 223c4762a1bSJed Brown } 224c4762a1bSJed Brown } 225c4762a1bSJed Brown } 2269566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&randCtx)); 2279566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 228b122ec5aSJacob Faibussowitsch return 0; 229c4762a1bSJed Brown } 230c4762a1bSJed Brown 231c4762a1bSJed Brown /*TEST 232c4762a1bSJed Brown 233c4762a1bSJed Brown test: 234c4762a1bSJed Brown suffix: 0 235c4762a1bSJed Brown args: -petscspace_degree 2 -tensor_petscspace_degree 2 2363886731fSPierre Jolivet output_file: output/empty.out 237c4762a1bSJed Brown 238c4762a1bSJed Brown TEST*/ 239