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