1c4762a1bSJed Brown static char help[] = "Demonstration of creating and viewing DMFields objects.\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscdmfield.h> 4c4762a1bSJed Brown #include <petscdmplex.h> 5c4762a1bSJed Brown #include <petscdmda.h> 6c4762a1bSJed Brown 7d71ae5a4SJacob Faibussowitsch static PetscErrorCode ViewResults(PetscViewer viewer, PetscInt N, PetscInt dim, PetscScalar *B, PetscScalar *D, PetscScalar *H, PetscReal *rB, PetscReal *rD, PetscReal *rH) 8d71ae5a4SJacob Faibussowitsch { 9c4762a1bSJed Brown PetscFunctionBegin; 109566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "B:\n")); 119566063dSJacob Faibussowitsch PetscCall(PetscScalarView(N, B, viewer)); 129566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "D:\n")); 139566063dSJacob Faibussowitsch PetscCall(PetscScalarView(N * dim, D, viewer)); 149566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "H:\n")); 159566063dSJacob Faibussowitsch PetscCall(PetscScalarView(N * dim * dim, H, viewer)); 16c4762a1bSJed Brown 179566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "rB:\n")); 189566063dSJacob Faibussowitsch PetscCall(PetscRealView(N, rB, viewer)); 199566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "rD:\n")); 209566063dSJacob Faibussowitsch PetscCall(PetscRealView(N * dim, rD, viewer)); 219566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "rH:\n")); 229566063dSJacob Faibussowitsch PetscCall(PetscRealView(N * dim * dim, rH, viewer)); 233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 24c4762a1bSJed Brown } 25c4762a1bSJed Brown 26d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestEvaluate(DMField field, PetscInt n, PetscRandom rand) 27d71ae5a4SJacob Faibussowitsch { 28c4762a1bSJed Brown DM dm; 29c4762a1bSJed Brown PetscInt dim, i, nc; 30c4762a1bSJed Brown PetscScalar *B, *D, *H; 31c4762a1bSJed Brown PetscReal *rB, *rD, *rH; 32c4762a1bSJed Brown Vec points; 33c4762a1bSJed Brown PetscScalar *array; 34c4762a1bSJed Brown PetscViewer viewer; 35c4762a1bSJed Brown MPI_Comm comm; 36c4762a1bSJed Brown 37c4762a1bSJed Brown PetscFunctionBegin; 38c4762a1bSJed Brown comm = PetscObjectComm((PetscObject)field); 399566063dSJacob Faibussowitsch PetscCall(DMFieldGetNumComponents(field, &nc)); 409566063dSJacob Faibussowitsch PetscCall(DMFieldGetDM(field, &dm)); 419566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 4277433607SBarry Smith PetscCall(VecCreateFromOptions(PetscObjectComm((PetscObject)field), NULL, 1, n * dim, PETSC_DETERMINE, &points)); 439566063dSJacob Faibussowitsch PetscCall(VecSetBlockSize(points, dim)); 449566063dSJacob Faibussowitsch PetscCall(VecGetArray(points, &array)); 459566063dSJacob Faibussowitsch for (i = 0; i < n * dim; i++) PetscCall(PetscRandomGetValue(rand, &array[i])); 469566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(points, &array)); 479566063dSJacob Faibussowitsch PetscCall(PetscMalloc6(n * nc, &B, n * nc, &rB, n * nc * dim, &D, n * nc * dim, &rD, n * nc * dim * dim, &H, n * nc * dim * dim, &rH)); 489566063dSJacob Faibussowitsch PetscCall(DMFieldEvaluate(field, points, PETSC_SCALAR, B, D, H)); 499566063dSJacob Faibussowitsch PetscCall(DMFieldEvaluate(field, points, PETSC_REAL, rB, rD, rH)); 50c4762a1bSJed Brown viewer = PETSC_VIEWER_STDOUT_(comm); 51c4762a1bSJed Brown 529566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)points, "Test Points")); 539566063dSJacob Faibussowitsch PetscCall(VecView(points, viewer)); 549566063dSJacob Faibussowitsch PetscCall(ViewResults(viewer, n * nc, dim, B, D, H, rB, rD, rH)); 55c4762a1bSJed Brown 569566063dSJacob Faibussowitsch PetscCall(PetscFree6(B, rB, D, rD, H, rH)); 579566063dSJacob Faibussowitsch PetscCall(VecDestroy(&points)); 583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 59c4762a1bSJed Brown } 60c4762a1bSJed Brown 61d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestEvaluateFE(DMField field, PetscInt n, PetscInt cStart, PetscInt cEnd, PetscQuadrature quad, PetscRandom rand) 62d71ae5a4SJacob Faibussowitsch { 63c4762a1bSJed Brown DM dm; 64c4762a1bSJed Brown PetscInt dim, i, nc, nq; 65c4762a1bSJed Brown PetscInt N; 66c4762a1bSJed Brown PetscScalar *B, *D, *H; 67c4762a1bSJed Brown PetscReal *rB, *rD, *rH; 68c4762a1bSJed Brown PetscInt *cells; 69c4762a1bSJed Brown IS cellIS; 70c4762a1bSJed Brown PetscViewer viewer; 71c4762a1bSJed Brown MPI_Comm comm; 72c4762a1bSJed Brown 73c4762a1bSJed Brown PetscFunctionBegin; 74c4762a1bSJed Brown comm = PetscObjectComm((PetscObject)field); 759566063dSJacob Faibussowitsch PetscCall(DMFieldGetNumComponents(field, &nc)); 769566063dSJacob Faibussowitsch PetscCall(DMFieldGetDM(field, &dm)); 779566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 789566063dSJacob Faibussowitsch PetscCall(PetscRandomSetInterval(rand, (PetscScalar)cStart, (PetscScalar)cEnd)); 799566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &cells)); 80c4762a1bSJed Brown for (i = 0; i < n; i++) { 81c4762a1bSJed Brown PetscReal rc; 82c4762a1bSJed Brown 839566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValueReal(rand, &rc)); 84c4762a1bSJed Brown cells[i] = PetscFloorReal(rc); 85c4762a1bSJed Brown } 869566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, n, cells, PETSC_OWN_POINTER, &cellIS)); 879566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)cellIS, "FE Test Cells")); 889566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(quad, NULL, NULL, &nq, NULL, NULL)); 89c4762a1bSJed Brown N = n * nq * nc; 909566063dSJacob Faibussowitsch PetscCall(PetscMalloc6(N, &B, N, &rB, N * dim, &D, N * dim, &rD, N * dim * dim, &H, N * dim * dim, &rH)); 919566063dSJacob Faibussowitsch PetscCall(DMFieldEvaluateFE(field, cellIS, quad, PETSC_SCALAR, B, D, H)); 929566063dSJacob Faibussowitsch PetscCall(DMFieldEvaluateFE(field, cellIS, quad, PETSC_REAL, rB, rD, rH)); 93c4762a1bSJed Brown viewer = PETSC_VIEWER_STDOUT_(comm); 94c4762a1bSJed Brown 959566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)quad, "Test quadrature")); 969566063dSJacob Faibussowitsch PetscCall(PetscQuadratureView(quad, viewer)); 979566063dSJacob Faibussowitsch PetscCall(ISView(cellIS, viewer)); 989566063dSJacob Faibussowitsch PetscCall(ViewResults(viewer, N, dim, B, D, H, rB, rD, rH)); 99c4762a1bSJed Brown 1009566063dSJacob Faibussowitsch PetscCall(PetscFree6(B, rB, D, rD, H, rH)); 1019566063dSJacob Faibussowitsch PetscCall(ISDestroy(&cellIS)); 1023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 103c4762a1bSJed Brown } 104c4762a1bSJed Brown 105d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestEvaluateFV(DMField field, PetscInt n, PetscInt cStart, PetscInt cEnd, PetscRandom rand) 106d71ae5a4SJacob Faibussowitsch { 107c4762a1bSJed Brown DM dm; 108c4762a1bSJed Brown PetscInt dim, i, nc; 109c4762a1bSJed Brown PetscInt N; 110c4762a1bSJed Brown PetscScalar *B, *D, *H; 111c4762a1bSJed Brown PetscReal *rB, *rD, *rH; 112c4762a1bSJed Brown PetscInt *cells; 113c4762a1bSJed Brown IS cellIS; 114c4762a1bSJed Brown PetscViewer viewer; 115c4762a1bSJed Brown MPI_Comm comm; 116c4762a1bSJed Brown 117c4762a1bSJed Brown PetscFunctionBegin; 118c4762a1bSJed Brown comm = PetscObjectComm((PetscObject)field); 1199566063dSJacob Faibussowitsch PetscCall(DMFieldGetNumComponents(field, &nc)); 1209566063dSJacob Faibussowitsch PetscCall(DMFieldGetDM(field, &dm)); 1219566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 1229566063dSJacob Faibussowitsch PetscCall(PetscRandomSetInterval(rand, (PetscScalar)cStart, (PetscScalar)cEnd)); 1239566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &cells)); 124c4762a1bSJed Brown for (i = 0; i < n; i++) { 125c4762a1bSJed Brown PetscReal rc; 126c4762a1bSJed Brown 1279566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValueReal(rand, &rc)); 128c4762a1bSJed Brown cells[i] = PetscFloorReal(rc); 129c4762a1bSJed Brown } 1309566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, n, cells, PETSC_OWN_POINTER, &cellIS)); 1319566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)cellIS, "FV Test Cells")); 132c4762a1bSJed Brown N = n * nc; 1339566063dSJacob Faibussowitsch PetscCall(PetscMalloc6(N, &B, N, &rB, N * dim, &D, N * dim, &rD, N * dim * dim, &H, N * dim * dim, &rH)); 1349566063dSJacob Faibussowitsch PetscCall(DMFieldEvaluateFV(field, cellIS, PETSC_SCALAR, B, D, H)); 1359566063dSJacob Faibussowitsch PetscCall(DMFieldEvaluateFV(field, cellIS, PETSC_REAL, rB, rD, rH)); 136c4762a1bSJed Brown viewer = PETSC_VIEWER_STDOUT_(comm); 137c4762a1bSJed Brown 1389566063dSJacob Faibussowitsch PetscCall(ISView(cellIS, viewer)); 1399566063dSJacob Faibussowitsch PetscCall(ViewResults(viewer, N, dim, B, D, H, rB, rD, rH)); 140c4762a1bSJed Brown 1419566063dSJacob Faibussowitsch PetscCall(PetscFree6(B, rB, D, rD, H, rH)); 1429566063dSJacob Faibussowitsch PetscCall(ISDestroy(&cellIS)); 1433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 144c4762a1bSJed Brown } 145c4762a1bSJed Brown 146d71ae5a4SJacob Faibussowitsch static PetscErrorCode radiusSquared(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar u[], void *ctx) 147d71ae5a4SJacob Faibussowitsch { 148c4762a1bSJed Brown PetscInt i; 149c4762a1bSJed Brown PetscReal r2 = 0.; 150c4762a1bSJed Brown 151c4762a1bSJed Brown PetscFunctionBegin; 152ad540459SPierre Jolivet for (i = 0; i < dim; i++) r2 += PetscSqr(x[i]); 153ad540459SPierre Jolivet for (i = 0; i < Nf; i++) u[i] = (i + 1) * r2; 1543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 155c4762a1bSJed Brown } 156c4762a1bSJed Brown 157d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestShellEvaluate(DMField field, Vec points, PetscDataType type, void *B, void *D, void *H) 158d71ae5a4SJacob Faibussowitsch { 159c4762a1bSJed Brown Vec ctxVec = NULL; 160c4762a1bSJed Brown const PetscScalar *mult; 161c4762a1bSJed Brown PetscInt dim; 162c4762a1bSJed Brown const PetscScalar *x; 163c4762a1bSJed Brown PetscInt Nc, n, i, j, k, l; 164c4762a1bSJed Brown 165c4762a1bSJed Brown PetscFunctionBegin; 1669566063dSJacob Faibussowitsch PetscCall(DMFieldGetNumComponents(field, &Nc)); 1679566063dSJacob Faibussowitsch PetscCall(DMFieldShellGetContext(field, &ctxVec)); 1689566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(points, &dim)); 1699566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(points, &n)); 170c4762a1bSJed Brown n /= Nc; 1719566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(ctxVec, &mult)); 1729566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(points, &x)); 173c4762a1bSJed Brown for (i = 0; i < n; i++) { 174c4762a1bSJed Brown PetscReal r2 = 0.; 175c4762a1bSJed Brown 176ad540459SPierre Jolivet for (j = 0; j < dim; j++) r2 += PetscSqr(PetscRealPart(x[i * dim + j])); 177c4762a1bSJed Brown for (j = 0; j < Nc; j++) { 178c4762a1bSJed Brown PetscReal m = PetscRealPart(mult[j]); 179c4762a1bSJed Brown if (B) { 180c4762a1bSJed Brown if (type == PETSC_SCALAR) { 181c4762a1bSJed Brown ((PetscScalar *)B)[i * Nc + j] = m * r2; 182c4762a1bSJed Brown } else { 183c4762a1bSJed Brown ((PetscReal *)B)[i * Nc + j] = m * r2; 184c4762a1bSJed Brown } 185c4762a1bSJed Brown } 186c4762a1bSJed Brown if (D) { 187c4762a1bSJed Brown if (type == PETSC_SCALAR) { 188c4762a1bSJed Brown for (k = 0; k < dim; k++) ((PetscScalar *)D)[(i * Nc + j) * dim + k] = 2. * m * x[i * dim + k]; 189c4762a1bSJed Brown } else { 190c4762a1bSJed Brown for (k = 0; k < dim; k++) ((PetscReal *)D)[(i * Nc + j) * dim + k] = 2. * m * PetscRealPart(x[i * dim + k]); 191c4762a1bSJed Brown } 192c4762a1bSJed Brown } 193c4762a1bSJed Brown if (H) { 194c4762a1bSJed Brown if (type == PETSC_SCALAR) { 1959371c9d4SSatish Balay for (k = 0; k < dim; k++) 1969371c9d4SSatish Balay for (l = 0; l < dim; l++) ((PetscScalar *)H)[((i * Nc + j) * dim + k) * dim + l] = (k == l) ? 2. * m : 0.; 197c4762a1bSJed Brown } else { 1989371c9d4SSatish Balay for (k = 0; k < dim; k++) 1999371c9d4SSatish Balay for (l = 0; l < dim; l++) ((PetscReal *)H)[((i * Nc + j) * dim + k) * dim + l] = (k == l) ? 2. * m : 0.; 200c4762a1bSJed Brown } 201c4762a1bSJed Brown } 202c4762a1bSJed Brown } 203c4762a1bSJed Brown } 2049566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(points, &x)); 2059566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(ctxVec, &mult)); 2063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 207c4762a1bSJed Brown } 208c4762a1bSJed Brown 209d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestShellDestroy(DMField field) 210d71ae5a4SJacob Faibussowitsch { 211c4762a1bSJed Brown Vec ctxVec = NULL; 212c4762a1bSJed Brown 213c4762a1bSJed Brown PetscFunctionBegin; 2149566063dSJacob Faibussowitsch PetscCall(DMFieldShellGetContext(field, &ctxVec)); 2159566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctxVec)); 2163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 217c4762a1bSJed Brown } 218c4762a1bSJed Brown 219d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 220d71ae5a4SJacob Faibussowitsch { 221c4762a1bSJed Brown DM dm = NULL; 222c4762a1bSJed Brown MPI_Comm comm; 223c4762a1bSJed Brown char type[256] = DMPLEX; 224c4762a1bSJed Brown PetscBool isda, isplex; 225c4762a1bSJed Brown PetscInt dim = 2; 226c4762a1bSJed Brown DMField field = NULL; 227c4762a1bSJed Brown PetscInt nc = 1; 228c4762a1bSJed Brown PetscInt cStart = -1, cEnd = -1; 229c4762a1bSJed Brown PetscRandom rand; 230c4762a1bSJed Brown PetscQuadrature quad = NULL; 231c4762a1bSJed Brown PetscInt pointsPerEdge = 2; 232c4762a1bSJed Brown PetscInt numPoint = 0, numFE = 0, numFV = 0; 233c4762a1bSJed Brown PetscBool testShell = PETSC_FALSE; 234c4762a1bSJed Brown 235327415f7SBarry Smith PetscFunctionBeginUser; 2369566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 237c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 238d0609cedSBarry Smith PetscOptionsBegin(comm, "", "DMField Tutorial Options", "DM"); 2399566063dSJacob Faibussowitsch PetscCall(PetscOptionsFList("-dm_type", "DM implementation on which to define field", "ex1.c", DMList, type, type, 256, NULL)); 2409566063dSJacob Faibussowitsch PetscCall(PetscOptionsRangeInt("-dim", "DM intrinsic dimension", "ex1.c", dim, &dim, NULL, 1, 3)); 2419566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-num_components", "Number of components in field", "ex1.c", nc, &nc, NULL, 1)); 2429566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-num_quad_points", "Number of quadrature points per dimension", "ex1.c", pointsPerEdge, &pointsPerEdge, NULL, 1)); 2439566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-num_point_tests", "Number of test points for DMFieldEvaluate()", "ex1.c", numPoint, &numPoint, NULL, 0)); 2449566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-num_fe_tests", "Number of test cells for DMFieldEvaluateFE()", "ex1.c", numFE, &numFE, NULL, 0)); 2459566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-num_fv_tests", "Number of test cells for DMFieldEvaluateFV()", "ex1.c", numFV, &numFV, NULL, 0)); 2469566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-test_shell", "Test the DMFIELDSHELL implementation of DMField", "ex1.c", testShell, &testShell, NULL)); 247d0609cedSBarry Smith PetscOptionsEnd(); 248c4762a1bSJed Brown 24963a3b9bcSJacob Faibussowitsch PetscCheck(dim <= 3, comm, PETSC_ERR_ARG_OUTOFRANGE, "This examples works for dim <= 3, not %" PetscInt_FMT, dim); 2509566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(type, DMPLEX, 256, &isplex)); 2519566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(type, DMDA, 256, &isda)); 252c4762a1bSJed Brown 2539566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rand)); 2549566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rand)); 255c4762a1bSJed Brown if (isplex) { 256c4762a1bSJed Brown PetscInt overlap = 0; 257c4762a1bSJed Brown Vec fieldvec; 258c4762a1bSJed Brown PetscInt cells[3] = {3, 3, 3}; 25930602db0SMatthew G. Knepley PetscBool simplex; 260c4762a1bSJed Brown PetscFE fe; 261c4762a1bSJed Brown 262d0609cedSBarry Smith PetscOptionsBegin(comm, "", "DMField DMPlex Options", "DM"); 2639566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-overlap", "DMPlex parallel overlap", "ex1.c", overlap, &overlap, NULL, 0)); 264d0609cedSBarry Smith PetscOptionsEnd(); 26530602db0SMatthew G. Knepley if (0) { 266*42108689Sksagiyam PetscCall(DMPlexCreateBoxMesh(comm, 2, PETSC_TRUE, cells, NULL, NULL, NULL, PETSC_TRUE, 0, PETSC_TRUE, &dm)); 26730602db0SMatthew G. Knepley } else { 2689566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, &dm)); 2699566063dSJacob Faibussowitsch PetscCall(DMSetType(dm, DMPLEX)); 2709566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dm)); 27130602db0SMatthew G. Knepley CHKMEMQ; 272c4762a1bSJed Brown } 2739566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 2749566063dSJacob Faibussowitsch PetscCall(DMPlexIsSimplex(dm, &simplex)); 2751baa6e33SBarry Smith if (simplex) PetscCall(PetscDTStroudConicalQuadrature(dim, 1, pointsPerEdge, -1.0, 1.0, &quad)); 2761baa6e33SBarry Smith else PetscCall(PetscDTGaussTensorQuadrature(dim, 1, pointsPerEdge, -1.0, 1.0, &quad)); 2779566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 278c4762a1bSJed Brown if (testShell) { 279c4762a1bSJed Brown Vec ctxVec; 280c4762a1bSJed Brown PetscInt i; 281c4762a1bSJed Brown PetscScalar *array; 282c4762a1bSJed Brown 2839566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, nc, &ctxVec)); 2849566063dSJacob Faibussowitsch PetscCall(VecSetUp(ctxVec)); 2859566063dSJacob Faibussowitsch PetscCall(VecGetArray(ctxVec, &array)); 286c4762a1bSJed Brown for (i = 0; i < nc; i++) array[i] = i + 1.; 2879566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(ctxVec, &array)); 2889566063dSJacob Faibussowitsch PetscCall(DMFieldCreateShell(dm, nc, DMFIELD_VERTEX, (void *)ctxVec, &field)); 2899566063dSJacob Faibussowitsch PetscCall(DMFieldShellSetEvaluate(field, TestShellEvaluate)); 2909566063dSJacob Faibussowitsch PetscCall(DMFieldShellSetDestroy(field, TestShellDestroy)); 291c4762a1bSJed Brown } else { 2929566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, nc, simplex, NULL, PETSC_DEFAULT, &fe)); 2939566063dSJacob Faibussowitsch PetscCall(PetscFESetName(fe, "MyPetscFE")); 2949566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe)); 2959566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 2969566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 2979566063dSJacob Faibussowitsch PetscCall(DMCreateLocalVector(dm, &fieldvec)); 298c4762a1bSJed Brown { 299c4762a1bSJed Brown PetscErrorCode (*func[1])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *); 300c4762a1bSJed Brown void *ctxs[1]; 301c4762a1bSJed Brown 302c4762a1bSJed Brown func[0] = radiusSquared; 303c4762a1bSJed Brown ctxs[0] = NULL; 304c4762a1bSJed Brown 3059566063dSJacob Faibussowitsch PetscCall(DMProjectFunctionLocal(dm, 0.0, func, ctxs, INSERT_ALL_VALUES, fieldvec)); 306c4762a1bSJed Brown } 3079566063dSJacob Faibussowitsch PetscCall(DMFieldCreateDS(dm, 0, fieldvec, &field)); 3089566063dSJacob Faibussowitsch PetscCall(VecDestroy(&fieldvec)); 309c4762a1bSJed Brown } 310c4762a1bSJed Brown } else if (isda) { 311c4762a1bSJed Brown PetscInt i; 312c4762a1bSJed Brown PetscScalar *cv; 313c4762a1bSJed Brown 314c4762a1bSJed Brown switch (dim) { 315d71ae5a4SJacob Faibussowitsch case 1: 316d71ae5a4SJacob Faibussowitsch PetscCall(DMDACreate1d(comm, DM_BOUNDARY_NONE, 3, 1, 1, NULL, &dm)); 317d71ae5a4SJacob Faibussowitsch break; 318d71ae5a4SJacob Faibussowitsch case 2: 319d71ae5a4SJacob Faibussowitsch PetscCall(DMDACreate2d(comm, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, 3, 3, PETSC_DETERMINE, PETSC_DETERMINE, 1, 1, NULL, NULL, &dm)); 320d71ae5a4SJacob Faibussowitsch break; 321d71ae5a4SJacob Faibussowitsch default: 322d71ae5a4SJacob Faibussowitsch PetscCall(DMDACreate3d(comm, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, 3, 3, 3, PETSC_DETERMINE, PETSC_DETERMINE, PETSC_DETERMINE, 1, 1, NULL, NULL, NULL, &dm)); 323d71ae5a4SJacob Faibussowitsch break; 324c4762a1bSJed Brown } 3259566063dSJacob Faibussowitsch PetscCall(DMSetUp(dm)); 3269566063dSJacob Faibussowitsch PetscCall(DMDAGetHeightStratum(dm, 0, &cStart, &cEnd)); 3279566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nc * (1 << dim), &cv)); 328c4762a1bSJed Brown for (i = 0; i < nc * (1 << dim); i++) { 329c4762a1bSJed Brown PetscReal rv; 330c4762a1bSJed Brown 3319566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValueReal(rand, &rv)); 332c4762a1bSJed Brown cv[i] = rv; 333c4762a1bSJed Brown } 3349566063dSJacob Faibussowitsch PetscCall(DMFieldCreateDA(dm, nc, cv, &field)); 3359566063dSJacob Faibussowitsch PetscCall(PetscFree(cv)); 3369566063dSJacob Faibussowitsch PetscCall(PetscDTGaussTensorQuadrature(dim, 1, pointsPerEdge, -1.0, 1.0, &quad)); 33798921bdaSJacob Faibussowitsch } else SETERRQ(comm, PETSC_ERR_SUP, "This test does not run for DM type %s", type); 338c4762a1bSJed Brown 3399566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)dm, "mesh")); 3409566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 3419566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 3429566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)field, "field")); 3439566063dSJacob Faibussowitsch PetscCall(PetscObjectViewFromOptions((PetscObject)field, NULL, "-dmfield_view")); 3449566063dSJacob Faibussowitsch if (numPoint) PetscCall(TestEvaluate(field, numPoint, rand)); 3459566063dSJacob Faibussowitsch if (numFE) PetscCall(TestEvaluateFE(field, numFE, cStart, cEnd, quad, rand)); 3469566063dSJacob Faibussowitsch if (numFV) PetscCall(TestEvaluateFV(field, numFV, cStart, cEnd, rand)); 3479566063dSJacob Faibussowitsch PetscCall(DMFieldDestroy(&field)); 3489566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(&quad)); 3499566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rand)); 3509566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 351b122ec5aSJacob Faibussowitsch return 0; 352c4762a1bSJed Brown } 353c4762a1bSJed Brown 354c4762a1bSJed Brown /*TEST 355c4762a1bSJed Brown 356c4762a1bSJed Brown test: 357c4762a1bSJed Brown suffix: da 358c4762a1bSJed Brown requires: !complex 359c4762a1bSJed Brown args: -dm_type da -dim 2 -num_components 2 -num_point_tests 2 -num_fe_tests 2 -num_fv_tests 2 -dmfield_view 360c4762a1bSJed Brown 361c4762a1bSJed Brown test: 362c4762a1bSJed Brown suffix: da_1 363c4762a1bSJed Brown requires: !complex 364c4762a1bSJed Brown args: -dm_type da -dim 1 -num_fe_tests 2 365c4762a1bSJed Brown 366c4762a1bSJed Brown test: 367c4762a1bSJed Brown suffix: da_2 368c4762a1bSJed Brown requires: !complex 369c4762a1bSJed Brown args: -dm_type da -dim 2 -num_fe_tests 2 370c4762a1bSJed Brown 371c4762a1bSJed Brown test: 372c4762a1bSJed Brown suffix: da_3 373c4762a1bSJed Brown requires: !complex 374c4762a1bSJed Brown args: -dm_type da -dim 3 -num_fe_tests 2 375c4762a1bSJed Brown 376c4762a1bSJed Brown test: 377c4762a1bSJed Brown suffix: ds 378c4762a1bSJed Brown requires: !complex triangle 37930602db0SMatthew G. Knepley args: -dm_type plex -dm_plex_dim 2 -dm_plex_box_faces 3,3 -num_components 2 -num_point_tests 2 -num_fe_tests 2 -num_fv_tests 2 -dmfield_view -petscspace_degree 2 -num_quad_points 1 380c4762a1bSJed Brown 381c4762a1bSJed Brown test: 382c4762a1bSJed Brown suffix: ds_simplex_0 383c4762a1bSJed Brown requires: !complex triangle 38430602db0SMatthew G. Knepley args: -dm_type plex -dm_plex_dim 2 -dm_plex_box_faces 3,3 -num_fe_tests 2 -petscspace_degree 0 385c4762a1bSJed Brown 386c4762a1bSJed Brown test: 387c4762a1bSJed Brown suffix: ds_simplex_1 388c4762a1bSJed Brown requires: !complex triangle 38930602db0SMatthew G. Knepley args: -dm_type plex -dm_plex_dim 2 -dm_plex_box_faces 3,3 -num_fe_tests 2 -petscspace_degree 1 390c4762a1bSJed Brown 391c4762a1bSJed Brown test: 392c4762a1bSJed Brown suffix: ds_simplex_2 393c4762a1bSJed Brown requires: !complex triangle 39430602db0SMatthew G. Knepley args: -dm_type plex -dm_plex_dim 2 -dm_plex_box_faces 3,3 -num_fe_tests 2 -petscspace_degree 2 395c4762a1bSJed Brown 396c4762a1bSJed Brown test: 397c4762a1bSJed Brown suffix: ds_tensor_2_0 398c4762a1bSJed Brown requires: !complex 39930602db0SMatthew G. Knepley args: -dm_type plex -dm_plex_dim 2 -dm_plex_box_faces 3,3 -num_fe_tests 2 -petscspace_poly_tensor 1 -petscspace_degree 0 -dm_plex_simplex 0 400c4762a1bSJed Brown 401c4762a1bSJed Brown test: 402c4762a1bSJed Brown suffix: ds_tensor_2_1 403c4762a1bSJed Brown requires: !complex 40430602db0SMatthew G. Knepley args: -dm_type plex -dm_plex_dim 2 -dm_plex_box_faces 3,3 -num_fe_tests 2 -petscspace_poly_tensor 1 -petscspace_degree 1 -dm_plex_simplex 0 405c4762a1bSJed Brown 406c4762a1bSJed Brown test: 407c4762a1bSJed Brown suffix: ds_tensor_2_2 408c4762a1bSJed Brown requires: !complex 40930602db0SMatthew G. Knepley args: -dm_type plex -dm_plex_dim 2 -dm_plex_box_faces 3,3 -num_fe_tests 2 -petscspace_poly_tensor 1 -petscspace_degree 2 -dm_plex_simplex 0 410c4762a1bSJed Brown 411c4762a1bSJed Brown test: 412c4762a1bSJed Brown suffix: ds_tensor_3_0 413c4762a1bSJed Brown requires: !complex 41430602db0SMatthew G. Knepley args: -dm_type plex -dm_plex_dim 3 -dm_plex_box_faces 3,3,3 -num_fe_tests 2 -petscspace_poly_tensor 1 -petscspace_degree 0 -dm_plex_simplex 0 415c4762a1bSJed Brown 416c4762a1bSJed Brown test: 417c4762a1bSJed Brown suffix: ds_tensor_3_1 418c4762a1bSJed Brown requires: !complex 41930602db0SMatthew G. Knepley args: -dm_type plex -dm_plex_dim 3 -dm_plex_box_faces 3,3,3 -num_fe_tests 2 -petscspace_poly_tensor 1 -petscspace_degree 1 -dm_plex_simplex 0 420c4762a1bSJed Brown 421c4762a1bSJed Brown test: 422c4762a1bSJed Brown suffix: ds_tensor_3_2 423c4762a1bSJed Brown requires: !complex 42430602db0SMatthew G. Knepley args: -dm_type plex -dm_plex_dim 3 -dm_plex_box_faces 3,3,3 -num_fe_tests 2 -petscspace_poly_tensor 1 -petscspace_degree 2 -dm_plex_simplex 0 425c4762a1bSJed Brown 426c4762a1bSJed Brown test: 427c4762a1bSJed Brown suffix: shell 428c4762a1bSJed Brown requires: !complex triangle 42930602db0SMatthew G. Knepley args: -dm_coord_space 0 -dm_type plex -dm_plex_dim 2 -dm_plex_box_faces 3,3 -num_components 2 -num_point_tests 2 -num_fe_tests 2 -num_fv_tests 2 -dmfield_view -num_quad_points 1 -test_shell 430c4762a1bSJed Brown 431c4762a1bSJed Brown TEST*/ 432