xref: /petsc/src/dm/impls/plex/tests/ex23.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
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   PetscFunctionBegin;
63c4762a1bSJed Brown   options->multifield  = PETSC_FALSE;
64c4762a1bSJed Brown   options->subdomain   = PETSC_FALSE;
65c4762a1bSJed Brown   options->submesh     = PETSC_FALSE;
66c4762a1bSJed Brown   options->auxfield    = PETSC_FALSE;
67c4762a1bSJed Brown 
68d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_SELF, "", "Meshing Problem Options", "DMPLEX");
699566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-multifield", "Flag for trying different numbers of input/output fields", "ex23.c", options->multifield, &options->multifield, NULL));
709566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-subdomain", "Flag for trying volumetric submesh", "ex23.c", options->subdomain, &options->subdomain, NULL));
719566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-submesh", "Flag for trying boundary submesh", "ex23.c", options->submesh, &options->submesh, NULL));
729566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-auxfield", "Flag for trying auxiliary fields", "ex23.c", options->auxfield, &options->auxfield, NULL));
73d0609cedSBarry Smith   PetscOptionsEnd();
74c4762a1bSJed Brown   PetscFunctionReturn(0);
75c4762a1bSJed Brown }
76c4762a1bSJed Brown 
77c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
78c4762a1bSJed Brown {
79c4762a1bSJed Brown   PetscFunctionBegin;
809566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
819566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
829566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
839566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-orig_dm_view"));
84c4762a1bSJed Brown   PetscFunctionReturn(0);
85c4762a1bSJed Brown }
86c4762a1bSJed Brown 
87c4762a1bSJed Brown static PetscErrorCode SetupDiscretization(DM dm, PetscInt dim, PetscBool simplex, AppCtx *user)
88c4762a1bSJed Brown {
89c4762a1bSJed Brown   PetscFE        fe;
90c4762a1bSJed Brown   MPI_Comm       comm;
91c4762a1bSJed Brown 
92c4762a1bSJed Brown   PetscFunctionBeginUser;
939566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject) dm, &comm));
949566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(comm, dim, dim, simplex, "velocity_", -1, &fe));
959566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) fe, "velocity"));
969566063dSJacob Faibussowitsch   PetscCall(DMSetField(dm, 0, NULL, (PetscObject) fe));
979566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe));
989566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, "pressure_", -1, &fe));
999566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) fe, "pressure"));
1009566063dSJacob Faibussowitsch   PetscCall(DMSetField(dm, 1, NULL, (PetscObject) fe));
1019566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe));
1029566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dm));
103c4762a1bSJed Brown   PetscFunctionReturn(0);
104c4762a1bSJed Brown }
105c4762a1bSJed Brown 
106c4762a1bSJed Brown static PetscErrorCode SetupOutputDiscretization(DM dm, PetscInt dim, PetscBool simplex, AppCtx *user)
107c4762a1bSJed Brown {
108c4762a1bSJed Brown   PetscFE        fe;
109c4762a1bSJed Brown   MPI_Comm       comm;
110c4762a1bSJed Brown 
111c4762a1bSJed Brown   PetscFunctionBeginUser;
1129566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject) dm, &comm));
1139566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(comm, dim, dim, simplex, "output_", -1, &fe));
1149566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) fe, "output"));
1159566063dSJacob Faibussowitsch   PetscCall(DMSetField(dm, 0, NULL, (PetscObject) fe));
1169566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe));
1179566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dm));
118c4762a1bSJed Brown   PetscFunctionReturn(0);
119c4762a1bSJed Brown }
120c4762a1bSJed Brown 
121c4762a1bSJed Brown static PetscErrorCode CreateSubdomainMesh(DM dm, DMLabel *domLabel, DM *subdm, AppCtx *user)
122c4762a1bSJed Brown {
123c4762a1bSJed Brown   DMLabel        label;
12430602db0SMatthew G. Knepley   PetscBool      simplex;
125c4762a1bSJed Brown   PetscInt       dim, cStart, cEnd, c;
126c4762a1bSJed Brown 
127c4762a1bSJed Brown   PetscFunctionBeginUser;
1289566063dSJacob Faibussowitsch   PetscCall(DMPlexIsSimplex(dm, &simplex));
1299566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
1309566063dSJacob Faibussowitsch   PetscCall(DMLabelCreate(PETSC_COMM_SELF, "subdomain", &label));
1319566063dSJacob Faibussowitsch   for (c = cStart + (cEnd-cStart)/2; c < cEnd; ++c) PetscCall(DMLabelSetValue(label, c, 1));
1329566063dSJacob Faibussowitsch   PetscCall(DMPlexFilter(dm, label, 1, subdm));
1339566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(*subdm, &dim));
1349566063dSJacob Faibussowitsch   PetscCall(SetupDiscretization(*subdm, dim, simplex, user));
1359566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) *subdm, "subdomain"));
1369566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*subdm, NULL, "-sub_dm_view"));
137c4762a1bSJed Brown   if (domLabel) *domLabel = label;
1389566063dSJacob Faibussowitsch   else          PetscCall(DMLabelDestroy(&label));
139c4762a1bSJed Brown   PetscFunctionReturn(0);
140c4762a1bSJed Brown }
141c4762a1bSJed Brown 
142c4762a1bSJed Brown static PetscErrorCode CreateBoundaryMesh(DM dm, DMLabel *bdLabel, DM *subdm, AppCtx *user)
143c4762a1bSJed Brown {
144c4762a1bSJed Brown   DMLabel        label;
14530602db0SMatthew G. Knepley   PetscBool      simplex;
146c4762a1bSJed Brown   PetscInt       dim;
147c4762a1bSJed Brown 
148c4762a1bSJed Brown   PetscFunctionBeginUser;
1499566063dSJacob Faibussowitsch   PetscCall(DMPlexIsSimplex(dm, &simplex));
1509566063dSJacob Faibussowitsch   PetscCall(DMLabelCreate(PETSC_COMM_SELF, "sub", &label));
1519566063dSJacob Faibussowitsch   PetscCall(DMPlexMarkBoundaryFaces(dm, 1, label));
1529566063dSJacob Faibussowitsch   PetscCall(DMPlexLabelComplete(dm, label));
1539566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateSubmesh(dm, label, 1, PETSC_TRUE, subdm));
1549566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(*subdm, &dim));
1559566063dSJacob Faibussowitsch   PetscCall(SetupDiscretization(*subdm, dim, simplex, user));
1569566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) *subdm, "boundary"));
1579566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*subdm, NULL, "-sub_dm_view"));
158c4762a1bSJed Brown   if (bdLabel) *bdLabel = label;
1599566063dSJacob Faibussowitsch   else         PetscCall(DMLabelDestroy(&label));
160c4762a1bSJed Brown   PetscFunctionReturn(0);
161c4762a1bSJed Brown }
162c4762a1bSJed Brown 
1639a2a23afSMatthew G. Knepley static PetscErrorCode CreateAuxiliaryVec(DM dm, DM *auxdm, Vec *la, AppCtx *user)
164c4762a1bSJed Brown {
165c4762a1bSJed Brown   PetscErrorCode (**afuncs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *);
16630602db0SMatthew G. Knepley   PetscBool         simplex;
167c4762a1bSJed Brown   PetscInt          dim, Nf, f;
168c4762a1bSJed Brown 
169c4762a1bSJed Brown   PetscFunctionBeginUser;
1709566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
1719566063dSJacob Faibussowitsch   PetscCall(DMPlexIsSimplex(dm, &simplex));
1729566063dSJacob Faibussowitsch   PetscCall(DMGetNumFields(dm, &Nf));
1739566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(Nf, &afuncs));
174c4762a1bSJed Brown   for (f = 0; f < Nf; ++f) afuncs[f]  = linear;
1759566063dSJacob Faibussowitsch   PetscCall(DMClone(dm, auxdm));
1769566063dSJacob Faibussowitsch   PetscCall(SetupDiscretization(*auxdm, dim, simplex, user));
1779566063dSJacob Faibussowitsch   PetscCall(DMCreateLocalVector(*auxdm, la));
1789566063dSJacob Faibussowitsch   PetscCall(DMProjectFunctionLocal(dm, 0.0, afuncs, NULL, INSERT_VALUES, *la));
1799566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(*la, NULL, "-local_aux_view"));
1809566063dSJacob Faibussowitsch   PetscCall(PetscFree(afuncs));
181c4762a1bSJed Brown   PetscFunctionReturn(0);
182c4762a1bSJed Brown }
183c4762a1bSJed Brown 
184c4762a1bSJed Brown static PetscErrorCode TestFunctionProjection(DM dm, DM dmAux, DMLabel label, Vec la, const char name[], AppCtx *user)
185c4762a1bSJed Brown {
186c4762a1bSJed Brown   PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *);
187c4762a1bSJed Brown   Vec               x, lx;
188c4762a1bSJed Brown   PetscInt          Nf, f;
189c4762a1bSJed Brown   PetscInt          val[1] = {1};
190c4762a1bSJed Brown   char              lname[PETSC_MAX_PATH_LEN];
191c4762a1bSJed Brown 
192c4762a1bSJed Brown   PetscFunctionBeginUser;
1939566063dSJacob Faibussowitsch   if (dmAux) PetscCall(DMSetAuxiliaryVec(dm, NULL, 0, 0, la));
1949566063dSJacob Faibussowitsch   PetscCall(DMGetNumFields(dm, &Nf));
1959566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(Nf, &funcs));
196c4762a1bSJed Brown   for (f = 0; f < Nf; ++f) funcs[f] = linear;
1979566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dm, &x));
1989566063dSJacob Faibussowitsch   PetscCall(PetscStrcpy(lname, "Function "));
1999566063dSJacob Faibussowitsch   PetscCall(PetscStrcat(lname, name));
2009566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) x, lname));
2019566063dSJacob Faibussowitsch   if (!label) PetscCall(DMProjectFunction(dm, 0.0, funcs, NULL, INSERT_VALUES, x));
2029566063dSJacob Faibussowitsch   else        PetscCall(DMProjectFunctionLabel(dm, 0.0, label, 1, val, 0, NULL, funcs, NULL, INSERT_VALUES, x));
2039566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(x, NULL, "-func_view"));
2049566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dm, &x));
2059566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &lx));
2069566063dSJacob Faibussowitsch   PetscCall(PetscStrcpy(lname, "Local Function "));
2079566063dSJacob Faibussowitsch   PetscCall(PetscStrcat(lname, name));
2089566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) lx, lname));
2099566063dSJacob Faibussowitsch   if (!label) PetscCall(DMProjectFunctionLocal(dm, 0.0, funcs, NULL, INSERT_VALUES, lx));
2109566063dSJacob Faibussowitsch   else        PetscCall(DMProjectFunctionLabelLocal(dm, 0.0, label, 1, val, 0, NULL, funcs, NULL, INSERT_VALUES, lx));
2119566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(lx, NULL, "-local_func_view"));
2129566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &lx));
2139566063dSJacob Faibussowitsch   PetscCall(PetscFree(funcs));
2149566063dSJacob Faibussowitsch   if (dmAux) PetscCall(DMSetAuxiliaryVec(dm, NULL, 0, 0, NULL));
215c4762a1bSJed Brown   PetscFunctionReturn(0);
216c4762a1bSJed Brown }
217c4762a1bSJed Brown 
218c4762a1bSJed Brown static PetscErrorCode TestFieldProjection(DM dm, DM dmAux, DMLabel label, Vec la, const char name[], AppCtx *user)
219c4762a1bSJed Brown {
220c4762a1bSJed Brown   PetscErrorCode (**afuncs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *);
221c4762a1bSJed Brown   void           (**funcs)(PetscInt, PetscInt, PetscInt,
222c4762a1bSJed Brown                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
223c4762a1bSJed Brown                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
224c4762a1bSJed Brown                            PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]);
225c4762a1bSJed Brown   Vec               lx, lu;
226c4762a1bSJed Brown   PetscInt          Nf, f;
227c4762a1bSJed Brown   PetscInt          val[1] = {1};
228c4762a1bSJed Brown   char              lname[PETSC_MAX_PATH_LEN];
229c4762a1bSJed Brown 
230c4762a1bSJed Brown   PetscFunctionBeginUser;
2319566063dSJacob Faibussowitsch   if (dmAux) PetscCall(DMSetAuxiliaryVec(dm, NULL, 0, 0, la));
2329566063dSJacob Faibussowitsch   PetscCall(DMGetNumFields(dm, &Nf));
2339566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(Nf, &funcs, Nf, &afuncs));
234c4762a1bSJed Brown   for (f = 0; f < Nf; ++f) afuncs[f]  = linear;
235c4762a1bSJed Brown   funcs[0] = linear_vector;
236c4762a1bSJed Brown   funcs[1] = linear_scalar;
2379566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &lu));
2389566063dSJacob Faibussowitsch   PetscCall(PetscStrcpy(lname, "Local Field Input "));
2399566063dSJacob Faibussowitsch   PetscCall(PetscStrcat(lname, name));
2409566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) lu, lname));
2419566063dSJacob Faibussowitsch   if (!label) PetscCall(DMProjectFunctionLocal(dm, 0.0, afuncs, NULL, INSERT_VALUES, lu));
2429566063dSJacob Faibussowitsch   else        PetscCall(DMProjectFunctionLabelLocal(dm, 0.0, label, 1, val, 0, NULL, afuncs, NULL, INSERT_VALUES, lu));
2439566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(lu, NULL, "-local_input_view"));
2449566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &lx));
2459566063dSJacob Faibussowitsch   PetscCall(PetscStrcpy(lname, "Local Field "));
2469566063dSJacob Faibussowitsch   PetscCall(PetscStrcat(lname, name));
2479566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) lx, lname));
2489566063dSJacob Faibussowitsch   if (!label) PetscCall(DMProjectFieldLocal(dm, 0.0, lu, funcs, INSERT_VALUES, lx));
2499566063dSJacob Faibussowitsch   else        PetscCall(DMProjectFieldLabelLocal(dm, 0.0, label, 1, val, 0, NULL, lu, funcs, INSERT_VALUES, lx));
2509566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(lx, NULL, "-local_field_view"));
2519566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &lx));
2529566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &lu));
2539566063dSJacob Faibussowitsch   PetscCall(PetscFree2(funcs, afuncs));
2549566063dSJacob Faibussowitsch   if (dmAux) PetscCall(DMSetAuxiliaryVec(dm, NULL, 0, 0, NULL));
255c4762a1bSJed Brown   PetscFunctionReturn(0);
256c4762a1bSJed Brown }
257c4762a1bSJed Brown 
258c4762a1bSJed Brown static PetscErrorCode TestFieldProjectionMultiple(DM dm, DM dmIn, DM dmAux, DMLabel label, Vec la, const char name[], AppCtx *user)
259c4762a1bSJed Brown {
260c4762a1bSJed Brown   PetscErrorCode (**afuncs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *);
261c4762a1bSJed Brown   void           (**funcs)(PetscInt, PetscInt, PetscInt,
262c4762a1bSJed Brown                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
263c4762a1bSJed Brown                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
264c4762a1bSJed Brown                            PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]);
265c4762a1bSJed Brown   Vec               lx, lu;
266c4762a1bSJed Brown   PetscInt          Nf, NfIn;
267c4762a1bSJed Brown   PetscInt          val[1] = {1};
268c4762a1bSJed Brown   char              lname[PETSC_MAX_PATH_LEN];
269c4762a1bSJed Brown 
270c4762a1bSJed Brown   PetscFunctionBeginUser;
2719566063dSJacob Faibussowitsch   if (dmAux) PetscCall(DMSetAuxiliaryVec(dm, NULL, 0, 0, la));
2729566063dSJacob Faibussowitsch   PetscCall(DMGetNumFields(dm, &Nf));
2739566063dSJacob Faibussowitsch   PetscCall(DMGetNumFields(dmIn, &NfIn));
2749566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(Nf, &funcs, NfIn, &afuncs));
275c4762a1bSJed Brown   funcs[0]  = divergence_sq;
276c4762a1bSJed Brown   afuncs[0] = linear2;
277c4762a1bSJed Brown   afuncs[1] = linear;
2789566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dmIn, &lu));
2799566063dSJacob Faibussowitsch   PetscCall(PetscStrcpy(lname, "Local MultiField Input "));
2809566063dSJacob Faibussowitsch   PetscCall(PetscStrcat(lname, name));
2819566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) lu, lname));
2829566063dSJacob Faibussowitsch   if (!label) PetscCall(DMProjectFunctionLocal(dmIn, 0.0, afuncs, NULL, INSERT_VALUES, lu));
2839566063dSJacob Faibussowitsch   else        PetscCall(DMProjectFunctionLabelLocal(dmIn, 0.0, label, 1, val, 0, NULL, afuncs, NULL, INSERT_VALUES, lu));
2849566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(lu, NULL, "-local_input_view"));
2859566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &lx));
2869566063dSJacob Faibussowitsch   PetscCall(PetscStrcpy(lname, "Local MultiField "));
2879566063dSJacob Faibussowitsch   PetscCall(PetscStrcat(lname, name));
2889566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) lx, lname));
2899566063dSJacob Faibussowitsch   if (!label) PetscCall(DMProjectFieldLocal(dm, 0.0, lu, funcs, INSERT_VALUES, lx));
2909566063dSJacob Faibussowitsch   else        PetscCall(DMProjectFieldLabelLocal(dm, 0.0, label, 1, val, 0, NULL, lu, funcs, INSERT_VALUES, lx));
2919566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(lx, NULL, "-local_field_view"));
2929566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &lx));
2939566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dmIn, &lu));
2949566063dSJacob Faibussowitsch   PetscCall(PetscFree2(funcs, afuncs));
2959566063dSJacob Faibussowitsch   if (dmAux) PetscCall(DMSetAuxiliaryVec(dm, NULL, 0, 0, NULL));
296c4762a1bSJed Brown   PetscFunctionReturn(0);
297c4762a1bSJed Brown }
298c4762a1bSJed Brown 
299c4762a1bSJed Brown int main(int argc, char **argv)
300c4762a1bSJed Brown {
301c4762a1bSJed Brown   DM             dm, subdm, auxdm;
302c4762a1bSJed Brown   Vec            la;
30330602db0SMatthew G. Knepley   PetscInt       dim;
30430602db0SMatthew G. Knepley   PetscBool      simplex;
305c4762a1bSJed Brown   AppCtx         user;
306c4762a1bSJed Brown 
307*327415f7SBarry Smith   PetscFunctionBeginUser;
3089566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL,help));
3099566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(&user));
3109566063dSJacob Faibussowitsch   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
3119566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
3129566063dSJacob Faibussowitsch   PetscCall(DMPlexIsSimplex(dm, &simplex));
3139566063dSJacob Faibussowitsch   PetscCall(SetupDiscretization(dm, dim, simplex, &user));
314c4762a1bSJed Brown   /* Volumetric Mesh Projection */
315c4762a1bSJed Brown   if (!user.multifield) {
3169566063dSJacob Faibussowitsch     PetscCall(TestFunctionProjection(dm, NULL, NULL, NULL, "Volumetric Primary", &user));
3179566063dSJacob Faibussowitsch     PetscCall(TestFieldProjection(dm, NULL, NULL, NULL, "Volumetric Primary", &user));
318c4762a1bSJed Brown   } else {
319c4762a1bSJed Brown     DM dmOut;
320c4762a1bSJed Brown 
3219566063dSJacob Faibussowitsch     PetscCall(DMClone(dm, &dmOut));
3229566063dSJacob Faibussowitsch     PetscCall(SetupOutputDiscretization(dmOut, dim, simplex, &user));
3239566063dSJacob Faibussowitsch     PetscCall(TestFieldProjectionMultiple(dmOut, dm, NULL, NULL, NULL, "Volumetric Primary", &user));
3249566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&dmOut));
325c4762a1bSJed Brown   }
326c4762a1bSJed Brown   if (user.auxfield) {
327c4762a1bSJed Brown     /* Volumetric Mesh Projection with Volumetric Data */
3289566063dSJacob Faibussowitsch     PetscCall(CreateAuxiliaryVec(dm, &auxdm, &la, &user));
3299566063dSJacob Faibussowitsch     PetscCall(TestFunctionProjection(dm, auxdm, NULL, la, "Volumetric Primary and Volumetric Auxiliary", &user));
3309566063dSJacob Faibussowitsch     PetscCall(TestFieldProjection(dm, auxdm, NULL, la, "Volumetric Primary and Volumetric Auxiliary", &user));
3319566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&la));
332c4762a1bSJed Brown     /* Update of Volumetric Auxiliary Data with primary Volumetric Data */
3339566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(dm, &la));
3349566063dSJacob Faibussowitsch     PetscCall(VecSet(la, 1.0));
3359566063dSJacob Faibussowitsch     PetscCall(TestFieldProjection(auxdm, dm, NULL, la, "Volumetric Auxiliary Update with Volumetric Primary", &user));
3369566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(dm, &la));
3379566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&auxdm));
338c4762a1bSJed Brown   }
339c4762a1bSJed Brown   if (user.subdomain) {
340c4762a1bSJed Brown     DMLabel domLabel;
341c4762a1bSJed Brown 
342c4762a1bSJed Brown     /* Subdomain Mesh Projection */
3439566063dSJacob Faibussowitsch     PetscCall(CreateSubdomainMesh(dm, &domLabel, &subdm, &user));
3449566063dSJacob Faibussowitsch     PetscCall(TestFunctionProjection(subdm, NULL, NULL, NULL, "Subdomain Primary", &user));
3459566063dSJacob Faibussowitsch     PetscCall(TestFieldProjection(subdm, NULL, NULL, NULL, "Subdomain Primary", &user));
346c4762a1bSJed Brown     if (user.auxfield) {
347c4762a1bSJed Brown       /* Subdomain Mesh Projection with Subdomain Data */
3489566063dSJacob Faibussowitsch       PetscCall(CreateAuxiliaryVec(subdm, &auxdm, &la, &user));
3499566063dSJacob Faibussowitsch       PetscCall(TestFunctionProjection(subdm, auxdm, NULL, la, "Subdomain Primary and Subdomain Auxiliary", &user));
3509566063dSJacob Faibussowitsch       PetscCall(TestFieldProjection(subdm, auxdm, NULL, la, "Subdomain Primary and Subdomain Auxiliary", &user));
3519566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&la));
3529566063dSJacob Faibussowitsch       PetscCall(DMDestroy(&auxdm));
353c4762a1bSJed Brown       /* Subdomain Mesh Projection with Volumetric Data */
3549566063dSJacob Faibussowitsch       PetscCall(CreateAuxiliaryVec(dm, &auxdm, &la, &user));
3559566063dSJacob Faibussowitsch       PetscCall(TestFunctionProjection(subdm, auxdm, NULL, la, "Subdomain Primary and Volumetric Auxiliary", &user));
3569566063dSJacob Faibussowitsch       PetscCall(TestFieldProjection(subdm, auxdm, NULL, la, "Subdomain Primary and Volumetric Auxiliary", &user));
3579566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&la));
3589566063dSJacob Faibussowitsch       PetscCall(DMDestroy(&auxdm));
359c4762a1bSJed Brown       /* Volumetric Mesh Projection with Subdomain Data */
3609566063dSJacob Faibussowitsch       PetscCall(CreateAuxiliaryVec(subdm, &auxdm, &la, &user));
3619566063dSJacob Faibussowitsch       PetscCall(TestFunctionProjection(subdm, auxdm, domLabel, la, "Volumetric Primary and Subdomain Auxiliary", &user));
3629566063dSJacob Faibussowitsch       PetscCall(TestFieldProjection(subdm, auxdm, domLabel, la, "Volumetric Primary and Subdomain Auxiliary", &user));
3639566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&la));
3649566063dSJacob Faibussowitsch       PetscCall(DMDestroy(&auxdm));
365c4762a1bSJed Brown     }
3669566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&subdm));
3679566063dSJacob Faibussowitsch     PetscCall(DMLabelDestroy(&domLabel));
368c4762a1bSJed Brown   }
369c4762a1bSJed Brown   if (user.submesh) {
370c4762a1bSJed Brown     DMLabel bdLabel;
371c4762a1bSJed Brown 
372c4762a1bSJed Brown     /* Boundary Mesh Projection */
3739566063dSJacob Faibussowitsch     PetscCall(CreateBoundaryMesh(dm, &bdLabel, &subdm, &user));
3749566063dSJacob Faibussowitsch     PetscCall(TestFunctionProjection(subdm, NULL, NULL, NULL, "Boundary Primary", &user));
3759566063dSJacob Faibussowitsch     PetscCall(TestFieldProjection(subdm, NULL, NULL, NULL, "Boundary Primary", &user));
376c4762a1bSJed Brown     if (user.auxfield) {
377c4762a1bSJed Brown       /* Boundary Mesh Projection with Boundary Data */
3789566063dSJacob Faibussowitsch       PetscCall(CreateAuxiliaryVec(subdm, &auxdm, &la, &user));
3799566063dSJacob Faibussowitsch       PetscCall(TestFunctionProjection(subdm, auxdm, NULL, la, "Boundary Primary and Boundary Auxiliary", &user));
3809566063dSJacob Faibussowitsch       PetscCall(TestFieldProjection(subdm, auxdm, NULL, la, "Boundary Primary and Boundary Auxiliary", &user));
3819566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&la));
3829566063dSJacob Faibussowitsch       PetscCall(DMDestroy(&auxdm));
383c4762a1bSJed Brown       /* Volumetric Mesh Projection with Boundary Data */
3849566063dSJacob Faibussowitsch       PetscCall(CreateAuxiliaryVec(subdm, &auxdm, &la, &user));
3859566063dSJacob Faibussowitsch       PetscCall(TestFunctionProjection(dm, auxdm, bdLabel, la, "Volumetric Primary and Boundary Auxiliary", &user));
3869566063dSJacob Faibussowitsch       PetscCall(TestFieldProjection(dm, auxdm, bdLabel, la, "Volumetric Primary and Boundary Auxiliary", &user));
3879566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&la));
3889566063dSJacob Faibussowitsch       PetscCall(DMDestroy(&auxdm));
389c4762a1bSJed Brown     }
3909566063dSJacob Faibussowitsch     PetscCall(DMLabelDestroy(&bdLabel));
3919566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&subdm));
392c4762a1bSJed Brown   }
3939566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
3949566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
395b122ec5aSJacob Faibussowitsch   return 0;
396c4762a1bSJed Brown }
397c4762a1bSJed Brown 
398c4762a1bSJed Brown /*TEST
399c4762a1bSJed Brown 
400c4762a1bSJed Brown   test:
401c4762a1bSJed Brown     suffix: 0
402c4762a1bSJed Brown     requires: triangle
40330602db0SMatthew G. Knepley     args: -dm_plex_box_faces 1,1 -func_view -local_func_view -local_input_view -local_field_view
404c4762a1bSJed Brown   test:
405c4762a1bSJed Brown     suffix: mf_0
406c4762a1bSJed Brown     requires: triangle
40730602db0SMatthew G. Knepley     args: -dm_plex_box_faces 1,1 -velocity_petscspace_degree 1 -velocity_petscfe_default_quadrature_order 2 \
408c4762a1bSJed Brown          -pressure_petscspace_degree 2 -pressure_petscfe_default_quadrature_order 2 \
409c4762a1bSJed Brown          -multifield -output_petscspace_degree 1 -output_petscfe_default_quadrature_order 2 \
410c4762a1bSJed Brown          -local_input_view -local_field_view
411c4762a1bSJed Brown   test:
412c4762a1bSJed Brown     suffix: 1
413c4762a1bSJed Brown     requires: triangle
41430602db0SMatthew 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
415c4762a1bSJed Brown   test:
416c4762a1bSJed Brown     suffix: 2
417c4762a1bSJed Brown     requires: triangle
41830602db0SMatthew 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
419c4762a1bSJed Brown 
420c4762a1bSJed Brown TEST*/
421c4762a1bSJed Brown 
422c4762a1bSJed Brown /*
423c4762a1bSJed Brown   Post-processing wants to project a function of the fields into some FE space
424c4762a1bSJed Brown   - This is DMProjectField()
425c4762a1bSJed Brown   - What about changing the number of components of the output, like displacement to stress? Aux vars
426c4762a1bSJed Brown 
427c4762a1bSJed Brown   Update of state variables
428c4762a1bSJed Brown   - This is DMProjectField(), but solution must be the aux var
429c4762a1bSJed Brown */
430