xref: /petsc/src/dm/dt/tests/ex10.c (revision 45480ffeba491aca5d3f3b8c78679069fb317d32)
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 {
106d092c84bSBrandon Whitchurch   PetscBool simplex;
107d092c84bSBrandon Whitchurch   PetscInt  dim;
108d092c84bSBrandon Whitchurch } UserCtx;
109d092c84bSBrandon Whitchurch 
110d092c84bSBrandon Whitchurch /* Process command line options and initialize the UserCtx struct */
111d092c84bSBrandon Whitchurch static PetscErrorCode ProcessOptions(MPI_Comm comm,UserCtx *user)
112d092c84bSBrandon Whitchurch {
113d092c84bSBrandon Whitchurch   PetscErrorCode ierr;
114d092c84bSBrandon Whitchurch 
115d092c84bSBrandon Whitchurch   PetscFunctionBegin;
116d092c84bSBrandon Whitchurch   /* Default to  2D, triangle mesh.*/
117d092c84bSBrandon Whitchurch   user->simplex = PETSC_TRUE;
118d092c84bSBrandon Whitchurch   user->dim     = 2;
119d092c84bSBrandon Whitchurch 
120d092c84bSBrandon Whitchurch   ierr = PetscOptionsBegin(comm,"","PetscSpaceSum Tests","PetscSpace");CHKERRQ(ierr);
121d092c84bSBrandon Whitchurch   ierr = PetscOptionsBool("-simplex","Whether to use simplices (true) or tensor-product (false) cells in " "the mesh","ex8.c",user->simplex,
122d092c84bSBrandon Whitchurch                           &user->simplex,NULL);CHKERRQ(ierr);
123d092c84bSBrandon Whitchurch   ierr = PetscOptionsInt("-dim","Number of solution dimensions","ex8.c",user->dim,&user->dim,NULL);CHKERRQ(ierr);
124d092c84bSBrandon Whitchurch   ierr = PetscOptionsEnd();CHKERRQ(ierr);
125d092c84bSBrandon Whitchurch   PetscFunctionReturn(0);
126d092c84bSBrandon Whitchurch }
127d092c84bSBrandon Whitchurch 
128d092c84bSBrandon Whitchurch static PetscErrorCode CreateMesh(MPI_Comm comm,UserCtx *user,DM *mesh)
129d092c84bSBrandon Whitchurch {
130d092c84bSBrandon Whitchurch   PetscErrorCode   ierr;
131d092c84bSBrandon Whitchurch   DMLabel          label;
132d092c84bSBrandon Whitchurch   const char       *name  = "marker";
133d092c84bSBrandon Whitchurch   DM               dmDist = NULL;
134d092c84bSBrandon Whitchurch   PetscPartitioner part;
135d092c84bSBrandon Whitchurch 
136d092c84bSBrandon Whitchurch   PetscFunctionBegin;
137d092c84bSBrandon Whitchurch   /* Create box mesh from user parameters */
138d092c84bSBrandon Whitchurch   ierr = DMPlexCreateBoxMesh(comm,user->dim,user->simplex,NULL,NULL,NULL,NULL,PETSC_TRUE,mesh);CHKERRQ(ierr);
139d092c84bSBrandon Whitchurch 
140d092c84bSBrandon Whitchurch   /* Make sure the mesh gets properly distributed if running in parallel */
141d092c84bSBrandon Whitchurch   ierr = DMPlexGetPartitioner(*mesh,&part);CHKERRQ(ierr);
142d092c84bSBrandon Whitchurch   ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr);
143d092c84bSBrandon Whitchurch   ierr = DMPlexDistribute(*mesh,0,NULL,&dmDist);CHKERRQ(ierr);
144d092c84bSBrandon Whitchurch   if (dmDist) {
145d092c84bSBrandon Whitchurch     ierr  = DMDestroy(mesh);CHKERRQ(ierr);
146d092c84bSBrandon Whitchurch     *mesh = dmDist;
147d092c84bSBrandon Whitchurch   }
148d092c84bSBrandon Whitchurch 
149d092c84bSBrandon Whitchurch   /* Mark the boundaries, we will need this later when setting up the system of
150d092c84bSBrandon Whitchurch    * equations */
151d092c84bSBrandon Whitchurch   ierr = DMCreateLabel(*mesh,name);CHKERRQ(ierr);
152d092c84bSBrandon Whitchurch   ierr = DMGetLabel(*mesh,name,&label);CHKERRQ(ierr);
153d092c84bSBrandon Whitchurch   ierr = DMPlexMarkBoundaryFaces(*mesh,1,label);CHKERRQ(ierr);
154d092c84bSBrandon Whitchurch   ierr = DMPlexLabelComplete(*mesh,label);CHKERRQ(ierr);
155d092c84bSBrandon Whitchurch   ierr = DMLocalizeCoordinates(*mesh);CHKERRQ(ierr);
156d092c84bSBrandon Whitchurch   ierr = PetscObjectSetName((PetscObject) *mesh,"Mesh");CHKERRQ(ierr);
157d092c84bSBrandon Whitchurch 
158d092c84bSBrandon Whitchurch   /* Get any other mesh options from the command line */
159d092c84bSBrandon Whitchurch   ierr = DMSetApplicationContext(*mesh,user);CHKERRQ(ierr);
160d092c84bSBrandon Whitchurch   ierr = DMSetFromOptions(*mesh);CHKERRQ(ierr);
161d092c84bSBrandon Whitchurch   ierr = DMViewFromOptions(*mesh,NULL,"-dm_view");CHKERRQ(ierr);
162d092c84bSBrandon Whitchurch 
163d092c84bSBrandon Whitchurch   ierr = DMDestroy(&dmDist);CHKERRQ(ierr);
164d092c84bSBrandon Whitchurch   PetscFunctionReturn(0);
165d092c84bSBrandon Whitchurch }
166d092c84bSBrandon Whitchurch 
167d092c84bSBrandon Whitchurch /* Setup the system of equations that we wish to solve */
168d092c84bSBrandon Whitchurch static PetscErrorCode SetupProblem(DM dm,UserCtx *user)
169d092c84bSBrandon Whitchurch {
170*45480ffeSMatthew G. Knepley   PetscDS        ds;
171*45480ffeSMatthew G. Knepley   DMLabel        label;
172*45480ffeSMatthew G. Knepley   PetscWeakForm  wf;
173d092c84bSBrandon Whitchurch   const PetscInt id = 1;
174*45480ffeSMatthew G. Knepley   PetscInt       bd;
175*45480ffeSMatthew G. Knepley   PetscErrorCode ierr;
176d092c84bSBrandon Whitchurch 
177d092c84bSBrandon Whitchurch   PetscFunctionBegin;
178*45480ffeSMatthew G. Knepley   ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
179d092c84bSBrandon Whitchurch   /* All of these are independent of the user's choice of solution */
180*45480ffeSMatthew G. Knepley   ierr = PetscDSSetResidual(ds,0,f0_v,f1_v);CHKERRQ(ierr);
181*45480ffeSMatthew G. Knepley   ierr = PetscDSSetResidual(ds,1,f0_q_linear,NULL);CHKERRQ(ierr);
182*45480ffeSMatthew G. Knepley   ierr = PetscDSSetJacobian(ds,0,0,g0_vu,NULL,NULL,NULL);CHKERRQ(ierr);
183*45480ffeSMatthew G. Knepley   ierr = PetscDSSetJacobian(ds,0,1,NULL,NULL,g2_vp,NULL);CHKERRQ(ierr);
184*45480ffeSMatthew G. Knepley   ierr = PetscDSSetJacobian(ds,1,0,NULL,g1_qu,NULL,NULL);CHKERRQ(ierr);
185d092c84bSBrandon Whitchurch 
186*45480ffeSMatthew G. Knepley   ierr = DMGetLabel(dm, "marker", &label);CHKERRQ(ierr);
187*45480ffeSMatthew G. Knepley   ierr = PetscDSAddBoundary(ds,DM_BC_NATURAL,"Boundary Integral",label,1,&id,0,0,NULL,(void (*)(void))NULL,NULL,user,&bd);CHKERRQ(ierr);
188*45480ffeSMatthew G. Knepley   ierr = PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
189*45480ffeSMatthew G. Knepley   ierr = PetscWeakFormSetIndexBdResidual(wf, label, 1, 0, 0, f0_bd_u_linear, 0, NULL);CHKERRQ(ierr);
190*45480ffeSMatthew G. Knepley 
191*45480ffeSMatthew G. Knepley   ierr = PetscDSSetExactSolution(ds,0,linear_u,NULL);CHKERRQ(ierr);
192*45480ffeSMatthew G. Knepley   ierr = PetscDSSetExactSolution(ds,1,linear_divu,NULL);CHKERRQ(ierr);
193d092c84bSBrandon Whitchurch   PetscFunctionReturn(0);
194d092c84bSBrandon Whitchurch }
195d092c84bSBrandon Whitchurch 
196d092c84bSBrandon Whitchurch /* Create the finite element spaces we will use for this system */
197d092c84bSBrandon Whitchurch static PetscErrorCode SetupDiscretization(DM mesh,DM mesh_sum,PetscErrorCode (*setup)(DM,UserCtx*),UserCtx *user)
198d092c84bSBrandon Whitchurch {
199d092c84bSBrandon Whitchurch   DM             cdm = mesh,cdm_sum = mesh_sum;
200d092c84bSBrandon Whitchurch   PetscFE        u,divu,u_sum,divu_sum;
201d092c84bSBrandon Whitchurch   const PetscInt dim = user->dim;
202d092c84bSBrandon Whitchurch   PetscErrorCode ierr;
203d092c84bSBrandon Whitchurch 
204d092c84bSBrandon Whitchurch   PetscFunctionBegin;
205d092c84bSBrandon Whitchurch   /* Create FE objects and give them names so that options can be set from
206d092c84bSBrandon Whitchurch    * command line. Each field gets 2 instances (i.e. velocity and velocity_sum)created twice so that we can compare between approaches. */
207d092c84bSBrandon Whitchurch   ierr = PetscFECreateDefault(PetscObjectComm((PetscObject)mesh),dim,dim,user->simplex,"velocity_",-1,&u);CHKERRQ(ierr);
208d092c84bSBrandon Whitchurch   ierr = PetscObjectSetName((PetscObject)u,"velocity");CHKERRQ(ierr);
209d092c84bSBrandon Whitchurch   ierr = PetscFECreateDefault(PetscObjectComm((PetscObject)mesh_sum),dim,dim,user->simplex,"velocity_sum_",-1,&u_sum);CHKERRQ(ierr);
210d092c84bSBrandon Whitchurch   ierr = PetscObjectSetName((PetscObject)u_sum,"velocity_sum");CHKERRQ(ierr);
211d092c84bSBrandon Whitchurch   ierr = PetscFECreateDefault(PetscObjectComm((PetscObject)mesh),dim,1,user->simplex,"divu_",-1,&divu);CHKERRQ(ierr);
212d092c84bSBrandon Whitchurch   ierr = PetscObjectSetName((PetscObject)divu,"divu");CHKERRQ(ierr);
213d092c84bSBrandon Whitchurch   ierr = PetscFECreateDefault(PetscObjectComm((PetscObject)mesh_sum),dim,1,user->simplex,"divu_sum_",-1,&divu_sum);CHKERRQ(ierr);
214d092c84bSBrandon Whitchurch   ierr = PetscObjectSetName((PetscObject)divu_sum,"divu_sum");CHKERRQ(ierr);
215d092c84bSBrandon Whitchurch 
216d092c84bSBrandon Whitchurch   ierr = PetscFECopyQuadrature(u,divu);CHKERRQ(ierr);
217d092c84bSBrandon Whitchurch   ierr = PetscFECopyQuadrature(u_sum,divu_sum);CHKERRQ(ierr);
218d092c84bSBrandon Whitchurch 
219d092c84bSBrandon Whitchurch   /* Associate the FE objects with the mesh and setup the system */
220d092c84bSBrandon Whitchurch   ierr = DMSetField(mesh,0,NULL,(PetscObject)u);CHKERRQ(ierr);
221d092c84bSBrandon Whitchurch   ierr = DMSetField(mesh,1,NULL,(PetscObject)divu);CHKERRQ(ierr);
222d092c84bSBrandon Whitchurch   ierr = DMCreateDS(mesh);CHKERRQ(ierr);
223d092c84bSBrandon Whitchurch   ierr = (*setup)(mesh,user);CHKERRQ(ierr);
224d092c84bSBrandon Whitchurch 
225d092c84bSBrandon Whitchurch   ierr = DMSetField(mesh_sum,0,NULL,(PetscObject)u_sum);CHKERRQ(ierr);
226d092c84bSBrandon Whitchurch   ierr = DMSetField(mesh_sum,1,NULL,(PetscObject)divu_sum);CHKERRQ(ierr);
227d092c84bSBrandon Whitchurch   ierr = DMCreateDS(mesh_sum);CHKERRQ(ierr);
228d092c84bSBrandon Whitchurch   ierr = (*setup)(mesh_sum,user);CHKERRQ(ierr);
229d092c84bSBrandon Whitchurch 
230d092c84bSBrandon Whitchurch   while (cdm) {
231d092c84bSBrandon Whitchurch     ierr = DMCopyDisc(mesh,cdm);CHKERRQ(ierr);
232d092c84bSBrandon Whitchurch     ierr = DMGetCoarseDM(cdm,&cdm);CHKERRQ(ierr);
233d092c84bSBrandon Whitchurch   }
234d092c84bSBrandon Whitchurch 
235d092c84bSBrandon Whitchurch   while (cdm_sum) {
236d092c84bSBrandon Whitchurch     ierr = DMCopyDisc(mesh_sum,cdm_sum);CHKERRQ(ierr);
237d092c84bSBrandon Whitchurch     ierr = DMGetCoarseDM(cdm_sum,&cdm_sum);CHKERRQ(ierr);
238d092c84bSBrandon Whitchurch   }
239d092c84bSBrandon Whitchurch 
240d092c84bSBrandon Whitchurch   /* The Mesh now owns the fields, so we can destroy the FEs created in this
241d092c84bSBrandon Whitchurch    * function */
242d092c84bSBrandon Whitchurch   ierr = PetscFEDestroy(&u);CHKERRQ(ierr);
243d092c84bSBrandon Whitchurch   ierr = PetscFEDestroy(&divu);CHKERRQ(ierr);
244d092c84bSBrandon Whitchurch   ierr = PetscFEDestroy(&u_sum);CHKERRQ(ierr);
245d092c84bSBrandon Whitchurch   ierr = PetscFEDestroy(&divu_sum);CHKERRQ(ierr);
246d092c84bSBrandon Whitchurch   ierr = DMDestroy(&cdm);CHKERRQ(ierr);
247d092c84bSBrandon Whitchurch   ierr = DMDestroy(&cdm_sum);CHKERRQ(ierr);
248d092c84bSBrandon Whitchurch   PetscFunctionReturn(0);
249d092c84bSBrandon Whitchurch }
250d092c84bSBrandon Whitchurch 
251d092c84bSBrandon Whitchurch int main(int argc,char **argv)
252d092c84bSBrandon Whitchurch {
253d092c84bSBrandon Whitchurch   UserCtx         user;
254d092c84bSBrandon Whitchurch   DM              dm,dm_sum;
255d092c84bSBrandon Whitchurch   SNES            snes,snes_sum;
256d092c84bSBrandon Whitchurch   Vec             u,u_sum;
257d092c84bSBrandon Whitchurch   PetscReal       errNorm;
258d092c84bSBrandon Whitchurch   const PetscReal errTol = PETSC_SMALL;
259d092c84bSBrandon Whitchurch   PetscErrorCode  ierr;
260d092c84bSBrandon Whitchurch 
261d092c84bSBrandon Whitchurch   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
262d092c84bSBrandon Whitchurch   ierr = ProcessOptions(PETSC_COMM_WORLD,&user);CHKERRQ(ierr);
263d092c84bSBrandon Whitchurch 
264d092c84bSBrandon Whitchurch   /* Set up a snes for the standard approach, one space with 2 components */
265d092c84bSBrandon Whitchurch   ierr = SNESCreate(PETSC_COMM_WORLD,&snes);CHKERRQ(ierr);
266d092c84bSBrandon Whitchurch   ierr = CreateMesh(PETSC_COMM_WORLD,&user,&dm);CHKERRQ(ierr);
267d092c84bSBrandon Whitchurch   ierr = SNESSetDM(snes,dm);CHKERRQ(ierr);
268d092c84bSBrandon Whitchurch 
269d092c84bSBrandon Whitchurch   /* Set up a snes for the sum space approach, where each subspace of the sum space represents one component */
270d092c84bSBrandon Whitchurch   ierr = SNESCreate(PETSC_COMM_WORLD,&snes_sum);CHKERRQ(ierr);
271d092c84bSBrandon Whitchurch   ierr = CreateMesh(PETSC_COMM_WORLD,&user,&dm_sum);CHKERRQ(ierr);
272d092c84bSBrandon Whitchurch   ierr = SNESSetDM(snes_sum,dm_sum);CHKERRQ(ierr);
273d092c84bSBrandon Whitchurch   ierr = SetupDiscretization(dm,dm_sum,SetupProblem,&user);CHKERRQ(ierr);
274d092c84bSBrandon Whitchurch 
275d092c84bSBrandon Whitchurch   /* Set up and solve the system using standard approach. */
276d092c84bSBrandon Whitchurch   ierr = DMCreateGlobalVector(dm,&u);CHKERRQ(ierr);
277d092c84bSBrandon Whitchurch   ierr = VecSet(u,0.0);CHKERRQ(ierr);
278d092c84bSBrandon Whitchurch   ierr = PetscObjectSetName((PetscObject)u,"solution");CHKERRQ(ierr);
279d092c84bSBrandon Whitchurch   ierr = DMPlexSetSNESLocalFEM(dm,&user,&user,&user);CHKERRQ(ierr);
280d092c84bSBrandon Whitchurch   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
281348a1646SMatthew G. Knepley   ierr = DMSNESCheckFromOptions(snes,u);CHKERRQ(ierr);
282d092c84bSBrandon Whitchurch   ierr = SNESSolve(snes,NULL,u);CHKERRQ(ierr);
283d092c84bSBrandon Whitchurch   ierr = SNESGetSolution(snes,&u);CHKERRQ(ierr);
284d092c84bSBrandon Whitchurch   ierr = VecViewFromOptions(u,NULL,"-solution_view");CHKERRQ(ierr);
285d092c84bSBrandon Whitchurch 
286d092c84bSBrandon Whitchurch   /* Set up and solve the sum space system */
287d092c84bSBrandon Whitchurch   ierr = DMCreateGlobalVector(dm_sum,&u_sum);CHKERRQ(ierr);
288d092c84bSBrandon Whitchurch   ierr = VecSet(u_sum,0.0);CHKERRQ(ierr);
289d092c84bSBrandon Whitchurch   ierr = PetscObjectSetName((PetscObject)u_sum,"solution_sum");CHKERRQ(ierr);
290d092c84bSBrandon Whitchurch   ierr = DMPlexSetSNESLocalFEM(dm_sum,&user,&user,&user);CHKERRQ(ierr);
291d092c84bSBrandon Whitchurch   ierr = SNESSetFromOptions(snes_sum);CHKERRQ(ierr);
292348a1646SMatthew G. Knepley   ierr = DMSNESCheckFromOptions(snes_sum,u_sum);CHKERRQ(ierr);
293d092c84bSBrandon Whitchurch   ierr = SNESSolve(snes_sum,NULL,u_sum);CHKERRQ(ierr);
294d092c84bSBrandon Whitchurch   ierr = SNESGetSolution(snes_sum,&u_sum);CHKERRQ(ierr);
295d092c84bSBrandon Whitchurch   ierr = VecViewFromOptions(u_sum,NULL,"-solution_sum_view");CHKERRQ(ierr);
296d092c84bSBrandon Whitchurch 
297d092c84bSBrandon Whitchurch   /* Check if standard solution and sum space solution match to machine precision */
298d092c84bSBrandon Whitchurch   ierr = VecAXPY(u_sum,-1,u);CHKERRQ(ierr);
299d092c84bSBrandon Whitchurch   ierr = VecNorm(u_sum,NORM_2,&errNorm);CHKERRQ(ierr);
300d092c84bSBrandon Whitchurch   ierr = PetscPrintf(PETSC_COMM_WORLD,"Sum space provides the same solution as a regular space: %s",(errNorm < errTol) ? "true" : "false");CHKERRQ(
301d092c84bSBrandon Whitchurch     ierr);
302d092c84bSBrandon Whitchurch 
303d092c84bSBrandon Whitchurch   /* Cleanup */
304d092c84bSBrandon Whitchurch   ierr = VecDestroy(&u_sum);CHKERRQ(ierr);
305d092c84bSBrandon Whitchurch   ierr = VecDestroy(&u);CHKERRQ(ierr);
306d092c84bSBrandon Whitchurch   ierr = SNESDestroy(&snes_sum);CHKERRQ(ierr);
307d092c84bSBrandon Whitchurch   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
308d092c84bSBrandon Whitchurch   ierr = DMDestroy(&dm_sum);CHKERRQ(ierr);
309d092c84bSBrandon Whitchurch   ierr = DMDestroy(&dm);CHKERRQ(ierr);
310d092c84bSBrandon Whitchurch   ierr = PetscFinalize();
311d092c84bSBrandon Whitchurch   return ierr;
312d092c84bSBrandon Whitchurch }
313d092c84bSBrandon Whitchurch 
314d092c84bSBrandon Whitchurch /*TEST
315d092c84bSBrandon Whitchurch   test:
316d092c84bSBrandon Whitchurch     suffix: 2d_lagrange
317d092c84bSBrandon Whitchurch     requires: triangle
318d092c84bSBrandon Whitchurch     args: -dim 2 \
319d092c84bSBrandon Whitchurch       -simplex true \
320d092c84bSBrandon Whitchurch       -velocity_petscspace_degree 1 \
321d092c84bSBrandon Whitchurch       -velocity_petscspace_type poly \
322d092c84bSBrandon Whitchurch       -velocity_petscspace_components 2\
323d092c84bSBrandon Whitchurch       -velocity_petscdualspace_type lagrange \
324d092c84bSBrandon Whitchurch       -divu_petscspace_degree 0 \
325d092c84bSBrandon Whitchurch       -divu_petscspace_type poly \
326d092c84bSBrandon Whitchurch       -divu_petscdualspace_lagrange_continuity false \
327d092c84bSBrandon Whitchurch       -velocity_sum_petscfe_default_quadrature_order 1 \
328d092c84bSBrandon Whitchurch       -velocity_sum_petscspace_degree 1 \
329d092c84bSBrandon Whitchurch       -velocity_sum_petscspace_type sum \
330d092c84bSBrandon Whitchurch       -velocity_sum_petscspace_variables 2 \
331d092c84bSBrandon Whitchurch       -velocity_sum_petscspace_components 2 \
332d092c84bSBrandon Whitchurch       -velocity_sum_petscspace_sum_spaces 2 \
333d092c84bSBrandon Whitchurch       -velocity_sum_petscspace_sum_concatenate true \
334d092c84bSBrandon Whitchurch       -velocity_sum_petscdualspace_type lagrange \
335d092c84bSBrandon Whitchurch       -velocity_sum_subspace0_petscspace_type poly \
336d092c84bSBrandon Whitchurch       -velocity_sum_subspace0_petscspace_degree 1 \
337d092c84bSBrandon Whitchurch       -velocity_sum_subspace0_petscspace_variables 2 \
338d092c84bSBrandon Whitchurch       -velocity_sum_subspace0_petscspace_components 1 \
339d092c84bSBrandon Whitchurch       -velocity_sum_subspace1_petscspace_type poly \
340d092c84bSBrandon Whitchurch       -velocity_sum_subspace1_petscspace_degree 1 \
341d092c84bSBrandon Whitchurch       -velocity_sum_subspace1_petscspace_variables 2 \
342d092c84bSBrandon Whitchurch       -velocity_sum_subspace1_petscspace_components 1 \
343d092c84bSBrandon Whitchurch       -divu_sum_petscspace_degree 0 \
344d092c84bSBrandon Whitchurch       -divu_sum_petscspace_type sum \
345d092c84bSBrandon Whitchurch       -divu_sum_petscspace_variables 2 \
346d092c84bSBrandon Whitchurch       -divu_sum_petscspace_components 1 \
347d092c84bSBrandon Whitchurch       -divu_sum_petscspace_sum_spaces 1 \
348d092c84bSBrandon Whitchurch       -divu_sum_petscspace_sum_concatenate true \
349d092c84bSBrandon Whitchurch       -divu_sum_petscdualspace_lagrange_continuity false \
350d092c84bSBrandon Whitchurch       -divu_sum_subspace0_petscspace_type poly \
351d092c84bSBrandon Whitchurch       -divu_sum_subspace0_petscspace_degree 0 \
352d092c84bSBrandon Whitchurch       -divu_sum_subspace0_petscspace_variables 2 \
353d092c84bSBrandon Whitchurch       -divu_sum_subspace0_petscspace_components 1 \
354d092c84bSBrandon Whitchurch       -dm_refine 0 \
355d092c84bSBrandon Whitchurch       -snes_error_if_not_converged \
356d092c84bSBrandon Whitchurch       -ksp_rtol 1e-10 \
357d092c84bSBrandon Whitchurch       -ksp_error_if_not_converged \
358d092c84bSBrandon Whitchurch       -pc_type fieldsplit\
359d092c84bSBrandon Whitchurch       -pc_fieldsplit_type schur\
360d092c84bSBrandon Whitchurch       -divu_sum_petscdualspace_lagrange_continuity false \
361d092c84bSBrandon Whitchurch       -pc_fieldsplit_schur_precondition full
362d092c84bSBrandon Whitchurch TEST*/
363