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