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); 71*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-multifield", "Flag for trying different numbers of input/output fields", "ex23.c", options->multifield, &options->multifield, NULL)); 72*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-subdomain", "Flag for trying volumetric submesh", "ex23.c", options->subdomain, &options->subdomain, NULL)); 73*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-submesh", "Flag for trying boundary submesh", "ex23.c", options->submesh, &options->submesh, NULL)); 74*5f80ce2aSJacob 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; 82*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreate(comm, dm)); 83*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetType(*dm, DMPLEX)); 84*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(*dm)); 85*5f80ce2aSJacob 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; 95*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetComm((PetscObject) dm, &comm)); 96*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFECreateDefault(comm, dim, dim, simplex, "velocity_", -1, &fe)); 97*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) fe, "velocity")); 98*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetField(dm, 0, NULL, (PetscObject) fe)); 99*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFEDestroy(&fe)); 100*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFECreateDefault(comm, dim, 1, simplex, "pressure_", -1, &fe)); 101*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) fe, "pressure")); 102*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetField(dm, 1, NULL, (PetscObject) fe)); 103*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFEDestroy(&fe)); 104*5f80ce2aSJacob 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; 114*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetComm((PetscObject) dm, &comm)); 115*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFECreateDefault(comm, dim, dim, simplex, "output_", -1, &fe)); 116*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) fe, "output")); 117*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetField(dm, 0, NULL, (PetscObject) fe)); 118*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFEDestroy(&fe)); 119*5f80ce2aSJacob 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; 130*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexIsSimplex(dm, &simplex)); 131*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 132*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelCreate(PETSC_COMM_SELF, "subdomain", &label)); 133*5f80ce2aSJacob Faibussowitsch for (c = cStart + (cEnd-cStart)/2; c < cEnd; ++c) CHKERRQ(DMLabelSetValue(label, c, 1)); 134*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexFilter(dm, label, 1, subdm)); 135*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(*subdm, &dim)); 136*5f80ce2aSJacob Faibussowitsch CHKERRQ(SetupDiscretization(*subdm, dim, simplex, user)); 137*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) *subdm, "subdomain")); 138*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(*subdm, NULL, "-sub_dm_view")); 139c4762a1bSJed Brown if (domLabel) *domLabel = label; 140*5f80ce2aSJacob 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; 151*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexIsSimplex(dm, &simplex)); 152*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelCreate(PETSC_COMM_SELF, "sub", &label)); 153*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexMarkBoundaryFaces(dm, 1, label)); 154*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexLabelComplete(dm, label)); 155*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateSubmesh(dm, label, 1, PETSC_TRUE, subdm)); 156*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(*subdm, &dim)); 157*5f80ce2aSJacob Faibussowitsch CHKERRQ(SetupDiscretization(*subdm, dim, simplex, user)); 158*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) *subdm, "boundary")); 159*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(*subdm, NULL, "-sub_dm_view")); 160c4762a1bSJed Brown if (bdLabel) *bdLabel = label; 161*5f80ce2aSJacob 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; 172*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(dm, &dim)); 173*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexIsSimplex(dm, &simplex)); 174*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetNumFields(dm, &Nf)); 175*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(Nf, &afuncs)); 176c4762a1bSJed Brown for (f = 0; f < Nf; ++f) afuncs[f] = linear; 177*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMClone(dm, auxdm)); 178*5f80ce2aSJacob Faibussowitsch CHKERRQ(SetupDiscretization(*auxdm, dim, simplex, user)); 179*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateLocalVector(*auxdm, la)); 180*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMProjectFunctionLocal(dm, 0.0, afuncs, NULL, INSERT_VALUES, *la)); 181*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecViewFromOptions(*la, NULL, "-local_aux_view")); 182*5f80ce2aSJacob 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; 195*5f80ce2aSJacob Faibussowitsch if (dmAux) CHKERRQ(DMSetAuxiliaryVec(dm, NULL, 0, 0, la)); 196*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetNumFields(dm, &Nf)); 197*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(Nf, &funcs)); 198c4762a1bSJed Brown for (f = 0; f < Nf; ++f) funcs[f] = linear; 199*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetGlobalVector(dm, &x)); 200*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrcpy(lname, "Function ")); 201*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrcat(lname, name)); 202*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) x, lname)); 203*5f80ce2aSJacob Faibussowitsch if (!label) CHKERRQ(DMProjectFunction(dm, 0.0, funcs, NULL, INSERT_VALUES, x)); 204*5f80ce2aSJacob Faibussowitsch else CHKERRQ(DMProjectFunctionLabel(dm, 0.0, label, 1, val, 0, NULL, funcs, NULL, INSERT_VALUES, x)); 205*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecViewFromOptions(x, NULL, "-func_view")); 206*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreGlobalVector(dm, &x)); 207*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(dm, &lx)); 208*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrcpy(lname, "Local Function ")); 209*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrcat(lname, name)); 210*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) lx, lname)); 211*5f80ce2aSJacob Faibussowitsch if (!label) CHKERRQ(DMProjectFunctionLocal(dm, 0.0, funcs, NULL, INSERT_VALUES, lx)); 212*5f80ce2aSJacob Faibussowitsch else CHKERRQ(DMProjectFunctionLabelLocal(dm, 0.0, label, 1, val, 0, NULL, funcs, NULL, INSERT_VALUES, lx)); 213*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecViewFromOptions(lx, NULL, "-local_func_view")); 214*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(dm, &lx)); 215*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(funcs)); 216*5f80ce2aSJacob 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; 233*5f80ce2aSJacob Faibussowitsch if (dmAux) CHKERRQ(DMSetAuxiliaryVec(dm, NULL, 0, 0, la)); 234*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetNumFields(dm, &Nf)); 235*5f80ce2aSJacob 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; 239*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(dm, &lu)); 240*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrcpy(lname, "Local Field Input ")); 241*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrcat(lname, name)); 242*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) lu, lname)); 243*5f80ce2aSJacob Faibussowitsch if (!label) CHKERRQ(DMProjectFunctionLocal(dm, 0.0, afuncs, NULL, INSERT_VALUES, lu)); 244*5f80ce2aSJacob Faibussowitsch else CHKERRQ(DMProjectFunctionLabelLocal(dm, 0.0, label, 1, val, 0, NULL, afuncs, NULL, INSERT_VALUES, lu)); 245*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecViewFromOptions(lu, NULL, "-local_input_view")); 246*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(dm, &lx)); 247*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrcpy(lname, "Local Field ")); 248*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrcat(lname, name)); 249*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) lx, lname)); 250*5f80ce2aSJacob Faibussowitsch if (!label) CHKERRQ(DMProjectFieldLocal(dm, 0.0, lu, funcs, INSERT_VALUES, lx)); 251*5f80ce2aSJacob Faibussowitsch else CHKERRQ(DMProjectFieldLabelLocal(dm, 0.0, label, 1, val, 0, NULL, lu, funcs, INSERT_VALUES, lx)); 252*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecViewFromOptions(lx, NULL, "-local_field_view")); 253*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(dm, &lx)); 254*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(dm, &lu)); 255*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree2(funcs, afuncs)); 256*5f80ce2aSJacob 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; 273*5f80ce2aSJacob Faibussowitsch if (dmAux) CHKERRQ(DMSetAuxiliaryVec(dm, NULL, 0, 0, la)); 274*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetNumFields(dm, &Nf)); 275*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetNumFields(dmIn, &NfIn)); 276*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc2(Nf, &funcs, NfIn, &afuncs)); 277c4762a1bSJed Brown funcs[0] = divergence_sq; 278c4762a1bSJed Brown afuncs[0] = linear2; 279c4762a1bSJed Brown afuncs[1] = linear; 280*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(dmIn, &lu)); 281*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrcpy(lname, "Local MultiField Input ")); 282*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrcat(lname, name)); 283*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) lu, lname)); 284*5f80ce2aSJacob Faibussowitsch if (!label) CHKERRQ(DMProjectFunctionLocal(dmIn, 0.0, afuncs, NULL, INSERT_VALUES, lu)); 285*5f80ce2aSJacob Faibussowitsch else CHKERRQ(DMProjectFunctionLabelLocal(dmIn, 0.0, label, 1, val, 0, NULL, afuncs, NULL, INSERT_VALUES, lu)); 286*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecViewFromOptions(lu, NULL, "-local_input_view")); 287*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(dm, &lx)); 288*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrcpy(lname, "Local MultiField ")); 289*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrcat(lname, name)); 290*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) lx, lname)); 291*5f80ce2aSJacob Faibussowitsch if (!label) CHKERRQ(DMProjectFieldLocal(dm, 0.0, lu, funcs, INSERT_VALUES, lx)); 292*5f80ce2aSJacob Faibussowitsch else CHKERRQ(DMProjectFieldLabelLocal(dm, 0.0, label, 1, val, 0, NULL, lu, funcs, INSERT_VALUES, lx)); 293*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecViewFromOptions(lx, NULL, "-local_field_view")); 294*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(dm, &lx)); 295*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(dmIn, &lu)); 296*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree2(funcs, afuncs)); 297*5f80ce2aSJacob 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 PetscErrorCode ierr; 309c4762a1bSJed Brown 310c4762a1bSJed Brown ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr; 311*5f80ce2aSJacob Faibussowitsch CHKERRQ(ProcessOptions(&user)); 312*5f80ce2aSJacob Faibussowitsch CHKERRQ(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 313*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(dm, &dim)); 314*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexIsSimplex(dm, &simplex)); 315*5f80ce2aSJacob Faibussowitsch CHKERRQ(SetupDiscretization(dm, dim, simplex, &user)); 316c4762a1bSJed Brown /* Volumetric Mesh Projection */ 317c4762a1bSJed Brown if (!user.multifield) { 318*5f80ce2aSJacob Faibussowitsch CHKERRQ(TestFunctionProjection(dm, NULL, NULL, NULL, "Volumetric Primary", &user)); 319*5f80ce2aSJacob Faibussowitsch CHKERRQ(TestFieldProjection(dm, NULL, NULL, NULL, "Volumetric Primary", &user)); 320c4762a1bSJed Brown } else { 321c4762a1bSJed Brown DM dmOut; 322c4762a1bSJed Brown 323*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMClone(dm, &dmOut)); 324*5f80ce2aSJacob Faibussowitsch CHKERRQ(SetupOutputDiscretization(dmOut, dim, simplex, &user)); 325*5f80ce2aSJacob Faibussowitsch CHKERRQ(TestFieldProjectionMultiple(dmOut, dm, NULL, NULL, NULL, "Volumetric Primary", &user)); 326*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dmOut)); 327c4762a1bSJed Brown } 328c4762a1bSJed Brown if (user.auxfield) { 329c4762a1bSJed Brown /* Volumetric Mesh Projection with Volumetric Data */ 330*5f80ce2aSJacob Faibussowitsch CHKERRQ(CreateAuxiliaryVec(dm, &auxdm, &la, &user)); 331*5f80ce2aSJacob Faibussowitsch CHKERRQ(TestFunctionProjection(dm, auxdm, NULL, la, "Volumetric Primary and Volumetric Auxiliary", &user)); 332*5f80ce2aSJacob Faibussowitsch CHKERRQ(TestFieldProjection(dm, auxdm, NULL, la, "Volumetric Primary and Volumetric Auxiliary", &user)); 333*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&la)); 334c4762a1bSJed Brown /* Update of Volumetric Auxiliary Data with primary Volumetric Data */ 335*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(dm, &la)); 336*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(la, 1.0)); 337*5f80ce2aSJacob Faibussowitsch CHKERRQ(TestFieldProjection(auxdm, dm, NULL, la, "Volumetric Auxiliary Update with Volumetric Primary", &user)); 338*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(dm, &la)); 339*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&auxdm)); 340c4762a1bSJed Brown } 341c4762a1bSJed Brown if (user.subdomain) { 342c4762a1bSJed Brown DMLabel domLabel; 343c4762a1bSJed Brown 344c4762a1bSJed Brown /* Subdomain Mesh Projection */ 345*5f80ce2aSJacob Faibussowitsch CHKERRQ(CreateSubdomainMesh(dm, &domLabel, &subdm, &user)); 346*5f80ce2aSJacob Faibussowitsch CHKERRQ(TestFunctionProjection(subdm, NULL, NULL, NULL, "Subdomain Primary", &user)); 347*5f80ce2aSJacob Faibussowitsch CHKERRQ(TestFieldProjection(subdm, NULL, NULL, NULL, "Subdomain Primary", &user)); 348c4762a1bSJed Brown if (user.auxfield) { 349c4762a1bSJed Brown /* Subdomain Mesh Projection with Subdomain Data */ 350*5f80ce2aSJacob Faibussowitsch CHKERRQ(CreateAuxiliaryVec(subdm, &auxdm, &la, &user)); 351*5f80ce2aSJacob Faibussowitsch CHKERRQ(TestFunctionProjection(subdm, auxdm, NULL, la, "Subdomain Primary and Subdomain Auxiliary", &user)); 352*5f80ce2aSJacob Faibussowitsch CHKERRQ(TestFieldProjection(subdm, auxdm, NULL, la, "Subdomain Primary and Subdomain Auxiliary", &user)); 353*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&la)); 354*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&auxdm)); 355c4762a1bSJed Brown /* Subdomain Mesh Projection with Volumetric Data */ 356*5f80ce2aSJacob Faibussowitsch CHKERRQ(CreateAuxiliaryVec(dm, &auxdm, &la, &user)); 357*5f80ce2aSJacob Faibussowitsch CHKERRQ(TestFunctionProjection(subdm, auxdm, NULL, la, "Subdomain Primary and Volumetric Auxiliary", &user)); 358*5f80ce2aSJacob Faibussowitsch CHKERRQ(TestFieldProjection(subdm, auxdm, NULL, la, "Subdomain Primary and Volumetric Auxiliary", &user)); 359*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&la)); 360*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&auxdm)); 361c4762a1bSJed Brown /* Volumetric Mesh Projection with Subdomain Data */ 362*5f80ce2aSJacob Faibussowitsch CHKERRQ(CreateAuxiliaryVec(subdm, &auxdm, &la, &user)); 363*5f80ce2aSJacob Faibussowitsch CHKERRQ(TestFunctionProjection(subdm, auxdm, domLabel, la, "Volumetric Primary and Subdomain Auxiliary", &user)); 364*5f80ce2aSJacob Faibussowitsch CHKERRQ(TestFieldProjection(subdm, auxdm, domLabel, la, "Volumetric Primary and Subdomain Auxiliary", &user)); 365*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&la)); 366*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&auxdm)); 367c4762a1bSJed Brown } 368*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&subdm)); 369*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelDestroy(&domLabel)); 370c4762a1bSJed Brown } 371c4762a1bSJed Brown if (user.submesh) { 372c4762a1bSJed Brown DMLabel bdLabel; 373c4762a1bSJed Brown 374c4762a1bSJed Brown /* Boundary Mesh Projection */ 375*5f80ce2aSJacob Faibussowitsch CHKERRQ(CreateBoundaryMesh(dm, &bdLabel, &subdm, &user)); 376*5f80ce2aSJacob Faibussowitsch CHKERRQ(TestFunctionProjection(subdm, NULL, NULL, NULL, "Boundary Primary", &user)); 377*5f80ce2aSJacob Faibussowitsch CHKERRQ(TestFieldProjection(subdm, NULL, NULL, NULL, "Boundary Primary", &user)); 378c4762a1bSJed Brown if (user.auxfield) { 379c4762a1bSJed Brown /* Boundary Mesh Projection with Boundary Data */ 380*5f80ce2aSJacob Faibussowitsch CHKERRQ(CreateAuxiliaryVec(subdm, &auxdm, &la, &user)); 381*5f80ce2aSJacob Faibussowitsch CHKERRQ(TestFunctionProjection(subdm, auxdm, NULL, la, "Boundary Primary and Boundary Auxiliary", &user)); 382*5f80ce2aSJacob Faibussowitsch CHKERRQ(TestFieldProjection(subdm, auxdm, NULL, la, "Boundary Primary and Boundary Auxiliary", &user)); 383*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&la)); 384*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&auxdm)); 385c4762a1bSJed Brown /* Volumetric Mesh Projection with Boundary Data */ 386*5f80ce2aSJacob Faibussowitsch CHKERRQ(CreateAuxiliaryVec(subdm, &auxdm, &la, &user)); 387*5f80ce2aSJacob Faibussowitsch CHKERRQ(TestFunctionProjection(dm, auxdm, bdLabel, la, "Volumetric Primary and Boundary Auxiliary", &user)); 388*5f80ce2aSJacob Faibussowitsch CHKERRQ(TestFieldProjection(dm, auxdm, bdLabel, la, "Volumetric Primary and Boundary Auxiliary", &user)); 389*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&la)); 390*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&auxdm)); 391c4762a1bSJed Brown } 392*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelDestroy(&bdLabel)); 393*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&subdm)); 394c4762a1bSJed Brown } 395*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dm)); 396c4762a1bSJed Brown ierr = PetscFinalize(); 397c4762a1bSJed Brown return ierr; 398c4762a1bSJed Brown } 399c4762a1bSJed Brown 400c4762a1bSJed Brown /*TEST 401c4762a1bSJed Brown 402c4762a1bSJed Brown test: 403c4762a1bSJed Brown suffix: 0 404c4762a1bSJed Brown requires: triangle 40530602db0SMatthew G. Knepley args: -dm_plex_box_faces 1,1 -func_view -local_func_view -local_input_view -local_field_view 406c4762a1bSJed Brown test: 407c4762a1bSJed Brown suffix: mf_0 408c4762a1bSJed Brown requires: triangle 40930602db0SMatthew G. Knepley args: -dm_plex_box_faces 1,1 -velocity_petscspace_degree 1 -velocity_petscfe_default_quadrature_order 2 \ 410c4762a1bSJed Brown -pressure_petscspace_degree 2 -pressure_petscfe_default_quadrature_order 2 \ 411c4762a1bSJed Brown -multifield -output_petscspace_degree 1 -output_petscfe_default_quadrature_order 2 \ 412c4762a1bSJed Brown -local_input_view -local_field_view 413c4762a1bSJed Brown test: 414c4762a1bSJed Brown suffix: 1 415c4762a1bSJed Brown requires: triangle 41630602db0SMatthew 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 417c4762a1bSJed Brown test: 418c4762a1bSJed Brown suffix: 2 419c4762a1bSJed Brown requires: triangle 42030602db0SMatthew 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 421c4762a1bSJed Brown 422c4762a1bSJed Brown TEST*/ 423c4762a1bSJed Brown 424c4762a1bSJed Brown /* 425c4762a1bSJed Brown Post-processing wants to project a function of the fields into some FE space 426c4762a1bSJed Brown - This is DMProjectField() 427c4762a1bSJed Brown - What about changing the number of components of the output, like displacement to stress? Aux vars 428c4762a1bSJed Brown 429c4762a1bSJed Brown Update of state variables 430c4762a1bSJed Brown - This is DMProjectField(), but solution must be the aux var 431c4762a1bSJed Brown */ 432