xref: /petsc/src/dm/dt/fe/tests/ex3.c (revision 503c0ea9b45bcfbcebbb1ea5341243bbc69f0bea)
1 static const char help[] = "Tests for determining whether a new finite element works";
2 
3 /*
4   Use -interpolation_view and -l2_projection_view to look at the interpolants.
5 */
6 
7 #include <petscdmplex.h>
8 #include <petscfe.h>
9 #include <petscds.h>
10 #include <petscsnes.h>
11 
12 static void constant(PetscInt dim, PetscInt Nf, PetscInt NfAux,
13                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
14                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
15                      PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
16 {
17   const PetscInt Nc = uOff[1] - uOff[0];
18   for (PetscInt c = 0; c < Nc; ++c) f0[c] += 5.;
19 }
20 
21 static void linear(PetscInt dim, PetscInt Nf, PetscInt NfAux,
22                    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
23                    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
24                    PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
25 {
26   const PetscInt Nc = uOff[1] - uOff[0];
27   for (PetscInt c = 0; c < Nc; ++c) f0[c] += 5.*x[c];
28 }
29 
30 static void quadratic(PetscInt dim, PetscInt Nf, PetscInt NfAux,
31                       const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
32                       const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
33                       PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
34 {
35   const PetscInt Nc = uOff[1] - uOff[0];
36   for (PetscInt c = 0; c < Nc; ++c) f0[c] += 5.*x[c]*x[c];
37 }
38 
39 static void trig(PetscInt dim, PetscInt Nf, PetscInt NfAux,
40                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
41                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
42                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
43 {
44   const PetscInt Nc = uOff[1] - uOff[0];
45   for (PetscInt c = 0; c < Nc; ++c) f0[c] += PetscCosReal(2.*PETSC_PI*x[c]);
46 }
47 
48 /*
49  The prime basis for the Wheeler-Yotov-Xue prism.
50  */
51 static void prime(PetscInt dim, PetscInt Nf, PetscInt NfAux,
52                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
53                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
54                   PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
55 {
56   PetscReal x = X[0], y = X[1], z = X[2], b = 1 + x + y + z;
57   f0[0] += b + 2.0*x*z + 2.0*y*z + x*y + x*x;
58   f0[1] += b + 2.0*x*z + 2.0*y*z + x*y + y*y;
59   f0[2] += b - 3.0*x*z - 3.0*y*z -   2.0*z*z;
60 }
61 
62 static const char    *names[]     = {"constant", "linear", "quadratic", "trig", "prime"};
63 static PetscPointFunc functions[] = { constant,   linear,   quadratic,   trig,   prime };
64 
65 typedef struct {
66   PetscPointFunc exactSol;
67   PetscReal shear,flatten;
68 } AppCtx;
69 
70 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
71 {
72   char           name[PETSC_MAX_PATH_LEN] = "constant";
73   PetscInt       Nfunc = sizeof(names)/sizeof(char *), i;
74   PetscErrorCode ierr;
75 
76   PetscFunctionBeginUser;
77   options->exactSol = NULL;
78   options->shear    = 0.;
79   options->flatten  = 1.;
80 
81   ierr = PetscOptionsBegin(comm, "", "FE Test Options", "PETSCFE");PetscCall(ierr);
82   PetscCall(PetscOptionsString("-func", "Function to project into space", "", name, name, PETSC_MAX_PATH_LEN, NULL));
83   PetscCall(PetscOptionsReal("-shear", "Factor by which to shear along the x-direction", "", options->shear, &(options->shear), NULL));
84   PetscCall(PetscOptionsReal("-flatten", "Factor by which to flatten", "", options->flatten, &(options->flatten), NULL));
85   ierr = PetscOptionsEnd();PetscCall(ierr);
86 
87   for (i = 0; i < Nfunc; ++i) {
88     PetscBool flg;
89 
90     PetscCall(PetscStrcmp(name, names[i], &flg));
91     if (flg) {options->exactSol = functions[i]; break;}
92   }
93   PetscCheck(options->exactSol, comm, PETSC_ERR_ARG_WRONG, "Invalid test function %s", name);
94   PetscFunctionReturn(0);
95 }
96 
97 /* The exact solution is the negative of the f0 contribution */
98 static PetscErrorCode exactSolution(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
99 {
100   AppCtx  *user    = (AppCtx *) ctx;
101   PetscInt uOff[2] = {0, Nc};
102 
103   user->exactSol(dim, 1, 0, uOff, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, time, x, 0, NULL, u);
104   for (PetscInt c = 0; c < Nc; ++c) u[c] *= -1.;
105   return 0;
106 }
107 
108 static void f0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
109                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
110                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
111                PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
112 {
113   const PetscInt Nc = uOff[1] - uOff[0];
114   for (PetscInt c = 0; c < Nc; ++c) f0[c] += u[c];
115 }
116 
117 static void g0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
118                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
119                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
120                PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
121 {
122   const PetscInt Nc = uOff[1] - uOff[0];
123   for (PetscInt c = 0; c < Nc; ++c) g0[c*Nc+c] = 1.0;
124 }
125 
126 static PetscErrorCode SetupProblem(DM dm, AppCtx *user)
127 {
128   PetscDS        ds;
129   PetscWeakForm  wf;
130 
131   PetscFunctionBeginUser;
132   PetscCall(DMGetDS(dm, &ds));
133   PetscCall(PetscDSGetWeakForm(ds, &wf));
134   PetscCall(PetscWeakFormSetIndexResidual(wf, NULL, 0, 0, 0, 0, f0, 0, NULL));
135   PetscCall(PetscWeakFormSetIndexResidual(wf, NULL, 0, 0, 0, 1, user->exactSol, 0, NULL));
136   PetscCall(PetscWeakFormSetIndexJacobian(wf, NULL, 0, 0, 0, 0, 0, g0, 0, NULL, 0, NULL, 0, NULL));
137   PetscCall(PetscDSSetExactSolution(ds, 0, exactSolution, user));
138   PetscFunctionReturn(0);
139 }
140 
141 static PetscErrorCode SetupDiscretization(DM dm, const char name[], AppCtx *user)
142 {
143   DM             cdm = dm;
144   PetscFE        fe;
145   char           prefix[PETSC_MAX_PATH_LEN];
146 
147   PetscFunctionBeginUser;
148   PetscCall(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name));
149   PetscCall(DMCreateFEDefault(dm, 1, name ? prefix : NULL, -1, &fe));
150   PetscCall(PetscObjectSetName((PetscObject) fe, name ? name : "Solution"));
151   /* Set discretization and boundary conditions for each mesh */
152   PetscCall(DMSetField(dm, 0, NULL, (PetscObject) fe));
153   PetscCall(DMCreateDS(dm));
154   PetscCall(SetupProblem(dm, user));
155   while (cdm) {
156     PetscCall(DMCopyDisc(dm, cdm));
157     PetscCall(DMGetCoarseDM(cdm, &cdm));
158   }
159   PetscCall(PetscFEDestroy(&fe));
160   PetscFunctionReturn(0);
161 }
162 
163 /* This test tells us whether the given function is contained in the approximation space */
164 static PetscErrorCode CheckInterpolation(DM dm, AppCtx *user)
165 {
166   PetscSimplePointFunc exactSol[1];
167   void                *exactCtx[1];
168   PetscDS              ds;
169   Vec                  u;
170   PetscReal            error, tol = PETSC_SMALL;
171   MPI_Comm             comm;
172 
173   PetscFunctionBeginUser;
174   PetscCall(PetscObjectGetComm((PetscObject) dm, &comm));
175   PetscCall(DMGetDS(dm, &ds));
176   PetscCall(DMGetGlobalVector(dm, &u));
177   PetscCall(PetscDSGetExactSolution(ds, 0, &exactSol[0], &exactCtx[0]));
178   PetscCall(DMProjectFunction(dm, 0.0, exactSol, exactCtx, INSERT_ALL_VALUES, u));
179   PetscCall(VecViewFromOptions(u, NULL, "-interpolation_view"));
180   PetscCall(DMComputeL2Diff(dm, 0.0, exactSol, exactCtx, u, &error));
181   PetscCall(DMRestoreGlobalVector(dm, &u));
182   if (error > tol) PetscCall(PetscPrintf(comm, "Interpolation tests FAIL at tolerance %g error %g\n", (double) tol, (double) error));
183   else             PetscCall(PetscPrintf(comm, "Interpolation tests pass at tolerance %g\n", (double) tol));
184   PetscFunctionReturn(0);
185 }
186 
187 /* This test tells us whether the element is unisolvent (the mass matrix has full rank), and what rate of convergence we achieve */
188 static PetscErrorCode CheckL2Projection(DM dm, AppCtx *user)
189 {
190   PetscSimplePointFunc exactSol[1];
191   void                *exactCtx[1];
192   SNES                 snes;
193   PetscDS              ds;
194   Vec                  u;
195   PetscReal            error, tol = PETSC_SMALL;
196   MPI_Comm             comm;
197 
198   PetscFunctionBeginUser;
199   PetscCall(PetscObjectGetComm((PetscObject) dm, &comm));
200   PetscCall(DMGetDS(dm, &ds));
201   PetscCall(DMGetGlobalVector(dm, &u));
202   PetscCall(PetscDSGetExactSolution(ds, 0, &exactSol[0], &exactCtx[0]));
203   PetscCall(SNESCreate(comm, &snes));
204   PetscCall(SNESSetDM(snes, dm));
205   PetscCall(VecSet(u, 0.0));
206   PetscCall(PetscObjectSetName((PetscObject) u, "solution"));
207   PetscCall(DMPlexSetSNESLocalFEM(dm, user, user, user));
208   PetscCall(SNESSetFromOptions(snes));
209   PetscCall(DMSNESCheckFromOptions(snes, u));
210   PetscCall(SNESSolve(snes, NULL, u));
211   PetscCall(SNESDestroy(&snes));
212   PetscCall(VecViewFromOptions(u, NULL, "-l2_projection_view"));
213   PetscCall(DMComputeL2Diff(dm, 0.0, exactSol, exactCtx, u, &error));
214   PetscCall(DMRestoreGlobalVector(dm, &u));
215   if (error > tol) PetscCall(PetscPrintf(comm, "L2 projection tests FAIL at tolerance %g error %g\n", (double) tol, (double) error));
216   else             PetscCall(PetscPrintf(comm, "L2 projection tests pass at tolerance %g\n", (double) tol));
217   PetscFunctionReturn(0);
218 }
219 
220 /* Distorts the mesh by shearing in the x-direction and flattening, factors provided in the options. */
221 static PetscErrorCode DistortMesh(DM dm, AppCtx *user)
222 {
223   Vec            coordinates;
224   PetscScalar   *ca;
225   PetscInt       dE, n, i;
226 
227   PetscFunctionBeginUser;
228   PetscCall(DMGetCoordinateDim(dm, &dE));
229   PetscCall(DMGetCoordinates(dm, &coordinates));
230   PetscCall(VecGetLocalSize(coordinates, &n));
231   PetscCall(VecGetArray(coordinates, &ca));
232   for (i = 0; i < (n/dE); ++i) {
233     ca[i*dE+0] += user->shear*ca[i*dE+0];
234     ca[i*dE+1] *= user->flatten;
235   }
236   PetscCall(VecRestoreArray(coordinates, &ca));
237   PetscFunctionReturn(0);
238 }
239 
240 int main(int argc, char **argv)
241 {
242   DM             dm;
243   AppCtx         user;
244   PetscMPIInt    size;
245 
246   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
247   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
248   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_SUP, "This is a uniprocessor example only.");
249   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
250   PetscCall(DMCreate(PETSC_COMM_WORLD, &dm));
251   PetscCall(DMSetType(dm, DMPLEX));
252   PetscCall(DMSetFromOptions(dm));
253   PetscCall(DistortMesh(dm,&user));
254   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
255   PetscCall(SetupDiscretization(dm, NULL, &user));
256 
257   PetscCall(CheckInterpolation(dm, &user));
258   PetscCall(CheckL2Projection(dm, &user));
259 
260   PetscCall(DMDestroy(&dm));
261   PetscCall(PetscFinalize());
262   return 0;
263 }
264 
265 /*TEST
266 
267   testset:
268     args: -dm_plex_reference_cell_domain -dm_plex_cell triangle -petscspace_degree 1\
269           -snes_error_if_not_converged -ksp_error_if_not_converged -pc_type lu
270 
271     test:
272       suffix: p1_0
273       args: -func {{constant linear}}
274 
275     # Using -dm_refine 2 -convest_num_refine 4 gives convergence rate 2.0
276     test:
277       suffix: p1_1
278       args: -func {{quadratic trig}} \
279             -snes_convergence_estimate -convest_num_refine 2
280 
281   testset:
282     requires: !complex double
283     args: -dm_plex_reference_cell_domain -dm_plex_cell triangular_prism \
284             -petscspace_type sum \
285             -petscspace_variables 3 \
286             -petscspace_components 3 \
287             -petscspace_sum_spaces 2 \
288             -petscspace_sum_concatenate false \
289               -sumcomp_0_petscspace_variables 3 \
290               -sumcomp_0_petscspace_components 3 \
291               -sumcomp_0_petscspace_degree 1 \
292               -sumcomp_1_petscspace_variables 3 \
293               -sumcomp_1_petscspace_components 3 \
294               -sumcomp_1_petscspace_type wxy \
295             -petscdualspace_form_degree 0 \
296             -petscdualspace_order 1 \
297             -petscdualspace_components 3 \
298           -snes_error_if_not_converged -ksp_error_if_not_converged -pc_type lu
299 
300     test:
301       suffix: wxy_0
302       args: -func constant
303 
304     test:
305       suffix: wxy_1
306       args: -func linear
307 
308     test:
309       suffix: wxy_2
310       args: -func prime
311 
312     test:
313       suffix: wxy_3
314       args: -func linear -shear 1 -flatten 1e-5
315 
316 TEST*/
317