xref: /petsc/src/dm/dt/fv/interface/fv.c (revision bcda9346efad4e5ba2d553af84eb238771ba1e25)
1c0ce001eSMatthew G. Knepley #include <petsc/private/petscfvimpl.h> /*I "petscfv.h" I*/
2412e9a14SMatthew G. Knepley #include <petscdmplex.h>
3012bc364SMatthew G. Knepley #include <petscdmplextransform.h>
4c0ce001eSMatthew G. Knepley #include <petscds.h>
5c0ce001eSMatthew G. Knepley 
6c0ce001eSMatthew G. Knepley PetscClassId PETSCLIMITER_CLASSID = 0;
7c0ce001eSMatthew G. Knepley 
8c0ce001eSMatthew G. Knepley PetscFunctionList PetscLimiterList              = NULL;
9c0ce001eSMatthew G. Knepley PetscBool         PetscLimiterRegisterAllCalled = PETSC_FALSE;
10c0ce001eSMatthew G. Knepley 
11c0ce001eSMatthew G. Knepley PetscBool  Limitercite       = PETSC_FALSE;
12c0ce001eSMatthew G. Knepley const char LimiterCitation[] = "@article{BergerAftosmisMurman2005,\n"
13c0ce001eSMatthew G. Knepley                                "  title   = {Analysis of slope limiters on irregular grids},\n"
14c0ce001eSMatthew G. Knepley                                "  journal = {AIAA paper},\n"
15c0ce001eSMatthew G. Knepley                                "  author  = {Marsha Berger and Michael J. Aftosmis and Scott M. Murman},\n"
16c0ce001eSMatthew G. Knepley                                "  volume  = {490},\n"
17c0ce001eSMatthew G. Knepley                                "  year    = {2005}\n}\n";
18c0ce001eSMatthew G. Knepley 
19c0ce001eSMatthew G. Knepley /*@C
20dce8aebaSBarry Smith   PetscLimiterRegister - Adds a new `PetscLimiter` implementation
21c0ce001eSMatthew G. Knepley 
22cc4c1da9SBarry Smith   Not Collective, No Fortran Support
23c0ce001eSMatthew G. Knepley 
24c0ce001eSMatthew G. Knepley   Input Parameters:
252fe279fdSBarry Smith + sname    - The name of a new user-defined creation routine
262fe279fdSBarry Smith - function - The creation routine
27c0ce001eSMatthew G. Knepley 
2860225df5SJacob Faibussowitsch   Example Usage:
29c0ce001eSMatthew G. Knepley .vb
30c0ce001eSMatthew G. Knepley     PetscLimiterRegister("my_lim", MyPetscLimiterCreate);
31c0ce001eSMatthew G. Knepley .ve
32c0ce001eSMatthew G. Knepley 
33dce8aebaSBarry Smith   Then, your `PetscLimiter` type can be chosen with the procedural interface via
34c0ce001eSMatthew G. Knepley .vb
35c0ce001eSMatthew G. Knepley     PetscLimiterCreate(MPI_Comm, PetscLimiter *);
36c0ce001eSMatthew G. Knepley     PetscLimiterSetType(PetscLimiter, "my_lim");
37c0ce001eSMatthew G. Knepley .ve
38c0ce001eSMatthew G. Knepley   or at runtime via the option
39c0ce001eSMatthew G. Knepley .vb
40c0ce001eSMatthew G. Knepley     -petsclimiter_type my_lim
41c0ce001eSMatthew G. Knepley .ve
42c0ce001eSMatthew G. Knepley 
43c0ce001eSMatthew G. Knepley   Level: advanced
44c0ce001eSMatthew G. Knepley 
45dce8aebaSBarry Smith   Note:
46dce8aebaSBarry Smith   `PetscLimiterRegister()` may be called multiple times to add several user-defined PetscLimiters
47c0ce001eSMatthew G. Knepley 
48dce8aebaSBarry Smith .seealso: `PetscLimiter`, `PetscLimiterType`, `PetscLimiterRegisterAll()`, `PetscLimiterRegisterDestroy()`
49c0ce001eSMatthew G. Knepley @*/
PetscLimiterRegister(const char sname[],PetscErrorCode (* function)(PetscLimiter))50d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscLimiterRegister(const char sname[], PetscErrorCode (*function)(PetscLimiter))
51d71ae5a4SJacob Faibussowitsch {
52c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
539566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&PetscLimiterList, sname, function));
543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
55c0ce001eSMatthew G. Knepley }
56c0ce001eSMatthew G. Knepley 
57cc4c1da9SBarry Smith /*@
58dce8aebaSBarry Smith   PetscLimiterSetType - Builds a `PetscLimiter` for a given `PetscLimiterType`
59c0ce001eSMatthew G. Knepley 
6020f4b53cSBarry Smith   Collective
61c0ce001eSMatthew G. Knepley 
62c0ce001eSMatthew G. Knepley   Input Parameters:
63dce8aebaSBarry Smith + lim  - The `PetscLimiter` object
64c0ce001eSMatthew G. Knepley - name - The kind of limiter
65c0ce001eSMatthew G. Knepley 
66c0ce001eSMatthew G. Knepley   Options Database Key:
67c0ce001eSMatthew G. Knepley . -petsclimiter_type <type> - Sets the PetscLimiter type; use -help for a list of available types
68c0ce001eSMatthew G. Knepley 
69c0ce001eSMatthew G. Knepley   Level: intermediate
70c0ce001eSMatthew G. Knepley 
71dce8aebaSBarry Smith .seealso: `PetscLimiter`, `PetscLimiterType`, `PetscLimiterGetType()`, `PetscLimiterCreate()`
72c0ce001eSMatthew G. Knepley @*/
PetscLimiterSetType(PetscLimiter lim,PetscLimiterType name)73d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscLimiterSetType(PetscLimiter lim, PetscLimiterType name)
74d71ae5a4SJacob Faibussowitsch {
75c0ce001eSMatthew G. Knepley   PetscErrorCode (*r)(PetscLimiter);
76c0ce001eSMatthew G. Knepley   PetscBool match;
77c0ce001eSMatthew G. Knepley 
78c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
79c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1);
809566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)lim, name, &match));
813ba16761SJacob Faibussowitsch   if (match) PetscFunctionReturn(PETSC_SUCCESS);
82c0ce001eSMatthew G. Knepley 
839566063dSJacob Faibussowitsch   PetscCall(PetscLimiterRegisterAll());
849566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListFind(PetscLimiterList, name, &r));
8528b400f6SJacob Faibussowitsch   PetscCheck(r, PetscObjectComm((PetscObject)lim), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscLimiter type: %s", name);
86c0ce001eSMatthew G. Knepley 
879927e4dfSBarry Smith   PetscTryTypeMethod(lim, destroy);
88c0ce001eSMatthew G. Knepley   lim->ops->destroy = NULL;
899927e4dfSBarry Smith 
909566063dSJacob Faibussowitsch   PetscCall((*r)(lim));
919566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)lim, name));
923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
93c0ce001eSMatthew G. Knepley }
94c0ce001eSMatthew G. Knepley 
95cc4c1da9SBarry Smith /*@
96dce8aebaSBarry Smith   PetscLimiterGetType - Gets the `PetscLimiterType` name (as a string) from the `PetscLimiter`.
97c0ce001eSMatthew G. Knepley 
98c0ce001eSMatthew G. Knepley   Not Collective
99c0ce001eSMatthew G. Knepley 
100c0ce001eSMatthew G. Knepley   Input Parameter:
101dce8aebaSBarry Smith . lim - The `PetscLimiter`
102c0ce001eSMatthew G. Knepley 
103c0ce001eSMatthew G. Knepley   Output Parameter:
104dce8aebaSBarry Smith . name - The `PetscLimiterType`
105c0ce001eSMatthew G. Knepley 
106c0ce001eSMatthew G. Knepley   Level: intermediate
107c0ce001eSMatthew G. Knepley 
108dce8aebaSBarry Smith .seealso: `PetscLimiter`, `PetscLimiterType`, `PetscLimiterSetType()`, `PetscLimiterCreate()`
109c0ce001eSMatthew G. Knepley @*/
PetscLimiterGetType(PetscLimiter lim,PetscLimiterType * name)110d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscLimiterGetType(PetscLimiter lim, PetscLimiterType *name)
111d71ae5a4SJacob Faibussowitsch {
112c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
113c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1);
1144f572ea9SToby Isaac   PetscAssertPointer(name, 2);
1159566063dSJacob Faibussowitsch   PetscCall(PetscLimiterRegisterAll());
116c0ce001eSMatthew G. Knepley   *name = ((PetscObject)lim)->type_name;
1173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
118c0ce001eSMatthew G. Knepley }
119c0ce001eSMatthew G. Knepley 
120ffeef943SBarry Smith /*@
121dce8aebaSBarry Smith   PetscLimiterViewFromOptions - View a `PetscLimiter` based on values in the options database
122fe2efc57SMark 
12320f4b53cSBarry Smith   Collective
124fe2efc57SMark 
125fe2efc57SMark   Input Parameters:
126dce8aebaSBarry Smith + A    - the `PetscLimiter` object to view
127dce8aebaSBarry Smith . obj  - Optional object that provides the options prefix to use
128dce8aebaSBarry Smith - name - command line option name
129fe2efc57SMark 
130fe2efc57SMark   Level: intermediate
131dce8aebaSBarry Smith 
132dce8aebaSBarry Smith .seealso: `PetscLimiter`, `PetscLimiterView()`, `PetscObjectViewFromOptions()`, `PetscLimiterCreate()`
133fe2efc57SMark @*/
PetscLimiterViewFromOptions(PetscLimiter A,PetscObject obj,const char name[])134d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscLimiterViewFromOptions(PetscLimiter A, PetscObject obj, const char name[])
135d71ae5a4SJacob Faibussowitsch {
136fe2efc57SMark   PetscFunctionBegin;
137fe2efc57SMark   PetscValidHeaderSpecific(A, PETSCLIMITER_CLASSID, 1);
1389566063dSJacob Faibussowitsch   PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name));
1393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
140fe2efc57SMark }
141fe2efc57SMark 
142ffeef943SBarry Smith /*@
143dce8aebaSBarry Smith   PetscLimiterView - Views a `PetscLimiter`
144c0ce001eSMatthew G. Knepley 
14520f4b53cSBarry Smith   Collective
146c0ce001eSMatthew G. Knepley 
147d8d19677SJose E. Roman   Input Parameters:
148dce8aebaSBarry Smith + lim - the `PetscLimiter` object to view
149c0ce001eSMatthew G. Knepley - v   - the viewer
150c0ce001eSMatthew G. Knepley 
15188f5f89eSMatthew G. Knepley   Level: beginner
152c0ce001eSMatthew G. Knepley 
153dce8aebaSBarry Smith .seealso: `PetscLimiter`, `PetscViewer`, `PetscLimiterDestroy()`, `PetscLimiterViewFromOptions()`
154c0ce001eSMatthew G. Knepley @*/
PetscLimiterView(PetscLimiter lim,PetscViewer v)155d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscLimiterView(PetscLimiter lim, PetscViewer v)
156d71ae5a4SJacob Faibussowitsch {
157c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
158c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1);
1599566063dSJacob Faibussowitsch   if (!v) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)lim), &v));
160dbbe0bcdSBarry Smith   PetscTryTypeMethod(lim, view, v);
1613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
162c0ce001eSMatthew G. Knepley }
163c0ce001eSMatthew G. Knepley 
164c0ce001eSMatthew G. Knepley /*@
165dce8aebaSBarry Smith   PetscLimiterSetFromOptions - sets parameters in a `PetscLimiter` from the options database
166c0ce001eSMatthew G. Knepley 
16720f4b53cSBarry Smith   Collective
168c0ce001eSMatthew G. Knepley 
169c0ce001eSMatthew G. Knepley   Input Parameter:
170dce8aebaSBarry Smith . lim - the `PetscLimiter` object to set options for
171c0ce001eSMatthew G. Knepley 
17288f5f89eSMatthew G. Knepley   Level: intermediate
173c0ce001eSMatthew G. Knepley 
174a94f484eSPierre Jolivet .seealso: `PetscLimiter`, `PetscLimiterView()`
175c0ce001eSMatthew G. Knepley @*/
PetscLimiterSetFromOptions(PetscLimiter lim)176d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscLimiterSetFromOptions(PetscLimiter lim)
177d71ae5a4SJacob Faibussowitsch {
178c0ce001eSMatthew G. Knepley   const char *defaultType;
179c0ce001eSMatthew G. Knepley   char        name[256];
180c0ce001eSMatthew G. Knepley   PetscBool   flg;
181c0ce001eSMatthew G. Knepley 
182c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
183c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1);
184c0ce001eSMatthew G. Knepley   if (!((PetscObject)lim)->type_name) defaultType = PETSCLIMITERSIN;
185c0ce001eSMatthew G. Knepley   else defaultType = ((PetscObject)lim)->type_name;
1869566063dSJacob Faibussowitsch   PetscCall(PetscLimiterRegisterAll());
187c0ce001eSMatthew G. Knepley 
188d0609cedSBarry Smith   PetscObjectOptionsBegin((PetscObject)lim);
1899566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFList("-petsclimiter_type", "Finite volume slope limiter", "PetscLimiterSetType", PetscLimiterList, defaultType, name, 256, &flg));
190c0ce001eSMatthew G. Knepley   if (flg) {
1919566063dSJacob Faibussowitsch     PetscCall(PetscLimiterSetType(lim, name));
192c0ce001eSMatthew G. Knepley   } else if (!((PetscObject)lim)->type_name) {
1939566063dSJacob Faibussowitsch     PetscCall(PetscLimiterSetType(lim, defaultType));
194c0ce001eSMatthew G. Knepley   }
195dbbe0bcdSBarry Smith   PetscTryTypeMethod(lim, setfromoptions);
196c0ce001eSMatthew G. Knepley   /* process any options handlers added with PetscObjectAddOptionsHandler() */
197dbbe0bcdSBarry Smith   PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)lim, PetscOptionsObject));
198d0609cedSBarry Smith   PetscOptionsEnd();
1999566063dSJacob Faibussowitsch   PetscCall(PetscLimiterViewFromOptions(lim, NULL, "-petsclimiter_view"));
2003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
201c0ce001eSMatthew G. Knepley }
202c0ce001eSMatthew G. Knepley 
203cc4c1da9SBarry Smith /*@
204dce8aebaSBarry Smith   PetscLimiterSetUp - Construct data structures for the `PetscLimiter`
205c0ce001eSMatthew G. Knepley 
20620f4b53cSBarry Smith   Collective
207c0ce001eSMatthew G. Knepley 
208c0ce001eSMatthew G. Knepley   Input Parameter:
209dce8aebaSBarry Smith . lim - the `PetscLimiter` object to setup
210c0ce001eSMatthew G. Knepley 
21188f5f89eSMatthew G. Knepley   Level: intermediate
212c0ce001eSMatthew G. Knepley 
213a94f484eSPierre Jolivet .seealso: `PetscLimiter`, `PetscLimiterView()`, `PetscLimiterDestroy()`
214c0ce001eSMatthew G. Knepley @*/
PetscLimiterSetUp(PetscLimiter lim)215d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscLimiterSetUp(PetscLimiter lim)
216d71ae5a4SJacob Faibussowitsch {
217c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
218c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1);
219dbbe0bcdSBarry Smith   PetscTryTypeMethod(lim, setup);
2203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
221c0ce001eSMatthew G. Knepley }
222c0ce001eSMatthew G. Knepley 
223c0ce001eSMatthew G. Knepley /*@
224dce8aebaSBarry Smith   PetscLimiterDestroy - Destroys a `PetscLimiter` object
225c0ce001eSMatthew G. Knepley 
22620f4b53cSBarry Smith   Collective
227c0ce001eSMatthew G. Knepley 
228c0ce001eSMatthew G. Knepley   Input Parameter:
229dce8aebaSBarry Smith . lim - the `PetscLimiter` object to destroy
230c0ce001eSMatthew G. Knepley 
23188f5f89eSMatthew G. Knepley   Level: beginner
232c0ce001eSMatthew G. Knepley 
233dce8aebaSBarry Smith .seealso: `PetscLimiter`, `PetscLimiterView()`
234c0ce001eSMatthew G. Knepley @*/
PetscLimiterDestroy(PetscLimiter * lim)235d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscLimiterDestroy(PetscLimiter *lim)
236d71ae5a4SJacob Faibussowitsch {
237c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
2383ba16761SJacob Faibussowitsch   if (!*lim) PetscFunctionReturn(PETSC_SUCCESS);
239f4f49eeaSPierre Jolivet   PetscValidHeaderSpecific(*lim, PETSCLIMITER_CLASSID, 1);
240c0ce001eSMatthew G. Knepley 
241f4f49eeaSPierre Jolivet   if (--((PetscObject)*lim)->refct > 0) {
2429371c9d4SSatish Balay     *lim = NULL;
2433ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
2449371c9d4SSatish Balay   }
245f4f49eeaSPierre Jolivet   ((PetscObject)*lim)->refct = 0;
246c0ce001eSMatthew G. Knepley 
247f4f49eeaSPierre Jolivet   PetscTryTypeMethod(*lim, destroy);
2489566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(lim));
2493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
250c0ce001eSMatthew G. Knepley }
251c0ce001eSMatthew G. Knepley 
252c0ce001eSMatthew G. Knepley /*@
253dce8aebaSBarry Smith   PetscLimiterCreate - Creates an empty `PetscLimiter` object. The type can then be set with `PetscLimiterSetType()`.
254c0ce001eSMatthew G. Knepley 
255c0ce001eSMatthew G. Knepley   Collective
256c0ce001eSMatthew G. Knepley 
257c0ce001eSMatthew G. Knepley   Input Parameter:
258dce8aebaSBarry Smith . comm - The communicator for the `PetscLimiter` object
259c0ce001eSMatthew G. Knepley 
260c0ce001eSMatthew G. Knepley   Output Parameter:
261dce8aebaSBarry Smith . lim - The `PetscLimiter` object
262c0ce001eSMatthew G. Knepley 
263c0ce001eSMatthew G. Knepley   Level: beginner
264c0ce001eSMatthew G. Knepley 
26560225df5SJacob Faibussowitsch .seealso: `PetscLimiter`, `PetscLimiterType`, `PetscLimiterSetType()`, `PETSCLIMITERSIN`
266c0ce001eSMatthew G. Knepley @*/
PetscLimiterCreate(MPI_Comm comm,PetscLimiter * lim)267d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscLimiterCreate(MPI_Comm comm, PetscLimiter *lim)
268d71ae5a4SJacob Faibussowitsch {
269c0ce001eSMatthew G. Knepley   PetscLimiter l;
270c0ce001eSMatthew G. Knepley 
271c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
2724f572ea9SToby Isaac   PetscAssertPointer(lim, 2);
2739566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(LimiterCitation, &Limitercite));
2749566063dSJacob Faibussowitsch   PetscCall(PetscFVInitializePackage());
275c0ce001eSMatthew G. Knepley 
2769566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(l, PETSCLIMITER_CLASSID, "PetscLimiter", "Finite Volume Slope Limiter", "PetscLimiter", comm, PetscLimiterDestroy, PetscLimiterView));
277c0ce001eSMatthew G. Knepley 
278c0ce001eSMatthew G. Knepley   *lim = l;
2793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
280c0ce001eSMatthew G. Knepley }
281c0ce001eSMatthew G. Knepley 
28288f5f89eSMatthew G. Knepley /*@
28388f5f89eSMatthew G. Knepley   PetscLimiterLimit - Limit the flux
28488f5f89eSMatthew G. Knepley 
28588f5f89eSMatthew G. Knepley   Input Parameters:
286dce8aebaSBarry Smith + lim  - The `PetscLimiter`
28788f5f89eSMatthew G. Knepley - flim - The input field
28888f5f89eSMatthew G. Knepley 
28988f5f89eSMatthew G. Knepley   Output Parameter:
29088f5f89eSMatthew G. Knepley . phi - The limited field
29188f5f89eSMatthew G. Knepley 
29288f5f89eSMatthew G. Knepley   Level: beginner
29388f5f89eSMatthew G. Knepley 
294dce8aebaSBarry Smith   Note:
295dce8aebaSBarry Smith   Limiters given in symmetric form following Berger, Aftosmis, and Murman 2005
296dce8aebaSBarry Smith .vb
297dce8aebaSBarry Smith  The classical flux-limited formulation is psi(r) where
298dce8aebaSBarry Smith 
299dce8aebaSBarry Smith  r = (u[0] - u[-1]) / (u[1] - u[0])
300dce8aebaSBarry Smith 
301dce8aebaSBarry Smith  The second order TVD region is bounded by
302dce8aebaSBarry Smith 
303dce8aebaSBarry Smith  psi_minmod(r) = min(r,1)      and        psi_superbee(r) = min(2, 2r, max(1,r))
304dce8aebaSBarry Smith 
305dce8aebaSBarry Smith  where all limiters are implicitly clipped to be non-negative. A more convenient slope-limited form is psi(r) =
306dce8aebaSBarry Smith  phi(r)(r+1)/2 in which the reconstructed interface values are
307dce8aebaSBarry Smith 
308dce8aebaSBarry Smith  u(v) = u[0] + phi(r) (grad u)[0] v
309dce8aebaSBarry Smith 
310dce8aebaSBarry Smith  where v is the vector from centroid to quadrature point. In these variables, the usual limiters become
311dce8aebaSBarry Smith 
312dce8aebaSBarry Smith  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))
313dce8aebaSBarry Smith 
314dce8aebaSBarry Smith  For a nicer symmetric formulation, rewrite in terms of
315dce8aebaSBarry Smith 
316dce8aebaSBarry Smith  f = (u[0] - u[-1]) / (u[1] - u[-1])
317dce8aebaSBarry Smith 
318dce8aebaSBarry Smith  where r(f) = f/(1-f). Not that r(1-f) = (1-f)/f = 1/r(f) so the symmetry condition
319dce8aebaSBarry Smith 
320dce8aebaSBarry Smith  phi(r) = phi(1/r)
321dce8aebaSBarry Smith 
322dce8aebaSBarry Smith  becomes
323dce8aebaSBarry Smith 
324dce8aebaSBarry Smith  w(f) = w(1-f).
325dce8aebaSBarry Smith 
326dce8aebaSBarry Smith  The limiters below implement this final form w(f). The reference methods are
327dce8aebaSBarry Smith 
328dce8aebaSBarry Smith  w_minmod(f) = 2 min(f,(1-f))             w_superbee(r) = 4 min((1-f), f)
329dce8aebaSBarry Smith .ve
330dce8aebaSBarry Smith 
331dce8aebaSBarry Smith .seealso: `PetscLimiter`, `PetscLimiterType`, `PetscLimiterSetType()`, `PetscLimiterCreate()`
33288f5f89eSMatthew G. Knepley @*/
PetscLimiterLimit(PetscLimiter lim,PetscReal flim,PetscReal * phi)333d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscLimiterLimit(PetscLimiter lim, PetscReal flim, PetscReal *phi)
334d71ae5a4SJacob Faibussowitsch {
335c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
336c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1);
3374f572ea9SToby Isaac   PetscAssertPointer(phi, 3);
338dbbe0bcdSBarry Smith   PetscUseTypeMethod(lim, limit, flim, phi);
3393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
340c0ce001eSMatthew G. Knepley }
341c0ce001eSMatthew G. Knepley 
PetscLimiterDestroy_Sin(PetscLimiter lim)342d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscLimiterDestroy_Sin(PetscLimiter lim)
343d71ae5a4SJacob Faibussowitsch {
344c0ce001eSMatthew G. Knepley   PetscLimiter_Sin *l = (PetscLimiter_Sin *)lim->data;
345c0ce001eSMatthew G. Knepley 
346c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
3479566063dSJacob Faibussowitsch   PetscCall(PetscFree(l));
3483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
349c0ce001eSMatthew G. Knepley }
350c0ce001eSMatthew G. Knepley 
PetscLimiterView_Sin_Ascii(PetscLimiter lim,PetscViewer viewer)351d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscLimiterView_Sin_Ascii(PetscLimiter lim, PetscViewer viewer)
352d71ae5a4SJacob Faibussowitsch {
353c0ce001eSMatthew G. Knepley   PetscViewerFormat format;
354c0ce001eSMatthew G. Knepley 
355c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
3569566063dSJacob Faibussowitsch   PetscCall(PetscViewerGetFormat(viewer, &format));
3579566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Sin Slope Limiter:\n"));
3583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
359c0ce001eSMatthew G. Knepley }
360c0ce001eSMatthew G. Knepley 
PetscLimiterView_Sin(PetscLimiter lim,PetscViewer viewer)361d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscLimiterView_Sin(PetscLimiter lim, PetscViewer viewer)
362d71ae5a4SJacob Faibussowitsch {
363*9f196a02SMartin Diehl   PetscBool isascii;
364c0ce001eSMatthew G. Knepley 
365c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
366c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1);
367c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
368*9f196a02SMartin Diehl   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
369*9f196a02SMartin Diehl   if (isascii) PetscCall(PetscLimiterView_Sin_Ascii(lim, viewer));
3703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
371c0ce001eSMatthew G. Knepley }
372c0ce001eSMatthew G. Knepley 
PetscLimiterLimit_Sin(PetscLimiter lim,PetscReal f,PetscReal * phi)373d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscLimiterLimit_Sin(PetscLimiter lim, PetscReal f, PetscReal *phi)
374d71ae5a4SJacob Faibussowitsch {
375c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
376c0ce001eSMatthew G. Knepley   *phi = PetscSinReal(PETSC_PI * PetscMax(0, PetscMin(f, 1)));
3773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
378c0ce001eSMatthew G. Knepley }
379c0ce001eSMatthew G. Knepley 
PetscLimiterInitialize_Sin(PetscLimiter lim)380d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscLimiterInitialize_Sin(PetscLimiter lim)
381d71ae5a4SJacob Faibussowitsch {
382c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
383c0ce001eSMatthew G. Knepley   lim->ops->view    = PetscLimiterView_Sin;
384c0ce001eSMatthew G. Knepley   lim->ops->destroy = PetscLimiterDestroy_Sin;
385c0ce001eSMatthew G. Knepley   lim->ops->limit   = PetscLimiterLimit_Sin;
3863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
387c0ce001eSMatthew G. Knepley }
388c0ce001eSMatthew G. Knepley 
389c0ce001eSMatthew G. Knepley /*MC
390dce8aebaSBarry Smith   PETSCLIMITERSIN = "sin" - A `PetscLimiter` implementation
391c0ce001eSMatthew G. Knepley 
392c0ce001eSMatthew G. Knepley   Level: intermediate
393c0ce001eSMatthew G. Knepley 
394dce8aebaSBarry Smith .seealso: `PetscLimiter`, `PetscLimiterType`, `PetscLimiterCreate()`, `PetscLimiterSetType()`
395c0ce001eSMatthew G. Knepley M*/
396c0ce001eSMatthew G. Knepley 
PetscLimiterCreate_Sin(PetscLimiter lim)397d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscLimiterCreate_Sin(PetscLimiter lim)
398d71ae5a4SJacob Faibussowitsch {
399c0ce001eSMatthew G. Knepley   PetscLimiter_Sin *l;
400c0ce001eSMatthew G. Knepley 
401c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
402c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1);
4034dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&l));
404c0ce001eSMatthew G. Knepley   lim->data = l;
405c0ce001eSMatthew G. Knepley 
4069566063dSJacob Faibussowitsch   PetscCall(PetscLimiterInitialize_Sin(lim));
4073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
408c0ce001eSMatthew G. Knepley }
409c0ce001eSMatthew G. Knepley 
PetscLimiterDestroy_Zero(PetscLimiter lim)410d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscLimiterDestroy_Zero(PetscLimiter lim)
411d71ae5a4SJacob Faibussowitsch {
412c0ce001eSMatthew G. Knepley   PetscLimiter_Zero *l = (PetscLimiter_Zero *)lim->data;
413c0ce001eSMatthew G. Knepley 
414c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
4159566063dSJacob Faibussowitsch   PetscCall(PetscFree(l));
4163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
417c0ce001eSMatthew G. Knepley }
418c0ce001eSMatthew G. Knepley 
PetscLimiterView_Zero_Ascii(PetscLimiter lim,PetscViewer viewer)419d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscLimiterView_Zero_Ascii(PetscLimiter lim, PetscViewer viewer)
420d71ae5a4SJacob Faibussowitsch {
421c0ce001eSMatthew G. Knepley   PetscViewerFormat format;
422c0ce001eSMatthew G. Knepley 
423c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
4249566063dSJacob Faibussowitsch   PetscCall(PetscViewerGetFormat(viewer, &format));
4259566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Zero Slope Limiter:\n"));
4263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
427c0ce001eSMatthew G. Knepley }
428c0ce001eSMatthew G. Knepley 
PetscLimiterView_Zero(PetscLimiter lim,PetscViewer viewer)429d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscLimiterView_Zero(PetscLimiter lim, PetscViewer viewer)
430d71ae5a4SJacob Faibussowitsch {
431*9f196a02SMartin Diehl   PetscBool isascii;
432c0ce001eSMatthew G. Knepley 
433c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
434c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1);
435c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
436*9f196a02SMartin Diehl   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
437*9f196a02SMartin Diehl   if (isascii) PetscCall(PetscLimiterView_Zero_Ascii(lim, viewer));
4383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
439c0ce001eSMatthew G. Knepley }
440c0ce001eSMatthew G. Knepley 
PetscLimiterLimit_Zero(PetscLimiter lim,PetscReal f,PetscReal * phi)441d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscLimiterLimit_Zero(PetscLimiter lim, PetscReal f, PetscReal *phi)
442d71ae5a4SJacob Faibussowitsch {
443c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
444c0ce001eSMatthew G. Knepley   *phi = 0.0;
4453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
446c0ce001eSMatthew G. Knepley }
447c0ce001eSMatthew G. Knepley 
PetscLimiterInitialize_Zero(PetscLimiter lim)448d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscLimiterInitialize_Zero(PetscLimiter lim)
449d71ae5a4SJacob Faibussowitsch {
450c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
451c0ce001eSMatthew G. Knepley   lim->ops->view    = PetscLimiterView_Zero;
452c0ce001eSMatthew G. Knepley   lim->ops->destroy = PetscLimiterDestroy_Zero;
453c0ce001eSMatthew G. Knepley   lim->ops->limit   = PetscLimiterLimit_Zero;
4543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
455c0ce001eSMatthew G. Knepley }
456c0ce001eSMatthew G. Knepley 
457c0ce001eSMatthew G. Knepley /*MC
458dce8aebaSBarry Smith   PETSCLIMITERZERO = "zero" - A simple `PetscLimiter` implementation
459c0ce001eSMatthew G. Knepley 
460c0ce001eSMatthew G. Knepley   Level: intermediate
461c0ce001eSMatthew G. Knepley 
462dce8aebaSBarry Smith .seealso: `PetscLimiter`, `PetscLimiterType`, `PetscLimiterCreate()`, `PetscLimiterSetType()`
463c0ce001eSMatthew G. Knepley M*/
464c0ce001eSMatthew G. Knepley 
PetscLimiterCreate_Zero(PetscLimiter lim)465d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscLimiterCreate_Zero(PetscLimiter lim)
466d71ae5a4SJacob Faibussowitsch {
467c0ce001eSMatthew G. Knepley   PetscLimiter_Zero *l;
468c0ce001eSMatthew G. Knepley 
469c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
470c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1);
4714dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&l));
472c0ce001eSMatthew G. Knepley   lim->data = l;
473c0ce001eSMatthew G. Knepley 
4749566063dSJacob Faibussowitsch   PetscCall(PetscLimiterInitialize_Zero(lim));
4753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
476c0ce001eSMatthew G. Knepley }
477c0ce001eSMatthew G. Knepley 
PetscLimiterDestroy_None(PetscLimiter lim)478d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscLimiterDestroy_None(PetscLimiter lim)
479d71ae5a4SJacob Faibussowitsch {
480c0ce001eSMatthew G. Knepley   PetscLimiter_None *l = (PetscLimiter_None *)lim->data;
481c0ce001eSMatthew G. Knepley 
482c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
4839566063dSJacob Faibussowitsch   PetscCall(PetscFree(l));
4843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
485c0ce001eSMatthew G. Knepley }
486c0ce001eSMatthew G. Knepley 
PetscLimiterView_None_Ascii(PetscLimiter lim,PetscViewer viewer)487d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscLimiterView_None_Ascii(PetscLimiter lim, PetscViewer viewer)
488d71ae5a4SJacob Faibussowitsch {
489c0ce001eSMatthew G. Knepley   PetscViewerFormat format;
490c0ce001eSMatthew G. Knepley 
491c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
4929566063dSJacob Faibussowitsch   PetscCall(PetscViewerGetFormat(viewer, &format));
4939566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "None Slope Limiter:\n"));
4943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
495c0ce001eSMatthew G. Knepley }
496c0ce001eSMatthew G. Knepley 
PetscLimiterView_None(PetscLimiter lim,PetscViewer viewer)497d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscLimiterView_None(PetscLimiter lim, PetscViewer viewer)
498d71ae5a4SJacob Faibussowitsch {
499*9f196a02SMartin Diehl   PetscBool isascii;
500c0ce001eSMatthew G. Knepley 
501c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
502c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1);
503c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
504*9f196a02SMartin Diehl   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
505*9f196a02SMartin Diehl   if (isascii) PetscCall(PetscLimiterView_None_Ascii(lim, viewer));
5063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
507c0ce001eSMatthew G. Knepley }
508c0ce001eSMatthew G. Knepley 
PetscLimiterLimit_None(PetscLimiter lim,PetscReal f,PetscReal * phi)509d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscLimiterLimit_None(PetscLimiter lim, PetscReal f, PetscReal *phi)
510d71ae5a4SJacob Faibussowitsch {
511c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
512c0ce001eSMatthew G. Knepley   *phi = 1.0;
5133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
514c0ce001eSMatthew G. Knepley }
515c0ce001eSMatthew G. Knepley 
PetscLimiterInitialize_None(PetscLimiter lim)516d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscLimiterInitialize_None(PetscLimiter lim)
517d71ae5a4SJacob Faibussowitsch {
518c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
519c0ce001eSMatthew G. Knepley   lim->ops->view    = PetscLimiterView_None;
520c0ce001eSMatthew G. Knepley   lim->ops->destroy = PetscLimiterDestroy_None;
521c0ce001eSMatthew G. Knepley   lim->ops->limit   = PetscLimiterLimit_None;
5223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
523c0ce001eSMatthew G. Knepley }
524c0ce001eSMatthew G. Knepley 
525c0ce001eSMatthew G. Knepley /*MC
526dce8aebaSBarry Smith   PETSCLIMITERNONE = "none" - A trivial `PetscLimiter` implementation
527c0ce001eSMatthew G. Knepley 
528c0ce001eSMatthew G. Knepley   Level: intermediate
529c0ce001eSMatthew G. Knepley 
530dce8aebaSBarry Smith .seealso: `PetscLimiter`, `PetscLimiterType`, `PetscLimiterCreate()`, `PetscLimiterSetType()`
531c0ce001eSMatthew G. Knepley M*/
532c0ce001eSMatthew G. Knepley 
PetscLimiterCreate_None(PetscLimiter lim)533d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscLimiterCreate_None(PetscLimiter lim)
534d71ae5a4SJacob Faibussowitsch {
535c0ce001eSMatthew G. Knepley   PetscLimiter_None *l;
536c0ce001eSMatthew G. Knepley 
537c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
538c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1);
5394dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&l));
540c0ce001eSMatthew G. Knepley   lim->data = l;
541c0ce001eSMatthew G. Knepley 
5429566063dSJacob Faibussowitsch   PetscCall(PetscLimiterInitialize_None(lim));
5433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
544c0ce001eSMatthew G. Knepley }
545c0ce001eSMatthew G. Knepley 
PetscLimiterDestroy_Minmod(PetscLimiter lim)546d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscLimiterDestroy_Minmod(PetscLimiter lim)
547d71ae5a4SJacob Faibussowitsch {
548c0ce001eSMatthew G. Knepley   PetscLimiter_Minmod *l = (PetscLimiter_Minmod *)lim->data;
549c0ce001eSMatthew G. Knepley 
550c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
5519566063dSJacob Faibussowitsch   PetscCall(PetscFree(l));
5523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
553c0ce001eSMatthew G. Knepley }
554c0ce001eSMatthew G. Knepley 
PetscLimiterView_Minmod_Ascii(PetscLimiter lim,PetscViewer viewer)555d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscLimiterView_Minmod_Ascii(PetscLimiter lim, PetscViewer viewer)
556d71ae5a4SJacob Faibussowitsch {
557c0ce001eSMatthew G. Knepley   PetscViewerFormat format;
558c0ce001eSMatthew G. Knepley 
559c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
5609566063dSJacob Faibussowitsch   PetscCall(PetscViewerGetFormat(viewer, &format));
5619566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Minmod Slope Limiter:\n"));
5623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
563c0ce001eSMatthew G. Knepley }
564c0ce001eSMatthew G. Knepley 
PetscLimiterView_Minmod(PetscLimiter lim,PetscViewer viewer)565d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscLimiterView_Minmod(PetscLimiter lim, PetscViewer viewer)
566d71ae5a4SJacob Faibussowitsch {
567*9f196a02SMartin Diehl   PetscBool isascii;
568c0ce001eSMatthew G. Knepley 
569c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
570c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1);
571c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
572*9f196a02SMartin Diehl   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
573*9f196a02SMartin Diehl   if (isascii) PetscCall(PetscLimiterView_Minmod_Ascii(lim, viewer));
5743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
575c0ce001eSMatthew G. Knepley }
576c0ce001eSMatthew G. Knepley 
PetscLimiterLimit_Minmod(PetscLimiter lim,PetscReal f,PetscReal * phi)577d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscLimiterLimit_Minmod(PetscLimiter lim, PetscReal f, PetscReal *phi)
578d71ae5a4SJacob Faibussowitsch {
579c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
580c0ce001eSMatthew G. Knepley   *phi = 2 * PetscMax(0, PetscMin(f, 1 - f));
5813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
582c0ce001eSMatthew G. Knepley }
583c0ce001eSMatthew G. Knepley 
PetscLimiterInitialize_Minmod(PetscLimiter lim)584d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscLimiterInitialize_Minmod(PetscLimiter lim)
585d71ae5a4SJacob Faibussowitsch {
586c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
587c0ce001eSMatthew G. Knepley   lim->ops->view    = PetscLimiterView_Minmod;
588c0ce001eSMatthew G. Knepley   lim->ops->destroy = PetscLimiterDestroy_Minmod;
589c0ce001eSMatthew G. Knepley   lim->ops->limit   = PetscLimiterLimit_Minmod;
5903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
591c0ce001eSMatthew G. Knepley }
592c0ce001eSMatthew G. Knepley 
593c0ce001eSMatthew G. Knepley /*MC
594dce8aebaSBarry Smith   PETSCLIMITERMINMOD = "minmod" - A `PetscLimiter` implementation
595c0ce001eSMatthew G. Knepley 
596c0ce001eSMatthew G. Knepley   Level: intermediate
597c0ce001eSMatthew G. Knepley 
598dce8aebaSBarry Smith .seealso: `PetscLimiter`, `PetscLimiterType`, `PetscLimiterCreate()`, `PetscLimiterSetType()`
599c0ce001eSMatthew G. Knepley M*/
600c0ce001eSMatthew G. Knepley 
PetscLimiterCreate_Minmod(PetscLimiter lim)601d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscLimiterCreate_Minmod(PetscLimiter lim)
602d71ae5a4SJacob Faibussowitsch {
603c0ce001eSMatthew G. Knepley   PetscLimiter_Minmod *l;
604c0ce001eSMatthew G. Knepley 
605c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
606c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1);
6074dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&l));
608c0ce001eSMatthew G. Knepley   lim->data = l;
609c0ce001eSMatthew G. Knepley 
6109566063dSJacob Faibussowitsch   PetscCall(PetscLimiterInitialize_Minmod(lim));
6113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
612c0ce001eSMatthew G. Knepley }
613c0ce001eSMatthew G. Knepley 
PetscLimiterDestroy_VanLeer(PetscLimiter lim)614d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscLimiterDestroy_VanLeer(PetscLimiter lim)
615d71ae5a4SJacob Faibussowitsch {
616c0ce001eSMatthew G. Knepley   PetscLimiter_VanLeer *l = (PetscLimiter_VanLeer *)lim->data;
617c0ce001eSMatthew G. Knepley 
618c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
6199566063dSJacob Faibussowitsch   PetscCall(PetscFree(l));
6203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
621c0ce001eSMatthew G. Knepley }
622c0ce001eSMatthew G. Knepley 
PetscLimiterView_VanLeer_Ascii(PetscLimiter lim,PetscViewer viewer)623d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscLimiterView_VanLeer_Ascii(PetscLimiter lim, PetscViewer viewer)
624d71ae5a4SJacob Faibussowitsch {
625c0ce001eSMatthew G. Knepley   PetscViewerFormat format;
626c0ce001eSMatthew G. Knepley 
627c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
6289566063dSJacob Faibussowitsch   PetscCall(PetscViewerGetFormat(viewer, &format));
6299566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Van Leer Slope Limiter:\n"));
6303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
631c0ce001eSMatthew G. Knepley }
632c0ce001eSMatthew G. Knepley 
PetscLimiterView_VanLeer(PetscLimiter lim,PetscViewer viewer)633d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscLimiterView_VanLeer(PetscLimiter lim, PetscViewer viewer)
634d71ae5a4SJacob Faibussowitsch {
635*9f196a02SMartin Diehl   PetscBool isascii;
636c0ce001eSMatthew G. Knepley 
637c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
638c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1);
639c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
640*9f196a02SMartin Diehl   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
641*9f196a02SMartin Diehl   if (isascii) PetscCall(PetscLimiterView_VanLeer_Ascii(lim, viewer));
6423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
643c0ce001eSMatthew G. Knepley }
644c0ce001eSMatthew G. Knepley 
PetscLimiterLimit_VanLeer(PetscLimiter lim,PetscReal f,PetscReal * phi)645d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscLimiterLimit_VanLeer(PetscLimiter lim, PetscReal f, PetscReal *phi)
646d71ae5a4SJacob Faibussowitsch {
647c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
648c0ce001eSMatthew G. Knepley   *phi = PetscMax(0, 4 * f * (1 - f));
6493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
650c0ce001eSMatthew G. Knepley }
651c0ce001eSMatthew G. Knepley 
PetscLimiterInitialize_VanLeer(PetscLimiter lim)652d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscLimiterInitialize_VanLeer(PetscLimiter lim)
653d71ae5a4SJacob Faibussowitsch {
654c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
655c0ce001eSMatthew G. Knepley   lim->ops->view    = PetscLimiterView_VanLeer;
656c0ce001eSMatthew G. Knepley   lim->ops->destroy = PetscLimiterDestroy_VanLeer;
657c0ce001eSMatthew G. Knepley   lim->ops->limit   = PetscLimiterLimit_VanLeer;
6583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
659c0ce001eSMatthew G. Knepley }
660c0ce001eSMatthew G. Knepley 
661c0ce001eSMatthew G. Knepley /*MC
662dce8aebaSBarry Smith   PETSCLIMITERVANLEER = "vanleer" - A `PetscLimiter` implementation
663c0ce001eSMatthew G. Knepley 
664c0ce001eSMatthew G. Knepley   Level: intermediate
665c0ce001eSMatthew G. Knepley 
666dce8aebaSBarry Smith .seealso: `PetscLimiter`, `PetscLimiterType`, `PetscLimiterCreate()`, `PetscLimiterSetType()`
667c0ce001eSMatthew G. Knepley M*/
668c0ce001eSMatthew G. Knepley 
PetscLimiterCreate_VanLeer(PetscLimiter lim)669d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscLimiterCreate_VanLeer(PetscLimiter lim)
670d71ae5a4SJacob Faibussowitsch {
671c0ce001eSMatthew G. Knepley   PetscLimiter_VanLeer *l;
672c0ce001eSMatthew G. Knepley 
673c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
674c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1);
6754dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&l));
676c0ce001eSMatthew G. Knepley   lim->data = l;
677c0ce001eSMatthew G. Knepley 
6789566063dSJacob Faibussowitsch   PetscCall(PetscLimiterInitialize_VanLeer(lim));
6793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
680c0ce001eSMatthew G. Knepley }
681c0ce001eSMatthew G. Knepley 
PetscLimiterDestroy_VanAlbada(PetscLimiter lim)682d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscLimiterDestroy_VanAlbada(PetscLimiter lim)
683d71ae5a4SJacob Faibussowitsch {
684c0ce001eSMatthew G. Knepley   PetscLimiter_VanAlbada *l = (PetscLimiter_VanAlbada *)lim->data;
685c0ce001eSMatthew G. Knepley 
686c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
6879566063dSJacob Faibussowitsch   PetscCall(PetscFree(l));
6883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
689c0ce001eSMatthew G. Knepley }
690c0ce001eSMatthew G. Knepley 
PetscLimiterView_VanAlbada_Ascii(PetscLimiter lim,PetscViewer viewer)691d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscLimiterView_VanAlbada_Ascii(PetscLimiter lim, PetscViewer viewer)
692d71ae5a4SJacob Faibussowitsch {
693c0ce001eSMatthew G. Knepley   PetscViewerFormat format;
694c0ce001eSMatthew G. Knepley 
695c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
6969566063dSJacob Faibussowitsch   PetscCall(PetscViewerGetFormat(viewer, &format));
6979566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Van Albada Slope Limiter:\n"));
6983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
699c0ce001eSMatthew G. Knepley }
700c0ce001eSMatthew G. Knepley 
PetscLimiterView_VanAlbada(PetscLimiter lim,PetscViewer viewer)701d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscLimiterView_VanAlbada(PetscLimiter lim, PetscViewer viewer)
702d71ae5a4SJacob Faibussowitsch {
703*9f196a02SMartin Diehl   PetscBool isascii;
704c0ce001eSMatthew G. Knepley 
705c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
706c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1);
707c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
708*9f196a02SMartin Diehl   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
709*9f196a02SMartin Diehl   if (isascii) PetscCall(PetscLimiterView_VanAlbada_Ascii(lim, viewer));
7103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
711c0ce001eSMatthew G. Knepley }
712c0ce001eSMatthew G. Knepley 
PetscLimiterLimit_VanAlbada(PetscLimiter lim,PetscReal f,PetscReal * phi)713d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscLimiterLimit_VanAlbada(PetscLimiter lim, PetscReal f, PetscReal *phi)
714d71ae5a4SJacob Faibussowitsch {
715c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
716c0ce001eSMatthew G. Knepley   *phi = PetscMax(0, 2 * f * (1 - f) / (PetscSqr(f) + PetscSqr(1 - f)));
7173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
718c0ce001eSMatthew G. Knepley }
719c0ce001eSMatthew G. Knepley 
PetscLimiterInitialize_VanAlbada(PetscLimiter lim)720d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscLimiterInitialize_VanAlbada(PetscLimiter lim)
721d71ae5a4SJacob Faibussowitsch {
722c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
723c0ce001eSMatthew G. Knepley   lim->ops->view    = PetscLimiterView_VanAlbada;
724c0ce001eSMatthew G. Knepley   lim->ops->destroy = PetscLimiterDestroy_VanAlbada;
725c0ce001eSMatthew G. Knepley   lim->ops->limit   = PetscLimiterLimit_VanAlbada;
7263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
727c0ce001eSMatthew G. Knepley }
728c0ce001eSMatthew G. Knepley 
729c0ce001eSMatthew G. Knepley /*MC
730dce8aebaSBarry Smith   PETSCLIMITERVANALBADA = "vanalbada" - A PetscLimiter implementation
731c0ce001eSMatthew G. Knepley 
732c0ce001eSMatthew G. Knepley   Level: intermediate
733c0ce001eSMatthew G. Knepley 
734dce8aebaSBarry Smith .seealso: `PetscLimiter`, `PetscLimiterType`, `PetscLimiterCreate()`, `PetscLimiterSetType()`
735c0ce001eSMatthew G. Knepley M*/
736c0ce001eSMatthew G. Knepley 
PetscLimiterCreate_VanAlbada(PetscLimiter lim)737d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscLimiterCreate_VanAlbada(PetscLimiter lim)
738d71ae5a4SJacob Faibussowitsch {
739c0ce001eSMatthew G. Knepley   PetscLimiter_VanAlbada *l;
740c0ce001eSMatthew G. Knepley 
741c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
742c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1);
7434dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&l));
744c0ce001eSMatthew G. Knepley   lim->data = l;
745c0ce001eSMatthew G. Knepley 
7469566063dSJacob Faibussowitsch   PetscCall(PetscLimiterInitialize_VanAlbada(lim));
7473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
748c0ce001eSMatthew G. Knepley }
749c0ce001eSMatthew G. Knepley 
PetscLimiterDestroy_Superbee(PetscLimiter lim)750d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscLimiterDestroy_Superbee(PetscLimiter lim)
751d71ae5a4SJacob Faibussowitsch {
752c0ce001eSMatthew G. Knepley   PetscLimiter_Superbee *l = (PetscLimiter_Superbee *)lim->data;
753c0ce001eSMatthew G. Knepley 
754c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
7559566063dSJacob Faibussowitsch   PetscCall(PetscFree(l));
7563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
757c0ce001eSMatthew G. Knepley }
758c0ce001eSMatthew G. Knepley 
PetscLimiterView_Superbee_Ascii(PetscLimiter lim,PetscViewer viewer)759d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscLimiterView_Superbee_Ascii(PetscLimiter lim, PetscViewer viewer)
760d71ae5a4SJacob Faibussowitsch {
761c0ce001eSMatthew G. Knepley   PetscViewerFormat format;
762c0ce001eSMatthew G. Knepley 
763c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
7649566063dSJacob Faibussowitsch   PetscCall(PetscViewerGetFormat(viewer, &format));
7659566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Superbee Slope Limiter:\n"));
7663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
767c0ce001eSMatthew G. Knepley }
768c0ce001eSMatthew G. Knepley 
PetscLimiterView_Superbee(PetscLimiter lim,PetscViewer viewer)769d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscLimiterView_Superbee(PetscLimiter lim, PetscViewer viewer)
770d71ae5a4SJacob Faibussowitsch {
771*9f196a02SMartin Diehl   PetscBool isascii;
772c0ce001eSMatthew G. Knepley 
773c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
774c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1);
775c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
776*9f196a02SMartin Diehl   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
777*9f196a02SMartin Diehl   if (isascii) PetscCall(PetscLimiterView_Superbee_Ascii(lim, viewer));
7783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
779c0ce001eSMatthew G. Knepley }
780c0ce001eSMatthew G. Knepley 
PetscLimiterLimit_Superbee(PetscLimiter lim,PetscReal f,PetscReal * phi)781d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscLimiterLimit_Superbee(PetscLimiter lim, PetscReal f, PetscReal *phi)
782d71ae5a4SJacob Faibussowitsch {
783c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
784c0ce001eSMatthew G. Knepley   *phi = 4 * PetscMax(0, PetscMin(f, 1 - f));
7853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
786c0ce001eSMatthew G. Knepley }
787c0ce001eSMatthew G. Knepley 
PetscLimiterInitialize_Superbee(PetscLimiter lim)788d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscLimiterInitialize_Superbee(PetscLimiter lim)
789d71ae5a4SJacob Faibussowitsch {
790c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
791c0ce001eSMatthew G. Knepley   lim->ops->view    = PetscLimiterView_Superbee;
792c0ce001eSMatthew G. Knepley   lim->ops->destroy = PetscLimiterDestroy_Superbee;
793c0ce001eSMatthew G. Knepley   lim->ops->limit   = PetscLimiterLimit_Superbee;
7943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
795c0ce001eSMatthew G. Knepley }
796c0ce001eSMatthew G. Knepley 
797c0ce001eSMatthew G. Knepley /*MC
798dce8aebaSBarry Smith   PETSCLIMITERSUPERBEE = "superbee" - A `PetscLimiter` implementation
799c0ce001eSMatthew G. Knepley 
800c0ce001eSMatthew G. Knepley   Level: intermediate
801c0ce001eSMatthew G. Knepley 
802dce8aebaSBarry Smith .seealso: `PetscLimiter`, `PetscLimiterType`, `PetscLimiterCreate()`, `PetscLimiterSetType()`
803c0ce001eSMatthew G. Knepley M*/
804c0ce001eSMatthew G. Knepley 
PetscLimiterCreate_Superbee(PetscLimiter lim)805d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscLimiterCreate_Superbee(PetscLimiter lim)
806d71ae5a4SJacob Faibussowitsch {
807c0ce001eSMatthew G. Knepley   PetscLimiter_Superbee *l;
808c0ce001eSMatthew G. Knepley 
809c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
810c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1);
8114dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&l));
812c0ce001eSMatthew G. Knepley   lim->data = l;
813c0ce001eSMatthew G. Knepley 
8149566063dSJacob Faibussowitsch   PetscCall(PetscLimiterInitialize_Superbee(lim));
8153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
816c0ce001eSMatthew G. Knepley }
817c0ce001eSMatthew G. Knepley 
PetscLimiterDestroy_MC(PetscLimiter lim)818d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscLimiterDestroy_MC(PetscLimiter lim)
819d71ae5a4SJacob Faibussowitsch {
820c0ce001eSMatthew G. Knepley   PetscLimiter_MC *l = (PetscLimiter_MC *)lim->data;
821c0ce001eSMatthew G. Knepley 
822c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
8239566063dSJacob Faibussowitsch   PetscCall(PetscFree(l));
8243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
825c0ce001eSMatthew G. Knepley }
826c0ce001eSMatthew G. Knepley 
PetscLimiterView_MC_Ascii(PetscLimiter lim,PetscViewer viewer)827d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscLimiterView_MC_Ascii(PetscLimiter lim, PetscViewer viewer)
828d71ae5a4SJacob Faibussowitsch {
829c0ce001eSMatthew G. Knepley   PetscViewerFormat format;
830c0ce001eSMatthew G. Knepley 
831c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
8329566063dSJacob Faibussowitsch   PetscCall(PetscViewerGetFormat(viewer, &format));
8339566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "MC Slope Limiter:\n"));
8343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
835c0ce001eSMatthew G. Knepley }
836c0ce001eSMatthew G. Knepley 
PetscLimiterView_MC(PetscLimiter lim,PetscViewer viewer)837d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscLimiterView_MC(PetscLimiter lim, PetscViewer viewer)
838d71ae5a4SJacob Faibussowitsch {
839*9f196a02SMartin Diehl   PetscBool isascii;
840c0ce001eSMatthew G. Knepley 
841c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
842c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1);
843c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
844*9f196a02SMartin Diehl   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
845*9f196a02SMartin Diehl   if (isascii) PetscCall(PetscLimiterView_MC_Ascii(lim, viewer));
8463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
847c0ce001eSMatthew G. Knepley }
848c0ce001eSMatthew G. Knepley 
849c0ce001eSMatthew G. Knepley /* aka Barth-Jespersen */
PetscLimiterLimit_MC(PetscLimiter lim,PetscReal f,PetscReal * phi)850d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscLimiterLimit_MC(PetscLimiter lim, PetscReal f, PetscReal *phi)
851d71ae5a4SJacob Faibussowitsch {
852c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
853c0ce001eSMatthew G. Knepley   *phi = PetscMin(1, 4 * PetscMax(0, PetscMin(f, 1 - f)));
8543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
855c0ce001eSMatthew G. Knepley }
856c0ce001eSMatthew G. Knepley 
PetscLimiterInitialize_MC(PetscLimiter lim)857d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscLimiterInitialize_MC(PetscLimiter lim)
858d71ae5a4SJacob Faibussowitsch {
859c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
860c0ce001eSMatthew G. Knepley   lim->ops->view    = PetscLimiterView_MC;
861c0ce001eSMatthew G. Knepley   lim->ops->destroy = PetscLimiterDestroy_MC;
862c0ce001eSMatthew G. Knepley   lim->ops->limit   = PetscLimiterLimit_MC;
8633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
864c0ce001eSMatthew G. Knepley }
865c0ce001eSMatthew G. Knepley 
866c0ce001eSMatthew G. Knepley /*MC
867dce8aebaSBarry Smith   PETSCLIMITERMC = "mc" - A `PetscLimiter` implementation
868c0ce001eSMatthew G. Knepley 
869c0ce001eSMatthew G. Knepley   Level: intermediate
870c0ce001eSMatthew G. Knepley 
871dce8aebaSBarry Smith .seealso: `PetscLimiter`, `PetscLimiterType`, `PetscLimiterCreate()`, `PetscLimiterSetType()`
872c0ce001eSMatthew G. Knepley M*/
873c0ce001eSMatthew G. Knepley 
PetscLimiterCreate_MC(PetscLimiter lim)874d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscLimiterCreate_MC(PetscLimiter lim)
875d71ae5a4SJacob Faibussowitsch {
876c0ce001eSMatthew G. Knepley   PetscLimiter_MC *l;
877c0ce001eSMatthew G. Knepley 
878c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
879c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1);
8804dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&l));
881c0ce001eSMatthew G. Knepley   lim->data = l;
882c0ce001eSMatthew G. Knepley 
8839566063dSJacob Faibussowitsch   PetscCall(PetscLimiterInitialize_MC(lim));
8843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
885c0ce001eSMatthew G. Knepley }
886c0ce001eSMatthew G. Knepley 
887c0ce001eSMatthew G. Knepley PetscClassId PETSCFV_CLASSID = 0;
888c0ce001eSMatthew G. Knepley 
889c0ce001eSMatthew G. Knepley PetscFunctionList PetscFVList              = NULL;
890c0ce001eSMatthew G. Knepley PetscBool         PetscFVRegisterAllCalled = PETSC_FALSE;
891c0ce001eSMatthew G. Knepley 
892c0ce001eSMatthew G. Knepley /*@C
893dce8aebaSBarry Smith   PetscFVRegister - Adds a new `PetscFV` implementation
894c0ce001eSMatthew G. Knepley 
895cc4c1da9SBarry Smith   Not Collective, No Fortran Support
896c0ce001eSMatthew G. Knepley 
897c0ce001eSMatthew G. Knepley   Input Parameters:
8982fe279fdSBarry Smith + sname    - The name of a new user-defined creation routine
8992fe279fdSBarry Smith - function - The creation routine itself
900c0ce001eSMatthew G. Knepley 
90160225df5SJacob Faibussowitsch   Example Usage:
902c0ce001eSMatthew G. Knepley .vb
903c0ce001eSMatthew G. Knepley     PetscFVRegister("my_fv", MyPetscFVCreate);
904c0ce001eSMatthew G. Knepley .ve
905c0ce001eSMatthew G. Knepley 
906c0ce001eSMatthew G. Knepley   Then, your PetscFV type can be chosen with the procedural interface via
907c0ce001eSMatthew G. Knepley .vb
908c0ce001eSMatthew G. Knepley     PetscFVCreate(MPI_Comm, PetscFV *);
909c0ce001eSMatthew G. Knepley     PetscFVSetType(PetscFV, "my_fv");
910c0ce001eSMatthew G. Knepley .ve
911c0ce001eSMatthew G. Knepley   or at runtime via the option
912c0ce001eSMatthew G. Knepley .vb
913c0ce001eSMatthew G. Knepley     -petscfv_type my_fv
914c0ce001eSMatthew G. Knepley .ve
915c0ce001eSMatthew G. Knepley 
916c0ce001eSMatthew G. Knepley   Level: advanced
917c0ce001eSMatthew G. Knepley 
918dce8aebaSBarry Smith   Note:
919dce8aebaSBarry Smith   `PetscFVRegister()` may be called multiple times to add several user-defined PetscFVs
920c0ce001eSMatthew G. Knepley 
921dce8aebaSBarry Smith .seealso: `PetscFV`, `PetscFVType`, `PetscFVRegisterAll()`, `PetscFVRegisterDestroy()`
922c0ce001eSMatthew G. Knepley @*/
PetscFVRegister(const char sname[],PetscErrorCode (* function)(PetscFV))923d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFVRegister(const char sname[], PetscErrorCode (*function)(PetscFV))
924d71ae5a4SJacob Faibussowitsch {
925c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
9269566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&PetscFVList, sname, function));
9273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
928c0ce001eSMatthew G. Knepley }
929c0ce001eSMatthew G. Knepley 
930cc4c1da9SBarry Smith /*@
931dce8aebaSBarry Smith   PetscFVSetType - Builds a particular `PetscFV`
932c0ce001eSMatthew G. Knepley 
93320f4b53cSBarry Smith   Collective
934c0ce001eSMatthew G. Knepley 
935c0ce001eSMatthew G. Knepley   Input Parameters:
936dce8aebaSBarry Smith + fvm  - The `PetscFV` object
937dce8aebaSBarry Smith - name - The type of FVM space
938c0ce001eSMatthew G. Knepley 
939c0ce001eSMatthew G. Knepley   Options Database Key:
940dce8aebaSBarry Smith . -petscfv_type <type> - Sets the `PetscFVType`; use -help for a list of available types
941c0ce001eSMatthew G. Knepley 
942c0ce001eSMatthew G. Knepley   Level: intermediate
943c0ce001eSMatthew G. Knepley 
944dce8aebaSBarry Smith .seealso: `PetscFV`, `PetscFVType`, `PetscFVGetType()`, `PetscFVCreate()`
945c0ce001eSMatthew G. Knepley @*/
PetscFVSetType(PetscFV fvm,PetscFVType name)946d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFVSetType(PetscFV fvm, PetscFVType name)
947d71ae5a4SJacob Faibussowitsch {
948c0ce001eSMatthew G. Knepley   PetscErrorCode (*r)(PetscFV);
949c0ce001eSMatthew G. Knepley   PetscBool match;
950c0ce001eSMatthew G. Knepley 
951c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
952c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1);
9539566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)fvm, name, &match));
9543ba16761SJacob Faibussowitsch   if (match) PetscFunctionReturn(PETSC_SUCCESS);
955c0ce001eSMatthew G. Knepley 
9569566063dSJacob Faibussowitsch   PetscCall(PetscFVRegisterAll());
9579566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListFind(PetscFVList, name, &r));
95828b400f6SJacob Faibussowitsch   PetscCheck(r, PetscObjectComm((PetscObject)fvm), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscFV type: %s", name);
959c0ce001eSMatthew G. Knepley 
960dbbe0bcdSBarry Smith   PetscTryTypeMethod(fvm, destroy);
961c0ce001eSMatthew G. Knepley   fvm->ops->destroy = NULL;
962dbbe0bcdSBarry Smith 
9639566063dSJacob Faibussowitsch   PetscCall((*r)(fvm));
9649566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)fvm, name));
9653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
966c0ce001eSMatthew G. Knepley }
967c0ce001eSMatthew G. Knepley 
968cc4c1da9SBarry Smith /*@
969dce8aebaSBarry Smith   PetscFVGetType - Gets the `PetscFVType` (as a string) from a `PetscFV`.
970c0ce001eSMatthew G. Knepley 
971c0ce001eSMatthew G. Knepley   Not Collective
972c0ce001eSMatthew G. Knepley 
973c0ce001eSMatthew G. Knepley   Input Parameter:
974dce8aebaSBarry Smith . fvm - The `PetscFV`
975c0ce001eSMatthew G. Knepley 
976c0ce001eSMatthew G. Knepley   Output Parameter:
977dce8aebaSBarry Smith . name - The `PetscFVType` name
978c0ce001eSMatthew G. Knepley 
979c0ce001eSMatthew G. Knepley   Level: intermediate
980c0ce001eSMatthew G. Knepley 
981dce8aebaSBarry Smith .seealso: `PetscFV`, `PetscFVType`, `PetscFVSetType()`, `PetscFVCreate()`
982c0ce001eSMatthew G. Knepley @*/
PetscFVGetType(PetscFV fvm,PetscFVType * name)983d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFVGetType(PetscFV fvm, PetscFVType *name)
984d71ae5a4SJacob Faibussowitsch {
985c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
986c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1);
9874f572ea9SToby Isaac   PetscAssertPointer(name, 2);
9889566063dSJacob Faibussowitsch   PetscCall(PetscFVRegisterAll());
989c0ce001eSMatthew G. Knepley   *name = ((PetscObject)fvm)->type_name;
9903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
991c0ce001eSMatthew G. Knepley }
992c0ce001eSMatthew G. Knepley 
993ffeef943SBarry Smith /*@
994dce8aebaSBarry Smith   PetscFVViewFromOptions - View a `PetscFV` based on values in the options database
995fe2efc57SMark 
99620f4b53cSBarry Smith   Collective
997fe2efc57SMark 
998fe2efc57SMark   Input Parameters:
999dce8aebaSBarry Smith + A    - the `PetscFV` object
1000dce8aebaSBarry Smith . obj  - Optional object that provides the options prefix
1001dce8aebaSBarry Smith - name - command line option name
1002fe2efc57SMark 
1003fe2efc57SMark   Level: intermediate
1004dce8aebaSBarry Smith 
1005dce8aebaSBarry Smith .seealso: `PetscFV`, `PetscFVView()`, `PetscObjectViewFromOptions()`, `PetscFVCreate()`
1006fe2efc57SMark @*/
PetscFVViewFromOptions(PetscFV A,PetscObject obj,const char name[])1007d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFVViewFromOptions(PetscFV A, PetscObject obj, const char name[])
1008d71ae5a4SJacob Faibussowitsch {
1009fe2efc57SMark   PetscFunctionBegin;
1010fe2efc57SMark   PetscValidHeaderSpecific(A, PETSCFV_CLASSID, 1);
10119566063dSJacob Faibussowitsch   PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name));
10123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1013fe2efc57SMark }
1014fe2efc57SMark 
1015ffeef943SBarry Smith /*@
1016dce8aebaSBarry Smith   PetscFVView - Views a `PetscFV`
1017c0ce001eSMatthew G. Knepley 
101820f4b53cSBarry Smith   Collective
1019c0ce001eSMatthew G. Knepley 
1020d8d19677SJose E. Roman   Input Parameters:
1021dce8aebaSBarry Smith + fvm - the `PetscFV` object to view
1022c0ce001eSMatthew G. Knepley - v   - the viewer
1023c0ce001eSMatthew G. Knepley 
102488f5f89eSMatthew G. Knepley   Level: beginner
1025c0ce001eSMatthew G. Knepley 
1026dce8aebaSBarry Smith .seealso: `PetscFV`, `PetscViewer`, `PetscFVDestroy()`
1027c0ce001eSMatthew G. Knepley @*/
PetscFVView(PetscFV fvm,PetscViewer v)1028d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFVView(PetscFV fvm, PetscViewer v)
1029d71ae5a4SJacob Faibussowitsch {
1030c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1031c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1);
10329566063dSJacob Faibussowitsch   if (!v) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)fvm), &v));
1033dbbe0bcdSBarry Smith   PetscTryTypeMethod(fvm, view, v);
10343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1035c0ce001eSMatthew G. Knepley }
1036c0ce001eSMatthew G. Knepley 
1037c0ce001eSMatthew G. Knepley /*@
1038dce8aebaSBarry Smith   PetscFVSetFromOptions - sets parameters in a `PetscFV` from the options database
1039c0ce001eSMatthew G. Knepley 
104020f4b53cSBarry Smith   Collective
1041c0ce001eSMatthew G. Knepley 
1042c0ce001eSMatthew G. Knepley   Input Parameter:
1043dce8aebaSBarry Smith . fvm - the `PetscFV` object to set options for
1044c0ce001eSMatthew G. Knepley 
1045c0ce001eSMatthew G. Knepley   Options Database Key:
1046c0ce001eSMatthew G. Knepley . -petscfv_compute_gradients <bool> - Determines whether cell gradients are calculated
1047c0ce001eSMatthew G. Knepley 
104888f5f89eSMatthew G. Knepley   Level: intermediate
1049c0ce001eSMatthew G. Knepley 
1050dce8aebaSBarry Smith .seealso: `PetscFV`, `PetscFVView()`
1051c0ce001eSMatthew G. Knepley @*/
PetscFVSetFromOptions(PetscFV fvm)1052d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFVSetFromOptions(PetscFV fvm)
1053d71ae5a4SJacob Faibussowitsch {
1054c0ce001eSMatthew G. Knepley   const char *defaultType;
1055c0ce001eSMatthew G. Knepley   char        name[256];
1056c0ce001eSMatthew G. Knepley   PetscBool   flg;
1057c0ce001eSMatthew G. Knepley 
1058c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1059c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1);
1060c0ce001eSMatthew G. Knepley   if (!((PetscObject)fvm)->type_name) defaultType = PETSCFVUPWIND;
1061c0ce001eSMatthew G. Knepley   else defaultType = ((PetscObject)fvm)->type_name;
10629566063dSJacob Faibussowitsch   PetscCall(PetscFVRegisterAll());
1063c0ce001eSMatthew G. Knepley 
1064d0609cedSBarry Smith   PetscObjectOptionsBegin((PetscObject)fvm);
10659566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFList("-petscfv_type", "Finite volume discretization", "PetscFVSetType", PetscFVList, defaultType, name, 256, &flg));
1066c0ce001eSMatthew G. Knepley   if (flg) {
10679566063dSJacob Faibussowitsch     PetscCall(PetscFVSetType(fvm, name));
1068c0ce001eSMatthew G. Knepley   } else if (!((PetscObject)fvm)->type_name) {
10699566063dSJacob Faibussowitsch     PetscCall(PetscFVSetType(fvm, defaultType));
1070c0ce001eSMatthew G. Knepley   }
10719566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-petscfv_compute_gradients", "Compute cell gradients", "PetscFVSetComputeGradients", fvm->computeGradients, &fvm->computeGradients, NULL));
1072dbbe0bcdSBarry Smith   PetscTryTypeMethod(fvm, setfromoptions);
1073c0ce001eSMatthew G. Knepley   /* process any options handlers added with PetscObjectAddOptionsHandler() */
1074dbbe0bcdSBarry Smith   PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)fvm, PetscOptionsObject));
10759566063dSJacob Faibussowitsch   PetscCall(PetscLimiterSetFromOptions(fvm->limiter));
1076d0609cedSBarry Smith   PetscOptionsEnd();
10779566063dSJacob Faibussowitsch   PetscCall(PetscFVViewFromOptions(fvm, NULL, "-petscfv_view"));
10783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1079c0ce001eSMatthew G. Knepley }
1080c0ce001eSMatthew G. Knepley 
1081c0ce001eSMatthew G. Knepley /*@
1082dce8aebaSBarry Smith   PetscFVSetUp - Setup the data structures for the `PetscFV` based on the `PetscFVType` provided by `PetscFVSetType()`
1083c0ce001eSMatthew G. Knepley 
108420f4b53cSBarry Smith   Collective
1085c0ce001eSMatthew G. Knepley 
1086c0ce001eSMatthew G. Knepley   Input Parameter:
1087dce8aebaSBarry Smith . fvm - the `PetscFV` object to setup
1088c0ce001eSMatthew G. Knepley 
108988f5f89eSMatthew G. Knepley   Level: intermediate
1090c0ce001eSMatthew G. Knepley 
1091dce8aebaSBarry Smith .seealso: `PetscFV`, `PetscFVView()`, `PetscFVDestroy()`
1092c0ce001eSMatthew G. Knepley @*/
PetscFVSetUp(PetscFV fvm)1093d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFVSetUp(PetscFV fvm)
1094d71ae5a4SJacob Faibussowitsch {
1095c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1096c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1);
10979566063dSJacob Faibussowitsch   PetscCall(PetscLimiterSetUp(fvm->limiter));
1098dbbe0bcdSBarry Smith   PetscTryTypeMethod(fvm, setup);
10993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1100c0ce001eSMatthew G. Knepley }
1101c0ce001eSMatthew G. Knepley 
1102c0ce001eSMatthew G. Knepley /*@
1103dce8aebaSBarry Smith   PetscFVDestroy - Destroys a `PetscFV` object
1104c0ce001eSMatthew G. Knepley 
110520f4b53cSBarry Smith   Collective
1106c0ce001eSMatthew G. Knepley 
1107c0ce001eSMatthew G. Knepley   Input Parameter:
1108dce8aebaSBarry Smith . fvm - the `PetscFV` object to destroy
1109c0ce001eSMatthew G. Knepley 
111088f5f89eSMatthew G. Knepley   Level: beginner
1111c0ce001eSMatthew G. Knepley 
1112dce8aebaSBarry Smith .seealso: `PetscFV`, `PetscFVCreate()`, `PetscFVView()`
1113c0ce001eSMatthew G. Knepley @*/
PetscFVDestroy(PetscFV * fvm)1114d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFVDestroy(PetscFV *fvm)
1115d71ae5a4SJacob Faibussowitsch {
1116c0ce001eSMatthew G. Knepley   PetscInt i;
1117c0ce001eSMatthew G. Knepley 
1118c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
11193ba16761SJacob Faibussowitsch   if (!*fvm) PetscFunctionReturn(PETSC_SUCCESS);
1120f4f49eeaSPierre Jolivet   PetscValidHeaderSpecific(*fvm, PETSCFV_CLASSID, 1);
1121c0ce001eSMatthew G. Knepley 
1122f4f49eeaSPierre Jolivet   if (--((PetscObject)*fvm)->refct > 0) {
11239371c9d4SSatish Balay     *fvm = NULL;
11243ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
11259371c9d4SSatish Balay   }
1126f4f49eeaSPierre Jolivet   ((PetscObject)*fvm)->refct = 0;
1127c0ce001eSMatthew G. Knepley 
112848a46eb9SPierre Jolivet   for (i = 0; i < (*fvm)->numComponents; i++) PetscCall(PetscFree((*fvm)->componentNames[i]));
11299566063dSJacob Faibussowitsch   PetscCall(PetscFree((*fvm)->componentNames));
11309566063dSJacob Faibussowitsch   PetscCall(PetscLimiterDestroy(&(*fvm)->limiter));
11319566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceDestroy(&(*fvm)->dualSpace));
11329566063dSJacob Faibussowitsch   PetscCall(PetscFree((*fvm)->fluxWork));
11339566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureDestroy(&(*fvm)->quadrature));
11349566063dSJacob Faibussowitsch   PetscCall(PetscTabulationDestroy(&(*fvm)->T));
1135c0ce001eSMatthew G. Knepley 
1136f4f49eeaSPierre Jolivet   PetscTryTypeMethod(*fvm, destroy);
11379566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(fvm));
11383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1139c0ce001eSMatthew G. Knepley }
1140c0ce001eSMatthew G. Knepley 
1141c0ce001eSMatthew G. Knepley /*@
1142dce8aebaSBarry Smith   PetscFVCreate - Creates an empty `PetscFV` object. The type can then be set with `PetscFVSetType()`.
1143c0ce001eSMatthew G. Knepley 
1144c0ce001eSMatthew G. Knepley   Collective
1145c0ce001eSMatthew G. Knepley 
1146c0ce001eSMatthew G. Knepley   Input Parameter:
1147dce8aebaSBarry Smith . comm - The communicator for the `PetscFV` object
1148c0ce001eSMatthew G. Knepley 
1149c0ce001eSMatthew G. Knepley   Output Parameter:
1150dce8aebaSBarry Smith . fvm - The `PetscFV` object
1151c0ce001eSMatthew G. Knepley 
1152c0ce001eSMatthew G. Knepley   Level: beginner
1153c0ce001eSMatthew G. Knepley 
11541d27aa22SBarry Smith .seealso: `PetscFVSetUp()`, `PetscFVSetType()`, `PETSCFVUPWIND`, `PetscFVDestroy()`
1155c0ce001eSMatthew G. Knepley @*/
PetscFVCreate(MPI_Comm comm,PetscFV * fvm)1156d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFVCreate(MPI_Comm comm, PetscFV *fvm)
1157d71ae5a4SJacob Faibussowitsch {
1158c0ce001eSMatthew G. Knepley   PetscFV f;
1159c0ce001eSMatthew G. Knepley 
1160c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
11614f572ea9SToby Isaac   PetscAssertPointer(fvm, 2);
11629566063dSJacob Faibussowitsch   PetscCall(PetscFVInitializePackage());
1163c0ce001eSMatthew G. Knepley 
11649566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(f, PETSCFV_CLASSID, "PetscFV", "Finite Volume", "PetscFV", comm, PetscFVDestroy, PetscFVView));
11659566063dSJacob Faibussowitsch   PetscCall(PetscMemzero(f->ops, sizeof(struct _PetscFVOps)));
11669566063dSJacob Faibussowitsch   PetscCall(PetscLimiterCreate(comm, &f->limiter));
1167c0ce001eSMatthew G. Knepley   f->numComponents    = 1;
1168c0ce001eSMatthew G. Knepley   f->dim              = 0;
1169c0ce001eSMatthew G. Knepley   f->computeGradients = PETSC_FALSE;
1170c0ce001eSMatthew G. Knepley   f->fluxWork         = NULL;
11719566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(f->numComponents, &f->componentNames));
1172c0ce001eSMatthew G. Knepley 
1173c0ce001eSMatthew G. Knepley   *fvm = f;
11743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1175c0ce001eSMatthew G. Knepley }
1176c0ce001eSMatthew G. Knepley 
1177c0ce001eSMatthew G. Knepley /*@
1178dce8aebaSBarry Smith   PetscFVSetLimiter - Set the `PetscLimiter` to the `PetscFV`
1179c0ce001eSMatthew G. Knepley 
118020f4b53cSBarry Smith   Logically Collective
1181c0ce001eSMatthew G. Knepley 
1182c0ce001eSMatthew G. Knepley   Input Parameters:
1183dce8aebaSBarry Smith + fvm - the `PetscFV` object
1184dce8aebaSBarry Smith - lim - The `PetscLimiter`
1185c0ce001eSMatthew G. Knepley 
118688f5f89eSMatthew G. Knepley   Level: intermediate
1187c0ce001eSMatthew G. Knepley 
1188dce8aebaSBarry Smith .seealso: `PetscFV`, `PetscLimiter`, `PetscFVGetLimiter()`
1189c0ce001eSMatthew G. Knepley @*/
PetscFVSetLimiter(PetscFV fvm,PetscLimiter lim)1190d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFVSetLimiter(PetscFV fvm, PetscLimiter lim)
1191d71ae5a4SJacob Faibussowitsch {
1192c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1193c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1);
1194c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 2);
11959566063dSJacob Faibussowitsch   PetscCall(PetscLimiterDestroy(&fvm->limiter));
11969566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)lim));
1197c0ce001eSMatthew G. Knepley   fvm->limiter = lim;
11983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1199c0ce001eSMatthew G. Knepley }
1200c0ce001eSMatthew G. Knepley 
1201c0ce001eSMatthew G. Knepley /*@
1202dce8aebaSBarry Smith   PetscFVGetLimiter - Get the `PetscLimiter` object from the `PetscFV`
1203c0ce001eSMatthew G. Knepley 
120420f4b53cSBarry Smith   Not Collective
1205c0ce001eSMatthew G. Knepley 
1206c0ce001eSMatthew G. Knepley   Input Parameter:
1207dce8aebaSBarry Smith . fvm - the `PetscFV` object
1208c0ce001eSMatthew G. Knepley 
1209c0ce001eSMatthew G. Knepley   Output Parameter:
1210dce8aebaSBarry Smith . lim - The `PetscLimiter`
1211c0ce001eSMatthew G. Knepley 
121288f5f89eSMatthew G. Knepley   Level: intermediate
1213c0ce001eSMatthew G. Knepley 
1214dce8aebaSBarry Smith .seealso: `PetscFV`, `PetscLimiter`, `PetscFVSetLimiter()`
1215c0ce001eSMatthew G. Knepley @*/
PetscFVGetLimiter(PetscFV fvm,PetscLimiter * lim)1216d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFVGetLimiter(PetscFV fvm, PetscLimiter *lim)
1217d71ae5a4SJacob Faibussowitsch {
1218c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1219c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1);
12204f572ea9SToby Isaac   PetscAssertPointer(lim, 2);
1221c0ce001eSMatthew G. Knepley   *lim = fvm->limiter;
12223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1223c0ce001eSMatthew G. Knepley }
1224c0ce001eSMatthew G. Knepley 
1225c0ce001eSMatthew G. Knepley /*@
1226dce8aebaSBarry Smith   PetscFVSetNumComponents - Set the number of field components in a `PetscFV`
1227c0ce001eSMatthew G. Knepley 
122820f4b53cSBarry Smith   Logically Collective
1229c0ce001eSMatthew G. Knepley 
1230c0ce001eSMatthew G. Knepley   Input Parameters:
1231dce8aebaSBarry Smith + fvm  - the `PetscFV` object
1232c0ce001eSMatthew G. Knepley - comp - The number of components
1233c0ce001eSMatthew G. Knepley 
123488f5f89eSMatthew G. Knepley   Level: intermediate
1235c0ce001eSMatthew G. Knepley 
1236dce8aebaSBarry Smith .seealso: `PetscFV`, `PetscFVGetNumComponents()`
1237c0ce001eSMatthew G. Knepley @*/
PetscFVSetNumComponents(PetscFV fvm,PetscInt comp)1238d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFVSetNumComponents(PetscFV fvm, PetscInt comp)
1239d71ae5a4SJacob Faibussowitsch {
1240c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1241c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1);
1242c0ce001eSMatthew G. Knepley   if (fvm->numComponents != comp) {
1243c0ce001eSMatthew G. Knepley     PetscInt i;
1244c0ce001eSMatthew G. Knepley 
124548a46eb9SPierre Jolivet     for (i = 0; i < fvm->numComponents; i++) PetscCall(PetscFree(fvm->componentNames[i]));
12469566063dSJacob Faibussowitsch     PetscCall(PetscFree(fvm->componentNames));
12479566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(comp, &fvm->componentNames));
1248c0ce001eSMatthew G. Knepley   }
1249c0ce001eSMatthew G. Knepley   fvm->numComponents = comp;
12509566063dSJacob Faibussowitsch   PetscCall(PetscFree(fvm->fluxWork));
12519566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(comp, &fvm->fluxWork));
12523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1253c0ce001eSMatthew G. Knepley }
1254c0ce001eSMatthew G. Knepley 
1255c0ce001eSMatthew G. Knepley /*@
1256dce8aebaSBarry Smith   PetscFVGetNumComponents - Get the number of field components in a `PetscFV`
1257c0ce001eSMatthew G. Knepley 
125820f4b53cSBarry Smith   Not Collective
1259c0ce001eSMatthew G. Knepley 
1260c0ce001eSMatthew G. Knepley   Input Parameter:
1261dce8aebaSBarry Smith . fvm - the `PetscFV` object
1262c0ce001eSMatthew G. Knepley 
1263c0ce001eSMatthew G. Knepley   Output Parameter:
1264a4e35b19SJacob Faibussowitsch . comp - The number of components
1265c0ce001eSMatthew G. Knepley 
126688f5f89eSMatthew G. Knepley   Level: intermediate
1267c0ce001eSMatthew G. Knepley 
1268dce8aebaSBarry Smith .seealso: `PetscFV`, `PetscFVSetNumComponents()`, `PetscFVSetComponentName()`
1269c0ce001eSMatthew G. Knepley @*/
PetscFVGetNumComponents(PetscFV fvm,PetscInt * comp)1270d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFVGetNumComponents(PetscFV fvm, PetscInt *comp)
1271d71ae5a4SJacob Faibussowitsch {
1272c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1273c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1);
12744f572ea9SToby Isaac   PetscAssertPointer(comp, 2);
1275c0ce001eSMatthew G. Knepley   *comp = fvm->numComponents;
12763ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1277c0ce001eSMatthew G. Knepley }
1278c0ce001eSMatthew G. Knepley 
12795d83a8b1SBarry Smith /*@
1280dce8aebaSBarry Smith   PetscFVSetComponentName - Set the name of a component (used in output and viewing) in a `PetscFV`
1281c0ce001eSMatthew G. Knepley 
128220f4b53cSBarry Smith   Logically Collective
1283dce8aebaSBarry Smith 
1284c0ce001eSMatthew G. Knepley   Input Parameters:
1285dce8aebaSBarry Smith + fvm  - the `PetscFV` object
1286c0ce001eSMatthew G. Knepley . comp - the component number
1287c0ce001eSMatthew G. Knepley - name - the component name
1288c0ce001eSMatthew G. Knepley 
128988f5f89eSMatthew G. Knepley   Level: intermediate
1290c0ce001eSMatthew G. Knepley 
1291dce8aebaSBarry Smith .seealso: `PetscFV`, `PetscFVGetComponentName()`
1292c0ce001eSMatthew G. Knepley @*/
PetscFVSetComponentName(PetscFV fvm,PetscInt comp,const char * name)1293d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFVSetComponentName(PetscFV fvm, PetscInt comp, const char *name)
1294d71ae5a4SJacob Faibussowitsch {
1295c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
12969566063dSJacob Faibussowitsch   PetscCall(PetscFree(fvm->componentNames[comp]));
12979566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(name, &fvm->componentNames[comp]));
12983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1299c0ce001eSMatthew G. Knepley }
1300c0ce001eSMatthew G. Knepley 
1301cc4c1da9SBarry Smith /*@
1302dce8aebaSBarry Smith   PetscFVGetComponentName - Get the name of a component (used in output and viewing) in a `PetscFV`
1303c0ce001eSMatthew G. Knepley 
130420f4b53cSBarry Smith   Logically Collective
130560225df5SJacob Faibussowitsch 
1306c0ce001eSMatthew G. Knepley   Input Parameters:
1307dce8aebaSBarry Smith + fvm  - the `PetscFV` object
1308c0ce001eSMatthew G. Knepley - comp - the component number
1309c0ce001eSMatthew G. Knepley 
1310c0ce001eSMatthew G. Knepley   Output Parameter:
1311c0ce001eSMatthew G. Knepley . name - the component name
1312c0ce001eSMatthew G. Knepley 
131388f5f89eSMatthew G. Knepley   Level: intermediate
1314c0ce001eSMatthew G. Knepley 
1315dce8aebaSBarry Smith .seealso: `PetscFV`, `PetscFVSetComponentName()`
1316c0ce001eSMatthew G. Knepley @*/
PetscFVGetComponentName(PetscFV fvm,PetscInt comp,const char * name[])1317cc4c1da9SBarry Smith PetscErrorCode PetscFVGetComponentName(PetscFV fvm, PetscInt comp, const char *name[])
1318d71ae5a4SJacob Faibussowitsch {
1319c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1320c0ce001eSMatthew G. Knepley   *name = fvm->componentNames[comp];
13213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1322c0ce001eSMatthew G. Knepley }
1323c0ce001eSMatthew G. Knepley 
1324c0ce001eSMatthew G. Knepley /*@
1325dce8aebaSBarry Smith   PetscFVSetSpatialDimension - Set the spatial dimension of a `PetscFV`
1326c0ce001eSMatthew G. Knepley 
132720f4b53cSBarry Smith   Logically Collective
1328c0ce001eSMatthew G. Knepley 
1329c0ce001eSMatthew G. Knepley   Input Parameters:
1330dce8aebaSBarry Smith + fvm - the `PetscFV` object
1331c0ce001eSMatthew G. Knepley - dim - The spatial dimension
1332c0ce001eSMatthew G. Knepley 
133388f5f89eSMatthew G. Knepley   Level: intermediate
1334c0ce001eSMatthew G. Knepley 
1335a94f484eSPierre Jolivet .seealso: `PetscFV`, `PetscFVGetSpatialDimension()`
1336c0ce001eSMatthew G. Knepley @*/
PetscFVSetSpatialDimension(PetscFV fvm,PetscInt dim)1337d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFVSetSpatialDimension(PetscFV fvm, PetscInt dim)
1338d71ae5a4SJacob Faibussowitsch {
1339c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1340c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1);
1341c0ce001eSMatthew G. Knepley   fvm->dim = dim;
13423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1343c0ce001eSMatthew G. Knepley }
1344c0ce001eSMatthew G. Knepley 
1345c0ce001eSMatthew G. Knepley /*@
1346dce8aebaSBarry Smith   PetscFVGetSpatialDimension - Get the spatial dimension of a `PetscFV`
1347c0ce001eSMatthew G. Knepley 
134820f4b53cSBarry Smith   Not Collective
1349c0ce001eSMatthew G. Knepley 
1350c0ce001eSMatthew G. Knepley   Input Parameter:
1351dce8aebaSBarry Smith . fvm - the `PetscFV` object
1352c0ce001eSMatthew G. Knepley 
1353c0ce001eSMatthew G. Knepley   Output Parameter:
1354c0ce001eSMatthew G. Knepley . dim - The spatial dimension
1355c0ce001eSMatthew G. Knepley 
135688f5f89eSMatthew G. Knepley   Level: intermediate
1357c0ce001eSMatthew G. Knepley 
1358dce8aebaSBarry Smith .seealso: `PetscFV`, `PetscFVSetSpatialDimension()`
1359c0ce001eSMatthew G. Knepley @*/
PetscFVGetSpatialDimension(PetscFV fvm,PetscInt * dim)1360d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFVGetSpatialDimension(PetscFV fvm, PetscInt *dim)
1361d71ae5a4SJacob Faibussowitsch {
1362c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1363c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1);
13644f572ea9SToby Isaac   PetscAssertPointer(dim, 2);
1365c0ce001eSMatthew G. Knepley   *dim = fvm->dim;
13663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1367c0ce001eSMatthew G. Knepley }
1368c0ce001eSMatthew G. Knepley 
1369c0ce001eSMatthew G. Knepley /*@
1370dce8aebaSBarry Smith   PetscFVSetComputeGradients - Toggle computation of cell gradients on a `PetscFV`
1371c0ce001eSMatthew G. Knepley 
137220f4b53cSBarry Smith   Logically Collective
1373c0ce001eSMatthew G. Knepley 
1374c0ce001eSMatthew G. Knepley   Input Parameters:
1375dce8aebaSBarry Smith + fvm              - the `PetscFV` object
1376c0ce001eSMatthew G. Knepley - computeGradients - Flag to compute cell gradients
1377c0ce001eSMatthew G. Knepley 
137888f5f89eSMatthew G. Knepley   Level: intermediate
1379c0ce001eSMatthew G. Knepley 
1380dce8aebaSBarry Smith .seealso: `PetscFV`, `PetscFVGetComputeGradients()`
1381c0ce001eSMatthew G. Knepley @*/
PetscFVSetComputeGradients(PetscFV fvm,PetscBool computeGradients)1382d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFVSetComputeGradients(PetscFV fvm, PetscBool computeGradients)
1383d71ae5a4SJacob Faibussowitsch {
1384c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1385c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1);
1386c0ce001eSMatthew G. Knepley   fvm->computeGradients = computeGradients;
13873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1388c0ce001eSMatthew G. Knepley }
1389c0ce001eSMatthew G. Knepley 
1390c0ce001eSMatthew G. Knepley /*@
1391dce8aebaSBarry Smith   PetscFVGetComputeGradients - Return flag for computation of cell gradients on a `PetscFV`
1392c0ce001eSMatthew G. Knepley 
139320f4b53cSBarry Smith   Not Collective
1394c0ce001eSMatthew G. Knepley 
1395c0ce001eSMatthew G. Knepley   Input Parameter:
1396dce8aebaSBarry Smith . fvm - the `PetscFV` object
1397c0ce001eSMatthew G. Knepley 
1398c0ce001eSMatthew G. Knepley   Output Parameter:
1399c0ce001eSMatthew G. Knepley . computeGradients - Flag to compute cell gradients
1400c0ce001eSMatthew G. Knepley 
140188f5f89eSMatthew G. Knepley   Level: intermediate
1402c0ce001eSMatthew G. Knepley 
1403dce8aebaSBarry Smith .seealso: `PetscFV`, `PetscFVSetComputeGradients()`
1404c0ce001eSMatthew G. Knepley @*/
PetscFVGetComputeGradients(PetscFV fvm,PetscBool * computeGradients)1405d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFVGetComputeGradients(PetscFV fvm, PetscBool *computeGradients)
1406d71ae5a4SJacob Faibussowitsch {
1407c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1408c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1);
14094f572ea9SToby Isaac   PetscAssertPointer(computeGradients, 2);
1410c0ce001eSMatthew G. Knepley   *computeGradients = fvm->computeGradients;
14113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1412c0ce001eSMatthew G. Knepley }
1413c0ce001eSMatthew G. Knepley 
1414c0ce001eSMatthew G. Knepley /*@
1415dce8aebaSBarry Smith   PetscFVSetQuadrature - Set the `PetscQuadrature` object for a `PetscFV`
1416c0ce001eSMatthew G. Knepley 
141720f4b53cSBarry Smith   Logically Collective
1418c0ce001eSMatthew G. Knepley 
1419c0ce001eSMatthew G. Knepley   Input Parameters:
1420dce8aebaSBarry Smith + fvm - the `PetscFV` object
1421dce8aebaSBarry Smith - q   - The `PetscQuadrature`
1422c0ce001eSMatthew G. Knepley 
142388f5f89eSMatthew G. Knepley   Level: intermediate
1424c0ce001eSMatthew G. Knepley 
1425dce8aebaSBarry Smith .seealso: `PetscQuadrature`, `PetscFV`, `PetscFVGetQuadrature()`
1426c0ce001eSMatthew G. Knepley @*/
PetscFVSetQuadrature(PetscFV fvm,PetscQuadrature q)1427d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFVSetQuadrature(PetscFV fvm, PetscQuadrature q)
1428d71ae5a4SJacob Faibussowitsch {
1429c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1430c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1);
14319566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)q));
1432f6feae9bSMatthew G. Knepley   PetscCall(PetscQuadratureDestroy(&fvm->quadrature));
1433c0ce001eSMatthew G. Knepley   fvm->quadrature = q;
14343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1435c0ce001eSMatthew G. Knepley }
1436c0ce001eSMatthew G. Knepley 
1437c0ce001eSMatthew G. Knepley /*@
1438dce8aebaSBarry Smith   PetscFVGetQuadrature - Get the `PetscQuadrature` from a `PetscFV`
1439c0ce001eSMatthew G. Knepley 
144020f4b53cSBarry Smith   Not Collective
1441c0ce001eSMatthew G. Knepley 
1442c0ce001eSMatthew G. Knepley   Input Parameter:
1443dce8aebaSBarry Smith . fvm - the `PetscFV` object
1444c0ce001eSMatthew G. Knepley 
1445c0ce001eSMatthew G. Knepley   Output Parameter:
144660225df5SJacob Faibussowitsch . q - The `PetscQuadrature`
1447c0ce001eSMatthew G. Knepley 
144888f5f89eSMatthew G. Knepley   Level: intermediate
1449c0ce001eSMatthew G. Knepley 
1450dce8aebaSBarry Smith .seealso: `PetscQuadrature`, `PetscFV`, `PetscFVSetQuadrature()`
1451c0ce001eSMatthew G. Knepley @*/
PetscFVGetQuadrature(PetscFV fvm,PetscQuadrature * q)1452d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFVGetQuadrature(PetscFV fvm, PetscQuadrature *q)
1453d71ae5a4SJacob Faibussowitsch {
1454c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1455c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1);
14564f572ea9SToby Isaac   PetscAssertPointer(q, 2);
1457c0ce001eSMatthew G. Knepley   if (!fvm->quadrature) {
1458c0ce001eSMatthew G. Knepley     /* Create default 1-point quadrature */
1459c0ce001eSMatthew G. Knepley     PetscReal *points, *weights;
1460c0ce001eSMatthew G. Knepley 
14619566063dSJacob Faibussowitsch     PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, &fvm->quadrature));
14629566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(fvm->dim, &points));
14639566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(1, &weights));
1464c0ce001eSMatthew G. Knepley     weights[0] = 1.0;
14659566063dSJacob Faibussowitsch     PetscCall(PetscQuadratureSetData(fvm->quadrature, fvm->dim, 1, 1, points, weights));
1466c0ce001eSMatthew G. Knepley   }
1467c0ce001eSMatthew G. Knepley   *q = fvm->quadrature;
14683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1469c0ce001eSMatthew G. Knepley }
1470c0ce001eSMatthew G. Knepley 
1471c0ce001eSMatthew G. Knepley /*@
14722f86f8c5SMatthew G. Knepley   PetscFVCreateDualSpace - Creates a `PetscDualSpace` appropriate for the `PetscFV`
14732f86f8c5SMatthew G. Knepley 
14742f86f8c5SMatthew G. Knepley   Not Collective
14752f86f8c5SMatthew G. Knepley 
14769cde84edSJose E. Roman   Input Parameters:
14772f86f8c5SMatthew G. Knepley + fvm - The `PetscFV` object
14782f86f8c5SMatthew G. Knepley - ct  - The `DMPolytopeType` for the cell
14792f86f8c5SMatthew G. Knepley 
14802f86f8c5SMatthew G. Knepley   Level: intermediate
14812f86f8c5SMatthew G. Knepley 
14822f86f8c5SMatthew G. Knepley .seealso: `PetscFVGetDualSpace()`, `PetscFVSetDualSpace()`, `PetscDualSpace`, `PetscFV`, `PetscFVCreate()`
14832f86f8c5SMatthew G. Knepley @*/
PetscFVCreateDualSpace(PetscFV fvm,DMPolytopeType ct)14842f86f8c5SMatthew G. Knepley PetscErrorCode PetscFVCreateDualSpace(PetscFV fvm, DMPolytopeType ct)
14852f86f8c5SMatthew G. Knepley {
14862f86f8c5SMatthew G. Knepley   DM       K;
14872f86f8c5SMatthew G. Knepley   PetscInt dim, Nc;
14882f86f8c5SMatthew G. Knepley 
14892f86f8c5SMatthew G. Knepley   PetscFunctionBegin;
14902f86f8c5SMatthew G. Knepley   PetscCall(PetscFVGetSpatialDimension(fvm, &dim));
14912f86f8c5SMatthew G. Knepley   PetscCall(PetscFVGetNumComponents(fvm, &Nc));
14922f86f8c5SMatthew G. Knepley   PetscCall(PetscDualSpaceCreate(PetscObjectComm((PetscObject)fvm), &fvm->dualSpace));
14932f86f8c5SMatthew G. Knepley   PetscCall(PetscDualSpaceSetType(fvm->dualSpace, PETSCDUALSPACESIMPLE));
14942f86f8c5SMatthew G. Knepley   PetscCall(DMPlexCreateReferenceCell(PETSC_COMM_SELF, ct, &K));
14952f86f8c5SMatthew G. Knepley   PetscCall(PetscDualSpaceSetNumComponents(fvm->dualSpace, Nc));
14962f86f8c5SMatthew G. Knepley   PetscCall(PetscDualSpaceSetDM(fvm->dualSpace, K));
14972f86f8c5SMatthew G. Knepley   PetscCall(DMDestroy(&K));
14982f86f8c5SMatthew G. Knepley   PetscCall(PetscDualSpaceSimpleSetDimension(fvm->dualSpace, Nc));
14992f86f8c5SMatthew G. Knepley   // Should we be using PetscFVGetQuadrature() here?
15002f86f8c5SMatthew G. Knepley   for (PetscInt c = 0; c < Nc; ++c) {
15012f86f8c5SMatthew G. Knepley     PetscQuadrature qc;
15022f86f8c5SMatthew G. Knepley     PetscReal      *points, *weights;
15032f86f8c5SMatthew G. Knepley 
15042f86f8c5SMatthew G. Knepley     PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, &qc));
15052f86f8c5SMatthew G. Knepley     PetscCall(PetscCalloc1(dim, &points));
15062f86f8c5SMatthew G. Knepley     PetscCall(PetscCalloc1(Nc, &weights));
15072f86f8c5SMatthew G. Knepley     weights[c] = 1.0;
15082f86f8c5SMatthew G. Knepley     PetscCall(PetscQuadratureSetData(qc, dim, Nc, 1, points, weights));
15092f86f8c5SMatthew G. Knepley     PetscCall(PetscDualSpaceSimpleSetFunctional(fvm->dualSpace, c, qc));
15102f86f8c5SMatthew G. Knepley     PetscCall(PetscQuadratureDestroy(&qc));
15112f86f8c5SMatthew G. Knepley   }
15122f86f8c5SMatthew G. Knepley   PetscCall(PetscDualSpaceSetUp(fvm->dualSpace));
15132f86f8c5SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
15142f86f8c5SMatthew G. Knepley }
15152f86f8c5SMatthew G. Knepley 
15162f86f8c5SMatthew G. Knepley /*@
1517dce8aebaSBarry Smith   PetscFVGetDualSpace - Returns the `PetscDualSpace` used to define the inner product on a `PetscFV`
1518c0ce001eSMatthew G. Knepley 
151920f4b53cSBarry Smith   Not Collective
1520c0ce001eSMatthew G. Knepley 
1521c0ce001eSMatthew G. Knepley   Input Parameter:
1522dce8aebaSBarry Smith . fvm - The `PetscFV` object
1523c0ce001eSMatthew G. Knepley 
1524c0ce001eSMatthew G. Knepley   Output Parameter:
152520f4b53cSBarry Smith . sp - The `PetscDualSpace` object
1526c0ce001eSMatthew G. Knepley 
152788f5f89eSMatthew G. Knepley   Level: intermediate
1528c0ce001eSMatthew G. Knepley 
152960225df5SJacob Faibussowitsch   Developer Notes:
1530dce8aebaSBarry Smith   There is overlap between the methods of `PetscFE` and `PetscFV`, they should probably share a common parent class
1531dce8aebaSBarry Smith 
15322f86f8c5SMatthew G. Knepley .seealso: `PetscFVSetDualSpace()`, `PetscFVCreateDualSpace()`, `PetscDualSpace`, `PetscFV`, `PetscFVCreate()`
1533c0ce001eSMatthew G. Knepley @*/
PetscFVGetDualSpace(PetscFV fvm,PetscDualSpace * sp)1534d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFVGetDualSpace(PetscFV fvm, PetscDualSpace *sp)
1535d71ae5a4SJacob Faibussowitsch {
1536c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1537c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1);
15384f572ea9SToby Isaac   PetscAssertPointer(sp, 2);
1539c0ce001eSMatthew G. Knepley   if (!fvm->dualSpace) {
15402f86f8c5SMatthew G. Knepley     PetscInt dim;
1541c0ce001eSMatthew G. Knepley 
15429566063dSJacob Faibussowitsch     PetscCall(PetscFVGetSpatialDimension(fvm, &dim));
15432f86f8c5SMatthew G. Knepley     PetscCall(PetscFVCreateDualSpace(fvm, DMPolytopeTypeSimpleShape(dim, PETSC_FALSE)));
1544c0ce001eSMatthew G. Knepley   }
1545c0ce001eSMatthew G. Knepley   *sp = fvm->dualSpace;
15463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1547c0ce001eSMatthew G. Knepley }
1548c0ce001eSMatthew G. Knepley 
1549c0ce001eSMatthew G. Knepley /*@
1550dce8aebaSBarry Smith   PetscFVSetDualSpace - Sets the `PetscDualSpace` used to define the inner product
1551c0ce001eSMatthew G. Knepley 
155220f4b53cSBarry Smith   Not Collective
1553c0ce001eSMatthew G. Knepley 
1554c0ce001eSMatthew G. Knepley   Input Parameters:
1555dce8aebaSBarry Smith + fvm - The `PetscFV` object
1556dce8aebaSBarry Smith - sp  - The `PetscDualSpace` object
1557c0ce001eSMatthew G. Knepley 
1558c0ce001eSMatthew G. Knepley   Level: intermediate
1559c0ce001eSMatthew G. Knepley 
1560dce8aebaSBarry Smith   Note:
1561dce8aebaSBarry Smith   A simple dual space is provided automatically, and the user typically will not need to override it.
1562c0ce001eSMatthew G. Knepley 
15632f86f8c5SMatthew G. Knepley .seealso: `PetscFVGetDualSpace()`, `PetscFVCreateDualSpace()`, `PetscDualSpace`, `PetscFV`, `PetscFVCreate()`
1564c0ce001eSMatthew G. Knepley @*/
PetscFVSetDualSpace(PetscFV fvm,PetscDualSpace sp)1565d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFVSetDualSpace(PetscFV fvm, PetscDualSpace sp)
1566d71ae5a4SJacob Faibussowitsch {
1567c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1568c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1);
1569c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 2);
15709566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceDestroy(&fvm->dualSpace));
1571c0ce001eSMatthew G. Knepley   fvm->dualSpace = sp;
15729566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)fvm->dualSpace));
15733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1574c0ce001eSMatthew G. Knepley }
1575c0ce001eSMatthew G. Knepley 
157688f5f89eSMatthew G. Knepley /*@C
1577ef0bb6c7SMatthew G. Knepley   PetscFVGetCellTabulation - Returns the tabulation of the basis functions at the quadrature points
157888f5f89eSMatthew G. Knepley 
157920f4b53cSBarry Smith   Not Collective
158088f5f89eSMatthew G. Knepley 
158188f5f89eSMatthew G. Knepley   Input Parameter:
1582dce8aebaSBarry Smith . fvm - The `PetscFV` object
158388f5f89eSMatthew G. Knepley 
1584ef0bb6c7SMatthew G. Knepley   Output Parameter:
1585a5b23f4aSJose E. Roman . T - The basis function values and derivatives at quadrature points
158688f5f89eSMatthew G. Knepley 
158788f5f89eSMatthew G. Knepley   Level: intermediate
158888f5f89eSMatthew G. Knepley 
1589dce8aebaSBarry Smith   Note:
1590dce8aebaSBarry Smith .vb
1591dce8aebaSBarry Smith   T->T[0] = B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c
1592dce8aebaSBarry Smith   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
1593dce8aebaSBarry Smith   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
1594dce8aebaSBarry Smith .ve
1595dce8aebaSBarry Smith 
1596dce8aebaSBarry Smith .seealso: `PetscFV`, `PetscTabulation`, `PetscFEGetCellTabulation()`, `PetscFVCreateTabulation()`, `PetscFVGetQuadrature()`, `PetscQuadratureGetData()`
159788f5f89eSMatthew G. Knepley @*/
PetscFVGetCellTabulation(PetscFV fvm,PetscTabulation * T)1598d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFVGetCellTabulation(PetscFV fvm, PetscTabulation *T)
1599d71ae5a4SJacob Faibussowitsch {
1600c0ce001eSMatthew G. Knepley   PetscInt         npoints;
1601c0ce001eSMatthew G. Knepley   const PetscReal *points;
1602c0ce001eSMatthew G. Knepley 
1603c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1604c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1);
16054f572ea9SToby Isaac   PetscAssertPointer(T, 2);
16069566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(fvm->quadrature, NULL, NULL, &npoints, &points, NULL));
16079566063dSJacob Faibussowitsch   if (!fvm->T) PetscCall(PetscFVCreateTabulation(fvm, 1, npoints, points, 1, &fvm->T));
1608ef0bb6c7SMatthew G. Knepley   *T = fvm->T;
16093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1610c0ce001eSMatthew G. Knepley }
1611c0ce001eSMatthew G. Knepley 
161288f5f89eSMatthew G. Knepley /*@C
1613ef0bb6c7SMatthew G. Knepley   PetscFVCreateTabulation - Tabulates the basis functions, and perhaps derivatives, at the points provided.
161488f5f89eSMatthew G. Knepley 
161520f4b53cSBarry Smith   Not Collective
161688f5f89eSMatthew G. Knepley 
161788f5f89eSMatthew G. Knepley   Input Parameters:
1618dce8aebaSBarry Smith + fvm     - The `PetscFV` object
1619ef0bb6c7SMatthew G. Knepley . nrepl   - The number of replicas
1620ef0bb6c7SMatthew G. Knepley . npoints - The number of tabulation points in a replica
1621ef0bb6c7SMatthew G. Knepley . points  - The tabulation point coordinates
1622ef0bb6c7SMatthew G. Knepley - K       - The order of derivative to tabulate
162388f5f89eSMatthew G. Knepley 
1624ef0bb6c7SMatthew G. Knepley   Output Parameter:
1625a5b23f4aSJose E. Roman . T - The basis function values and derivative at tabulation points
162688f5f89eSMatthew G. Knepley 
162788f5f89eSMatthew G. Knepley   Level: intermediate
162888f5f89eSMatthew G. Knepley 
1629dce8aebaSBarry Smith   Note:
1630dce8aebaSBarry Smith .vb
1631dce8aebaSBarry Smith   T->T[0] = B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c
1632dce8aebaSBarry Smith   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
1633dce8aebaSBarry Smith   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
1634dce8aebaSBarry Smith .ve
1635dce8aebaSBarry Smith 
1636dce8aebaSBarry Smith .seealso: `PetscFV`, `PetscTabulation`, `PetscFECreateTabulation()`, `PetscTabulationDestroy()`, `PetscFEGetCellTabulation()`
163788f5f89eSMatthew G. Knepley @*/
PetscFVCreateTabulation(PetscFV fvm,PetscInt nrepl,PetscInt npoints,const PetscReal points[],PetscInt K,PetscTabulation * T)1638d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFVCreateTabulation(PetscFV fvm, PetscInt nrepl, PetscInt npoints, const PetscReal points[], PetscInt K, PetscTabulation *T)
1639d71ae5a4SJacob Faibussowitsch {
16402f86f8c5SMatthew G. Knepley   PetscInt pdim; // Dimension of approximation space P
16412f86f8c5SMatthew G. Knepley   PetscInt cdim; // Spatial dimension
16422f86f8c5SMatthew G. Knepley   PetscInt Nc;   // Field components
1643ef0bb6c7SMatthew G. Knepley   PetscInt k, p, d, c, e;
1644c0ce001eSMatthew G. Knepley 
1645c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1646ef0bb6c7SMatthew G. Knepley   if (!npoints || K < 0) {
1647ef0bb6c7SMatthew G. Knepley     *T = NULL;
16483ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
1649c0ce001eSMatthew G. Knepley   }
1650c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1);
16514f572ea9SToby Isaac   PetscAssertPointer(points, 4);
16524f572ea9SToby Isaac   PetscAssertPointer(T, 6);
16539566063dSJacob Faibussowitsch   PetscCall(PetscFVGetSpatialDimension(fvm, &cdim));
16549566063dSJacob Faibussowitsch   PetscCall(PetscFVGetNumComponents(fvm, &Nc));
16552f86f8c5SMatthew G. Knepley   pdim = Nc;
16569566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(1, T));
1657ef0bb6c7SMatthew G. Knepley   (*T)->K    = !cdim ? 0 : K;
1658ef0bb6c7SMatthew G. Knepley   (*T)->Nr   = nrepl;
1659ef0bb6c7SMatthew G. Knepley   (*T)->Np   = npoints;
1660ef0bb6c7SMatthew G. Knepley   (*T)->Nb   = pdim;
1661ef0bb6c7SMatthew G. Knepley   (*T)->Nc   = Nc;
1662ef0bb6c7SMatthew G. Knepley   (*T)->cdim = cdim;
16639566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1((*T)->K + 1, &(*T)->T));
166448a46eb9SPierre Jolivet   for (k = 0; k <= (*T)->K; ++k) PetscCall(PetscMalloc1(nrepl * npoints * pdim * Nc * PetscPowInt(cdim, k), &(*T)->T[k]));
16659371c9d4SSatish Balay   if (K >= 0) {
16669371c9d4SSatish Balay     for (p = 0; p < nrepl * npoints; ++p)
16679371c9d4SSatish Balay       for (d = 0; d < pdim; ++d)
16682f86f8c5SMatthew G. Knepley         for (c = 0; c < Nc; ++c) (*T)->T[0][(p * pdim + d) * Nc + c] = 1.;
1669ef0bb6c7SMatthew G. Knepley   }
16709371c9d4SSatish Balay   if (K >= 1) {
16719371c9d4SSatish Balay     for (p = 0; p < nrepl * npoints; ++p)
16729371c9d4SSatish Balay       for (d = 0; d < pdim; ++d)
16739371c9d4SSatish Balay         for (c = 0; c < Nc; ++c)
16749371c9d4SSatish Balay           for (e = 0; e < cdim; ++e) (*T)->T[1][((p * pdim + d) * Nc + c) * cdim + e] = 0.0;
16759371c9d4SSatish Balay   }
16769371c9d4SSatish Balay   if (K >= 2) {
16779371c9d4SSatish Balay     for (p = 0; p < nrepl * npoints; ++p)
16789371c9d4SSatish Balay       for (d = 0; d < pdim; ++d)
16799371c9d4SSatish Balay         for (c = 0; c < Nc; ++c)
16809371c9d4SSatish Balay           for (e = 0; e < cdim * cdim; ++e) (*T)->T[2][((p * pdim + d) * Nc + c) * cdim * cdim + e] = 0.0;
16819371c9d4SSatish Balay   }
16823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1683c0ce001eSMatthew G. Knepley }
1684c0ce001eSMatthew G. Knepley 
1685cc4c1da9SBarry Smith /*@
1686c0ce001eSMatthew G. Knepley   PetscFVComputeGradient - Compute the gradient reconstruction matrix for a given cell
1687c0ce001eSMatthew G. Knepley 
1688c0ce001eSMatthew G. Knepley   Input Parameters:
1689dce8aebaSBarry Smith + fvm      - The `PetscFV` object
1690c0ce001eSMatthew G. Knepley . numFaces - The number of cell faces which are not constrained
1691c0ce001eSMatthew G. Knepley - dx       - The vector from the cell centroid to the neighboring cell centroid for each face
1692c0ce001eSMatthew G. Knepley 
1693a4e35b19SJacob Faibussowitsch   Output Parameter:
1694a4e35b19SJacob Faibussowitsch . grad - the gradient
1695a4e35b19SJacob Faibussowitsch 
169688f5f89eSMatthew G. Knepley   Level: advanced
1697c0ce001eSMatthew G. Knepley 
1698dce8aebaSBarry Smith .seealso: `PetscFV`, `PetscFVCreate()`
1699c0ce001eSMatthew G. Knepley @*/
PetscFVComputeGradient(PetscFV fvm,PetscInt numFaces,PetscScalar dx[],PetscScalar grad[])1700d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFVComputeGradient(PetscFV fvm, PetscInt numFaces, PetscScalar dx[], PetscScalar grad[])
1701d71ae5a4SJacob Faibussowitsch {
1702c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1703c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1);
1704dbbe0bcdSBarry Smith   PetscTryTypeMethod(fvm, computegradient, numFaces, dx, grad);
17053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1706c0ce001eSMatthew G. Knepley }
1707c0ce001eSMatthew G. Knepley 
170888f5f89eSMatthew G. Knepley /*@C
1709c0ce001eSMatthew G. Knepley   PetscFVIntegrateRHSFunction - Produce the cell residual vector for a chunk of elements by quadrature integration
1710c0ce001eSMatthew G. Knepley 
171120f4b53cSBarry Smith   Not Collective
1712c0ce001eSMatthew G. Knepley 
1713c0ce001eSMatthew G. Knepley   Input Parameters:
1714dce8aebaSBarry Smith + fvm         - The `PetscFV` object for the field being integrated
1715da81f932SPierre Jolivet . prob        - The `PetscDS` specifying the discretizations and continuum functions
1716c0ce001eSMatthew G. Knepley . field       - The field being integrated
1717c0ce001eSMatthew G. Knepley . Nf          - The number of faces in the chunk
1718c0ce001eSMatthew G. Knepley . fgeom       - The face geometry for each face in the chunk
1719c0ce001eSMatthew G. Knepley . neighborVol - The volume for each pair of cells in the chunk
1720c0ce001eSMatthew G. Knepley . uL          - The state from the cell on the left
1721c0ce001eSMatthew G. Knepley - uR          - The state from the cell on the right
1722c0ce001eSMatthew G. Knepley 
1723d8d19677SJose E. Roman   Output Parameters:
1724c0ce001eSMatthew G. Knepley + fluxL - the left fluxes for each face
1725c0ce001eSMatthew G. Knepley - fluxR - the right fluxes for each face
172688f5f89eSMatthew G. Knepley 
172788f5f89eSMatthew G. Knepley   Level: developer
172888f5f89eSMatthew G. Knepley 
1729dce8aebaSBarry Smith .seealso: `PetscFV`, `PetscDS`, `PetscFVFaceGeom`, `PetscFVCreate()`
173088f5f89eSMatthew G. Knepley @*/
PetscFVIntegrateRHSFunction(PetscFV fvm,PetscDS prob,PetscInt field,PetscInt Nf,PetscFVFaceGeom * fgeom,PetscReal * neighborVol,PetscScalar uL[],PetscScalar uR[],PetscScalar fluxL[],PetscScalar fluxR[])1731d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFVIntegrateRHSFunction(PetscFV fvm, PetscDS prob, PetscInt field, PetscInt Nf, PetscFVFaceGeom *fgeom, PetscReal *neighborVol, PetscScalar uL[], PetscScalar uR[], PetscScalar fluxL[], PetscScalar fluxR[])
1732d71ae5a4SJacob Faibussowitsch {
1733c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1734c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1);
1735dbbe0bcdSBarry Smith   PetscTryTypeMethod(fvm, integraterhsfunction, prob, field, Nf, fgeom, neighborVol, uL, uR, fluxL, fluxR);
17363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1737c0ce001eSMatthew G. Knepley }
1738c0ce001eSMatthew G. Knepley 
1739c0ce001eSMatthew G. Knepley /*@
1740d8b4a066SPierre Jolivet   PetscFVClone - Create a shallow copy of a `PetscFV` object that just references the internal objects.
1741f6feae9bSMatthew G. Knepley 
1742f6feae9bSMatthew G. Knepley   Input Parameter:
1743f6feae9bSMatthew G. Knepley . fv - The initial `PetscFV`
1744f6feae9bSMatthew G. Knepley 
1745f6feae9bSMatthew G. Knepley   Output Parameter:
1746f6feae9bSMatthew G. Knepley . fvNew - A clone of the `PetscFV`
1747f6feae9bSMatthew G. Knepley 
1748f6feae9bSMatthew G. Knepley   Level: advanced
1749f6feae9bSMatthew G. Knepley 
1750f6feae9bSMatthew G. Knepley   Notes:
1751f6feae9bSMatthew G. Knepley   This is typically used to change the number of components.
1752f6feae9bSMatthew G. Knepley 
1753f6feae9bSMatthew G. Knepley .seealso: `PetscFV`, `PetscFVType`, `PetscFVCreate()`, `PetscFVSetType()`
1754f6feae9bSMatthew G. Knepley @*/
PetscFVClone(PetscFV fv,PetscFV * fvNew)1755f6feae9bSMatthew G. Knepley PetscErrorCode PetscFVClone(PetscFV fv, PetscFV *fvNew)
1756f6feae9bSMatthew G. Knepley {
1757f6feae9bSMatthew G. Knepley   PetscDualSpace  Q;
1758f6feae9bSMatthew G. Knepley   DM              K;
1759f6feae9bSMatthew G. Knepley   PetscQuadrature q;
1760f6feae9bSMatthew G. Knepley   PetscInt        Nc, cdim;
1761f6feae9bSMatthew G. Knepley 
1762f6feae9bSMatthew G. Knepley   PetscFunctionBegin;
1763f6feae9bSMatthew G. Knepley   PetscCall(PetscFVGetDualSpace(fv, &Q));
1764f6feae9bSMatthew G. Knepley   PetscCall(PetscFVGetQuadrature(fv, &q));
1765f6feae9bSMatthew G. Knepley   PetscCall(PetscDualSpaceGetDM(Q, &K));
1766f6feae9bSMatthew G. Knepley 
1767f6feae9bSMatthew G. Knepley   PetscCall(PetscFVCreate(PetscObjectComm((PetscObject)fv), fvNew));
1768f6feae9bSMatthew G. Knepley   PetscCall(PetscFVSetDualSpace(*fvNew, Q));
1769f6feae9bSMatthew G. Knepley   PetscCall(PetscFVGetNumComponents(fv, &Nc));
1770f6feae9bSMatthew G. Knepley   PetscCall(PetscFVSetNumComponents(*fvNew, Nc));
1771f6feae9bSMatthew G. Knepley   PetscCall(PetscFVGetSpatialDimension(fv, &cdim));
1772f6feae9bSMatthew G. Knepley   PetscCall(PetscFVSetSpatialDimension(*fvNew, cdim));
1773f6feae9bSMatthew G. Knepley   PetscCall(PetscFVSetQuadrature(*fvNew, q));
1774f6feae9bSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
1775f6feae9bSMatthew G. Knepley }
1776f6feae9bSMatthew G. Knepley 
1777f6feae9bSMatthew G. Knepley /*@
1778a4e35b19SJacob Faibussowitsch   PetscFVRefine - Create a "refined" `PetscFV` object that refines the reference cell into
1779a4e35b19SJacob Faibussowitsch   smaller copies.
1780c0ce001eSMatthew G. Knepley 
1781c0ce001eSMatthew G. Knepley   Input Parameter:
1782dce8aebaSBarry Smith . fv - The initial `PetscFV`
1783c0ce001eSMatthew G. Knepley 
1784c0ce001eSMatthew G. Knepley   Output Parameter:
1785dce8aebaSBarry Smith . fvRef - The refined `PetscFV`
1786c0ce001eSMatthew G. Knepley 
178788f5f89eSMatthew G. Knepley   Level: advanced
1788c0ce001eSMatthew G. Knepley 
1789a4e35b19SJacob Faibussowitsch   Notes:
1790a4e35b19SJacob Faibussowitsch   This is typically used to generate a preconditioner for a high order method from a lower order method on a
1791a4e35b19SJacob Faibussowitsch   refined mesh having the same number of dofs (but more sparsity). It is also used to create an
1792a4e35b19SJacob Faibussowitsch   interpolation between regularly refined meshes.
1793a4e35b19SJacob Faibussowitsch 
1794dce8aebaSBarry Smith .seealso: `PetscFV`, `PetscFVType`, `PetscFVCreate()`, `PetscFVSetType()`
1795c0ce001eSMatthew G. Knepley @*/
PetscFVRefine(PetscFV fv,PetscFV * fvRef)1796d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFVRefine(PetscFV fv, PetscFV *fvRef)
1797d71ae5a4SJacob Faibussowitsch {
1798c0ce001eSMatthew G. Knepley   PetscDualSpace  Q, Qref;
1799c0ce001eSMatthew G. Knepley   DM              K, Kref;
1800c0ce001eSMatthew G. Knepley   PetscQuadrature q, qref;
1801412e9a14SMatthew G. Knepley   DMPolytopeType  ct;
1802012bc364SMatthew G. Knepley   DMPlexTransform tr;
1803c0ce001eSMatthew G. Knepley   PetscReal      *v0;
1804c0ce001eSMatthew G. Knepley   PetscReal      *jac, *invjac;
1805c0ce001eSMatthew G. Knepley   PetscInt        numComp, numSubelements, s;
1806c0ce001eSMatthew G. Knepley 
1807c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
18089566063dSJacob Faibussowitsch   PetscCall(PetscFVGetDualSpace(fv, &Q));
18099566063dSJacob Faibussowitsch   PetscCall(PetscFVGetQuadrature(fv, &q));
18109566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetDM(Q, &K));
1811c0ce001eSMatthew G. Knepley   /* Create dual space */
18129566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceDuplicate(Q, &Qref));
18139566063dSJacob Faibussowitsch   PetscCall(DMRefine(K, PetscObjectComm((PetscObject)fv), &Kref));
18149566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceSetDM(Qref, Kref));
18159566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&Kref));
18169566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceSetUp(Qref));
1817c0ce001eSMatthew G. Knepley   /* Create volume */
18189566063dSJacob Faibussowitsch   PetscCall(PetscFVCreate(PetscObjectComm((PetscObject)fv), fvRef));
18199566063dSJacob Faibussowitsch   PetscCall(PetscFVSetDualSpace(*fvRef, Qref));
18209566063dSJacob Faibussowitsch   PetscCall(PetscFVGetNumComponents(fv, &numComp));
18219566063dSJacob Faibussowitsch   PetscCall(PetscFVSetNumComponents(*fvRef, numComp));
18229566063dSJacob Faibussowitsch   PetscCall(PetscFVSetUp(*fvRef));
1823c0ce001eSMatthew G. Knepley   /* Create quadrature */
18249566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCellType(K, 0, &ct));
18259566063dSJacob Faibussowitsch   PetscCall(DMPlexTransformCreate(PETSC_COMM_SELF, &tr));
18269566063dSJacob Faibussowitsch   PetscCall(DMPlexTransformSetType(tr, DMPLEXREFINEREGULAR));
18279566063dSJacob Faibussowitsch   PetscCall(DMPlexRefineRegularGetAffineTransforms(tr, ct, &numSubelements, &v0, &jac, &invjac));
18289566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureExpandComposite(q, numSubelements, v0, jac, &qref));
18299566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceSimpleSetDimension(Qref, numSubelements));
1830c0ce001eSMatthew G. Knepley   for (s = 0; s < numSubelements; ++s) {
1831c0ce001eSMatthew G. Knepley     PetscQuadrature  qs;
1832c0ce001eSMatthew G. Knepley     const PetscReal *points, *weights;
1833c0ce001eSMatthew G. Knepley     PetscReal       *p, *w;
1834c0ce001eSMatthew G. Knepley     PetscInt         dim, Nc, npoints, np;
1835c0ce001eSMatthew G. Knepley 
18369566063dSJacob Faibussowitsch     PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, &qs));
18379566063dSJacob Faibussowitsch     PetscCall(PetscQuadratureGetData(q, &dim, &Nc, &npoints, &points, &weights));
1838c0ce001eSMatthew G. Knepley     np = npoints / numSubelements;
18399566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(np * dim, &p));
18409566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(np * Nc, &w));
18419566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(p, &points[s * np * dim], np * dim));
18429566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(w, &weights[s * np * Nc], np * Nc));
18439566063dSJacob Faibussowitsch     PetscCall(PetscQuadratureSetData(qs, dim, Nc, np, p, w));
18449566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceSimpleSetFunctional(Qref, s, qs));
18459566063dSJacob Faibussowitsch     PetscCall(PetscQuadratureDestroy(&qs));
1846c0ce001eSMatthew G. Knepley   }
18479566063dSJacob Faibussowitsch   PetscCall(PetscFVSetQuadrature(*fvRef, qref));
18489566063dSJacob Faibussowitsch   PetscCall(DMPlexTransformDestroy(&tr));
18499566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureDestroy(&qref));
18509566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceDestroy(&Qref));
18513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1852c0ce001eSMatthew G. Knepley }
1853c0ce001eSMatthew G. Knepley 
PetscFVDestroy_Upwind(PetscFV fvm)1854d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscFVDestroy_Upwind(PetscFV fvm)
1855d71ae5a4SJacob Faibussowitsch {
1856c0ce001eSMatthew G. Knepley   PetscFV_Upwind *b = (PetscFV_Upwind *)fvm->data;
1857c0ce001eSMatthew G. Knepley 
1858c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
18599566063dSJacob Faibussowitsch   PetscCall(PetscFree(b));
18603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1861c0ce001eSMatthew G. Knepley }
1862c0ce001eSMatthew G. Knepley 
PetscFVView_Upwind_Ascii(PetscFV fv,PetscViewer viewer)1863d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscFVView_Upwind_Ascii(PetscFV fv, PetscViewer viewer)
1864d71ae5a4SJacob Faibussowitsch {
1865c0ce001eSMatthew G. Knepley   PetscInt          Nc, c;
1866c0ce001eSMatthew G. Knepley   PetscViewerFormat format;
1867c0ce001eSMatthew G. Knepley 
1868c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
18699566063dSJacob Faibussowitsch   PetscCall(PetscFVGetNumComponents(fv, &Nc));
18709566063dSJacob Faibussowitsch   PetscCall(PetscViewerGetFormat(viewer, &format));
18719566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Upwind Finite Volume:\n"));
187263a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "  num components: %" PetscInt_FMT "\n", Nc));
1873c0ce001eSMatthew G. Knepley   for (c = 0; c < Nc; c++) {
187448a46eb9SPierre Jolivet     if (fv->componentNames[c]) PetscCall(PetscViewerASCIIPrintf(viewer, "    component %" PetscInt_FMT ": %s\n", c, fv->componentNames[c]));
1875c0ce001eSMatthew G. Knepley   }
18763ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1877c0ce001eSMatthew G. Knepley }
1878c0ce001eSMatthew G. Knepley 
PetscFVView_Upwind(PetscFV fv,PetscViewer viewer)1879d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscFVView_Upwind(PetscFV fv, PetscViewer viewer)
1880d71ae5a4SJacob Faibussowitsch {
1881*9f196a02SMartin Diehl   PetscBool isascii;
1882c0ce001eSMatthew G. Knepley 
1883c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1884c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(fv, PETSCFV_CLASSID, 1);
1885c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
1886*9f196a02SMartin Diehl   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
1887*9f196a02SMartin Diehl   if (isascii) PetscCall(PetscFVView_Upwind_Ascii(fv, viewer));
18883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1889c0ce001eSMatthew G. Knepley }
1890c0ce001eSMatthew G. Knepley 
PetscFVComputeGradient_Upwind(PetscFV fv,PetscInt numFaces,const PetscScalar dx[],PetscScalar grad[])18916f6f020cSMatthew G. Knepley static PetscErrorCode PetscFVComputeGradient_Upwind(PetscFV fv, PetscInt numFaces, const PetscScalar dx[], PetscScalar grad[])
18926f6f020cSMatthew G. Knepley {
18936f6f020cSMatthew G. Knepley   PetscInt dim;
18946f6f020cSMatthew G. Knepley 
18956f6f020cSMatthew G. Knepley   PetscFunctionBegin;
18966f6f020cSMatthew G. Knepley   PetscCall(PetscFVGetSpatialDimension(fv, &dim));
18976f6f020cSMatthew G. Knepley   for (PetscInt f = 0; f < numFaces; ++f) {
18986f6f020cSMatthew G. Knepley     for (PetscInt d = 0; d < dim; ++d) grad[f * dim + d] = 0.;
18996f6f020cSMatthew G. Knepley   }
19006f6f020cSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
19016f6f020cSMatthew G. Knepley }
19026f6f020cSMatthew G. Knepley 
1903c0ce001eSMatthew G. Knepley /*
1904c0ce001eSMatthew G. Knepley   neighborVol[f*2+0] contains the left  geom
1905c0ce001eSMatthew G. Knepley   neighborVol[f*2+1] contains the right geom
1906c0ce001eSMatthew G. Knepley */
PetscFVIntegrateRHSFunction_Upwind(PetscFV fvm,PetscDS prob,PetscInt field,PetscInt Nf,PetscFVFaceGeom * fgeom,PetscReal * neighborVol,PetscScalar uL[],PetscScalar uR[],PetscScalar fluxL[],PetscScalar fluxR[])1907d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscFVIntegrateRHSFunction_Upwind(PetscFV fvm, PetscDS prob, PetscInt field, PetscInt Nf, PetscFVFaceGeom *fgeom, PetscReal *neighborVol, PetscScalar uL[], PetscScalar uR[], PetscScalar fluxL[], PetscScalar fluxR[])
1908d71ae5a4SJacob Faibussowitsch {
1909c0ce001eSMatthew G. Knepley   void (*riemann)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscInt, const PetscScalar[], PetscScalar[], void *);
1910c0ce001eSMatthew G. Knepley   void              *rctx;
1911c0ce001eSMatthew G. Knepley   PetscScalar       *flux = fvm->fluxWork;
1912c0ce001eSMatthew G. Knepley   const PetscScalar *constants;
1913c0ce001eSMatthew G. Knepley   PetscInt           dim, numConstants, pdim, totDim, Nc, off, f, d;
1914c0ce001eSMatthew G. Knepley 
1915c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
19169566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalComponents(prob, &Nc));
19179566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(prob, &totDim));
19189566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffset(prob, field, &off));
19199566063dSJacob Faibussowitsch   PetscCall(PetscDSGetRiemannSolver(prob, field, &riemann));
19209566063dSJacob Faibussowitsch   PetscCall(PetscDSGetContext(prob, field, &rctx));
19219566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(prob, &numConstants, &constants));
19229566063dSJacob Faibussowitsch   PetscCall(PetscFVGetSpatialDimension(fvm, &dim));
19239566063dSJacob Faibussowitsch   PetscCall(PetscFVGetNumComponents(fvm, &pdim));
1924c0ce001eSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
1925c0ce001eSMatthew G. Knepley     (*riemann)(dim, pdim, fgeom[f].centroid, fgeom[f].normal, &uL[f * Nc], &uR[f * Nc], numConstants, constants, flux, rctx);
1926c0ce001eSMatthew G. Knepley     for (d = 0; d < pdim; ++d) {
1927c0ce001eSMatthew G. Knepley       fluxL[f * totDim + off + d] = flux[d] / neighborVol[f * 2 + 0];
1928c0ce001eSMatthew G. Knepley       fluxR[f * totDim + off + d] = flux[d] / neighborVol[f * 2 + 1];
1929c0ce001eSMatthew G. Knepley     }
1930c0ce001eSMatthew G. Knepley   }
19313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1932c0ce001eSMatthew G. Knepley }
1933c0ce001eSMatthew G. Knepley 
PetscFVInitialize_Upwind(PetscFV fvm)1934d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscFVInitialize_Upwind(PetscFV fvm)
1935d71ae5a4SJacob Faibussowitsch {
1936c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1937c0ce001eSMatthew G. Knepley   fvm->ops->setfromoptions       = NULL;
1938c0ce001eSMatthew G. Knepley   fvm->ops->view                 = PetscFVView_Upwind;
1939c0ce001eSMatthew G. Knepley   fvm->ops->destroy              = PetscFVDestroy_Upwind;
19406f6f020cSMatthew G. Knepley   fvm->ops->computegradient      = PetscFVComputeGradient_Upwind;
1941c0ce001eSMatthew G. Knepley   fvm->ops->integraterhsfunction = PetscFVIntegrateRHSFunction_Upwind;
19423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1943c0ce001eSMatthew G. Knepley }
1944c0ce001eSMatthew G. Knepley 
1945c0ce001eSMatthew G. Knepley /*MC
1946dce8aebaSBarry Smith   PETSCFVUPWIND = "upwind" - A `PetscFV` implementation
1947c0ce001eSMatthew G. Knepley 
1948c0ce001eSMatthew G. Knepley   Level: intermediate
1949c0ce001eSMatthew G. Knepley 
1950dce8aebaSBarry Smith .seealso: `PetscFV`, `PetscFVType`, `PetscFVCreate()`, `PetscFVSetType()`
1951c0ce001eSMatthew G. Knepley M*/
1952c0ce001eSMatthew G. Knepley 
PetscFVCreate_Upwind(PetscFV fvm)1953d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscFVCreate_Upwind(PetscFV fvm)
1954d71ae5a4SJacob Faibussowitsch {
1955c0ce001eSMatthew G. Knepley   PetscFV_Upwind *b;
1956c0ce001eSMatthew G. Knepley 
1957c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
1958c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1);
19594dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&b));
1960c0ce001eSMatthew G. Knepley   fvm->data = b;
1961c0ce001eSMatthew G. Knepley 
19629566063dSJacob Faibussowitsch   PetscCall(PetscFVInitialize_Upwind(fvm));
19633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1964c0ce001eSMatthew G. Knepley }
1965c0ce001eSMatthew G. Knepley 
1966c0ce001eSMatthew G. Knepley #include <petscblaslapack.h>
1967c0ce001eSMatthew G. Knepley 
PetscFVDestroy_LeastSquares(PetscFV fvm)1968d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscFVDestroy_LeastSquares(PetscFV fvm)
1969d71ae5a4SJacob Faibussowitsch {
1970c0ce001eSMatthew G. Knepley   PetscFV_LeastSquares *ls = (PetscFV_LeastSquares *)fvm->data;
1971c0ce001eSMatthew G. Knepley 
1972c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
19739566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)fvm, "PetscFVLeastSquaresSetMaxFaces_C", NULL));
19749566063dSJacob Faibussowitsch   PetscCall(PetscFree4(ls->B, ls->Binv, ls->tau, ls->work));
19759566063dSJacob Faibussowitsch   PetscCall(PetscFree(ls));
19763ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1977c0ce001eSMatthew G. Knepley }
1978c0ce001eSMatthew G. Knepley 
PetscFVView_LeastSquares_Ascii(PetscFV fv,PetscViewer viewer)1979d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscFVView_LeastSquares_Ascii(PetscFV fv, PetscViewer viewer)
1980d71ae5a4SJacob Faibussowitsch {
1981c0ce001eSMatthew G. Knepley   PetscInt          Nc, c;
1982c0ce001eSMatthew G. Knepley   PetscViewerFormat format;
1983c0ce001eSMatthew G. Knepley 
1984c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
19859566063dSJacob Faibussowitsch   PetscCall(PetscFVGetNumComponents(fv, &Nc));
19869566063dSJacob Faibussowitsch   PetscCall(PetscViewerGetFormat(viewer, &format));
19879566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Finite Volume with Least Squares Reconstruction:\n"));
198863a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "  num components: %" PetscInt_FMT "\n", Nc));
1989c0ce001eSMatthew G. Knepley   for (c = 0; c < Nc; c++) {
199048a46eb9SPierre Jolivet     if (fv->componentNames[c]) PetscCall(PetscViewerASCIIPrintf(viewer, "    component %" PetscInt_FMT ": %s\n", c, fv->componentNames[c]));
1991c0ce001eSMatthew G. Knepley   }
19923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1993c0ce001eSMatthew G. Knepley }
1994c0ce001eSMatthew G. Knepley 
PetscFVView_LeastSquares(PetscFV fv,PetscViewer viewer)1995d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscFVView_LeastSquares(PetscFV fv, PetscViewer viewer)
1996d71ae5a4SJacob Faibussowitsch {
1997*9f196a02SMartin Diehl   PetscBool isascii;
1998c0ce001eSMatthew G. Knepley 
1999c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
2000c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(fv, PETSCFV_CLASSID, 1);
2001c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
2002*9f196a02SMartin Diehl   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
2003*9f196a02SMartin Diehl   if (isascii) PetscCall(PetscFVView_LeastSquares_Ascii(fv, viewer));
20043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2005c0ce001eSMatthew G. Knepley }
2006c0ce001eSMatthew G. Knepley 
2007c0ce001eSMatthew G. Knepley /* Overwrites A. Can only handle full-rank problems with m>=n */
PetscFVLeastSquaresPseudoInverse_Static(PetscInt m,PetscInt mstride,PetscInt n,PetscScalar * A,PetscScalar * Ainv,PetscScalar * tau,PetscInt worksize,PetscScalar * work)2008d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscFVLeastSquaresPseudoInverse_Static(PetscInt m, PetscInt mstride, PetscInt n, PetscScalar *A, PetscScalar *Ainv, PetscScalar *tau, PetscInt worksize, PetscScalar *work)
2009d71ae5a4SJacob Faibussowitsch {
2010c0ce001eSMatthew G. Knepley   PetscBool    debug = PETSC_FALSE;
2011c0ce001eSMatthew G. Knepley   PetscBLASInt M, N, K, lda, ldb, ldwork, info;
2012c0ce001eSMatthew G. Knepley   PetscScalar *R, *Q, *Aback, Alpha;
2013c0ce001eSMatthew G. Knepley 
2014c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
2015c0ce001eSMatthew G. Knepley   if (debug) {
20169566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(m * n, &Aback));
20179566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(Aback, A, m * n));
2018c0ce001eSMatthew G. Knepley   }
2019c0ce001eSMatthew G. Knepley 
20209566063dSJacob Faibussowitsch   PetscCall(PetscBLASIntCast(m, &M));
20219566063dSJacob Faibussowitsch   PetscCall(PetscBLASIntCast(n, &N));
20229566063dSJacob Faibussowitsch   PetscCall(PetscBLASIntCast(mstride, &lda));
20239566063dSJacob Faibussowitsch   PetscCall(PetscBLASIntCast(worksize, &ldwork));
20249566063dSJacob Faibussowitsch   PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
2025792fecdfSBarry Smith   PetscCallBLAS("LAPACKgeqrf", LAPACKgeqrf_(&M, &N, A, &lda, tau, work, &ldwork, &info));
20269566063dSJacob Faibussowitsch   PetscCall(PetscFPTrapPop());
202728b400f6SJacob Faibussowitsch   PetscCheck(!info, PETSC_COMM_SELF, PETSC_ERR_LIB, "xGEQRF error");
2028c0ce001eSMatthew G. Knepley   R = A; /* Upper triangular part of A now contains R, the rest contains the elementary reflectors */
2029c0ce001eSMatthew G. Knepley 
2030c0ce001eSMatthew G. Knepley   /* Extract an explicit representation of Q */
2031c0ce001eSMatthew G. Knepley   Q = Ainv;
20329566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(Q, A, mstride * n));
2033c0ce001eSMatthew G. Knepley   K = N; /* full rank */
2034792fecdfSBarry Smith   PetscCallBLAS("LAPACKorgqr", LAPACKorgqr_(&M, &N, &K, Q, &lda, tau, work, &ldwork, &info));
203528b400f6SJacob Faibussowitsch   PetscCheck(!info, PETSC_COMM_SELF, PETSC_ERR_LIB, "xORGQR/xUNGQR error");
2036c0ce001eSMatthew G. Knepley 
2037c0ce001eSMatthew G. Knepley   /* Compute A^{-T} = (R^{-1} Q^T)^T = Q R^{-T} */
2038c0ce001eSMatthew G. Knepley   Alpha = 1.0;
2039c0ce001eSMatthew G. Knepley   ldb   = lda;
2040c0ce001eSMatthew G. Knepley   BLAStrsm_("Right", "Upper", "ConjugateTranspose", "NotUnitTriangular", &M, &N, &Alpha, R, &lda, Q, &ldb);
2041c0ce001eSMatthew G. Knepley   /* Ainv is Q, overwritten with inverse */
2042c0ce001eSMatthew G. Knepley 
2043c0ce001eSMatthew G. Knepley   if (debug) { /* Check that pseudo-inverse worked */
2044c0ce001eSMatthew G. Knepley     PetscScalar  Beta = 0.0;
2045c0ce001eSMatthew G. Knepley     PetscBLASInt ldc;
2046c0ce001eSMatthew G. Knepley     K   = N;
2047c0ce001eSMatthew G. Knepley     ldc = N;
2048c0ce001eSMatthew G. Knepley     BLASgemm_("ConjugateTranspose", "Normal", &N, &K, &M, &Alpha, Ainv, &lda, Aback, &ldb, &Beta, work, &ldc);
20499566063dSJacob Faibussowitsch     PetscCall(PetscScalarView(n * n, work, PETSC_VIEWER_STDOUT_SELF));
20509566063dSJacob Faibussowitsch     PetscCall(PetscFree(Aback));
2051c0ce001eSMatthew G. Knepley   }
20523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2053c0ce001eSMatthew G. Knepley }
2054c0ce001eSMatthew G. Knepley 
2055c0ce001eSMatthew G. Knepley /* Overwrites A. Can handle degenerate problems and m<n. */
PetscFVLeastSquaresPseudoInverseSVD_Static(PetscInt m,PetscInt mstride,PetscInt n,PetscScalar * A,PetscScalar * Ainv,PetscScalar * tau,PetscInt worksize,PetscScalar * work)2056d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscFVLeastSquaresPseudoInverseSVD_Static(PetscInt m, PetscInt mstride, PetscInt n, PetscScalar *A, PetscScalar *Ainv, PetscScalar *tau, PetscInt worksize, PetscScalar *work)
2057d71ae5a4SJacob Faibussowitsch {
20586bb27163SBarry Smith   PetscScalar *Brhs;
2059c0ce001eSMatthew G. Knepley   PetscScalar *tmpwork;
2060c0ce001eSMatthew G. Knepley   PetscReal    rcond;
2061c0ce001eSMatthew G. Knepley #if defined(PETSC_USE_COMPLEX)
2062c0ce001eSMatthew G. Knepley   PetscInt   rworkSize;
2063835f2295SStefano Zampini   PetscReal *rwork, *rtau;
2064c0ce001eSMatthew G. Knepley #endif
2065c0ce001eSMatthew G. Knepley   PetscInt     i, j, maxmn;
2066c0ce001eSMatthew G. Knepley   PetscBLASInt M, N, lda, ldb, ldwork;
2067c0ce001eSMatthew G. Knepley   PetscBLASInt nrhs, irank, info;
2068c0ce001eSMatthew G. Knepley 
2069c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
2070c0ce001eSMatthew G. Knepley   /* initialize to identity */
2071736d4f42SpierreXVI   tmpwork = work;
2072736d4f42SpierreXVI   Brhs    = Ainv;
2073c0ce001eSMatthew G. Knepley   maxmn   = PetscMax(m, n);
2074c0ce001eSMatthew G. Knepley   for (j = 0; j < maxmn; j++) {
2075c0ce001eSMatthew G. Knepley     for (i = 0; i < maxmn; i++) Brhs[i + j * maxmn] = 1.0 * (i == j);
2076c0ce001eSMatthew G. Knepley   }
2077c0ce001eSMatthew G. Knepley 
20789566063dSJacob Faibussowitsch   PetscCall(PetscBLASIntCast(m, &M));
20799566063dSJacob Faibussowitsch   PetscCall(PetscBLASIntCast(n, &N));
20809566063dSJacob Faibussowitsch   PetscCall(PetscBLASIntCast(mstride, &lda));
20819566063dSJacob Faibussowitsch   PetscCall(PetscBLASIntCast(maxmn, &ldb));
20829566063dSJacob Faibussowitsch   PetscCall(PetscBLASIntCast(worksize, &ldwork));
2083c0ce001eSMatthew G. Knepley   rcond = -1;
2084c0ce001eSMatthew G. Knepley   nrhs  = M;
2085c0ce001eSMatthew G. Knepley #if defined(PETSC_USE_COMPLEX)
2086c0ce001eSMatthew G. Knepley   rworkSize = 5 * PetscMin(M, N);
20879566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(rworkSize, &rwork));
2088835f2295SStefano Zampini   PetscCall(PetscMalloc1(PetscMin(M, N), &rtau));
20896bb27163SBarry Smith   PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
2090835f2295SStefano Zampini   PetscCallBLAS("LAPACKgelss", LAPACKgelss_(&M, &N, &nrhs, A, &lda, Brhs, &ldb, rtau, &rcond, &irank, tmpwork, &ldwork, rwork, &info));
20919566063dSJacob Faibussowitsch   PetscCall(PetscFPTrapPop());
20929566063dSJacob Faibussowitsch   PetscCall(PetscFree(rwork));
2093835f2295SStefano Zampini   for (i = 0; i < PetscMin(M, N); i++) tau[i] = rtau[i];
2094835f2295SStefano Zampini   PetscCall(PetscFree(rtau));
2095c0ce001eSMatthew G. Knepley #else
2096c0ce001eSMatthew G. Knepley   nrhs = M;
20976bb27163SBarry Smith   PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
2098835f2295SStefano Zampini   PetscCallBLAS("LAPACKgelss", LAPACKgelss_(&M, &N, &nrhs, A, &lda, Brhs, &ldb, tau, &rcond, &irank, tmpwork, &ldwork, &info));
20999566063dSJacob Faibussowitsch   PetscCall(PetscFPTrapPop());
2100c0ce001eSMatthew G. Knepley #endif
210128b400f6SJacob Faibussowitsch   PetscCheck(!info, PETSC_COMM_SELF, PETSC_ERR_LIB, "xGELSS error");
2102c0ce001eSMatthew G. Knepley   /* The following check should be turned into a diagnostic as soon as someone wants to do this intentionally */
210346091a0eSPierre Jolivet   PetscCheck(irank >= PetscMin(M, N), PETSC_COMM_SELF, PETSC_ERR_USER, "Rank deficient least squares fit, indicates an isolated cell with two collinear points");
21043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2105c0ce001eSMatthew G. Knepley }
2106c0ce001eSMatthew G. Knepley 
2107c0ce001eSMatthew G. Knepley #if 0
2108c0ce001eSMatthew G. Knepley static PetscErrorCode PetscFVLeastSquaresDebugCell_Static(PetscFV fvm, PetscInt cell, DM dm, DM dmFace, PetscScalar *fgeom, DM dmCell, PetscScalar *cgeom)
2109c0ce001eSMatthew G. Knepley {
2110c0ce001eSMatthew G. Knepley   PetscReal       grad[2] = {0, 0};
2111c0ce001eSMatthew G. Knepley   const PetscInt *faces;
2112c0ce001eSMatthew G. Knepley   PetscInt        numFaces, f;
2113c0ce001eSMatthew G. Knepley 
2114c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
21159566063dSJacob Faibussowitsch   PetscCall(DMPlexGetConeSize(dm, cell, &numFaces));
21169566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCone(dm, cell, &faces));
2117c0ce001eSMatthew G. Knepley   for (f = 0; f < numFaces; ++f) {
2118c0ce001eSMatthew G. Knepley     const PetscInt *fcells;
2119c0ce001eSMatthew G. Knepley     const CellGeom *cg1;
2120c0ce001eSMatthew G. Knepley     const FaceGeom *fg;
2121c0ce001eSMatthew G. Knepley 
21229566063dSJacob Faibussowitsch     PetscCall(DMPlexGetSupport(dm, faces[f], &fcells));
21239566063dSJacob Faibussowitsch     PetscCall(DMPlexPointLocalRead(dmFace, faces[f], fgeom, &fg));
2124c0ce001eSMatthew G. Knepley     for (i = 0; i < 2; ++i) {
2125c0ce001eSMatthew G. Knepley       PetscScalar du;
2126c0ce001eSMatthew G. Knepley 
2127c0ce001eSMatthew G. Knepley       if (fcells[i] == c) continue;
21289566063dSJacob Faibussowitsch       PetscCall(DMPlexPointLocalRead(dmCell, fcells[i], cgeom, &cg1));
2129c0ce001eSMatthew G. Knepley       du   = cg1->centroid[0] + 3*cg1->centroid[1] - (cg->centroid[0] + 3*cg->centroid[1]);
2130c0ce001eSMatthew G. Knepley       grad[0] += fg->grad[!i][0] * du;
2131c0ce001eSMatthew G. Knepley       grad[1] += fg->grad[!i][1] * du;
2132c0ce001eSMatthew G. Knepley     }
2133c0ce001eSMatthew G. Knepley   }
21349566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF, "cell[%d] grad (%g, %g)\n", cell, grad[0], grad[1]));
21353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2136c0ce001eSMatthew G. Knepley }
2137c0ce001eSMatthew G. Knepley #endif
2138c0ce001eSMatthew G. Knepley 
2139c0ce001eSMatthew G. Knepley /*
21406f6f020cSMatthew G. Knepley   PetscFVComputeGradient_LeastSquares - Compute the gradient reconstruction matrix for a given cell
2141c0ce001eSMatthew G. Knepley 
2142c0ce001eSMatthew G. Knepley   Input Parameters:
2143dce8aebaSBarry Smith + fvm      - The `PetscFV` object
2144c0ce001eSMatthew G. Knepley . numFaces - The number of cell faces which are not constrained
2145c0ce001eSMatthew G. Knepley . dx       - The vector from the cell centroid to the neighboring cell centroid for each face
2146c0ce001eSMatthew G. Knepley 
2147c0ce001eSMatthew G. Knepley   Level: developer
2148c0ce001eSMatthew G. Knepley 
2149dce8aebaSBarry Smith .seealso: `PetscFV`, `PetscFVCreate()`
2150c0ce001eSMatthew G. Knepley */
PetscFVComputeGradient_LeastSquares(PetscFV fvm,PetscInt numFaces,const PetscScalar dx[],PetscScalar grad[])2151d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscFVComputeGradient_LeastSquares(PetscFV fvm, PetscInt numFaces, const PetscScalar dx[], PetscScalar grad[])
2152d71ae5a4SJacob Faibussowitsch {
2153c0ce001eSMatthew G. Knepley   PetscFV_LeastSquares *ls       = (PetscFV_LeastSquares *)fvm->data;
2154c0ce001eSMatthew G. Knepley   const PetscBool       useSVD   = PETSC_TRUE;
2155c0ce001eSMatthew G. Knepley   const PetscInt        maxFaces = ls->maxFaces;
2156c0ce001eSMatthew G. Knepley   PetscInt              dim, f, d;
2157c0ce001eSMatthew G. Knepley 
2158c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
2159c0ce001eSMatthew G. Knepley   if (numFaces > maxFaces) {
216008401ef6SPierre Jolivet     PetscCheck(maxFaces >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Reconstruction has not been initialized, call PetscFVLeastSquaresSetMaxFaces()");
216163a3b9bcSJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of input faces %" PetscInt_FMT " > %" PetscInt_FMT " maxfaces", numFaces, maxFaces);
2162c0ce001eSMatthew G. Knepley   }
21639566063dSJacob Faibussowitsch   PetscCall(PetscFVGetSpatialDimension(fvm, &dim));
2164c0ce001eSMatthew G. Knepley   for (f = 0; f < numFaces; ++f) {
2165c0ce001eSMatthew G. Knepley     for (d = 0; d < dim; ++d) ls->B[d * maxFaces + f] = dx[f * dim + d];
2166c0ce001eSMatthew G. Knepley   }
2167c0ce001eSMatthew G. Knepley   /* Overwrites B with garbage, returns Binv in row-major format */
2168736d4f42SpierreXVI   if (useSVD) {
2169736d4f42SpierreXVI     PetscInt maxmn = PetscMax(numFaces, dim);
21709566063dSJacob Faibussowitsch     PetscCall(PetscFVLeastSquaresPseudoInverseSVD_Static(numFaces, maxFaces, dim, ls->B, ls->Binv, ls->tau, ls->workSize, ls->work));
2171736d4f42SpierreXVI     /* Binv shaped in column-major, coldim=maxmn.*/
2172736d4f42SpierreXVI     for (f = 0; f < numFaces; ++f) {
2173736d4f42SpierreXVI       for (d = 0; d < dim; ++d) grad[f * dim + d] = ls->Binv[d + maxmn * f];
2174736d4f42SpierreXVI     }
2175736d4f42SpierreXVI   } else {
21769566063dSJacob Faibussowitsch     PetscCall(PetscFVLeastSquaresPseudoInverse_Static(numFaces, maxFaces, dim, ls->B, ls->Binv, ls->tau, ls->workSize, ls->work));
2177736d4f42SpierreXVI     /* Binv shaped in row-major, rowdim=maxFaces.*/
2178c0ce001eSMatthew G. Knepley     for (f = 0; f < numFaces; ++f) {
2179c0ce001eSMatthew G. Knepley       for (d = 0; d < dim; ++d) grad[f * dim + d] = ls->Binv[d * maxFaces + f];
2180c0ce001eSMatthew G. Knepley     }
2181736d4f42SpierreXVI   }
21823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2183c0ce001eSMatthew G. Knepley }
2184c0ce001eSMatthew G. Knepley 
2185c0ce001eSMatthew G. Knepley /*
2186c0ce001eSMatthew G. Knepley   neighborVol[f*2+0] contains the left  geom
2187c0ce001eSMatthew G. Knepley   neighborVol[f*2+1] contains the right geom
2188c0ce001eSMatthew G. Knepley */
PetscFVIntegrateRHSFunction_LeastSquares(PetscFV fvm,PetscDS prob,PetscInt field,PetscInt Nf,PetscFVFaceGeom * fgeom,PetscReal * neighborVol,PetscScalar uL[],PetscScalar uR[],PetscScalar fluxL[],PetscScalar fluxR[])2189d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscFVIntegrateRHSFunction_LeastSquares(PetscFV fvm, PetscDS prob, PetscInt field, PetscInt Nf, PetscFVFaceGeom *fgeom, PetscReal *neighborVol, PetscScalar uL[], PetscScalar uR[], PetscScalar fluxL[], PetscScalar fluxR[])
2190d71ae5a4SJacob Faibussowitsch {
2191c0ce001eSMatthew G. Knepley   void (*riemann)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscInt, const PetscScalar[], PetscScalar[], void *);
2192c0ce001eSMatthew G. Knepley   void              *rctx;
2193c0ce001eSMatthew G. Knepley   PetscScalar       *flux = fvm->fluxWork;
2194c0ce001eSMatthew G. Knepley   const PetscScalar *constants;
2195c0ce001eSMatthew G. Knepley   PetscInt           dim, numConstants, pdim, Nc, totDim, off, f, d;
2196c0ce001eSMatthew G. Knepley 
2197c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
21989566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalComponents(prob, &Nc));
21999566063dSJacob Faibussowitsch   PetscCall(PetscDSGetTotalDimension(prob, &totDim));
22009566063dSJacob Faibussowitsch   PetscCall(PetscDSGetFieldOffset(prob, field, &off));
22019566063dSJacob Faibussowitsch   PetscCall(PetscDSGetRiemannSolver(prob, field, &riemann));
22029566063dSJacob Faibussowitsch   PetscCall(PetscDSGetContext(prob, field, &rctx));
22039566063dSJacob Faibussowitsch   PetscCall(PetscDSGetConstants(prob, &numConstants, &constants));
22049566063dSJacob Faibussowitsch   PetscCall(PetscFVGetSpatialDimension(fvm, &dim));
22059566063dSJacob Faibussowitsch   PetscCall(PetscFVGetNumComponents(fvm, &pdim));
2206c0ce001eSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
2207c0ce001eSMatthew G. Knepley     (*riemann)(dim, pdim, fgeom[f].centroid, fgeom[f].normal, &uL[f * Nc], &uR[f * Nc], numConstants, constants, flux, rctx);
2208c0ce001eSMatthew G. Knepley     for (d = 0; d < pdim; ++d) {
2209c0ce001eSMatthew G. Knepley       fluxL[f * totDim + off + d] = flux[d] / neighborVol[f * 2 + 0];
2210c0ce001eSMatthew G. Knepley       fluxR[f * totDim + off + d] = flux[d] / neighborVol[f * 2 + 1];
2211c0ce001eSMatthew G. Knepley     }
2212c0ce001eSMatthew G. Knepley   }
22133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2214c0ce001eSMatthew G. Knepley }
2215c0ce001eSMatthew G. Knepley 
PetscFVLeastSquaresSetMaxFaces_LS(PetscFV fvm,PetscInt maxFaces)2216d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscFVLeastSquaresSetMaxFaces_LS(PetscFV fvm, PetscInt maxFaces)
2217d71ae5a4SJacob Faibussowitsch {
2218c0ce001eSMatthew G. Knepley   PetscFV_LeastSquares *ls = (PetscFV_LeastSquares *)fvm->data;
2219736d4f42SpierreXVI   PetscInt              dim, m, n, nrhs, minmn, maxmn;
2220c0ce001eSMatthew G. Knepley 
2221c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
2222c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1);
22239566063dSJacob Faibussowitsch   PetscCall(PetscFVGetSpatialDimension(fvm, &dim));
22249566063dSJacob Faibussowitsch   PetscCall(PetscFree4(ls->B, ls->Binv, ls->tau, ls->work));
2225c0ce001eSMatthew G. Knepley   ls->maxFaces = maxFaces;
2226c0ce001eSMatthew G. Knepley   m            = ls->maxFaces;
2227c0ce001eSMatthew G. Knepley   n            = dim;
2228c0ce001eSMatthew G. Knepley   nrhs         = ls->maxFaces;
2229736d4f42SpierreXVI   minmn        = PetscMin(m, n);
2230736d4f42SpierreXVI   maxmn        = PetscMax(m, n);
2231736d4f42SpierreXVI   ls->workSize = 3 * minmn + PetscMax(2 * minmn, PetscMax(maxmn, nrhs)); /* required by LAPACK */
22329566063dSJacob Faibussowitsch   PetscCall(PetscMalloc4(m * n, &ls->B, maxmn * maxmn, &ls->Binv, minmn, &ls->tau, ls->workSize, &ls->work));
22333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2234c0ce001eSMatthew G. Knepley }
2235c0ce001eSMatthew G. Knepley 
PetscFVInitialize_LeastSquares(PetscFV fvm)223666976f2fSJacob Faibussowitsch static PetscErrorCode PetscFVInitialize_LeastSquares(PetscFV fvm)
2237d71ae5a4SJacob Faibussowitsch {
2238c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
2239c0ce001eSMatthew G. Knepley   fvm->ops->setfromoptions       = NULL;
2240c0ce001eSMatthew G. Knepley   fvm->ops->view                 = PetscFVView_LeastSquares;
2241c0ce001eSMatthew G. Knepley   fvm->ops->destroy              = PetscFVDestroy_LeastSquares;
2242c0ce001eSMatthew G. Knepley   fvm->ops->computegradient      = PetscFVComputeGradient_LeastSquares;
2243c0ce001eSMatthew G. Knepley   fvm->ops->integraterhsfunction = PetscFVIntegrateRHSFunction_LeastSquares;
22443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2245c0ce001eSMatthew G. Knepley }
2246c0ce001eSMatthew G. Knepley 
2247c0ce001eSMatthew G. Knepley /*MC
2248dce8aebaSBarry Smith   PETSCFVLEASTSQUARES = "leastsquares" - A `PetscFV` implementation
2249c0ce001eSMatthew G. Knepley 
2250c0ce001eSMatthew G. Knepley   Level: intermediate
2251c0ce001eSMatthew G. Knepley 
2252dce8aebaSBarry Smith .seealso: `PetscFV`, `PetscFVType`, `PetscFVCreate()`, `PetscFVSetType()`
2253c0ce001eSMatthew G. Knepley M*/
2254c0ce001eSMatthew G. Knepley 
PetscFVCreate_LeastSquares(PetscFV fvm)2255d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscFVCreate_LeastSquares(PetscFV fvm)
2256d71ae5a4SJacob Faibussowitsch {
2257c0ce001eSMatthew G. Knepley   PetscFV_LeastSquares *ls;
2258c0ce001eSMatthew G. Knepley 
2259c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
2260c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1);
22614dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&ls));
2262c0ce001eSMatthew G. Knepley   fvm->data = ls;
2263c0ce001eSMatthew G. Knepley 
2264c0ce001eSMatthew G. Knepley   ls->maxFaces = -1;
2265c0ce001eSMatthew G. Knepley   ls->workSize = -1;
2266c0ce001eSMatthew G. Knepley   ls->B        = NULL;
2267c0ce001eSMatthew G. Knepley   ls->Binv     = NULL;
2268c0ce001eSMatthew G. Knepley   ls->tau      = NULL;
2269c0ce001eSMatthew G. Knepley   ls->work     = NULL;
2270c0ce001eSMatthew G. Knepley 
22719566063dSJacob Faibussowitsch   PetscCall(PetscFVSetComputeGradients(fvm, PETSC_TRUE));
22729566063dSJacob Faibussowitsch   PetscCall(PetscFVInitialize_LeastSquares(fvm));
22739566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)fvm, "PetscFVLeastSquaresSetMaxFaces_C", PetscFVLeastSquaresSetMaxFaces_LS));
22743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2275c0ce001eSMatthew G. Knepley }
2276c0ce001eSMatthew G. Knepley 
2277c0ce001eSMatthew G. Knepley /*@
2278c0ce001eSMatthew G. Knepley   PetscFVLeastSquaresSetMaxFaces - Set the maximum number of cell faces for gradient reconstruction
2279c0ce001eSMatthew G. Knepley 
228020f4b53cSBarry Smith   Not Collective
2281c0ce001eSMatthew G. Knepley 
228260225df5SJacob Faibussowitsch   Input Parameters:
2283dce8aebaSBarry Smith + fvm      - The `PetscFV` object
2284c0ce001eSMatthew G. Knepley - maxFaces - The maximum number of cell faces
2285c0ce001eSMatthew G. Knepley 
2286c0ce001eSMatthew G. Knepley   Level: intermediate
2287c0ce001eSMatthew G. Knepley 
2288dce8aebaSBarry Smith .seealso: `PetscFV`, `PetscFVCreate()`, `PETSCFVLEASTSQUARES`, `PetscFVComputeGradient()`
2289c0ce001eSMatthew G. Knepley @*/
PetscFVLeastSquaresSetMaxFaces(PetscFV fvm,PetscInt maxFaces)2290d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFVLeastSquaresSetMaxFaces(PetscFV fvm, PetscInt maxFaces)
2291d71ae5a4SJacob Faibussowitsch {
2292c0ce001eSMatthew G. Knepley   PetscFunctionBegin;
2293c0ce001eSMatthew G. Knepley   PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1);
2294cac4c232SBarry Smith   PetscTryMethod(fvm, "PetscFVLeastSquaresSetMaxFaces_C", (PetscFV, PetscInt), (fvm, maxFaces));
22953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2296c0ce001eSMatthew G. Knepley }
2297