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