130602db0SMatthew G. Knepley const char help[] = "A test of H-div conforming discretizations on different cell types.\n";
2c4762a1bSJed Brown
3c4762a1bSJed Brown #include <petscdmplex.h>
4c4762a1bSJed Brown #include <petscds.h>
5c4762a1bSJed Brown #include <petscsnes.h>
6c4762a1bSJed Brown #include <petscconvest.h>
7c4762a1bSJed Brown #include <petscfe.h>
8c4762a1bSJed Brown #include <petsc/private/petscfeimpl.h>
9c4762a1bSJed Brown
10c4762a1bSJed Brown /*
110e3d61c9SBarry Smith We are using the system
120e3d61c9SBarry Smith
130e3d61c9SBarry Smith \vec{u} = \vec{\hat{u}}
140e3d61c9SBarry Smith p = \div{\vec{u}} in low degree approximation space
150e3d61c9SBarry Smith d = \div{\vec{u}} - p == 0 in higher degree approximation space
160e3d61c9SBarry Smith
170e3d61c9SBarry Smith That is, we are using the field d to compute the error between \div{\vec{u}}
180e3d61c9SBarry Smith computed in a space 1 degree higher than p and the value of p which is
190e3d61c9SBarry Smith \div{u} computed in the low degree space. If H-div
200e3d61c9SBarry Smith elements are implemented correctly then this should be identically zero since
210e3d61c9SBarry Smith the divergence of a function in H(div) should be exactly representable in L_2
220e3d61c9SBarry Smith by definition.
23c4762a1bSJed Brown */
zero_func(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)24*2a8381b2SBarry Smith static PetscErrorCode zero_func(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
25d71ae5a4SJacob Faibussowitsch {
26c4762a1bSJed Brown PetscInt c;
27c4762a1bSJed Brown for (c = 0; c < Nc; ++c) u[c] = 0;
283ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
29c4762a1bSJed Brown }
30c4762a1bSJed Brown /* Linear Exact Functions
31c4762a1bSJed Brown \vec{u} = \vec{x};
32c4762a1bSJed Brown p = dim;
33c4762a1bSJed Brown */
linear_u(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)34*2a8381b2SBarry Smith static PetscErrorCode linear_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
35d71ae5a4SJacob Faibussowitsch {
36c4762a1bSJed Brown PetscInt c;
37c4762a1bSJed Brown for (c = 0; c < Nc; ++c) u[c] = x[c];
383ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
39c4762a1bSJed Brown }
linear_p(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)40*2a8381b2SBarry Smith static PetscErrorCode linear_p(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
41d71ae5a4SJacob Faibussowitsch {
42c4762a1bSJed Brown u[0] = dim;
433ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
44c4762a1bSJed Brown }
45c4762a1bSJed Brown
46c4762a1bSJed Brown /* Sinusoidal Exact Functions
47c4762a1bSJed Brown * u_i = \sin{2*\pi*x_i}
48c4762a1bSJed Brown * p = \Sum_{i=1}^{dim} 2*\pi*cos{2*\pi*x_i}
49c4762a1bSJed Brown * */
50c4762a1bSJed Brown
sinusoid_u(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)51*2a8381b2SBarry Smith static PetscErrorCode sinusoid_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
52d71ae5a4SJacob Faibussowitsch {
53c4762a1bSJed Brown PetscInt c;
54c4762a1bSJed Brown for (c = 0; c < Nc; ++c) u[c] = PetscSinReal(2 * PETSC_PI * x[c]);
553ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
56c4762a1bSJed Brown }
sinusoid_p(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)57*2a8381b2SBarry Smith static PetscErrorCode sinusoid_p(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
58d71ae5a4SJacob Faibussowitsch {
59c4762a1bSJed Brown PetscInt d;
60c4762a1bSJed Brown u[0] = 0;
61c4762a1bSJed Brown for (d = 0; d < dim; ++d) u[0] += 2 * PETSC_PI * PetscCosReal(2 * PETSC_PI * x[d]);
623ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
63c4762a1bSJed Brown }
64c4762a1bSJed Brown
65c4762a1bSJed Brown /* Pointwise residual for u = u*. Need one of these for each possible u* */
f0_v_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[])66d71ae5a4SJacob Faibussowitsch static void f0_v_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[])
67d71ae5a4SJacob Faibussowitsch {
68c4762a1bSJed Brown PetscInt i;
69c4762a1bSJed Brown PetscScalar *u_rhs;
70c4762a1bSJed Brown
713ba16761SJacob Faibussowitsch PetscCallAbort(PETSC_COMM_SELF, PetscCalloc1(dim, &u_rhs));
723ba16761SJacob Faibussowitsch PetscCallAbort(PETSC_COMM_SELF, linear_u(dim, t, x, dim, u_rhs, NULL));
73c4762a1bSJed Brown for (i = 0; i < dim; ++i) f0[i] = u[uOff[0] + i] - u_rhs[i];
743ba16761SJacob Faibussowitsch PetscCallAbort(PETSC_COMM_SELF, PetscFree(u_rhs));
75c4762a1bSJed Brown }
76c4762a1bSJed Brown
f0_v_sinusoid(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[])77d71ae5a4SJacob Faibussowitsch static void f0_v_sinusoid(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[])
78d71ae5a4SJacob Faibussowitsch {
79c4762a1bSJed Brown PetscInt i;
80c4762a1bSJed Brown PetscScalar *u_rhs;
81c4762a1bSJed Brown
823ba16761SJacob Faibussowitsch PetscCallAbort(PETSC_COMM_SELF, PetscCalloc1(dim, &u_rhs));
833ba16761SJacob Faibussowitsch PetscCallAbort(PETSC_COMM_SELF, sinusoid_u(dim, t, x, dim, u_rhs, NULL));
84c4762a1bSJed Brown for (i = 0; i < dim; ++i) f0[i] = u[uOff[0] + i] - u_rhs[i];
853ba16761SJacob Faibussowitsch PetscCallAbort(PETSC_COMM_SELF, PetscFree(u_rhs));
86c4762a1bSJed Brown }
87c4762a1bSJed Brown
88c4762a1bSJed Brown /* Residual function for enforcing p = \div{u}. */
f0_q(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[])89d71ae5a4SJacob Faibussowitsch static void f0_q(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[])
90d71ae5a4SJacob Faibussowitsch {
91c4762a1bSJed Brown PetscInt i;
92c4762a1bSJed Brown PetscScalar divu;
93c4762a1bSJed Brown
94c4762a1bSJed Brown divu = 0.;
95c4762a1bSJed Brown for (i = 0; i < dim; ++i) divu += u_x[uOff_x[0] + i * dim + i];
96c4762a1bSJed Brown f0[0] = u[uOff[1]] - divu;
97c4762a1bSJed Brown }
98c4762a1bSJed Brown
99c4762a1bSJed Brown /* Residual function for p_err = \div{u} - p. */
f0_w(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[])100d71ae5a4SJacob Faibussowitsch static void f0_w(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[])
101d71ae5a4SJacob Faibussowitsch {
102c4762a1bSJed Brown PetscInt i;
103c4762a1bSJed Brown PetscScalar divu;
104c4762a1bSJed Brown
105c4762a1bSJed Brown divu = 0.;
106c4762a1bSJed Brown for (i = 0; i < dim; ++i) divu += u_x[uOff_x[0] + i * dim + i];
107c4762a1bSJed Brown f0[0] = u[uOff[2]] - u[uOff[1]] + divu;
108c4762a1bSJed Brown }
109c4762a1bSJed Brown
110c4762a1bSJed Brown /* Boundary residual for the embedding system. Need one for each form of
111c4762a1bSJed Brown * solution. These enforce u = \hat{u} at the boundary. */
f0_bd_u_sinusoid(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[])112d71ae5a4SJacob Faibussowitsch static void f0_bd_u_sinusoid(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[])
113d71ae5a4SJacob Faibussowitsch {
114c4762a1bSJed Brown PetscInt d;
115c4762a1bSJed Brown PetscScalar *u_rhs;
116c4762a1bSJed Brown
1173ba16761SJacob Faibussowitsch PetscCallAbort(PETSC_COMM_SELF, PetscCalloc1(dim, &u_rhs));
1183ba16761SJacob Faibussowitsch PetscCallAbort(PETSC_COMM_SELF, sinusoid_u(dim, t, x, dim, u_rhs, NULL));
119c4762a1bSJed Brown for (d = 0; d < dim; ++d) f0[d] = u_rhs[d];
1203ba16761SJacob Faibussowitsch PetscCallAbort(PETSC_COMM_SELF, PetscFree(u_rhs));
121c4762a1bSJed Brown }
122c4762a1bSJed Brown
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[])123d71ae5a4SJacob Faibussowitsch 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[])
124d71ae5a4SJacob Faibussowitsch {
125c4762a1bSJed Brown PetscInt d;
126c4762a1bSJed Brown PetscScalar *u_rhs;
127c4762a1bSJed Brown
1283ba16761SJacob Faibussowitsch PetscCallAbort(PETSC_COMM_SELF, PetscCalloc1(dim, &u_rhs));
1293ba16761SJacob Faibussowitsch PetscCallAbort(PETSC_COMM_SELF, linear_u(dim, t, x, dim, u_rhs, NULL));
130c4762a1bSJed Brown for (d = 0; d < dim; ++d) f0[d] = u_rhs[d];
1313ba16761SJacob Faibussowitsch PetscCallAbort(PETSC_COMM_SELF, PetscFree(u_rhs));
132c4762a1bSJed Brown }
133c4762a1bSJed Brown /* Jacobian functions. For the following, v is the test function associated with
134c4762a1bSJed Brown * u, q the test function associated with p, and w the test function associated
135c4762a1bSJed Brown * with d. */
136c4762a1bSJed Brown /* <v, u> */
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[])137d71ae5a4SJacob Faibussowitsch 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[])
138d71ae5a4SJacob Faibussowitsch {
139c4762a1bSJed Brown PetscInt c;
140c4762a1bSJed Brown for (c = 0; c < dim; ++c) g0[c * dim + c] = 1.0;
141c4762a1bSJed Brown }
142c4762a1bSJed Brown
143c4762a1bSJed Brown /* <q, p> */
g0_qp(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[])144d71ae5a4SJacob Faibussowitsch static void g0_qp(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[])
145d71ae5a4SJacob Faibussowitsch {
146c4762a1bSJed Brown PetscInt d;
147c4762a1bSJed Brown for (d = 0; d < dim; ++d) g0[d * dim + d] = 1.0;
148c4762a1bSJed Brown }
149c4762a1bSJed Brown
150c4762a1bSJed Brown /* -<q, \div{u}> For the embedded system. This is different from the method of
151c4762a1bSJed Brown * manufactured solution because instead of computing <q,\div{u}> - <q,f> we
152c4762a1bSJed Brown * need <q,p> - <q,\div{u}.*/
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[])153d71ae5a4SJacob Faibussowitsch 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[])
154d71ae5a4SJacob Faibussowitsch {
155c4762a1bSJed Brown PetscInt d;
156c4762a1bSJed Brown for (d = 0; d < dim; ++d) g1[d * dim + d] = -1.0;
157c4762a1bSJed Brown }
158c4762a1bSJed Brown
159c4762a1bSJed Brown /* <w, p> This is only used by the embedded system. Where we need to compute
160c4762a1bSJed Brown * <w,d> - <w,p> + <w, \div{u}>*/
g0_wp(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[])161d71ae5a4SJacob Faibussowitsch static void g0_wp(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[])
162d71ae5a4SJacob Faibussowitsch {
163c4762a1bSJed Brown PetscInt d;
164c4762a1bSJed Brown for (d = 0; d < dim; ++d) g0[d * dim + d] = -1.0;
165c4762a1bSJed Brown }
166c4762a1bSJed Brown
167c4762a1bSJed Brown /* <w, d> */
g0_wd(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[])168d71ae5a4SJacob Faibussowitsch static void g0_wd(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[])
169d71ae5a4SJacob Faibussowitsch {
170c4762a1bSJed Brown PetscInt c;
171c4762a1bSJed Brown for (c = 0; c < dim; ++c) g0[c * dim + c] = 1.0;
172c4762a1bSJed Brown }
173c4762a1bSJed Brown
174c4762a1bSJed Brown /* <w, \div{u}> for the embedded system. */
g1_wu(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[])175d71ae5a4SJacob Faibussowitsch static void g1_wu(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[])
176d71ae5a4SJacob Faibussowitsch {
177c4762a1bSJed Brown PetscInt d;
178c4762a1bSJed Brown for (d = 0; d < dim; ++d) g1[d * dim + d] = 1.0;
179c4762a1bSJed Brown }
180c4762a1bSJed Brown
181c4762a1bSJed Brown /* Enum and string array for selecting mesh perturbation options */
1829371c9d4SSatish Balay typedef enum {
1839371c9d4SSatish Balay NONE = 0,
1849371c9d4SSatish Balay PERTURB = 1,
1859371c9d4SSatish Balay SKEW = 2,
1869371c9d4SSatish Balay SKEW_PERTURB = 3
1879371c9d4SSatish Balay } Transform;
188c4762a1bSJed Brown const char *const TransformTypes[] = {"none", "perturb", "skew", "skew_perturb", "Perturbation", "", NULL};
189c4762a1bSJed Brown
190c4762a1bSJed Brown /* Enum and string array for selecting the form of the exact solution*/
1919371c9d4SSatish Balay typedef enum {
1929371c9d4SSatish Balay LINEAR = 0,
1939371c9d4SSatish Balay SINUSOIDAL = 1
1949371c9d4SSatish Balay } Solution;
195c4762a1bSJed Brown const char *const SolutionTypes[] = {"linear", "sinusoidal", "Solution", "", NULL};
196c4762a1bSJed Brown
1979371c9d4SSatish Balay typedef struct {
198c4762a1bSJed Brown Transform mesh_transform;
199c4762a1bSJed Brown Solution sol_form;
200c4762a1bSJed Brown } UserCtx;
201c4762a1bSJed Brown
202c4762a1bSJed Brown /* Process command line options and initialize the UserCtx struct */
ProcessOptions(MPI_Comm comm,UserCtx * user)203d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, UserCtx *user)
204d71ae5a4SJacob Faibussowitsch {
205c4762a1bSJed Brown PetscFunctionBegin;
206c4762a1bSJed Brown /* Default to 2D, unperturbed triangle mesh and Linear solution.*/
207c4762a1bSJed Brown user->mesh_transform = NONE;
208c4762a1bSJed Brown user->sol_form = LINEAR;
209c4762a1bSJed Brown
210d0609cedSBarry Smith PetscOptionsBegin(comm, "", "H-div Test Options", "DMPLEX");
2119566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-mesh_transform", "Method used to perturb the mesh vertices. Options are skew, perturb, skew_perturb,or none", "ex39.c", TransformTypes, (PetscEnum)user->mesh_transform, (PetscEnum *)&user->mesh_transform, NULL));
2129566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-sol_form", "Form of the exact solution. Options are Linear or Sinusoidal", "ex39.c", SolutionTypes, (PetscEnum)user->sol_form, (PetscEnum *)&user->sol_form, NULL));
213d0609cedSBarry Smith PetscOptionsEnd();
2143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
215c4762a1bSJed Brown }
216c4762a1bSJed Brown
217c4762a1bSJed Brown /* Perturb the position of each mesh vertex by a small amount.*/
PerturbMesh(DM * mesh,PetscScalar * coordVals,PetscInt npoints,PetscInt dim)218d71ae5a4SJacob Faibussowitsch static PetscErrorCode PerturbMesh(DM *mesh, PetscScalar *coordVals, PetscInt npoints, PetscInt dim)
219d71ae5a4SJacob Faibussowitsch {
220c4762a1bSJed Brown PetscInt i, j, k;
221d092c84bSBrandon Whitchurch PetscReal minCoords[3], maxCoords[3], maxPert[3], randVal, amp;
222c4762a1bSJed Brown PetscRandom ran;
223c4762a1bSJed Brown
224c4762a1bSJed Brown PetscFunctionBegin;
2259566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(*mesh, &dim));
2269566063dSJacob Faibussowitsch PetscCall(DMGetLocalBoundingBox(*mesh, minCoords, maxCoords));
2279566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &ran));
228c4762a1bSJed Brown
229c4762a1bSJed Brown /* Compute something approximately equal to half an edge length. This is the
230a5b23f4aSJose E. Roman * most we can perturb points and guarantee that there won't be any topology
231c4762a1bSJed Brown * issues. */
232d092c84bSBrandon Whitchurch for (k = 0; k < dim; ++k) maxPert[k] = 0.025 * (maxCoords[k] - minCoords[k]) / (PetscPowReal(npoints, 1. / dim) - 1);
233c4762a1bSJed Brown /* For each mesh vertex */
234c4762a1bSJed Brown for (i = 0; i < npoints; ++i) {
235c4762a1bSJed Brown /* For each coordinate of the vertex */
236c4762a1bSJed Brown for (j = 0; j < dim; ++j) {
237c4762a1bSJed Brown /* Generate a random amplitude in [-0.5*maxPert, 0.5*maxPert] */
2389566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValueReal(ran, &randVal));
239c4762a1bSJed Brown amp = maxPert[j] * (randVal - 0.5);
240c4762a1bSJed Brown /* Add the perturbation to the vertex*/
241d092c84bSBrandon Whitchurch coordVals[dim * i + j] += amp;
242c4762a1bSJed Brown }
243c4762a1bSJed Brown }
244c4762a1bSJed Brown
2453ba16761SJacob Faibussowitsch PetscCall(PetscRandomDestroy(&ran));
2463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
247c4762a1bSJed Brown }
248c4762a1bSJed Brown
249c4762a1bSJed Brown /* Apply a global skew transformation to the mesh. */
SkewMesh(DM * mesh,PetscScalar * coordVals,PetscInt npoints,PetscInt dim)250d71ae5a4SJacob Faibussowitsch static PetscErrorCode SkewMesh(DM *mesh, PetscScalar *coordVals, PetscInt npoints, PetscInt dim)
251d71ae5a4SJacob Faibussowitsch {
252c4762a1bSJed Brown PetscInt i, j, k, l;
253c4762a1bSJed Brown PetscScalar *transMat;
254c4762a1bSJed Brown PetscScalar tmpcoord[3];
255c4762a1bSJed Brown PetscRandom ran;
256c4762a1bSJed Brown PetscReal randVal;
257c4762a1bSJed Brown
258c4762a1bSJed Brown PetscFunctionBegin;
2599566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(dim * dim, &transMat));
2609566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &ran));
261c4762a1bSJed Brown
262c4762a1bSJed Brown /* Make a matrix representing a skew transformation */
263c4762a1bSJed Brown for (i = 0; i < dim; ++i) {
264c4762a1bSJed Brown for (j = 0; j < dim; ++j) {
2659566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValueReal(ran, &randVal));
266d092c84bSBrandon Whitchurch if (i == j) transMat[i * dim + j] = 1.;
267c4762a1bSJed Brown else if (j < i) transMat[i * dim + j] = 2 * (j + i) * randVal;
268c4762a1bSJed Brown else transMat[i * dim + j] = 0;
269c4762a1bSJed Brown }
270c4762a1bSJed Brown }
271c4762a1bSJed Brown
272da81f932SPierre Jolivet /* Multiply each coordinate vector by our transformation.*/
273c4762a1bSJed Brown for (i = 0; i < npoints; ++i) {
274c4762a1bSJed Brown for (j = 0; j < dim; ++j) {
275c4762a1bSJed Brown tmpcoord[j] = 0;
276c4762a1bSJed Brown for (k = 0; k < dim; ++k) tmpcoord[j] += coordVals[dim * i + k] * transMat[dim * k + j];
277c4762a1bSJed Brown }
278c4762a1bSJed Brown for (l = 0; l < dim; ++l) coordVals[dim * i + l] = tmpcoord[l];
279c4762a1bSJed Brown }
2809566063dSJacob Faibussowitsch PetscCall(PetscFree(transMat));
2819566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&ran));
2823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
283c4762a1bSJed Brown }
284c4762a1bSJed Brown
285c4762a1bSJed Brown /* Accesses the mesh coordinate array and performs the transformation operations
286c4762a1bSJed Brown * specified by the user options */
TransformMesh(UserCtx * user,DM * mesh)287d71ae5a4SJacob Faibussowitsch static PetscErrorCode TransformMesh(UserCtx *user, DM *mesh)
288d71ae5a4SJacob Faibussowitsch {
289c4762a1bSJed Brown PetscInt dim, npoints;
290c4762a1bSJed Brown PetscScalar *coordVals;
291c4762a1bSJed Brown Vec coords;
292c4762a1bSJed Brown
293c4762a1bSJed Brown PetscFunctionBegin;
2949566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(*mesh, &coords));
2959566063dSJacob Faibussowitsch PetscCall(VecGetArray(coords, &coordVals));
2969566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(coords, &npoints));
2979566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(*mesh, &dim));
298c4762a1bSJed Brown npoints = npoints / dim;
299c4762a1bSJed Brown
300c4762a1bSJed Brown switch (user->mesh_transform) {
301d71ae5a4SJacob Faibussowitsch case PERTURB:
302d71ae5a4SJacob Faibussowitsch PetscCall(PerturbMesh(mesh, coordVals, npoints, dim));
303d71ae5a4SJacob Faibussowitsch break;
304d71ae5a4SJacob Faibussowitsch case SKEW:
305d71ae5a4SJacob Faibussowitsch PetscCall(SkewMesh(mesh, coordVals, npoints, dim));
306d71ae5a4SJacob Faibussowitsch break;
307c4762a1bSJed Brown case SKEW_PERTURB:
3089566063dSJacob Faibussowitsch PetscCall(SkewMesh(mesh, coordVals, npoints, dim));
3099566063dSJacob Faibussowitsch PetscCall(PerturbMesh(mesh, coordVals, npoints, dim));
310c4762a1bSJed Brown break;
311d71ae5a4SJacob Faibussowitsch default:
312d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "invalid mesh transformation");
313c4762a1bSJed Brown }
3149566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(coords, &coordVals));
3159566063dSJacob Faibussowitsch PetscCall(DMSetCoordinates(*mesh, coords));
3163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
317c4762a1bSJed Brown }
318c4762a1bSJed Brown
CreateMesh(MPI_Comm comm,UserCtx * user,DM * mesh)319d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, UserCtx *user, DM *mesh)
320d71ae5a4SJacob Faibussowitsch {
321c4762a1bSJed Brown PetscFunctionBegin;
3229566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, mesh));
3239566063dSJacob Faibussowitsch PetscCall(DMSetType(*mesh, DMPLEX));
3249566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*mesh));
325c4762a1bSJed Brown
326c4762a1bSJed Brown /* Perform any mesh transformations if specified by user */
32748a46eb9SPierre Jolivet if (user->mesh_transform != NONE) PetscCall(TransformMesh(user, mesh));
328c4762a1bSJed Brown
329c4762a1bSJed Brown /* Get any other mesh options from the command line */
3309566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(*mesh, user));
3319566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*mesh, NULL, "-dm_view"));
3323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
333c4762a1bSJed Brown }
334c4762a1bSJed Brown
335c4762a1bSJed Brown /* Setup the system of equations that we wish to solve */
SetupProblem(DM dm,UserCtx * user)336d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupProblem(DM dm, UserCtx *user)
337d71ae5a4SJacob Faibussowitsch {
338c4762a1bSJed Brown PetscDS prob;
33945480ffeSMatthew G. Knepley DMLabel label;
340c4762a1bSJed Brown const PetscInt id = 1;
341c4762a1bSJed Brown
342c4762a1bSJed Brown PetscFunctionBegin;
3439566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &prob));
344c4762a1bSJed Brown /* All of these are independent of the user's choice of solution */
3459566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(prob, 1, f0_q, NULL));
3469566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(prob, 2, f0_w, NULL));
3479566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(prob, 0, 0, g0_vu, NULL, NULL, NULL));
3489566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(prob, 1, 0, NULL, g1_qu, NULL, NULL));
3499566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(prob, 1, 1, g0_qp, NULL, NULL, NULL));
3509566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(prob, 2, 0, NULL, g1_wu, NULL, NULL));
3519566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(prob, 2, 1, g0_wp, NULL, NULL, NULL));
3529566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(prob, 2, 2, g0_wd, NULL, NULL, NULL));
353c4762a1bSJed Brown
354c4762a1bSJed Brown /* Field 2 is the error between \div{u} and pressure in a higher dimensional
355c4762a1bSJed Brown * space. If all is right this should be machine zero. */
3569566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(prob, 2, zero_func, NULL));
357c4762a1bSJed Brown
358c4762a1bSJed Brown switch (user->sol_form) {
359c4762a1bSJed Brown case LINEAR:
3609566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(prob, 0, f0_v_linear, NULL));
3619566063dSJacob Faibussowitsch PetscCall(PetscDSSetBdResidual(prob, 0, f0_bd_u_linear, NULL));
3629566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(prob, 0, linear_u, NULL));
3639566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(prob, 1, linear_p, NULL));
364c4762a1bSJed Brown break;
365c4762a1bSJed Brown case SINUSOIDAL:
3669566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(prob, 0, f0_v_sinusoid, NULL));
3679566063dSJacob Faibussowitsch PetscCall(PetscDSSetBdResidual(prob, 0, f0_bd_u_sinusoid, NULL));
3689566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(prob, 0, sinusoid_u, NULL));
3699566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(prob, 1, sinusoid_p, NULL));
370c4762a1bSJed Brown break;
371d71ae5a4SJacob Faibussowitsch default:
372d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "invalid solution form");
373c4762a1bSJed Brown }
374c4762a1bSJed Brown
3759566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label));
37657d50842SBarry Smith PetscCall(PetscDSAddBoundary(prob, DM_BC_NATURAL, "Boundary Integral", label, 1, &id, 0, 0, NULL, (PetscVoidFn *)NULL, NULL, user, NULL));
3773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
378c4762a1bSJed Brown }
379c4762a1bSJed Brown
380c4762a1bSJed Brown /* Create the finite element spaces we will use for this system */
SetupDiscretization(DM mesh,PetscErrorCode (* setup)(DM,UserCtx *),UserCtx * user)381d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupDiscretization(DM mesh, PetscErrorCode (*setup)(DM, UserCtx *), UserCtx *user)
382d71ae5a4SJacob Faibussowitsch {
383c4762a1bSJed Brown DM cdm = mesh;
384c4762a1bSJed Brown PetscFE fevel, fepres, fedivErr;
38530602db0SMatthew G. Knepley PetscInt dim;
38630602db0SMatthew G. Knepley PetscBool simplex;
387c4762a1bSJed Brown
388c4762a1bSJed Brown PetscFunctionBegin;
3899566063dSJacob Faibussowitsch PetscCall(DMGetDimension(mesh, &dim));
3909566063dSJacob Faibussowitsch PetscCall(DMPlexIsSimplex(mesh, &simplex));
391c4762a1bSJed Brown /* Create FE objects and give them names so that options can be set from
392c4762a1bSJed Brown * command line */
3939566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)mesh), dim, dim, simplex, "velocity_", -1, &fevel));
3949566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fevel, "velocity"));
395c4762a1bSJed Brown
3969566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)mesh), dim, 1, simplex, "pressure_", -1, &fepres));
3979566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fepres, "pressure"));
398c4762a1bSJed Brown
399d0609cedSBarry Smith PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)mesh), dim, 1, simplex, "divErr_", -1, &fedivErr));
4009566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fedivErr, "divErr"));
401c4762a1bSJed Brown
4029566063dSJacob Faibussowitsch PetscCall(PetscFECopyQuadrature(fevel, fepres));
4039566063dSJacob Faibussowitsch PetscCall(PetscFECopyQuadrature(fevel, fedivErr));
404c4762a1bSJed Brown
405c4762a1bSJed Brown /* Associate the FE objects with the mesh and setup the system */
4069566063dSJacob Faibussowitsch PetscCall(DMSetField(mesh, 0, NULL, (PetscObject)fevel));
4079566063dSJacob Faibussowitsch PetscCall(DMSetField(mesh, 1, NULL, (PetscObject)fepres));
4089566063dSJacob Faibussowitsch PetscCall(DMSetField(mesh, 2, NULL, (PetscObject)fedivErr));
4099566063dSJacob Faibussowitsch PetscCall(DMCreateDS(mesh));
4109566063dSJacob Faibussowitsch PetscCall((*setup)(mesh, user));
411c4762a1bSJed Brown
412c4762a1bSJed Brown while (cdm) {
4139566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(mesh, cdm));
4149566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(cdm, &cdm));
415c4762a1bSJed Brown }
416c4762a1bSJed Brown
417c4762a1bSJed Brown /* The Mesh now owns the fields, so we can destroy the FEs created in this
418c4762a1bSJed Brown * function */
4199566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fevel));
4209566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fepres));
4219566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fedivErr));
4229566063dSJacob Faibussowitsch PetscCall(DMDestroy(&cdm));
4233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
424c4762a1bSJed Brown }
425c4762a1bSJed Brown
main(int argc,char ** argv)426d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
427d71ae5a4SJacob Faibussowitsch {
428c4762a1bSJed Brown PetscInt i;
429c4762a1bSJed Brown UserCtx user;
430c4762a1bSJed Brown DM mesh;
431c4762a1bSJed Brown SNES snes;
432c4762a1bSJed Brown Vec computed, divErr;
433c4762a1bSJed Brown PetscReal divErrNorm;
434c4762a1bSJed Brown IS *fieldIS;
435c4762a1bSJed Brown PetscBool exampleSuccess = PETSC_FALSE;
436d092c84bSBrandon Whitchurch const PetscReal errTol = 10. * PETSC_SMALL;
437c4762a1bSJed Brown
438c4762a1bSJed Brown char stdFormat[] = "L2 Norm of the Divergence Error is: %g\n H(div) elements working correctly: %s\n";
439c4762a1bSJed Brown
440c4762a1bSJed Brown /* Initialize PETSc */
441327415f7SBarry Smith PetscFunctionBeginUser;
4429566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help));
4439566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
444c4762a1bSJed Brown
445c4762a1bSJed Brown /* Set up the system, we need to create a solver and a mesh and then assign
446c4762a1bSJed Brown * the correct spaces into the mesh*/
4479566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
4489566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &mesh));
4499566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snes, mesh));
4509566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(mesh, SetupProblem, &user));
4516493148fSStefano Zampini PetscCall(DMPlexSetSNESLocalFEM(mesh, PETSC_FALSE, &user));
4529566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes));
453c4762a1bSJed Brown
454c4762a1bSJed Brown /* Grab field IS so that we can view the solution by field */
4559566063dSJacob Faibussowitsch PetscCall(DMCreateFieldIS(mesh, NULL, NULL, &fieldIS));
456c4762a1bSJed Brown
457c4762a1bSJed Brown /* Create a vector to store the SNES solution, solve the system and grab the
458c4762a1bSJed Brown * solution from SNES */
4599566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(mesh, &computed));
4609566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)computed, "computedSolution"));
4619566063dSJacob Faibussowitsch PetscCall(VecSet(computed, 0.0));
4629566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, computed));
4639566063dSJacob Faibussowitsch PetscCall(SNESGetSolution(snes, &computed));
4649566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(computed, NULL, "-computedSolution_view"));
465c4762a1bSJed Brown
466c4762a1bSJed Brown /* Now we pull out the portion of the vector corresponding to the 3rd field
467c4762a1bSJed Brown * which is the error between \div{u} computed in a higher dimensional space
468c4762a1bSJed Brown * and p = \div{u} computed in a low dimension space. We report the L2 norm of
469c4762a1bSJed Brown * this vector which should be zero if the H(div) spaces are implemented
470c4762a1bSJed Brown * correctly. */
4719566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(computed, fieldIS[2], &divErr));
4729566063dSJacob Faibussowitsch PetscCall(VecNorm(divErr, NORM_2, &divErrNorm));
4739566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(computed, fieldIS[2], &divErr));
474c4762a1bSJed Brown exampleSuccess = (PetscBool)(divErrNorm <= errTol);
475c4762a1bSJed Brown
4769566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, stdFormat, divErrNorm, exampleSuccess ? "true" : "false"));
477c4762a1bSJed Brown
478c4762a1bSJed Brown /* Tear down */
4799566063dSJacob Faibussowitsch PetscCall(VecDestroy(&divErr));
4809566063dSJacob Faibussowitsch PetscCall(VecDestroy(&computed));
48148a46eb9SPierre Jolivet for (i = 0; i < 3; ++i) PetscCall(ISDestroy(&fieldIS[i]));
4829566063dSJacob Faibussowitsch PetscCall(PetscFree(fieldIS));
4839566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes));
4849566063dSJacob Faibussowitsch PetscCall(DMDestroy(&mesh));
4859566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
486b122ec5aSJacob Faibussowitsch return 0;
487c4762a1bSJed Brown }
488c4762a1bSJed Brown
489c4762a1bSJed Brown /*TEST
490c4762a1bSJed Brown testset:
491c4762a1bSJed Brown suffix: 2d_bdm
492c4762a1bSJed Brown requires: triangle
49330602db0SMatthew G. Knepley args: -velocity_petscfe_default_quadrature_order 1 \
494c4762a1bSJed Brown -velocity_petscspace_degree 1 \
495c4762a1bSJed Brown -velocity_petscdualspace_type bdm \
496c4762a1bSJed Brown -divErr_petscspace_degree 1 \
497c4762a1bSJed Brown -divErr_petscdualspace_lagrange_continuity false \
498c4762a1bSJed Brown -snes_error_if_not_converged \
499c4762a1bSJed Brown -ksp_rtol 1e-10 \
500c4762a1bSJed Brown -ksp_error_if_not_converged \
501c4762a1bSJed Brown -pc_type fieldsplit\
502c4762a1bSJed Brown -pc_fieldsplit_detect_saddle_point\
503c4762a1bSJed Brown -pc_fieldsplit_type schur\
504c4762a1bSJed Brown -pc_fieldsplit_schur_precondition full
505c4762a1bSJed Brown test:
506c4762a1bSJed Brown suffix: linear
507c4762a1bSJed Brown args: -sol_form linear -mesh_transform none
508c4762a1bSJed Brown test:
509c4762a1bSJed Brown suffix: sinusoidal
510c4762a1bSJed Brown args: -sol_form sinusoidal -mesh_transform none
511c4762a1bSJed Brown test:
512c4762a1bSJed Brown suffix: sinusoidal_skew
513c4762a1bSJed Brown args: -sol_form sinusoidal -mesh_transform skew
514c4762a1bSJed Brown test:
515c4762a1bSJed Brown suffix: sinusoidal_perturb
516c4762a1bSJed Brown args: -sol_form sinusoidal -mesh_transform perturb
517c4762a1bSJed Brown test:
518c4762a1bSJed Brown suffix: sinusoidal_skew_perturb
519c4762a1bSJed Brown args: -sol_form sinusoidal -mesh_transform skew_perturb
520c4762a1bSJed Brown
521c4762a1bSJed Brown testset:
522c4762a1bSJed Brown TODO: broken
523c4762a1bSJed Brown suffix: 2d_bdmq
5243886731fSPierre Jolivet output_file: output/empty.out
52530602db0SMatthew G. Knepley args: -dm_plex_simplex false \
526c4762a1bSJed Brown -velocity_petscspace_degree 1 \
527c4762a1bSJed Brown -velocity_petscdualspace_type bdm \
528d092c84bSBrandon Whitchurch -velocity_petscdualspace_lagrange_tensor 1 \
529c4762a1bSJed Brown -divErr_petscspace_degree 1 \
530c4762a1bSJed Brown -divErr_petscdualspace_lagrange_continuity false \
531c4762a1bSJed Brown -snes_error_if_not_converged \
532c4762a1bSJed Brown -ksp_rtol 1e-10 \
533c4762a1bSJed Brown -ksp_error_if_not_converged \
534c4762a1bSJed Brown -pc_type fieldsplit\
535c4762a1bSJed Brown -pc_fieldsplit_detect_saddle_point\
536c4762a1bSJed Brown -pc_fieldsplit_type schur\
537c4762a1bSJed Brown -pc_fieldsplit_schur_precondition full
538c4762a1bSJed Brown test:
539c4762a1bSJed Brown suffix: linear
540c4762a1bSJed Brown args: -sol_form linear -mesh_transform none
541c4762a1bSJed Brown test:
542c4762a1bSJed Brown suffix: sinusoidal
543c4762a1bSJed Brown args: -sol_form sinusoidal -mesh_transform none
544c4762a1bSJed Brown test:
545c4762a1bSJed Brown suffix: sinusoidal_skew
546c4762a1bSJed Brown args: -sol_form sinusoidal -mesh_transform skew
547c4762a1bSJed Brown test:
548c4762a1bSJed Brown suffix: sinusoidal_perturb
549c4762a1bSJed Brown args: -sol_form sinusoidal -mesh_transform perturb
550c4762a1bSJed Brown test:
551c4762a1bSJed Brown suffix: sinusoidal_skew_perturb
552c4762a1bSJed Brown args: -sol_form sinusoidal -mesh_transform skew_perturb
553c4762a1bSJed Brown
554c4762a1bSJed Brown testset:
555c4762a1bSJed Brown suffix: 3d_bdm
556c4762a1bSJed Brown requires: ctetgen
55730602db0SMatthew G. Knepley args: -dm_plex_dim 3 \
558c4762a1bSJed Brown -velocity_petscspace_degree 1 \
559c4762a1bSJed Brown -velocity_petscdualspace_type bdm \
560c4762a1bSJed Brown -divErr_petscspace_degree 1 \
561c4762a1bSJed Brown -divErr_petscdualspace_lagrange_continuity false \
562c4762a1bSJed Brown -snes_error_if_not_converged \
563c4762a1bSJed Brown -ksp_rtol 1e-10 \
564c4762a1bSJed Brown -ksp_error_if_not_converged \
565c4762a1bSJed Brown -pc_type fieldsplit \
566c4762a1bSJed Brown -pc_fieldsplit_detect_saddle_point \
567c4762a1bSJed Brown -pc_fieldsplit_type schur \
568c4762a1bSJed Brown -pc_fieldsplit_schur_precondition full
569c4762a1bSJed Brown test:
570c4762a1bSJed Brown suffix: linear
571c4762a1bSJed Brown args: -sol_form linear -mesh_transform none
572c4762a1bSJed Brown test:
573c4762a1bSJed Brown suffix: sinusoidal
574c4762a1bSJed Brown args: -sol_form sinusoidal -mesh_transform none
575c4762a1bSJed Brown test:
576c4762a1bSJed Brown suffix: sinusoidal_skew
577c4762a1bSJed Brown args: -sol_form sinusoidal -mesh_transform skew
578c4762a1bSJed Brown test:
579c4762a1bSJed Brown suffix: sinusoidal_perturb
580c4762a1bSJed Brown args: -sol_form sinusoidal -mesh_transform perturb
581c4762a1bSJed Brown test:
582c4762a1bSJed Brown suffix: sinusoidal_skew_perturb
583c4762a1bSJed Brown args: -sol_form sinusoidal -mesh_transform skew_perturb
584c4762a1bSJed Brown
585c4762a1bSJed Brown testset:
586c4762a1bSJed Brown TODO: broken
587c4762a1bSJed Brown suffix: 3d_bdmq
5883886731fSPierre Jolivet output_file: output/empty.out
589c4762a1bSJed Brown requires: ctetgen
59030602db0SMatthew G. Knepley args: -dm_plex_dim 3 \
59130602db0SMatthew G. Knepley -dm_plex_simplex false \
592c4762a1bSJed Brown -velocity_petscspace_degree 1 \
593c4762a1bSJed Brown -velocity_petscdualspace_type bdm \
594d092c84bSBrandon Whitchurch -velocity_petscdualspace_lagrange_tensor 1 \
595c4762a1bSJed Brown -divErr_petscspace_degree 1 \
596c4762a1bSJed Brown -divErr_petscdualspace_lagrange_continuity false \
597c4762a1bSJed Brown -snes_error_if_not_converged \
598c4762a1bSJed Brown -ksp_rtol 1e-10 \
599c4762a1bSJed Brown -ksp_error_if_not_converged \
600c4762a1bSJed Brown -pc_type fieldsplit \
601c4762a1bSJed Brown -pc_fieldsplit_detect_saddle_point \
602c4762a1bSJed Brown -pc_fieldsplit_type schur \
603c4762a1bSJed Brown -pc_fieldsplit_schur_precondition full
604c4762a1bSJed Brown test:
605c4762a1bSJed Brown suffix: linear
606c4762a1bSJed Brown args: -sol_form linear -mesh_transform none
607c4762a1bSJed Brown test:
608c4762a1bSJed Brown suffix: sinusoidal
609c4762a1bSJed Brown args: -sol_form sinusoidal -mesh_transform none
610c4762a1bSJed Brown test:
611c4762a1bSJed Brown suffix: sinusoidal_skew
612c4762a1bSJed Brown args: -sol_form sinusoidal -mesh_transform skew
613c4762a1bSJed Brown test:
614c4762a1bSJed Brown suffix: sinusoidal_perturb
615c4762a1bSJed Brown args: -sol_form sinusoidal -mesh_transform perturb
616c4762a1bSJed Brown test:
617c4762a1bSJed Brown suffix: sinusoidal_skew_perturb
618c4762a1bSJed Brown args: -sol_form sinusoidal -mesh_transform skew_perturb
619d092c84bSBrandon Whitchurch
620d092c84bSBrandon Whitchurch test:
621d092c84bSBrandon Whitchurch suffix: quad_rt_0
62230602db0SMatthew G. Knepley args: -dm_plex_simplex false -mesh_transform skew \
623d092c84bSBrandon Whitchurch -divErr_petscspace_degree 1 \
624d092c84bSBrandon Whitchurch -divErr_petscdualspace_lagrange_continuity false \
625d092c84bSBrandon Whitchurch -snes_error_if_not_converged \
626d092c84bSBrandon Whitchurch -ksp_rtol 1e-10 \
627d092c84bSBrandon Whitchurch -ksp_error_if_not_converged \
628d092c84bSBrandon Whitchurch -pc_type fieldsplit\
629d092c84bSBrandon Whitchurch -pc_fieldsplit_detect_saddle_point\
630d092c84bSBrandon Whitchurch -pc_fieldsplit_type schur\
631d092c84bSBrandon Whitchurch -pc_fieldsplit_schur_precondition full \
632d092c84bSBrandon Whitchurch -velocity_petscfe_default_quadrature_order 1 \
633d092c84bSBrandon Whitchurch -velocity_petscspace_type sum \
634d092c84bSBrandon Whitchurch -velocity_petscspace_variables 2 \
635d092c84bSBrandon Whitchurch -velocity_petscspace_components 2 \
636d092c84bSBrandon Whitchurch -velocity_petscspace_sum_spaces 2 \
637d092c84bSBrandon Whitchurch -velocity_petscspace_sum_concatenate true \
638417c287bSToby Isaac -velocity_sumcomp_0_petscspace_variables 2 \
639417c287bSToby Isaac -velocity_sumcomp_0_petscspace_type tensor \
640417c287bSToby Isaac -velocity_sumcomp_0_petscspace_tensor_spaces 2 \
641417c287bSToby Isaac -velocity_sumcomp_0_petscspace_tensor_uniform false \
642417c287bSToby Isaac -velocity_sumcomp_0_tensorcomp_0_petscspace_degree 1 \
643417c287bSToby Isaac -velocity_sumcomp_0_tensorcomp_1_petscspace_degree 0 \
644417c287bSToby Isaac -velocity_sumcomp_1_petscspace_variables 2 \
645417c287bSToby Isaac -velocity_sumcomp_1_petscspace_type tensor \
646417c287bSToby Isaac -velocity_sumcomp_1_petscspace_tensor_spaces 2 \
647417c287bSToby Isaac -velocity_sumcomp_1_petscspace_tensor_uniform false \
648417c287bSToby Isaac -velocity_sumcomp_1_tensorcomp_0_petscspace_degree 0 \
649417c287bSToby Isaac -velocity_sumcomp_1_tensorcomp_1_petscspace_degree 1 \
650d092c84bSBrandon Whitchurch -velocity_petscdualspace_form_degree -1 \
651d092c84bSBrandon Whitchurch -velocity_petscdualspace_order 1 \
652d092c84bSBrandon Whitchurch -velocity_petscdualspace_lagrange_trimmed true
653c4762a1bSJed Brown TEST*/
654