xref: /petsc/doc/manual/fe.md (revision 7f296bb328fcd4c99f2da7bfe8ba7ed8a4ebceee)
1*7f296bb3SBarry Smith(ch_fe)=
2*7f296bb3SBarry Smith
3*7f296bb3SBarry Smith# PetscFE: Finite Element Infrastructure in PETSc
4*7f296bb3SBarry Smith
5*7f296bb3SBarry SmithThis chapter introduces the `PetscFE` class, and related subclasses `PetscSpace` and `PetscDualSpace`, which are used to represent finite element discretizations. It details there interaction with the `DMPLEX` class to assemble functions and operators over computational meshes, and produce optimal solvers by constructing multilevel iterations, for example using `PCPATCH`. The idea behind these classes is not to encompass all of computational finite elements, but rather to establish an interface and infrastructure that will allow PETSc to leverage the excellent work done in packages such as Firedrake, FEniCS, LibMesh, and Deal.II.
6*7f296bb3SBarry Smith
7*7f296bb3SBarry Smith## Using Pointwise Functions to Specify Finite Element Problems
8*7f296bb3SBarry Smith
9*7f296bb3SBarry SmithSee the paper about [Unified Residual Evaluation](https://arxiv.org/abs/1309.1204), which explains the use of pointwise evaluation functions to describe weak forms.
10*7f296bb3SBarry Smith
11*7f296bb3SBarry Smith## Describing a particular finite element problem to PETSc
12*7f296bb3SBarry Smith
13*7f296bb3SBarry SmithA finite element problem is presented to PETSc in a series of steps. This is both to facilitate automation, and to allow multiple entry points for user code and external packages since so much finite element software already exists. First, we tell the `DM`, usually a `DMPLEX` or `DMFOREST`, that we have a set of finite element fields which we intended to solve for in our problem, using
14*7f296bb3SBarry Smith
15*7f296bb3SBarry Smith```
16*7f296bb3SBarry SmithDMAddField(dm, NULL, presDisc);
17*7f296bb3SBarry SmithDMAddField(dm, channelLabel, velDisc);
18*7f296bb3SBarry Smith```
19*7f296bb3SBarry Smith
20*7f296bb3SBarry SmithThe second argument is a `DMLabel` object indicating the support of the field on the mesh, with `NULL` indicating the entire domain. Once we have a set of fields, we calls
21*7f296bb3SBarry Smith
22*7f296bb3SBarry Smith```
23*7f296bb3SBarry SmithDMCreateDS(dm);
24*7f296bb3SBarry Smith```
25*7f296bb3SBarry Smith
26*7f296bb3SBarry SmithThis divides the computational domain into subdomains, called *regions* in PETSc, each with a unique set of fields supported on it. These subdomain are identified by labels, and each one has a `PetscDS` object describing the discrete system on that subdomain. There are query functions to get the set of `PetscDS` objects for the `DM`, but it is usually easiest to get the proper `PetscDS` for a given cell using
27*7f296bb3SBarry Smith
28*7f296bb3SBarry Smith```
29*7f296bb3SBarry SmithDMGetCellDS(dm, cell, &ds, NULL);
30*7f296bb3SBarry Smith```
31*7f296bb3SBarry Smith
32*7f296bb3SBarry SmithEach `PetscDS` object has a set of fields, each with a `PetscFE` or `PetscFV` discretization. This allows it to calculate the size of the local discrete approximation, as well as allocate scratch space for all the associated computations. The final thing needed is to specify the actual equations to be enforced on each region. The `PetscDS` contains a `PetscWeakForm` object that holds callback function pointers that define the equations. A simplified, top-level interface through `PetscDS` allows users to quickly define problems for a single region. For example, in
33*7f296bb3SBarry Smith<a href="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/snes/tutorials/ex13.c.html">SNES Tutorial ex13</a>, we define the Poisson problem using
34*7f296bb3SBarry Smith
35*7f296bb3SBarry Smith```
36*7f296bb3SBarry SmithDMLabel  label;
37*7f296bb3SBarry SmithPetscInt f = 0, id = 1;
38*7f296bb3SBarry Smith
39*7f296bb3SBarry SmithPetscDSSetResidual(ds, f, f0_trig_inhomogeneous_u, f1_u);
40*7f296bb3SBarry SmithPetscDSSetJacobian(ds, f, f, NULL, NULL, NULL, g3_uu);
41*7f296bb3SBarry SmithPetscDSSetExactSolution(ds, f, trig_inhomogeneous_u, user);
42*7f296bb3SBarry SmithDMGetLabel(dm, "marker", &label);
43*7f296bb3SBarry SmithDMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, f, 0, NULL, (void (*)(void)) ex, NULL, user, NULL);
44*7f296bb3SBarry Smith```
45*7f296bb3SBarry Smith
46*7f296bb3SBarry Smithwhere the pointwise functions are
47*7f296bb3SBarry Smith
48*7f296bb3SBarry Smith```
49*7f296bb3SBarry Smithstatic PetscErrorCode trig_inhomogeneous_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
50*7f296bb3SBarry Smith{
51*7f296bb3SBarry Smith  PetscInt d;
52*7f296bb3SBarry Smith  *u = 0.0;
53*7f296bb3SBarry Smith  for (d = 0; d < dim; ++d) *u += PetscSinReal(2.0*PETSC_PI*x[d]);
54*7f296bb3SBarry Smith  return 0;
55*7f296bb3SBarry Smith}
56*7f296bb3SBarry Smith
57*7f296bb3SBarry Smithstatic void f0_trig_inhomogeneous_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
58*7f296bb3SBarry Smith                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
59*7f296bb3SBarry Smith                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
60*7f296bb3SBarry Smith                      PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
61*7f296bb3SBarry Smith{
62*7f296bb3SBarry Smith  PetscInt d;
63*7f296bb3SBarry Smith  for (d = 0; d < dim; ++d) f0[0] += -4.0*PetscSqr(PETSC_PI)*PetscSinReal(2.0*PETSC_PI*x[d]);
64*7f296bb3SBarry Smith}
65*7f296bb3SBarry Smith
66*7f296bb3SBarry Smithstatic void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
67*7f296bb3SBarry Smith                 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
68*7f296bb3SBarry Smith                 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
69*7f296bb3SBarry Smith                 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
70*7f296bb3SBarry Smith{
71*7f296bb3SBarry Smith  PetscInt d;
72*7f296bb3SBarry Smith  for (d = 0; d < dim; ++d) f1[d] = u_x[d];
73*7f296bb3SBarry Smith}
74*7f296bb3SBarry Smith
75*7f296bb3SBarry Smithstatic void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
76*7f296bb3SBarry Smith                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
77*7f296bb3SBarry Smith                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
78*7f296bb3SBarry Smith                  PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
79*7f296bb3SBarry Smith{
80*7f296bb3SBarry Smith  PetscInt d;
81*7f296bb3SBarry Smith  for (d = 0; d < dim; ++d) g3[d*dim+d] = 1.0;
82*7f296bb3SBarry Smith}
83*7f296bb3SBarry Smith```
84*7f296bb3SBarry Smith
85*7f296bb3SBarry SmithNotice that we set boundary conditions using `DMAddBoundary`, which will be described later in this chapter. Also we set an exact solution for the field. This can be used to automatically calculate mesh convergence using the `PetscConvEst` object described later in this chapter.
86*7f296bb3SBarry Smith
87*7f296bb3SBarry SmithFor more complex cases with multiple regions, we need to use the `PetscWeakForm` interface directly. The weak form object allows you to set any number of functions for a given field, and also allows functions to be associated with particular subsets of the mesh using labels and label values. We can reproduce the above problem using the *SetIndex* variants which only set a single function at the specified index, rather than a list of functions. We use a `NULL` label and value, meaning that the entire domain is used.
88*7f296bb3SBarry Smith
89*7f296bb3SBarry Smith```
90*7f296bb3SBarry SmithPetscInt f = 0, val = 0;
91*7f296bb3SBarry Smith
92*7f296bb3SBarry SmithPetscDSGetWeakForm(ds, &wf);
93*7f296bb3SBarry SmithPetscWeakFormSetIndexResidual(ds, NULL, val, f, 0, 0, f0_trig_inhomogeneous_u, 0, f1_u);
94*7f296bb3SBarry SmithPetscWeakFormSetIndexJacobian(ds, NULL, val, f, f, 0, 0, NULL, 0, NULL, 0, NULL, 0, g3_uu);
95*7f296bb3SBarry Smith```
96*7f296bb3SBarry Smith
97*7f296bb3SBarry SmithIn <a href="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/snes/tutorials/ex23.c.html">SNES Tutorial ex23</a>, we define the Poisson problem over the entire domain, but in the top half we also define a pressure. The entire problem can be specified as follows
98*7f296bb3SBarry Smith
99*7f296bb3SBarry Smith```
100*7f296bb3SBarry SmithDMGetRegionNumDS(dm, 0, &label, NULL, &ds, NULL);
101*7f296bb3SBarry SmithPetscDSGetWeakForm(ds, &wf);
102*7f296bb3SBarry SmithPetscWeakFormSetIndexResidual(wf, label, 1, 0, 0, 0, f0_quad_u, 0, f1_u);
103*7f296bb3SBarry SmithPetscWeakFormSetIndexJacobian(wf, label, 1, 0, 0, 0, 0, NULL, 0, NULL, 0, NULL, 0, g3_uu);
104*7f296bb3SBarry SmithPetscDSSetExactSolution(ds, 0, quad_u, user);
105*7f296bb3SBarry SmithDMGetRegionNumDS(dm, 1, &label, NULL, &ds, NULL);
106*7f296bb3SBarry SmithPetscDSGetWeakForm(ds, &wf);
107*7f296bb3SBarry SmithPetscWeakFormSetIndexResidual(wf, label, 1, 0, 0, 0, f0_quad_u, 0, f1_u);
108*7f296bb3SBarry SmithPetscWeakFormSetIndexJacobian(wf, label, 1, 0, 0, 0, 0, NULL, 0, NULL, 0, NULL, 0, g3_uu);
109*7f296bb3SBarry SmithPetscWeakFormSetIndexResidual(wf, label, 1, 1, 0, 0, f0_quad_p, 0, NULL);
110*7f296bb3SBarry SmithPetscWeakFormSetIndexJacobian(wf, label, 1, 1, 1, 0, 0, g0_pp, 0, NULL, 0, NULL, 0, NULL);
111*7f296bb3SBarry SmithPetscDSSetExactSolution(ds, 0, quad_u, user);
112*7f296bb3SBarry SmithPetscDSSetExactSolution(ds, 1, quad_p, user);
113*7f296bb3SBarry SmithDMGetLabel(dm, "marker", &label);
114*7f296bb3SBarry SmithDMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) quad_u, NULL, user, NULL);
115*7f296bb3SBarry Smith```
116*7f296bb3SBarry Smith
117*7f296bb3SBarry SmithIn the [PyLith software](https://geodynamics.org/cig/software/pylith/) we use this capability to combine bulk elasticity with a fault constitutive model integrated over the embedded manifolds corresponding to earthquake faults.
118*7f296bb3SBarry Smith
119*7f296bb3SBarry Smith## Assembling finite element residuals and Jacobians
120*7f296bb3SBarry Smith
121*7f296bb3SBarry SmithOnce the pointwise functions are set in each `PetscDS`, mesh traversals can be automatically determined from the `DMLabel` and value specifications in the keys. This default traversal strategy can be activated by attaching the `DM` and default callbacks to a solver
122*7f296bb3SBarry Smith
123*7f296bb3SBarry Smith```
124*7f296bb3SBarry SmithSNESSetDM(snes, dm);
125*7f296bb3SBarry SmithDMPlexSetSNESLocalFEM(dm, &user, &user, &user);
126*7f296bb3SBarry Smith
127*7f296bb3SBarry SmithTSSetDM(ts, dm);
128*7f296bb3SBarry SmithDMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &user);
129*7f296bb3SBarry SmithDMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &user);
130*7f296bb3SBarry SmithDMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &user);
131*7f296bb3SBarry Smith```
132