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