xref: /petsc/src/dm/dt/tests/ex10.c (revision 8fb5bd83c3955fefcf33a54e3bb66920a9fa884b)
1d092c84bSBrandon Whitchurch static char help[] =
2d092c84bSBrandon Whitchurch   "Tests implementation of PetscSpace_Sum by solving the Poisson equations using a PetscSpace_Poly and a PetscSpace_Sum and checking that \
3d092c84bSBrandon Whitchurch   solutions agree up to machine precision.\n\n";
4d092c84bSBrandon Whitchurch 
5d092c84bSBrandon Whitchurch #include <petscdmplex.h>
6d092c84bSBrandon Whitchurch #include <petscds.h>
7d092c84bSBrandon Whitchurch #include <petscfe.h>
8d092c84bSBrandon Whitchurch #include <petscsnes.h>
9d092c84bSBrandon Whitchurch /* We are solving the system of equations:
10d092c84bSBrandon Whitchurch  * \vec{u} = -\grad{p}
11d092c84bSBrandon Whitchurch  * \div{u} = f
12d092c84bSBrandon Whitchurch  */
13d092c84bSBrandon Whitchurch 
14d092c84bSBrandon Whitchurch /* Exact solutions for linear velocity
15d092c84bSBrandon Whitchurch    \vec{u} = \vec{x};
16d092c84bSBrandon Whitchurch    p = -0.5*(\vec{x} \cdot \vec{x});
17d092c84bSBrandon Whitchurch    */
18d092c84bSBrandon Whitchurch static PetscErrorCode linear_u(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar *u,void *ctx)
19d092c84bSBrandon Whitchurch {
20d092c84bSBrandon Whitchurch   PetscInt c;
21d092c84bSBrandon Whitchurch 
22d092c84bSBrandon Whitchurch   for (c = 0; c < Nc; ++c) u[c] = x[c];
23d092c84bSBrandon Whitchurch   return 0;
24d092c84bSBrandon Whitchurch }
25d092c84bSBrandon Whitchurch 
26d092c84bSBrandon Whitchurch static PetscErrorCode linear_p(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar *u,void *ctx)
27d092c84bSBrandon Whitchurch {
28d092c84bSBrandon Whitchurch   PetscInt d;
29d092c84bSBrandon Whitchurch 
30d092c84bSBrandon Whitchurch   u[0] = 0.;
31d092c84bSBrandon Whitchurch   for (d=0; d<dim; ++d) u[0] += -0.5*x[d]*x[d];
32d092c84bSBrandon Whitchurch   return 0;
33d092c84bSBrandon Whitchurch }
34d092c84bSBrandon Whitchurch 
35d092c84bSBrandon Whitchurch static PetscErrorCode linear_divu(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar *u,void *ctx)
36d092c84bSBrandon Whitchurch {
37d092c84bSBrandon Whitchurch   u[0] = dim;
38d092c84bSBrandon Whitchurch   return 0;
39d092c84bSBrandon Whitchurch }
40d092c84bSBrandon Whitchurch 
41d092c84bSBrandon Whitchurch /* fx_v are the residual functions for the equation \vec{u} = \grad{p}. f0_v is the term <v,u>.*/
42d092c84bSBrandon Whitchurch 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[])
43d092c84bSBrandon Whitchurch {
44d092c84bSBrandon Whitchurch   PetscInt i;
45d092c84bSBrandon Whitchurch 
46d092c84bSBrandon Whitchurch   for (i=0; i<dim; ++i) f0[i] = u[uOff[0] + i];
47d092c84bSBrandon Whitchurch }
48d092c84bSBrandon Whitchurch 
49d092c84bSBrandon Whitchurch /* f1_v is the term <v,-\grad{p}> but we integrate by parts to get <\grad{v}, -p*I> */
50d092c84bSBrandon Whitchurch 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[])
51d092c84bSBrandon Whitchurch {
52d092c84bSBrandon Whitchurch   PetscInt c;
53d092c84bSBrandon Whitchurch 
54d092c84bSBrandon Whitchurch   for (c=0; c<dim; ++c) {
55d092c84bSBrandon Whitchurch     PetscInt d;
56d092c84bSBrandon Whitchurch 
57d092c84bSBrandon Whitchurch     for (d=0; d<dim; ++d) f1[c*dim + d] = (c==d) ? -u[uOff[1]] : 0;
58d092c84bSBrandon Whitchurch   }
59d092c84bSBrandon Whitchurch }
60d092c84bSBrandon Whitchurch 
61d092c84bSBrandon Whitchurch /* Residual function for enforcing \div{u} = f. */
62d092c84bSBrandon Whitchurch 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[])
63d092c84bSBrandon Whitchurch {
64d092c84bSBrandon Whitchurch   PetscScalar rhs,divu=0;
65d092c84bSBrandon Whitchurch   PetscInt    i;
66d092c84bSBrandon Whitchurch 
672da392ccSBarry Smith   (void)linear_divu(dim,t,x,dim,&rhs,NULL);
68d092c84bSBrandon Whitchurch   for (i=0; i< dim; ++i) divu += u_x[uOff_x[0]+i*dim+i];
69d092c84bSBrandon Whitchurch   f0[0] = divu-rhs;
70d092c84bSBrandon Whitchurch }
71d092c84bSBrandon Whitchurch 
72d092c84bSBrandon Whitchurch /* Boundary residual. Dirichlet boundary for u means u_bdy=p*n */
73d092c84bSBrandon Whitchurch 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[])
74d092c84bSBrandon Whitchurch {
75d092c84bSBrandon Whitchurch   PetscScalar pressure;
76d092c84bSBrandon Whitchurch   PetscInt    d;
77d092c84bSBrandon Whitchurch 
78d092c84bSBrandon Whitchurch   (void)linear_p(dim,t,x,dim,&pressure,NULL);
79d092c84bSBrandon Whitchurch   for (d=0; d<dim; ++d) f0[d] = pressure*n[d];
80d092c84bSBrandon Whitchurch }
81d092c84bSBrandon Whitchurch 
82d092c84bSBrandon Whitchurch /* gx_yz are the jacobian functions obtained by taking the derivative of the y residual w.r.t z*/
83d092c84bSBrandon Whitchurch 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[])
84d092c84bSBrandon Whitchurch {
85d092c84bSBrandon Whitchurch   PetscInt c;
86d092c84bSBrandon Whitchurch 
87d092c84bSBrandon Whitchurch   for (c=0; c<dim; ++c) g0[c*dim + c] = 1.0;
88d092c84bSBrandon Whitchurch }
89d092c84bSBrandon Whitchurch 
90d092c84bSBrandon Whitchurch 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[])
91d092c84bSBrandon Whitchurch {
92d092c84bSBrandon Whitchurch   PetscInt c;
93d092c84bSBrandon Whitchurch 
94d092c84bSBrandon Whitchurch   for (c=0; c<dim; ++c) g1[c*dim + c] = 1.0;
95d092c84bSBrandon Whitchurch }
96d092c84bSBrandon Whitchurch 
97d092c84bSBrandon Whitchurch 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[])
98d092c84bSBrandon Whitchurch {
99d092c84bSBrandon Whitchurch   PetscInt c;
100d092c84bSBrandon Whitchurch 
101d092c84bSBrandon Whitchurch   for (c=0; c<dim; ++c) g2[c*dim + c] = -1.0;
102d092c84bSBrandon Whitchurch }
103d092c84bSBrandon Whitchurch 
104d092c84bSBrandon Whitchurch typedef struct
105d092c84bSBrandon Whitchurch {
10630602db0SMatthew G. Knepley   PetscInt dummy;
107d092c84bSBrandon Whitchurch } UserCtx;
108d092c84bSBrandon Whitchurch 
109d092c84bSBrandon Whitchurch static PetscErrorCode CreateMesh(MPI_Comm comm,UserCtx *user,DM *mesh)
110d092c84bSBrandon Whitchurch {
111d092c84bSBrandon Whitchurch   PetscFunctionBegin;
1129566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, mesh));
1139566063dSJacob Faibussowitsch   PetscCall(DMSetType(*mesh, DMPLEX));
1149566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*mesh));
1159566063dSJacob Faibussowitsch   PetscCall(DMSetApplicationContext(*mesh,user));
1169566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*mesh,NULL,"-dm_view"));
117d092c84bSBrandon Whitchurch   PetscFunctionReturn(0);
118d092c84bSBrandon Whitchurch }
119d092c84bSBrandon Whitchurch 
120d092c84bSBrandon Whitchurch /* Setup the system of equations that we wish to solve */
121d092c84bSBrandon Whitchurch static PetscErrorCode SetupProblem(DM dm,UserCtx *user)
122d092c84bSBrandon Whitchurch {
12345480ffeSMatthew G. Knepley   PetscDS        ds;
12445480ffeSMatthew G. Knepley   DMLabel        label;
12545480ffeSMatthew G. Knepley   PetscWeakForm  wf;
126d092c84bSBrandon Whitchurch   const PetscInt id = 1;
12745480ffeSMatthew G. Knepley   PetscInt       bd;
128d092c84bSBrandon Whitchurch 
129d092c84bSBrandon Whitchurch   PetscFunctionBegin;
1309566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &ds));
131d092c84bSBrandon Whitchurch   /* All of these are independent of the user's choice of solution */
1329566063dSJacob Faibussowitsch   PetscCall(PetscDSSetResidual(ds,0,f0_v,f1_v));
1339566063dSJacob Faibussowitsch   PetscCall(PetscDSSetResidual(ds,1,f0_q_linear,NULL));
1349566063dSJacob Faibussowitsch   PetscCall(PetscDSSetJacobian(ds,0,0,g0_vu,NULL,NULL,NULL));
1359566063dSJacob Faibussowitsch   PetscCall(PetscDSSetJacobian(ds,0,1,NULL,NULL,g2_vp,NULL));
1369566063dSJacob Faibussowitsch   PetscCall(PetscDSSetJacobian(ds,1,0,NULL,g1_qu,NULL,NULL));
137d092c84bSBrandon Whitchurch 
1389566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "marker", &label));
1399566063dSJacob Faibussowitsch   PetscCall(PetscDSAddBoundary(ds,DM_BC_NATURAL,"Boundary Integral",label,1,&id,0,0,NULL,(void (*)(void))NULL,NULL,user,&bd));
1409566063dSJacob Faibussowitsch   PetscCall(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
1419566063dSJacob Faibussowitsch   PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, 1, 0, 0, 0, f0_bd_u_linear, 0, NULL));
14245480ffeSMatthew G. Knepley 
1439566063dSJacob Faibussowitsch   PetscCall(PetscDSSetExactSolution(ds,0,linear_u,NULL));
144*8fb5bd83SMatthew G. Knepley   PetscCall(PetscDSSetExactSolution(ds,1,linear_p,NULL));
145d092c84bSBrandon Whitchurch   PetscFunctionReturn(0);
146d092c84bSBrandon Whitchurch }
147d092c84bSBrandon Whitchurch 
148d092c84bSBrandon Whitchurch /* Create the finite element spaces we will use for this system */
149d092c84bSBrandon Whitchurch static PetscErrorCode SetupDiscretization(DM mesh,DM mesh_sum,PetscErrorCode (*setup)(DM,UserCtx*),UserCtx *user)
150d092c84bSBrandon Whitchurch {
151d092c84bSBrandon Whitchurch   DM             cdm = mesh,cdm_sum = mesh_sum;
152d092c84bSBrandon Whitchurch   PetscFE        u,divu,u_sum,divu_sum;
15330602db0SMatthew G. Knepley   PetscInt       dim;
15430602db0SMatthew G. Knepley   PetscBool      simplex;
155d092c84bSBrandon Whitchurch 
156d092c84bSBrandon Whitchurch   PetscFunctionBegin;
1579566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(mesh, &dim));
1589566063dSJacob Faibussowitsch   PetscCall(DMPlexIsSimplex(mesh, &simplex));
159d092c84bSBrandon Whitchurch   /* Create FE objects and give them names so that options can be set from
160d092c84bSBrandon Whitchurch    * command line. Each field gets 2 instances (i.e. velocity and velocity_sum)created twice so that we can compare between approaches. */
1619566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)mesh),dim,dim,simplex,"velocity_",-1,&u));
1629566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)u,"velocity"));
1639566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)mesh_sum),dim,dim,simplex,"velocity_sum_",-1,&u_sum));
1649566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)u_sum,"velocity_sum"));
1659566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)mesh),dim,1,simplex,"divu_",-1,&divu));
1669566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)divu,"divu"));
1679566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)mesh_sum),dim,1,simplex,"divu_sum_",-1,&divu_sum));
1689566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)divu_sum,"divu_sum"));
169d092c84bSBrandon Whitchurch 
1709566063dSJacob Faibussowitsch   PetscCall(PetscFECopyQuadrature(u,divu));
1719566063dSJacob Faibussowitsch   PetscCall(PetscFECopyQuadrature(u_sum,divu_sum));
172d092c84bSBrandon Whitchurch 
173d092c84bSBrandon Whitchurch   /* Associate the FE objects with the mesh and setup the system */
1749566063dSJacob Faibussowitsch   PetscCall(DMSetField(mesh,0,NULL,(PetscObject)u));
1759566063dSJacob Faibussowitsch   PetscCall(DMSetField(mesh,1,NULL,(PetscObject)divu));
1769566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(mesh));
1779566063dSJacob Faibussowitsch   PetscCall((*setup)(mesh,user));
178d092c84bSBrandon Whitchurch 
1799566063dSJacob Faibussowitsch   PetscCall(DMSetField(mesh_sum,0,NULL,(PetscObject)u_sum));
1809566063dSJacob Faibussowitsch   PetscCall(DMSetField(mesh_sum,1,NULL,(PetscObject)divu_sum));
1819566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(mesh_sum));
1829566063dSJacob Faibussowitsch   PetscCall((*setup)(mesh_sum,user));
183d092c84bSBrandon Whitchurch 
184d092c84bSBrandon Whitchurch   while (cdm) {
1859566063dSJacob Faibussowitsch     PetscCall(DMCopyDisc(mesh,cdm));
1869566063dSJacob Faibussowitsch     PetscCall(DMGetCoarseDM(cdm,&cdm));
187d092c84bSBrandon Whitchurch   }
188d092c84bSBrandon Whitchurch 
189d092c84bSBrandon Whitchurch   while (cdm_sum) {
1909566063dSJacob Faibussowitsch     PetscCall(DMCopyDisc(mesh_sum,cdm_sum));
1919566063dSJacob Faibussowitsch     PetscCall(DMGetCoarseDM(cdm_sum,&cdm_sum));
192d092c84bSBrandon Whitchurch   }
193d092c84bSBrandon Whitchurch 
194d092c84bSBrandon Whitchurch   /* The Mesh now owns the fields, so we can destroy the FEs created in this
195d092c84bSBrandon Whitchurch    * function */
1969566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&u));
1979566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&divu));
1989566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&u_sum));
1999566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&divu_sum));
2009566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&cdm));
2019566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&cdm_sum));
202d092c84bSBrandon Whitchurch   PetscFunctionReturn(0);
203d092c84bSBrandon Whitchurch }
204d092c84bSBrandon Whitchurch 
205d092c84bSBrandon Whitchurch int main(int argc,char **argv)
206d092c84bSBrandon Whitchurch {
207d092c84bSBrandon Whitchurch   UserCtx         user;
208d092c84bSBrandon Whitchurch   DM              dm,dm_sum;
209d092c84bSBrandon Whitchurch   SNES            snes,snes_sum;
210d092c84bSBrandon Whitchurch   Vec             u,u_sum;
211d092c84bSBrandon Whitchurch   PetscReal       errNorm;
212d092c84bSBrandon Whitchurch   const PetscReal errTol = PETSC_SMALL;
213d092c84bSBrandon Whitchurch 
2149566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
215d092c84bSBrandon Whitchurch 
216d092c84bSBrandon Whitchurch   /* Set up a snes for the standard approach, one space with 2 components */
2179566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_WORLD,&snes));
2189566063dSJacob Faibussowitsch   PetscCall(CreateMesh(PETSC_COMM_WORLD,&user,&dm));
2199566063dSJacob Faibussowitsch   PetscCall(SNESSetDM(snes,dm));
220d092c84bSBrandon Whitchurch 
221d092c84bSBrandon Whitchurch   /* Set up a snes for the sum space approach, where each subspace of the sum space represents one component */
2229566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_WORLD,&snes_sum));
2239566063dSJacob Faibussowitsch   PetscCall(CreateMesh(PETSC_COMM_WORLD,&user,&dm_sum));
2249566063dSJacob Faibussowitsch   PetscCall(SNESSetDM(snes_sum,dm_sum));
2259566063dSJacob Faibussowitsch   PetscCall(SetupDiscretization(dm,dm_sum,SetupProblem,&user));
226d092c84bSBrandon Whitchurch 
227d092c84bSBrandon Whitchurch   /* Set up and solve the system using standard approach. */
2289566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(dm,&u));
2299566063dSJacob Faibussowitsch   PetscCall(VecSet(u,0.0));
2309566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)u,"solution"));
2319566063dSJacob Faibussowitsch   PetscCall(DMPlexSetSNESLocalFEM(dm,&user,&user,&user));
2329566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
2339566063dSJacob Faibussowitsch   PetscCall(DMSNESCheckFromOptions(snes,u));
2349566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes,NULL,u));
2359566063dSJacob Faibussowitsch   PetscCall(SNESGetSolution(snes,&u));
2369566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(u,NULL,"-solution_view"));
237d092c84bSBrandon Whitchurch 
238d092c84bSBrandon Whitchurch   /* Set up and solve the sum space system */
2399566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(dm_sum,&u_sum));
2409566063dSJacob Faibussowitsch   PetscCall(VecSet(u_sum,0.0));
2419566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)u_sum,"solution_sum"));
2429566063dSJacob Faibussowitsch   PetscCall(DMPlexSetSNESLocalFEM(dm_sum,&user,&user,&user));
2439566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes_sum));
2449566063dSJacob Faibussowitsch   PetscCall(DMSNESCheckFromOptions(snes_sum,u_sum));
2459566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes_sum,NULL,u_sum));
2469566063dSJacob Faibussowitsch   PetscCall(SNESGetSolution(snes_sum,&u_sum));
2479566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(u_sum,NULL,"-solution_sum_view"));
248d092c84bSBrandon Whitchurch 
249d092c84bSBrandon Whitchurch   /* Check if standard solution and sum space solution match to machine precision */
2509566063dSJacob Faibussowitsch   PetscCall(VecAXPY(u_sum,-1,u));
2519566063dSJacob Faibussowitsch   PetscCall(VecNorm(u_sum,NORM_2,&errNorm));
252d0609cedSBarry Smith   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Sum space provides the same solution as a regular space: %s",(errNorm < errTol) ? "true" : "false"));
253d092c84bSBrandon Whitchurch 
254d092c84bSBrandon Whitchurch   /* Cleanup */
2559566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u_sum));
2569566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
2579566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snes_sum));
2589566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snes));
2599566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm_sum));
2609566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
2619566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
262b122ec5aSJacob Faibussowitsch   return 0;
263d092c84bSBrandon Whitchurch }
264d092c84bSBrandon Whitchurch 
265d092c84bSBrandon Whitchurch /*TEST
266d092c84bSBrandon Whitchurch   test:
267d092c84bSBrandon Whitchurch     suffix: 2d_lagrange
268d092c84bSBrandon Whitchurch     requires: triangle
26930602db0SMatthew G. Knepley     args: -velocity_petscspace_degree 1 \
270d092c84bSBrandon Whitchurch       -velocity_petscspace_type poly \
271d092c84bSBrandon Whitchurch       -velocity_petscspace_components 2\
272d092c84bSBrandon Whitchurch       -velocity_petscdualspace_type lagrange \
273d092c84bSBrandon Whitchurch       -divu_petscspace_degree 0 \
274d092c84bSBrandon Whitchurch       -divu_petscspace_type poly \
275d092c84bSBrandon Whitchurch       -divu_petscdualspace_lagrange_continuity false \
276d092c84bSBrandon Whitchurch       -velocity_sum_petscfe_default_quadrature_order 1 \
277d092c84bSBrandon Whitchurch       -velocity_sum_petscspace_degree 1 \
278d092c84bSBrandon Whitchurch       -velocity_sum_petscspace_type sum \
279d092c84bSBrandon Whitchurch       -velocity_sum_petscspace_variables 2 \
280d092c84bSBrandon Whitchurch       -velocity_sum_petscspace_components 2 \
281d092c84bSBrandon Whitchurch       -velocity_sum_petscspace_sum_spaces 2 \
282d092c84bSBrandon Whitchurch       -velocity_sum_petscspace_sum_concatenate true \
283d092c84bSBrandon Whitchurch       -velocity_sum_petscdualspace_type lagrange \
284417c287bSToby Isaac       -velocity_sum_sumcomp_0_petscspace_type poly \
285417c287bSToby Isaac       -velocity_sum_sumcomp_0_petscspace_degree 1 \
286417c287bSToby Isaac       -velocity_sum_sumcomp_0_petscspace_variables 2 \
287417c287bSToby Isaac       -velocity_sum_sumcomp_0_petscspace_components 1 \
288417c287bSToby Isaac       -velocity_sum_sumcomp_1_petscspace_type poly \
289417c287bSToby Isaac       -velocity_sum_sumcomp_1_petscspace_degree 1 \
290417c287bSToby Isaac       -velocity_sum_sumcomp_1_petscspace_variables 2 \
291417c287bSToby Isaac       -velocity_sum_sumcomp_1_petscspace_components 1 \
292d092c84bSBrandon Whitchurch       -divu_sum_petscspace_degree 0 \
293d092c84bSBrandon Whitchurch       -divu_sum_petscspace_type sum \
294d092c84bSBrandon Whitchurch       -divu_sum_petscspace_variables 2 \
295d092c84bSBrandon Whitchurch       -divu_sum_petscspace_components 1 \
296d092c84bSBrandon Whitchurch       -divu_sum_petscspace_sum_spaces 1 \
297d092c84bSBrandon Whitchurch       -divu_sum_petscspace_sum_concatenate true \
298d092c84bSBrandon Whitchurch       -divu_sum_petscdualspace_lagrange_continuity false \
299417c287bSToby Isaac       -divu_sum_sumcomp_0_petscspace_type poly \
300417c287bSToby Isaac       -divu_sum_sumcomp_0_petscspace_degree 0 \
301417c287bSToby Isaac       -divu_sum_sumcomp_0_petscspace_variables 2 \
302417c287bSToby Isaac       -divu_sum_sumcomp_0_petscspace_components 1 \
303d092c84bSBrandon Whitchurch       -dm_refine 0 \
304d092c84bSBrandon Whitchurch       -snes_error_if_not_converged \
305d092c84bSBrandon Whitchurch       -ksp_rtol 1e-10 \
306d092c84bSBrandon Whitchurch       -ksp_error_if_not_converged \
307d092c84bSBrandon Whitchurch       -pc_type fieldsplit\
308d092c84bSBrandon Whitchurch       -pc_fieldsplit_type schur\
309d092c84bSBrandon Whitchurch       -divu_sum_petscdualspace_lagrange_continuity false \
310d092c84bSBrandon Whitchurch       -pc_fieldsplit_schur_precondition full
311d092c84bSBrandon Whitchurch TEST*/
312