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 */ 14*9371c9d4SSatish Balay static PetscErrorCode linear(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) { 15c4762a1bSJed Brown PetscInt c; 16c4762a1bSJed Brown for (c = 0; c < Nc; ++c) u[c] = (x[0] + x[1]) * Nc + c; 17c4762a1bSJed Brown return 0; 18c4762a1bSJed Brown } 19c4762a1bSJed Brown 20c4762a1bSJed Brown /* {x, y, z} */ 21*9371c9d4SSatish Balay static PetscErrorCode linear2(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) { 22c4762a1bSJed Brown PetscInt c; 23c4762a1bSJed Brown for (c = 0; c < Nc; ++c) u[c] = x[c]; 24c4762a1bSJed Brown return 0; 25c4762a1bSJed Brown } 26c4762a1bSJed Brown 27c4762a1bSJed Brown /* {u_x, u_y, u_z} */ 28*9371c9d4SSatish Balay static void linear_vector(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f[]) { 29c4762a1bSJed Brown PetscInt d; 30c4762a1bSJed Brown for (d = 0; d < uOff[1] - uOff[0]; ++d) f[d] = u[d + uOff[0]]; 31c4762a1bSJed Brown } 32c4762a1bSJed Brown 33c4762a1bSJed Brown /* p */ 34*9371c9d4SSatish Balay static void linear_scalar(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f[]) { 35c4762a1bSJed Brown f[0] = u[uOff[1]]; 36c4762a1bSJed Brown } 37c4762a1bSJed Brown 38c4762a1bSJed Brown /* {div u, p^2} */ 39*9371c9d4SSatish Balay static void divergence_sq(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f[]) { 40c4762a1bSJed Brown PetscInt d; 41c4762a1bSJed Brown f[0] = 0.0; 42c4762a1bSJed Brown for (d = 0; d < dim; ++d) f[0] += u_x[uOff_x[0] + d * dim + d]; 43c4762a1bSJed Brown f[1] = PetscSqr(u[uOff[1]]); 44c4762a1bSJed Brown } 45c4762a1bSJed Brown 46*9371c9d4SSatish Balay static PetscErrorCode ProcessOptions(AppCtx *options) { 47c4762a1bSJed Brown PetscFunctionBegin; 48c4762a1bSJed Brown options->multifield = PETSC_FALSE; 49c4762a1bSJed Brown options->subdomain = PETSC_FALSE; 50c4762a1bSJed Brown options->submesh = PETSC_FALSE; 51c4762a1bSJed Brown options->auxfield = PETSC_FALSE; 52c4762a1bSJed Brown 53d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_SELF, "", "Meshing Problem Options", "DMPLEX"); 549566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-multifield", "Flag for trying different numbers of input/output fields", "ex23.c", options->multifield, &options->multifield, NULL)); 559566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-subdomain", "Flag for trying volumetric submesh", "ex23.c", options->subdomain, &options->subdomain, NULL)); 569566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-submesh", "Flag for trying boundary submesh", "ex23.c", options->submesh, &options->submesh, NULL)); 579566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-auxfield", "Flag for trying auxiliary fields", "ex23.c", options->auxfield, &options->auxfield, NULL)); 58d0609cedSBarry Smith PetscOptionsEnd(); 59c4762a1bSJed Brown PetscFunctionReturn(0); 60c4762a1bSJed Brown } 61c4762a1bSJed Brown 62*9371c9d4SSatish Balay static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) { 63c4762a1bSJed Brown PetscFunctionBegin; 649566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 659566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 669566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 679566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-orig_dm_view")); 68c4762a1bSJed Brown PetscFunctionReturn(0); 69c4762a1bSJed Brown } 70c4762a1bSJed Brown 71*9371c9d4SSatish Balay static PetscErrorCode SetupDiscretization(DM dm, PetscInt dim, PetscBool simplex, AppCtx *user) { 72c4762a1bSJed Brown PetscFE fe; 73c4762a1bSJed Brown MPI_Comm comm; 74c4762a1bSJed Brown 75c4762a1bSJed Brown PetscFunctionBeginUser; 769566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 779566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, dim, simplex, "velocity_", -1, &fe)); 789566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe, "velocity")); 799566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe)); 809566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 819566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, "pressure_", -1, &fe)); 829566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe, "pressure")); 839566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 1, NULL, (PetscObject)fe)); 849566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 859566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 86c4762a1bSJed Brown PetscFunctionReturn(0); 87c4762a1bSJed Brown } 88c4762a1bSJed Brown 89*9371c9d4SSatish Balay static PetscErrorCode SetupOutputDiscretization(DM dm, PetscInt dim, PetscBool simplex, AppCtx *user) { 90c4762a1bSJed Brown PetscFE fe; 91c4762a1bSJed Brown MPI_Comm comm; 92c4762a1bSJed Brown 93c4762a1bSJed Brown PetscFunctionBeginUser; 949566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 959566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, dim, simplex, "output_", -1, &fe)); 969566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe, "output")); 979566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe)); 989566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 999566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 100c4762a1bSJed Brown PetscFunctionReturn(0); 101c4762a1bSJed Brown } 102c4762a1bSJed Brown 103*9371c9d4SSatish Balay static PetscErrorCode CreateSubdomainMesh(DM dm, DMLabel *domLabel, DM *subdm, AppCtx *user) { 104c4762a1bSJed Brown DMLabel label; 10530602db0SMatthew G. Knepley PetscBool simplex; 106c4762a1bSJed Brown PetscInt dim, cStart, cEnd, c; 107c4762a1bSJed Brown 108c4762a1bSJed Brown PetscFunctionBeginUser; 1099566063dSJacob Faibussowitsch PetscCall(DMPlexIsSimplex(dm, &simplex)); 1109566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 1119566063dSJacob Faibussowitsch PetscCall(DMLabelCreate(PETSC_COMM_SELF, "subdomain", &label)); 1129566063dSJacob Faibussowitsch for (c = cStart + (cEnd - cStart) / 2; c < cEnd; ++c) PetscCall(DMLabelSetValue(label, c, 1)); 1139566063dSJacob Faibussowitsch PetscCall(DMPlexFilter(dm, label, 1, subdm)); 1149566063dSJacob Faibussowitsch PetscCall(DMGetDimension(*subdm, &dim)); 1159566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(*subdm, dim, simplex, user)); 1169566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)*subdm, "subdomain")); 1179566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*subdm, NULL, "-sub_dm_view")); 118c4762a1bSJed Brown if (domLabel) *domLabel = label; 1199566063dSJacob Faibussowitsch else PetscCall(DMLabelDestroy(&label)); 120c4762a1bSJed Brown PetscFunctionReturn(0); 121c4762a1bSJed Brown } 122c4762a1bSJed Brown 123*9371c9d4SSatish Balay static PetscErrorCode CreateBoundaryMesh(DM dm, DMLabel *bdLabel, DM *subdm, AppCtx *user) { 124c4762a1bSJed Brown DMLabel label; 12530602db0SMatthew G. Knepley PetscBool simplex; 126c4762a1bSJed Brown PetscInt dim; 127c4762a1bSJed Brown 128c4762a1bSJed Brown PetscFunctionBeginUser; 1299566063dSJacob Faibussowitsch PetscCall(DMPlexIsSimplex(dm, &simplex)); 1309566063dSJacob Faibussowitsch PetscCall(DMLabelCreate(PETSC_COMM_SELF, "sub", &label)); 1319566063dSJacob Faibussowitsch PetscCall(DMPlexMarkBoundaryFaces(dm, 1, label)); 1329566063dSJacob Faibussowitsch PetscCall(DMPlexLabelComplete(dm, label)); 1339566063dSJacob Faibussowitsch PetscCall(DMPlexCreateSubmesh(dm, label, 1, PETSC_TRUE, subdm)); 1349566063dSJacob Faibussowitsch PetscCall(DMGetDimension(*subdm, &dim)); 1359566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(*subdm, dim, simplex, user)); 1369566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)*subdm, "boundary")); 1379566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*subdm, NULL, "-sub_dm_view")); 138c4762a1bSJed Brown if (bdLabel) *bdLabel = label; 1399566063dSJacob Faibussowitsch else PetscCall(DMLabelDestroy(&label)); 140c4762a1bSJed Brown PetscFunctionReturn(0); 141c4762a1bSJed Brown } 142c4762a1bSJed Brown 143*9371c9d4SSatish Balay static PetscErrorCode CreateAuxiliaryVec(DM dm, DM *auxdm, Vec *la, AppCtx *user) { 144c4762a1bSJed Brown PetscErrorCode (**afuncs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *); 14530602db0SMatthew G. Knepley PetscBool simplex; 146c4762a1bSJed Brown PetscInt dim, Nf, f; 147c4762a1bSJed Brown 148c4762a1bSJed Brown PetscFunctionBeginUser; 1499566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 1509566063dSJacob Faibussowitsch PetscCall(DMPlexIsSimplex(dm, &simplex)); 1519566063dSJacob Faibussowitsch PetscCall(DMGetNumFields(dm, &Nf)); 1529566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Nf, &afuncs)); 153c4762a1bSJed Brown for (f = 0; f < Nf; ++f) afuncs[f] = linear; 1549566063dSJacob Faibussowitsch PetscCall(DMClone(dm, auxdm)); 1559566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(*auxdm, dim, simplex, user)); 1569566063dSJacob Faibussowitsch PetscCall(DMCreateLocalVector(*auxdm, la)); 1579566063dSJacob Faibussowitsch PetscCall(DMProjectFunctionLocal(dm, 0.0, afuncs, NULL, INSERT_VALUES, *la)); 1589566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(*la, NULL, "-local_aux_view")); 1599566063dSJacob Faibussowitsch PetscCall(PetscFree(afuncs)); 160c4762a1bSJed Brown PetscFunctionReturn(0); 161c4762a1bSJed Brown } 162c4762a1bSJed Brown 163*9371c9d4SSatish Balay static PetscErrorCode TestFunctionProjection(DM dm, DM dmAux, DMLabel label, Vec la, const char name[], AppCtx *user) { 164c4762a1bSJed Brown PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *); 165c4762a1bSJed Brown Vec x, lx; 166c4762a1bSJed Brown PetscInt Nf, f; 167c4762a1bSJed Brown PetscInt val[1] = {1}; 168c4762a1bSJed Brown char lname[PETSC_MAX_PATH_LEN]; 169c4762a1bSJed Brown 170c4762a1bSJed Brown PetscFunctionBeginUser; 1719566063dSJacob Faibussowitsch if (dmAux) PetscCall(DMSetAuxiliaryVec(dm, NULL, 0, 0, la)); 1729566063dSJacob Faibussowitsch PetscCall(DMGetNumFields(dm, &Nf)); 1739566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Nf, &funcs)); 174c4762a1bSJed Brown for (f = 0; f < Nf; ++f) funcs[f] = linear; 1759566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &x)); 1769566063dSJacob Faibussowitsch PetscCall(PetscStrcpy(lname, "Function ")); 1779566063dSJacob Faibussowitsch PetscCall(PetscStrcat(lname, name)); 1789566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)x, lname)); 1799566063dSJacob Faibussowitsch if (!label) PetscCall(DMProjectFunction(dm, 0.0, funcs, NULL, INSERT_VALUES, x)); 1809566063dSJacob Faibussowitsch else PetscCall(DMProjectFunctionLabel(dm, 0.0, label, 1, val, 0, NULL, funcs, NULL, INSERT_VALUES, x)); 1819566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(x, NULL, "-func_view")); 1829566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &x)); 1839566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &lx)); 1849566063dSJacob Faibussowitsch PetscCall(PetscStrcpy(lname, "Local Function ")); 1859566063dSJacob Faibussowitsch PetscCall(PetscStrcat(lname, name)); 1869566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)lx, lname)); 1879566063dSJacob Faibussowitsch if (!label) PetscCall(DMProjectFunctionLocal(dm, 0.0, funcs, NULL, INSERT_VALUES, lx)); 1889566063dSJacob Faibussowitsch else PetscCall(DMProjectFunctionLabelLocal(dm, 0.0, label, 1, val, 0, NULL, funcs, NULL, INSERT_VALUES, lx)); 1899566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(lx, NULL, "-local_func_view")); 1909566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &lx)); 1919566063dSJacob Faibussowitsch PetscCall(PetscFree(funcs)); 1929566063dSJacob Faibussowitsch if (dmAux) PetscCall(DMSetAuxiliaryVec(dm, NULL, 0, 0, NULL)); 193c4762a1bSJed Brown PetscFunctionReturn(0); 194c4762a1bSJed Brown } 195c4762a1bSJed Brown 196*9371c9d4SSatish Balay static PetscErrorCode TestFieldProjection(DM dm, DM dmAux, DMLabel label, Vec la, const char name[], AppCtx *user) { 197c4762a1bSJed Brown PetscErrorCode (**afuncs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *); 198*9371c9d4SSatish Balay void (**funcs)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]); 199c4762a1bSJed Brown Vec lx, lu; 200c4762a1bSJed Brown PetscInt Nf, f; 201c4762a1bSJed Brown PetscInt val[1] = {1}; 202c4762a1bSJed Brown char lname[PETSC_MAX_PATH_LEN]; 203c4762a1bSJed Brown 204c4762a1bSJed Brown PetscFunctionBeginUser; 2059566063dSJacob Faibussowitsch if (dmAux) PetscCall(DMSetAuxiliaryVec(dm, NULL, 0, 0, la)); 2069566063dSJacob Faibussowitsch PetscCall(DMGetNumFields(dm, &Nf)); 2079566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(Nf, &funcs, Nf, &afuncs)); 208c4762a1bSJed Brown for (f = 0; f < Nf; ++f) afuncs[f] = linear; 209c4762a1bSJed Brown funcs[0] = linear_vector; 210c4762a1bSJed Brown funcs[1] = linear_scalar; 2119566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &lu)); 2129566063dSJacob Faibussowitsch PetscCall(PetscStrcpy(lname, "Local Field Input ")); 2139566063dSJacob Faibussowitsch PetscCall(PetscStrcat(lname, name)); 2149566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)lu, lname)); 2159566063dSJacob Faibussowitsch if (!label) PetscCall(DMProjectFunctionLocal(dm, 0.0, afuncs, NULL, INSERT_VALUES, lu)); 2169566063dSJacob Faibussowitsch else PetscCall(DMProjectFunctionLabelLocal(dm, 0.0, label, 1, val, 0, NULL, afuncs, NULL, INSERT_VALUES, lu)); 2179566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(lu, NULL, "-local_input_view")); 2189566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &lx)); 2199566063dSJacob Faibussowitsch PetscCall(PetscStrcpy(lname, "Local Field ")); 2209566063dSJacob Faibussowitsch PetscCall(PetscStrcat(lname, name)); 2219566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)lx, lname)); 2229566063dSJacob Faibussowitsch if (!label) PetscCall(DMProjectFieldLocal(dm, 0.0, lu, funcs, INSERT_VALUES, lx)); 2239566063dSJacob Faibussowitsch else PetscCall(DMProjectFieldLabelLocal(dm, 0.0, label, 1, val, 0, NULL, lu, funcs, INSERT_VALUES, lx)); 2249566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(lx, NULL, "-local_field_view")); 2259566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &lx)); 2269566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &lu)); 2279566063dSJacob Faibussowitsch PetscCall(PetscFree2(funcs, afuncs)); 2289566063dSJacob Faibussowitsch if (dmAux) PetscCall(DMSetAuxiliaryVec(dm, NULL, 0, 0, NULL)); 229c4762a1bSJed Brown PetscFunctionReturn(0); 230c4762a1bSJed Brown } 231c4762a1bSJed Brown 232*9371c9d4SSatish Balay static PetscErrorCode TestFieldProjectionMultiple(DM dm, DM dmIn, DM dmAux, DMLabel label, Vec la, const char name[], AppCtx *user) { 233c4762a1bSJed Brown PetscErrorCode (**afuncs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *); 234*9371c9d4SSatish Balay void (**funcs)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]); 235c4762a1bSJed Brown Vec lx, lu; 236c4762a1bSJed Brown PetscInt Nf, NfIn; 237c4762a1bSJed Brown PetscInt val[1] = {1}; 238c4762a1bSJed Brown char lname[PETSC_MAX_PATH_LEN]; 239c4762a1bSJed Brown 240c4762a1bSJed Brown PetscFunctionBeginUser; 2419566063dSJacob Faibussowitsch if (dmAux) PetscCall(DMSetAuxiliaryVec(dm, NULL, 0, 0, la)); 2429566063dSJacob Faibussowitsch PetscCall(DMGetNumFields(dm, &Nf)); 2439566063dSJacob Faibussowitsch PetscCall(DMGetNumFields(dmIn, &NfIn)); 2449566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(Nf, &funcs, NfIn, &afuncs)); 245c4762a1bSJed Brown funcs[0] = divergence_sq; 246c4762a1bSJed Brown afuncs[0] = linear2; 247c4762a1bSJed Brown afuncs[1] = linear; 2489566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dmIn, &lu)); 2499566063dSJacob Faibussowitsch PetscCall(PetscStrcpy(lname, "Local MultiField Input ")); 2509566063dSJacob Faibussowitsch PetscCall(PetscStrcat(lname, name)); 2519566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)lu, lname)); 2529566063dSJacob Faibussowitsch if (!label) PetscCall(DMProjectFunctionLocal(dmIn, 0.0, afuncs, NULL, INSERT_VALUES, lu)); 2539566063dSJacob Faibussowitsch else PetscCall(DMProjectFunctionLabelLocal(dmIn, 0.0, label, 1, val, 0, NULL, afuncs, NULL, INSERT_VALUES, lu)); 2549566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(lu, NULL, "-local_input_view")); 2559566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &lx)); 2569566063dSJacob Faibussowitsch PetscCall(PetscStrcpy(lname, "Local MultiField ")); 2579566063dSJacob Faibussowitsch PetscCall(PetscStrcat(lname, name)); 2589566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)lx, lname)); 2599566063dSJacob Faibussowitsch if (!label) PetscCall(DMProjectFieldLocal(dm, 0.0, lu, funcs, INSERT_VALUES, lx)); 2609566063dSJacob Faibussowitsch else PetscCall(DMProjectFieldLabelLocal(dm, 0.0, label, 1, val, 0, NULL, lu, funcs, INSERT_VALUES, lx)); 2619566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(lx, NULL, "-local_field_view")); 2629566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &lx)); 2639566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dmIn, &lu)); 2649566063dSJacob Faibussowitsch PetscCall(PetscFree2(funcs, afuncs)); 2659566063dSJacob Faibussowitsch if (dmAux) PetscCall(DMSetAuxiliaryVec(dm, NULL, 0, 0, NULL)); 266c4762a1bSJed Brown PetscFunctionReturn(0); 267c4762a1bSJed Brown } 268c4762a1bSJed Brown 269*9371c9d4SSatish Balay int main(int argc, char **argv) { 270c4762a1bSJed Brown DM dm, subdm, auxdm; 271c4762a1bSJed Brown Vec la; 27230602db0SMatthew G. Knepley PetscInt dim; 27330602db0SMatthew G. Knepley PetscBool simplex; 274c4762a1bSJed Brown AppCtx user; 275c4762a1bSJed Brown 276327415f7SBarry Smith PetscFunctionBeginUser; 2779566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 2789566063dSJacob Faibussowitsch PetscCall(ProcessOptions(&user)); 2799566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 2809566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 2819566063dSJacob Faibussowitsch PetscCall(DMPlexIsSimplex(dm, &simplex)); 2829566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(dm, dim, simplex, &user)); 283c4762a1bSJed Brown /* Volumetric Mesh Projection */ 284c4762a1bSJed Brown if (!user.multifield) { 2859566063dSJacob Faibussowitsch PetscCall(TestFunctionProjection(dm, NULL, NULL, NULL, "Volumetric Primary", &user)); 2869566063dSJacob Faibussowitsch PetscCall(TestFieldProjection(dm, NULL, NULL, NULL, "Volumetric Primary", &user)); 287c4762a1bSJed Brown } else { 288c4762a1bSJed Brown DM dmOut; 289c4762a1bSJed Brown 2909566063dSJacob Faibussowitsch PetscCall(DMClone(dm, &dmOut)); 2919566063dSJacob Faibussowitsch PetscCall(SetupOutputDiscretization(dmOut, dim, simplex, &user)); 2929566063dSJacob Faibussowitsch PetscCall(TestFieldProjectionMultiple(dmOut, dm, NULL, NULL, NULL, "Volumetric Primary", &user)); 2939566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmOut)); 294c4762a1bSJed Brown } 295c4762a1bSJed Brown if (user.auxfield) { 296c4762a1bSJed Brown /* Volumetric Mesh Projection with Volumetric Data */ 2979566063dSJacob Faibussowitsch PetscCall(CreateAuxiliaryVec(dm, &auxdm, &la, &user)); 2989566063dSJacob Faibussowitsch PetscCall(TestFunctionProjection(dm, auxdm, NULL, la, "Volumetric Primary and Volumetric Auxiliary", &user)); 2999566063dSJacob Faibussowitsch PetscCall(TestFieldProjection(dm, auxdm, NULL, la, "Volumetric Primary and Volumetric Auxiliary", &user)); 3009566063dSJacob Faibussowitsch PetscCall(VecDestroy(&la)); 301c4762a1bSJed Brown /* Update of Volumetric Auxiliary Data with primary Volumetric Data */ 3029566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &la)); 3039566063dSJacob Faibussowitsch PetscCall(VecSet(la, 1.0)); 3049566063dSJacob Faibussowitsch PetscCall(TestFieldProjection(auxdm, dm, NULL, la, "Volumetric Auxiliary Update with Volumetric Primary", &user)); 3059566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &la)); 3069566063dSJacob Faibussowitsch PetscCall(DMDestroy(&auxdm)); 307c4762a1bSJed Brown } 308c4762a1bSJed Brown if (user.subdomain) { 309c4762a1bSJed Brown DMLabel domLabel; 310c4762a1bSJed Brown 311c4762a1bSJed Brown /* Subdomain Mesh Projection */ 3129566063dSJacob Faibussowitsch PetscCall(CreateSubdomainMesh(dm, &domLabel, &subdm, &user)); 3139566063dSJacob Faibussowitsch PetscCall(TestFunctionProjection(subdm, NULL, NULL, NULL, "Subdomain Primary", &user)); 3149566063dSJacob Faibussowitsch PetscCall(TestFieldProjection(subdm, NULL, NULL, NULL, "Subdomain Primary", &user)); 315c4762a1bSJed Brown if (user.auxfield) { 316c4762a1bSJed Brown /* Subdomain Mesh Projection with Subdomain Data */ 3179566063dSJacob Faibussowitsch PetscCall(CreateAuxiliaryVec(subdm, &auxdm, &la, &user)); 3189566063dSJacob Faibussowitsch PetscCall(TestFunctionProjection(subdm, auxdm, NULL, la, "Subdomain Primary and Subdomain Auxiliary", &user)); 3199566063dSJacob Faibussowitsch PetscCall(TestFieldProjection(subdm, auxdm, NULL, la, "Subdomain Primary and Subdomain Auxiliary", &user)); 3209566063dSJacob Faibussowitsch PetscCall(VecDestroy(&la)); 3219566063dSJacob Faibussowitsch PetscCall(DMDestroy(&auxdm)); 322c4762a1bSJed Brown /* Subdomain Mesh Projection with Volumetric Data */ 3239566063dSJacob Faibussowitsch PetscCall(CreateAuxiliaryVec(dm, &auxdm, &la, &user)); 3249566063dSJacob Faibussowitsch PetscCall(TestFunctionProjection(subdm, auxdm, NULL, la, "Subdomain Primary and Volumetric Auxiliary", &user)); 3259566063dSJacob Faibussowitsch PetscCall(TestFieldProjection(subdm, auxdm, NULL, la, "Subdomain Primary and Volumetric Auxiliary", &user)); 3269566063dSJacob Faibussowitsch PetscCall(VecDestroy(&la)); 3279566063dSJacob Faibussowitsch PetscCall(DMDestroy(&auxdm)); 328c4762a1bSJed Brown /* Volumetric Mesh Projection with Subdomain Data */ 3299566063dSJacob Faibussowitsch PetscCall(CreateAuxiliaryVec(subdm, &auxdm, &la, &user)); 3309566063dSJacob Faibussowitsch PetscCall(TestFunctionProjection(subdm, auxdm, domLabel, la, "Volumetric Primary and Subdomain Auxiliary", &user)); 3319566063dSJacob Faibussowitsch PetscCall(TestFieldProjection(subdm, auxdm, domLabel, la, "Volumetric Primary and Subdomain Auxiliary", &user)); 3329566063dSJacob Faibussowitsch PetscCall(VecDestroy(&la)); 3339566063dSJacob Faibussowitsch PetscCall(DMDestroy(&auxdm)); 334c4762a1bSJed Brown } 3359566063dSJacob Faibussowitsch PetscCall(DMDestroy(&subdm)); 3369566063dSJacob Faibussowitsch PetscCall(DMLabelDestroy(&domLabel)); 337c4762a1bSJed Brown } 338c4762a1bSJed Brown if (user.submesh) { 339c4762a1bSJed Brown DMLabel bdLabel; 340c4762a1bSJed Brown 341c4762a1bSJed Brown /* Boundary Mesh Projection */ 3429566063dSJacob Faibussowitsch PetscCall(CreateBoundaryMesh(dm, &bdLabel, &subdm, &user)); 3439566063dSJacob Faibussowitsch PetscCall(TestFunctionProjection(subdm, NULL, NULL, NULL, "Boundary Primary", &user)); 3449566063dSJacob Faibussowitsch PetscCall(TestFieldProjection(subdm, NULL, NULL, NULL, "Boundary Primary", &user)); 345c4762a1bSJed Brown if (user.auxfield) { 346c4762a1bSJed Brown /* Boundary Mesh Projection with Boundary Data */ 3479566063dSJacob Faibussowitsch PetscCall(CreateAuxiliaryVec(subdm, &auxdm, &la, &user)); 3489566063dSJacob Faibussowitsch PetscCall(TestFunctionProjection(subdm, auxdm, NULL, la, "Boundary Primary and Boundary Auxiliary", &user)); 3499566063dSJacob Faibussowitsch PetscCall(TestFieldProjection(subdm, auxdm, NULL, la, "Boundary Primary and Boundary Auxiliary", &user)); 3509566063dSJacob Faibussowitsch PetscCall(VecDestroy(&la)); 3519566063dSJacob Faibussowitsch PetscCall(DMDestroy(&auxdm)); 352c4762a1bSJed Brown /* Volumetric Mesh Projection with Boundary Data */ 3539566063dSJacob Faibussowitsch PetscCall(CreateAuxiliaryVec(subdm, &auxdm, &la, &user)); 3549566063dSJacob Faibussowitsch PetscCall(TestFunctionProjection(dm, auxdm, bdLabel, la, "Volumetric Primary and Boundary Auxiliary", &user)); 3559566063dSJacob Faibussowitsch PetscCall(TestFieldProjection(dm, auxdm, bdLabel, la, "Volumetric Primary and Boundary Auxiliary", &user)); 3569566063dSJacob Faibussowitsch PetscCall(VecDestroy(&la)); 3579566063dSJacob Faibussowitsch PetscCall(DMDestroy(&auxdm)); 358c4762a1bSJed Brown } 3599566063dSJacob Faibussowitsch PetscCall(DMLabelDestroy(&bdLabel)); 3609566063dSJacob Faibussowitsch PetscCall(DMDestroy(&subdm)); 361c4762a1bSJed Brown } 3629566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 3639566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 364b122ec5aSJacob Faibussowitsch return 0; 365c4762a1bSJed Brown } 366c4762a1bSJed Brown 367c4762a1bSJed Brown /*TEST 368c4762a1bSJed Brown 369c4762a1bSJed Brown test: 370c4762a1bSJed Brown suffix: 0 371c4762a1bSJed Brown requires: triangle 37230602db0SMatthew G. Knepley args: -dm_plex_box_faces 1,1 -func_view -local_func_view -local_input_view -local_field_view 373c4762a1bSJed Brown test: 374c4762a1bSJed Brown suffix: mf_0 375c4762a1bSJed Brown requires: triangle 37630602db0SMatthew G. Knepley args: -dm_plex_box_faces 1,1 -velocity_petscspace_degree 1 -velocity_petscfe_default_quadrature_order 2 \ 377c4762a1bSJed Brown -pressure_petscspace_degree 2 -pressure_petscfe_default_quadrature_order 2 \ 378c4762a1bSJed Brown -multifield -output_petscspace_degree 1 -output_petscfe_default_quadrature_order 2 \ 379c4762a1bSJed Brown -local_input_view -local_field_view 380c4762a1bSJed Brown test: 381c4762a1bSJed Brown suffix: 1 382c4762a1bSJed Brown requires: triangle 38330602db0SMatthew 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 384c4762a1bSJed Brown test: 385c4762a1bSJed Brown suffix: 2 386c4762a1bSJed Brown requires: triangle 38730602db0SMatthew 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 388c4762a1bSJed Brown 389c4762a1bSJed Brown TEST*/ 390c4762a1bSJed Brown 391c4762a1bSJed Brown /* 392c4762a1bSJed Brown Post-processing wants to project a function of the fields into some FE space 393c4762a1bSJed Brown - This is DMProjectField() 394c4762a1bSJed Brown - What about changing the number of components of the output, like displacement to stress? Aux vars 395c4762a1bSJed Brown 396c4762a1bSJed Brown Update of state variables 397c4762a1bSJed Brown - This is DMProjectField(), but solution must be the aux var 398c4762a1bSJed Brown */ 399