1*c4762a1bSJed Brown static char help[] = "Test for function and field projection\n\n"; 2*c4762a1bSJed Brown 3*c4762a1bSJed Brown #include <petscdmplex.h> 4*c4762a1bSJed Brown #include <petscds.h> 5*c4762a1bSJed Brown 6*c4762a1bSJed Brown typedef struct { 7*c4762a1bSJed Brown PetscInt dim; /* The topological mesh dimension */ 8*c4762a1bSJed Brown PetscBool cellSimplex; /* Flag for simplices */ 9*c4762a1bSJed Brown PetscBool multifield; /* Different numbers of input and output fields */ 10*c4762a1bSJed Brown PetscBool subdomain; /* Try with a volumetric submesh */ 11*c4762a1bSJed Brown PetscBool submesh; /* Try with a boundary submesh */ 12*c4762a1bSJed Brown PetscBool auxfield; /* Try with auxiliary fields */ 13*c4762a1bSJed Brown } AppCtx; 14*c4762a1bSJed Brown 15*c4762a1bSJed Brown /* (x + y)*dim + d */ 16*c4762a1bSJed Brown static PetscErrorCode linear(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 17*c4762a1bSJed Brown { 18*c4762a1bSJed Brown PetscInt c; 19*c4762a1bSJed Brown for (c = 0; c < Nc; ++c) u[c] = (x[0] + x[1])*Nc + c; 20*c4762a1bSJed Brown return 0; 21*c4762a1bSJed Brown } 22*c4762a1bSJed Brown 23*c4762a1bSJed Brown /* {x, y, z} */ 24*c4762a1bSJed Brown static PetscErrorCode linear2(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 25*c4762a1bSJed Brown { 26*c4762a1bSJed Brown PetscInt c; 27*c4762a1bSJed Brown for (c = 0; c < Nc; ++c) u[c] = x[c]; 28*c4762a1bSJed Brown return 0; 29*c4762a1bSJed Brown } 30*c4762a1bSJed Brown 31*c4762a1bSJed Brown /* {u_x, u_y, u_z} */ 32*c4762a1bSJed Brown static void linear_vector(PetscInt dim, PetscInt Nf, PetscInt NfAux, 33*c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 34*c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 35*c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f[]) 36*c4762a1bSJed Brown { 37*c4762a1bSJed Brown PetscInt d; 38*c4762a1bSJed Brown for (d = 0; d < uOff[1]-uOff[0]; ++d) f[d] = u[d+uOff[0]]; 39*c4762a1bSJed Brown } 40*c4762a1bSJed Brown 41*c4762a1bSJed Brown /* p */ 42*c4762a1bSJed Brown static void linear_scalar(PetscInt dim, PetscInt Nf, PetscInt NfAux, 43*c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 44*c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 45*c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f[]) 46*c4762a1bSJed Brown { 47*c4762a1bSJed Brown f[0] = u[uOff[1]]; 48*c4762a1bSJed Brown } 49*c4762a1bSJed Brown 50*c4762a1bSJed Brown /* {div u, p^2} */ 51*c4762a1bSJed Brown static void divergence_sq(PetscInt dim, PetscInt Nf, PetscInt NfAux, 52*c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 53*c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 54*c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f[]) 55*c4762a1bSJed Brown { 56*c4762a1bSJed Brown PetscInt d; 57*c4762a1bSJed Brown f[0] = 0.0; 58*c4762a1bSJed Brown for (d = 0; d < dim; ++d) f[0] += u_x[uOff_x[0]+d*dim+d]; 59*c4762a1bSJed Brown f[1] = PetscSqr(u[uOff[1]]); 60*c4762a1bSJed Brown } 61*c4762a1bSJed Brown 62*c4762a1bSJed Brown static PetscErrorCode ProcessOptions(AppCtx *options) 63*c4762a1bSJed Brown { 64*c4762a1bSJed Brown PetscErrorCode ierr; 65*c4762a1bSJed Brown 66*c4762a1bSJed Brown PetscFunctionBegin; 67*c4762a1bSJed Brown options->dim = 2; 68*c4762a1bSJed Brown options->cellSimplex = PETSC_TRUE; 69*c4762a1bSJed Brown options->multifield = PETSC_FALSE; 70*c4762a1bSJed Brown options->subdomain = PETSC_FALSE; 71*c4762a1bSJed Brown options->submesh = PETSC_FALSE; 72*c4762a1bSJed Brown options->auxfield = PETSC_FALSE; 73*c4762a1bSJed Brown 74*c4762a1bSJed Brown ierr = PetscOptionsBegin(PETSC_COMM_SELF, "", "Meshing Problem Options", "DMPLEX");CHKERRQ(ierr); 75*c4762a1bSJed Brown ierr = PetscOptionsRangeInt("-dim", "The topological mesh dimension", "ex23.c", options->dim, &options->dim, NULL,1,3);CHKERRQ(ierr); 76*c4762a1bSJed Brown ierr = PetscOptionsBool("-cellSimplex", "Flag for simplices", "ex23.c", options->cellSimplex, &options->cellSimplex, NULL);CHKERRQ(ierr); 77*c4762a1bSJed Brown ierr = PetscOptionsBool("-multifield", "Flag for trying different numbers of input/output fields", "ex23.c", options->multifield, &options->multifield, NULL);CHKERRQ(ierr); 78*c4762a1bSJed Brown ierr = PetscOptionsBool("-subdomain", "Flag for trying volumetric submesh", "ex23.c", options->subdomain, &options->subdomain, NULL);CHKERRQ(ierr); 79*c4762a1bSJed Brown ierr = PetscOptionsBool("-submesh", "Flag for trying boundary submesh", "ex23.c", options->submesh, &options->submesh, NULL);CHKERRQ(ierr); 80*c4762a1bSJed Brown ierr = PetscOptionsBool("-auxfield", "Flag for trying auxiliary fields", "ex23.c", options->auxfield, &options->auxfield, NULL);CHKERRQ(ierr); 81*c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 82*c4762a1bSJed Brown PetscFunctionReturn(0); 83*c4762a1bSJed Brown } 84*c4762a1bSJed Brown 85*c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 86*c4762a1bSJed Brown { 87*c4762a1bSJed Brown PetscInt dim = user->dim; 88*c4762a1bSJed Brown PetscBool simplex = user->cellSimplex; 89*c4762a1bSJed Brown PetscInt cells[3] = {1, 1, 1}; 90*c4762a1bSJed Brown PetscErrorCode ierr; 91*c4762a1bSJed Brown 92*c4762a1bSJed Brown PetscFunctionBegin; 93*c4762a1bSJed Brown ierr = DMPlexCreateBoxMesh(comm, dim, simplex, cells, NULL, NULL, NULL, PETSC_TRUE, dm);CHKERRQ(ierr); 94*c4762a1bSJed Brown { 95*c4762a1bSJed Brown DM ddm = NULL; 96*c4762a1bSJed Brown 97*c4762a1bSJed Brown ierr = DMPlexDistribute(*dm, 0, NULL, &ddm);CHKERRQ(ierr); 98*c4762a1bSJed Brown if (ddm) { 99*c4762a1bSJed Brown ierr = DMDestroy(dm);CHKERRQ(ierr); 100*c4762a1bSJed Brown *dm = ddm; 101*c4762a1bSJed Brown } 102*c4762a1bSJed Brown } 103*c4762a1bSJed Brown ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 104*c4762a1bSJed Brown ierr = DMViewFromOptions(*dm, NULL, "-orig_dm_view");CHKERRQ(ierr); 105*c4762a1bSJed Brown PetscFunctionReturn(0); 106*c4762a1bSJed Brown } 107*c4762a1bSJed Brown 108*c4762a1bSJed Brown static PetscErrorCode SetupDiscretization(DM dm, PetscInt dim, PetscBool simplex, AppCtx *user) 109*c4762a1bSJed Brown { 110*c4762a1bSJed Brown PetscFE fe; 111*c4762a1bSJed Brown MPI_Comm comm; 112*c4762a1bSJed Brown PetscErrorCode ierr; 113*c4762a1bSJed Brown 114*c4762a1bSJed Brown PetscFunctionBeginUser; 115*c4762a1bSJed Brown ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 116*c4762a1bSJed Brown ierr = PetscFECreateDefault(comm, dim, dim, simplex, "velocity_", -1, &fe);CHKERRQ(ierr); 117*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) fe, "velocity");CHKERRQ(ierr); 118*c4762a1bSJed Brown ierr = DMSetField(dm, 0, NULL, (PetscObject) fe);CHKERRQ(ierr); 119*c4762a1bSJed Brown ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 120*c4762a1bSJed Brown ierr = PetscFECreateDefault(comm, dim, 1, simplex, "pressure_", -1, &fe);CHKERRQ(ierr); 121*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) fe, "pressure");CHKERRQ(ierr); 122*c4762a1bSJed Brown ierr = DMSetField(dm, 1, NULL, (PetscObject) fe);CHKERRQ(ierr); 123*c4762a1bSJed Brown ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 124*c4762a1bSJed Brown ierr = DMCreateDS(dm);CHKERRQ(ierr); 125*c4762a1bSJed Brown PetscFunctionReturn(0); 126*c4762a1bSJed Brown } 127*c4762a1bSJed Brown 128*c4762a1bSJed Brown static PetscErrorCode SetupOutputDiscretization(DM dm, PetscInt dim, PetscBool simplex, AppCtx *user) 129*c4762a1bSJed Brown { 130*c4762a1bSJed Brown PetscFE fe; 131*c4762a1bSJed Brown MPI_Comm comm; 132*c4762a1bSJed Brown PetscErrorCode ierr; 133*c4762a1bSJed Brown 134*c4762a1bSJed Brown PetscFunctionBeginUser; 135*c4762a1bSJed Brown ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 136*c4762a1bSJed Brown ierr = PetscFECreateDefault(comm, dim, dim, simplex, "output_", -1, &fe);CHKERRQ(ierr); 137*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) fe, "output");CHKERRQ(ierr); 138*c4762a1bSJed Brown ierr = DMSetField(dm, 0, NULL, (PetscObject) fe);CHKERRQ(ierr); 139*c4762a1bSJed Brown ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 140*c4762a1bSJed Brown ierr = DMCreateDS(dm);CHKERRQ(ierr); 141*c4762a1bSJed Brown PetscFunctionReturn(0); 142*c4762a1bSJed Brown } 143*c4762a1bSJed Brown 144*c4762a1bSJed Brown static PetscErrorCode CreateSubdomainMesh(DM dm, DMLabel *domLabel, DM *subdm, AppCtx *user) 145*c4762a1bSJed Brown { 146*c4762a1bSJed Brown DMLabel label; 147*c4762a1bSJed Brown PetscInt dim, cStart, cEnd, c; 148*c4762a1bSJed Brown PetscErrorCode ierr; 149*c4762a1bSJed Brown 150*c4762a1bSJed Brown PetscFunctionBeginUser; 151*c4762a1bSJed Brown ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 152*c4762a1bSJed Brown ierr = DMLabelCreate(PETSC_COMM_SELF, "subdomain", &label);CHKERRQ(ierr); 153*c4762a1bSJed Brown for (c = cStart + (cEnd-cStart)/2; c < cEnd; ++c) {ierr = DMLabelSetValue(label, c, 1);CHKERRQ(ierr);} 154*c4762a1bSJed Brown ierr = DMPlexFilter(dm, label, 1, subdm);CHKERRQ(ierr); 155*c4762a1bSJed Brown ierr = DMGetDimension(*subdm, &dim);CHKERRQ(ierr); 156*c4762a1bSJed Brown ierr = SetupDiscretization(*subdm, dim, user->cellSimplex, user);CHKERRQ(ierr); 157*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) *subdm, "subdomain");CHKERRQ(ierr); 158*c4762a1bSJed Brown ierr = DMViewFromOptions(*subdm, NULL, "-sub_dm_view");CHKERRQ(ierr); 159*c4762a1bSJed Brown if (domLabel) *domLabel = label; 160*c4762a1bSJed Brown else {ierr = DMLabelDestroy(&label);CHKERRQ(ierr);} 161*c4762a1bSJed Brown PetscFunctionReturn(0); 162*c4762a1bSJed Brown } 163*c4762a1bSJed Brown 164*c4762a1bSJed Brown static PetscErrorCode CreateBoundaryMesh(DM dm, DMLabel *bdLabel, DM *subdm, AppCtx *user) 165*c4762a1bSJed Brown { 166*c4762a1bSJed Brown DMLabel label; 167*c4762a1bSJed Brown PetscInt dim; 168*c4762a1bSJed Brown PetscErrorCode ierr; 169*c4762a1bSJed Brown 170*c4762a1bSJed Brown PetscFunctionBeginUser; 171*c4762a1bSJed Brown ierr = DMLabelCreate(PETSC_COMM_SELF, "sub", &label);CHKERRQ(ierr); 172*c4762a1bSJed Brown ierr = DMPlexMarkBoundaryFaces(dm, 1, label);CHKERRQ(ierr); 173*c4762a1bSJed Brown ierr = DMPlexLabelComplete(dm, label);CHKERRQ(ierr); 174*c4762a1bSJed Brown ierr = DMPlexCreateSubmesh(dm, label, 1, PETSC_TRUE, subdm);CHKERRQ(ierr); 175*c4762a1bSJed Brown ierr = DMGetDimension(*subdm, &dim);CHKERRQ(ierr); 176*c4762a1bSJed Brown ierr = SetupDiscretization(*subdm, dim, user->cellSimplex, user);CHKERRQ(ierr); 177*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) *subdm, "boundary");CHKERRQ(ierr); 178*c4762a1bSJed Brown ierr = DMViewFromOptions(*subdm, NULL, "-sub_dm_view");CHKERRQ(ierr); 179*c4762a1bSJed Brown if (bdLabel) *bdLabel = label; 180*c4762a1bSJed Brown else {ierr = DMLabelDestroy(&label);CHKERRQ(ierr);} 181*c4762a1bSJed Brown PetscFunctionReturn(0); 182*c4762a1bSJed Brown } 183*c4762a1bSJed Brown 184*c4762a1bSJed Brown static PetscErrorCode CreateAuxiliaryData(DM dm, DM *auxdm, Vec *la, AppCtx *user) 185*c4762a1bSJed Brown { 186*c4762a1bSJed Brown PetscErrorCode (**afuncs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *); 187*c4762a1bSJed Brown PetscInt dim, Nf, f; 188*c4762a1bSJed Brown PetscErrorCode ierr; 189*c4762a1bSJed Brown 190*c4762a1bSJed Brown PetscFunctionBeginUser; 191*c4762a1bSJed Brown ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 192*c4762a1bSJed Brown ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr); 193*c4762a1bSJed Brown ierr = PetscMalloc1(Nf, &afuncs);CHKERRQ(ierr); 194*c4762a1bSJed Brown for (f = 0; f < Nf; ++f) afuncs[f] = linear; 195*c4762a1bSJed Brown ierr = DMClone(dm, auxdm);CHKERRQ(ierr); 196*c4762a1bSJed Brown ierr = SetupDiscretization(*auxdm, dim, user->cellSimplex, user);CHKERRQ(ierr); 197*c4762a1bSJed Brown ierr = DMCreateLocalVector(*auxdm, la);CHKERRQ(ierr); 198*c4762a1bSJed Brown ierr = DMProjectFunctionLocal(dm, 0.0, afuncs, NULL, INSERT_VALUES, *la);CHKERRQ(ierr); 199*c4762a1bSJed Brown ierr = VecViewFromOptions(*la, NULL, "-local_aux_view");CHKERRQ(ierr); 200*c4762a1bSJed Brown ierr = PetscFree(afuncs);CHKERRQ(ierr); 201*c4762a1bSJed Brown PetscFunctionReturn(0); 202*c4762a1bSJed Brown } 203*c4762a1bSJed Brown 204*c4762a1bSJed Brown static PetscErrorCode TestFunctionProjection(DM dm, DM dmAux, DMLabel label, Vec la, const char name[], AppCtx *user) 205*c4762a1bSJed Brown { 206*c4762a1bSJed Brown PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *); 207*c4762a1bSJed Brown Vec x, lx; 208*c4762a1bSJed Brown PetscInt Nf, f; 209*c4762a1bSJed Brown PetscInt val[1] = {1}; 210*c4762a1bSJed Brown char lname[PETSC_MAX_PATH_LEN]; 211*c4762a1bSJed Brown PetscErrorCode ierr; 212*c4762a1bSJed Brown 213*c4762a1bSJed Brown PetscFunctionBeginUser; 214*c4762a1bSJed Brown if (dmAux) { 215*c4762a1bSJed Brown ierr = PetscObjectCompose((PetscObject) dm, "dmAux", (PetscObject) dmAux);CHKERRQ(ierr); 216*c4762a1bSJed Brown ierr = PetscObjectCompose((PetscObject) dm, "A", (PetscObject) la);CHKERRQ(ierr); 217*c4762a1bSJed Brown } 218*c4762a1bSJed Brown ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr); 219*c4762a1bSJed Brown ierr = PetscMalloc1(Nf, &funcs);CHKERRQ(ierr); 220*c4762a1bSJed Brown for (f = 0; f < Nf; ++f) funcs[f] = linear; 221*c4762a1bSJed Brown ierr = DMGetGlobalVector(dm, &x);CHKERRQ(ierr); 222*c4762a1bSJed Brown ierr = PetscStrcpy(lname, "Function ");CHKERRQ(ierr); 223*c4762a1bSJed Brown ierr = PetscStrcat(lname, name);CHKERRQ(ierr); 224*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) x, lname);CHKERRQ(ierr); 225*c4762a1bSJed Brown if (!label) {ierr = DMProjectFunction(dm, 0.0, funcs, NULL, INSERT_VALUES, x);CHKERRQ(ierr);} 226*c4762a1bSJed Brown else {ierr = DMProjectFunctionLabel(dm, 0.0, label, 1, val, 0, NULL, funcs, NULL, INSERT_VALUES, x);CHKERRQ(ierr);} 227*c4762a1bSJed Brown ierr = VecViewFromOptions(x, NULL, "-func_view");CHKERRQ(ierr); 228*c4762a1bSJed Brown ierr = DMRestoreGlobalVector(dm, &x);CHKERRQ(ierr); 229*c4762a1bSJed Brown ierr = DMGetLocalVector(dm, &lx);CHKERRQ(ierr); 230*c4762a1bSJed Brown ierr = PetscStrcpy(lname, "Local Function ");CHKERRQ(ierr); 231*c4762a1bSJed Brown ierr = PetscStrcat(lname, name);CHKERRQ(ierr); 232*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) lx, lname);CHKERRQ(ierr); 233*c4762a1bSJed Brown if (!label) {ierr = DMProjectFunctionLocal(dm, 0.0, funcs, NULL, INSERT_VALUES, lx);CHKERRQ(ierr);} 234*c4762a1bSJed Brown else {ierr = DMProjectFunctionLabelLocal(dm, 0.0, label, 1, val, 0, NULL, funcs, NULL, INSERT_VALUES, lx);CHKERRQ(ierr);} 235*c4762a1bSJed Brown ierr = VecViewFromOptions(lx, NULL, "-local_func_view");CHKERRQ(ierr); 236*c4762a1bSJed Brown ierr = DMRestoreLocalVector(dm, &lx);CHKERRQ(ierr); 237*c4762a1bSJed Brown ierr = PetscFree(funcs);CHKERRQ(ierr); 238*c4762a1bSJed Brown if (dmAux) { 239*c4762a1bSJed Brown ierr = PetscObjectCompose((PetscObject) dm, "dmAux", NULL);CHKERRQ(ierr); 240*c4762a1bSJed Brown ierr = PetscObjectCompose((PetscObject) dm, "A", NULL);CHKERRQ(ierr); 241*c4762a1bSJed Brown } 242*c4762a1bSJed Brown PetscFunctionReturn(0); 243*c4762a1bSJed Brown } 244*c4762a1bSJed Brown 245*c4762a1bSJed Brown static PetscErrorCode TestFieldProjection(DM dm, DM dmAux, DMLabel label, Vec la, const char name[], AppCtx *user) 246*c4762a1bSJed Brown { 247*c4762a1bSJed Brown PetscErrorCode (**afuncs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *); 248*c4762a1bSJed Brown void (**funcs)(PetscInt, PetscInt, PetscInt, 249*c4762a1bSJed Brown const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 250*c4762a1bSJed Brown const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 251*c4762a1bSJed Brown PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]); 252*c4762a1bSJed Brown Vec lx, lu; 253*c4762a1bSJed Brown PetscInt Nf, f; 254*c4762a1bSJed Brown PetscInt val[1] = {1}; 255*c4762a1bSJed Brown char lname[PETSC_MAX_PATH_LEN]; 256*c4762a1bSJed Brown PetscErrorCode ierr; 257*c4762a1bSJed Brown 258*c4762a1bSJed Brown PetscFunctionBeginUser; 259*c4762a1bSJed Brown if (dmAux) { 260*c4762a1bSJed Brown ierr = PetscObjectCompose((PetscObject) dm, "dmAux", (PetscObject) dmAux);CHKERRQ(ierr); 261*c4762a1bSJed Brown ierr = PetscObjectCompose((PetscObject) dm, "A", (PetscObject) la);CHKERRQ(ierr); 262*c4762a1bSJed Brown } 263*c4762a1bSJed Brown ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr); 264*c4762a1bSJed Brown ierr = PetscMalloc2(Nf, &funcs, Nf, &afuncs);CHKERRQ(ierr); 265*c4762a1bSJed Brown for (f = 0; f < Nf; ++f) afuncs[f] = linear; 266*c4762a1bSJed Brown funcs[0] = linear_vector; 267*c4762a1bSJed Brown funcs[1] = linear_scalar; 268*c4762a1bSJed Brown ierr = DMGetLocalVector(dm, &lu);CHKERRQ(ierr); 269*c4762a1bSJed Brown ierr = PetscStrcpy(lname, "Local Field Input ");CHKERRQ(ierr); 270*c4762a1bSJed Brown ierr = PetscStrcat(lname, name);CHKERRQ(ierr); 271*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) lu, lname);CHKERRQ(ierr); 272*c4762a1bSJed Brown if (!label) {ierr = DMProjectFunctionLocal(dm, 0.0, afuncs, NULL, INSERT_VALUES, lu);CHKERRQ(ierr);} 273*c4762a1bSJed Brown else {ierr = DMProjectFunctionLabelLocal(dm, 0.0, label, 1, val, 0, NULL, afuncs, NULL, INSERT_VALUES, lu);CHKERRQ(ierr);} 274*c4762a1bSJed Brown ierr = VecViewFromOptions(lu, NULL, "-local_input_view");CHKERRQ(ierr); 275*c4762a1bSJed Brown ierr = DMGetLocalVector(dm, &lx);CHKERRQ(ierr); 276*c4762a1bSJed Brown ierr = PetscStrcpy(lname, "Local Field ");CHKERRQ(ierr); 277*c4762a1bSJed Brown ierr = PetscStrcat(lname, name);CHKERRQ(ierr); 278*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) lx, lname);CHKERRQ(ierr); 279*c4762a1bSJed Brown if (!label) {ierr = DMProjectFieldLocal(dm, 0.0, lu, funcs, INSERT_VALUES, lx);CHKERRQ(ierr);} 280*c4762a1bSJed Brown else {ierr = DMProjectFieldLabelLocal(dm, 0.0, label, 1, val, 0, NULL, lu, funcs, INSERT_VALUES, lx);CHKERRQ(ierr);} 281*c4762a1bSJed Brown ierr = VecViewFromOptions(lx, NULL, "-local_field_view");CHKERRQ(ierr); 282*c4762a1bSJed Brown ierr = DMRestoreLocalVector(dm, &lx);CHKERRQ(ierr); 283*c4762a1bSJed Brown ierr = DMRestoreLocalVector(dm, &lu);CHKERRQ(ierr); 284*c4762a1bSJed Brown ierr = PetscFree2(funcs, afuncs);CHKERRQ(ierr); 285*c4762a1bSJed Brown if (dmAux) { 286*c4762a1bSJed Brown ierr = PetscObjectCompose((PetscObject) dm, "dmAux", NULL);CHKERRQ(ierr); 287*c4762a1bSJed Brown ierr = PetscObjectCompose((PetscObject) dm, "A", NULL);CHKERRQ(ierr); 288*c4762a1bSJed Brown } 289*c4762a1bSJed Brown PetscFunctionReturn(0); 290*c4762a1bSJed Brown } 291*c4762a1bSJed Brown 292*c4762a1bSJed Brown static PetscErrorCode TestFieldProjectionMultiple(DM dm, DM dmIn, DM dmAux, DMLabel label, Vec la, const char name[], AppCtx *user) 293*c4762a1bSJed Brown { 294*c4762a1bSJed Brown PetscErrorCode (**afuncs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *); 295*c4762a1bSJed Brown void (**funcs)(PetscInt, PetscInt, PetscInt, 296*c4762a1bSJed Brown const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 297*c4762a1bSJed Brown const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 298*c4762a1bSJed Brown PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]); 299*c4762a1bSJed Brown Vec lx, lu; 300*c4762a1bSJed Brown PetscInt Nf, NfIn; 301*c4762a1bSJed Brown PetscInt val[1] = {1}; 302*c4762a1bSJed Brown char lname[PETSC_MAX_PATH_LEN]; 303*c4762a1bSJed Brown PetscErrorCode ierr; 304*c4762a1bSJed Brown 305*c4762a1bSJed Brown PetscFunctionBeginUser; 306*c4762a1bSJed Brown if (dmAux) { 307*c4762a1bSJed Brown ierr = PetscObjectCompose((PetscObject) dm, "dmAux", (PetscObject) dmAux);CHKERRQ(ierr); 308*c4762a1bSJed Brown ierr = PetscObjectCompose((PetscObject) dm, "A", (PetscObject) la);CHKERRQ(ierr); 309*c4762a1bSJed Brown } 310*c4762a1bSJed Brown ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr); 311*c4762a1bSJed Brown ierr = DMGetNumFields(dmIn, &NfIn);CHKERRQ(ierr); 312*c4762a1bSJed Brown ierr = PetscMalloc2(Nf, &funcs, NfIn, &afuncs);CHKERRQ(ierr); 313*c4762a1bSJed Brown funcs[0] = divergence_sq; 314*c4762a1bSJed Brown afuncs[0] = linear2; 315*c4762a1bSJed Brown afuncs[1] = linear; 316*c4762a1bSJed Brown ierr = DMGetLocalVector(dmIn, &lu);CHKERRQ(ierr); 317*c4762a1bSJed Brown ierr = PetscStrcpy(lname, "Local MultiField Input ");CHKERRQ(ierr); 318*c4762a1bSJed Brown ierr = PetscStrcat(lname, name);CHKERRQ(ierr); 319*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) lu, lname);CHKERRQ(ierr); 320*c4762a1bSJed Brown if (!label) {ierr = DMProjectFunctionLocal(dmIn, 0.0, afuncs, NULL, INSERT_VALUES, lu);CHKERRQ(ierr);} 321*c4762a1bSJed Brown else {ierr = DMProjectFunctionLabelLocal(dmIn, 0.0, label, 1, val, 0, NULL, afuncs, NULL, INSERT_VALUES, lu);CHKERRQ(ierr);} 322*c4762a1bSJed Brown ierr = VecViewFromOptions(lu, NULL, "-local_input_view");CHKERRQ(ierr); 323*c4762a1bSJed Brown ierr = DMGetLocalVector(dm, &lx);CHKERRQ(ierr); 324*c4762a1bSJed Brown ierr = PetscStrcpy(lname, "Local MultiField ");CHKERRQ(ierr); 325*c4762a1bSJed Brown ierr = PetscStrcat(lname, name);CHKERRQ(ierr); 326*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) lx, lname);CHKERRQ(ierr); 327*c4762a1bSJed Brown if (!label) {ierr = DMProjectFieldLocal(dm, 0.0, lu, funcs, INSERT_VALUES, lx);CHKERRQ(ierr);} 328*c4762a1bSJed Brown else {ierr = DMProjectFieldLabelLocal(dm, 0.0, label, 1, val, 0, NULL, lu, funcs, INSERT_VALUES, lx);CHKERRQ(ierr);} 329*c4762a1bSJed Brown ierr = VecViewFromOptions(lx, NULL, "-local_field_view");CHKERRQ(ierr); 330*c4762a1bSJed Brown ierr = DMRestoreLocalVector(dm, &lx);CHKERRQ(ierr); 331*c4762a1bSJed Brown ierr = DMRestoreLocalVector(dmIn, &lu);CHKERRQ(ierr); 332*c4762a1bSJed Brown ierr = PetscFree2(funcs, afuncs);CHKERRQ(ierr); 333*c4762a1bSJed Brown if (dmAux) { 334*c4762a1bSJed Brown ierr = PetscObjectCompose((PetscObject) dm, "dmAux", NULL);CHKERRQ(ierr); 335*c4762a1bSJed Brown ierr = PetscObjectCompose((PetscObject) dm, "A", NULL);CHKERRQ(ierr); 336*c4762a1bSJed Brown } 337*c4762a1bSJed Brown PetscFunctionReturn(0); 338*c4762a1bSJed Brown } 339*c4762a1bSJed Brown 340*c4762a1bSJed Brown int main(int argc, char **argv) 341*c4762a1bSJed Brown { 342*c4762a1bSJed Brown DM dm, subdm, auxdm; 343*c4762a1bSJed Brown Vec la; 344*c4762a1bSJed Brown AppCtx user; 345*c4762a1bSJed Brown PetscErrorCode ierr; 346*c4762a1bSJed Brown 347*c4762a1bSJed Brown ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr; 348*c4762a1bSJed Brown ierr = ProcessOptions(&user);CHKERRQ(ierr); 349*c4762a1bSJed Brown ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr); 350*c4762a1bSJed Brown ierr = SetupDiscretization(dm, user.dim, user.cellSimplex, &user);CHKERRQ(ierr); 351*c4762a1bSJed Brown /* Volumetric Mesh Projection */ 352*c4762a1bSJed Brown if (!user.multifield) { 353*c4762a1bSJed Brown ierr = TestFunctionProjection(dm, NULL, NULL, NULL, "Volumetric Primary", &user);CHKERRQ(ierr); 354*c4762a1bSJed Brown ierr = TestFieldProjection(dm, NULL, NULL, NULL, "Volumetric Primary", &user);CHKERRQ(ierr); 355*c4762a1bSJed Brown } else { 356*c4762a1bSJed Brown DM dmOut; 357*c4762a1bSJed Brown 358*c4762a1bSJed Brown ierr = DMClone(dm, &dmOut);CHKERRQ(ierr); 359*c4762a1bSJed Brown ierr = SetupOutputDiscretization(dmOut, user.dim, user.cellSimplex, &user);CHKERRQ(ierr); 360*c4762a1bSJed Brown ierr = TestFieldProjectionMultiple(dmOut, dm, NULL, NULL, NULL, "Volumetric Primary", &user);CHKERRQ(ierr); 361*c4762a1bSJed Brown ierr = DMDestroy(&dmOut);CHKERRQ(ierr); 362*c4762a1bSJed Brown } 363*c4762a1bSJed Brown if (user.auxfield) { 364*c4762a1bSJed Brown /* Volumetric Mesh Projection with Volumetric Data */ 365*c4762a1bSJed Brown ierr = CreateAuxiliaryData(dm, &auxdm, &la, &user);CHKERRQ(ierr); 366*c4762a1bSJed Brown ierr = TestFunctionProjection(dm, auxdm, NULL, la, "Volumetric Primary and Volumetric Auxiliary", &user);CHKERRQ(ierr); 367*c4762a1bSJed Brown ierr = TestFieldProjection(dm, auxdm, NULL, la, "Volumetric Primary and Volumetric Auxiliary", &user);CHKERRQ(ierr); 368*c4762a1bSJed Brown ierr = VecDestroy(&la);CHKERRQ(ierr); 369*c4762a1bSJed Brown /* Update of Volumetric Auxiliary Data with primary Volumetric Data */ 370*c4762a1bSJed Brown ierr = DMGetLocalVector(dm, &la);CHKERRQ(ierr); 371*c4762a1bSJed Brown ierr = VecSet(la, 1.0);CHKERRQ(ierr); 372*c4762a1bSJed Brown ierr = TestFieldProjection(auxdm, dm, NULL, la, "Volumetric Auxiliary Update with Volumetric Primary", &user);CHKERRQ(ierr); 373*c4762a1bSJed Brown ierr = DMRestoreLocalVector(dm, &la);CHKERRQ(ierr); 374*c4762a1bSJed Brown ierr = DMDestroy(&auxdm);CHKERRQ(ierr); 375*c4762a1bSJed Brown } 376*c4762a1bSJed Brown if (user.subdomain) { 377*c4762a1bSJed Brown DMLabel domLabel; 378*c4762a1bSJed Brown 379*c4762a1bSJed Brown /* Subdomain Mesh Projection */ 380*c4762a1bSJed Brown ierr = CreateSubdomainMesh(dm, &domLabel, &subdm, &user);CHKERRQ(ierr); 381*c4762a1bSJed Brown ierr = TestFunctionProjection(subdm, NULL, NULL, NULL, "Subdomain Primary", &user);CHKERRQ(ierr); 382*c4762a1bSJed Brown ierr = TestFieldProjection(subdm, NULL, NULL, NULL, "Subdomain Primary", &user);CHKERRQ(ierr); 383*c4762a1bSJed Brown if (user.auxfield) { 384*c4762a1bSJed Brown /* Subdomain Mesh Projection with Subdomain Data */ 385*c4762a1bSJed Brown ierr = CreateAuxiliaryData(subdm, &auxdm, &la, &user);CHKERRQ(ierr); 386*c4762a1bSJed Brown ierr = TestFunctionProjection(subdm, auxdm, NULL, la, "Subdomain Primary and Subdomain Auxiliary", &user);CHKERRQ(ierr); 387*c4762a1bSJed Brown ierr = TestFieldProjection(subdm, auxdm, NULL, la, "Subdomain Primary and Subdomain Auxiliary", &user);CHKERRQ(ierr); 388*c4762a1bSJed Brown ierr = VecDestroy(&la);CHKERRQ(ierr); 389*c4762a1bSJed Brown ierr = DMDestroy(&auxdm);CHKERRQ(ierr); 390*c4762a1bSJed Brown /* Subdomain Mesh Projection with Volumetric Data */ 391*c4762a1bSJed Brown ierr = CreateAuxiliaryData(dm, &auxdm, &la, &user);CHKERRQ(ierr); 392*c4762a1bSJed Brown ierr = TestFunctionProjection(subdm, auxdm, NULL, la, "Subdomain Primary and Volumetric Auxiliary", &user);CHKERRQ(ierr); 393*c4762a1bSJed Brown ierr = TestFieldProjection(subdm, auxdm, NULL, la, "Subdomain Primary and Volumetric Auxiliary", &user);CHKERRQ(ierr); 394*c4762a1bSJed Brown ierr = VecDestroy(&la);CHKERRQ(ierr); 395*c4762a1bSJed Brown ierr = DMDestroy(&auxdm);CHKERRQ(ierr); 396*c4762a1bSJed Brown /* Volumetric Mesh Projection with Subdomain Data */ 397*c4762a1bSJed Brown ierr = CreateAuxiliaryData(subdm, &auxdm, &la, &user);CHKERRQ(ierr); 398*c4762a1bSJed Brown ierr = TestFunctionProjection(subdm, auxdm, domLabel, la, "Volumetric Primary and Subdomain Auxiliary", &user);CHKERRQ(ierr); 399*c4762a1bSJed Brown ierr = TestFieldProjection(subdm, auxdm, domLabel, la, "Volumetric Primary and Subdomain Auxiliary", &user);CHKERRQ(ierr); 400*c4762a1bSJed Brown ierr = VecDestroy(&la);CHKERRQ(ierr); 401*c4762a1bSJed Brown ierr = DMDestroy(&auxdm);CHKERRQ(ierr); 402*c4762a1bSJed Brown } 403*c4762a1bSJed Brown ierr = DMDestroy(&subdm);CHKERRQ(ierr); 404*c4762a1bSJed Brown ierr = DMLabelDestroy(&domLabel);CHKERRQ(ierr); 405*c4762a1bSJed Brown } 406*c4762a1bSJed Brown if (user.submesh) { 407*c4762a1bSJed Brown DMLabel bdLabel; 408*c4762a1bSJed Brown 409*c4762a1bSJed Brown /* Boundary Mesh Projection */ 410*c4762a1bSJed Brown ierr = CreateBoundaryMesh(dm, &bdLabel, &subdm, &user);CHKERRQ(ierr); 411*c4762a1bSJed Brown ierr = TestFunctionProjection(subdm, NULL, NULL, NULL, "Boundary Primary", &user);CHKERRQ(ierr); 412*c4762a1bSJed Brown ierr = TestFieldProjection(subdm, NULL, NULL, NULL, "Boundary Primary", &user);CHKERRQ(ierr); 413*c4762a1bSJed Brown if (user.auxfield) { 414*c4762a1bSJed Brown /* Boundary Mesh Projection with Boundary Data */ 415*c4762a1bSJed Brown ierr = CreateAuxiliaryData(subdm, &auxdm, &la, &user);CHKERRQ(ierr); 416*c4762a1bSJed Brown ierr = TestFunctionProjection(subdm, auxdm, NULL, la, "Boundary Primary and Boundary Auxiliary", &user);CHKERRQ(ierr); 417*c4762a1bSJed Brown ierr = TestFieldProjection(subdm, auxdm, NULL, la, "Boundary Primary and Boundary Auxiliary", &user);CHKERRQ(ierr); 418*c4762a1bSJed Brown ierr = VecDestroy(&la);CHKERRQ(ierr); 419*c4762a1bSJed Brown ierr = DMDestroy(&auxdm);CHKERRQ(ierr); 420*c4762a1bSJed Brown /* Volumetric Mesh Projection with Boundary Data */ 421*c4762a1bSJed Brown ierr = CreateAuxiliaryData(subdm, &auxdm, &la, &user);CHKERRQ(ierr); 422*c4762a1bSJed Brown ierr = TestFunctionProjection(dm, auxdm, bdLabel, la, "Volumetric Primary and Boundary Auxiliary", &user);CHKERRQ(ierr); 423*c4762a1bSJed Brown ierr = TestFieldProjection(dm, auxdm, bdLabel, la, "Volumetric Primary and Boundary Auxiliary", &user);CHKERRQ(ierr); 424*c4762a1bSJed Brown ierr = VecDestroy(&la);CHKERRQ(ierr); 425*c4762a1bSJed Brown ierr = DMDestroy(&auxdm);CHKERRQ(ierr); 426*c4762a1bSJed Brown } 427*c4762a1bSJed Brown ierr = DMLabelDestroy(&bdLabel);CHKERRQ(ierr); 428*c4762a1bSJed Brown ierr = DMDestroy(&subdm);CHKERRQ(ierr); 429*c4762a1bSJed Brown } 430*c4762a1bSJed Brown ierr = DMDestroy(&dm);CHKERRQ(ierr); 431*c4762a1bSJed Brown ierr = PetscFinalize(); 432*c4762a1bSJed Brown return ierr; 433*c4762a1bSJed Brown } 434*c4762a1bSJed Brown 435*c4762a1bSJed Brown /*TEST 436*c4762a1bSJed Brown 437*c4762a1bSJed Brown test: 438*c4762a1bSJed Brown suffix: 0 439*c4762a1bSJed Brown requires: triangle 440*c4762a1bSJed Brown args: -dim 2 -func_view -local_func_view -local_input_view -local_field_view 441*c4762a1bSJed Brown test: 442*c4762a1bSJed Brown suffix: mf_0 443*c4762a1bSJed Brown requires: triangle 444*c4762a1bSJed Brown args: -dim 2 -velocity_petscspace_degree 1 -velocity_petscfe_default_quadrature_order 2 \ 445*c4762a1bSJed Brown -pressure_petscspace_degree 2 -pressure_petscfe_default_quadrature_order 2 \ 446*c4762a1bSJed Brown -multifield -output_petscspace_degree 1 -output_petscfe_default_quadrature_order 2 \ 447*c4762a1bSJed Brown -local_input_view -local_field_view 448*c4762a1bSJed Brown test: 449*c4762a1bSJed Brown suffix: 1 450*c4762a1bSJed Brown requires: triangle 451*c4762a1bSJed Brown args: -dim 2 -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 452*c4762a1bSJed Brown test: 453*c4762a1bSJed Brown suffix: 2 454*c4762a1bSJed Brown requires: triangle 455*c4762a1bSJed Brown args: -dim 2 -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 456*c4762a1bSJed Brown 457*c4762a1bSJed Brown TEST*/ 458*c4762a1bSJed Brown 459*c4762a1bSJed Brown /* 460*c4762a1bSJed Brown Post-processing wants to project a function of the fields into some FE space 461*c4762a1bSJed Brown - This is DMProjectField() 462*c4762a1bSJed Brown - What about changing the number of components of the output, like displacement to stress? Aux vars 463*c4762a1bSJed Brown 464*c4762a1bSJed Brown Update of state variables 465*c4762a1bSJed Brown - This is DMProjectField(), but solution must be the aux var 466*c4762a1bSJed Brown */ 467