convectionScheme.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Description
25  Abstract base class for finite volume calculus convection schemes.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "fv.H"
30 #include "HashTable.H"
31 #include "linear.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 
38 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39 
40 namespace fv
41 {
42 
43 
44 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45 
46 template<class Type>
48 :
49  refCount(),
50  mesh_(cs.mesh_)
51 {}
52 
53 
54 // * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
55 
56 template<class Type>
58 (
59  const fvMesh& mesh,
60  const surfaceScalarField& faceFlux,
61  Istream& schemeData
62 )
63 {
64  if (fv::debug)
65  {
66  Info<< "convectionScheme<Type>::New"
67  "(const fvMesh&, const surfaceScalarField&, Istream&) : "
68  "constructing convectionScheme<Type>"
69  << endl;
70  }
71 
72  if (schemeData.eof())
73  {
75  (
76  "convectionScheme<Type>::New"
77  "(const fvMesh&, const surfaceScalarField&, Istream&)",
78  schemeData
79  ) << "Convection scheme not specified" << endl << endl
80  << "Valid convection schemes are :" << endl
81  << IstreamConstructorTablePtr_->sortedToc()
82  << exit(FatalIOError);
83  }
84 
85  const word schemeName(schemeData);
86 
87  typename IstreamConstructorTable::iterator cstrIter =
88  IstreamConstructorTablePtr_->find(schemeName);
89 
90  if (cstrIter == IstreamConstructorTablePtr_->end())
91  {
93  (
94  "convectionScheme<Type>::New"
95  "(const fvMesh&, const surfaceScalarField&, Istream&)",
96  schemeData
97  ) << "Unknown convection scheme " << schemeName << nl << nl
98  << "Valid convection schemes are :" << endl
99  << IstreamConstructorTablePtr_->sortedToc()
100  << exit(FatalIOError);
101  }
102 
103  return cstrIter()(mesh, faceFlux, schemeData);
104 }
105 
106 
107 template<class Type>
109 (
110  const fvMesh& mesh,
112  fieldTable& fields,
113  const surfaceScalarField& faceFlux,
114  Istream& schemeData
115 )
116 {
117  if (fv::debug)
118  {
119  Info<< "convectionScheme<Type>::New"
120  "(const fvMesh&, "
121  "const typename multivariateSurfaceInterpolationScheme<Type>"
122  "::fieldTable&, const surfaceScalarField&, Istream&) : "
123  "constructing convectionScheme<Type>"
124  << endl;
125  }
126 
127  if (schemeData.eof())
128  {
130  (
131  "convectionScheme<Type>::New"
132  "(const fvMesh&, "
133  "const typename multivariateSurfaceInterpolationScheme<Type>"
134  "::fieldTable&, const surfaceScalarField&, Istream&)",
135  schemeData
136  ) << "Convection scheme not specified" << endl << endl
137  << "Valid convection schemes are :" << endl
138  << MultivariateConstructorTablePtr_->sortedToc()
139  << exit(FatalIOError);
140  }
141 
142  const word schemeName(schemeData);
143 
144  typename MultivariateConstructorTable::iterator cstrIter =
145  MultivariateConstructorTablePtr_->find(schemeName);
146 
147  if (cstrIter == MultivariateConstructorTablePtr_->end())
148  {
150  (
151  "convectionScheme<Type>::New"
152  "(const fvMesh&, "
153  "const typename multivariateSurfaceInterpolationScheme<Type>"
154  "::fieldTable&, const surfaceScalarField&, Istream&)",
155  schemeData
156  ) << "Unknown convection scheme " << schemeName << nl << nl
157  << "Valid convection schemes are :" << endl
158  << MultivariateConstructorTablePtr_->sortedToc()
159  << exit(FatalIOError);
160  }
161 
162  return cstrIter()(mesh, fields, faceFlux, schemeData);
163 }
164 
165 
166 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
167 
168 template<class Type>
170 {}
171 
172 
173 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
174 
175 template<class Type>
177 {
178  if (this == &cs)
179  {
181  (
182  "convectionScheme<Type>::operator=(const convectionScheme<Type>&)"
183  ) << "attempted assignment to self"
184  << abort(FatalError);
185  }
186 }
187 
188 
189 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
190 
191 } // End namespace fv
192 
193 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
194 
195 } // End namespace Foam
196 
197 // ************************************************************************* //
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Info<< "Creating field dpdt\n"<< endl;volScalarField dpdt(IOobject("dpdt", runTime.timeName(), mesh), mesh, dimensionedScalar("dpdt", p.dimensions()/dimTime, 0));Info<< "Creating field kinetic energy K\n"<< endl;volScalarField K("K", 0.5 *magSqr(U));volScalarField p_rgh(IOobject("p_rgh", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);p_rgh=p-rho *gh;mesh.setFluxRequired(p_rgh.name());multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:127
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
A class for handling words, derived from string.
Definition: word.H:59
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
const fvMesh & mesh() const
Return mesh reference.
messageStream Info
Namespace for OpenFOAM.
bool eof() const
Return true if end of input seen.
Definition: IOstream.H:339
convectionScheme(const convectionScheme &)
Copy construct.
static const char nl
Definition: Ostream.H:260
static tmp< convectionScheme< Type > > New(const fvMesh &mesh, const surfaceScalarField &faceFlux, Istream &schemeData)
Return a pointer to a new convectionScheme created on freestore.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
IOerror FatalIOError
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Reference counter for various OpenFOAM components.
Definition: refCount.H:45
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:314
Abstract base class for convection schemes.
error FatalError
labelList fv(nPoints)
Abstract base class for multi-variate surface interpolation schemes.
virtual ~convectionScheme()
Destructor.
void operator=(const convectionScheme< Type > &)
#define FatalIOErrorIn(functionName, ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:325
A class for managing temporary objects.
Definition: PtrList.H:118