1c4762a1bSJed Brown static char help[] = "Test for function and field projection\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscdmplex.h> 4c4762a1bSJed Brown #include <petscds.h> 5c4762a1bSJed Brown 6c4762a1bSJed Brown typedef struct { 7c4762a1bSJed Brown PetscBool multifield; /* Different numbers of input and output fields */ 8c4762a1bSJed Brown PetscBool subdomain; /* Try with a volumetric submesh */ 9c4762a1bSJed Brown PetscBool submesh; /* Try with a boundary submesh */ 10c4762a1bSJed Brown PetscBool auxfield; /* Try with auxiliary fields */ 11c4762a1bSJed Brown } AppCtx; 12c4762a1bSJed Brown 13c4762a1bSJed Brown /* (x + y)*dim + d */ 14c4762a1bSJed Brown static PetscErrorCode linear(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 15c4762a1bSJed Brown { 16c4762a1bSJed Brown PetscInt c; 17c4762a1bSJed Brown for (c = 0; c < Nc; ++c) u[c] = (x[0] + x[1])*Nc + c; 18c4762a1bSJed Brown return 0; 19c4762a1bSJed Brown } 20c4762a1bSJed Brown 21c4762a1bSJed Brown /* {x, y, z} */ 22c4762a1bSJed Brown static PetscErrorCode linear2(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 23c4762a1bSJed Brown { 24c4762a1bSJed Brown PetscInt c; 25c4762a1bSJed Brown for (c = 0; c < Nc; ++c) u[c] = x[c]; 26c4762a1bSJed Brown return 0; 27c4762a1bSJed Brown } 28c4762a1bSJed Brown 29c4762a1bSJed Brown /* {u_x, u_y, u_z} */ 30c4762a1bSJed Brown static void linear_vector(PetscInt dim, PetscInt Nf, PetscInt NfAux, 31c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 32c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 33c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f[]) 34c4762a1bSJed Brown { 35c4762a1bSJed Brown PetscInt d; 36c4762a1bSJed Brown for (d = 0; d < uOff[1]-uOff[0]; ++d) f[d] = u[d+uOff[0]]; 37c4762a1bSJed Brown } 38c4762a1bSJed Brown 39c4762a1bSJed Brown /* p */ 40c4762a1bSJed Brown static void linear_scalar(PetscInt dim, PetscInt Nf, PetscInt NfAux, 41c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 42c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 43c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f[]) 44c4762a1bSJed Brown { 45c4762a1bSJed Brown f[0] = u[uOff[1]]; 46c4762a1bSJed Brown } 47c4762a1bSJed Brown 48c4762a1bSJed Brown /* {div u, p^2} */ 49c4762a1bSJed Brown static void divergence_sq(PetscInt dim, PetscInt Nf, PetscInt NfAux, 50c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 51c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 52c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f[]) 53c4762a1bSJed Brown { 54c4762a1bSJed Brown PetscInt d; 55c4762a1bSJed Brown f[0] = 0.0; 56c4762a1bSJed Brown for (d = 0; d < dim; ++d) f[0] += u_x[uOff_x[0]+d*dim+d]; 57c4762a1bSJed Brown f[1] = PetscSqr(u[uOff[1]]); 58c4762a1bSJed Brown } 59c4762a1bSJed Brown 60c4762a1bSJed Brown static PetscErrorCode ProcessOptions(AppCtx *options) 61c4762a1bSJed Brown { 62c4762a1bSJed Brown PetscErrorCode ierr; 63c4762a1bSJed Brown 64c4762a1bSJed Brown PetscFunctionBegin; 65c4762a1bSJed Brown options->multifield = PETSC_FALSE; 66c4762a1bSJed Brown options->subdomain = PETSC_FALSE; 67c4762a1bSJed Brown options->submesh = PETSC_FALSE; 68c4762a1bSJed Brown options->auxfield = PETSC_FALSE; 69c4762a1bSJed Brown 70c4762a1bSJed Brown ierr = PetscOptionsBegin(PETSC_COMM_SELF, "", "Meshing Problem Options", "DMPLEX");CHKERRQ(ierr); 715f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-multifield", "Flag for trying different numbers of input/output fields", "ex23.c", options->multifield, &options->multifield, NULL)); 725f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-subdomain", "Flag for trying volumetric submesh", "ex23.c", options->subdomain, &options->subdomain, NULL)); 735f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-submesh", "Flag for trying boundary submesh", "ex23.c", options->submesh, &options->submesh, NULL)); 745f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-auxfield", "Flag for trying auxiliary fields", "ex23.c", options->auxfield, &options->auxfield, NULL)); 75c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 76c4762a1bSJed Brown PetscFunctionReturn(0); 77c4762a1bSJed Brown } 78c4762a1bSJed Brown 79c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 80c4762a1bSJed Brown { 81c4762a1bSJed Brown PetscFunctionBegin; 825f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreate(comm, dm)); 835f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetType(*dm, DMPLEX)); 845f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(*dm)); 855f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(*dm, NULL, "-orig_dm_view")); 86c4762a1bSJed Brown PetscFunctionReturn(0); 87c4762a1bSJed Brown } 88c4762a1bSJed Brown 89c4762a1bSJed Brown static PetscErrorCode SetupDiscretization(DM dm, PetscInt dim, PetscBool simplex, AppCtx *user) 90c4762a1bSJed Brown { 91c4762a1bSJed Brown PetscFE fe; 92c4762a1bSJed Brown MPI_Comm comm; 93c4762a1bSJed Brown 94c4762a1bSJed Brown PetscFunctionBeginUser; 955f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetComm((PetscObject) dm, &comm)); 965f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFECreateDefault(comm, dim, dim, simplex, "velocity_", -1, &fe)); 975f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) fe, "velocity")); 985f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetField(dm, 0, NULL, (PetscObject) fe)); 995f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFEDestroy(&fe)); 1005f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFECreateDefault(comm, dim, 1, simplex, "pressure_", -1, &fe)); 1015f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) fe, "pressure")); 1025f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetField(dm, 1, NULL, (PetscObject) fe)); 1035f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFEDestroy(&fe)); 1045f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateDS(dm)); 105c4762a1bSJed Brown PetscFunctionReturn(0); 106c4762a1bSJed Brown } 107c4762a1bSJed Brown 108c4762a1bSJed Brown static PetscErrorCode SetupOutputDiscretization(DM dm, PetscInt dim, PetscBool simplex, AppCtx *user) 109c4762a1bSJed Brown { 110c4762a1bSJed Brown PetscFE fe; 111c4762a1bSJed Brown MPI_Comm comm; 112c4762a1bSJed Brown 113c4762a1bSJed Brown PetscFunctionBeginUser; 1145f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetComm((PetscObject) dm, &comm)); 1155f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFECreateDefault(comm, dim, dim, simplex, "output_", -1, &fe)); 1165f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) fe, "output")); 1175f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetField(dm, 0, NULL, (PetscObject) fe)); 1185f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFEDestroy(&fe)); 1195f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateDS(dm)); 120c4762a1bSJed Brown PetscFunctionReturn(0); 121c4762a1bSJed Brown } 122c4762a1bSJed Brown 123c4762a1bSJed Brown static PetscErrorCode CreateSubdomainMesh(DM dm, DMLabel *domLabel, DM *subdm, AppCtx *user) 124c4762a1bSJed Brown { 125c4762a1bSJed Brown DMLabel label; 12630602db0SMatthew G. Knepley PetscBool simplex; 127c4762a1bSJed Brown PetscInt dim, cStart, cEnd, c; 128c4762a1bSJed Brown 129c4762a1bSJed Brown PetscFunctionBeginUser; 1305f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexIsSimplex(dm, &simplex)); 1315f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 1325f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelCreate(PETSC_COMM_SELF, "subdomain", &label)); 1335f80ce2aSJacob Faibussowitsch for (c = cStart + (cEnd-cStart)/2; c < cEnd; ++c) CHKERRQ(DMLabelSetValue(label, c, 1)); 1345f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexFilter(dm, label, 1, subdm)); 1355f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(*subdm, &dim)); 1365f80ce2aSJacob Faibussowitsch CHKERRQ(SetupDiscretization(*subdm, dim, simplex, user)); 1375f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) *subdm, "subdomain")); 1385f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(*subdm, NULL, "-sub_dm_view")); 139c4762a1bSJed Brown if (domLabel) *domLabel = label; 1405f80ce2aSJacob Faibussowitsch else CHKERRQ(DMLabelDestroy(&label)); 141c4762a1bSJed Brown PetscFunctionReturn(0); 142c4762a1bSJed Brown } 143c4762a1bSJed Brown 144c4762a1bSJed Brown static PetscErrorCode CreateBoundaryMesh(DM dm, DMLabel *bdLabel, DM *subdm, AppCtx *user) 145c4762a1bSJed Brown { 146c4762a1bSJed Brown DMLabel label; 14730602db0SMatthew G. Knepley PetscBool simplex; 148c4762a1bSJed Brown PetscInt dim; 149c4762a1bSJed Brown 150c4762a1bSJed Brown PetscFunctionBeginUser; 1515f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexIsSimplex(dm, &simplex)); 1525f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelCreate(PETSC_COMM_SELF, "sub", &label)); 1535f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexMarkBoundaryFaces(dm, 1, label)); 1545f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexLabelComplete(dm, label)); 1555f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateSubmesh(dm, label, 1, PETSC_TRUE, subdm)); 1565f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(*subdm, &dim)); 1575f80ce2aSJacob Faibussowitsch CHKERRQ(SetupDiscretization(*subdm, dim, simplex, user)); 1585f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) *subdm, "boundary")); 1595f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(*subdm, NULL, "-sub_dm_view")); 160c4762a1bSJed Brown if (bdLabel) *bdLabel = label; 1615f80ce2aSJacob Faibussowitsch else CHKERRQ(DMLabelDestroy(&label)); 162c4762a1bSJed Brown PetscFunctionReturn(0); 163c4762a1bSJed Brown } 164c4762a1bSJed Brown 1659a2a23afSMatthew G. Knepley static PetscErrorCode CreateAuxiliaryVec(DM dm, DM *auxdm, Vec *la, AppCtx *user) 166c4762a1bSJed Brown { 167c4762a1bSJed Brown PetscErrorCode (**afuncs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *); 16830602db0SMatthew G. Knepley PetscBool simplex; 169c4762a1bSJed Brown PetscInt dim, Nf, f; 170c4762a1bSJed Brown 171c4762a1bSJed Brown PetscFunctionBeginUser; 1725f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(dm, &dim)); 1735f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexIsSimplex(dm, &simplex)); 1745f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetNumFields(dm, &Nf)); 1755f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(Nf, &afuncs)); 176c4762a1bSJed Brown for (f = 0; f < Nf; ++f) afuncs[f] = linear; 1775f80ce2aSJacob Faibussowitsch CHKERRQ(DMClone(dm, auxdm)); 1785f80ce2aSJacob Faibussowitsch CHKERRQ(SetupDiscretization(*auxdm, dim, simplex, user)); 1795f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateLocalVector(*auxdm, la)); 1805f80ce2aSJacob Faibussowitsch CHKERRQ(DMProjectFunctionLocal(dm, 0.0, afuncs, NULL, INSERT_VALUES, *la)); 1815f80ce2aSJacob Faibussowitsch CHKERRQ(VecViewFromOptions(*la, NULL, "-local_aux_view")); 1825f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(afuncs)); 183c4762a1bSJed Brown PetscFunctionReturn(0); 184c4762a1bSJed Brown } 185c4762a1bSJed Brown 186c4762a1bSJed Brown static PetscErrorCode TestFunctionProjection(DM dm, DM dmAux, DMLabel label, Vec la, const char name[], AppCtx *user) 187c4762a1bSJed Brown { 188c4762a1bSJed Brown PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *); 189c4762a1bSJed Brown Vec x, lx; 190c4762a1bSJed Brown PetscInt Nf, f; 191c4762a1bSJed Brown PetscInt val[1] = {1}; 192c4762a1bSJed Brown char lname[PETSC_MAX_PATH_LEN]; 193c4762a1bSJed Brown 194c4762a1bSJed Brown PetscFunctionBeginUser; 1955f80ce2aSJacob Faibussowitsch if (dmAux) CHKERRQ(DMSetAuxiliaryVec(dm, NULL, 0, 0, la)); 1965f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetNumFields(dm, &Nf)); 1975f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(Nf, &funcs)); 198c4762a1bSJed Brown for (f = 0; f < Nf; ++f) funcs[f] = linear; 1995f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetGlobalVector(dm, &x)); 2005f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrcpy(lname, "Function ")); 2015f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrcat(lname, name)); 2025f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) x, lname)); 2035f80ce2aSJacob Faibussowitsch if (!label) CHKERRQ(DMProjectFunction(dm, 0.0, funcs, NULL, INSERT_VALUES, x)); 2045f80ce2aSJacob Faibussowitsch else CHKERRQ(DMProjectFunctionLabel(dm, 0.0, label, 1, val, 0, NULL, funcs, NULL, INSERT_VALUES, x)); 2055f80ce2aSJacob Faibussowitsch CHKERRQ(VecViewFromOptions(x, NULL, "-func_view")); 2065f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreGlobalVector(dm, &x)); 2075f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(dm, &lx)); 2085f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrcpy(lname, "Local Function ")); 2095f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrcat(lname, name)); 2105f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) lx, lname)); 2115f80ce2aSJacob Faibussowitsch if (!label) CHKERRQ(DMProjectFunctionLocal(dm, 0.0, funcs, NULL, INSERT_VALUES, lx)); 2125f80ce2aSJacob Faibussowitsch else CHKERRQ(DMProjectFunctionLabelLocal(dm, 0.0, label, 1, val, 0, NULL, funcs, NULL, INSERT_VALUES, lx)); 2135f80ce2aSJacob Faibussowitsch CHKERRQ(VecViewFromOptions(lx, NULL, "-local_func_view")); 2145f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(dm, &lx)); 2155f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(funcs)); 2165f80ce2aSJacob Faibussowitsch if (dmAux) CHKERRQ(DMSetAuxiliaryVec(dm, NULL, 0, 0, NULL)); 217c4762a1bSJed Brown PetscFunctionReturn(0); 218c4762a1bSJed Brown } 219c4762a1bSJed Brown 220c4762a1bSJed Brown static PetscErrorCode TestFieldProjection(DM dm, DM dmAux, DMLabel label, Vec la, const char name[], AppCtx *user) 221c4762a1bSJed Brown { 222c4762a1bSJed Brown PetscErrorCode (**afuncs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *); 223c4762a1bSJed Brown void (**funcs)(PetscInt, PetscInt, PetscInt, 224c4762a1bSJed Brown const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 225c4762a1bSJed Brown const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 226c4762a1bSJed Brown PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]); 227c4762a1bSJed Brown Vec lx, lu; 228c4762a1bSJed Brown PetscInt Nf, f; 229c4762a1bSJed Brown PetscInt val[1] = {1}; 230c4762a1bSJed Brown char lname[PETSC_MAX_PATH_LEN]; 231c4762a1bSJed Brown 232c4762a1bSJed Brown PetscFunctionBeginUser; 2335f80ce2aSJacob Faibussowitsch if (dmAux) CHKERRQ(DMSetAuxiliaryVec(dm, NULL, 0, 0, la)); 2345f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetNumFields(dm, &Nf)); 2355f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc2(Nf, &funcs, Nf, &afuncs)); 236c4762a1bSJed Brown for (f = 0; f < Nf; ++f) afuncs[f] = linear; 237c4762a1bSJed Brown funcs[0] = linear_vector; 238c4762a1bSJed Brown funcs[1] = linear_scalar; 2395f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(dm, &lu)); 2405f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrcpy(lname, "Local Field Input ")); 2415f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrcat(lname, name)); 2425f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) lu, lname)); 2435f80ce2aSJacob Faibussowitsch if (!label) CHKERRQ(DMProjectFunctionLocal(dm, 0.0, afuncs, NULL, INSERT_VALUES, lu)); 2445f80ce2aSJacob Faibussowitsch else CHKERRQ(DMProjectFunctionLabelLocal(dm, 0.0, label, 1, val, 0, NULL, afuncs, NULL, INSERT_VALUES, lu)); 2455f80ce2aSJacob Faibussowitsch CHKERRQ(VecViewFromOptions(lu, NULL, "-local_input_view")); 2465f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(dm, &lx)); 2475f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrcpy(lname, "Local Field ")); 2485f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrcat(lname, name)); 2495f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) lx, lname)); 2505f80ce2aSJacob Faibussowitsch if (!label) CHKERRQ(DMProjectFieldLocal(dm, 0.0, lu, funcs, INSERT_VALUES, lx)); 2515f80ce2aSJacob Faibussowitsch else CHKERRQ(DMProjectFieldLabelLocal(dm, 0.0, label, 1, val, 0, NULL, lu, funcs, INSERT_VALUES, lx)); 2525f80ce2aSJacob Faibussowitsch CHKERRQ(VecViewFromOptions(lx, NULL, "-local_field_view")); 2535f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(dm, &lx)); 2545f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(dm, &lu)); 2555f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree2(funcs, afuncs)); 2565f80ce2aSJacob Faibussowitsch if (dmAux) CHKERRQ(DMSetAuxiliaryVec(dm, NULL, 0, 0, NULL)); 257c4762a1bSJed Brown PetscFunctionReturn(0); 258c4762a1bSJed Brown } 259c4762a1bSJed Brown 260c4762a1bSJed Brown static PetscErrorCode TestFieldProjectionMultiple(DM dm, DM dmIn, DM dmAux, DMLabel label, Vec la, const char name[], AppCtx *user) 261c4762a1bSJed Brown { 262c4762a1bSJed Brown PetscErrorCode (**afuncs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *); 263c4762a1bSJed Brown void (**funcs)(PetscInt, PetscInt, PetscInt, 264c4762a1bSJed Brown const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 265c4762a1bSJed Brown const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 266c4762a1bSJed Brown PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]); 267c4762a1bSJed Brown Vec lx, lu; 268c4762a1bSJed Brown PetscInt Nf, NfIn; 269c4762a1bSJed Brown PetscInt val[1] = {1}; 270c4762a1bSJed Brown char lname[PETSC_MAX_PATH_LEN]; 271c4762a1bSJed Brown 272c4762a1bSJed Brown PetscFunctionBeginUser; 2735f80ce2aSJacob Faibussowitsch if (dmAux) CHKERRQ(DMSetAuxiliaryVec(dm, NULL, 0, 0, la)); 2745f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetNumFields(dm, &Nf)); 2755f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetNumFields(dmIn, &NfIn)); 2765f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc2(Nf, &funcs, NfIn, &afuncs)); 277c4762a1bSJed Brown funcs[0] = divergence_sq; 278c4762a1bSJed Brown afuncs[0] = linear2; 279c4762a1bSJed Brown afuncs[1] = linear; 2805f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(dmIn, &lu)); 2815f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrcpy(lname, "Local MultiField Input ")); 2825f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrcat(lname, name)); 2835f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) lu, lname)); 2845f80ce2aSJacob Faibussowitsch if (!label) CHKERRQ(DMProjectFunctionLocal(dmIn, 0.0, afuncs, NULL, INSERT_VALUES, lu)); 2855f80ce2aSJacob Faibussowitsch else CHKERRQ(DMProjectFunctionLabelLocal(dmIn, 0.0, label, 1, val, 0, NULL, afuncs, NULL, INSERT_VALUES, lu)); 2865f80ce2aSJacob Faibussowitsch CHKERRQ(VecViewFromOptions(lu, NULL, "-local_input_view")); 2875f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(dm, &lx)); 2885f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrcpy(lname, "Local MultiField ")); 2895f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrcat(lname, name)); 2905f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) lx, lname)); 2915f80ce2aSJacob Faibussowitsch if (!label) CHKERRQ(DMProjectFieldLocal(dm, 0.0, lu, funcs, INSERT_VALUES, lx)); 2925f80ce2aSJacob Faibussowitsch else CHKERRQ(DMProjectFieldLabelLocal(dm, 0.0, label, 1, val, 0, NULL, lu, funcs, INSERT_VALUES, lx)); 2935f80ce2aSJacob Faibussowitsch CHKERRQ(VecViewFromOptions(lx, NULL, "-local_field_view")); 2945f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(dm, &lx)); 2955f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(dmIn, &lu)); 2965f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree2(funcs, afuncs)); 2975f80ce2aSJacob Faibussowitsch if (dmAux) CHKERRQ(DMSetAuxiliaryVec(dm, NULL, 0, 0, NULL)); 298c4762a1bSJed Brown PetscFunctionReturn(0); 299c4762a1bSJed Brown } 300c4762a1bSJed Brown 301c4762a1bSJed Brown int main(int argc, char **argv) 302c4762a1bSJed Brown { 303c4762a1bSJed Brown DM dm, subdm, auxdm; 304c4762a1bSJed Brown Vec la; 30530602db0SMatthew G. Knepley PetscInt dim; 30630602db0SMatthew G. Knepley PetscBool simplex; 307c4762a1bSJed Brown AppCtx user; 308c4762a1bSJed Brown 309*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc, &argv, NULL,help)); 3105f80ce2aSJacob Faibussowitsch CHKERRQ(ProcessOptions(&user)); 3115f80ce2aSJacob Faibussowitsch CHKERRQ(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 3125f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(dm, &dim)); 3135f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexIsSimplex(dm, &simplex)); 3145f80ce2aSJacob Faibussowitsch CHKERRQ(SetupDiscretization(dm, dim, simplex, &user)); 315c4762a1bSJed Brown /* Volumetric Mesh Projection */ 316c4762a1bSJed Brown if (!user.multifield) { 3175f80ce2aSJacob Faibussowitsch CHKERRQ(TestFunctionProjection(dm, NULL, NULL, NULL, "Volumetric Primary", &user)); 3185f80ce2aSJacob Faibussowitsch CHKERRQ(TestFieldProjection(dm, NULL, NULL, NULL, "Volumetric Primary", &user)); 319c4762a1bSJed Brown } else { 320c4762a1bSJed Brown DM dmOut; 321c4762a1bSJed Brown 3225f80ce2aSJacob Faibussowitsch CHKERRQ(DMClone(dm, &dmOut)); 3235f80ce2aSJacob Faibussowitsch CHKERRQ(SetupOutputDiscretization(dmOut, dim, simplex, &user)); 3245f80ce2aSJacob Faibussowitsch CHKERRQ(TestFieldProjectionMultiple(dmOut, dm, NULL, NULL, NULL, "Volumetric Primary", &user)); 3255f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dmOut)); 326c4762a1bSJed Brown } 327c4762a1bSJed Brown if (user.auxfield) { 328c4762a1bSJed Brown /* Volumetric Mesh Projection with Volumetric Data */ 3295f80ce2aSJacob Faibussowitsch CHKERRQ(CreateAuxiliaryVec(dm, &auxdm, &la, &user)); 3305f80ce2aSJacob Faibussowitsch CHKERRQ(TestFunctionProjection(dm, auxdm, NULL, la, "Volumetric Primary and Volumetric Auxiliary", &user)); 3315f80ce2aSJacob Faibussowitsch CHKERRQ(TestFieldProjection(dm, auxdm, NULL, la, "Volumetric Primary and Volumetric Auxiliary", &user)); 3325f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&la)); 333c4762a1bSJed Brown /* Update of Volumetric Auxiliary Data with primary Volumetric Data */ 3345f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(dm, &la)); 3355f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(la, 1.0)); 3365f80ce2aSJacob Faibussowitsch CHKERRQ(TestFieldProjection(auxdm, dm, NULL, la, "Volumetric Auxiliary Update with Volumetric Primary", &user)); 3375f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(dm, &la)); 3385f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&auxdm)); 339c4762a1bSJed Brown } 340c4762a1bSJed Brown if (user.subdomain) { 341c4762a1bSJed Brown DMLabel domLabel; 342c4762a1bSJed Brown 343c4762a1bSJed Brown /* Subdomain Mesh Projection */ 3445f80ce2aSJacob Faibussowitsch CHKERRQ(CreateSubdomainMesh(dm, &domLabel, &subdm, &user)); 3455f80ce2aSJacob Faibussowitsch CHKERRQ(TestFunctionProjection(subdm, NULL, NULL, NULL, "Subdomain Primary", &user)); 3465f80ce2aSJacob Faibussowitsch CHKERRQ(TestFieldProjection(subdm, NULL, NULL, NULL, "Subdomain Primary", &user)); 347c4762a1bSJed Brown if (user.auxfield) { 348c4762a1bSJed Brown /* Subdomain Mesh Projection with Subdomain Data */ 3495f80ce2aSJacob Faibussowitsch CHKERRQ(CreateAuxiliaryVec(subdm, &auxdm, &la, &user)); 3505f80ce2aSJacob Faibussowitsch CHKERRQ(TestFunctionProjection(subdm, auxdm, NULL, la, "Subdomain Primary and Subdomain Auxiliary", &user)); 3515f80ce2aSJacob Faibussowitsch CHKERRQ(TestFieldProjection(subdm, auxdm, NULL, la, "Subdomain Primary and Subdomain Auxiliary", &user)); 3525f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&la)); 3535f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&auxdm)); 354c4762a1bSJed Brown /* Subdomain Mesh Projection with Volumetric Data */ 3555f80ce2aSJacob Faibussowitsch CHKERRQ(CreateAuxiliaryVec(dm, &auxdm, &la, &user)); 3565f80ce2aSJacob Faibussowitsch CHKERRQ(TestFunctionProjection(subdm, auxdm, NULL, la, "Subdomain Primary and Volumetric Auxiliary", &user)); 3575f80ce2aSJacob Faibussowitsch CHKERRQ(TestFieldProjection(subdm, auxdm, NULL, la, "Subdomain Primary and Volumetric Auxiliary", &user)); 3585f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&la)); 3595f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&auxdm)); 360c4762a1bSJed Brown /* Volumetric Mesh Projection with Subdomain Data */ 3615f80ce2aSJacob Faibussowitsch CHKERRQ(CreateAuxiliaryVec(subdm, &auxdm, &la, &user)); 3625f80ce2aSJacob Faibussowitsch CHKERRQ(TestFunctionProjection(subdm, auxdm, domLabel, la, "Volumetric Primary and Subdomain Auxiliary", &user)); 3635f80ce2aSJacob Faibussowitsch CHKERRQ(TestFieldProjection(subdm, auxdm, domLabel, la, "Volumetric Primary and Subdomain Auxiliary", &user)); 3645f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&la)); 3655f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&auxdm)); 366c4762a1bSJed Brown } 3675f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&subdm)); 3685f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelDestroy(&domLabel)); 369c4762a1bSJed Brown } 370c4762a1bSJed Brown if (user.submesh) { 371c4762a1bSJed Brown DMLabel bdLabel; 372c4762a1bSJed Brown 373c4762a1bSJed Brown /* Boundary Mesh Projection */ 3745f80ce2aSJacob Faibussowitsch CHKERRQ(CreateBoundaryMesh(dm, &bdLabel, &subdm, &user)); 3755f80ce2aSJacob Faibussowitsch CHKERRQ(TestFunctionProjection(subdm, NULL, NULL, NULL, "Boundary Primary", &user)); 3765f80ce2aSJacob Faibussowitsch CHKERRQ(TestFieldProjection(subdm, NULL, NULL, NULL, "Boundary Primary", &user)); 377c4762a1bSJed Brown if (user.auxfield) { 378c4762a1bSJed Brown /* Boundary Mesh Projection with Boundary Data */ 3795f80ce2aSJacob Faibussowitsch CHKERRQ(CreateAuxiliaryVec(subdm, &auxdm, &la, &user)); 3805f80ce2aSJacob Faibussowitsch CHKERRQ(TestFunctionProjection(subdm, auxdm, NULL, la, "Boundary Primary and Boundary Auxiliary", &user)); 3815f80ce2aSJacob Faibussowitsch CHKERRQ(TestFieldProjection(subdm, auxdm, NULL, la, "Boundary Primary and Boundary Auxiliary", &user)); 3825f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&la)); 3835f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&auxdm)); 384c4762a1bSJed Brown /* Volumetric Mesh Projection with Boundary Data */ 3855f80ce2aSJacob Faibussowitsch CHKERRQ(CreateAuxiliaryVec(subdm, &auxdm, &la, &user)); 3865f80ce2aSJacob Faibussowitsch CHKERRQ(TestFunctionProjection(dm, auxdm, bdLabel, la, "Volumetric Primary and Boundary Auxiliary", &user)); 3875f80ce2aSJacob Faibussowitsch CHKERRQ(TestFieldProjection(dm, auxdm, bdLabel, la, "Volumetric Primary and Boundary Auxiliary", &user)); 3885f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&la)); 3895f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&auxdm)); 390c4762a1bSJed Brown } 3915f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelDestroy(&bdLabel)); 3925f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&subdm)); 393c4762a1bSJed Brown } 3945f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dm)); 395*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 396*b122ec5aSJacob Faibussowitsch return 0; 397c4762a1bSJed Brown } 398c4762a1bSJed Brown 399c4762a1bSJed Brown /*TEST 400c4762a1bSJed Brown 401c4762a1bSJed Brown test: 402c4762a1bSJed Brown suffix: 0 403c4762a1bSJed Brown requires: triangle 40430602db0SMatthew G. Knepley args: -dm_plex_box_faces 1,1 -func_view -local_func_view -local_input_view -local_field_view 405c4762a1bSJed Brown test: 406c4762a1bSJed Brown suffix: mf_0 407c4762a1bSJed Brown requires: triangle 40830602db0SMatthew G. Knepley args: -dm_plex_box_faces 1,1 -velocity_petscspace_degree 1 -velocity_petscfe_default_quadrature_order 2 \ 409c4762a1bSJed Brown -pressure_petscspace_degree 2 -pressure_petscfe_default_quadrature_order 2 \ 410c4762a1bSJed Brown -multifield -output_petscspace_degree 1 -output_petscfe_default_quadrature_order 2 \ 411c4762a1bSJed Brown -local_input_view -local_field_view 412c4762a1bSJed Brown test: 413c4762a1bSJed Brown suffix: 1 414c4762a1bSJed Brown requires: triangle 41530602db0SMatthew G. Knepley args: -dm_plex_box_faces 1,1 -velocity_petscspace_degree 1 -velocity_petscfe_default_quadrature_order 2 -pressure_petscspace_degree 2 -pressure_petscfe_default_quadrature_order 2 -func_view -local_func_view -local_input_view -local_field_view -submesh -auxfield 416c4762a1bSJed Brown test: 417c4762a1bSJed Brown suffix: 2 418c4762a1bSJed Brown requires: triangle 41930602db0SMatthew G. Knepley args: -dm_plex_box_faces 1,1 -velocity_petscspace_degree 1 -velocity_petscfe_default_quadrature_order 2 -pressure_petscspace_degree 2 -pressure_petscfe_default_quadrature_order 2 -func_view -local_func_view -local_input_view -local_field_view -subdomain -auxfield 420c4762a1bSJed Brown 421c4762a1bSJed Brown TEST*/ 422c4762a1bSJed Brown 423c4762a1bSJed Brown /* 424c4762a1bSJed Brown Post-processing wants to project a function of the fields into some FE space 425c4762a1bSJed Brown - This is DMProjectField() 426c4762a1bSJed Brown - What about changing the number of components of the output, like displacement to stress? Aux vars 427c4762a1bSJed Brown 428c4762a1bSJed Brown Update of state variables 429c4762a1bSJed Brown - This is DMProjectField(), but solution must be the aux var 430c4762a1bSJed Brown */ 431