Actual source code: ex10.c
1: static char help[] = "Tests implementation of PetscSpace_Sum by solving the Poisson equations using a PetscSpace_Poly and a PetscSpace_Sum and checking that \
2: solutions agree up to machine precision.\n\n";
4: #include <petscdmplex.h>
5: #include <petscds.h>
6: #include <petscfe.h>
7: #include <petscsnes.h>
8: /* We are solving the system of equations:
9: * \vec{u} = -\grad{p}
10: * \div{u} = f
11: */
13: /* Exact solutions for linear velocity
14: \vec{u} = \vec{x};
15: p = -0.5*(\vec{x} \cdot \vec{x});
16: */
17: static PetscErrorCode linear_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
18: {
19: PetscInt c;
21: for (c = 0; c < Nc; ++c) u[c] = x[c];
22: return PETSC_SUCCESS;
23: }
25: static PetscErrorCode linear_p(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
26: {
27: PetscInt d;
29: u[0] = 0.;
30: for (d = 0; d < dim; ++d) u[0] += -0.5 * x[d] * x[d];
31: return PETSC_SUCCESS;
32: }
34: static PetscErrorCode linear_divu(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
35: {
36: u[0] = dim;
37: return PETSC_SUCCESS;
38: }
40: /* fx_v are the residual functions for the equation \vec{u} = \grad{p}. f0_v is the term <v,u>.*/
41: static void f0_v(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 f0[])
42: {
43: PetscInt i;
45: for (i = 0; i < dim; ++i) f0[i] = u[uOff[0] + i];
46: }
48: /* f1_v is the term <v,-\grad{p}> but we integrate by parts to get <\grad{v}, -p*I> */
49: static void f1_v(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 f1[])
50: {
51: PetscInt c;
53: for (c = 0; c < dim; ++c) {
54: PetscInt d;
56: for (d = 0; d < dim; ++d) f1[c * dim + d] = (c == d) ? -u[uOff[1]] : 0;
57: }
58: }
60: /* Residual function for enforcing \div{u} = f. */
61: static void f0_q_linear(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 f0[])
62: {
63: PetscScalar rhs, divu = 0;
64: PetscInt i;
66: (void)linear_divu(dim, t, x, dim, &rhs, NULL);
67: for (i = 0; i < dim; ++i) divu += u_x[uOff_x[0] + i * dim + i];
68: f0[0] = divu - rhs;
69: }
71: /* Boundary residual. Dirichlet boundary for u means u_bdy=p*n */
72: static void f0_bd_u_linear(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[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
73: {
74: PetscScalar pressure;
75: PetscInt d;
77: (void)linear_p(dim, t, x, dim, &pressure, NULL);
78: for (d = 0; d < dim; ++d) f0[d] = pressure * n[d];
79: }
81: /* gx_yz are the jacobian functions obtained by taking the derivative of the y residual w.r.t z*/
82: static void g0_vu(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, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
83: {
84: PetscInt c;
86: for (c = 0; c < dim; ++c) g0[c * dim + c] = 1.0;
87: }
89: static void g1_qu(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, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
90: {
91: PetscInt c;
93: for (c = 0; c < dim; ++c) g1[c * dim + c] = 1.0;
94: }
96: static void g2_vp(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, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[])
97: {
98: PetscInt c;
100: for (c = 0; c < dim; ++c) g2[c * dim + c] = -1.0;
101: }
103: typedef struct {
104: PetscInt dummy;
105: } UserCtx;
107: static PetscErrorCode CreateMesh(MPI_Comm comm, UserCtx *user, DM *mesh)
108: {
109: PetscFunctionBegin;
110: PetscCall(DMCreate(comm, mesh));
111: PetscCall(DMSetType(*mesh, DMPLEX));
112: PetscCall(DMSetFromOptions(*mesh));
113: PetscCall(DMSetApplicationContext(*mesh, user));
114: PetscCall(DMViewFromOptions(*mesh, NULL, "-dm_view"));
115: PetscFunctionReturn(PETSC_SUCCESS);
116: }
118: /* Setup the system of equations that we wish to solve */
119: static PetscErrorCode SetupProblem(DM dm, UserCtx *user)
120: {
121: PetscDS ds;
122: DMLabel label;
123: PetscWeakForm wf;
124: const PetscInt id = 1;
125: PetscInt bd;
127: PetscFunctionBegin;
128: PetscCall(DMGetDS(dm, &ds));
129: /* All of these are independent of the user's choice of solution */
130: PetscCall(PetscDSSetResidual(ds, 0, f0_v, f1_v));
131: PetscCall(PetscDSSetResidual(ds, 1, f0_q_linear, NULL));
132: PetscCall(PetscDSSetJacobian(ds, 0, 0, g0_vu, NULL, NULL, NULL));
133: PetscCall(PetscDSSetJacobian(ds, 0, 1, NULL, NULL, g2_vp, NULL));
134: PetscCall(PetscDSSetJacobian(ds, 1, 0, NULL, g1_qu, NULL, NULL));
136: PetscCall(DMGetLabel(dm, "marker", &label));
137: PetscCall(PetscDSAddBoundary(ds, DM_BC_NATURAL, "Boundary Integral", label, 1, &id, 0, 0, NULL, (void (*)(void))NULL, NULL, user, &bd));
138: PetscCall(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
139: PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, 1, 0, 0, 0, f0_bd_u_linear, 0, NULL));
141: PetscCall(PetscDSSetExactSolution(ds, 0, linear_u, NULL));
142: PetscCall(PetscDSSetExactSolution(ds, 1, linear_p, NULL));
143: PetscFunctionReturn(PETSC_SUCCESS);
144: }
146: /* Create the finite element spaces we will use for this system */
147: static PetscErrorCode SetupDiscretization(DM mesh, DM mesh_sum, PetscErrorCode (*setup)(DM, UserCtx *), UserCtx *user)
148: {
149: DM cdm = mesh, cdm_sum = mesh_sum;
150: PetscFE u, divu, u_sum, divu_sum;
151: PetscInt dim;
152: PetscBool simplex;
154: PetscFunctionBegin;
155: PetscCall(DMGetDimension(mesh, &dim));
156: PetscCall(DMPlexIsSimplex(mesh, &simplex));
157: /* Create FE objects and give them names so that options can be set from
158: * command line. Each field gets 2 instances (i.e. velocity and velocity_sum)created twice so that we can compare between approaches. */
159: PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)mesh), dim, dim, simplex, "velocity_", -1, &u));
160: PetscCall(PetscObjectSetName((PetscObject)u, "velocity"));
161: PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)mesh_sum), dim, dim, simplex, "velocity_sum_", -1, &u_sum));
162: PetscCall(PetscObjectSetName((PetscObject)u_sum, "velocity_sum"));
163: PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)mesh), dim, 1, simplex, "divu_", -1, &divu));
164: PetscCall(PetscObjectSetName((PetscObject)divu, "divu"));
165: PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)mesh_sum), dim, 1, simplex, "divu_sum_", -1, &divu_sum));
166: PetscCall(PetscObjectSetName((PetscObject)divu_sum, "divu_sum"));
168: PetscCall(PetscFECopyQuadrature(u, divu));
169: PetscCall(PetscFECopyQuadrature(u_sum, divu_sum));
171: /* Associate the FE objects with the mesh and setup the system */
172: PetscCall(DMSetField(mesh, 0, NULL, (PetscObject)u));
173: PetscCall(DMSetField(mesh, 1, NULL, (PetscObject)divu));
174: PetscCall(DMCreateDS(mesh));
175: PetscCall((*setup)(mesh, user));
177: PetscCall(DMSetField(mesh_sum, 0, NULL, (PetscObject)u_sum));
178: PetscCall(DMSetField(mesh_sum, 1, NULL, (PetscObject)divu_sum));
179: PetscCall(DMCreateDS(mesh_sum));
180: PetscCall((*setup)(mesh_sum, user));
182: while (cdm) {
183: PetscCall(DMCopyDisc(mesh, cdm));
184: PetscCall(DMGetCoarseDM(cdm, &cdm));
185: }
187: while (cdm_sum) {
188: PetscCall(DMCopyDisc(mesh_sum, cdm_sum));
189: PetscCall(DMGetCoarseDM(cdm_sum, &cdm_sum));
190: }
192: /* The Mesh now owns the fields, so we can destroy the FEs created in this
193: * function */
194: PetscCall(PetscFEDestroy(&u));
195: PetscCall(PetscFEDestroy(&divu));
196: PetscCall(PetscFEDestroy(&u_sum));
197: PetscCall(PetscFEDestroy(&divu_sum));
198: PetscCall(DMDestroy(&cdm));
199: PetscCall(DMDestroy(&cdm_sum));
200: PetscFunctionReturn(PETSC_SUCCESS);
201: }
203: int main(int argc, char **argv)
204: {
205: UserCtx user;
206: DM dm, dm_sum;
207: SNES snes, snes_sum;
208: Vec u, u_sum;
209: PetscReal errNorm;
210: const PetscReal errTol = PETSC_SMALL;
212: PetscFunctionBeginUser;
213: PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
215: /* Set up a snes for the standard approach, one space with 2 components */
216: PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
217: PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
218: PetscCall(SNESSetDM(snes, dm));
220: /* Set up a snes for the sum space approach, where each subspace of the sum space represents one component */
221: PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes_sum));
222: PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm_sum));
223: PetscCall(SNESSetDM(snes_sum, dm_sum));
224: PetscCall(SetupDiscretization(dm, dm_sum, SetupProblem, &user));
226: /* Set up and solve the system using standard approach. */
227: PetscCall(DMCreateGlobalVector(dm, &u));
228: PetscCall(VecSet(u, 0.0));
229: PetscCall(PetscObjectSetName((PetscObject)u, "solution"));
230: PetscCall(DMPlexSetSNESLocalFEM(dm, &user, &user, &user));
231: PetscCall(SNESSetFromOptions(snes));
232: PetscCall(DMSNESCheckFromOptions(snes, u));
233: PetscCall(SNESSolve(snes, NULL, u));
234: PetscCall(SNESGetSolution(snes, &u));
235: PetscCall(VecViewFromOptions(u, NULL, "-solution_view"));
237: /* Set up and solve the sum space system */
238: PetscCall(DMCreateGlobalVector(dm_sum, &u_sum));
239: PetscCall(VecSet(u_sum, 0.0));
240: PetscCall(PetscObjectSetName((PetscObject)u_sum, "solution_sum"));
241: PetscCall(DMPlexSetSNESLocalFEM(dm_sum, &user, &user, &user));
242: PetscCall(SNESSetFromOptions(snes_sum));
243: PetscCall(DMSNESCheckFromOptions(snes_sum, u_sum));
244: PetscCall(SNESSolve(snes_sum, NULL, u_sum));
245: PetscCall(SNESGetSolution(snes_sum, &u_sum));
246: PetscCall(VecViewFromOptions(u_sum, NULL, "-solution_sum_view"));
248: /* Check if standard solution and sum space solution match to machine precision */
249: PetscCall(VecAXPY(u_sum, -1, u));
250: PetscCall(VecNorm(u_sum, NORM_2, &errNorm));
251: PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Sum space provides the same solution as a regular space: %s", (errNorm < errTol) ? "true" : "false"));
253: /* Cleanup */
254: PetscCall(VecDestroy(&u_sum));
255: PetscCall(VecDestroy(&u));
256: PetscCall(SNESDestroy(&snes_sum));
257: PetscCall(SNESDestroy(&snes));
258: PetscCall(DMDestroy(&dm_sum));
259: PetscCall(DMDestroy(&dm));
260: PetscCall(PetscFinalize());
261: return 0;
262: }
264: /*TEST
265: test:
266: suffix: 2d_lagrange
267: requires: triangle
268: args: -velocity_petscspace_degree 1 \
269: -velocity_petscspace_type poly \
270: -velocity_petscspace_components 2\
271: -velocity_petscdualspace_type lagrange \
272: -divu_petscspace_degree 0 \
273: -divu_petscspace_type poly \
274: -divu_petscdualspace_lagrange_continuity false \
275: -velocity_sum_petscfe_default_quadrature_order 1 \
276: -velocity_sum_petscspace_degree 1 \
277: -velocity_sum_petscspace_type sum \
278: -velocity_sum_petscspace_variables 2 \
279: -velocity_sum_petscspace_components 2 \
280: -velocity_sum_petscspace_sum_spaces 2 \
281: -velocity_sum_petscspace_sum_concatenate true \
282: -velocity_sum_petscdualspace_type lagrange \
283: -velocity_sum_sumcomp_0_petscspace_type poly \
284: -velocity_sum_sumcomp_0_petscspace_degree 1 \
285: -velocity_sum_sumcomp_0_petscspace_variables 2 \
286: -velocity_sum_sumcomp_0_petscspace_components 1 \
287: -velocity_sum_sumcomp_1_petscspace_type poly \
288: -velocity_sum_sumcomp_1_petscspace_degree 1 \
289: -velocity_sum_sumcomp_1_petscspace_variables 2 \
290: -velocity_sum_sumcomp_1_petscspace_components 1 \
291: -divu_sum_petscspace_degree 0 \
292: -divu_sum_petscspace_type sum \
293: -divu_sum_petscspace_variables 2 \
294: -divu_sum_petscspace_components 1 \
295: -divu_sum_petscspace_sum_spaces 1 \
296: -divu_sum_petscspace_sum_concatenate true \
297: -divu_sum_petscdualspace_lagrange_continuity false \
298: -divu_sum_sumcomp_0_petscspace_type poly \
299: -divu_sum_sumcomp_0_petscspace_degree 0 \
300: -divu_sum_sumcomp_0_petscspace_variables 2 \
301: -divu_sum_sumcomp_0_petscspace_components 1 \
302: -dm_refine 0 \
303: -snes_error_if_not_converged \
304: -ksp_rtol 1e-10 \
305: -ksp_error_if_not_converged \
306: -pc_type fieldsplit\
307: -pc_fieldsplit_type schur\
308: -divu_sum_petscdualspace_lagrange_continuity false \
309: -pc_fieldsplit_schur_precondition full
310: TEST*/