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