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 7c4762a1bSJed Brown static PetscErrorCode ViewResults(PetscViewer viewer, PetscInt N, PetscInt dim, PetscScalar *B, PetscScalar *D, PetscScalar *H, PetscReal *rB, PetscReal *rD, PetscReal *rH) 8c4762a1bSJed Brown { 9c4762a1bSJed Brown PetscFunctionBegin; 105f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer,"B:\n")); 115f80ce2aSJacob Faibussowitsch CHKERRQ(PetscScalarView(N,B,viewer)); 125f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer,"D:\n")); 135f80ce2aSJacob Faibussowitsch CHKERRQ(PetscScalarView(N*dim,D,viewer)); 145f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer,"H:\n")); 155f80ce2aSJacob Faibussowitsch CHKERRQ(PetscScalarView(N*dim*dim,H,viewer)); 16c4762a1bSJed Brown 175f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer,"rB:\n")); 185f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRealView(N,rB,viewer)); 195f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer,"rD:\n")); 205f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRealView(N*dim,rD,viewer)); 215f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer,"rH:\n")); 225f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRealView(N*dim*dim,rH,viewer)); 23c4762a1bSJed Brown PetscFunctionReturn(0); 24c4762a1bSJed Brown } 25c4762a1bSJed Brown 26c4762a1bSJed Brown static PetscErrorCode TestEvaluate(DMField field, PetscInt n, PetscRandom rand) 27c4762a1bSJed Brown { 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); 395f80ce2aSJacob Faibussowitsch CHKERRQ(DMFieldGetNumComponents(field,&nc)); 405f80ce2aSJacob Faibussowitsch CHKERRQ(DMFieldGetDM(field,&dm)); 415f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(dm,&dim)); 425f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateMPI(PetscObjectComm((PetscObject)field),n * dim,PETSC_DETERMINE,&points)); 435f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetBlockSize(points,dim)); 445f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(points,&array)); 455f80ce2aSJacob Faibussowitsch for (i = 0; i < n * dim; i++) CHKERRQ(PetscRandomGetValue(rand,&array[i])); 465f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(points,&array)); 475f80ce2aSJacob Faibussowitsch CHKERRQ(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)); 485f80ce2aSJacob Faibussowitsch CHKERRQ(DMFieldEvaluate(field,points,PETSC_SCALAR,B,D,H)); 495f80ce2aSJacob Faibussowitsch CHKERRQ(DMFieldEvaluate(field,points,PETSC_REAL,rB,rD,rH)); 50c4762a1bSJed Brown viewer = PETSC_VIEWER_STDOUT_(comm); 51c4762a1bSJed Brown 525f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)points,"Test Points")); 535f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(points,viewer)); 545f80ce2aSJacob Faibussowitsch CHKERRQ(ViewResults(viewer,n*nc,dim,B,D,H,rB,rD,rH)); 55c4762a1bSJed Brown 565f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree6(B,rB,D,rD,H,rH)); 575f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&points)); 58c4762a1bSJed Brown PetscFunctionReturn(0); 59c4762a1bSJed Brown } 60c4762a1bSJed Brown 61c4762a1bSJed Brown static PetscErrorCode TestEvaluateFE(DMField field, PetscInt n, PetscInt cStart, PetscInt cEnd, PetscQuadrature quad, PetscRandom rand) 62c4762a1bSJed Brown { 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); 755f80ce2aSJacob Faibussowitsch CHKERRQ(DMFieldGetNumComponents(field,&nc)); 765f80ce2aSJacob Faibussowitsch CHKERRQ(DMFieldGetDM(field,&dm)); 775f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(dm,&dim)); 785f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetInterval(rand,(PetscScalar) cStart, (PetscScalar) cEnd)); 795f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(n,&cells)); 80c4762a1bSJed Brown for (i = 0; i < n; i++) { 81c4762a1bSJed Brown PetscReal rc; 82c4762a1bSJed Brown 835f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomGetValueReal(rand,&rc)); 84c4762a1bSJed Brown cells[i] = PetscFloorReal(rc); 85c4762a1bSJed Brown } 865f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateGeneral(comm,n,cells,PETSC_OWN_POINTER,&cellIS)); 875f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)cellIS,"FE Test Cells")); 885f80ce2aSJacob Faibussowitsch CHKERRQ(PetscQuadratureGetData(quad,NULL,NULL,&nq,NULL,NULL)); 89c4762a1bSJed Brown N = n * nq * nc; 905f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc6(N,&B,N,&rB,N*dim,&D,N*dim,&rD,N*dim*dim,&H,N*dim*dim,&rH)); 915f80ce2aSJacob Faibussowitsch CHKERRQ(DMFieldEvaluateFE(field,cellIS,quad,PETSC_SCALAR,B,D,H)); 925f80ce2aSJacob Faibussowitsch CHKERRQ(DMFieldEvaluateFE(field,cellIS,quad,PETSC_REAL,rB,rD,rH)); 93c4762a1bSJed Brown viewer = PETSC_VIEWER_STDOUT_(comm); 94c4762a1bSJed Brown 955f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)quad,"Test quadrature")); 965f80ce2aSJacob Faibussowitsch CHKERRQ(PetscQuadratureView(quad,viewer)); 975f80ce2aSJacob Faibussowitsch CHKERRQ(ISView(cellIS,viewer)); 985f80ce2aSJacob Faibussowitsch CHKERRQ(ViewResults(viewer,N,dim,B,D,H,rB,rD,rH)); 99c4762a1bSJed Brown 1005f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree6(B,rB,D,rD,H,rH)); 1015f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&cellIS)); 102c4762a1bSJed Brown PetscFunctionReturn(0); 103c4762a1bSJed Brown } 104c4762a1bSJed Brown 105c4762a1bSJed Brown static PetscErrorCode TestEvaluateFV(DMField field, PetscInt n, PetscInt cStart, PetscInt cEnd, PetscRandom rand) 106c4762a1bSJed Brown { 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); 1195f80ce2aSJacob Faibussowitsch CHKERRQ(DMFieldGetNumComponents(field,&nc)); 1205f80ce2aSJacob Faibussowitsch CHKERRQ(DMFieldGetDM(field,&dm)); 1215f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(dm,&dim)); 1225f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetInterval(rand,(PetscScalar) cStart, (PetscScalar) cEnd)); 1235f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(n,&cells)); 124c4762a1bSJed Brown for (i = 0; i < n; i++) { 125c4762a1bSJed Brown PetscReal rc; 126c4762a1bSJed Brown 1275f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomGetValueReal(rand,&rc)); 128c4762a1bSJed Brown cells[i] = PetscFloorReal(rc); 129c4762a1bSJed Brown } 1305f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateGeneral(comm,n,cells,PETSC_OWN_POINTER,&cellIS)); 1315f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)cellIS,"FV Test Cells")); 132c4762a1bSJed Brown N = n * nc; 1335f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc6(N,&B,N,&rB,N*dim,&D,N*dim,&rD,N*dim*dim,&H,N*dim*dim,&rH)); 1345f80ce2aSJacob Faibussowitsch CHKERRQ(DMFieldEvaluateFV(field,cellIS,PETSC_SCALAR,B,D,H)); 1355f80ce2aSJacob Faibussowitsch CHKERRQ(DMFieldEvaluateFV(field,cellIS,PETSC_REAL,rB,rD,rH)); 136c4762a1bSJed Brown viewer = PETSC_VIEWER_STDOUT_(comm); 137c4762a1bSJed Brown 1385f80ce2aSJacob Faibussowitsch CHKERRQ(ISView(cellIS,viewer)); 1395f80ce2aSJacob Faibussowitsch CHKERRQ(ViewResults(viewer,N,dim,B,D,H,rB,rD,rH)); 140c4762a1bSJed Brown 1415f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree6(B,rB,D,rD,H,rH)); 1425f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&cellIS)); 143c4762a1bSJed Brown PetscFunctionReturn(0); 144c4762a1bSJed Brown } 145c4762a1bSJed Brown 146c4762a1bSJed Brown static PetscErrorCode radiusSquared(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar u[], void *ctx) 147c4762a1bSJed Brown { 148c4762a1bSJed Brown PetscInt i; 149c4762a1bSJed Brown PetscReal r2 = 0.; 150c4762a1bSJed Brown 151c4762a1bSJed Brown PetscFunctionBegin; 152c4762a1bSJed Brown for (i = 0; i < dim; i++) {r2 += PetscSqr(x[i]);} 153c4762a1bSJed Brown for (i = 0; i < Nf; i++) { 154c4762a1bSJed Brown u[i] = (i + 1) * r2; 155c4762a1bSJed Brown } 156c4762a1bSJed Brown PetscFunctionReturn(0); 157c4762a1bSJed Brown } 158c4762a1bSJed Brown 159c4762a1bSJed Brown static PetscErrorCode TestShellEvaluate(DMField field, Vec points, PetscDataType type, void *B, void *D, void *H) 160c4762a1bSJed Brown { 161c4762a1bSJed Brown Vec ctxVec = NULL; 162c4762a1bSJed Brown const PetscScalar *mult; 163c4762a1bSJed Brown PetscInt dim; 164c4762a1bSJed Brown const PetscScalar *x; 165c4762a1bSJed Brown PetscInt Nc, n, i, j, k, l; 166c4762a1bSJed Brown 167c4762a1bSJed Brown PetscFunctionBegin; 1685f80ce2aSJacob Faibussowitsch CHKERRQ(DMFieldGetNumComponents(field, &Nc)); 1695f80ce2aSJacob Faibussowitsch CHKERRQ(DMFieldShellGetContext(field, &ctxVec)); 1705f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetBlockSize(points, &dim)); 1715f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetLocalSize(points, &n)); 172c4762a1bSJed Brown n /= Nc; 1735f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(ctxVec, &mult)); 1745f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(points, &x)); 175c4762a1bSJed Brown for (i = 0; i < n; i++) { 176c4762a1bSJed Brown PetscReal r2 = 0.; 177c4762a1bSJed Brown 178c4762a1bSJed Brown for (j = 0; j < dim; j++) {r2 += PetscSqr(PetscRealPart(x[i * dim + j]));} 179c4762a1bSJed Brown for (j = 0; j < Nc; j++) { 180c4762a1bSJed Brown PetscReal m = PetscRealPart(mult[j]); 181c4762a1bSJed Brown if (B) { 182c4762a1bSJed Brown if (type == PETSC_SCALAR) { 183c4762a1bSJed Brown ((PetscScalar *)B)[i * Nc + j] = m * r2; 184c4762a1bSJed Brown } else { 185c4762a1bSJed Brown ((PetscReal *)B)[i * Nc + j] = m * r2; 186c4762a1bSJed Brown } 187c4762a1bSJed Brown } 188c4762a1bSJed Brown if (D) { 189c4762a1bSJed Brown if (type == PETSC_SCALAR) { 190c4762a1bSJed Brown for (k = 0; k < dim; k++) ((PetscScalar *)D)[(i * Nc + j) * dim + k] = 2. * m * x[i * dim + k]; 191c4762a1bSJed Brown } else { 192c4762a1bSJed Brown for (k = 0; k < dim; k++) ((PetscReal *)D)[(i * Nc + j) * dim + k] = 2. * m * PetscRealPart(x[i * dim + k]); 193c4762a1bSJed Brown } 194c4762a1bSJed Brown } 195c4762a1bSJed Brown if (H) { 196c4762a1bSJed Brown if (type == PETSC_SCALAR) { 197c4762a1bSJed Brown for (k = 0; k < dim; k++) for (l = 0; l < dim; l++) ((PetscScalar *)H)[((i * Nc + j) * dim + k) * dim + l] = (k == l) ? 2. * m : 0.; 198c4762a1bSJed Brown } else { 199c4762a1bSJed Brown for (k = 0; k < dim; k++) 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 } 2045f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(points, &x)); 2055f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(ctxVec, &mult)); 206c4762a1bSJed Brown PetscFunctionReturn(0); 207c4762a1bSJed Brown } 208c4762a1bSJed Brown 209c4762a1bSJed Brown static PetscErrorCode TestShellDestroy(DMField field) 210c4762a1bSJed Brown { 211c4762a1bSJed Brown Vec ctxVec = NULL; 212c4762a1bSJed Brown 213c4762a1bSJed Brown PetscFunctionBegin; 2145f80ce2aSJacob Faibussowitsch CHKERRQ(DMFieldShellGetContext(field, &ctxVec)); 2155f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&ctxVec)); 216c4762a1bSJed Brown PetscFunctionReturn(0); 217c4762a1bSJed Brown } 218c4762a1bSJed Brown 219c4762a1bSJed Brown int main(int argc, char **argv) 220c4762a1bSJed Brown { 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 PetscErrorCode ierr; 235c4762a1bSJed Brown 236*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc, &argv, NULL, help)); 237c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 238c4762a1bSJed Brown ierr = PetscOptionsBegin(comm, "", "DMField Tutorial Options", "DM");CHKERRQ(ierr); 2395f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsFList("-dm_type","DM implementation on which to define field","ex1.c",DMList,type,type,256,NULL)); 2405f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsRangeInt("-dim","DM intrinsic dimension", "ex1.c", dim, &dim, NULL,1,3)); 2415f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBoundedInt("-num_components","Number of components in field", "ex1.c", nc, &nc, NULL,1)); 2425f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBoundedInt("-num_quad_points","Number of quadrature points per dimension", "ex1.c", pointsPerEdge, &pointsPerEdge, NULL,1)); 2435f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBoundedInt("-num_point_tests", "Number of test points for DMFieldEvaluate()", "ex1.c", numPoint, &numPoint, NULL,0)); 2445f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBoundedInt("-num_fe_tests", "Number of test cells for DMFieldEvaluateFE()", "ex1.c", numFE, &numFE, NULL,0)); 2455f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBoundedInt("-num_fv_tests", "Number of test cells for DMFieldEvaluateFV()", "ex1.c", numFV, &numFV, NULL,0)); 2465f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-test_shell", "Test the DMFIELDSHELL implementation of DMField", "ex1.c", testShell, &testShell, NULL)); 247c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 248c4762a1bSJed Brown 2492c71b3e2SJacob Faibussowitsch PetscCheckFalse(dim > 3,comm,PETSC_ERR_ARG_OUTOFRANGE,"This examples works for dim <= 3, not %D",dim); 2505f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrncmp(type,DMPLEX,256,&isplex)); 2515f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrncmp(type,DMDA,256,&isda)); 252c4762a1bSJed Brown 2535f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomCreate(PETSC_COMM_SELF,&rand)); 2545f80ce2aSJacob Faibussowitsch CHKERRQ(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 262c4762a1bSJed Brown ierr = PetscOptionsBegin(comm, "", "DMField DMPlex Options", "DM");CHKERRQ(ierr); 2635f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBoundedInt("-overlap","DMPlex parallel overlap","ex1.c",overlap,&overlap,NULL,0)); 264c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 26530602db0SMatthew G. Knepley if (0) { 2665f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateBoxMesh(comm,2,PETSC_TRUE,cells,NULL,NULL,NULL,PETSC_TRUE,&dm)); 26730602db0SMatthew G. Knepley } else { 2685f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreate(comm, &dm)); 2695f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetType(dm, DMPLEX)); 2705f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(dm)); 27130602db0SMatthew G. Knepley CHKMEMQ; 272c4762a1bSJed Brown } 2735f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(dm, &dim)); 2745f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexIsSimplex(dm, &simplex)); 275c4762a1bSJed Brown if (simplex) { 2765f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDTStroudConicalQuadrature(dim, 1, pointsPerEdge, -1.0, 1.0, &quad)); 277c4762a1bSJed Brown } else { 2785f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDTGaussTensorQuadrature(dim, 1, pointsPerEdge, -1.0, 1.0, &quad)); 279c4762a1bSJed Brown } 2805f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetHeightStratum(dm,0,&cStart,&cEnd)); 281c4762a1bSJed Brown if (testShell) { 282c4762a1bSJed Brown Vec ctxVec; 283c4762a1bSJed Brown PetscInt i; 284c4762a1bSJed Brown PetscScalar *array; 285c4762a1bSJed Brown 2865f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeq(PETSC_COMM_SELF, nc, &ctxVec)); 2875f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetUp(ctxVec)); 2885f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(ctxVec,&array)); 289c4762a1bSJed Brown for (i = 0; i < nc; i++) array[i] = i + 1.; 2905f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(ctxVec,&array)); 2915f80ce2aSJacob Faibussowitsch CHKERRQ(DMFieldCreateShell(dm, nc, DMFIELD_VERTEX, (void *) ctxVec, &field)); 2925f80ce2aSJacob Faibussowitsch CHKERRQ(DMFieldShellSetEvaluate(field, TestShellEvaluate)); 2935f80ce2aSJacob Faibussowitsch CHKERRQ(DMFieldShellSetDestroy(field, TestShellDestroy)); 294c4762a1bSJed Brown } else { 2955f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFECreateDefault(PETSC_COMM_SELF,dim,nc,simplex,NULL,PETSC_DEFAULT,&fe)); 2965f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFESetName(fe,"MyPetscFE")); 2975f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetField(dm,0,NULL,(PetscObject)fe)); 2985f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFEDestroy(&fe)); 2995f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateDS(dm)); 3005f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateLocalVector(dm,&fieldvec)); 301c4762a1bSJed Brown { 302c4762a1bSJed Brown PetscErrorCode (*func[1]) (PetscInt,PetscReal,const PetscReal [],PetscInt, PetscScalar *,void *); 303c4762a1bSJed Brown void *ctxs[1]; 304c4762a1bSJed Brown 305c4762a1bSJed Brown func[0] = radiusSquared; 306c4762a1bSJed Brown ctxs[0] = NULL; 307c4762a1bSJed Brown 3085f80ce2aSJacob Faibussowitsch CHKERRQ(DMProjectFunctionLocal(dm,0.0,func,ctxs,INSERT_ALL_VALUES,fieldvec)); 309c4762a1bSJed Brown } 3105f80ce2aSJacob Faibussowitsch CHKERRQ(DMFieldCreateDS(dm,0,fieldvec,&field)); 3115f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&fieldvec)); 312c4762a1bSJed Brown } 313c4762a1bSJed Brown } else if (isda) { 314c4762a1bSJed Brown PetscInt i; 315c4762a1bSJed Brown PetscScalar *cv; 316c4762a1bSJed Brown 317c4762a1bSJed Brown switch (dim) { 318c4762a1bSJed Brown case 1: 3195f80ce2aSJacob Faibussowitsch CHKERRQ(DMDACreate1d(comm, DM_BOUNDARY_NONE, 3, 1, 1, NULL, &dm)); 320c4762a1bSJed Brown break; 321c4762a1bSJed Brown case 2: 3225f80ce2aSJacob Faibussowitsch CHKERRQ(DMDACreate2d(comm, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, 3, 3, PETSC_DETERMINE, PETSC_DETERMINE, 1, 1, NULL, NULL, &dm)); 323c4762a1bSJed Brown break; 324c4762a1bSJed Brown default: 3255f80ce2aSJacob Faibussowitsch CHKERRQ(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)); 326c4762a1bSJed Brown break; 327c4762a1bSJed Brown } 3285f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(dm)); 3295f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetHeightStratum(dm,0,&cStart,&cEnd)); 3305f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(nc * (1 << dim),&cv)); 331c4762a1bSJed Brown for (i = 0; i < nc * (1 << dim); i++) { 332c4762a1bSJed Brown PetscReal rv; 333c4762a1bSJed Brown 3345f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomGetValueReal(rand,&rv)); 335c4762a1bSJed Brown cv[i] = rv; 336c4762a1bSJed Brown } 3375f80ce2aSJacob Faibussowitsch CHKERRQ(DMFieldCreateDA(dm,nc,cv,&field)); 3385f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(cv)); 3395f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDTGaussTensorQuadrature(dim, 1, pointsPerEdge, -1.0, 1.0, &quad)); 34098921bdaSJacob Faibussowitsch } else SETERRQ(comm,PETSC_ERR_SUP,"This test does not run for DM type %s",type); 341c4762a1bSJed Brown 3425f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)dm,"mesh")); 3435f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(dm,NULL,"-dm_view")); 3445f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dm)); 3455f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)field,"field")); 3465f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectViewFromOptions((PetscObject)field,NULL,"-dmfield_view")); 3475f80ce2aSJacob Faibussowitsch if (numPoint) CHKERRQ(TestEvaluate(field,numPoint,rand)); 3485f80ce2aSJacob Faibussowitsch if (numFE) CHKERRQ(TestEvaluateFE(field,numFE,cStart,cEnd,quad,rand)); 3495f80ce2aSJacob Faibussowitsch if (numFV) CHKERRQ(TestEvaluateFV(field,numFV,cStart,cEnd,rand)); 3505f80ce2aSJacob Faibussowitsch CHKERRQ(DMFieldDestroy(&field)); 3515f80ce2aSJacob Faibussowitsch CHKERRQ(PetscQuadratureDestroy(&quad)); 3525f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomDestroy(&rand)); 353*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 354*b122ec5aSJacob Faibussowitsch return 0; 355c4762a1bSJed Brown } 356c4762a1bSJed Brown 357c4762a1bSJed Brown /*TEST 358c4762a1bSJed Brown 359c4762a1bSJed Brown test: 360c4762a1bSJed Brown suffix: da 361c4762a1bSJed Brown requires: !complex 362c4762a1bSJed Brown args: -dm_type da -dim 2 -num_components 2 -num_point_tests 2 -num_fe_tests 2 -num_fv_tests 2 -dmfield_view 363c4762a1bSJed Brown 364c4762a1bSJed Brown test: 365c4762a1bSJed Brown suffix: da_1 366c4762a1bSJed Brown requires: !complex 367c4762a1bSJed Brown args: -dm_type da -dim 1 -num_fe_tests 2 368c4762a1bSJed Brown 369c4762a1bSJed Brown test: 370c4762a1bSJed Brown suffix: da_2 371c4762a1bSJed Brown requires: !complex 372c4762a1bSJed Brown args: -dm_type da -dim 2 -num_fe_tests 2 373c4762a1bSJed Brown 374c4762a1bSJed Brown test: 375c4762a1bSJed Brown suffix: da_3 376c4762a1bSJed Brown requires: !complex 377c4762a1bSJed Brown args: -dm_type da -dim 3 -num_fe_tests 2 378c4762a1bSJed Brown 379c4762a1bSJed Brown test: 380c4762a1bSJed Brown suffix: ds 381c4762a1bSJed Brown requires: !complex triangle 38230602db0SMatthew 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 383c4762a1bSJed Brown 384c4762a1bSJed Brown test: 385c4762a1bSJed Brown suffix: ds_simplex_0 386c4762a1bSJed Brown requires: !complex triangle 38730602db0SMatthew G. Knepley args: -dm_type plex -dm_plex_dim 2 -dm_plex_box_faces 3,3 -num_fe_tests 2 -petscspace_degree 0 388c4762a1bSJed Brown 389c4762a1bSJed Brown test: 390c4762a1bSJed Brown suffix: ds_simplex_1 391c4762a1bSJed Brown requires: !complex triangle 39230602db0SMatthew G. Knepley args: -dm_type plex -dm_plex_dim 2 -dm_plex_box_faces 3,3 -num_fe_tests 2 -petscspace_degree 1 393c4762a1bSJed Brown 394c4762a1bSJed Brown test: 395c4762a1bSJed Brown suffix: ds_simplex_2 396c4762a1bSJed Brown requires: !complex triangle 39730602db0SMatthew G. Knepley args: -dm_type plex -dm_plex_dim 2 -dm_plex_box_faces 3,3 -num_fe_tests 2 -petscspace_degree 2 398c4762a1bSJed Brown 399c4762a1bSJed Brown test: 400c4762a1bSJed Brown suffix: ds_tensor_2_0 401c4762a1bSJed Brown requires: !complex 40230602db0SMatthew 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 403c4762a1bSJed Brown 404c4762a1bSJed Brown test: 405c4762a1bSJed Brown suffix: ds_tensor_2_1 406c4762a1bSJed Brown requires: !complex 40730602db0SMatthew 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 408c4762a1bSJed Brown 409c4762a1bSJed Brown test: 410c4762a1bSJed Brown suffix: ds_tensor_2_2 411c4762a1bSJed Brown requires: !complex 41230602db0SMatthew 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 413c4762a1bSJed Brown 414c4762a1bSJed Brown test: 415c4762a1bSJed Brown suffix: ds_tensor_3_0 416c4762a1bSJed Brown requires: !complex 41730602db0SMatthew 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 418c4762a1bSJed Brown 419c4762a1bSJed Brown test: 420c4762a1bSJed Brown suffix: ds_tensor_3_1 421c4762a1bSJed Brown requires: !complex 42230602db0SMatthew 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 423c4762a1bSJed Brown 424c4762a1bSJed Brown test: 425c4762a1bSJed Brown suffix: ds_tensor_3_2 426c4762a1bSJed Brown requires: !complex 42730602db0SMatthew 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 428c4762a1bSJed Brown 429c4762a1bSJed Brown test: 430c4762a1bSJed Brown suffix: shell 431c4762a1bSJed Brown requires: !complex triangle 43230602db0SMatthew 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 433c4762a1bSJed Brown 434c4762a1bSJed Brown TEST*/ 435