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