xref: /petsc/src/dm/dt/fv/interface/fv.c (revision c0ce001ee2722011d4ad0f8c03a06ea86b62b0c2)
1*c0ce001eSMatthew G. Knepley #include <petsc/private/petscfvimpl.h> /*I "petscfv.h" I*/
2*c0ce001eSMatthew G. Knepley #include <petsc/private/dmpleximpl.h> /* For CellRefiner */
3*c0ce001eSMatthew G. Knepley #include <petscds.h>
4*c0ce001eSMatthew G. Knepley 
5*c0ce001eSMatthew G. Knepley PetscClassId PETSCLIMITER_CLASSID = 0;
6*c0ce001eSMatthew G. Knepley 
7*c0ce001eSMatthew G. Knepley PetscFunctionList PetscLimiterList              = NULL;
8*c0ce001eSMatthew G. Knepley PetscBool         PetscLimiterRegisterAllCalled = PETSC_FALSE;
9*c0ce001eSMatthew G. Knepley 
10*c0ce001eSMatthew G. Knepley PetscBool Limitercite = PETSC_FALSE;
11*c0ce001eSMatthew G. Knepley const char LimiterCitation[] = "@article{BergerAftosmisMurman2005,\n"
12*c0ce001eSMatthew G. Knepley                                "  title   = {Analysis of slope limiters on irregular grids},\n"
13*c0ce001eSMatthew G. Knepley                                "  journal = {AIAA paper},\n"
14*c0ce001eSMatthew G. Knepley                                "  author  = {Marsha Berger and Michael J. Aftosmis and Scott M. Murman},\n"
15*c0ce001eSMatthew G. Knepley                                "  volume  = {490},\n"
16*c0ce001eSMatthew G. Knepley                                "  year    = {2005}\n}\n";
17*c0ce001eSMatthew G. Knepley 
18*c0ce001eSMatthew G. Knepley /*@C
19*c0ce001eSMatthew G. Knepley   PetscLimiterRegister - Adds a new PetscLimiter implementation
20*c0ce001eSMatthew G. Knepley 
21*c0ce001eSMatthew G. Knepley   Not Collective
22*c0ce001eSMatthew G. Knepley 
23*c0ce001eSMatthew G. Knepley   Input Parameters:
24*c0ce001eSMatthew G. Knepley + name        - The name of a new user-defined creation routine
25*c0ce001eSMatthew G. Knepley - create_func - The creation routine itself
26*c0ce001eSMatthew G. Knepley 
27*c0ce001eSMatthew G. Knepley   Notes:
28*c0ce001eSMatthew G. Knepley   PetscLimiterRegister() may be called multiple times to add several user-defined PetscLimiters
29*c0ce001eSMatthew G. Knepley 
30*c0ce001eSMatthew G. Knepley   Sample usage:
31*c0ce001eSMatthew G. Knepley .vb
32*c0ce001eSMatthew G. Knepley     PetscLimiterRegister("my_lim", MyPetscLimiterCreate);
33*c0ce001eSMatthew G. Knepley .ve
34*c0ce001eSMatthew G. Knepley 
35*c0ce001eSMatthew G. Knepley   Then, your PetscLimiter type can be chosen with the procedural interface via
36*c0ce001eSMatthew G. Knepley .vb
37*c0ce001eSMatthew G. Knepley     PetscLimiterCreate(MPI_Comm, PetscLimiter *);
38*c0ce001eSMatthew G. Knepley     PetscLimiterSetType(PetscLimiter, "my_lim");
39*c0ce001eSMatthew G. Knepley .ve
40*c0ce001eSMatthew G. Knepley    or at runtime via the option
41*c0ce001eSMatthew G. Knepley .vb
42*c0ce001eSMatthew G. Knepley     -petsclimiter_type my_lim
43*c0ce001eSMatthew G. Knepley .ve
44*c0ce001eSMatthew G. Knepley 
45*c0ce001eSMatthew G. Knepley   Level: advanced
46*c0ce001eSMatthew G. Knepley 
47*c0ce001eSMatthew G. Knepley .seealso: PetscLimiterRegisterAll(), PetscLimiterRegisterDestroy()
48*c0ce001eSMatthew G. Knepley 
49*c0ce001eSMatthew G. Knepley @*/
50*c0ce001eSMatthew G. Knepley PetscErrorCode PetscLimiterRegister(const char sname[], PetscErrorCode (*function)(PetscLimiter))
51*c0ce001eSMatthew G. Knepley {
52*c0ce001eSMatthew G. Knepley   PetscErrorCode ierr;
53*c0ce001eSMatthew G. Knepley 
54*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
55*c0ce001eSMatthew G. Knepley   ierr = PetscFunctionListAdd(&PetscLimiterList, sname, function);CHKERRQ(ierr);
56*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
57*c0ce001eSMatthew G. Knepley }
58*c0ce001eSMatthew G. Knepley 
59*c0ce001eSMatthew G. Knepley /*@C
60*c0ce001eSMatthew G. Knepley   PetscLimiterSetType - Builds a particular PetscLimiter
61*c0ce001eSMatthew G. Knepley 
62*c0ce001eSMatthew G. Knepley   Collective on lim
63*c0ce001eSMatthew G. Knepley 
64*c0ce001eSMatthew G. Knepley   Input Parameters:
65*c0ce001eSMatthew G. Knepley + lim  - The PetscLimiter object
66*c0ce001eSMatthew G. Knepley - name - The kind of limiter
67*c0ce001eSMatthew G. Knepley 
68*c0ce001eSMatthew G. Knepley   Options Database Key:
69*c0ce001eSMatthew G. Knepley . -petsclimiter_type <type> - Sets the PetscLimiter type; use -help for a list of available types
70*c0ce001eSMatthew G. Knepley 
71*c0ce001eSMatthew G. Knepley   Level: intermediate
72*c0ce001eSMatthew G. Knepley 
73*c0ce001eSMatthew G. Knepley .seealso: PetscLimiterGetType(), PetscLimiterCreate()
74*c0ce001eSMatthew G. Knepley @*/
75*c0ce001eSMatthew G. Knepley PetscErrorCode PetscLimiterSetType(PetscLimiter lim, PetscLimiterType name)
76*c0ce001eSMatthew G. Knepley {
77*c0ce001eSMatthew G. Knepley   PetscErrorCode (*r)(PetscLimiter);
78*c0ce001eSMatthew G. Knepley   PetscBool      match;
79*c0ce001eSMatthew G. Knepley   PetscErrorCode ierr;
80*c0ce001eSMatthew G. Knepley 
81*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
82*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1);
83*c0ce001eSMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject) lim, name, &match);CHKERRQ(ierr);
84*c0ce001eSMatthew G. Knepley   if (match) PetscFunctionReturn(0);
85*c0ce001eSMatthew G. Knepley 
86*c0ce001eSMatthew G. Knepley   ierr = PetscLimiterRegisterAll();CHKERRQ(ierr);
87*c0ce001eSMatthew G. Knepley   ierr = PetscFunctionListFind(PetscLimiterList, name, &r);CHKERRQ(ierr);
88*c0ce001eSMatthew G. Knepley   if (!r) SETERRQ1(PetscObjectComm((PetscObject) lim), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscLimiter type: %s", name);
89*c0ce001eSMatthew G. Knepley 
90*c0ce001eSMatthew G. Knepley   if (lim->ops->destroy) {
91*c0ce001eSMatthew G. Knepley     ierr              = (*lim->ops->destroy)(lim);CHKERRQ(ierr);
92*c0ce001eSMatthew G. Knepley     lim->ops->destroy = NULL;
93*c0ce001eSMatthew G. Knepley   }
94*c0ce001eSMatthew G. Knepley   ierr = (*r)(lim);CHKERRQ(ierr);
95*c0ce001eSMatthew G. Knepley   ierr = PetscObjectChangeTypeName((PetscObject) lim, name);CHKERRQ(ierr);
96*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
97*c0ce001eSMatthew G. Knepley }
98*c0ce001eSMatthew G. Knepley 
99*c0ce001eSMatthew G. Knepley /*@C
100*c0ce001eSMatthew G. Knepley   PetscLimiterGetType - Gets the PetscLimiter type name (as a string) from the object.
101*c0ce001eSMatthew G. Knepley 
102*c0ce001eSMatthew G. Knepley   Not Collective
103*c0ce001eSMatthew G. Knepley 
104*c0ce001eSMatthew G. Knepley   Input Parameter:
105*c0ce001eSMatthew G. Knepley . lim  - The PetscLimiter
106*c0ce001eSMatthew G. Knepley 
107*c0ce001eSMatthew G. Knepley   Output Parameter:
108*c0ce001eSMatthew G. Knepley . name - The PetscLimiter type name
109*c0ce001eSMatthew G. Knepley 
110*c0ce001eSMatthew G. Knepley   Level: intermediate
111*c0ce001eSMatthew G. Knepley 
112*c0ce001eSMatthew G. Knepley .seealso: PetscLimiterSetType(), PetscLimiterCreate()
113*c0ce001eSMatthew G. Knepley @*/
114*c0ce001eSMatthew G. Knepley PetscErrorCode PetscLimiterGetType(PetscLimiter lim, PetscLimiterType *name)
115*c0ce001eSMatthew G. Knepley {
116*c0ce001eSMatthew G. Knepley   PetscErrorCode ierr;
117*c0ce001eSMatthew G. Knepley 
118*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
119*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1);
120*c0ce001eSMatthew G. Knepley   PetscValidPointer(name, 2);
121*c0ce001eSMatthew G. Knepley   ierr = PetscLimiterRegisterAll();CHKERRQ(ierr);
122*c0ce001eSMatthew G. Knepley   *name = ((PetscObject) lim)->type_name;
123*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
124*c0ce001eSMatthew G. Knepley }
125*c0ce001eSMatthew G. Knepley 
126*c0ce001eSMatthew G. Knepley /*@C
127*c0ce001eSMatthew G. Knepley   PetscLimiterView - Views a PetscLimiter
128*c0ce001eSMatthew G. Knepley 
129*c0ce001eSMatthew G. Knepley   Collective on lim
130*c0ce001eSMatthew G. Knepley 
131*c0ce001eSMatthew G. Knepley   Input Parameter:
132*c0ce001eSMatthew G. Knepley + lim - the PetscLimiter object to view
133*c0ce001eSMatthew G. Knepley - v   - the viewer
134*c0ce001eSMatthew G. Knepley 
135*c0ce001eSMatthew G. Knepley   Level: developer
136*c0ce001eSMatthew G. Knepley 
137*c0ce001eSMatthew G. Knepley .seealso: PetscLimiterDestroy()
138*c0ce001eSMatthew G. Knepley @*/
139*c0ce001eSMatthew G. Knepley PetscErrorCode PetscLimiterView(PetscLimiter lim, PetscViewer v)
140*c0ce001eSMatthew G. Knepley {
141*c0ce001eSMatthew G. Knepley   PetscErrorCode ierr;
142*c0ce001eSMatthew G. Knepley 
143*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
144*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1);
145*c0ce001eSMatthew G. Knepley   if (!v) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) lim), &v);CHKERRQ(ierr);}
146*c0ce001eSMatthew G. Knepley   if (lim->ops->view) {ierr = (*lim->ops->view)(lim, v);CHKERRQ(ierr);}
147*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
148*c0ce001eSMatthew G. Knepley }
149*c0ce001eSMatthew G. Knepley 
150*c0ce001eSMatthew G. Knepley /*@
151*c0ce001eSMatthew G. Knepley   PetscLimiterSetFromOptions - sets parameters in a PetscLimiter from the options database
152*c0ce001eSMatthew G. Knepley 
153*c0ce001eSMatthew G. Knepley   Collective on lim
154*c0ce001eSMatthew G. Knepley 
155*c0ce001eSMatthew G. Knepley   Input Parameter:
156*c0ce001eSMatthew G. Knepley . lim - the PetscLimiter object to set options for
157*c0ce001eSMatthew G. Knepley 
158*c0ce001eSMatthew G. Knepley   Level: developer
159*c0ce001eSMatthew G. Knepley 
160*c0ce001eSMatthew G. Knepley .seealso: PetscLimiterView()
161*c0ce001eSMatthew G. Knepley @*/
162*c0ce001eSMatthew G. Knepley PetscErrorCode PetscLimiterSetFromOptions(PetscLimiter lim)
163*c0ce001eSMatthew G. Knepley {
164*c0ce001eSMatthew G. Knepley   const char    *defaultType;
165*c0ce001eSMatthew G. Knepley   char           name[256];
166*c0ce001eSMatthew G. Knepley   PetscBool      flg;
167*c0ce001eSMatthew G. Knepley   PetscErrorCode ierr;
168*c0ce001eSMatthew G. Knepley 
169*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
170*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1);
171*c0ce001eSMatthew G. Knepley   if (!((PetscObject) lim)->type_name) defaultType = PETSCLIMITERSIN;
172*c0ce001eSMatthew G. Knepley   else                                 defaultType = ((PetscObject) lim)->type_name;
173*c0ce001eSMatthew G. Knepley   ierr = PetscLimiterRegisterAll();CHKERRQ(ierr);
174*c0ce001eSMatthew G. Knepley 
175*c0ce001eSMatthew G. Knepley   ierr = PetscObjectOptionsBegin((PetscObject) lim);CHKERRQ(ierr);
176*c0ce001eSMatthew G. Knepley   ierr = PetscOptionsFList("-petsclimiter_type", "Finite volume slope limiter", "PetscLimiterSetType", PetscLimiterList, defaultType, name, 256, &flg);CHKERRQ(ierr);
177*c0ce001eSMatthew G. Knepley   if (flg) {
178*c0ce001eSMatthew G. Knepley     ierr = PetscLimiterSetType(lim, name);CHKERRQ(ierr);
179*c0ce001eSMatthew G. Knepley   } else if (!((PetscObject) lim)->type_name) {
180*c0ce001eSMatthew G. Knepley     ierr = PetscLimiterSetType(lim, defaultType);CHKERRQ(ierr);
181*c0ce001eSMatthew G. Knepley   }
182*c0ce001eSMatthew G. Knepley   if (lim->ops->setfromoptions) {ierr = (*lim->ops->setfromoptions)(lim);CHKERRQ(ierr);}
183*c0ce001eSMatthew G. Knepley   /* process any options handlers added with PetscObjectAddOptionsHandler() */
184*c0ce001eSMatthew G. Knepley   ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject) lim);CHKERRQ(ierr);
185*c0ce001eSMatthew G. Knepley   ierr = PetscOptionsEnd();CHKERRQ(ierr);
186*c0ce001eSMatthew G. Knepley   ierr = PetscLimiterViewFromOptions(lim, NULL, "-petsclimiter_view");CHKERRQ(ierr);
187*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
188*c0ce001eSMatthew G. Knepley }
189*c0ce001eSMatthew G. Knepley 
190*c0ce001eSMatthew G. Knepley /*@C
191*c0ce001eSMatthew G. Knepley   PetscLimiterSetUp - Construct data structures for the PetscLimiter
192*c0ce001eSMatthew G. Knepley 
193*c0ce001eSMatthew G. Knepley   Collective on lim
194*c0ce001eSMatthew G. Knepley 
195*c0ce001eSMatthew G. Knepley   Input Parameter:
196*c0ce001eSMatthew G. Knepley . lim - the PetscLimiter object to setup
197*c0ce001eSMatthew G. Knepley 
198*c0ce001eSMatthew G. Knepley   Level: developer
199*c0ce001eSMatthew G. Knepley 
200*c0ce001eSMatthew G. Knepley .seealso: PetscLimiterView(), PetscLimiterDestroy()
201*c0ce001eSMatthew G. Knepley @*/
202*c0ce001eSMatthew G. Knepley PetscErrorCode PetscLimiterSetUp(PetscLimiter lim)
203*c0ce001eSMatthew G. Knepley {
204*c0ce001eSMatthew G. Knepley   PetscErrorCode ierr;
205*c0ce001eSMatthew G. Knepley 
206*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
207*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1);
208*c0ce001eSMatthew G. Knepley   if (lim->ops->setup) {ierr = (*lim->ops->setup)(lim);CHKERRQ(ierr);}
209*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
210*c0ce001eSMatthew G. Knepley }
211*c0ce001eSMatthew G. Knepley 
212*c0ce001eSMatthew G. Knepley /*@
213*c0ce001eSMatthew G. Knepley   PetscLimiterDestroy - Destroys a PetscLimiter object
214*c0ce001eSMatthew G. Knepley 
215*c0ce001eSMatthew G. Knepley   Collective on lim
216*c0ce001eSMatthew G. Knepley 
217*c0ce001eSMatthew G. Knepley   Input Parameter:
218*c0ce001eSMatthew G. Knepley . lim - the PetscLimiter object to destroy
219*c0ce001eSMatthew G. Knepley 
220*c0ce001eSMatthew G. Knepley   Level: developer
221*c0ce001eSMatthew G. Knepley 
222*c0ce001eSMatthew G. Knepley .seealso: PetscLimiterView()
223*c0ce001eSMatthew G. Knepley @*/
224*c0ce001eSMatthew G. Knepley PetscErrorCode PetscLimiterDestroy(PetscLimiter *lim)
225*c0ce001eSMatthew G. Knepley {
226*c0ce001eSMatthew G. Knepley   PetscErrorCode ierr;
227*c0ce001eSMatthew G. Knepley 
228*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
229*c0ce001eSMatthew G. Knepley   if (!*lim) PetscFunctionReturn(0);
230*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific((*lim), PETSCLIMITER_CLASSID, 1);
231*c0ce001eSMatthew G. Knepley 
232*c0ce001eSMatthew G. Knepley   if (--((PetscObject)(*lim))->refct > 0) {*lim = 0; PetscFunctionReturn(0);}
233*c0ce001eSMatthew G. Knepley   ((PetscObject) (*lim))->refct = 0;
234*c0ce001eSMatthew G. Knepley 
235*c0ce001eSMatthew G. Knepley   if ((*lim)->ops->destroy) {ierr = (*(*lim)->ops->destroy)(*lim);CHKERRQ(ierr);}
236*c0ce001eSMatthew G. Knepley   ierr = PetscHeaderDestroy(lim);CHKERRQ(ierr);
237*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
238*c0ce001eSMatthew G. Knepley }
239*c0ce001eSMatthew G. Knepley 
240*c0ce001eSMatthew G. Knepley /*@
241*c0ce001eSMatthew G. Knepley   PetscLimiterCreate - Creates an empty PetscLimiter object. The type can then be set with PetscLimiterSetType().
242*c0ce001eSMatthew G. Knepley 
243*c0ce001eSMatthew G. Knepley   Collective
244*c0ce001eSMatthew G. Knepley 
245*c0ce001eSMatthew G. Knepley   Input Parameter:
246*c0ce001eSMatthew G. Knepley . comm - The communicator for the PetscLimiter object
247*c0ce001eSMatthew G. Knepley 
248*c0ce001eSMatthew G. Knepley   Output Parameter:
249*c0ce001eSMatthew G. Knepley . lim - The PetscLimiter object
250*c0ce001eSMatthew G. Knepley 
251*c0ce001eSMatthew G. Knepley   Level: beginner
252*c0ce001eSMatthew G. Knepley 
253*c0ce001eSMatthew G. Knepley .seealso: PetscLimiterSetType(), PETSCLIMITERSIN
254*c0ce001eSMatthew G. Knepley @*/
255*c0ce001eSMatthew G. Knepley PetscErrorCode PetscLimiterCreate(MPI_Comm comm, PetscLimiter *lim)
256*c0ce001eSMatthew G. Knepley {
257*c0ce001eSMatthew G. Knepley   PetscLimiter   l;
258*c0ce001eSMatthew G. Knepley   PetscErrorCode ierr;
259*c0ce001eSMatthew G. Knepley 
260*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
261*c0ce001eSMatthew G. Knepley   PetscValidPointer(lim, 2);
262*c0ce001eSMatthew G. Knepley   ierr = PetscCitationsRegister(LimiterCitation,&Limitercite);CHKERRQ(ierr);
263*c0ce001eSMatthew G. Knepley   *lim = NULL;
264*c0ce001eSMatthew G. Knepley   ierr = PetscFVInitializePackage();CHKERRQ(ierr);
265*c0ce001eSMatthew G. Knepley 
266*c0ce001eSMatthew G. Knepley   ierr = PetscHeaderCreate(l, PETSCLIMITER_CLASSID, "PetscLimiter", "Finite Volume Slope Limiter", "PetscLimiter", comm, PetscLimiterDestroy, PetscLimiterView);CHKERRQ(ierr);
267*c0ce001eSMatthew G. Knepley 
268*c0ce001eSMatthew G. Knepley   *lim = l;
269*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
270*c0ce001eSMatthew G. Knepley }
271*c0ce001eSMatthew G. Knepley 
272*c0ce001eSMatthew G. Knepley /* Limiters given in symmetric form following Berger, Aftosmis, and Murman 2005
273*c0ce001eSMatthew G. Knepley  *
274*c0ce001eSMatthew G. Knepley  * The classical flux-limited formulation is psi(r) where
275*c0ce001eSMatthew G. Knepley  *
276*c0ce001eSMatthew G. Knepley  * r = (u[0] - u[-1]) / (u[1] - u[0])
277*c0ce001eSMatthew G. Knepley  *
278*c0ce001eSMatthew G. Knepley  * The second order TVD region is bounded by
279*c0ce001eSMatthew G. Knepley  *
280*c0ce001eSMatthew G. Knepley  * psi_minmod(r) = min(r,1)      and        psi_superbee(r) = min(2, 2r, max(1,r))
281*c0ce001eSMatthew G. Knepley  *
282*c0ce001eSMatthew G. Knepley  * where all limiters are implicitly clipped to be non-negative. A more convenient slope-limited form is psi(r) =
283*c0ce001eSMatthew G. Knepley  * phi(r)(r+1)/2 in which the reconstructed interface values are
284*c0ce001eSMatthew G. Knepley  *
285*c0ce001eSMatthew G. Knepley  * u(v) = u[0] + phi(r) (grad u)[0] v
286*c0ce001eSMatthew G. Knepley  *
287*c0ce001eSMatthew G. Knepley  * where v is the vector from centroid to quadrature point. In these variables, the usual limiters become
288*c0ce001eSMatthew G. Knepley  *
289*c0ce001eSMatthew G. Knepley  * phi_minmod(r) = 2 min(1/(1+r),r/(1+r))   phi_superbee(r) = 2 min(2/(1+r), 2r/(1+r), max(1,r)/(1+r))
290*c0ce001eSMatthew G. Knepley  *
291*c0ce001eSMatthew G. Knepley  * For a nicer symmetric formulation, rewrite in terms of
292*c0ce001eSMatthew G. Knepley  *
293*c0ce001eSMatthew G. Knepley  * f = (u[0] - u[-1]) / (u[1] - u[-1])
294*c0ce001eSMatthew G. Knepley  *
295*c0ce001eSMatthew G. Knepley  * where r(f) = f/(1-f). Not that r(1-f) = (1-f)/f = 1/r(f) so the symmetry condition
296*c0ce001eSMatthew G. Knepley  *
297*c0ce001eSMatthew G. Knepley  * phi(r) = phi(1/r)
298*c0ce001eSMatthew G. Knepley  *
299*c0ce001eSMatthew G. Knepley  * becomes
300*c0ce001eSMatthew G. Knepley  *
301*c0ce001eSMatthew G. Knepley  * w(f) = w(1-f).
302*c0ce001eSMatthew G. Knepley  *
303*c0ce001eSMatthew G. Knepley  * The limiters below implement this final form w(f). The reference methods are
304*c0ce001eSMatthew G. Knepley  *
305*c0ce001eSMatthew G. Knepley  * w_minmod(f) = 2 min(f,(1-f))             w_superbee(r) = 4 min((1-f), f)
306*c0ce001eSMatthew G. Knepley  * */
307*c0ce001eSMatthew G. Knepley PetscErrorCode PetscLimiterLimit(PetscLimiter lim, PetscReal flim, PetscReal *phi)
308*c0ce001eSMatthew G. Knepley {
309*c0ce001eSMatthew G. Knepley   PetscErrorCode ierr;
310*c0ce001eSMatthew G. Knepley 
311*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
312*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1);
313*c0ce001eSMatthew G. Knepley   PetscValidPointer(phi, 3);
314*c0ce001eSMatthew G. Knepley   ierr = (*lim->ops->limit)(lim, flim, phi);CHKERRQ(ierr);
315*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
316*c0ce001eSMatthew G. Knepley }
317*c0ce001eSMatthew G. Knepley 
318*c0ce001eSMatthew G. Knepley PetscErrorCode PetscLimiterDestroy_Sin(PetscLimiter lim)
319*c0ce001eSMatthew G. Knepley {
320*c0ce001eSMatthew G. Knepley   PetscLimiter_Sin *l = (PetscLimiter_Sin *) lim->data;
321*c0ce001eSMatthew G. Knepley   PetscErrorCode    ierr;
322*c0ce001eSMatthew G. Knepley 
323*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
324*c0ce001eSMatthew G. Knepley   ierr = PetscFree(l);CHKERRQ(ierr);
325*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
326*c0ce001eSMatthew G. Knepley }
327*c0ce001eSMatthew G. Knepley 
328*c0ce001eSMatthew G. Knepley PetscErrorCode PetscLimiterView_Sin_Ascii(PetscLimiter lim, PetscViewer viewer)
329*c0ce001eSMatthew G. Knepley {
330*c0ce001eSMatthew G. Knepley   PetscViewerFormat format;
331*c0ce001eSMatthew G. Knepley   PetscErrorCode    ierr;
332*c0ce001eSMatthew G. Knepley 
333*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
334*c0ce001eSMatthew G. Knepley   ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr);
335*c0ce001eSMatthew G. Knepley   ierr = PetscViewerASCIIPrintf(viewer, "Sin Slope Limiter:\n");CHKERRQ(ierr);
336*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
337*c0ce001eSMatthew G. Knepley }
338*c0ce001eSMatthew G. Knepley 
339*c0ce001eSMatthew G. Knepley PetscErrorCode PetscLimiterView_Sin(PetscLimiter lim, PetscViewer viewer)
340*c0ce001eSMatthew G. Knepley {
341*c0ce001eSMatthew G. Knepley   PetscBool      iascii;
342*c0ce001eSMatthew G. Knepley   PetscErrorCode ierr;
343*c0ce001eSMatthew G. Knepley 
344*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
345*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1);
346*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
347*c0ce001eSMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
348*c0ce001eSMatthew G. Knepley   if (iascii) {ierr = PetscLimiterView_Sin_Ascii(lim, viewer);CHKERRQ(ierr);}
349*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
350*c0ce001eSMatthew G. Knepley }
351*c0ce001eSMatthew G. Knepley 
352*c0ce001eSMatthew G. Knepley PetscErrorCode PetscLimiterLimit_Sin(PetscLimiter lim, PetscReal f, PetscReal *phi)
353*c0ce001eSMatthew G. Knepley {
354*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
355*c0ce001eSMatthew G. Knepley   *phi = PetscSinReal(PETSC_PI*PetscMax(0, PetscMin(f, 1)));
356*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
357*c0ce001eSMatthew G. Knepley }
358*c0ce001eSMatthew G. Knepley 
359*c0ce001eSMatthew G. Knepley PetscErrorCode PetscLimiterInitialize_Sin(PetscLimiter lim)
360*c0ce001eSMatthew G. Knepley {
361*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
362*c0ce001eSMatthew G. Knepley   lim->ops->view    = PetscLimiterView_Sin;
363*c0ce001eSMatthew G. Knepley   lim->ops->destroy = PetscLimiterDestroy_Sin;
364*c0ce001eSMatthew G. Knepley   lim->ops->limit   = PetscLimiterLimit_Sin;
365*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
366*c0ce001eSMatthew G. Knepley }
367*c0ce001eSMatthew G. Knepley 
368*c0ce001eSMatthew G. Knepley /*MC
369*c0ce001eSMatthew G. Knepley   PETSCLIMITERSIN = "sin" - A PetscLimiter object
370*c0ce001eSMatthew G. Knepley 
371*c0ce001eSMatthew G. Knepley   Level: intermediate
372*c0ce001eSMatthew G. Knepley 
373*c0ce001eSMatthew G. Knepley .seealso: PetscLimiterType, PetscLimiterCreate(), PetscLimiterSetType()
374*c0ce001eSMatthew G. Knepley M*/
375*c0ce001eSMatthew G. Knepley 
376*c0ce001eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscLimiterCreate_Sin(PetscLimiter lim)
377*c0ce001eSMatthew G. Knepley {
378*c0ce001eSMatthew G. Knepley   PetscLimiter_Sin *l;
379*c0ce001eSMatthew G. Knepley   PetscErrorCode    ierr;
380*c0ce001eSMatthew G. Knepley 
381*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
382*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1);
383*c0ce001eSMatthew G. Knepley   ierr      = PetscNewLog(lim, &l);CHKERRQ(ierr);
384*c0ce001eSMatthew G. Knepley   lim->data = l;
385*c0ce001eSMatthew G. Knepley 
386*c0ce001eSMatthew G. Knepley   ierr = PetscLimiterInitialize_Sin(lim);CHKERRQ(ierr);
387*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
388*c0ce001eSMatthew G. Knepley }
389*c0ce001eSMatthew G. Knepley 
390*c0ce001eSMatthew G. Knepley PetscErrorCode PetscLimiterDestroy_Zero(PetscLimiter lim)
391*c0ce001eSMatthew G. Knepley {
392*c0ce001eSMatthew G. Knepley   PetscLimiter_Zero *l = (PetscLimiter_Zero *) lim->data;
393*c0ce001eSMatthew G. Knepley   PetscErrorCode    ierr;
394*c0ce001eSMatthew G. Knepley 
395*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
396*c0ce001eSMatthew G. Knepley   ierr = PetscFree(l);CHKERRQ(ierr);
397*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
398*c0ce001eSMatthew G. Knepley }
399*c0ce001eSMatthew G. Knepley 
400*c0ce001eSMatthew G. Knepley PetscErrorCode PetscLimiterView_Zero_Ascii(PetscLimiter lim, PetscViewer viewer)
401*c0ce001eSMatthew G. Knepley {
402*c0ce001eSMatthew G. Knepley   PetscViewerFormat format;
403*c0ce001eSMatthew G. Knepley   PetscErrorCode    ierr;
404*c0ce001eSMatthew G. Knepley 
405*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
406*c0ce001eSMatthew G. Knepley   ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr);
407*c0ce001eSMatthew G. Knepley   ierr = PetscViewerASCIIPrintf(viewer, "Zero Slope Limiter:\n");CHKERRQ(ierr);
408*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
409*c0ce001eSMatthew G. Knepley }
410*c0ce001eSMatthew G. Knepley 
411*c0ce001eSMatthew G. Knepley PetscErrorCode PetscLimiterView_Zero(PetscLimiter lim, PetscViewer viewer)
412*c0ce001eSMatthew G. Knepley {
413*c0ce001eSMatthew G. Knepley   PetscBool      iascii;
414*c0ce001eSMatthew G. Knepley   PetscErrorCode ierr;
415*c0ce001eSMatthew G. Knepley 
416*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
417*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1);
418*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
419*c0ce001eSMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
420*c0ce001eSMatthew G. Knepley   if (iascii) {ierr = PetscLimiterView_Zero_Ascii(lim, viewer);CHKERRQ(ierr);}
421*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
422*c0ce001eSMatthew G. Knepley }
423*c0ce001eSMatthew G. Knepley 
424*c0ce001eSMatthew G. Knepley PetscErrorCode PetscLimiterLimit_Zero(PetscLimiter lim, PetscReal f, PetscReal *phi)
425*c0ce001eSMatthew G. Knepley {
426*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
427*c0ce001eSMatthew G. Knepley   *phi = 0.0;
428*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
429*c0ce001eSMatthew G. Knepley }
430*c0ce001eSMatthew G. Knepley 
431*c0ce001eSMatthew G. Knepley PetscErrorCode PetscLimiterInitialize_Zero(PetscLimiter lim)
432*c0ce001eSMatthew G. Knepley {
433*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
434*c0ce001eSMatthew G. Knepley   lim->ops->view    = PetscLimiterView_Zero;
435*c0ce001eSMatthew G. Knepley   lim->ops->destroy = PetscLimiterDestroy_Zero;
436*c0ce001eSMatthew G. Knepley   lim->ops->limit   = PetscLimiterLimit_Zero;
437*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
438*c0ce001eSMatthew G. Knepley }
439*c0ce001eSMatthew G. Knepley 
440*c0ce001eSMatthew G. Knepley /*MC
441*c0ce001eSMatthew G. Knepley   PETSCLIMITERZERO = "zero" - A PetscLimiter object
442*c0ce001eSMatthew G. Knepley 
443*c0ce001eSMatthew G. Knepley   Level: intermediate
444*c0ce001eSMatthew G. Knepley 
445*c0ce001eSMatthew G. Knepley .seealso: PetscLimiterType, PetscLimiterCreate(), PetscLimiterSetType()
446*c0ce001eSMatthew G. Knepley M*/
447*c0ce001eSMatthew G. Knepley 
448*c0ce001eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscLimiterCreate_Zero(PetscLimiter lim)
449*c0ce001eSMatthew G. Knepley {
450*c0ce001eSMatthew G. Knepley   PetscLimiter_Zero *l;
451*c0ce001eSMatthew G. Knepley   PetscErrorCode     ierr;
452*c0ce001eSMatthew G. Knepley 
453*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
454*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1);
455*c0ce001eSMatthew G. Knepley   ierr      = PetscNewLog(lim, &l);CHKERRQ(ierr);
456*c0ce001eSMatthew G. Knepley   lim->data = l;
457*c0ce001eSMatthew G. Knepley 
458*c0ce001eSMatthew G. Knepley   ierr = PetscLimiterInitialize_Zero(lim);CHKERRQ(ierr);
459*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
460*c0ce001eSMatthew G. Knepley }
461*c0ce001eSMatthew G. Knepley 
462*c0ce001eSMatthew G. Knepley PetscErrorCode PetscLimiterDestroy_None(PetscLimiter lim)
463*c0ce001eSMatthew G. Knepley {
464*c0ce001eSMatthew G. Knepley   PetscLimiter_None *l = (PetscLimiter_None *) lim->data;
465*c0ce001eSMatthew G. Knepley   PetscErrorCode    ierr;
466*c0ce001eSMatthew G. Knepley 
467*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
468*c0ce001eSMatthew G. Knepley   ierr = PetscFree(l);CHKERRQ(ierr);
469*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
470*c0ce001eSMatthew G. Knepley }
471*c0ce001eSMatthew G. Knepley 
472*c0ce001eSMatthew G. Knepley PetscErrorCode PetscLimiterView_None_Ascii(PetscLimiter lim, PetscViewer viewer)
473*c0ce001eSMatthew G. Knepley {
474*c0ce001eSMatthew G. Knepley   PetscViewerFormat format;
475*c0ce001eSMatthew G. Knepley   PetscErrorCode    ierr;
476*c0ce001eSMatthew G. Knepley 
477*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
478*c0ce001eSMatthew G. Knepley   ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr);
479*c0ce001eSMatthew G. Knepley   ierr = PetscViewerASCIIPrintf(viewer, "None Slope Limiter:\n");CHKERRQ(ierr);
480*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
481*c0ce001eSMatthew G. Knepley }
482*c0ce001eSMatthew G. Knepley 
483*c0ce001eSMatthew G. Knepley PetscErrorCode PetscLimiterView_None(PetscLimiter lim, PetscViewer viewer)
484*c0ce001eSMatthew G. Knepley {
485*c0ce001eSMatthew G. Knepley   PetscBool      iascii;
486*c0ce001eSMatthew G. Knepley   PetscErrorCode ierr;
487*c0ce001eSMatthew G. Knepley 
488*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
489*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1);
490*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
491*c0ce001eSMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
492*c0ce001eSMatthew G. Knepley   if (iascii) {ierr = PetscLimiterView_None_Ascii(lim, viewer);CHKERRQ(ierr);}
493*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
494*c0ce001eSMatthew G. Knepley }
495*c0ce001eSMatthew G. Knepley 
496*c0ce001eSMatthew G. Knepley PetscErrorCode PetscLimiterLimit_None(PetscLimiter lim, PetscReal f, PetscReal *phi)
497*c0ce001eSMatthew G. Knepley {
498*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
499*c0ce001eSMatthew G. Knepley   *phi = 1.0;
500*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
501*c0ce001eSMatthew G. Knepley }
502*c0ce001eSMatthew G. Knepley 
503*c0ce001eSMatthew G. Knepley PetscErrorCode PetscLimiterInitialize_None(PetscLimiter lim)
504*c0ce001eSMatthew G. Knepley {
505*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
506*c0ce001eSMatthew G. Knepley   lim->ops->view    = PetscLimiterView_None;
507*c0ce001eSMatthew G. Knepley   lim->ops->destroy = PetscLimiterDestroy_None;
508*c0ce001eSMatthew G. Knepley   lim->ops->limit   = PetscLimiterLimit_None;
509*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
510*c0ce001eSMatthew G. Knepley }
511*c0ce001eSMatthew G. Knepley 
512*c0ce001eSMatthew G. Knepley /*MC
513*c0ce001eSMatthew G. Knepley   PETSCLIMITERNONE = "none" - A PetscLimiter object
514*c0ce001eSMatthew G. Knepley 
515*c0ce001eSMatthew G. Knepley   Level: intermediate
516*c0ce001eSMatthew G. Knepley 
517*c0ce001eSMatthew G. Knepley .seealso: PetscLimiterType, PetscLimiterCreate(), PetscLimiterSetType()
518*c0ce001eSMatthew G. Knepley M*/
519*c0ce001eSMatthew G. Knepley 
520*c0ce001eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscLimiterCreate_None(PetscLimiter lim)
521*c0ce001eSMatthew G. Knepley {
522*c0ce001eSMatthew G. Knepley   PetscLimiter_None *l;
523*c0ce001eSMatthew G. Knepley   PetscErrorCode    ierr;
524*c0ce001eSMatthew G. Knepley 
525*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
526*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1);
527*c0ce001eSMatthew G. Knepley   ierr      = PetscNewLog(lim, &l);CHKERRQ(ierr);
528*c0ce001eSMatthew G. Knepley   lim->data = l;
529*c0ce001eSMatthew G. Knepley 
530*c0ce001eSMatthew G. Knepley   ierr = PetscLimiterInitialize_None(lim);CHKERRQ(ierr);
531*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
532*c0ce001eSMatthew G. Knepley }
533*c0ce001eSMatthew G. Knepley 
534*c0ce001eSMatthew G. Knepley PetscErrorCode PetscLimiterDestroy_Minmod(PetscLimiter lim)
535*c0ce001eSMatthew G. Knepley {
536*c0ce001eSMatthew G. Knepley   PetscLimiter_Minmod *l = (PetscLimiter_Minmod *) lim->data;
537*c0ce001eSMatthew G. Knepley   PetscErrorCode    ierr;
538*c0ce001eSMatthew G. Knepley 
539*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
540*c0ce001eSMatthew G. Knepley   ierr = PetscFree(l);CHKERRQ(ierr);
541*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
542*c0ce001eSMatthew G. Knepley }
543*c0ce001eSMatthew G. Knepley 
544*c0ce001eSMatthew G. Knepley PetscErrorCode PetscLimiterView_Minmod_Ascii(PetscLimiter lim, PetscViewer viewer)
545*c0ce001eSMatthew G. Knepley {
546*c0ce001eSMatthew G. Knepley   PetscViewerFormat format;
547*c0ce001eSMatthew G. Knepley   PetscErrorCode    ierr;
548*c0ce001eSMatthew G. Knepley 
549*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
550*c0ce001eSMatthew G. Knepley   ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr);
551*c0ce001eSMatthew G. Knepley   ierr = PetscViewerASCIIPrintf(viewer, "Minmod Slope Limiter:\n");CHKERRQ(ierr);
552*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
553*c0ce001eSMatthew G. Knepley }
554*c0ce001eSMatthew G. Knepley 
555*c0ce001eSMatthew G. Knepley PetscErrorCode PetscLimiterView_Minmod(PetscLimiter lim, PetscViewer viewer)
556*c0ce001eSMatthew G. Knepley {
557*c0ce001eSMatthew G. Knepley   PetscBool      iascii;
558*c0ce001eSMatthew G. Knepley   PetscErrorCode ierr;
559*c0ce001eSMatthew G. Knepley 
560*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
561*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1);
562*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
563*c0ce001eSMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
564*c0ce001eSMatthew G. Knepley   if (iascii) {ierr = PetscLimiterView_Minmod_Ascii(lim, viewer);CHKERRQ(ierr);}
565*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
566*c0ce001eSMatthew G. Knepley }
567*c0ce001eSMatthew G. Knepley 
568*c0ce001eSMatthew G. Knepley PetscErrorCode PetscLimiterLimit_Minmod(PetscLimiter lim, PetscReal f, PetscReal *phi)
569*c0ce001eSMatthew G. Knepley {
570*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
571*c0ce001eSMatthew G. Knepley   *phi = 2*PetscMax(0, PetscMin(f, 1-f));
572*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
573*c0ce001eSMatthew G. Knepley }
574*c0ce001eSMatthew G. Knepley 
575*c0ce001eSMatthew G. Knepley PetscErrorCode PetscLimiterInitialize_Minmod(PetscLimiter lim)
576*c0ce001eSMatthew G. Knepley {
577*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
578*c0ce001eSMatthew G. Knepley   lim->ops->view    = PetscLimiterView_Minmod;
579*c0ce001eSMatthew G. Knepley   lim->ops->destroy = PetscLimiterDestroy_Minmod;
580*c0ce001eSMatthew G. Knepley   lim->ops->limit   = PetscLimiterLimit_Minmod;
581*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
582*c0ce001eSMatthew G. Knepley }
583*c0ce001eSMatthew G. Knepley 
584*c0ce001eSMatthew G. Knepley /*MC
585*c0ce001eSMatthew G. Knepley   PETSCLIMITERMINMOD = "minmod" - A PetscLimiter object
586*c0ce001eSMatthew G. Knepley 
587*c0ce001eSMatthew G. Knepley   Level: intermediate
588*c0ce001eSMatthew G. Knepley 
589*c0ce001eSMatthew G. Knepley .seealso: PetscLimiterType, PetscLimiterCreate(), PetscLimiterSetType()
590*c0ce001eSMatthew G. Knepley M*/
591*c0ce001eSMatthew G. Knepley 
592*c0ce001eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscLimiterCreate_Minmod(PetscLimiter lim)
593*c0ce001eSMatthew G. Knepley {
594*c0ce001eSMatthew G. Knepley   PetscLimiter_Minmod *l;
595*c0ce001eSMatthew G. Knepley   PetscErrorCode    ierr;
596*c0ce001eSMatthew G. Knepley 
597*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
598*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1);
599*c0ce001eSMatthew G. Knepley   ierr      = PetscNewLog(lim, &l);CHKERRQ(ierr);
600*c0ce001eSMatthew G. Knepley   lim->data = l;
601*c0ce001eSMatthew G. Knepley 
602*c0ce001eSMatthew G. Knepley   ierr = PetscLimiterInitialize_Minmod(lim);CHKERRQ(ierr);
603*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
604*c0ce001eSMatthew G. Knepley }
605*c0ce001eSMatthew G. Knepley 
606*c0ce001eSMatthew G. Knepley PetscErrorCode PetscLimiterDestroy_VanLeer(PetscLimiter lim)
607*c0ce001eSMatthew G. Knepley {
608*c0ce001eSMatthew G. Knepley   PetscLimiter_VanLeer *l = (PetscLimiter_VanLeer *) lim->data;
609*c0ce001eSMatthew G. Knepley   PetscErrorCode    ierr;
610*c0ce001eSMatthew G. Knepley 
611*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
612*c0ce001eSMatthew G. Knepley   ierr = PetscFree(l);CHKERRQ(ierr);
613*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
614*c0ce001eSMatthew G. Knepley }
615*c0ce001eSMatthew G. Knepley 
616*c0ce001eSMatthew G. Knepley PetscErrorCode PetscLimiterView_VanLeer_Ascii(PetscLimiter lim, PetscViewer viewer)
617*c0ce001eSMatthew G. Knepley {
618*c0ce001eSMatthew G. Knepley   PetscViewerFormat format;
619*c0ce001eSMatthew G. Knepley   PetscErrorCode    ierr;
620*c0ce001eSMatthew G. Knepley 
621*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
622*c0ce001eSMatthew G. Knepley   ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr);
623*c0ce001eSMatthew G. Knepley   ierr = PetscViewerASCIIPrintf(viewer, "Van Leer Slope Limiter:\n");CHKERRQ(ierr);
624*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
625*c0ce001eSMatthew G. Knepley }
626*c0ce001eSMatthew G. Knepley 
627*c0ce001eSMatthew G. Knepley PetscErrorCode PetscLimiterView_VanLeer(PetscLimiter lim, PetscViewer viewer)
628*c0ce001eSMatthew G. Knepley {
629*c0ce001eSMatthew G. Knepley   PetscBool      iascii;
630*c0ce001eSMatthew G. Knepley   PetscErrorCode ierr;
631*c0ce001eSMatthew G. Knepley 
632*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
633*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1);
634*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
635*c0ce001eSMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
636*c0ce001eSMatthew G. Knepley   if (iascii) {ierr = PetscLimiterView_VanLeer_Ascii(lim, viewer);CHKERRQ(ierr);}
637*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
638*c0ce001eSMatthew G. Knepley }
639*c0ce001eSMatthew G. Knepley 
640*c0ce001eSMatthew G. Knepley PetscErrorCode PetscLimiterLimit_VanLeer(PetscLimiter lim, PetscReal f, PetscReal *phi)
641*c0ce001eSMatthew G. Knepley {
642*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
643*c0ce001eSMatthew G. Knepley   *phi = PetscMax(0, 4*f*(1-f));
644*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
645*c0ce001eSMatthew G. Knepley }
646*c0ce001eSMatthew G. Knepley 
647*c0ce001eSMatthew G. Knepley PetscErrorCode PetscLimiterInitialize_VanLeer(PetscLimiter lim)
648*c0ce001eSMatthew G. Knepley {
649*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
650*c0ce001eSMatthew G. Knepley   lim->ops->view    = PetscLimiterView_VanLeer;
651*c0ce001eSMatthew G. Knepley   lim->ops->destroy = PetscLimiterDestroy_VanLeer;
652*c0ce001eSMatthew G. Knepley   lim->ops->limit   = PetscLimiterLimit_VanLeer;
653*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
654*c0ce001eSMatthew G. Knepley }
655*c0ce001eSMatthew G. Knepley 
656*c0ce001eSMatthew G. Knepley /*MC
657*c0ce001eSMatthew G. Knepley   PETSCLIMITERVANLEER = "vanleer" - A PetscLimiter object
658*c0ce001eSMatthew G. Knepley 
659*c0ce001eSMatthew G. Knepley   Level: intermediate
660*c0ce001eSMatthew G. Knepley 
661*c0ce001eSMatthew G. Knepley .seealso: PetscLimiterType, PetscLimiterCreate(), PetscLimiterSetType()
662*c0ce001eSMatthew G. Knepley M*/
663*c0ce001eSMatthew G. Knepley 
664*c0ce001eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscLimiterCreate_VanLeer(PetscLimiter lim)
665*c0ce001eSMatthew G. Knepley {
666*c0ce001eSMatthew G. Knepley   PetscLimiter_VanLeer *l;
667*c0ce001eSMatthew G. Knepley   PetscErrorCode    ierr;
668*c0ce001eSMatthew G. Knepley 
669*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
670*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1);
671*c0ce001eSMatthew G. Knepley   ierr      = PetscNewLog(lim, &l);CHKERRQ(ierr);
672*c0ce001eSMatthew G. Knepley   lim->data = l;
673*c0ce001eSMatthew G. Knepley 
674*c0ce001eSMatthew G. Knepley   ierr = PetscLimiterInitialize_VanLeer(lim);CHKERRQ(ierr);
675*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
676*c0ce001eSMatthew G. Knepley }
677*c0ce001eSMatthew G. Knepley 
678*c0ce001eSMatthew G. Knepley PetscErrorCode PetscLimiterDestroy_VanAlbada(PetscLimiter lim)
679*c0ce001eSMatthew G. Knepley {
680*c0ce001eSMatthew G. Knepley   PetscLimiter_VanAlbada *l = (PetscLimiter_VanAlbada *) lim->data;
681*c0ce001eSMatthew G. Knepley   PetscErrorCode    ierr;
682*c0ce001eSMatthew G. Knepley 
683*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
684*c0ce001eSMatthew G. Knepley   ierr = PetscFree(l);CHKERRQ(ierr);
685*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
686*c0ce001eSMatthew G. Knepley }
687*c0ce001eSMatthew G. Knepley 
688*c0ce001eSMatthew G. Knepley PetscErrorCode PetscLimiterView_VanAlbada_Ascii(PetscLimiter lim, PetscViewer viewer)
689*c0ce001eSMatthew G. Knepley {
690*c0ce001eSMatthew G. Knepley   PetscViewerFormat format;
691*c0ce001eSMatthew G. Knepley   PetscErrorCode    ierr;
692*c0ce001eSMatthew G. Knepley 
693*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
694*c0ce001eSMatthew G. Knepley   ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr);
695*c0ce001eSMatthew G. Knepley   ierr = PetscViewerASCIIPrintf(viewer, "Van Albada Slope Limiter:\n");CHKERRQ(ierr);
696*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
697*c0ce001eSMatthew G. Knepley }
698*c0ce001eSMatthew G. Knepley 
699*c0ce001eSMatthew G. Knepley PetscErrorCode PetscLimiterView_VanAlbada(PetscLimiter lim, PetscViewer viewer)
700*c0ce001eSMatthew G. Knepley {
701*c0ce001eSMatthew G. Knepley   PetscBool      iascii;
702*c0ce001eSMatthew G. Knepley   PetscErrorCode ierr;
703*c0ce001eSMatthew G. Knepley 
704*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
705*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1);
706*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
707*c0ce001eSMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
708*c0ce001eSMatthew G. Knepley   if (iascii) {ierr = PetscLimiterView_VanAlbada_Ascii(lim, viewer);CHKERRQ(ierr);}
709*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
710*c0ce001eSMatthew G. Knepley }
711*c0ce001eSMatthew G. Knepley 
712*c0ce001eSMatthew G. Knepley PetscErrorCode PetscLimiterLimit_VanAlbada(PetscLimiter lim, PetscReal f, PetscReal *phi)
713*c0ce001eSMatthew G. Knepley {
714*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
715*c0ce001eSMatthew G. Knepley   *phi = PetscMax(0, 2*f*(1-f) / (PetscSqr(f) + PetscSqr(1-f)));
716*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
717*c0ce001eSMatthew G. Knepley }
718*c0ce001eSMatthew G. Knepley 
719*c0ce001eSMatthew G. Knepley PetscErrorCode PetscLimiterInitialize_VanAlbada(PetscLimiter lim)
720*c0ce001eSMatthew G. Knepley {
721*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
722*c0ce001eSMatthew G. Knepley   lim->ops->view    = PetscLimiterView_VanAlbada;
723*c0ce001eSMatthew G. Knepley   lim->ops->destroy = PetscLimiterDestroy_VanAlbada;
724*c0ce001eSMatthew G. Knepley   lim->ops->limit   = PetscLimiterLimit_VanAlbada;
725*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
726*c0ce001eSMatthew G. Knepley }
727*c0ce001eSMatthew G. Knepley 
728*c0ce001eSMatthew G. Knepley /*MC
729*c0ce001eSMatthew G. Knepley   PETSCLIMITERVANALBADA = "vanalbada" - A PetscLimiter object
730*c0ce001eSMatthew G. Knepley 
731*c0ce001eSMatthew G. Knepley   Level: intermediate
732*c0ce001eSMatthew G. Knepley 
733*c0ce001eSMatthew G. Knepley .seealso: PetscLimiterType, PetscLimiterCreate(), PetscLimiterSetType()
734*c0ce001eSMatthew G. Knepley M*/
735*c0ce001eSMatthew G. Knepley 
736*c0ce001eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscLimiterCreate_VanAlbada(PetscLimiter lim)
737*c0ce001eSMatthew G. Knepley {
738*c0ce001eSMatthew G. Knepley   PetscLimiter_VanAlbada *l;
739*c0ce001eSMatthew G. Knepley   PetscErrorCode    ierr;
740*c0ce001eSMatthew G. Knepley 
741*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
742*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1);
743*c0ce001eSMatthew G. Knepley   ierr      = PetscNewLog(lim, &l);CHKERRQ(ierr);
744*c0ce001eSMatthew G. Knepley   lim->data = l;
745*c0ce001eSMatthew G. Knepley 
746*c0ce001eSMatthew G. Knepley   ierr = PetscLimiterInitialize_VanAlbada(lim);CHKERRQ(ierr);
747*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
748*c0ce001eSMatthew G. Knepley }
749*c0ce001eSMatthew G. Knepley 
750*c0ce001eSMatthew G. Knepley PetscErrorCode PetscLimiterDestroy_Superbee(PetscLimiter lim)
751*c0ce001eSMatthew G. Knepley {
752*c0ce001eSMatthew G. Knepley   PetscLimiter_Superbee *l = (PetscLimiter_Superbee *) lim->data;
753*c0ce001eSMatthew G. Knepley   PetscErrorCode    ierr;
754*c0ce001eSMatthew G. Knepley 
755*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
756*c0ce001eSMatthew G. Knepley   ierr = PetscFree(l);CHKERRQ(ierr);
757*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
758*c0ce001eSMatthew G. Knepley }
759*c0ce001eSMatthew G. Knepley 
760*c0ce001eSMatthew G. Knepley PetscErrorCode PetscLimiterView_Superbee_Ascii(PetscLimiter lim, PetscViewer viewer)
761*c0ce001eSMatthew G. Knepley {
762*c0ce001eSMatthew G. Knepley   PetscViewerFormat format;
763*c0ce001eSMatthew G. Knepley   PetscErrorCode    ierr;
764*c0ce001eSMatthew G. Knepley 
765*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
766*c0ce001eSMatthew G. Knepley   ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr);
767*c0ce001eSMatthew G. Knepley   ierr = PetscViewerASCIIPrintf(viewer, "Superbee Slope Limiter:\n");CHKERRQ(ierr);
768*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
769*c0ce001eSMatthew G. Knepley }
770*c0ce001eSMatthew G. Knepley 
771*c0ce001eSMatthew G. Knepley PetscErrorCode PetscLimiterView_Superbee(PetscLimiter lim, PetscViewer viewer)
772*c0ce001eSMatthew G. Knepley {
773*c0ce001eSMatthew G. Knepley   PetscBool      iascii;
774*c0ce001eSMatthew G. Knepley   PetscErrorCode ierr;
775*c0ce001eSMatthew G. Knepley 
776*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
777*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1);
778*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
779*c0ce001eSMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
780*c0ce001eSMatthew G. Knepley   if (iascii) {ierr = PetscLimiterView_Superbee_Ascii(lim, viewer);CHKERRQ(ierr);}
781*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
782*c0ce001eSMatthew G. Knepley }
783*c0ce001eSMatthew G. Knepley 
784*c0ce001eSMatthew G. Knepley PetscErrorCode PetscLimiterLimit_Superbee(PetscLimiter lim, PetscReal f, PetscReal *phi)
785*c0ce001eSMatthew G. Knepley {
786*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
787*c0ce001eSMatthew G. Knepley   *phi = 4*PetscMax(0, PetscMin(f, 1-f));
788*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
789*c0ce001eSMatthew G. Knepley }
790*c0ce001eSMatthew G. Knepley 
791*c0ce001eSMatthew G. Knepley PetscErrorCode PetscLimiterInitialize_Superbee(PetscLimiter lim)
792*c0ce001eSMatthew G. Knepley {
793*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
794*c0ce001eSMatthew G. Knepley   lim->ops->view    = PetscLimiterView_Superbee;
795*c0ce001eSMatthew G. Knepley   lim->ops->destroy = PetscLimiterDestroy_Superbee;
796*c0ce001eSMatthew G. Knepley   lim->ops->limit   = PetscLimiterLimit_Superbee;
797*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
798*c0ce001eSMatthew G. Knepley }
799*c0ce001eSMatthew G. Knepley 
800*c0ce001eSMatthew G. Knepley /*MC
801*c0ce001eSMatthew G. Knepley   PETSCLIMITERSUPERBEE = "superbee" - A PetscLimiter object
802*c0ce001eSMatthew G. Knepley 
803*c0ce001eSMatthew G. Knepley   Level: intermediate
804*c0ce001eSMatthew G. Knepley 
805*c0ce001eSMatthew G. Knepley .seealso: PetscLimiterType, PetscLimiterCreate(), PetscLimiterSetType()
806*c0ce001eSMatthew G. Knepley M*/
807*c0ce001eSMatthew G. Knepley 
808*c0ce001eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscLimiterCreate_Superbee(PetscLimiter lim)
809*c0ce001eSMatthew G. Knepley {
810*c0ce001eSMatthew G. Knepley   PetscLimiter_Superbee *l;
811*c0ce001eSMatthew G. Knepley   PetscErrorCode    ierr;
812*c0ce001eSMatthew G. Knepley 
813*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
814*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1);
815*c0ce001eSMatthew G. Knepley   ierr      = PetscNewLog(lim, &l);CHKERRQ(ierr);
816*c0ce001eSMatthew G. Knepley   lim->data = l;
817*c0ce001eSMatthew G. Knepley 
818*c0ce001eSMatthew G. Knepley   ierr = PetscLimiterInitialize_Superbee(lim);CHKERRQ(ierr);
819*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
820*c0ce001eSMatthew G. Knepley }
821*c0ce001eSMatthew G. Knepley 
822*c0ce001eSMatthew G. Knepley PetscErrorCode PetscLimiterDestroy_MC(PetscLimiter lim)
823*c0ce001eSMatthew G. Knepley {
824*c0ce001eSMatthew G. Knepley   PetscLimiter_MC *l = (PetscLimiter_MC *) lim->data;
825*c0ce001eSMatthew G. Knepley   PetscErrorCode    ierr;
826*c0ce001eSMatthew G. Knepley 
827*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
828*c0ce001eSMatthew G. Knepley   ierr = PetscFree(l);CHKERRQ(ierr);
829*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
830*c0ce001eSMatthew G. Knepley }
831*c0ce001eSMatthew G. Knepley 
832*c0ce001eSMatthew G. Knepley PetscErrorCode PetscLimiterView_MC_Ascii(PetscLimiter lim, PetscViewer viewer)
833*c0ce001eSMatthew G. Knepley {
834*c0ce001eSMatthew G. Knepley   PetscViewerFormat format;
835*c0ce001eSMatthew G. Knepley   PetscErrorCode    ierr;
836*c0ce001eSMatthew G. Knepley 
837*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
838*c0ce001eSMatthew G. Knepley   ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr);
839*c0ce001eSMatthew G. Knepley   ierr = PetscViewerASCIIPrintf(viewer, "MC Slope Limiter:\n");CHKERRQ(ierr);
840*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
841*c0ce001eSMatthew G. Knepley }
842*c0ce001eSMatthew G. Knepley 
843*c0ce001eSMatthew G. Knepley PetscErrorCode PetscLimiterView_MC(PetscLimiter lim, PetscViewer viewer)
844*c0ce001eSMatthew G. Knepley {
845*c0ce001eSMatthew G. Knepley   PetscBool      iascii;
846*c0ce001eSMatthew G. Knepley   PetscErrorCode ierr;
847*c0ce001eSMatthew G. Knepley 
848*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
849*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1);
850*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
851*c0ce001eSMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
852*c0ce001eSMatthew G. Knepley   if (iascii) {ierr = PetscLimiterView_MC_Ascii(lim, viewer);CHKERRQ(ierr);}
853*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
854*c0ce001eSMatthew G. Knepley }
855*c0ce001eSMatthew G. Knepley 
856*c0ce001eSMatthew G. Knepley /* aka Barth-Jespersen */
857*c0ce001eSMatthew G. Knepley PetscErrorCode PetscLimiterLimit_MC(PetscLimiter lim, PetscReal f, PetscReal *phi)
858*c0ce001eSMatthew G. Knepley {
859*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
860*c0ce001eSMatthew G. Knepley   *phi = PetscMin(1, 4*PetscMax(0, PetscMin(f, 1-f)));
861*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
862*c0ce001eSMatthew G. Knepley }
863*c0ce001eSMatthew G. Knepley 
864*c0ce001eSMatthew G. Knepley PetscErrorCode PetscLimiterInitialize_MC(PetscLimiter lim)
865*c0ce001eSMatthew G. Knepley {
866*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
867*c0ce001eSMatthew G. Knepley   lim->ops->view    = PetscLimiterView_MC;
868*c0ce001eSMatthew G. Knepley   lim->ops->destroy = PetscLimiterDestroy_MC;
869*c0ce001eSMatthew G. Knepley   lim->ops->limit   = PetscLimiterLimit_MC;
870*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
871*c0ce001eSMatthew G. Knepley }
872*c0ce001eSMatthew G. Knepley 
873*c0ce001eSMatthew G. Knepley /*MC
874*c0ce001eSMatthew G. Knepley   PETSCLIMITERMC = "mc" - A PetscLimiter object
875*c0ce001eSMatthew G. Knepley 
876*c0ce001eSMatthew G. Knepley   Level: intermediate
877*c0ce001eSMatthew G. Knepley 
878*c0ce001eSMatthew G. Knepley .seealso: PetscLimiterType, PetscLimiterCreate(), PetscLimiterSetType()
879*c0ce001eSMatthew G. Knepley M*/
880*c0ce001eSMatthew G. Knepley 
881*c0ce001eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscLimiterCreate_MC(PetscLimiter lim)
882*c0ce001eSMatthew G. Knepley {
883*c0ce001eSMatthew G. Knepley   PetscLimiter_MC *l;
884*c0ce001eSMatthew G. Knepley   PetscErrorCode    ierr;
885*c0ce001eSMatthew G. Knepley 
886*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
887*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1);
888*c0ce001eSMatthew G. Knepley   ierr      = PetscNewLog(lim, &l);CHKERRQ(ierr);
889*c0ce001eSMatthew G. Knepley   lim->data = l;
890*c0ce001eSMatthew G. Knepley 
891*c0ce001eSMatthew G. Knepley   ierr = PetscLimiterInitialize_MC(lim);CHKERRQ(ierr);
892*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
893*c0ce001eSMatthew G. Knepley }
894*c0ce001eSMatthew G. Knepley 
895*c0ce001eSMatthew G. Knepley PetscClassId PETSCFV_CLASSID = 0;
896*c0ce001eSMatthew G. Knepley 
897*c0ce001eSMatthew G. Knepley PetscFunctionList PetscFVList              = NULL;
898*c0ce001eSMatthew G. Knepley PetscBool         PetscFVRegisterAllCalled = PETSC_FALSE;
899*c0ce001eSMatthew G. Knepley 
900*c0ce001eSMatthew G. Knepley /*@C
901*c0ce001eSMatthew G. Knepley   PetscFVRegister - Adds a new PetscFV implementation
902*c0ce001eSMatthew G. Knepley 
903*c0ce001eSMatthew G. Knepley   Not Collective
904*c0ce001eSMatthew G. Knepley 
905*c0ce001eSMatthew G. Knepley   Input Parameters:
906*c0ce001eSMatthew G. Knepley + name        - The name of a new user-defined creation routine
907*c0ce001eSMatthew G. Knepley - create_func - The creation routine itself
908*c0ce001eSMatthew G. Knepley 
909*c0ce001eSMatthew G. Knepley   Notes:
910*c0ce001eSMatthew G. Knepley   PetscFVRegister() may be called multiple times to add several user-defined PetscFVs
911*c0ce001eSMatthew G. Knepley 
912*c0ce001eSMatthew G. Knepley   Sample usage:
913*c0ce001eSMatthew G. Knepley .vb
914*c0ce001eSMatthew G. Knepley     PetscFVRegister("my_fv", MyPetscFVCreate);
915*c0ce001eSMatthew G. Knepley .ve
916*c0ce001eSMatthew G. Knepley 
917*c0ce001eSMatthew G. Knepley   Then, your PetscFV type can be chosen with the procedural interface via
918*c0ce001eSMatthew G. Knepley .vb
919*c0ce001eSMatthew G. Knepley     PetscFVCreate(MPI_Comm, PetscFV *);
920*c0ce001eSMatthew G. Knepley     PetscFVSetType(PetscFV, "my_fv");
921*c0ce001eSMatthew G. Knepley .ve
922*c0ce001eSMatthew G. Knepley    or at runtime via the option
923*c0ce001eSMatthew G. Knepley .vb
924*c0ce001eSMatthew G. Knepley     -petscfv_type my_fv
925*c0ce001eSMatthew G. Knepley .ve
926*c0ce001eSMatthew G. Knepley 
927*c0ce001eSMatthew G. Knepley   Level: advanced
928*c0ce001eSMatthew G. Knepley 
929*c0ce001eSMatthew G. Knepley .seealso: PetscFVRegisterAll(), PetscFVRegisterDestroy()
930*c0ce001eSMatthew G. Knepley 
931*c0ce001eSMatthew G. Knepley @*/
932*c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVRegister(const char sname[], PetscErrorCode (*function)(PetscFV))
933*c0ce001eSMatthew G. Knepley {
934*c0ce001eSMatthew G. Knepley   PetscErrorCode ierr;
935*c0ce001eSMatthew G. Knepley 
936*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
937*c0ce001eSMatthew G. Knepley   ierr = PetscFunctionListAdd(&PetscFVList, sname, function);CHKERRQ(ierr);
938*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
939*c0ce001eSMatthew G. Knepley }
940*c0ce001eSMatthew G. Knepley 
941*c0ce001eSMatthew G. Knepley /*@C
942*c0ce001eSMatthew G. Knepley   PetscFVSetType - Builds a particular PetscFV
943*c0ce001eSMatthew G. Knepley 
944*c0ce001eSMatthew G. Knepley   Collective on fvm
945*c0ce001eSMatthew G. Knepley 
946*c0ce001eSMatthew G. Knepley   Input Parameters:
947*c0ce001eSMatthew G. Knepley + fvm  - The PetscFV object
948*c0ce001eSMatthew G. Knepley - name - The kind of FVM space
949*c0ce001eSMatthew G. Knepley 
950*c0ce001eSMatthew G. Knepley   Options Database Key:
951*c0ce001eSMatthew G. Knepley . -petscfv_type <type> - Sets the PetscFV type; use -help for a list of available types
952*c0ce001eSMatthew G. Knepley 
953*c0ce001eSMatthew G. Knepley   Level: intermediate
954*c0ce001eSMatthew G. Knepley 
955*c0ce001eSMatthew G. Knepley .seealso: PetscFVGetType(), PetscFVCreate()
956*c0ce001eSMatthew G. Knepley @*/
957*c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVSetType(PetscFV fvm, PetscFVType name)
958*c0ce001eSMatthew G. Knepley {
959*c0ce001eSMatthew G. Knepley   PetscErrorCode (*r)(PetscFV);
960*c0ce001eSMatthew G. Knepley   PetscBool      match;
961*c0ce001eSMatthew G. Knepley   PetscErrorCode ierr;
962*c0ce001eSMatthew G. Knepley 
963*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
964*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1);
965*c0ce001eSMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject) fvm, name, &match);CHKERRQ(ierr);
966*c0ce001eSMatthew G. Knepley   if (match) PetscFunctionReturn(0);
967*c0ce001eSMatthew G. Knepley 
968*c0ce001eSMatthew G. Knepley   ierr = PetscFVRegisterAll();CHKERRQ(ierr);
969*c0ce001eSMatthew G. Knepley   ierr = PetscFunctionListFind(PetscFVList, name, &r);CHKERRQ(ierr);
970*c0ce001eSMatthew G. Knepley   if (!r) SETERRQ1(PetscObjectComm((PetscObject) fvm), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscFV type: %s", name);
971*c0ce001eSMatthew G. Knepley 
972*c0ce001eSMatthew G. Knepley   if (fvm->ops->destroy) {
973*c0ce001eSMatthew G. Knepley     ierr              = (*fvm->ops->destroy)(fvm);CHKERRQ(ierr);
974*c0ce001eSMatthew G. Knepley     fvm->ops->destroy = NULL;
975*c0ce001eSMatthew G. Knepley   }
976*c0ce001eSMatthew G. Knepley   ierr = (*r)(fvm);CHKERRQ(ierr);
977*c0ce001eSMatthew G. Knepley   ierr = PetscObjectChangeTypeName((PetscObject) fvm, name);CHKERRQ(ierr);
978*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
979*c0ce001eSMatthew G. Knepley }
980*c0ce001eSMatthew G. Knepley 
981*c0ce001eSMatthew G. Knepley /*@C
982*c0ce001eSMatthew G. Knepley   PetscFVGetType - Gets the PetscFV type name (as a string) from the object.
983*c0ce001eSMatthew G. Knepley 
984*c0ce001eSMatthew G. Knepley   Not Collective
985*c0ce001eSMatthew G. Knepley 
986*c0ce001eSMatthew G. Knepley   Input Parameter:
987*c0ce001eSMatthew G. Knepley . fvm  - The PetscFV
988*c0ce001eSMatthew G. Knepley 
989*c0ce001eSMatthew G. Knepley   Output Parameter:
990*c0ce001eSMatthew G. Knepley . name - The PetscFV type name
991*c0ce001eSMatthew G. Knepley 
992*c0ce001eSMatthew G. Knepley   Level: intermediate
993*c0ce001eSMatthew G. Knepley 
994*c0ce001eSMatthew G. Knepley .seealso: PetscFVSetType(), PetscFVCreate()
995*c0ce001eSMatthew G. Knepley @*/
996*c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVGetType(PetscFV fvm, PetscFVType *name)
997*c0ce001eSMatthew G. Knepley {
998*c0ce001eSMatthew G. Knepley   PetscErrorCode ierr;
999*c0ce001eSMatthew G. Knepley 
1000*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1001*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1);
1002*c0ce001eSMatthew G. Knepley   PetscValidPointer(name, 2);
1003*c0ce001eSMatthew G. Knepley   ierr = PetscFVRegisterAll();CHKERRQ(ierr);
1004*c0ce001eSMatthew G. Knepley   *name = ((PetscObject) fvm)->type_name;
1005*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
1006*c0ce001eSMatthew G. Knepley }
1007*c0ce001eSMatthew G. Knepley 
1008*c0ce001eSMatthew G. Knepley /*@C
1009*c0ce001eSMatthew G. Knepley   PetscFVView - Views a PetscFV
1010*c0ce001eSMatthew G. Knepley 
1011*c0ce001eSMatthew G. Knepley   Collective on fvm
1012*c0ce001eSMatthew G. Knepley 
1013*c0ce001eSMatthew G. Knepley   Input Parameter:
1014*c0ce001eSMatthew G. Knepley + fvm - the PetscFV object to view
1015*c0ce001eSMatthew G. Knepley - v   - the viewer
1016*c0ce001eSMatthew G. Knepley 
1017*c0ce001eSMatthew G. Knepley   Level: developer
1018*c0ce001eSMatthew G. Knepley 
1019*c0ce001eSMatthew G. Knepley .seealso: PetscFVDestroy()
1020*c0ce001eSMatthew G. Knepley @*/
1021*c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVView(PetscFV fvm, PetscViewer v)
1022*c0ce001eSMatthew G. Knepley {
1023*c0ce001eSMatthew G. Knepley   PetscErrorCode ierr;
1024*c0ce001eSMatthew G. Knepley 
1025*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1026*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1);
1027*c0ce001eSMatthew G. Knepley   if (!v) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) fvm), &v);CHKERRQ(ierr);}
1028*c0ce001eSMatthew G. Knepley   if (fvm->ops->view) {ierr = (*fvm->ops->view)(fvm, v);CHKERRQ(ierr);}
1029*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
1030*c0ce001eSMatthew G. Knepley }
1031*c0ce001eSMatthew G. Knepley 
1032*c0ce001eSMatthew G. Knepley /*@
1033*c0ce001eSMatthew G. Knepley   PetscFVSetFromOptions - sets parameters in a PetscFV from the options database
1034*c0ce001eSMatthew G. Knepley 
1035*c0ce001eSMatthew G. Knepley   Collective on fvm
1036*c0ce001eSMatthew G. Knepley 
1037*c0ce001eSMatthew G. Knepley   Input Parameter:
1038*c0ce001eSMatthew G. Knepley . fvm - the PetscFV object to set options for
1039*c0ce001eSMatthew G. Knepley 
1040*c0ce001eSMatthew G. Knepley   Options Database Key:
1041*c0ce001eSMatthew G. Knepley . -petscfv_compute_gradients <bool> - Determines whether cell gradients are calculated
1042*c0ce001eSMatthew G. Knepley 
1043*c0ce001eSMatthew G. Knepley   Level: developer
1044*c0ce001eSMatthew G. Knepley 
1045*c0ce001eSMatthew G. Knepley .seealso: PetscFVView()
1046*c0ce001eSMatthew G. Knepley @*/
1047*c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVSetFromOptions(PetscFV fvm)
1048*c0ce001eSMatthew G. Knepley {
1049*c0ce001eSMatthew G. Knepley   const char    *defaultType;
1050*c0ce001eSMatthew G. Knepley   char           name[256];
1051*c0ce001eSMatthew G. Knepley   PetscBool      flg;
1052*c0ce001eSMatthew G. Knepley   PetscErrorCode ierr;
1053*c0ce001eSMatthew G. Knepley 
1054*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1055*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1);
1056*c0ce001eSMatthew G. Knepley   if (!((PetscObject) fvm)->type_name) defaultType = PETSCFVUPWIND;
1057*c0ce001eSMatthew G. Knepley   else                                 defaultType = ((PetscObject) fvm)->type_name;
1058*c0ce001eSMatthew G. Knepley   ierr = PetscFVRegisterAll();CHKERRQ(ierr);
1059*c0ce001eSMatthew G. Knepley 
1060*c0ce001eSMatthew G. Knepley   ierr = PetscObjectOptionsBegin((PetscObject) fvm);CHKERRQ(ierr);
1061*c0ce001eSMatthew G. Knepley   ierr = PetscOptionsFList("-petscfv_type", "Finite volume discretization", "PetscFVSetType", PetscFVList, defaultType, name, 256, &flg);CHKERRQ(ierr);
1062*c0ce001eSMatthew G. Knepley   if (flg) {
1063*c0ce001eSMatthew G. Knepley     ierr = PetscFVSetType(fvm, name);CHKERRQ(ierr);
1064*c0ce001eSMatthew G. Knepley   } else if (!((PetscObject) fvm)->type_name) {
1065*c0ce001eSMatthew G. Knepley     ierr = PetscFVSetType(fvm, defaultType);CHKERRQ(ierr);
1066*c0ce001eSMatthew G. Knepley 
1067*c0ce001eSMatthew G. Knepley   }
1068*c0ce001eSMatthew G. Knepley   ierr = PetscOptionsBool("-petscfv_compute_gradients", "Compute cell gradients", "PetscFVSetComputeGradients", fvm->computeGradients, &fvm->computeGradients, NULL);CHKERRQ(ierr);
1069*c0ce001eSMatthew G. Knepley   if (fvm->ops->setfromoptions) {ierr = (*fvm->ops->setfromoptions)(fvm);CHKERRQ(ierr);}
1070*c0ce001eSMatthew G. Knepley   /* process any options handlers added with PetscObjectAddOptionsHandler() */
1071*c0ce001eSMatthew G. Knepley   ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject) fvm);CHKERRQ(ierr);
1072*c0ce001eSMatthew G. Knepley   ierr = PetscLimiterSetFromOptions(fvm->limiter);CHKERRQ(ierr);
1073*c0ce001eSMatthew G. Knepley   ierr = PetscOptionsEnd();CHKERRQ(ierr);
1074*c0ce001eSMatthew G. Knepley   ierr = PetscFVViewFromOptions(fvm, NULL, "-petscfv_view");CHKERRQ(ierr);
1075*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
1076*c0ce001eSMatthew G. Knepley }
1077*c0ce001eSMatthew G. Knepley 
1078*c0ce001eSMatthew G. Knepley /*@
1079*c0ce001eSMatthew G. Knepley   PetscFVSetUp - Construct data structures for the PetscFV
1080*c0ce001eSMatthew G. Knepley 
1081*c0ce001eSMatthew G. Knepley   Collective on fvm
1082*c0ce001eSMatthew G. Knepley 
1083*c0ce001eSMatthew G. Knepley   Input Parameter:
1084*c0ce001eSMatthew G. Knepley . fvm - the PetscFV object to setup
1085*c0ce001eSMatthew G. Knepley 
1086*c0ce001eSMatthew G. Knepley   Level: developer
1087*c0ce001eSMatthew G. Knepley 
1088*c0ce001eSMatthew G. Knepley .seealso: PetscFVView(), PetscFVDestroy()
1089*c0ce001eSMatthew G. Knepley @*/
1090*c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVSetUp(PetscFV fvm)
1091*c0ce001eSMatthew G. Knepley {
1092*c0ce001eSMatthew G. Knepley   PetscErrorCode ierr;
1093*c0ce001eSMatthew G. Knepley 
1094*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1095*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1);
1096*c0ce001eSMatthew G. Knepley   ierr = PetscLimiterSetUp(fvm->limiter);CHKERRQ(ierr);
1097*c0ce001eSMatthew G. Knepley   if (fvm->ops->setup) {ierr = (*fvm->ops->setup)(fvm);CHKERRQ(ierr);}
1098*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
1099*c0ce001eSMatthew G. Knepley }
1100*c0ce001eSMatthew G. Knepley 
1101*c0ce001eSMatthew G. Knepley /*@
1102*c0ce001eSMatthew G. Knepley   PetscFVDestroy - Destroys a PetscFV object
1103*c0ce001eSMatthew G. Knepley 
1104*c0ce001eSMatthew G. Knepley   Collective on fvm
1105*c0ce001eSMatthew G. Knepley 
1106*c0ce001eSMatthew G. Knepley   Input Parameter:
1107*c0ce001eSMatthew G. Knepley . fvm - the PetscFV object to destroy
1108*c0ce001eSMatthew G. Knepley 
1109*c0ce001eSMatthew G. Knepley   Level: developer
1110*c0ce001eSMatthew G. Knepley 
1111*c0ce001eSMatthew G. Knepley .seealso: PetscFVView()
1112*c0ce001eSMatthew G. Knepley @*/
1113*c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVDestroy(PetscFV *fvm)
1114*c0ce001eSMatthew G. Knepley {
1115*c0ce001eSMatthew G. Knepley   PetscInt       i;
1116*c0ce001eSMatthew G. Knepley   PetscErrorCode ierr;
1117*c0ce001eSMatthew G. Knepley 
1118*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1119*c0ce001eSMatthew G. Knepley   if (!*fvm) PetscFunctionReturn(0);
1120*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific((*fvm), PETSCFV_CLASSID, 1);
1121*c0ce001eSMatthew G. Knepley 
1122*c0ce001eSMatthew G. Knepley   if (--((PetscObject)(*fvm))->refct > 0) {*fvm = 0; PetscFunctionReturn(0);}
1123*c0ce001eSMatthew G. Knepley   ((PetscObject) (*fvm))->refct = 0;
1124*c0ce001eSMatthew G. Knepley 
1125*c0ce001eSMatthew G. Knepley   for (i = 0; i < (*fvm)->numComponents; i++) {
1126*c0ce001eSMatthew G. Knepley     ierr = PetscFree((*fvm)->componentNames[i]);CHKERRQ(ierr);
1127*c0ce001eSMatthew G. Knepley   }
1128*c0ce001eSMatthew G. Knepley   ierr = PetscFree((*fvm)->componentNames);CHKERRQ(ierr);
1129*c0ce001eSMatthew G. Knepley   ierr = PetscLimiterDestroy(&(*fvm)->limiter);CHKERRQ(ierr);
1130*c0ce001eSMatthew G. Knepley   ierr = PetscDualSpaceDestroy(&(*fvm)->dualSpace);CHKERRQ(ierr);
1131*c0ce001eSMatthew G. Knepley   ierr = PetscFree((*fvm)->fluxWork);CHKERRQ(ierr);
1132*c0ce001eSMatthew G. Knepley   ierr = PetscQuadratureDestroy(&(*fvm)->quadrature);CHKERRQ(ierr);
1133*c0ce001eSMatthew G. Knepley   ierr = PetscFVRestoreTabulation((*fvm), 0, NULL, &(*fvm)->B, &(*fvm)->D, NULL /*&(*fvm)->H*/);CHKERRQ(ierr);
1134*c0ce001eSMatthew G. Knepley 
1135*c0ce001eSMatthew G. Knepley   if ((*fvm)->ops->destroy) {ierr = (*(*fvm)->ops->destroy)(*fvm);CHKERRQ(ierr);}
1136*c0ce001eSMatthew G. Knepley   ierr = PetscHeaderDestroy(fvm);CHKERRQ(ierr);
1137*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
1138*c0ce001eSMatthew G. Knepley }
1139*c0ce001eSMatthew G. Knepley 
1140*c0ce001eSMatthew G. Knepley /*@
1141*c0ce001eSMatthew G. Knepley   PetscFVCreate - Creates an empty PetscFV object. The type can then be set with PetscFVSetType().
1142*c0ce001eSMatthew G. Knepley 
1143*c0ce001eSMatthew G. Knepley   Collective
1144*c0ce001eSMatthew G. Knepley 
1145*c0ce001eSMatthew G. Knepley   Input Parameter:
1146*c0ce001eSMatthew G. Knepley . comm - The communicator for the PetscFV object
1147*c0ce001eSMatthew G. Knepley 
1148*c0ce001eSMatthew G. Knepley   Output Parameter:
1149*c0ce001eSMatthew G. Knepley . fvm - The PetscFV object
1150*c0ce001eSMatthew G. Knepley 
1151*c0ce001eSMatthew G. Knepley   Level: beginner
1152*c0ce001eSMatthew G. Knepley 
1153*c0ce001eSMatthew G. Knepley .seealso: PetscFVSetType(), PETSCFVUPWIND
1154*c0ce001eSMatthew G. Knepley @*/
1155*c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVCreate(MPI_Comm comm, PetscFV *fvm)
1156*c0ce001eSMatthew G. Knepley {
1157*c0ce001eSMatthew G. Knepley   PetscFV        f;
1158*c0ce001eSMatthew G. Knepley   PetscErrorCode ierr;
1159*c0ce001eSMatthew G. Knepley 
1160*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1161*c0ce001eSMatthew G. Knepley   PetscValidPointer(fvm, 2);
1162*c0ce001eSMatthew G. Knepley   *fvm = NULL;
1163*c0ce001eSMatthew G. Knepley   ierr = PetscFVInitializePackage();CHKERRQ(ierr);
1164*c0ce001eSMatthew G. Knepley 
1165*c0ce001eSMatthew G. Knepley   ierr = PetscHeaderCreate(f, PETSCFV_CLASSID, "PetscFV", "Finite Volume", "PetscFV", comm, PetscFVDestroy, PetscFVView);CHKERRQ(ierr);
1166*c0ce001eSMatthew G. Knepley   ierr = PetscMemzero(f->ops, sizeof(struct _PetscFVOps));CHKERRQ(ierr);
1167*c0ce001eSMatthew G. Knepley 
1168*c0ce001eSMatthew G. Knepley   ierr = PetscLimiterCreate(comm, &f->limiter);CHKERRQ(ierr);
1169*c0ce001eSMatthew G. Knepley   f->numComponents    = 1;
1170*c0ce001eSMatthew G. Knepley   f->dim              = 0;
1171*c0ce001eSMatthew G. Knepley   f->computeGradients = PETSC_FALSE;
1172*c0ce001eSMatthew G. Knepley   f->fluxWork         = NULL;
1173*c0ce001eSMatthew G. Knepley   ierr = PetscCalloc1(f->numComponents,&f->componentNames);CHKERRQ(ierr);
1174*c0ce001eSMatthew G. Knepley 
1175*c0ce001eSMatthew G. Knepley   *fvm = f;
1176*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
1177*c0ce001eSMatthew G. Knepley }
1178*c0ce001eSMatthew G. Knepley 
1179*c0ce001eSMatthew G. Knepley /*@
1180*c0ce001eSMatthew G. Knepley   PetscFVSetLimiter - Set the limiter object
1181*c0ce001eSMatthew G. Knepley 
1182*c0ce001eSMatthew G. Knepley   Logically collective on fvm
1183*c0ce001eSMatthew G. Knepley 
1184*c0ce001eSMatthew G. Knepley   Input Parameters:
1185*c0ce001eSMatthew G. Knepley + fvm - the PetscFV object
1186*c0ce001eSMatthew G. Knepley - lim - The PetscLimiter
1187*c0ce001eSMatthew G. Knepley 
1188*c0ce001eSMatthew G. Knepley   Level: developer
1189*c0ce001eSMatthew G. Knepley 
1190*c0ce001eSMatthew G. Knepley .seealso: PetscFVGetLimiter()
1191*c0ce001eSMatthew G. Knepley @*/
1192*c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVSetLimiter(PetscFV fvm, PetscLimiter lim)
1193*c0ce001eSMatthew G. Knepley {
1194*c0ce001eSMatthew G. Knepley   PetscErrorCode ierr;
1195*c0ce001eSMatthew G. Knepley 
1196*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1197*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1);
1198*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 2);
1199*c0ce001eSMatthew G. Knepley   ierr = PetscLimiterDestroy(&fvm->limiter);CHKERRQ(ierr);
1200*c0ce001eSMatthew G. Knepley   ierr = PetscObjectReference((PetscObject) lim);CHKERRQ(ierr);
1201*c0ce001eSMatthew G. Knepley   fvm->limiter = lim;
1202*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
1203*c0ce001eSMatthew G. Knepley }
1204*c0ce001eSMatthew G. Knepley 
1205*c0ce001eSMatthew G. Knepley /*@
1206*c0ce001eSMatthew G. Knepley   PetscFVGetLimiter - Get the limiter object
1207*c0ce001eSMatthew G. Knepley 
1208*c0ce001eSMatthew G. Knepley   Not collective
1209*c0ce001eSMatthew G. Knepley 
1210*c0ce001eSMatthew G. Knepley   Input Parameter:
1211*c0ce001eSMatthew G. Knepley . fvm - the PetscFV object
1212*c0ce001eSMatthew G. Knepley 
1213*c0ce001eSMatthew G. Knepley   Output Parameter:
1214*c0ce001eSMatthew G. Knepley . lim - The PetscLimiter
1215*c0ce001eSMatthew G. Knepley 
1216*c0ce001eSMatthew G. Knepley   Level: developer
1217*c0ce001eSMatthew G. Knepley 
1218*c0ce001eSMatthew G. Knepley .seealso: PetscFVSetLimiter()
1219*c0ce001eSMatthew G. Knepley @*/
1220*c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVGetLimiter(PetscFV fvm, PetscLimiter *lim)
1221*c0ce001eSMatthew G. Knepley {
1222*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1223*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1);
1224*c0ce001eSMatthew G. Knepley   PetscValidPointer(lim, 2);
1225*c0ce001eSMatthew G. Knepley   *lim = fvm->limiter;
1226*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
1227*c0ce001eSMatthew G. Knepley }
1228*c0ce001eSMatthew G. Knepley 
1229*c0ce001eSMatthew G. Knepley /*@
1230*c0ce001eSMatthew G. Knepley   PetscFVSetNumComponents - Set the number of field components
1231*c0ce001eSMatthew G. Knepley 
1232*c0ce001eSMatthew G. Knepley   Logically collective on fvm
1233*c0ce001eSMatthew G. Knepley 
1234*c0ce001eSMatthew G. Knepley   Input Parameters:
1235*c0ce001eSMatthew G. Knepley + fvm - the PetscFV object
1236*c0ce001eSMatthew G. Knepley - comp - The number of components
1237*c0ce001eSMatthew G. Knepley 
1238*c0ce001eSMatthew G. Knepley   Level: developer
1239*c0ce001eSMatthew G. Knepley 
1240*c0ce001eSMatthew G. Knepley .seealso: PetscFVGetNumComponents()
1241*c0ce001eSMatthew G. Knepley @*/
1242*c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVSetNumComponents(PetscFV fvm, PetscInt comp)
1243*c0ce001eSMatthew G. Knepley {
1244*c0ce001eSMatthew G. Knepley   PetscErrorCode ierr;
1245*c0ce001eSMatthew G. Knepley 
1246*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1247*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1);
1248*c0ce001eSMatthew G. Knepley   if (fvm->numComponents != comp) {
1249*c0ce001eSMatthew G. Knepley     PetscInt i;
1250*c0ce001eSMatthew G. Knepley 
1251*c0ce001eSMatthew G. Knepley     for (i = 0; i < fvm->numComponents; i++) {
1252*c0ce001eSMatthew G. Knepley       ierr = PetscFree(fvm->componentNames[i]);CHKERRQ(ierr);
1253*c0ce001eSMatthew G. Knepley     }
1254*c0ce001eSMatthew G. Knepley     ierr = PetscFree(fvm->componentNames);CHKERRQ(ierr);
1255*c0ce001eSMatthew G. Knepley     ierr = PetscCalloc1(comp,&fvm->componentNames);CHKERRQ(ierr);
1256*c0ce001eSMatthew G. Knepley   }
1257*c0ce001eSMatthew G. Knepley   fvm->numComponents = comp;
1258*c0ce001eSMatthew G. Knepley   ierr = PetscFree(fvm->fluxWork);CHKERRQ(ierr);
1259*c0ce001eSMatthew G. Knepley   ierr = PetscMalloc1(comp, &fvm->fluxWork);CHKERRQ(ierr);
1260*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
1261*c0ce001eSMatthew G. Knepley }
1262*c0ce001eSMatthew G. Knepley 
1263*c0ce001eSMatthew G. Knepley /*@
1264*c0ce001eSMatthew G. Knepley   PetscFVGetNumComponents - Get the number of field components
1265*c0ce001eSMatthew G. Knepley 
1266*c0ce001eSMatthew G. Knepley   Not collective
1267*c0ce001eSMatthew G. Knepley 
1268*c0ce001eSMatthew G. Knepley   Input Parameter:
1269*c0ce001eSMatthew G. Knepley . fvm - the PetscFV object
1270*c0ce001eSMatthew G. Knepley 
1271*c0ce001eSMatthew G. Knepley   Output Parameter:
1272*c0ce001eSMatthew G. Knepley , comp - The number of components
1273*c0ce001eSMatthew G. Knepley 
1274*c0ce001eSMatthew G. Knepley   Level: developer
1275*c0ce001eSMatthew G. Knepley 
1276*c0ce001eSMatthew G. Knepley .seealso: PetscFVSetNumComponents()
1277*c0ce001eSMatthew G. Knepley @*/
1278*c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVGetNumComponents(PetscFV fvm, PetscInt *comp)
1279*c0ce001eSMatthew G. Knepley {
1280*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1281*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1);
1282*c0ce001eSMatthew G. Knepley   PetscValidPointer(comp, 2);
1283*c0ce001eSMatthew G. Knepley   *comp = fvm->numComponents;
1284*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
1285*c0ce001eSMatthew G. Knepley }
1286*c0ce001eSMatthew G. Knepley 
1287*c0ce001eSMatthew G. Knepley /*@C
1288*c0ce001eSMatthew G. Knepley   PetscFVSetComponentName - Set the name of a component (used in output and viewing)
1289*c0ce001eSMatthew G. Knepley 
1290*c0ce001eSMatthew G. Knepley   Logically collective on fvm
1291*c0ce001eSMatthew G. Knepley   Input Parameters:
1292*c0ce001eSMatthew G. Knepley + fvm - the PetscFV object
1293*c0ce001eSMatthew G. Knepley . comp - the component number
1294*c0ce001eSMatthew G. Knepley - name - the component name
1295*c0ce001eSMatthew G. Knepley 
1296*c0ce001eSMatthew G. Knepley   Level: developer
1297*c0ce001eSMatthew G. Knepley 
1298*c0ce001eSMatthew G. Knepley .seealso: PetscFVGetComponentName()
1299*c0ce001eSMatthew G. Knepley @*/
1300*c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVSetComponentName(PetscFV fvm, PetscInt comp, const char *name)
1301*c0ce001eSMatthew G. Knepley {
1302*c0ce001eSMatthew G. Knepley   PetscErrorCode ierr;
1303*c0ce001eSMatthew G. Knepley 
1304*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1305*c0ce001eSMatthew G. Knepley   ierr = PetscFree(fvm->componentNames[comp]);CHKERRQ(ierr);
1306*c0ce001eSMatthew G. Knepley   ierr = PetscStrallocpy(name,&fvm->componentNames[comp]);CHKERRQ(ierr);
1307*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
1308*c0ce001eSMatthew G. Knepley }
1309*c0ce001eSMatthew G. Knepley 
1310*c0ce001eSMatthew G. Knepley /*@C
1311*c0ce001eSMatthew G. Knepley   PetscFVGetComponentName - Get the name of a component (used in output and viewing)
1312*c0ce001eSMatthew G. Knepley 
1313*c0ce001eSMatthew G. Knepley   Logically collective on fvm
1314*c0ce001eSMatthew G. Knepley   Input Parameters:
1315*c0ce001eSMatthew G. Knepley + fvm - the PetscFV object
1316*c0ce001eSMatthew G. Knepley - comp - the component number
1317*c0ce001eSMatthew G. Knepley 
1318*c0ce001eSMatthew G. Knepley   Output Parameter:
1319*c0ce001eSMatthew G. Knepley . name - the component name
1320*c0ce001eSMatthew G. Knepley 
1321*c0ce001eSMatthew G. Knepley   Level: developer
1322*c0ce001eSMatthew G. Knepley 
1323*c0ce001eSMatthew G. Knepley .seealso: PetscFVSetComponentName()
1324*c0ce001eSMatthew G. Knepley @*/
1325*c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVGetComponentName(PetscFV fvm, PetscInt comp, const char **name)
1326*c0ce001eSMatthew G. Knepley {
1327*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1328*c0ce001eSMatthew G. Knepley   *name = fvm->componentNames[comp];
1329*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
1330*c0ce001eSMatthew G. Knepley }
1331*c0ce001eSMatthew G. Knepley 
1332*c0ce001eSMatthew G. Knepley /*@
1333*c0ce001eSMatthew G. Knepley   PetscFVSetSpatialDimension - Set the spatial dimension
1334*c0ce001eSMatthew G. Knepley 
1335*c0ce001eSMatthew G. Knepley   Logically collective on fvm
1336*c0ce001eSMatthew G. Knepley 
1337*c0ce001eSMatthew G. Knepley   Input Parameters:
1338*c0ce001eSMatthew G. Knepley + fvm - the PetscFV object
1339*c0ce001eSMatthew G. Knepley - dim - The spatial dimension
1340*c0ce001eSMatthew G. Knepley 
1341*c0ce001eSMatthew G. Knepley   Level: developer
1342*c0ce001eSMatthew G. Knepley 
1343*c0ce001eSMatthew G. Knepley .seealso: PetscFVGetSpatialDimension()
1344*c0ce001eSMatthew G. Knepley @*/
1345*c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVSetSpatialDimension(PetscFV fvm, PetscInt dim)
1346*c0ce001eSMatthew G. Knepley {
1347*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1348*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1);
1349*c0ce001eSMatthew G. Knepley   fvm->dim = dim;
1350*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
1351*c0ce001eSMatthew G. Knepley }
1352*c0ce001eSMatthew G. Knepley 
1353*c0ce001eSMatthew G. Knepley /*@
1354*c0ce001eSMatthew G. Knepley   PetscFVGetSpatialDimension - Get the spatial dimension
1355*c0ce001eSMatthew G. Knepley 
1356*c0ce001eSMatthew G. Knepley   Logically collective on fvm
1357*c0ce001eSMatthew G. Knepley 
1358*c0ce001eSMatthew G. Knepley   Input Parameter:
1359*c0ce001eSMatthew G. Knepley . fvm - the PetscFV object
1360*c0ce001eSMatthew G. Knepley 
1361*c0ce001eSMatthew G. Knepley   Output Parameter:
1362*c0ce001eSMatthew G. Knepley . dim - The spatial dimension
1363*c0ce001eSMatthew G. Knepley 
1364*c0ce001eSMatthew G. Knepley   Level: developer
1365*c0ce001eSMatthew G. Knepley 
1366*c0ce001eSMatthew G. Knepley .seealso: PetscFVSetSpatialDimension()
1367*c0ce001eSMatthew G. Knepley @*/
1368*c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVGetSpatialDimension(PetscFV fvm, PetscInt *dim)
1369*c0ce001eSMatthew G. Knepley {
1370*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1371*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1);
1372*c0ce001eSMatthew G. Knepley   PetscValidPointer(dim, 2);
1373*c0ce001eSMatthew G. Knepley   *dim = fvm->dim;
1374*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
1375*c0ce001eSMatthew G. Knepley }
1376*c0ce001eSMatthew G. Knepley 
1377*c0ce001eSMatthew G. Knepley /*@
1378*c0ce001eSMatthew G. Knepley   PetscFVSetComputeGradients - Toggle computation of cell gradients
1379*c0ce001eSMatthew G. Knepley 
1380*c0ce001eSMatthew G. Knepley   Logically collective on fvm
1381*c0ce001eSMatthew G. Knepley 
1382*c0ce001eSMatthew G. Knepley   Input Parameters:
1383*c0ce001eSMatthew G. Knepley + fvm - the PetscFV object
1384*c0ce001eSMatthew G. Knepley - computeGradients - Flag to compute cell gradients
1385*c0ce001eSMatthew G. Knepley 
1386*c0ce001eSMatthew G. Knepley   Level: developer
1387*c0ce001eSMatthew G. Knepley 
1388*c0ce001eSMatthew G. Knepley .seealso: PetscFVGetComputeGradients()
1389*c0ce001eSMatthew G. Knepley @*/
1390*c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVSetComputeGradients(PetscFV fvm, PetscBool computeGradients)
1391*c0ce001eSMatthew G. Knepley {
1392*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1393*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1);
1394*c0ce001eSMatthew G. Knepley   fvm->computeGradients = computeGradients;
1395*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
1396*c0ce001eSMatthew G. Knepley }
1397*c0ce001eSMatthew G. Knepley 
1398*c0ce001eSMatthew G. Knepley /*@
1399*c0ce001eSMatthew G. Knepley   PetscFVGetComputeGradients - Return flag for computation of cell gradients
1400*c0ce001eSMatthew G. Knepley 
1401*c0ce001eSMatthew G. Knepley   Not collective
1402*c0ce001eSMatthew G. Knepley 
1403*c0ce001eSMatthew G. Knepley   Input Parameter:
1404*c0ce001eSMatthew G. Knepley . fvm - the PetscFV object
1405*c0ce001eSMatthew G. Knepley 
1406*c0ce001eSMatthew G. Knepley   Output Parameter:
1407*c0ce001eSMatthew G. Knepley . computeGradients - Flag to compute cell gradients
1408*c0ce001eSMatthew G. Knepley 
1409*c0ce001eSMatthew G. Knepley   Level: developer
1410*c0ce001eSMatthew G. Knepley 
1411*c0ce001eSMatthew G. Knepley .seealso: PetscFVSetComputeGradients()
1412*c0ce001eSMatthew G. Knepley @*/
1413*c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVGetComputeGradients(PetscFV fvm, PetscBool *computeGradients)
1414*c0ce001eSMatthew G. Knepley {
1415*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1416*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1);
1417*c0ce001eSMatthew G. Knepley   PetscValidPointer(computeGradients, 2);
1418*c0ce001eSMatthew G. Knepley   *computeGradients = fvm->computeGradients;
1419*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
1420*c0ce001eSMatthew G. Knepley }
1421*c0ce001eSMatthew G. Knepley 
1422*c0ce001eSMatthew G. Knepley /*@
1423*c0ce001eSMatthew G. Knepley   PetscFVSetQuadrature - Set the quadrature object
1424*c0ce001eSMatthew G. Knepley 
1425*c0ce001eSMatthew G. Knepley   Logically collective on fvm
1426*c0ce001eSMatthew G. Knepley 
1427*c0ce001eSMatthew G. Knepley   Input Parameters:
1428*c0ce001eSMatthew G. Knepley + fvm - the PetscFV object
1429*c0ce001eSMatthew G. Knepley - q - The PetscQuadrature
1430*c0ce001eSMatthew G. Knepley 
1431*c0ce001eSMatthew G. Knepley   Level: developer
1432*c0ce001eSMatthew G. Knepley 
1433*c0ce001eSMatthew G. Knepley .seealso: PetscFVGetQuadrature()
1434*c0ce001eSMatthew G. Knepley @*/
1435*c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVSetQuadrature(PetscFV fvm, PetscQuadrature q)
1436*c0ce001eSMatthew G. Knepley {
1437*c0ce001eSMatthew G. Knepley   PetscErrorCode ierr;
1438*c0ce001eSMatthew G. Knepley 
1439*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1440*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1);
1441*c0ce001eSMatthew G. Knepley   ierr = PetscQuadratureDestroy(&fvm->quadrature);CHKERRQ(ierr);
1442*c0ce001eSMatthew G. Knepley   ierr = PetscObjectReference((PetscObject) q);CHKERRQ(ierr);
1443*c0ce001eSMatthew G. Knepley   fvm->quadrature = q;
1444*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
1445*c0ce001eSMatthew G. Knepley }
1446*c0ce001eSMatthew G. Knepley 
1447*c0ce001eSMatthew G. Knepley /*@
1448*c0ce001eSMatthew G. Knepley   PetscFVGetQuadrature - Get the quadrature object
1449*c0ce001eSMatthew G. Knepley 
1450*c0ce001eSMatthew G. Knepley   Not collective
1451*c0ce001eSMatthew G. Knepley 
1452*c0ce001eSMatthew G. Knepley   Input Parameter:
1453*c0ce001eSMatthew G. Knepley . fvm - the PetscFV object
1454*c0ce001eSMatthew G. Knepley 
1455*c0ce001eSMatthew G. Knepley   Output Parameter:
1456*c0ce001eSMatthew G. Knepley . lim - The PetscQuadrature
1457*c0ce001eSMatthew G. Knepley 
1458*c0ce001eSMatthew G. Knepley   Level: developer
1459*c0ce001eSMatthew G. Knepley 
1460*c0ce001eSMatthew G. Knepley .seealso: PetscFVSetQuadrature()
1461*c0ce001eSMatthew G. Knepley @*/
1462*c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVGetQuadrature(PetscFV fvm, PetscQuadrature *q)
1463*c0ce001eSMatthew G. Knepley {
1464*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1465*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1);
1466*c0ce001eSMatthew G. Knepley   PetscValidPointer(q, 2);
1467*c0ce001eSMatthew G. Knepley   if (!fvm->quadrature) {
1468*c0ce001eSMatthew G. Knepley     /* Create default 1-point quadrature */
1469*c0ce001eSMatthew G. Knepley     PetscReal     *points, *weights;
1470*c0ce001eSMatthew G. Knepley     PetscErrorCode ierr;
1471*c0ce001eSMatthew G. Knepley 
1472*c0ce001eSMatthew G. Knepley     ierr = PetscQuadratureCreate(PETSC_COMM_SELF, &fvm->quadrature);CHKERRQ(ierr);
1473*c0ce001eSMatthew G. Knepley     ierr = PetscCalloc1(fvm->dim, &points);CHKERRQ(ierr);
1474*c0ce001eSMatthew G. Knepley     ierr = PetscMalloc1(1, &weights);CHKERRQ(ierr);
1475*c0ce001eSMatthew G. Knepley     weights[0] = 1.0;
1476*c0ce001eSMatthew G. Knepley     ierr = PetscQuadratureSetData(fvm->quadrature, fvm->dim, 1, 1, points, weights);CHKERRQ(ierr);
1477*c0ce001eSMatthew G. Knepley   }
1478*c0ce001eSMatthew G. Knepley   *q = fvm->quadrature;
1479*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
1480*c0ce001eSMatthew G. Knepley }
1481*c0ce001eSMatthew G. Knepley 
1482*c0ce001eSMatthew G. Knepley /*@
1483*c0ce001eSMatthew G. Knepley   PetscFVGetDualSpace - Returns the PetscDualSpace used to define the inner product
1484*c0ce001eSMatthew G. Knepley 
1485*c0ce001eSMatthew G. Knepley   Not collective
1486*c0ce001eSMatthew G. Knepley 
1487*c0ce001eSMatthew G. Knepley   Input Parameter:
1488*c0ce001eSMatthew G. Knepley . fvm - The PetscFV object
1489*c0ce001eSMatthew G. Knepley 
1490*c0ce001eSMatthew G. Knepley   Output Parameter:
1491*c0ce001eSMatthew G. Knepley . sp - The PetscDualSpace object
1492*c0ce001eSMatthew G. Knepley 
1493*c0ce001eSMatthew G. Knepley   Note: A simple dual space is provided automatically, and the user typically will not need to override it.
1494*c0ce001eSMatthew G. Knepley 
1495*c0ce001eSMatthew G. Knepley   Level: developer
1496*c0ce001eSMatthew G. Knepley 
1497*c0ce001eSMatthew G. Knepley .seealso: PetscFVCreate()
1498*c0ce001eSMatthew G. Knepley @*/
1499*c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVGetDualSpace(PetscFV fvm, PetscDualSpace *sp)
1500*c0ce001eSMatthew G. Knepley {
1501*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1502*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1);
1503*c0ce001eSMatthew G. Knepley   PetscValidPointer(sp, 2);
1504*c0ce001eSMatthew G. Knepley   if (!fvm->dualSpace) {
1505*c0ce001eSMatthew G. Knepley     DM              K;
1506*c0ce001eSMatthew G. Knepley     PetscInt        dim, Nc, c;
1507*c0ce001eSMatthew G. Knepley     PetscErrorCode  ierr;
1508*c0ce001eSMatthew G. Knepley 
1509*c0ce001eSMatthew G. Knepley     ierr = PetscFVGetSpatialDimension(fvm, &dim);CHKERRQ(ierr);
1510*c0ce001eSMatthew G. Knepley     ierr = PetscFVGetNumComponents(fvm, &Nc);CHKERRQ(ierr);
1511*c0ce001eSMatthew G. Knepley     ierr = PetscDualSpaceCreate(PetscObjectComm((PetscObject) fvm), &fvm->dualSpace);CHKERRQ(ierr);
1512*c0ce001eSMatthew G. Knepley     ierr = PetscDualSpaceSetType(fvm->dualSpace, PETSCDUALSPACESIMPLE);CHKERRQ(ierr);
1513*c0ce001eSMatthew G. Knepley     ierr = PetscDualSpaceCreateReferenceCell(fvm->dualSpace, dim, PETSC_FALSE, &K);CHKERRQ(ierr); /* TODO: The reference cell type should be held by the discretization object */
1514*c0ce001eSMatthew G. Knepley     ierr = PetscDualSpaceSetNumComponents(fvm->dualSpace, Nc);CHKERRQ(ierr);
1515*c0ce001eSMatthew G. Knepley     ierr = PetscDualSpaceSetDM(fvm->dualSpace, K);CHKERRQ(ierr);
1516*c0ce001eSMatthew G. Knepley     ierr = DMDestroy(&K);CHKERRQ(ierr);
1517*c0ce001eSMatthew G. Knepley     ierr = PetscDualSpaceSimpleSetDimension(fvm->dualSpace, Nc);CHKERRQ(ierr);
1518*c0ce001eSMatthew G. Knepley     /* Should we be using PetscFVGetQuadrature() here? */
1519*c0ce001eSMatthew G. Knepley     for (c = 0; c < Nc; ++c) {
1520*c0ce001eSMatthew G. Knepley       PetscQuadrature qc;
1521*c0ce001eSMatthew G. Knepley       PetscReal      *points, *weights;
1522*c0ce001eSMatthew G. Knepley       PetscErrorCode  ierr;
1523*c0ce001eSMatthew G. Knepley 
1524*c0ce001eSMatthew G. Knepley       ierr = PetscQuadratureCreate(PETSC_COMM_SELF, &qc);CHKERRQ(ierr);
1525*c0ce001eSMatthew G. Knepley       ierr = PetscCalloc1(dim, &points);CHKERRQ(ierr);
1526*c0ce001eSMatthew G. Knepley       ierr = PetscCalloc1(Nc, &weights);CHKERRQ(ierr);
1527*c0ce001eSMatthew G. Knepley       weights[c] = 1.0;
1528*c0ce001eSMatthew G. Knepley       ierr = PetscQuadratureSetData(qc, dim, Nc, 1, points, weights);CHKERRQ(ierr);
1529*c0ce001eSMatthew G. Knepley       ierr = PetscDualSpaceSimpleSetFunctional(fvm->dualSpace, c, qc);CHKERRQ(ierr);
1530*c0ce001eSMatthew G. Knepley       ierr = PetscQuadratureDestroy(&qc);CHKERRQ(ierr);
1531*c0ce001eSMatthew G. Knepley     }
1532*c0ce001eSMatthew G. Knepley   }
1533*c0ce001eSMatthew G. Knepley   *sp = fvm->dualSpace;
1534*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
1535*c0ce001eSMatthew G. Knepley }
1536*c0ce001eSMatthew G. Knepley 
1537*c0ce001eSMatthew G. Knepley /*@
1538*c0ce001eSMatthew G. Knepley   PetscFVSetDualSpace - Sets the PetscDualSpace used to define the inner product
1539*c0ce001eSMatthew G. Knepley 
1540*c0ce001eSMatthew G. Knepley   Not collective
1541*c0ce001eSMatthew G. Knepley 
1542*c0ce001eSMatthew G. Knepley   Input Parameters:
1543*c0ce001eSMatthew G. Knepley + fvm - The PetscFV object
1544*c0ce001eSMatthew G. Knepley - sp  - The PetscDualSpace object
1545*c0ce001eSMatthew G. Knepley 
1546*c0ce001eSMatthew G. Knepley   Level: intermediate
1547*c0ce001eSMatthew G. Knepley 
1548*c0ce001eSMatthew G. Knepley   Note: A simple dual space is provided automatically, and the user typically will not need to override it.
1549*c0ce001eSMatthew G. Knepley 
1550*c0ce001eSMatthew G. Knepley .seealso: PetscFVCreate()
1551*c0ce001eSMatthew G. Knepley @*/
1552*c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVSetDualSpace(PetscFV fvm, PetscDualSpace sp)
1553*c0ce001eSMatthew G. Knepley {
1554*c0ce001eSMatthew G. Knepley   PetscErrorCode ierr;
1555*c0ce001eSMatthew G. Knepley 
1556*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1557*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1);
1558*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 2);
1559*c0ce001eSMatthew G. Knepley   ierr = PetscDualSpaceDestroy(&fvm->dualSpace);CHKERRQ(ierr);
1560*c0ce001eSMatthew G. Knepley   fvm->dualSpace = sp;
1561*c0ce001eSMatthew G. Knepley   ierr = PetscObjectReference((PetscObject) fvm->dualSpace);CHKERRQ(ierr);
1562*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
1563*c0ce001eSMatthew G. Knepley }
1564*c0ce001eSMatthew G. Knepley 
1565*c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVGetDefaultTabulation(PetscFV fvm, PetscReal **B, PetscReal **D, PetscReal **H)
1566*c0ce001eSMatthew G. Knepley {
1567*c0ce001eSMatthew G. Knepley   PetscInt         npoints;
1568*c0ce001eSMatthew G. Knepley   const PetscReal *points;
1569*c0ce001eSMatthew G. Knepley   PetscErrorCode   ierr;
1570*c0ce001eSMatthew G. Knepley 
1571*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1572*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1);
1573*c0ce001eSMatthew G. Knepley   if (B) PetscValidPointer(B, 2);
1574*c0ce001eSMatthew G. Knepley   if (D) PetscValidPointer(D, 3);
1575*c0ce001eSMatthew G. Knepley   if (H) PetscValidPointer(H, 4);
1576*c0ce001eSMatthew G. Knepley   ierr = PetscQuadratureGetData(fvm->quadrature, NULL, NULL, &npoints, &points, NULL);CHKERRQ(ierr);
1577*c0ce001eSMatthew G. Knepley   if (!fvm->B) {ierr = PetscFVGetTabulation(fvm, npoints, points, &fvm->B, &fvm->D, NULL/*&fvm->H*/);CHKERRQ(ierr);}
1578*c0ce001eSMatthew G. Knepley   if (B) *B = fvm->B;
1579*c0ce001eSMatthew G. Knepley   if (D) *D = fvm->D;
1580*c0ce001eSMatthew G. Knepley   if (H) *H = fvm->H;
1581*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
1582*c0ce001eSMatthew G. Knepley }
1583*c0ce001eSMatthew G. Knepley 
1584*c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVGetTabulation(PetscFV fvm, PetscInt npoints, const PetscReal points[], PetscReal **B, PetscReal **D, PetscReal **H)
1585*c0ce001eSMatthew G. Knepley {
1586*c0ce001eSMatthew G. Knepley   PetscInt         pdim = 1; /* Dimension of approximation space P */
1587*c0ce001eSMatthew G. Knepley   PetscInt         dim;      /* Spatial dimension */
1588*c0ce001eSMatthew G. Knepley   PetscInt         comp;     /* Field components */
1589*c0ce001eSMatthew G. Knepley   PetscInt         p, d, c, e;
1590*c0ce001eSMatthew G. Knepley   PetscErrorCode   ierr;
1591*c0ce001eSMatthew G. Knepley 
1592*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1593*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1);
1594*c0ce001eSMatthew G. Knepley   PetscValidPointer(points, 3);
1595*c0ce001eSMatthew G. Knepley   if (B) PetscValidPointer(B, 4);
1596*c0ce001eSMatthew G. Knepley   if (D) PetscValidPointer(D, 5);
1597*c0ce001eSMatthew G. Knepley   if (H) PetscValidPointer(H, 6);
1598*c0ce001eSMatthew G. Knepley   ierr = PetscFVGetSpatialDimension(fvm, &dim);CHKERRQ(ierr);
1599*c0ce001eSMatthew G. Knepley   ierr = PetscFVGetNumComponents(fvm, &comp);CHKERRQ(ierr);
1600*c0ce001eSMatthew G. Knepley   if (B) {ierr = PetscMalloc1(npoints*pdim*comp, B);CHKERRQ(ierr);}
1601*c0ce001eSMatthew G. Knepley   if (D) {ierr = PetscMalloc1(npoints*pdim*comp*dim, D);CHKERRQ(ierr);}
1602*c0ce001eSMatthew G. Knepley   if (H) {ierr = PetscMalloc1(npoints*pdim*comp*dim*dim, H);CHKERRQ(ierr);}
1603*c0ce001eSMatthew G. Knepley   if (B) {for (p = 0; p < npoints; ++p) for (d = 0; d < pdim; ++d) for (c = 0; c < comp; ++c) (*B)[(p*pdim + d)*comp + c] = 1.0;}
1604*c0ce001eSMatthew G. Knepley   if (D) {for (p = 0; p < npoints; ++p) for (d = 0; d < pdim; ++d) for (c = 0; c < comp; ++c) for (e = 0; e < dim; ++e) (*D)[((p*pdim + d)*comp + c)*dim + e] = 0.0;}
1605*c0ce001eSMatthew G. Knepley   if (H) {for (p = 0; p < npoints; ++p) for (d = 0; d < pdim; ++d) for (c = 0; c < comp; ++c) for (e = 0; e < dim*dim; ++e) (*H)[((p*pdim + d)*comp + c)*dim*dim + e] = 0.0;}
1606*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
1607*c0ce001eSMatthew G. Knepley }
1608*c0ce001eSMatthew G. Knepley 
1609*c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVRestoreTabulation(PetscFV fvm, PetscInt npoints, const PetscReal points[], PetscReal **B, PetscReal **D, PetscReal **H)
1610*c0ce001eSMatthew G. Knepley {
1611*c0ce001eSMatthew G. Knepley   PetscErrorCode ierr;
1612*c0ce001eSMatthew G. Knepley 
1613*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1614*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1);
1615*c0ce001eSMatthew G. Knepley   if (B && *B) {ierr = PetscFree(*B);CHKERRQ(ierr);}
1616*c0ce001eSMatthew G. Knepley   if (D && *D) {ierr = PetscFree(*D);CHKERRQ(ierr);}
1617*c0ce001eSMatthew G. Knepley   if (H && *H) {ierr = PetscFree(*H);CHKERRQ(ierr);}
1618*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
1619*c0ce001eSMatthew G. Knepley }
1620*c0ce001eSMatthew G. Knepley 
1621*c0ce001eSMatthew G. Knepley /*@C
1622*c0ce001eSMatthew G. Knepley   PetscFVComputeGradient - Compute the gradient reconstruction matrix for a given cell
1623*c0ce001eSMatthew G. Knepley 
1624*c0ce001eSMatthew G. Knepley   Input Parameters:
1625*c0ce001eSMatthew G. Knepley + fvm      - The PetscFV object
1626*c0ce001eSMatthew G. Knepley . numFaces - The number of cell faces which are not constrained
1627*c0ce001eSMatthew G. Knepley - dx       - The vector from the cell centroid to the neighboring cell centroid for each face
1628*c0ce001eSMatthew G. Knepley 
1629*c0ce001eSMatthew G. Knepley   Level: developer
1630*c0ce001eSMatthew G. Knepley 
1631*c0ce001eSMatthew G. Knepley .seealso: PetscFVCreate()
1632*c0ce001eSMatthew G. Knepley @*/
1633*c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVComputeGradient(PetscFV fvm, PetscInt numFaces, PetscScalar dx[], PetscScalar grad[])
1634*c0ce001eSMatthew G. Knepley {
1635*c0ce001eSMatthew G. Knepley   PetscErrorCode ierr;
1636*c0ce001eSMatthew G. Knepley 
1637*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1638*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1);
1639*c0ce001eSMatthew G. Knepley   if (fvm->ops->computegradient) {ierr = (*fvm->ops->computegradient)(fvm, numFaces, dx, grad);CHKERRQ(ierr);}
1640*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
1641*c0ce001eSMatthew G. Knepley }
1642*c0ce001eSMatthew G. Knepley 
1643*c0ce001eSMatthew G. Knepley /*C
1644*c0ce001eSMatthew G. Knepley   PetscFVIntegrateRHSFunction - Produce the cell residual vector for a chunk of elements by quadrature integration
1645*c0ce001eSMatthew G. Knepley 
1646*c0ce001eSMatthew G. Knepley   Not collective
1647*c0ce001eSMatthew G. Knepley 
1648*c0ce001eSMatthew G. Knepley   Input Parameters:
1649*c0ce001eSMatthew G. Knepley + fvm          - The PetscFV object for the field being integrated
1650*c0ce001eSMatthew G. Knepley . prob         - The PetscDS specifing the discretizations and continuum functions
1651*c0ce001eSMatthew G. Knepley . field        - The field being integrated
1652*c0ce001eSMatthew G. Knepley . Nf           - The number of faces in the chunk
1653*c0ce001eSMatthew G. Knepley . fgeom        - The face geometry for each face in the chunk
1654*c0ce001eSMatthew G. Knepley . neighborVol  - The volume for each pair of cells in the chunk
1655*c0ce001eSMatthew G. Knepley . uL           - The state from the cell on the left
1656*c0ce001eSMatthew G. Knepley - uR           - The state from the cell on the right
1657*c0ce001eSMatthew G. Knepley 
1658*c0ce001eSMatthew G. Knepley   Output Parameter
1659*c0ce001eSMatthew G. Knepley + fluxL        - the left fluxes for each face
1660*c0ce001eSMatthew G. Knepley - fluxR        - the right fluxes for each face
1661*c0ce001eSMatthew G. Knepley */
1662*c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVIntegrateRHSFunction(PetscFV fvm, PetscDS prob, PetscInt field, PetscInt Nf, PetscFVFaceGeom *fgeom, PetscReal *neighborVol,
1663*c0ce001eSMatthew G. Knepley                                            PetscScalar uL[], PetscScalar uR[], PetscScalar fluxL[], PetscScalar fluxR[])
1664*c0ce001eSMatthew G. Knepley {
1665*c0ce001eSMatthew G. Knepley   PetscErrorCode ierr;
1666*c0ce001eSMatthew G. Knepley 
1667*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1668*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1);
1669*c0ce001eSMatthew G. Knepley   if (fvm->ops->integraterhsfunction) {ierr = (*fvm->ops->integraterhsfunction)(fvm, prob, field, Nf, fgeom, neighborVol, uL, uR, fluxL, fluxR);CHKERRQ(ierr);}
1670*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
1671*c0ce001eSMatthew G. Knepley }
1672*c0ce001eSMatthew G. Knepley 
1673*c0ce001eSMatthew G. Knepley /*@
1674*c0ce001eSMatthew G. Knepley   PetscFVRefine - Create a "refined" PetscFV object that refines the reference cell into smaller copies. This is typically used
1675*c0ce001eSMatthew G. Knepley   to precondition a higher order method with a lower order method on a refined mesh having the same number of dofs (but more
1676*c0ce001eSMatthew G. Knepley   sparsity). It is also used to create an interpolation between regularly refined meshes.
1677*c0ce001eSMatthew G. Knepley 
1678*c0ce001eSMatthew G. Knepley   Input Parameter:
1679*c0ce001eSMatthew G. Knepley . fv - The initial PetscFV
1680*c0ce001eSMatthew G. Knepley 
1681*c0ce001eSMatthew G. Knepley   Output Parameter:
1682*c0ce001eSMatthew G. Knepley . fvRef - The refined PetscFV
1683*c0ce001eSMatthew G. Knepley 
1684*c0ce001eSMatthew G. Knepley   Level: developer
1685*c0ce001eSMatthew G. Knepley 
1686*c0ce001eSMatthew G. Knepley .seealso: PetscFVType, PetscFVCreate(), PetscFVSetType()
1687*c0ce001eSMatthew G. Knepley @*/
1688*c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVRefine(PetscFV fv, PetscFV *fvRef)
1689*c0ce001eSMatthew G. Knepley {
1690*c0ce001eSMatthew G. Knepley   PetscDualSpace   Q, Qref;
1691*c0ce001eSMatthew G. Knepley   DM               K, Kref;
1692*c0ce001eSMatthew G. Knepley   PetscQuadrature  q, qref;
1693*c0ce001eSMatthew G. Knepley   CellRefiner      cellRefiner;
1694*c0ce001eSMatthew G. Knepley   PetscReal       *v0;
1695*c0ce001eSMatthew G. Knepley   PetscReal       *jac, *invjac;
1696*c0ce001eSMatthew G. Knepley   PetscInt         numComp, numSubelements, s;
1697*c0ce001eSMatthew G. Knepley   PetscErrorCode   ierr;
1698*c0ce001eSMatthew G. Knepley 
1699*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1700*c0ce001eSMatthew G. Knepley   ierr = PetscFVGetDualSpace(fv, &Q);CHKERRQ(ierr);
1701*c0ce001eSMatthew G. Knepley   ierr = PetscFVGetQuadrature(fv, &q);CHKERRQ(ierr);
1702*c0ce001eSMatthew G. Knepley   ierr = PetscDualSpaceGetDM(Q, &K);CHKERRQ(ierr);
1703*c0ce001eSMatthew G. Knepley   /* Create dual space */
1704*c0ce001eSMatthew G. Knepley   ierr = PetscDualSpaceDuplicate(Q, &Qref);CHKERRQ(ierr);
1705*c0ce001eSMatthew G. Knepley   ierr = DMRefine(K, PetscObjectComm((PetscObject) fv), &Kref);CHKERRQ(ierr);
1706*c0ce001eSMatthew G. Knepley   ierr = PetscDualSpaceSetDM(Qref, Kref);CHKERRQ(ierr);
1707*c0ce001eSMatthew G. Knepley   ierr = DMDestroy(&Kref);CHKERRQ(ierr);
1708*c0ce001eSMatthew G. Knepley   ierr = PetscDualSpaceSetUp(Qref);CHKERRQ(ierr);
1709*c0ce001eSMatthew G. Knepley   /* Create volume */
1710*c0ce001eSMatthew G. Knepley   ierr = PetscFVCreate(PetscObjectComm((PetscObject) fv), fvRef);CHKERRQ(ierr);
1711*c0ce001eSMatthew G. Knepley   ierr = PetscFVSetDualSpace(*fvRef, Qref);CHKERRQ(ierr);
1712*c0ce001eSMatthew G. Knepley   ierr = PetscFVGetNumComponents(fv,    &numComp);CHKERRQ(ierr);
1713*c0ce001eSMatthew G. Knepley   ierr = PetscFVSetNumComponents(*fvRef, numComp);CHKERRQ(ierr);
1714*c0ce001eSMatthew G. Knepley   ierr = PetscFVSetUp(*fvRef);CHKERRQ(ierr);
1715*c0ce001eSMatthew G. Knepley   /* Create quadrature */
1716*c0ce001eSMatthew G. Knepley   ierr = DMPlexGetCellRefiner_Internal(K, &cellRefiner);CHKERRQ(ierr);
1717*c0ce001eSMatthew G. Knepley   ierr = CellRefinerGetAffineTransforms_Internal(cellRefiner, &numSubelements, &v0, &jac, &invjac);CHKERRQ(ierr);
1718*c0ce001eSMatthew G. Knepley   ierr = PetscQuadratureExpandComposite(q, numSubelements, v0, jac, &qref);CHKERRQ(ierr);
1719*c0ce001eSMatthew G. Knepley   ierr = PetscDualSpaceSimpleSetDimension(Qref, numSubelements);CHKERRQ(ierr);
1720*c0ce001eSMatthew G. Knepley   for (s = 0; s < numSubelements; ++s) {
1721*c0ce001eSMatthew G. Knepley     PetscQuadrature  qs;
1722*c0ce001eSMatthew G. Knepley     const PetscReal *points, *weights;
1723*c0ce001eSMatthew G. Knepley     PetscReal       *p, *w;
1724*c0ce001eSMatthew G. Knepley     PetscInt         dim, Nc, npoints, np;
1725*c0ce001eSMatthew G. Knepley 
1726*c0ce001eSMatthew G. Knepley     ierr = PetscQuadratureCreate(PETSC_COMM_SELF, &qs);CHKERRQ(ierr);
1727*c0ce001eSMatthew G. Knepley     ierr = PetscQuadratureGetData(q, &dim, &Nc, &npoints, &points, &weights);CHKERRQ(ierr);
1728*c0ce001eSMatthew G. Knepley     np   = npoints/numSubelements;
1729*c0ce001eSMatthew G. Knepley     ierr = PetscMalloc1(np*dim,&p);CHKERRQ(ierr);
1730*c0ce001eSMatthew G. Knepley     ierr = PetscMalloc1(np*Nc,&w);CHKERRQ(ierr);
1731*c0ce001eSMatthew G. Knepley     ierr = PetscArraycpy(p, &points[s*np*dim], np*dim);CHKERRQ(ierr);
1732*c0ce001eSMatthew G. Knepley     ierr = PetscArraycpy(w, &weights[s*np*Nc], np*Nc);CHKERRQ(ierr);
1733*c0ce001eSMatthew G. Knepley     ierr = PetscQuadratureSetData(qs, dim, Nc, np, p, w);CHKERRQ(ierr);
1734*c0ce001eSMatthew G. Knepley     ierr = PetscDualSpaceSimpleSetFunctional(Qref, s, qs);CHKERRQ(ierr);
1735*c0ce001eSMatthew G. Knepley     ierr = PetscQuadratureDestroy(&qs);CHKERRQ(ierr);
1736*c0ce001eSMatthew G. Knepley   }
1737*c0ce001eSMatthew G. Knepley   ierr = CellRefinerRestoreAffineTransforms_Internal(cellRefiner, &numSubelements, &v0, &jac, &invjac);CHKERRQ(ierr);
1738*c0ce001eSMatthew G. Knepley   ierr = PetscFVSetQuadrature(*fvRef, qref);CHKERRQ(ierr);
1739*c0ce001eSMatthew G. Knepley   ierr = PetscQuadratureDestroy(&qref);CHKERRQ(ierr);
1740*c0ce001eSMatthew G. Knepley   ierr = PetscDualSpaceDestroy(&Qref);CHKERRQ(ierr);
1741*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
1742*c0ce001eSMatthew G. Knepley }
1743*c0ce001eSMatthew G. Knepley 
1744*c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVDestroy_Upwind(PetscFV fvm)
1745*c0ce001eSMatthew G. Knepley {
1746*c0ce001eSMatthew G. Knepley   PetscFV_Upwind *b = (PetscFV_Upwind *) fvm->data;
1747*c0ce001eSMatthew G. Knepley   PetscErrorCode ierr;
1748*c0ce001eSMatthew G. Knepley 
1749*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1750*c0ce001eSMatthew G. Knepley   ierr = PetscFree(b);CHKERRQ(ierr);
1751*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
1752*c0ce001eSMatthew G. Knepley }
1753*c0ce001eSMatthew G. Knepley 
1754*c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVView_Upwind_Ascii(PetscFV fv, PetscViewer viewer)
1755*c0ce001eSMatthew G. Knepley {
1756*c0ce001eSMatthew G. Knepley   PetscInt          Nc, c;
1757*c0ce001eSMatthew G. Knepley   PetscViewerFormat format;
1758*c0ce001eSMatthew G. Knepley   PetscErrorCode    ierr;
1759*c0ce001eSMatthew G. Knepley 
1760*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1761*c0ce001eSMatthew G. Knepley   ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr);
1762*c0ce001eSMatthew G. Knepley   ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr);
1763*c0ce001eSMatthew G. Knepley   ierr = PetscViewerASCIIPrintf(viewer, "Upwind Finite Volume:\n");CHKERRQ(ierr);
1764*c0ce001eSMatthew G. Knepley   ierr = PetscViewerASCIIPrintf(viewer, "  num components: %d\n", Nc);CHKERRQ(ierr);
1765*c0ce001eSMatthew G. Knepley   for (c = 0; c < Nc; c++) {
1766*c0ce001eSMatthew G. Knepley     if (fv->componentNames[c]) {
1767*c0ce001eSMatthew G. Knepley       ierr = PetscViewerASCIIPrintf(viewer, "    component %d: %s\n", c, fv->componentNames[c]);CHKERRQ(ierr);
1768*c0ce001eSMatthew G. Knepley     }
1769*c0ce001eSMatthew G. Knepley   }
1770*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
1771*c0ce001eSMatthew G. Knepley }
1772*c0ce001eSMatthew G. Knepley 
1773*c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVView_Upwind(PetscFV fv, PetscViewer viewer)
1774*c0ce001eSMatthew G. Knepley {
1775*c0ce001eSMatthew G. Knepley   PetscBool      iascii;
1776*c0ce001eSMatthew G. Knepley   PetscErrorCode ierr;
1777*c0ce001eSMatthew G. Knepley 
1778*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1779*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(fv, PETSCFV_CLASSID, 1);
1780*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
1781*c0ce001eSMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
1782*c0ce001eSMatthew G. Knepley   if (iascii) {ierr = PetscFVView_Upwind_Ascii(fv, viewer);CHKERRQ(ierr);}
1783*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
1784*c0ce001eSMatthew G. Knepley }
1785*c0ce001eSMatthew G. Knepley 
1786*c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVSetUp_Upwind(PetscFV fvm)
1787*c0ce001eSMatthew G. Knepley {
1788*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1789*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
1790*c0ce001eSMatthew G. Knepley }
1791*c0ce001eSMatthew G. Knepley 
1792*c0ce001eSMatthew G. Knepley /*
1793*c0ce001eSMatthew G. Knepley   neighborVol[f*2+0] contains the left  geom
1794*c0ce001eSMatthew G. Knepley   neighborVol[f*2+1] contains the right geom
1795*c0ce001eSMatthew G. Knepley */
1796*c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVIntegrateRHSFunction_Upwind(PetscFV fvm, PetscDS prob, PetscInt field, PetscInt Nf, PetscFVFaceGeom *fgeom, PetscReal *neighborVol,
1797*c0ce001eSMatthew G. Knepley                                                   PetscScalar uL[], PetscScalar uR[], PetscScalar fluxL[], PetscScalar fluxR[])
1798*c0ce001eSMatthew G. Knepley {
1799*c0ce001eSMatthew G. Knepley   void             (*riemann)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscInt, const PetscScalar[], PetscScalar[], void *);
1800*c0ce001eSMatthew G. Knepley   void              *rctx;
1801*c0ce001eSMatthew G. Knepley   PetscScalar       *flux = fvm->fluxWork;
1802*c0ce001eSMatthew G. Knepley   const PetscScalar *constants;
1803*c0ce001eSMatthew G. Knepley   PetscInt           dim, numConstants, pdim, totDim, Nc, off, f, d;
1804*c0ce001eSMatthew G. Knepley   PetscErrorCode     ierr;
1805*c0ce001eSMatthew G. Knepley 
1806*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1807*c0ce001eSMatthew G. Knepley   ierr = PetscDSGetTotalComponents(prob, &Nc);CHKERRQ(ierr);
1808*c0ce001eSMatthew G. Knepley   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
1809*c0ce001eSMatthew G. Knepley   ierr = PetscDSGetFieldOffset(prob, field, &off);CHKERRQ(ierr);
1810*c0ce001eSMatthew G. Knepley   ierr = PetscDSGetRiemannSolver(prob, field, &riemann);CHKERRQ(ierr);
1811*c0ce001eSMatthew G. Knepley   ierr = PetscDSGetContext(prob, field, &rctx);CHKERRQ(ierr);
1812*c0ce001eSMatthew G. Knepley   ierr = PetscDSGetConstants(prob, &numConstants, &constants);CHKERRQ(ierr);
1813*c0ce001eSMatthew G. Knepley   ierr = PetscFVGetSpatialDimension(fvm, &dim);CHKERRQ(ierr);
1814*c0ce001eSMatthew G. Knepley   ierr = PetscFVGetNumComponents(fvm, &pdim);CHKERRQ(ierr);
1815*c0ce001eSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
1816*c0ce001eSMatthew G. Knepley     (*riemann)(dim, pdim, fgeom[f].centroid, fgeom[f].normal, &uL[f*Nc], &uR[f*Nc], numConstants, constants, flux, rctx);
1817*c0ce001eSMatthew G. Knepley     for (d = 0; d < pdim; ++d) {
1818*c0ce001eSMatthew G. Knepley       fluxL[f*totDim+off+d] = flux[d] / neighborVol[f*2+0];
1819*c0ce001eSMatthew G. Knepley       fluxR[f*totDim+off+d] = flux[d] / neighborVol[f*2+1];
1820*c0ce001eSMatthew G. Knepley     }
1821*c0ce001eSMatthew G. Knepley   }
1822*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
1823*c0ce001eSMatthew G. Knepley }
1824*c0ce001eSMatthew G. Knepley 
1825*c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVInitialize_Upwind(PetscFV fvm)
1826*c0ce001eSMatthew G. Knepley {
1827*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1828*c0ce001eSMatthew G. Knepley   fvm->ops->setfromoptions          = NULL;
1829*c0ce001eSMatthew G. Knepley   fvm->ops->setup                   = PetscFVSetUp_Upwind;
1830*c0ce001eSMatthew G. Knepley   fvm->ops->view                    = PetscFVView_Upwind;
1831*c0ce001eSMatthew G. Knepley   fvm->ops->destroy                 = PetscFVDestroy_Upwind;
1832*c0ce001eSMatthew G. Knepley   fvm->ops->integraterhsfunction    = PetscFVIntegrateRHSFunction_Upwind;
1833*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
1834*c0ce001eSMatthew G. Knepley }
1835*c0ce001eSMatthew G. Knepley 
1836*c0ce001eSMatthew G. Knepley /*MC
1837*c0ce001eSMatthew G. Knepley   PETSCFVUPWIND = "upwind" - A PetscFV object
1838*c0ce001eSMatthew G. Knepley 
1839*c0ce001eSMatthew G. Knepley   Level: intermediate
1840*c0ce001eSMatthew G. Knepley 
1841*c0ce001eSMatthew G. Knepley .seealso: PetscFVType, PetscFVCreate(), PetscFVSetType()
1842*c0ce001eSMatthew G. Knepley M*/
1843*c0ce001eSMatthew G. Knepley 
1844*c0ce001eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFVCreate_Upwind(PetscFV fvm)
1845*c0ce001eSMatthew G. Knepley {
1846*c0ce001eSMatthew G. Knepley   PetscFV_Upwind *b;
1847*c0ce001eSMatthew G. Knepley   PetscErrorCode  ierr;
1848*c0ce001eSMatthew G. Knepley 
1849*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1850*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1);
1851*c0ce001eSMatthew G. Knepley   ierr      = PetscNewLog(fvm,&b);CHKERRQ(ierr);
1852*c0ce001eSMatthew G. Knepley   fvm->data = b;
1853*c0ce001eSMatthew G. Knepley 
1854*c0ce001eSMatthew G. Knepley   ierr = PetscFVInitialize_Upwind(fvm);CHKERRQ(ierr);
1855*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
1856*c0ce001eSMatthew G. Knepley }
1857*c0ce001eSMatthew G. Knepley 
1858*c0ce001eSMatthew G. Knepley #include <petscblaslapack.h>
1859*c0ce001eSMatthew G. Knepley 
1860*c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVDestroy_LeastSquares(PetscFV fvm)
1861*c0ce001eSMatthew G. Knepley {
1862*c0ce001eSMatthew G. Knepley   PetscFV_LeastSquares *ls = (PetscFV_LeastSquares *) fvm->data;
1863*c0ce001eSMatthew G. Knepley   PetscErrorCode        ierr;
1864*c0ce001eSMatthew G. Knepley 
1865*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1866*c0ce001eSMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject) fvm, "PetscFVLeastSquaresSetMaxFaces_C", NULL);CHKERRQ(ierr);
1867*c0ce001eSMatthew G. Knepley   ierr = PetscFree4(ls->B, ls->Binv, ls->tau, ls->work);CHKERRQ(ierr);
1868*c0ce001eSMatthew G. Knepley   ierr = PetscFree(ls);CHKERRQ(ierr);
1869*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
1870*c0ce001eSMatthew G. Knepley }
1871*c0ce001eSMatthew G. Knepley 
1872*c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVView_LeastSquares_Ascii(PetscFV fv, PetscViewer viewer)
1873*c0ce001eSMatthew G. Knepley {
1874*c0ce001eSMatthew G. Knepley   PetscInt          Nc, c;
1875*c0ce001eSMatthew G. Knepley   PetscViewerFormat format;
1876*c0ce001eSMatthew G. Knepley   PetscErrorCode    ierr;
1877*c0ce001eSMatthew G. Knepley 
1878*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1879*c0ce001eSMatthew G. Knepley   ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr);
1880*c0ce001eSMatthew G. Knepley   ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr);
1881*c0ce001eSMatthew G. Knepley   ierr = PetscViewerASCIIPrintf(viewer, "Finite Volume with Least Squares Reconstruction:\n");CHKERRQ(ierr);
1882*c0ce001eSMatthew G. Knepley   ierr = PetscViewerASCIIPrintf(viewer, "  num components: %d\n", Nc);CHKERRQ(ierr);
1883*c0ce001eSMatthew G. Knepley   for (c = 0; c < Nc; c++) {
1884*c0ce001eSMatthew G. Knepley     if (fv->componentNames[c]) {
1885*c0ce001eSMatthew G. Knepley       ierr = PetscViewerASCIIPrintf(viewer, "    component %d: %s\n", c, fv->componentNames[c]);CHKERRQ(ierr);
1886*c0ce001eSMatthew G. Knepley     }
1887*c0ce001eSMatthew G. Knepley   }
1888*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
1889*c0ce001eSMatthew G. Knepley }
1890*c0ce001eSMatthew G. Knepley 
1891*c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVView_LeastSquares(PetscFV fv, PetscViewer viewer)
1892*c0ce001eSMatthew G. Knepley {
1893*c0ce001eSMatthew G. Knepley   PetscBool      iascii;
1894*c0ce001eSMatthew G. Knepley   PetscErrorCode ierr;
1895*c0ce001eSMatthew G. Knepley 
1896*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1897*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(fv, PETSCFV_CLASSID, 1);
1898*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
1899*c0ce001eSMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
1900*c0ce001eSMatthew G. Knepley   if (iascii) {ierr = PetscFVView_LeastSquares_Ascii(fv, viewer);CHKERRQ(ierr);}
1901*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
1902*c0ce001eSMatthew G. Knepley }
1903*c0ce001eSMatthew G. Knepley 
1904*c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVSetUp_LeastSquares(PetscFV fvm)
1905*c0ce001eSMatthew G. Knepley {
1906*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1907*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
1908*c0ce001eSMatthew G. Knepley }
1909*c0ce001eSMatthew G. Knepley 
1910*c0ce001eSMatthew G. Knepley /* Overwrites A. Can only handle full-rank problems with m>=n */
1911*c0ce001eSMatthew G. Knepley static PetscErrorCode PetscFVLeastSquaresPseudoInverse_Static(PetscInt m,PetscInt mstride,PetscInt n,PetscScalar *A,PetscScalar *Ainv,PetscScalar *tau,PetscInt worksize,PetscScalar *work)
1912*c0ce001eSMatthew G. Knepley {
1913*c0ce001eSMatthew G. Knepley   PetscBool      debug = PETSC_FALSE;
1914*c0ce001eSMatthew G. Knepley   PetscErrorCode ierr;
1915*c0ce001eSMatthew G. Knepley   PetscBLASInt   M,N,K,lda,ldb,ldwork,info;
1916*c0ce001eSMatthew G. Knepley   PetscScalar    *R,*Q,*Aback,Alpha;
1917*c0ce001eSMatthew G. Knepley 
1918*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1919*c0ce001eSMatthew G. Knepley   if (debug) {
1920*c0ce001eSMatthew G. Knepley     ierr = PetscMalloc1(m*n,&Aback);CHKERRQ(ierr);
1921*c0ce001eSMatthew G. Knepley     ierr = PetscArraycpy(Aback,A,m*n);CHKERRQ(ierr);
1922*c0ce001eSMatthew G. Knepley   }
1923*c0ce001eSMatthew G. Knepley 
1924*c0ce001eSMatthew G. Knepley   ierr = PetscBLASIntCast(m,&M);CHKERRQ(ierr);
1925*c0ce001eSMatthew G. Knepley   ierr = PetscBLASIntCast(n,&N);CHKERRQ(ierr);
1926*c0ce001eSMatthew G. Knepley   ierr = PetscBLASIntCast(mstride,&lda);CHKERRQ(ierr);
1927*c0ce001eSMatthew G. Knepley   ierr = PetscBLASIntCast(worksize,&ldwork);CHKERRQ(ierr);
1928*c0ce001eSMatthew G. Knepley   ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
1929*c0ce001eSMatthew G. Knepley   LAPACKgeqrf_(&M,&N,A,&lda,tau,work,&ldwork,&info);
1930*c0ce001eSMatthew G. Knepley   ierr = PetscFPTrapPop();CHKERRQ(ierr);
1931*c0ce001eSMatthew G. Knepley   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xGEQRF error");
1932*c0ce001eSMatthew G. Knepley   R = A; /* Upper triangular part of A now contains R, the rest contains the elementary reflectors */
1933*c0ce001eSMatthew G. Knepley 
1934*c0ce001eSMatthew G. Knepley   /* Extract an explicit representation of Q */
1935*c0ce001eSMatthew G. Knepley   Q    = Ainv;
1936*c0ce001eSMatthew G. Knepley   ierr = PetscArraycpy(Q,A,mstride*n);CHKERRQ(ierr);
1937*c0ce001eSMatthew G. Knepley   K    = N;                     /* full rank */
1938*c0ce001eSMatthew G. Knepley   PetscStackCallBLAS("LAPACKorgqr",LAPACKorgqr_(&M,&N,&K,Q,&lda,tau,work,&ldwork,&info));
1939*c0ce001eSMatthew G. Knepley   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xORGQR/xUNGQR error");
1940*c0ce001eSMatthew G. Knepley 
1941*c0ce001eSMatthew G. Knepley   /* Compute A^{-T} = (R^{-1} Q^T)^T = Q R^{-T} */
1942*c0ce001eSMatthew G. Knepley   Alpha = 1.0;
1943*c0ce001eSMatthew G. Knepley   ldb   = lda;
1944*c0ce001eSMatthew G. Knepley   BLAStrsm_("Right","Upper","ConjugateTranspose","NotUnitTriangular",&M,&N,&Alpha,R,&lda,Q,&ldb);
1945*c0ce001eSMatthew G. Knepley   /* Ainv is Q, overwritten with inverse */
1946*c0ce001eSMatthew G. Knepley 
1947*c0ce001eSMatthew G. Knepley   if (debug) {                      /* Check that pseudo-inverse worked */
1948*c0ce001eSMatthew G. Knepley     PetscScalar  Beta = 0.0;
1949*c0ce001eSMatthew G. Knepley     PetscBLASInt ldc;
1950*c0ce001eSMatthew G. Knepley     K   = N;
1951*c0ce001eSMatthew G. Knepley     ldc = N;
1952*c0ce001eSMatthew G. Knepley     BLASgemm_("ConjugateTranspose","Normal",&N,&K,&M,&Alpha,Ainv,&lda,Aback,&ldb,&Beta,work,&ldc);
1953*c0ce001eSMatthew G. Knepley     ierr = PetscScalarView(n*n,work,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
1954*c0ce001eSMatthew G. Knepley     ierr = PetscFree(Aback);CHKERRQ(ierr);
1955*c0ce001eSMatthew G. Knepley   }
1956*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
1957*c0ce001eSMatthew G. Knepley }
1958*c0ce001eSMatthew G. Knepley 
1959*c0ce001eSMatthew G. Knepley /* Overwrites A. Can handle degenerate problems and m<n. */
1960*c0ce001eSMatthew G. Knepley static PetscErrorCode PetscFVLeastSquaresPseudoInverseSVD_Static(PetscInt m,PetscInt mstride,PetscInt n,PetscScalar *A,PetscScalar *Ainv,PetscScalar *tau,PetscInt worksize,PetscScalar *work)
1961*c0ce001eSMatthew G. Knepley {
1962*c0ce001eSMatthew G. Knepley   PetscBool      debug = PETSC_FALSE;
1963*c0ce001eSMatthew G. Knepley   PetscScalar   *Brhs, *Aback;
1964*c0ce001eSMatthew G. Knepley   PetscScalar   *tmpwork;
1965*c0ce001eSMatthew G. Knepley   PetscReal      rcond;
1966*c0ce001eSMatthew G. Knepley #if defined (PETSC_USE_COMPLEX)
1967*c0ce001eSMatthew G. Knepley   PetscInt       rworkSize;
1968*c0ce001eSMatthew G. Knepley   PetscReal     *rwork;
1969*c0ce001eSMatthew G. Knepley #endif
1970*c0ce001eSMatthew G. Knepley   PetscInt       i, j, maxmn;
1971*c0ce001eSMatthew G. Knepley   PetscBLASInt   M, N, lda, ldb, ldwork;
1972*c0ce001eSMatthew G. Knepley   PetscBLASInt   nrhs, irank, info;
1973*c0ce001eSMatthew G. Knepley   PetscErrorCode ierr;
1974*c0ce001eSMatthew G. Knepley 
1975*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1976*c0ce001eSMatthew G. Knepley   if (debug) {
1977*c0ce001eSMatthew G. Knepley     ierr = PetscMalloc1(m*n,&Aback);CHKERRQ(ierr);
1978*c0ce001eSMatthew G. Knepley     ierr = PetscArraycpy(Aback,A,m*n);CHKERRQ(ierr);
1979*c0ce001eSMatthew G. Knepley   }
1980*c0ce001eSMatthew G. Knepley 
1981*c0ce001eSMatthew G. Knepley   /* initialize to identity */
1982*c0ce001eSMatthew G. Knepley   tmpwork = Ainv;
1983*c0ce001eSMatthew G. Knepley   Brhs = work;
1984*c0ce001eSMatthew G. Knepley   maxmn = PetscMax(m,n);
1985*c0ce001eSMatthew G. Knepley   for (j=0; j<maxmn; j++) {
1986*c0ce001eSMatthew G. Knepley     for (i=0; i<maxmn; i++) Brhs[i + j*maxmn] = 1.0*(i == j);
1987*c0ce001eSMatthew G. Knepley   }
1988*c0ce001eSMatthew G. Knepley 
1989*c0ce001eSMatthew G. Knepley   ierr  = PetscBLASIntCast(m,&M);CHKERRQ(ierr);
1990*c0ce001eSMatthew G. Knepley   ierr  = PetscBLASIntCast(n,&N);CHKERRQ(ierr);
1991*c0ce001eSMatthew G. Knepley   ierr  = PetscBLASIntCast(mstride,&lda);CHKERRQ(ierr);
1992*c0ce001eSMatthew G. Knepley   ierr  = PetscBLASIntCast(maxmn,&ldb);CHKERRQ(ierr);
1993*c0ce001eSMatthew G. Knepley   ierr  = PetscBLASIntCast(worksize,&ldwork);CHKERRQ(ierr);
1994*c0ce001eSMatthew G. Knepley   rcond = -1;
1995*c0ce001eSMatthew G. Knepley   ierr  = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
1996*c0ce001eSMatthew G. Knepley   nrhs  = M;
1997*c0ce001eSMatthew G. Knepley #if defined(PETSC_USE_COMPLEX)
1998*c0ce001eSMatthew G. Knepley   rworkSize = 5 * PetscMin(M,N);
1999*c0ce001eSMatthew G. Knepley   ierr  = PetscMalloc1(rworkSize,&rwork);CHKERRQ(ierr);
2000*c0ce001eSMatthew G. Knepley   LAPACKgelss_(&M,&N,&nrhs,A,&lda,Brhs,&ldb, (PetscReal *) tau,&rcond,&irank,tmpwork,&ldwork,rwork,&info);
2001*c0ce001eSMatthew G. Knepley   ierr = PetscFPTrapPop();CHKERRQ(ierr);
2002*c0ce001eSMatthew G. Knepley   ierr = PetscFree(rwork);CHKERRQ(ierr);
2003*c0ce001eSMatthew G. Knepley #else
2004*c0ce001eSMatthew G. Knepley   nrhs  = M;
2005*c0ce001eSMatthew G. Knepley   LAPACKgelss_(&M,&N,&nrhs,A,&lda,Brhs,&ldb, (PetscReal *) tau,&rcond,&irank,tmpwork,&ldwork,&info);
2006*c0ce001eSMatthew G. Knepley   ierr = PetscFPTrapPop();CHKERRQ(ierr);
2007*c0ce001eSMatthew G. Knepley #endif
2008*c0ce001eSMatthew G. Knepley   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xGELSS error");
2009*c0ce001eSMatthew G. Knepley   /* The following check should be turned into a diagnostic as soon as someone wants to do this intentionally */
2010*c0ce001eSMatthew G. Knepley   if (irank < PetscMin(M,N)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Rank deficient least squares fit, indicates an isolated cell with two colinear points");
2011*c0ce001eSMatthew G. Knepley 
2012*c0ce001eSMatthew G. Knepley   /* Brhs shaped (M,nrhs) column-major coldim=mstride was overwritten by Ainv shaped (N,nrhs) column-major coldim=maxmn.
2013*c0ce001eSMatthew G. Knepley    * Here we transpose to (N,nrhs) row-major rowdim=mstride. */
2014*c0ce001eSMatthew G. Knepley   for (i=0; i<n; i++) {
2015*c0ce001eSMatthew G. Knepley     for (j=0; j<nrhs; j++) Ainv[i*mstride+j] = Brhs[i + j*maxmn];
2016*c0ce001eSMatthew G. Knepley   }
2017*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
2018*c0ce001eSMatthew G. Knepley }
2019*c0ce001eSMatthew G. Knepley 
2020*c0ce001eSMatthew G. Knepley #if 0
2021*c0ce001eSMatthew G. Knepley static PetscErrorCode PetscFVLeastSquaresDebugCell_Static(PetscFV fvm, PetscInt cell, DM dm, DM dmFace, PetscScalar *fgeom, DM dmCell, PetscScalar *cgeom)
2022*c0ce001eSMatthew G. Knepley {
2023*c0ce001eSMatthew G. Knepley   PetscReal       grad[2] = {0, 0};
2024*c0ce001eSMatthew G. Knepley   const PetscInt *faces;
2025*c0ce001eSMatthew G. Knepley   PetscInt        numFaces, f;
2026*c0ce001eSMatthew G. Knepley   PetscErrorCode  ierr;
2027*c0ce001eSMatthew G. Knepley 
2028*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
2029*c0ce001eSMatthew G. Knepley   ierr = DMPlexGetConeSize(dm, cell, &numFaces);CHKERRQ(ierr);
2030*c0ce001eSMatthew G. Knepley   ierr = DMPlexGetCone(dm, cell, &faces);CHKERRQ(ierr);
2031*c0ce001eSMatthew G. Knepley   for (f = 0; f < numFaces; ++f) {
2032*c0ce001eSMatthew G. Knepley     const PetscInt *fcells;
2033*c0ce001eSMatthew G. Knepley     const CellGeom *cg1;
2034*c0ce001eSMatthew G. Knepley     const FaceGeom *fg;
2035*c0ce001eSMatthew G. Knepley 
2036*c0ce001eSMatthew G. Knepley     ierr = DMPlexGetSupport(dm, faces[f], &fcells);CHKERRQ(ierr);
2037*c0ce001eSMatthew G. Knepley     ierr = DMPlexPointLocalRead(dmFace, faces[f], fgeom, &fg);CHKERRQ(ierr);
2038*c0ce001eSMatthew G. Knepley     for (i = 0; i < 2; ++i) {
2039*c0ce001eSMatthew G. Knepley       PetscScalar du;
2040*c0ce001eSMatthew G. Knepley 
2041*c0ce001eSMatthew G. Knepley       if (fcells[i] == c) continue;
2042*c0ce001eSMatthew G. Knepley       ierr = DMPlexPointLocalRead(dmCell, fcells[i], cgeom, &cg1);CHKERRQ(ierr);
2043*c0ce001eSMatthew G. Knepley       du   = cg1->centroid[0] + 3*cg1->centroid[1] - (cg->centroid[0] + 3*cg->centroid[1]);
2044*c0ce001eSMatthew G. Knepley       grad[0] += fg->grad[!i][0] * du;
2045*c0ce001eSMatthew G. Knepley       grad[1] += fg->grad[!i][1] * du;
2046*c0ce001eSMatthew G. Knepley     }
2047*c0ce001eSMatthew G. Knepley   }
2048*c0ce001eSMatthew G. Knepley   PetscPrintf(PETSC_COMM_SELF, "cell[%d] grad (%g, %g)\n", cell, grad[0], grad[1]);CHKERRQ(ierr);
2049*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
2050*c0ce001eSMatthew G. Knepley }
2051*c0ce001eSMatthew G. Knepley #endif
2052*c0ce001eSMatthew G. Knepley 
2053*c0ce001eSMatthew G. Knepley /*
2054*c0ce001eSMatthew G. Knepley   PetscFVComputeGradient - Compute the gradient reconstruction matrix for a given cell
2055*c0ce001eSMatthew G. Knepley 
2056*c0ce001eSMatthew G. Knepley   Input Parameters:
2057*c0ce001eSMatthew G. Knepley + fvm      - The PetscFV object
2058*c0ce001eSMatthew G. Knepley . numFaces - The number of cell faces which are not constrained
2059*c0ce001eSMatthew G. Knepley . dx       - The vector from the cell centroid to the neighboring cell centroid for each face
2060*c0ce001eSMatthew G. Knepley 
2061*c0ce001eSMatthew G. Knepley   Level: developer
2062*c0ce001eSMatthew G. Knepley 
2063*c0ce001eSMatthew G. Knepley .seealso: PetscFVCreate()
2064*c0ce001eSMatthew G. Knepley */
2065*c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVComputeGradient_LeastSquares(PetscFV fvm, PetscInt numFaces, const PetscScalar dx[], PetscScalar grad[])
2066*c0ce001eSMatthew G. Knepley {
2067*c0ce001eSMatthew G. Knepley   PetscFV_LeastSquares *ls       = (PetscFV_LeastSquares *) fvm->data;
2068*c0ce001eSMatthew G. Knepley   const PetscBool       useSVD   = PETSC_TRUE;
2069*c0ce001eSMatthew G. Knepley   const PetscInt        maxFaces = ls->maxFaces;
2070*c0ce001eSMatthew G. Knepley   PetscInt              dim, f, d;
2071*c0ce001eSMatthew G. Knepley   PetscErrorCode        ierr;
2072*c0ce001eSMatthew G. Knepley 
2073*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
2074*c0ce001eSMatthew G. Knepley   if (numFaces > maxFaces) {
2075*c0ce001eSMatthew G. Knepley     if (maxFaces < 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Reconstruction has not been initialized, call PetscFVLeastSquaresSetMaxFaces()");
2076*c0ce001eSMatthew G. Knepley     SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of input faces %D > %D maxfaces", numFaces, maxFaces);
2077*c0ce001eSMatthew G. Knepley   }
2078*c0ce001eSMatthew G. Knepley   ierr = PetscFVGetSpatialDimension(fvm, &dim);CHKERRQ(ierr);
2079*c0ce001eSMatthew G. Knepley   for (f = 0; f < numFaces; ++f) {
2080*c0ce001eSMatthew G. Knepley     for (d = 0; d < dim; ++d) ls->B[d*maxFaces+f] = dx[f*dim+d];
2081*c0ce001eSMatthew G. Knepley   }
2082*c0ce001eSMatthew G. Knepley   /* Overwrites B with garbage, returns Binv in row-major format */
2083*c0ce001eSMatthew G. Knepley   if (useSVD) {ierr = PetscFVLeastSquaresPseudoInverseSVD_Static(numFaces, maxFaces, dim, ls->B, ls->Binv, ls->tau, ls->workSize, ls->work);CHKERRQ(ierr);}
2084*c0ce001eSMatthew G. Knepley   else        {ierr = PetscFVLeastSquaresPseudoInverse_Static(numFaces, maxFaces, dim, ls->B, ls->Binv, ls->tau, ls->workSize, ls->work);CHKERRQ(ierr);}
2085*c0ce001eSMatthew G. Knepley   for (f = 0; f < numFaces; ++f) {
2086*c0ce001eSMatthew G. Knepley     for (d = 0; d < dim; ++d) grad[f*dim+d] = ls->Binv[d*maxFaces+f];
2087*c0ce001eSMatthew G. Knepley   }
2088*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
2089*c0ce001eSMatthew G. Knepley }
2090*c0ce001eSMatthew G. Knepley 
2091*c0ce001eSMatthew G. Knepley /*
2092*c0ce001eSMatthew G. Knepley   neighborVol[f*2+0] contains the left  geom
2093*c0ce001eSMatthew G. Knepley   neighborVol[f*2+1] contains the right geom
2094*c0ce001eSMatthew G. Knepley */
2095*c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVIntegrateRHSFunction_LeastSquares(PetscFV fvm, PetscDS prob, PetscInt field, PetscInt Nf, PetscFVFaceGeom *fgeom, PetscReal *neighborVol,
2096*c0ce001eSMatthew G. Knepley                                                         PetscScalar uL[], PetscScalar uR[], PetscScalar fluxL[], PetscScalar fluxR[])
2097*c0ce001eSMatthew G. Knepley {
2098*c0ce001eSMatthew G. Knepley   void             (*riemann)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscInt, const PetscScalar[], PetscScalar[], void *);
2099*c0ce001eSMatthew G. Knepley   void              *rctx;
2100*c0ce001eSMatthew G. Knepley   PetscScalar       *flux = fvm->fluxWork;
2101*c0ce001eSMatthew G. Knepley   const PetscScalar *constants;
2102*c0ce001eSMatthew G. Knepley   PetscInt           dim, numConstants, pdim, Nc, totDim, off, f, d;
2103*c0ce001eSMatthew G. Knepley   PetscErrorCode     ierr;
2104*c0ce001eSMatthew G. Knepley 
2105*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
2106*c0ce001eSMatthew G. Knepley   ierr = PetscDSGetTotalComponents(prob, &Nc);CHKERRQ(ierr);
2107*c0ce001eSMatthew G. Knepley   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
2108*c0ce001eSMatthew G. Knepley   ierr = PetscDSGetFieldOffset(prob, field, &off);CHKERRQ(ierr);
2109*c0ce001eSMatthew G. Knepley   ierr = PetscDSGetRiemannSolver(prob, field, &riemann);CHKERRQ(ierr);
2110*c0ce001eSMatthew G. Knepley   ierr = PetscDSGetContext(prob, field, &rctx);CHKERRQ(ierr);
2111*c0ce001eSMatthew G. Knepley   ierr = PetscDSGetConstants(prob, &numConstants, &constants);CHKERRQ(ierr);
2112*c0ce001eSMatthew G. Knepley   ierr = PetscFVGetSpatialDimension(fvm, &dim);CHKERRQ(ierr);
2113*c0ce001eSMatthew G. Knepley   ierr = PetscFVGetNumComponents(fvm, &pdim);CHKERRQ(ierr);
2114*c0ce001eSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
2115*c0ce001eSMatthew G. Knepley     (*riemann)(dim, pdim, fgeom[f].centroid, fgeom[f].normal, &uL[f*Nc], &uR[f*Nc], numConstants, constants, flux, rctx);
2116*c0ce001eSMatthew G. Knepley     for (d = 0; d < pdim; ++d) {
2117*c0ce001eSMatthew G. Knepley       fluxL[f*totDim+off+d] = flux[d] / neighborVol[f*2+0];
2118*c0ce001eSMatthew G. Knepley       fluxR[f*totDim+off+d] = flux[d] / neighborVol[f*2+1];
2119*c0ce001eSMatthew G. Knepley     }
2120*c0ce001eSMatthew G. Knepley   }
2121*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
2122*c0ce001eSMatthew G. Knepley }
2123*c0ce001eSMatthew G. Knepley 
2124*c0ce001eSMatthew G. Knepley static PetscErrorCode PetscFVLeastSquaresSetMaxFaces_LS(PetscFV fvm, PetscInt maxFaces)
2125*c0ce001eSMatthew G. Knepley {
2126*c0ce001eSMatthew G. Knepley   PetscFV_LeastSquares *ls = (PetscFV_LeastSquares *) fvm->data;
2127*c0ce001eSMatthew G. Knepley   PetscInt              dim, m, n, nrhs, minwork;
2128*c0ce001eSMatthew G. Knepley   PetscErrorCode        ierr;
2129*c0ce001eSMatthew G. Knepley 
2130*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
2131*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1);
2132*c0ce001eSMatthew G. Knepley   ierr = PetscFVGetSpatialDimension(fvm, &dim);CHKERRQ(ierr);
2133*c0ce001eSMatthew G. Knepley   ierr = PetscFree4(ls->B, ls->Binv, ls->tau, ls->work);CHKERRQ(ierr);
2134*c0ce001eSMatthew G. Knepley   ls->maxFaces = maxFaces;
2135*c0ce001eSMatthew G. Knepley   m       = ls->maxFaces;
2136*c0ce001eSMatthew G. Knepley   n       = dim;
2137*c0ce001eSMatthew G. Knepley   nrhs    = ls->maxFaces;
2138*c0ce001eSMatthew G. Knepley   minwork = 3*PetscMin(m,n) + PetscMax(2*PetscMin(m,n), PetscMax(PetscMax(m,n), nrhs)); /* required by LAPACK */
2139*c0ce001eSMatthew G. Knepley   ls->workSize = 5*minwork; /* We can afford to be extra generous */
2140*c0ce001eSMatthew G. Knepley   ierr = PetscMalloc4(ls->maxFaces*dim,&ls->B,ls->workSize,&ls->Binv,ls->maxFaces,&ls->tau,ls->workSize,&ls->work);CHKERRQ(ierr);
2141*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
2142*c0ce001eSMatthew G. Knepley }
2143*c0ce001eSMatthew G. Knepley 
2144*c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVInitialize_LeastSquares(PetscFV fvm)
2145*c0ce001eSMatthew G. Knepley {
2146*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
2147*c0ce001eSMatthew G. Knepley   fvm->ops->setfromoptions          = NULL;
2148*c0ce001eSMatthew G. Knepley   fvm->ops->setup                   = PetscFVSetUp_LeastSquares;
2149*c0ce001eSMatthew G. Knepley   fvm->ops->view                    = PetscFVView_LeastSquares;
2150*c0ce001eSMatthew G. Knepley   fvm->ops->destroy                 = PetscFVDestroy_LeastSquares;
2151*c0ce001eSMatthew G. Knepley   fvm->ops->computegradient         = PetscFVComputeGradient_LeastSquares;
2152*c0ce001eSMatthew G. Knepley   fvm->ops->integraterhsfunction    = PetscFVIntegrateRHSFunction_LeastSquares;
2153*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
2154*c0ce001eSMatthew G. Knepley }
2155*c0ce001eSMatthew G. Knepley 
2156*c0ce001eSMatthew G. Knepley /*MC
2157*c0ce001eSMatthew G. Knepley   PETSCFVLEASTSQUARES = "leastsquares" - A PetscFV object
2158*c0ce001eSMatthew G. Knepley 
2159*c0ce001eSMatthew G. Knepley   Level: intermediate
2160*c0ce001eSMatthew G. Knepley 
2161*c0ce001eSMatthew G. Knepley .seealso: PetscFVType, PetscFVCreate(), PetscFVSetType()
2162*c0ce001eSMatthew G. Knepley M*/
2163*c0ce001eSMatthew G. Knepley 
2164*c0ce001eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFVCreate_LeastSquares(PetscFV fvm)
2165*c0ce001eSMatthew G. Knepley {
2166*c0ce001eSMatthew G. Knepley   PetscFV_LeastSquares *ls;
2167*c0ce001eSMatthew G. Knepley   PetscErrorCode        ierr;
2168*c0ce001eSMatthew G. Knepley 
2169*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
2170*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1);
2171*c0ce001eSMatthew G. Knepley   ierr      = PetscNewLog(fvm, &ls);CHKERRQ(ierr);
2172*c0ce001eSMatthew G. Knepley   fvm->data = ls;
2173*c0ce001eSMatthew G. Knepley 
2174*c0ce001eSMatthew G. Knepley   ls->maxFaces = -1;
2175*c0ce001eSMatthew G. Knepley   ls->workSize = -1;
2176*c0ce001eSMatthew G. Knepley   ls->B        = NULL;
2177*c0ce001eSMatthew G. Knepley   ls->Binv     = NULL;
2178*c0ce001eSMatthew G. Knepley   ls->tau      = NULL;
2179*c0ce001eSMatthew G. Knepley   ls->work     = NULL;
2180*c0ce001eSMatthew G. Knepley 
2181*c0ce001eSMatthew G. Knepley   ierr = PetscFVSetComputeGradients(fvm, PETSC_TRUE);CHKERRQ(ierr);
2182*c0ce001eSMatthew G. Knepley   ierr = PetscFVInitialize_LeastSquares(fvm);CHKERRQ(ierr);
2183*c0ce001eSMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject) fvm, "PetscFVLeastSquaresSetMaxFaces_C", PetscFVLeastSquaresSetMaxFaces_LS);CHKERRQ(ierr);
2184*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
2185*c0ce001eSMatthew G. Knepley }
2186*c0ce001eSMatthew G. Knepley 
2187*c0ce001eSMatthew G. Knepley /*@
2188*c0ce001eSMatthew G. Knepley   PetscFVLeastSquaresSetMaxFaces - Set the maximum number of cell faces for gradient reconstruction
2189*c0ce001eSMatthew G. Knepley 
2190*c0ce001eSMatthew G. Knepley   Not collective
2191*c0ce001eSMatthew G. Knepley 
2192*c0ce001eSMatthew G. Knepley   Input parameters:
2193*c0ce001eSMatthew G. Knepley + fvm      - The PetscFV object
2194*c0ce001eSMatthew G. Knepley - maxFaces - The maximum number of cell faces
2195*c0ce001eSMatthew G. Knepley 
2196*c0ce001eSMatthew G. Knepley   Level: intermediate
2197*c0ce001eSMatthew G. Knepley 
2198*c0ce001eSMatthew G. Knepley .seealso: PetscFVCreate(), PETSCFVLEASTSQUARES
2199*c0ce001eSMatthew G. Knepley @*/
2200*c0ce001eSMatthew G. Knepley PetscErrorCode PetscFVLeastSquaresSetMaxFaces(PetscFV fvm, PetscInt maxFaces)
2201*c0ce001eSMatthew G. Knepley {
2202*c0ce001eSMatthew G. Knepley   PetscErrorCode ierr;
2203*c0ce001eSMatthew G. Knepley 
2204*c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
2205*c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1);
2206*c0ce001eSMatthew G. Knepley   ierr = PetscTryMethod(fvm, "PetscFVLeastSquaresSetMaxFaces_C", (PetscFV,PetscInt), (fvm,maxFaces));CHKERRQ(ierr);
2207*c0ce001eSMatthew G. Knepley   PetscFunctionReturn(0);
2208*c0ce001eSMatthew G. Knepley }
2209