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); 71c4762a1bSJed Brown ierr = PetscOptionsBool("-multifield", "Flag for trying different numbers of input/output fields", "ex23.c", options->multifield, &options->multifield, NULL);CHKERRQ(ierr); 72c4762a1bSJed Brown ierr = PetscOptionsBool("-subdomain", "Flag for trying volumetric submesh", "ex23.c", options->subdomain, &options->subdomain, NULL);CHKERRQ(ierr); 73c4762a1bSJed Brown ierr = PetscOptionsBool("-submesh", "Flag for trying boundary submesh", "ex23.c", options->submesh, &options->submesh, NULL);CHKERRQ(ierr); 74c4762a1bSJed Brown ierr = PetscOptionsBool("-auxfield", "Flag for trying auxiliary fields", "ex23.c", options->auxfield, &options->auxfield, NULL);CHKERRQ(ierr); 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 PetscErrorCode ierr; 82c4762a1bSJed Brown 83c4762a1bSJed Brown PetscFunctionBegin; 84*30602db0SMatthew G. Knepley ierr = DMCreate(comm, dm);CHKERRQ(ierr); 85*30602db0SMatthew G. Knepley ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 86c4762a1bSJed Brown ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 87c4762a1bSJed Brown ierr = DMViewFromOptions(*dm, NULL, "-orig_dm_view");CHKERRQ(ierr); 88c4762a1bSJed Brown PetscFunctionReturn(0); 89c4762a1bSJed Brown } 90c4762a1bSJed Brown 91c4762a1bSJed Brown static PetscErrorCode SetupDiscretization(DM dm, PetscInt dim, PetscBool simplex, AppCtx *user) 92c4762a1bSJed Brown { 93c4762a1bSJed Brown PetscFE fe; 94c4762a1bSJed Brown MPI_Comm comm; 95c4762a1bSJed Brown PetscErrorCode ierr; 96c4762a1bSJed Brown 97c4762a1bSJed Brown PetscFunctionBeginUser; 98c4762a1bSJed Brown ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 99c4762a1bSJed Brown ierr = PetscFECreateDefault(comm, dim, dim, simplex, "velocity_", -1, &fe);CHKERRQ(ierr); 100c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) fe, "velocity");CHKERRQ(ierr); 101c4762a1bSJed Brown ierr = DMSetField(dm, 0, NULL, (PetscObject) fe);CHKERRQ(ierr); 102c4762a1bSJed Brown ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 103c4762a1bSJed Brown ierr = PetscFECreateDefault(comm, dim, 1, simplex, "pressure_", -1, &fe);CHKERRQ(ierr); 104c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) fe, "pressure");CHKERRQ(ierr); 105c4762a1bSJed Brown ierr = DMSetField(dm, 1, NULL, (PetscObject) fe);CHKERRQ(ierr); 106c4762a1bSJed Brown ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 107c4762a1bSJed Brown ierr = DMCreateDS(dm);CHKERRQ(ierr); 108c4762a1bSJed Brown PetscFunctionReturn(0); 109c4762a1bSJed Brown } 110c4762a1bSJed Brown 111c4762a1bSJed Brown static PetscErrorCode SetupOutputDiscretization(DM dm, PetscInt dim, PetscBool simplex, AppCtx *user) 112c4762a1bSJed Brown { 113c4762a1bSJed Brown PetscFE fe; 114c4762a1bSJed Brown MPI_Comm comm; 115c4762a1bSJed Brown PetscErrorCode ierr; 116c4762a1bSJed Brown 117c4762a1bSJed Brown PetscFunctionBeginUser; 118c4762a1bSJed Brown ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 119c4762a1bSJed Brown ierr = PetscFECreateDefault(comm, dim, dim, simplex, "output_", -1, &fe);CHKERRQ(ierr); 120c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) fe, "output");CHKERRQ(ierr); 121c4762a1bSJed Brown ierr = DMSetField(dm, 0, NULL, (PetscObject) fe);CHKERRQ(ierr); 122c4762a1bSJed Brown ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 123c4762a1bSJed Brown ierr = DMCreateDS(dm);CHKERRQ(ierr); 124c4762a1bSJed Brown PetscFunctionReturn(0); 125c4762a1bSJed Brown } 126c4762a1bSJed Brown 127c4762a1bSJed Brown static PetscErrorCode CreateSubdomainMesh(DM dm, DMLabel *domLabel, DM *subdm, AppCtx *user) 128c4762a1bSJed Brown { 129c4762a1bSJed Brown DMLabel label; 130*30602db0SMatthew G. Knepley PetscBool simplex; 131c4762a1bSJed Brown PetscInt dim, cStart, cEnd, c; 132c4762a1bSJed Brown PetscErrorCode ierr; 133c4762a1bSJed Brown 134c4762a1bSJed Brown PetscFunctionBeginUser; 135*30602db0SMatthew G. Knepley ierr = DMPlexIsSimplex(dm, &simplex);CHKERRQ(ierr); 136c4762a1bSJed Brown ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 137c4762a1bSJed Brown ierr = DMLabelCreate(PETSC_COMM_SELF, "subdomain", &label);CHKERRQ(ierr); 138c4762a1bSJed Brown for (c = cStart + (cEnd-cStart)/2; c < cEnd; ++c) {ierr = DMLabelSetValue(label, c, 1);CHKERRQ(ierr);} 139c4762a1bSJed Brown ierr = DMPlexFilter(dm, label, 1, subdm);CHKERRQ(ierr); 140c4762a1bSJed Brown ierr = DMGetDimension(*subdm, &dim);CHKERRQ(ierr); 141*30602db0SMatthew G. Knepley ierr = SetupDiscretization(*subdm, dim, simplex, user);CHKERRQ(ierr); 142c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) *subdm, "subdomain");CHKERRQ(ierr); 143c4762a1bSJed Brown ierr = DMViewFromOptions(*subdm, NULL, "-sub_dm_view");CHKERRQ(ierr); 144c4762a1bSJed Brown if (domLabel) *domLabel = label; 145c4762a1bSJed Brown else {ierr = DMLabelDestroy(&label);CHKERRQ(ierr);} 146c4762a1bSJed Brown PetscFunctionReturn(0); 147c4762a1bSJed Brown } 148c4762a1bSJed Brown 149c4762a1bSJed Brown static PetscErrorCode CreateBoundaryMesh(DM dm, DMLabel *bdLabel, DM *subdm, AppCtx *user) 150c4762a1bSJed Brown { 151c4762a1bSJed Brown DMLabel label; 152*30602db0SMatthew G. Knepley PetscBool simplex; 153c4762a1bSJed Brown PetscInt dim; 154c4762a1bSJed Brown PetscErrorCode ierr; 155c4762a1bSJed Brown 156c4762a1bSJed Brown PetscFunctionBeginUser; 157*30602db0SMatthew G. Knepley ierr = DMPlexIsSimplex(dm, &simplex);CHKERRQ(ierr); 158c4762a1bSJed Brown ierr = DMLabelCreate(PETSC_COMM_SELF, "sub", &label);CHKERRQ(ierr); 159c4762a1bSJed Brown ierr = DMPlexMarkBoundaryFaces(dm, 1, label);CHKERRQ(ierr); 160c4762a1bSJed Brown ierr = DMPlexLabelComplete(dm, label);CHKERRQ(ierr); 161c4762a1bSJed Brown ierr = DMPlexCreateSubmesh(dm, label, 1, PETSC_TRUE, subdm);CHKERRQ(ierr); 162c4762a1bSJed Brown ierr = DMGetDimension(*subdm, &dim);CHKERRQ(ierr); 163*30602db0SMatthew G. Knepley ierr = SetupDiscretization(*subdm, dim, simplex, user);CHKERRQ(ierr); 164c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) *subdm, "boundary");CHKERRQ(ierr); 165c4762a1bSJed Brown ierr = DMViewFromOptions(*subdm, NULL, "-sub_dm_view");CHKERRQ(ierr); 166c4762a1bSJed Brown if (bdLabel) *bdLabel = label; 167c4762a1bSJed Brown else {ierr = DMLabelDestroy(&label);CHKERRQ(ierr);} 168c4762a1bSJed Brown PetscFunctionReturn(0); 169c4762a1bSJed Brown } 170c4762a1bSJed Brown 1719a2a23afSMatthew G. Knepley static PetscErrorCode CreateAuxiliaryVec(DM dm, DM *auxdm, Vec *la, AppCtx *user) 172c4762a1bSJed Brown { 173c4762a1bSJed Brown PetscErrorCode (**afuncs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *); 174*30602db0SMatthew G. Knepley PetscBool simplex; 175c4762a1bSJed Brown PetscInt dim, Nf, f; 176c4762a1bSJed Brown PetscErrorCode ierr; 177c4762a1bSJed Brown 178c4762a1bSJed Brown PetscFunctionBeginUser; 179c4762a1bSJed Brown ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 180*30602db0SMatthew G. Knepley ierr = DMPlexIsSimplex(dm, &simplex);CHKERRQ(ierr); 181c4762a1bSJed Brown ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr); 182c4762a1bSJed Brown ierr = PetscMalloc1(Nf, &afuncs);CHKERRQ(ierr); 183c4762a1bSJed Brown for (f = 0; f < Nf; ++f) afuncs[f] = linear; 184c4762a1bSJed Brown ierr = DMClone(dm, auxdm);CHKERRQ(ierr); 185*30602db0SMatthew G. Knepley ierr = SetupDiscretization(*auxdm, dim, simplex, user);CHKERRQ(ierr); 186c4762a1bSJed Brown ierr = DMCreateLocalVector(*auxdm, la);CHKERRQ(ierr); 187c4762a1bSJed Brown ierr = DMProjectFunctionLocal(dm, 0.0, afuncs, NULL, INSERT_VALUES, *la);CHKERRQ(ierr); 188c4762a1bSJed Brown ierr = VecViewFromOptions(*la, NULL, "-local_aux_view");CHKERRQ(ierr); 189c4762a1bSJed Brown ierr = PetscFree(afuncs);CHKERRQ(ierr); 190c4762a1bSJed Brown PetscFunctionReturn(0); 191c4762a1bSJed Brown } 192c4762a1bSJed Brown 193c4762a1bSJed Brown static PetscErrorCode TestFunctionProjection(DM dm, DM dmAux, DMLabel label, Vec la, const char name[], AppCtx *user) 194c4762a1bSJed Brown { 195c4762a1bSJed Brown PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *); 196c4762a1bSJed Brown Vec x, lx; 197c4762a1bSJed Brown PetscInt Nf, f; 198c4762a1bSJed Brown PetscInt val[1] = {1}; 199c4762a1bSJed Brown char lname[PETSC_MAX_PATH_LEN]; 200c4762a1bSJed Brown PetscErrorCode ierr; 201c4762a1bSJed Brown 202c4762a1bSJed Brown PetscFunctionBeginUser; 2039a2a23afSMatthew G. Knepley if (dmAux) {ierr = DMSetAuxiliaryVec(dm, NULL, 0, la);CHKERRQ(ierr);} 204c4762a1bSJed Brown ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr); 205c4762a1bSJed Brown ierr = PetscMalloc1(Nf, &funcs);CHKERRQ(ierr); 206c4762a1bSJed Brown for (f = 0; f < Nf; ++f) funcs[f] = linear; 207c4762a1bSJed Brown ierr = DMGetGlobalVector(dm, &x);CHKERRQ(ierr); 208c4762a1bSJed Brown ierr = PetscStrcpy(lname, "Function ");CHKERRQ(ierr); 209c4762a1bSJed Brown ierr = PetscStrcat(lname, name);CHKERRQ(ierr); 210c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) x, lname);CHKERRQ(ierr); 211c4762a1bSJed Brown if (!label) {ierr = DMProjectFunction(dm, 0.0, funcs, NULL, INSERT_VALUES, x);CHKERRQ(ierr);} 212c4762a1bSJed Brown else {ierr = DMProjectFunctionLabel(dm, 0.0, label, 1, val, 0, NULL, funcs, NULL, INSERT_VALUES, x);CHKERRQ(ierr);} 213c4762a1bSJed Brown ierr = VecViewFromOptions(x, NULL, "-func_view");CHKERRQ(ierr); 214c4762a1bSJed Brown ierr = DMRestoreGlobalVector(dm, &x);CHKERRQ(ierr); 215c4762a1bSJed Brown ierr = DMGetLocalVector(dm, &lx);CHKERRQ(ierr); 216c4762a1bSJed Brown ierr = PetscStrcpy(lname, "Local Function ");CHKERRQ(ierr); 217c4762a1bSJed Brown ierr = PetscStrcat(lname, name);CHKERRQ(ierr); 218c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) lx, lname);CHKERRQ(ierr); 219c4762a1bSJed Brown if (!label) {ierr = DMProjectFunctionLocal(dm, 0.0, funcs, NULL, INSERT_VALUES, lx);CHKERRQ(ierr);} 220c4762a1bSJed Brown else {ierr = DMProjectFunctionLabelLocal(dm, 0.0, label, 1, val, 0, NULL, funcs, NULL, INSERT_VALUES, lx);CHKERRQ(ierr);} 221c4762a1bSJed Brown ierr = VecViewFromOptions(lx, NULL, "-local_func_view");CHKERRQ(ierr); 222c4762a1bSJed Brown ierr = DMRestoreLocalVector(dm, &lx);CHKERRQ(ierr); 223c4762a1bSJed Brown ierr = PetscFree(funcs);CHKERRQ(ierr); 2249a2a23afSMatthew G. Knepley if (dmAux) {ierr = DMSetAuxiliaryVec(dm, NULL, 0, NULL);CHKERRQ(ierr);} 225c4762a1bSJed Brown PetscFunctionReturn(0); 226c4762a1bSJed Brown } 227c4762a1bSJed Brown 228c4762a1bSJed Brown static PetscErrorCode TestFieldProjection(DM dm, DM dmAux, DMLabel label, Vec la, const char name[], AppCtx *user) 229c4762a1bSJed Brown { 230c4762a1bSJed Brown PetscErrorCode (**afuncs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *); 231c4762a1bSJed Brown void (**funcs)(PetscInt, PetscInt, PetscInt, 232c4762a1bSJed Brown const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 233c4762a1bSJed Brown const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 234c4762a1bSJed Brown PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]); 235c4762a1bSJed Brown Vec lx, lu; 236c4762a1bSJed Brown PetscInt Nf, f; 237c4762a1bSJed Brown PetscInt val[1] = {1}; 238c4762a1bSJed Brown char lname[PETSC_MAX_PATH_LEN]; 239c4762a1bSJed Brown PetscErrorCode ierr; 240c4762a1bSJed Brown 241c4762a1bSJed Brown PetscFunctionBeginUser; 2429a2a23afSMatthew G. Knepley if (dmAux) {ierr = DMSetAuxiliaryVec(dm, NULL, 0, la);CHKERRQ(ierr);} 243c4762a1bSJed Brown ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr); 244c4762a1bSJed Brown ierr = PetscMalloc2(Nf, &funcs, Nf, &afuncs);CHKERRQ(ierr); 245c4762a1bSJed Brown for (f = 0; f < Nf; ++f) afuncs[f] = linear; 246c4762a1bSJed Brown funcs[0] = linear_vector; 247c4762a1bSJed Brown funcs[1] = linear_scalar; 248c4762a1bSJed Brown ierr = DMGetLocalVector(dm, &lu);CHKERRQ(ierr); 249c4762a1bSJed Brown ierr = PetscStrcpy(lname, "Local Field Input ");CHKERRQ(ierr); 250c4762a1bSJed Brown ierr = PetscStrcat(lname, name);CHKERRQ(ierr); 251c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) lu, lname);CHKERRQ(ierr); 252c4762a1bSJed Brown if (!label) {ierr = DMProjectFunctionLocal(dm, 0.0, afuncs, NULL, INSERT_VALUES, lu);CHKERRQ(ierr);} 253c4762a1bSJed Brown else {ierr = DMProjectFunctionLabelLocal(dm, 0.0, label, 1, val, 0, NULL, afuncs, NULL, INSERT_VALUES, lu);CHKERRQ(ierr);} 254c4762a1bSJed Brown ierr = VecViewFromOptions(lu, NULL, "-local_input_view");CHKERRQ(ierr); 255c4762a1bSJed Brown ierr = DMGetLocalVector(dm, &lx);CHKERRQ(ierr); 256c4762a1bSJed Brown ierr = PetscStrcpy(lname, "Local Field ");CHKERRQ(ierr); 257c4762a1bSJed Brown ierr = PetscStrcat(lname, name);CHKERRQ(ierr); 258c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) lx, lname);CHKERRQ(ierr); 259c4762a1bSJed Brown if (!label) {ierr = DMProjectFieldLocal(dm, 0.0, lu, funcs, INSERT_VALUES, lx);CHKERRQ(ierr);} 260c4762a1bSJed Brown else {ierr = DMProjectFieldLabelLocal(dm, 0.0, label, 1, val, 0, NULL, lu, funcs, INSERT_VALUES, lx);CHKERRQ(ierr);} 261c4762a1bSJed Brown ierr = VecViewFromOptions(lx, NULL, "-local_field_view");CHKERRQ(ierr); 262c4762a1bSJed Brown ierr = DMRestoreLocalVector(dm, &lx);CHKERRQ(ierr); 263c4762a1bSJed Brown ierr = DMRestoreLocalVector(dm, &lu);CHKERRQ(ierr); 264c4762a1bSJed Brown ierr = PetscFree2(funcs, afuncs);CHKERRQ(ierr); 2659a2a23afSMatthew G. Knepley if (dmAux) {ierr = DMSetAuxiliaryVec(dm, NULL, 0, NULL);CHKERRQ(ierr);} 266c4762a1bSJed Brown PetscFunctionReturn(0); 267c4762a1bSJed Brown } 268c4762a1bSJed Brown 269c4762a1bSJed Brown static PetscErrorCode TestFieldProjectionMultiple(DM dm, DM dmIn, DM dmAux, DMLabel label, Vec la, const char name[], AppCtx *user) 270c4762a1bSJed Brown { 271c4762a1bSJed Brown PetscErrorCode (**afuncs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *); 272c4762a1bSJed Brown void (**funcs)(PetscInt, PetscInt, PetscInt, 273c4762a1bSJed Brown const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 274c4762a1bSJed Brown const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 275c4762a1bSJed Brown PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]); 276c4762a1bSJed Brown Vec lx, lu; 277c4762a1bSJed Brown PetscInt Nf, NfIn; 278c4762a1bSJed Brown PetscInt val[1] = {1}; 279c4762a1bSJed Brown char lname[PETSC_MAX_PATH_LEN]; 280c4762a1bSJed Brown PetscErrorCode ierr; 281c4762a1bSJed Brown 282c4762a1bSJed Brown PetscFunctionBeginUser; 2839a2a23afSMatthew G. Knepley if (dmAux) {ierr = DMSetAuxiliaryVec(dm, NULL, 0, la);CHKERRQ(ierr);} 284c4762a1bSJed Brown ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr); 285c4762a1bSJed Brown ierr = DMGetNumFields(dmIn, &NfIn);CHKERRQ(ierr); 286c4762a1bSJed Brown ierr = PetscMalloc2(Nf, &funcs, NfIn, &afuncs);CHKERRQ(ierr); 287c4762a1bSJed Brown funcs[0] = divergence_sq; 288c4762a1bSJed Brown afuncs[0] = linear2; 289c4762a1bSJed Brown afuncs[1] = linear; 290c4762a1bSJed Brown ierr = DMGetLocalVector(dmIn, &lu);CHKERRQ(ierr); 291c4762a1bSJed Brown ierr = PetscStrcpy(lname, "Local MultiField Input ");CHKERRQ(ierr); 292c4762a1bSJed Brown ierr = PetscStrcat(lname, name);CHKERRQ(ierr); 293c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) lu, lname);CHKERRQ(ierr); 294c4762a1bSJed Brown if (!label) {ierr = DMProjectFunctionLocal(dmIn, 0.0, afuncs, NULL, INSERT_VALUES, lu);CHKERRQ(ierr);} 295c4762a1bSJed Brown else {ierr = DMProjectFunctionLabelLocal(dmIn, 0.0, label, 1, val, 0, NULL, afuncs, NULL, INSERT_VALUES, lu);CHKERRQ(ierr);} 296c4762a1bSJed Brown ierr = VecViewFromOptions(lu, NULL, "-local_input_view");CHKERRQ(ierr); 297c4762a1bSJed Brown ierr = DMGetLocalVector(dm, &lx);CHKERRQ(ierr); 298c4762a1bSJed Brown ierr = PetscStrcpy(lname, "Local MultiField ");CHKERRQ(ierr); 299c4762a1bSJed Brown ierr = PetscStrcat(lname, name);CHKERRQ(ierr); 300c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) lx, lname);CHKERRQ(ierr); 301c4762a1bSJed Brown if (!label) {ierr = DMProjectFieldLocal(dm, 0.0, lu, funcs, INSERT_VALUES, lx);CHKERRQ(ierr);} 302c4762a1bSJed Brown else {ierr = DMProjectFieldLabelLocal(dm, 0.0, label, 1, val, 0, NULL, lu, funcs, INSERT_VALUES, lx);CHKERRQ(ierr);} 303c4762a1bSJed Brown ierr = VecViewFromOptions(lx, NULL, "-local_field_view");CHKERRQ(ierr); 304c4762a1bSJed Brown ierr = DMRestoreLocalVector(dm, &lx);CHKERRQ(ierr); 305c4762a1bSJed Brown ierr = DMRestoreLocalVector(dmIn, &lu);CHKERRQ(ierr); 306c4762a1bSJed Brown ierr = PetscFree2(funcs, afuncs);CHKERRQ(ierr); 3079a2a23afSMatthew G. Knepley if (dmAux) {ierr = DMSetAuxiliaryVec(dm, NULL, 0, NULL);CHKERRQ(ierr);} 308c4762a1bSJed Brown PetscFunctionReturn(0); 309c4762a1bSJed Brown } 310c4762a1bSJed Brown 311c4762a1bSJed Brown int main(int argc, char **argv) 312c4762a1bSJed Brown { 313c4762a1bSJed Brown DM dm, subdm, auxdm; 314c4762a1bSJed Brown Vec la; 315*30602db0SMatthew G. Knepley PetscInt dim; 316*30602db0SMatthew G. Knepley PetscBool simplex; 317c4762a1bSJed Brown AppCtx user; 318c4762a1bSJed Brown PetscErrorCode ierr; 319c4762a1bSJed Brown 320c4762a1bSJed Brown ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr; 321c4762a1bSJed Brown ierr = ProcessOptions(&user);CHKERRQ(ierr); 322c4762a1bSJed Brown ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr); 323*30602db0SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 324*30602db0SMatthew G. Knepley ierr = DMPlexIsSimplex(dm, &simplex);CHKERRQ(ierr); 325*30602db0SMatthew G. Knepley ierr = SetupDiscretization(dm, dim, simplex, &user);CHKERRQ(ierr); 326c4762a1bSJed Brown /* Volumetric Mesh Projection */ 327c4762a1bSJed Brown if (!user.multifield) { 328c4762a1bSJed Brown ierr = TestFunctionProjection(dm, NULL, NULL, NULL, "Volumetric Primary", &user);CHKERRQ(ierr); 329c4762a1bSJed Brown ierr = TestFieldProjection(dm, NULL, NULL, NULL, "Volumetric Primary", &user);CHKERRQ(ierr); 330c4762a1bSJed Brown } else { 331c4762a1bSJed Brown DM dmOut; 332c4762a1bSJed Brown 333c4762a1bSJed Brown ierr = DMClone(dm, &dmOut);CHKERRQ(ierr); 334*30602db0SMatthew G. Knepley ierr = SetupOutputDiscretization(dmOut, dim, simplex, &user);CHKERRQ(ierr); 335c4762a1bSJed Brown ierr = TestFieldProjectionMultiple(dmOut, dm, NULL, NULL, NULL, "Volumetric Primary", &user);CHKERRQ(ierr); 336c4762a1bSJed Brown ierr = DMDestroy(&dmOut);CHKERRQ(ierr); 337c4762a1bSJed Brown } 338c4762a1bSJed Brown if (user.auxfield) { 339c4762a1bSJed Brown /* Volumetric Mesh Projection with Volumetric Data */ 3409a2a23afSMatthew G. Knepley ierr = CreateAuxiliaryVec(dm, &auxdm, &la, &user);CHKERRQ(ierr); 341c4762a1bSJed Brown ierr = TestFunctionProjection(dm, auxdm, NULL, la, "Volumetric Primary and Volumetric Auxiliary", &user);CHKERRQ(ierr); 342c4762a1bSJed Brown ierr = TestFieldProjection(dm, auxdm, NULL, la, "Volumetric Primary and Volumetric Auxiliary", &user);CHKERRQ(ierr); 343c4762a1bSJed Brown ierr = VecDestroy(&la);CHKERRQ(ierr); 344c4762a1bSJed Brown /* Update of Volumetric Auxiliary Data with primary Volumetric Data */ 345c4762a1bSJed Brown ierr = DMGetLocalVector(dm, &la);CHKERRQ(ierr); 346c4762a1bSJed Brown ierr = VecSet(la, 1.0);CHKERRQ(ierr); 347c4762a1bSJed Brown ierr = TestFieldProjection(auxdm, dm, NULL, la, "Volumetric Auxiliary Update with Volumetric Primary", &user);CHKERRQ(ierr); 348c4762a1bSJed Brown ierr = DMRestoreLocalVector(dm, &la);CHKERRQ(ierr); 349c4762a1bSJed Brown ierr = DMDestroy(&auxdm);CHKERRQ(ierr); 350c4762a1bSJed Brown } 351c4762a1bSJed Brown if (user.subdomain) { 352c4762a1bSJed Brown DMLabel domLabel; 353c4762a1bSJed Brown 354c4762a1bSJed Brown /* Subdomain Mesh Projection */ 355c4762a1bSJed Brown ierr = CreateSubdomainMesh(dm, &domLabel, &subdm, &user);CHKERRQ(ierr); 356c4762a1bSJed Brown ierr = TestFunctionProjection(subdm, NULL, NULL, NULL, "Subdomain Primary", &user);CHKERRQ(ierr); 357c4762a1bSJed Brown ierr = TestFieldProjection(subdm, NULL, NULL, NULL, "Subdomain Primary", &user);CHKERRQ(ierr); 358c4762a1bSJed Brown if (user.auxfield) { 359c4762a1bSJed Brown /* Subdomain Mesh Projection with Subdomain Data */ 3609a2a23afSMatthew G. Knepley ierr = CreateAuxiliaryVec(subdm, &auxdm, &la, &user);CHKERRQ(ierr); 361c4762a1bSJed Brown ierr = TestFunctionProjection(subdm, auxdm, NULL, la, "Subdomain Primary and Subdomain Auxiliary", &user);CHKERRQ(ierr); 362c4762a1bSJed Brown ierr = TestFieldProjection(subdm, auxdm, NULL, la, "Subdomain Primary and Subdomain Auxiliary", &user);CHKERRQ(ierr); 363c4762a1bSJed Brown ierr = VecDestroy(&la);CHKERRQ(ierr); 364c4762a1bSJed Brown ierr = DMDestroy(&auxdm);CHKERRQ(ierr); 365c4762a1bSJed Brown /* Subdomain Mesh Projection with Volumetric Data */ 3669a2a23afSMatthew G. Knepley ierr = CreateAuxiliaryVec(dm, &auxdm, &la, &user);CHKERRQ(ierr); 367c4762a1bSJed Brown ierr = TestFunctionProjection(subdm, auxdm, NULL, la, "Subdomain Primary and Volumetric Auxiliary", &user);CHKERRQ(ierr); 368c4762a1bSJed Brown ierr = TestFieldProjection(subdm, auxdm, NULL, la, "Subdomain Primary and Volumetric Auxiliary", &user);CHKERRQ(ierr); 369c4762a1bSJed Brown ierr = VecDestroy(&la);CHKERRQ(ierr); 370c4762a1bSJed Brown ierr = DMDestroy(&auxdm);CHKERRQ(ierr); 371c4762a1bSJed Brown /* Volumetric Mesh Projection with Subdomain Data */ 3729a2a23afSMatthew G. Knepley ierr = CreateAuxiliaryVec(subdm, &auxdm, &la, &user);CHKERRQ(ierr); 373c4762a1bSJed Brown ierr = TestFunctionProjection(subdm, auxdm, domLabel, la, "Volumetric Primary and Subdomain Auxiliary", &user);CHKERRQ(ierr); 374c4762a1bSJed Brown ierr = TestFieldProjection(subdm, auxdm, domLabel, la, "Volumetric Primary and Subdomain Auxiliary", &user);CHKERRQ(ierr); 375c4762a1bSJed Brown ierr = VecDestroy(&la);CHKERRQ(ierr); 376c4762a1bSJed Brown ierr = DMDestroy(&auxdm);CHKERRQ(ierr); 377c4762a1bSJed Brown } 378c4762a1bSJed Brown ierr = DMDestroy(&subdm);CHKERRQ(ierr); 379c4762a1bSJed Brown ierr = DMLabelDestroy(&domLabel);CHKERRQ(ierr); 380c4762a1bSJed Brown } 381c4762a1bSJed Brown if (user.submesh) { 382c4762a1bSJed Brown DMLabel bdLabel; 383c4762a1bSJed Brown 384c4762a1bSJed Brown /* Boundary Mesh Projection */ 385c4762a1bSJed Brown ierr = CreateBoundaryMesh(dm, &bdLabel, &subdm, &user);CHKERRQ(ierr); 386c4762a1bSJed Brown ierr = TestFunctionProjection(subdm, NULL, NULL, NULL, "Boundary Primary", &user);CHKERRQ(ierr); 387c4762a1bSJed Brown ierr = TestFieldProjection(subdm, NULL, NULL, NULL, "Boundary Primary", &user);CHKERRQ(ierr); 388c4762a1bSJed Brown if (user.auxfield) { 389c4762a1bSJed Brown /* Boundary Mesh Projection with Boundary Data */ 3909a2a23afSMatthew G. Knepley ierr = CreateAuxiliaryVec(subdm, &auxdm, &la, &user);CHKERRQ(ierr); 391c4762a1bSJed Brown ierr = TestFunctionProjection(subdm, auxdm, NULL, la, "Boundary Primary and Boundary Auxiliary", &user);CHKERRQ(ierr); 392c4762a1bSJed Brown ierr = TestFieldProjection(subdm, auxdm, NULL, la, "Boundary Primary and Boundary Auxiliary", &user);CHKERRQ(ierr); 393c4762a1bSJed Brown ierr = VecDestroy(&la);CHKERRQ(ierr); 394c4762a1bSJed Brown ierr = DMDestroy(&auxdm);CHKERRQ(ierr); 395c4762a1bSJed Brown /* Volumetric Mesh Projection with Boundary Data */ 3969a2a23afSMatthew G. Knepley ierr = CreateAuxiliaryVec(subdm, &auxdm, &la, &user);CHKERRQ(ierr); 397c4762a1bSJed Brown ierr = TestFunctionProjection(dm, auxdm, bdLabel, la, "Volumetric Primary and Boundary Auxiliary", &user);CHKERRQ(ierr); 398c4762a1bSJed Brown ierr = TestFieldProjection(dm, auxdm, bdLabel, la, "Volumetric Primary and Boundary Auxiliary", &user);CHKERRQ(ierr); 399c4762a1bSJed Brown ierr = VecDestroy(&la);CHKERRQ(ierr); 400c4762a1bSJed Brown ierr = DMDestroy(&auxdm);CHKERRQ(ierr); 401c4762a1bSJed Brown } 402c4762a1bSJed Brown ierr = DMLabelDestroy(&bdLabel);CHKERRQ(ierr); 403c4762a1bSJed Brown ierr = DMDestroy(&subdm);CHKERRQ(ierr); 404c4762a1bSJed Brown } 405c4762a1bSJed Brown ierr = DMDestroy(&dm);CHKERRQ(ierr); 406c4762a1bSJed Brown ierr = PetscFinalize(); 407c4762a1bSJed Brown return ierr; 408c4762a1bSJed Brown } 409c4762a1bSJed Brown 410c4762a1bSJed Brown /*TEST 411c4762a1bSJed Brown 412c4762a1bSJed Brown test: 413c4762a1bSJed Brown suffix: 0 414c4762a1bSJed Brown requires: triangle 415*30602db0SMatthew G. Knepley args: -dm_plex_box_faces 1,1 -func_view -local_func_view -local_input_view -local_field_view 416c4762a1bSJed Brown test: 417c4762a1bSJed Brown suffix: mf_0 418c4762a1bSJed Brown requires: triangle 419*30602db0SMatthew G. Knepley args: -dm_plex_box_faces 1,1 -velocity_petscspace_degree 1 -velocity_petscfe_default_quadrature_order 2 \ 420c4762a1bSJed Brown -pressure_petscspace_degree 2 -pressure_petscfe_default_quadrature_order 2 \ 421c4762a1bSJed Brown -multifield -output_petscspace_degree 1 -output_petscfe_default_quadrature_order 2 \ 422c4762a1bSJed Brown -local_input_view -local_field_view 423c4762a1bSJed Brown test: 424c4762a1bSJed Brown suffix: 1 425c4762a1bSJed Brown requires: triangle 426*30602db0SMatthew 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 427c4762a1bSJed Brown test: 428c4762a1bSJed Brown suffix: 2 429c4762a1bSJed Brown requires: triangle 430*30602db0SMatthew 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 431c4762a1bSJed Brown 432c4762a1bSJed Brown TEST*/ 433c4762a1bSJed Brown 434c4762a1bSJed Brown /* 435c4762a1bSJed Brown Post-processing wants to project a function of the fields into some FE space 436c4762a1bSJed Brown - This is DMProjectField() 437c4762a1bSJed Brown - What about changing the number of components of the output, like displacement to stress? Aux vars 438c4762a1bSJed Brown 439c4762a1bSJed Brown Update of state variables 440c4762a1bSJed Brown - This is DMProjectField(), but solution must be the aux var 441c4762a1bSJed Brown */ 442