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