xref: /petsc/src/dm/impls/plex/tests/ex23.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
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 */
linear(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)14*2a8381b2SBarry Smith static PetscErrorCode linear(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
15d71ae5a4SJacob Faibussowitsch {
16c4762a1bSJed Brown   PetscInt c;
17c4762a1bSJed Brown   for (c = 0; c < Nc; ++c) u[c] = (x[0] + x[1]) * Nc + c;
183ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
19c4762a1bSJed Brown }
20c4762a1bSJed Brown 
21c4762a1bSJed Brown /* {x, y, z} */
linear2(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)22*2a8381b2SBarry Smith static PetscErrorCode linear2(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
23d71ae5a4SJacob Faibussowitsch {
24c4762a1bSJed Brown   PetscInt c;
25c4762a1bSJed Brown   for (c = 0; c < Nc; ++c) u[c] = x[c];
263ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
27c4762a1bSJed Brown }
28c4762a1bSJed Brown 
29c4762a1bSJed Brown /* {u_x, u_y, u_z} */
linear_vector(PetscInt dim,PetscInt Nf,PetscInt NfAux,const PetscInt uOff[],const PetscInt uOff_x[],const PetscScalar u[],const PetscScalar u_t[],const PetscScalar u_x[],const PetscInt aOff[],const PetscInt aOff_x[],const PetscScalar a[],const PetscScalar a_t[],const PetscScalar a_x[],PetscReal t,const PetscReal x[],PetscInt numConstants,const PetscScalar constants[],PetscScalar f[])30d71ae5a4SJacob Faibussowitsch static void linear_vector(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f[])
31d71ae5a4SJacob Faibussowitsch {
32c4762a1bSJed Brown   PetscInt d;
33c4762a1bSJed Brown   for (d = 0; d < uOff[1] - uOff[0]; ++d) f[d] = u[d + uOff[0]];
34c4762a1bSJed Brown }
35c4762a1bSJed Brown 
36c4762a1bSJed Brown /* p */
linear_scalar(PetscInt dim,PetscInt Nf,PetscInt NfAux,const PetscInt uOff[],const PetscInt uOff_x[],const PetscScalar u[],const PetscScalar u_t[],const PetscScalar u_x[],const PetscInt aOff[],const PetscInt aOff_x[],const PetscScalar a[],const PetscScalar a_t[],const PetscScalar a_x[],PetscReal t,const PetscReal x[],PetscInt numConstants,const PetscScalar constants[],PetscScalar f[])37d71ae5a4SJacob Faibussowitsch static void linear_scalar(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f[])
38d71ae5a4SJacob Faibussowitsch {
39c4762a1bSJed Brown   f[0] = u[uOff[1]];
40c4762a1bSJed Brown }
41c4762a1bSJed Brown 
42c4762a1bSJed Brown /* {div u, p^2} */
divergence_sq(PetscInt dim,PetscInt Nf,PetscInt NfAux,const PetscInt uOff[],const PetscInt uOff_x[],const PetscScalar u[],const PetscScalar u_t[],const PetscScalar u_x[],const PetscInt aOff[],const PetscInt aOff_x[],const PetscScalar a[],const PetscScalar a_t[],const PetscScalar a_x[],PetscReal t,const PetscReal x[],PetscInt numConstants,const PetscScalar constants[],PetscScalar f[])43d71ae5a4SJacob Faibussowitsch static void divergence_sq(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f[])
44d71ae5a4SJacob Faibussowitsch {
45c4762a1bSJed Brown   PetscInt d;
46c4762a1bSJed Brown   f[0] = 0.0;
47c4762a1bSJed Brown   for (d = 0; d < dim; ++d) f[0] += u_x[uOff_x[0] + d * dim + d];
48c4762a1bSJed Brown   f[1] = PetscSqr(u[uOff[1]]);
49c4762a1bSJed Brown }
50c4762a1bSJed Brown 
ProcessOptions(AppCtx * options)51d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(AppCtx *options)
52d71ae5a4SJacob Faibussowitsch {
53c4762a1bSJed Brown   PetscFunctionBegin;
54c4762a1bSJed Brown   options->multifield = PETSC_FALSE;
55c4762a1bSJed Brown   options->subdomain  = PETSC_FALSE;
56c4762a1bSJed Brown   options->submesh    = PETSC_FALSE;
57c4762a1bSJed Brown   options->auxfield   = PETSC_FALSE;
58c4762a1bSJed Brown 
59d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_SELF, "", "Meshing Problem Options", "DMPLEX");
609566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-multifield", "Flag for trying different numbers of input/output fields", "ex23.c", options->multifield, &options->multifield, NULL));
619566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-subdomain", "Flag for trying volumetric submesh", "ex23.c", options->subdomain, &options->subdomain, NULL));
629566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-submesh", "Flag for trying boundary submesh", "ex23.c", options->submesh, &options->submesh, NULL));
639566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-auxfield", "Flag for trying auxiliary fields", "ex23.c", options->auxfield, &options->auxfield, NULL));
64d0609cedSBarry Smith   PetscOptionsEnd();
653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
66c4762a1bSJed Brown }
67c4762a1bSJed Brown 
CreateMesh(MPI_Comm comm,AppCtx * user,DM * dm)68d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
69d71ae5a4SJacob Faibussowitsch {
70c4762a1bSJed Brown   PetscFunctionBegin;
719566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
729566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
739566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
749566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-orig_dm_view"));
753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
76c4762a1bSJed Brown }
77c4762a1bSJed Brown 
SetupDiscretization(DM dm,PetscInt dim,PetscBool simplex,AppCtx * user)78d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupDiscretization(DM dm, PetscInt dim, PetscBool simplex, AppCtx *user)
79d71ae5a4SJacob Faibussowitsch {
80c4762a1bSJed Brown   PetscFE  fe;
81c4762a1bSJed Brown   MPI_Comm comm;
82c4762a1bSJed Brown 
83c4762a1bSJed Brown   PetscFunctionBeginUser;
849566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
859566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(comm, dim, dim, simplex, "velocity_", -1, &fe));
869566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)fe, "velocity"));
879566063dSJacob Faibussowitsch   PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe));
889566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe));
899566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, "pressure_", -1, &fe));
909566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)fe, "pressure"));
919566063dSJacob Faibussowitsch   PetscCall(DMSetField(dm, 1, NULL, (PetscObject)fe));
929566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe));
939566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dm));
943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
95c4762a1bSJed Brown }
96c4762a1bSJed Brown 
SetupOutputDiscretization(DM dm,PetscInt dim,PetscBool simplex,AppCtx * user)97d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupOutputDiscretization(DM dm, PetscInt dim, PetscBool simplex, AppCtx *user)
98d71ae5a4SJacob Faibussowitsch {
99c4762a1bSJed Brown   PetscFE  fe;
100c4762a1bSJed Brown   MPI_Comm comm;
101c4762a1bSJed Brown 
102c4762a1bSJed Brown   PetscFunctionBeginUser;
1039566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
1049566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(comm, dim, dim, simplex, "output_", -1, &fe));
1059566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)fe, "output"));
1069566063dSJacob Faibussowitsch   PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe));
1079566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe));
1089566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dm));
1093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
110c4762a1bSJed Brown }
111c4762a1bSJed Brown 
CreateSubdomainMesh(DM dm,DMLabel * domLabel,DM * subdm,AppCtx * user)112d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateSubdomainMesh(DM dm, DMLabel *domLabel, DM *subdm, AppCtx *user)
113d71ae5a4SJacob Faibussowitsch {
114c4762a1bSJed Brown   DMLabel   label;
11530602db0SMatthew G. Knepley   PetscBool simplex;
116c4762a1bSJed Brown   PetscInt  dim, cStart, cEnd, c;
117c4762a1bSJed Brown 
118c4762a1bSJed Brown   PetscFunctionBeginUser;
1199566063dSJacob Faibussowitsch   PetscCall(DMPlexIsSimplex(dm, &simplex));
1209566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
1219566063dSJacob Faibussowitsch   PetscCall(DMLabelCreate(PETSC_COMM_SELF, "subdomain", &label));
1229566063dSJacob Faibussowitsch   for (c = cStart + (cEnd - cStart) / 2; c < cEnd; ++c) PetscCall(DMLabelSetValue(label, c, 1));
12371f1c950SStefano Zampini   PetscCall(DMPlexFilter(dm, label, 1, PETSC_FALSE, PETSC_FALSE, PetscObjectComm((PetscObject)dm), NULL, subdm));
1249566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(*subdm, &dim));
1259566063dSJacob Faibussowitsch   PetscCall(SetupDiscretization(*subdm, dim, simplex, user));
1269566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)*subdm, "subdomain"));
1279566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*subdm, NULL, "-sub_dm_view"));
128c4762a1bSJed Brown   if (domLabel) *domLabel = label;
1299566063dSJacob Faibussowitsch   else PetscCall(DMLabelDestroy(&label));
1303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
131c4762a1bSJed Brown }
132c4762a1bSJed Brown 
CreateBoundaryMesh(DM dm,DMLabel * bdLabel,DM * subdm,AppCtx * user)133d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateBoundaryMesh(DM dm, DMLabel *bdLabel, DM *subdm, AppCtx *user)
134d71ae5a4SJacob Faibussowitsch {
135c4762a1bSJed Brown   DMLabel   label;
13630602db0SMatthew G. Knepley   PetscBool simplex;
137c4762a1bSJed Brown   PetscInt  dim;
138c4762a1bSJed Brown 
139c4762a1bSJed Brown   PetscFunctionBeginUser;
1409566063dSJacob Faibussowitsch   PetscCall(DMPlexIsSimplex(dm, &simplex));
1419566063dSJacob Faibussowitsch   PetscCall(DMLabelCreate(PETSC_COMM_SELF, "sub", &label));
1429566063dSJacob Faibussowitsch   PetscCall(DMPlexMarkBoundaryFaces(dm, 1, label));
1439566063dSJacob Faibussowitsch   PetscCall(DMPlexLabelComplete(dm, label));
1449566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateSubmesh(dm, label, 1, PETSC_TRUE, subdm));
1459566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(*subdm, &dim));
1469566063dSJacob Faibussowitsch   PetscCall(SetupDiscretization(*subdm, dim, simplex, user));
1479566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)*subdm, "boundary"));
1489566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*subdm, NULL, "-sub_dm_view"));
149c4762a1bSJed Brown   if (bdLabel) *bdLabel = label;
1509566063dSJacob Faibussowitsch   else PetscCall(DMLabelDestroy(&label));
1513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
152c4762a1bSJed Brown }
153c4762a1bSJed Brown 
CreateAuxiliaryVec(DM dm,DM * auxdm,Vec * la,AppCtx * user)154d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateAuxiliaryVec(DM dm, DM *auxdm, Vec *la, AppCtx *user)
155d71ae5a4SJacob Faibussowitsch {
156c4762a1bSJed Brown   PetscErrorCode (**afuncs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *);
15730602db0SMatthew G. Knepley   PetscBool simplex;
158c4762a1bSJed Brown   PetscInt  dim, Nf, f;
159c4762a1bSJed Brown 
160c4762a1bSJed Brown   PetscFunctionBeginUser;
1619566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
1629566063dSJacob Faibussowitsch   PetscCall(DMPlexIsSimplex(dm, &simplex));
1639566063dSJacob Faibussowitsch   PetscCall(DMGetNumFields(dm, &Nf));
1649566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(Nf, &afuncs));
165c4762a1bSJed Brown   for (f = 0; f < Nf; ++f) afuncs[f] = linear;
1669566063dSJacob Faibussowitsch   PetscCall(DMClone(dm, auxdm));
1679566063dSJacob Faibussowitsch   PetscCall(SetupDiscretization(*auxdm, dim, simplex, user));
1689566063dSJacob Faibussowitsch   PetscCall(DMCreateLocalVector(*auxdm, la));
1699566063dSJacob Faibussowitsch   PetscCall(DMProjectFunctionLocal(dm, 0.0, afuncs, NULL, INSERT_VALUES, *la));
1709566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(*la, NULL, "-local_aux_view"));
1719566063dSJacob Faibussowitsch   PetscCall(PetscFree(afuncs));
1723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
173c4762a1bSJed Brown }
174c4762a1bSJed Brown 
TestFunctionProjection(DM dm,DM dmAux,DMLabel label,Vec la,const char name[],AppCtx * user)175d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestFunctionProjection(DM dm, DM dmAux, DMLabel label, Vec la, const char name[], AppCtx *user)
176d71ae5a4SJacob Faibussowitsch {
177c4762a1bSJed Brown   PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *);
178c4762a1bSJed Brown   Vec      x, lx;
179c4762a1bSJed Brown   PetscInt Nf, f;
180c4762a1bSJed Brown   PetscInt val[1] = {1};
181c4762a1bSJed Brown   char     lname[PETSC_MAX_PATH_LEN];
182c4762a1bSJed Brown 
183c4762a1bSJed Brown   PetscFunctionBeginUser;
1849566063dSJacob Faibussowitsch   if (dmAux) PetscCall(DMSetAuxiliaryVec(dm, NULL, 0, 0, la));
1859566063dSJacob Faibussowitsch   PetscCall(DMGetNumFields(dm, &Nf));
1869566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(Nf, &funcs));
187c4762a1bSJed Brown   for (f = 0; f < Nf; ++f) funcs[f] = linear;
1889566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dm, &x));
189c6a7a370SJeremy L Thompson   PetscCall(PetscStrncpy(lname, "Function ", sizeof(lname)));
190c6a7a370SJeremy L Thompson   PetscCall(PetscStrlcat(lname, name, sizeof(lname)));
1919566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)x, lname));
1929566063dSJacob Faibussowitsch   if (!label) PetscCall(DMProjectFunction(dm, 0.0, funcs, NULL, INSERT_VALUES, x));
1939566063dSJacob Faibussowitsch   else PetscCall(DMProjectFunctionLabel(dm, 0.0, label, 1, val, 0, NULL, funcs, NULL, INSERT_VALUES, x));
1949566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(x, NULL, "-func_view"));
1959566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dm, &x));
1969566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &lx));
197c6a7a370SJeremy L Thompson   PetscCall(PetscStrncpy(lname, "Local Function ", sizeof(lname)));
198c6a7a370SJeremy L Thompson   PetscCall(PetscStrlcat(lname, name, sizeof(lname)));
1999566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)lx, lname));
2009566063dSJacob Faibussowitsch   if (!label) PetscCall(DMProjectFunctionLocal(dm, 0.0, funcs, NULL, INSERT_VALUES, lx));
2019566063dSJacob Faibussowitsch   else PetscCall(DMProjectFunctionLabelLocal(dm, 0.0, label, 1, val, 0, NULL, funcs, NULL, INSERT_VALUES, lx));
2029566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(lx, NULL, "-local_func_view"));
2039566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &lx));
2049566063dSJacob Faibussowitsch   PetscCall(PetscFree(funcs));
2059566063dSJacob Faibussowitsch   if (dmAux) PetscCall(DMSetAuxiliaryVec(dm, NULL, 0, 0, NULL));
2063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
207c4762a1bSJed Brown }
208c4762a1bSJed Brown 
TestFieldProjection(DM dm,DM dmAux,DMLabel label,Vec la,const char name[],AppCtx * user)209d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestFieldProjection(DM dm, DM dmAux, DMLabel label, Vec la, const char name[], AppCtx *user)
210d71ae5a4SJacob Faibussowitsch {
211c4762a1bSJed Brown   PetscErrorCode (**afuncs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *);
2129371c9d4SSatish Balay   void (**funcs)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]);
213c4762a1bSJed Brown   Vec      lx, lu;
214c4762a1bSJed Brown   PetscInt Nf, f;
215c4762a1bSJed Brown   PetscInt val[1] = {1};
216c4762a1bSJed Brown   char     lname[PETSC_MAX_PATH_LEN];
217c4762a1bSJed Brown 
218c4762a1bSJed Brown   PetscFunctionBeginUser;
2199566063dSJacob Faibussowitsch   if (dmAux) PetscCall(DMSetAuxiliaryVec(dm, NULL, 0, 0, la));
2209566063dSJacob Faibussowitsch   PetscCall(DMGetNumFields(dm, &Nf));
2219566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(Nf, &funcs, Nf, &afuncs));
222c4762a1bSJed Brown   for (f = 0; f < Nf; ++f) afuncs[f] = linear;
223c4762a1bSJed Brown   funcs[0] = linear_vector;
224c4762a1bSJed Brown   funcs[1] = linear_scalar;
2259566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &lu));
226c6a7a370SJeremy L Thompson   PetscCall(PetscStrncpy(lname, "Local Field Input ", sizeof(lname)));
227c6a7a370SJeremy L Thompson   PetscCall(PetscStrlcat(lname, name, sizeof(lname)));
2289566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)lu, lname));
2299566063dSJacob Faibussowitsch   if (!label) PetscCall(DMProjectFunctionLocal(dm, 0.0, afuncs, NULL, INSERT_VALUES, lu));
2309566063dSJacob Faibussowitsch   else PetscCall(DMProjectFunctionLabelLocal(dm, 0.0, label, 1, val, 0, NULL, afuncs, NULL, INSERT_VALUES, lu));
2319566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(lu, NULL, "-local_input_view"));
2329566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &lx));
233c6a7a370SJeremy L Thompson   PetscCall(PetscStrncpy(lname, "Local Field ", sizeof(lname)));
234c6a7a370SJeremy L Thompson   PetscCall(PetscStrlcat(lname, name, sizeof(lname)));
2359566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)lx, lname));
2369566063dSJacob Faibussowitsch   if (!label) PetscCall(DMProjectFieldLocal(dm, 0.0, lu, funcs, INSERT_VALUES, lx));
2379566063dSJacob Faibussowitsch   else PetscCall(DMProjectFieldLabelLocal(dm, 0.0, label, 1, val, 0, NULL, lu, funcs, INSERT_VALUES, lx));
2389566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(lx, NULL, "-local_field_view"));
2399566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &lx));
2409566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &lu));
2419566063dSJacob Faibussowitsch   PetscCall(PetscFree2(funcs, afuncs));
2429566063dSJacob Faibussowitsch   if (dmAux) PetscCall(DMSetAuxiliaryVec(dm, NULL, 0, 0, NULL));
2433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
244c4762a1bSJed Brown }
245c4762a1bSJed Brown 
TestFieldProjectionMultiple(DM dm,DM dmIn,DM dmAux,DMLabel label,Vec la,const char name[],AppCtx * user)246d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestFieldProjectionMultiple(DM dm, DM dmIn, DM dmAux, DMLabel label, Vec la, const char name[], AppCtx *user)
247d71ae5a4SJacob Faibussowitsch {
248c4762a1bSJed Brown   PetscErrorCode (**afuncs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *);
2499371c9d4SSatish Balay   void (**funcs)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]);
250c4762a1bSJed Brown   Vec      lx, lu;
251c4762a1bSJed Brown   PetscInt Nf, NfIn;
252c4762a1bSJed Brown   PetscInt val[1] = {1};
253c4762a1bSJed Brown   char     lname[PETSC_MAX_PATH_LEN];
254c4762a1bSJed Brown 
255c4762a1bSJed Brown   PetscFunctionBeginUser;
2569566063dSJacob Faibussowitsch   if (dmAux) PetscCall(DMSetAuxiliaryVec(dm, NULL, 0, 0, la));
2579566063dSJacob Faibussowitsch   PetscCall(DMGetNumFields(dm, &Nf));
2589566063dSJacob Faibussowitsch   PetscCall(DMGetNumFields(dmIn, &NfIn));
2599566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(Nf, &funcs, NfIn, &afuncs));
260c4762a1bSJed Brown   funcs[0]  = divergence_sq;
261c4762a1bSJed Brown   afuncs[0] = linear2;
262c4762a1bSJed Brown   afuncs[1] = linear;
2639566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dmIn, &lu));
264c6a7a370SJeremy L Thompson   PetscCall(PetscStrncpy(lname, "Local MultiField Input ", sizeof(lname)));
265c6a7a370SJeremy L Thompson   PetscCall(PetscStrlcat(lname, name, sizeof(lname)));
2669566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)lu, lname));
2679566063dSJacob Faibussowitsch   if (!label) PetscCall(DMProjectFunctionLocal(dmIn, 0.0, afuncs, NULL, INSERT_VALUES, lu));
2689566063dSJacob Faibussowitsch   else PetscCall(DMProjectFunctionLabelLocal(dmIn, 0.0, label, 1, val, 0, NULL, afuncs, NULL, INSERT_VALUES, lu));
2699566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(lu, NULL, "-local_input_view"));
2709566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(dm, &lx));
271c6a7a370SJeremy L Thompson   PetscCall(PetscStrncpy(lname, "Local MultiField ", sizeof(lname)));
272c6a7a370SJeremy L Thompson   PetscCall(PetscStrlcat(lname, name, sizeof(lname)));
2739566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)lx, lname));
2749566063dSJacob Faibussowitsch   if (!label) PetscCall(DMProjectFieldLocal(dm, 0.0, lu, funcs, INSERT_VALUES, lx));
2759566063dSJacob Faibussowitsch   else PetscCall(DMProjectFieldLabelLocal(dm, 0.0, label, 1, val, 0, NULL, lu, funcs, INSERT_VALUES, lx));
2769566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(lx, NULL, "-local_field_view"));
2779566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dm, &lx));
2789566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(dmIn, &lu));
2799566063dSJacob Faibussowitsch   PetscCall(PetscFree2(funcs, afuncs));
2809566063dSJacob Faibussowitsch   if (dmAux) PetscCall(DMSetAuxiliaryVec(dm, NULL, 0, 0, NULL));
2813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
282c4762a1bSJed Brown }
283c4762a1bSJed Brown 
main(int argc,char ** argv)284d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
285d71ae5a4SJacob Faibussowitsch {
286c4762a1bSJed Brown   DM        dm, subdm, auxdm;
287c4762a1bSJed Brown   Vec       la;
28830602db0SMatthew G. Knepley   PetscInt  dim;
28930602db0SMatthew G. Knepley   PetscBool simplex;
290c4762a1bSJed Brown   AppCtx    user;
291c4762a1bSJed Brown 
292327415f7SBarry Smith   PetscFunctionBeginUser;
2939566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
2949566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(&user));
2959566063dSJacob Faibussowitsch   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
2969566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
2979566063dSJacob Faibussowitsch   PetscCall(DMPlexIsSimplex(dm, &simplex));
2989566063dSJacob Faibussowitsch   PetscCall(SetupDiscretization(dm, dim, simplex, &user));
299c4762a1bSJed Brown   /* Volumetric Mesh Projection */
300c4762a1bSJed Brown   if (!user.multifield) {
3019566063dSJacob Faibussowitsch     PetscCall(TestFunctionProjection(dm, NULL, NULL, NULL, "Volumetric Primary", &user));
3029566063dSJacob Faibussowitsch     PetscCall(TestFieldProjection(dm, NULL, NULL, NULL, "Volumetric Primary", &user));
303c4762a1bSJed Brown   } else {
304c4762a1bSJed Brown     DM dmOut;
305c4762a1bSJed Brown 
3069566063dSJacob Faibussowitsch     PetscCall(DMClone(dm, &dmOut));
3079566063dSJacob Faibussowitsch     PetscCall(SetupOutputDiscretization(dmOut, dim, simplex, &user));
3089566063dSJacob Faibussowitsch     PetscCall(TestFieldProjectionMultiple(dmOut, dm, NULL, NULL, NULL, "Volumetric Primary", &user));
3099566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&dmOut));
310c4762a1bSJed Brown   }
311c4762a1bSJed Brown   if (user.auxfield) {
312c4762a1bSJed Brown     /* Volumetric Mesh Projection with Volumetric Data */
3139566063dSJacob Faibussowitsch     PetscCall(CreateAuxiliaryVec(dm, &auxdm, &la, &user));
3149566063dSJacob Faibussowitsch     PetscCall(TestFunctionProjection(dm, auxdm, NULL, la, "Volumetric Primary and Volumetric Auxiliary", &user));
3159566063dSJacob Faibussowitsch     PetscCall(TestFieldProjection(dm, auxdm, NULL, la, "Volumetric Primary and Volumetric Auxiliary", &user));
3169566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&la));
317c4762a1bSJed Brown     /* Update of Volumetric Auxiliary Data with primary Volumetric Data */
3189566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(dm, &la));
3199566063dSJacob Faibussowitsch     PetscCall(VecSet(la, 1.0));
3209566063dSJacob Faibussowitsch     PetscCall(TestFieldProjection(auxdm, dm, NULL, la, "Volumetric Auxiliary Update with Volumetric Primary", &user));
3219566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(dm, &la));
3229566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&auxdm));
323c4762a1bSJed Brown   }
324c4762a1bSJed Brown   if (user.subdomain) {
325c4762a1bSJed Brown     DMLabel domLabel;
326c4762a1bSJed Brown 
327c4762a1bSJed Brown     /* Subdomain Mesh Projection */
3289566063dSJacob Faibussowitsch     PetscCall(CreateSubdomainMesh(dm, &domLabel, &subdm, &user));
3299566063dSJacob Faibussowitsch     PetscCall(TestFunctionProjection(subdm, NULL, NULL, NULL, "Subdomain Primary", &user));
3309566063dSJacob Faibussowitsch     PetscCall(TestFieldProjection(subdm, NULL, NULL, NULL, "Subdomain Primary", &user));
331c4762a1bSJed Brown     if (user.auxfield) {
332c4762a1bSJed Brown       /* Subdomain Mesh Projection with Subdomain Data */
3339566063dSJacob Faibussowitsch       PetscCall(CreateAuxiliaryVec(subdm, &auxdm, &la, &user));
3349566063dSJacob Faibussowitsch       PetscCall(TestFunctionProjection(subdm, auxdm, NULL, la, "Subdomain Primary and Subdomain Auxiliary", &user));
3359566063dSJacob Faibussowitsch       PetscCall(TestFieldProjection(subdm, auxdm, NULL, la, "Subdomain Primary and Subdomain Auxiliary", &user));
3369566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&la));
3379566063dSJacob Faibussowitsch       PetscCall(DMDestroy(&auxdm));
338c4762a1bSJed Brown       /* Subdomain Mesh Projection with Volumetric Data */
3399566063dSJacob Faibussowitsch       PetscCall(CreateAuxiliaryVec(dm, &auxdm, &la, &user));
3409566063dSJacob Faibussowitsch       PetscCall(TestFunctionProjection(subdm, auxdm, NULL, la, "Subdomain Primary and Volumetric Auxiliary", &user));
3419566063dSJacob Faibussowitsch       PetscCall(TestFieldProjection(subdm, auxdm, NULL, la, "Subdomain Primary and Volumetric Auxiliary", &user));
3429566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&la));
3439566063dSJacob Faibussowitsch       PetscCall(DMDestroy(&auxdm));
344c4762a1bSJed Brown       /* Volumetric Mesh Projection with Subdomain Data */
3459566063dSJacob Faibussowitsch       PetscCall(CreateAuxiliaryVec(subdm, &auxdm, &la, &user));
3469566063dSJacob Faibussowitsch       PetscCall(TestFunctionProjection(subdm, auxdm, domLabel, la, "Volumetric Primary and Subdomain Auxiliary", &user));
3479566063dSJacob Faibussowitsch       PetscCall(TestFieldProjection(subdm, auxdm, domLabel, la, "Volumetric Primary and Subdomain Auxiliary", &user));
3489566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&la));
3499566063dSJacob Faibussowitsch       PetscCall(DMDestroy(&auxdm));
350c4762a1bSJed Brown     }
3519566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&subdm));
3529566063dSJacob Faibussowitsch     PetscCall(DMLabelDestroy(&domLabel));
353c4762a1bSJed Brown   }
354c4762a1bSJed Brown   if (user.submesh) {
355c4762a1bSJed Brown     DMLabel bdLabel;
356c4762a1bSJed Brown 
357c4762a1bSJed Brown     /* Boundary Mesh Projection */
3589566063dSJacob Faibussowitsch     PetscCall(CreateBoundaryMesh(dm, &bdLabel, &subdm, &user));
3599566063dSJacob Faibussowitsch     PetscCall(TestFunctionProjection(subdm, NULL, NULL, NULL, "Boundary Primary", &user));
3609566063dSJacob Faibussowitsch     PetscCall(TestFieldProjection(subdm, NULL, NULL, NULL, "Boundary Primary", &user));
361c4762a1bSJed Brown     if (user.auxfield) {
362c4762a1bSJed Brown       /* Boundary Mesh Projection with Boundary Data */
3639566063dSJacob Faibussowitsch       PetscCall(CreateAuxiliaryVec(subdm, &auxdm, &la, &user));
3649566063dSJacob Faibussowitsch       PetscCall(TestFunctionProjection(subdm, auxdm, NULL, la, "Boundary Primary and Boundary Auxiliary", &user));
3659566063dSJacob Faibussowitsch       PetscCall(TestFieldProjection(subdm, auxdm, NULL, la, "Boundary Primary and Boundary Auxiliary", &user));
3669566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&la));
3679566063dSJacob Faibussowitsch       PetscCall(DMDestroy(&auxdm));
368c4762a1bSJed Brown       /* Volumetric Mesh Projection with Boundary Data */
3699566063dSJacob Faibussowitsch       PetscCall(CreateAuxiliaryVec(subdm, &auxdm, &la, &user));
3709566063dSJacob Faibussowitsch       PetscCall(TestFunctionProjection(dm, auxdm, bdLabel, la, "Volumetric Primary and Boundary Auxiliary", &user));
3719566063dSJacob Faibussowitsch       PetscCall(TestFieldProjection(dm, auxdm, bdLabel, la, "Volumetric Primary and Boundary Auxiliary", &user));
3729566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&la));
3739566063dSJacob Faibussowitsch       PetscCall(DMDestroy(&auxdm));
374c4762a1bSJed Brown     }
3759566063dSJacob Faibussowitsch     PetscCall(DMLabelDestroy(&bdLabel));
3769566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&subdm));
377c4762a1bSJed Brown   }
3789566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
3799566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
380b122ec5aSJacob Faibussowitsch   return 0;
381c4762a1bSJed Brown }
382c4762a1bSJed Brown 
383c4762a1bSJed Brown /*TEST
384c4762a1bSJed Brown 
385c4762a1bSJed Brown   test:
386c4762a1bSJed Brown     suffix: 0
387c4762a1bSJed Brown     requires: triangle
38830602db0SMatthew G. Knepley     args: -dm_plex_box_faces 1,1 -func_view -local_func_view -local_input_view -local_field_view
389c4762a1bSJed Brown   test:
390c4762a1bSJed Brown     suffix: mf_0
391c4762a1bSJed Brown     requires: triangle
39230602db0SMatthew G. Knepley     args: -dm_plex_box_faces 1,1 -velocity_petscspace_degree 1 -velocity_petscfe_default_quadrature_order 2 \
393c4762a1bSJed Brown          -pressure_petscspace_degree 2 -pressure_petscfe_default_quadrature_order 2 \
394c4762a1bSJed Brown          -multifield -output_petscspace_degree 1 -output_petscfe_default_quadrature_order 2 \
395c4762a1bSJed Brown          -local_input_view -local_field_view
396c4762a1bSJed Brown   test:
397c4762a1bSJed Brown     suffix: 1
398c4762a1bSJed Brown     requires: triangle
39930602db0SMatthew 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
400c4762a1bSJed Brown   test:
401c4762a1bSJed Brown     suffix: 2
402c4762a1bSJed Brown     requires: triangle
40330602db0SMatthew 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
404c4762a1bSJed Brown 
405c4762a1bSJed Brown TEST*/
406c4762a1bSJed Brown 
407c4762a1bSJed Brown /*
408c4762a1bSJed Brown   Post-processing wants to project a function of the fields into some FE space
409c4762a1bSJed Brown   - This is DMProjectField()
410c4762a1bSJed Brown   - What about changing the number of components of the output, like displacement to stress? Aux vars
411c4762a1bSJed Brown 
412c4762a1bSJed Brown   Update of state variables
413c4762a1bSJed Brown   - This is DMProjectField(), but solution must be the aux var
414c4762a1bSJed Brown */
415