convectionScheme.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 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  tmp<convectionScheme<Type>>::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  InfoInFunction << "Constructing convectionScheme<Type>" << endl;
67  }
68 
69  if (schemeData.eof())
70  {
72  (
73  schemeData
74  ) << "Convection scheme not specified" << endl << endl
75  << "Valid convection schemes are :" << endl
76  << IstreamConstructorTablePtr_->sortedToc()
77  << exit(FatalIOError);
78  }
79 
80  const word schemeName(schemeData);
81 
82  typename IstreamConstructorTable::iterator cstrIter =
83  IstreamConstructorTablePtr_->find(schemeName);
84 
85  if (cstrIter == IstreamConstructorTablePtr_->end())
86  {
88  (
89  schemeData
90  ) << "Unknown convection scheme " << schemeName << nl << nl
91  << "Valid convection schemes are :" << endl
92  << IstreamConstructorTablePtr_->sortedToc()
93  << exit(FatalIOError);
94  }
95 
96  return cstrIter()(mesh, faceFlux, schemeData);
97 }
98 
99 
100 template<class Type>
102 (
103  const fvMesh& mesh,
105  fieldTable& fields,
106  const surfaceScalarField& faceFlux,
107  Istream& schemeData
108 )
109 {
110  if (fv::debug)
111  {
112  InfoInFunction << "Constructing convectionScheme<Type>" << endl;
113  }
114 
115  if (schemeData.eof())
116  {
118  (
119  schemeData
120  ) << "Convection scheme not specified" << endl << endl
121  << "Valid convection schemes are :" << endl
122  << MultivariateConstructorTablePtr_->sortedToc()
123  << exit(FatalIOError);
124  }
125 
126  const word schemeName(schemeData);
127 
128  typename MultivariateConstructorTable::iterator cstrIter =
129  MultivariateConstructorTablePtr_->find(schemeName);
130 
131  if (cstrIter == MultivariateConstructorTablePtr_->end())
132  {
134  (
135  schemeData
136  ) << "Unknown convection scheme " << schemeName << nl << nl
137  << "Valid convection schemes are :" << endl
138  << MultivariateConstructorTablePtr_->sortedToc()
139  << exit(FatalIOError);
140  }
141 
142  return cstrIter()(mesh, fields, faceFlux, schemeData);
143 }
144 
145 
146 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
147 
148 template<class Type>
150 {}
151 
152 
153 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
154 
155 template<class Type>
157 {
158  if (this == &cs)
159  {
161  << "attempted assignment to self"
162  << abort(FatalError);
163  }
164 }
165 
166 
167 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
168 
169 } // End namespace fv
170 
171 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
172 
173 } // End namespace Foam
174 
175 // ************************************************************************* //
convectionScheme(const convectionScheme &)
Copy construct.
void operator=(const convectionScheme< Type > &)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Reference counter for various OpenFOAM components.
Definition: refCount.H:49
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
virtual ~convectionScheme()
Destructor.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Abstract base class for convection schemes.
A class for handling words, derived from string.
Definition: word.H:59
labelList fv(nPoints)
errorManip< error > abort(error &err)
Definition: errorManip.H:131
bool eof() const
Return true if end of input seen.
Definition: IOstream.H:336
static const char nl
Definition: Ostream.H:260
Info<< "Reading field p_rgh\"<< endl;volScalarField p_rgh(IOobject("p_rgh", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);pressureReference pressureReference(p, p_rgh, pimple.dict(), thermo.incompressible());mesh.schemes().setFluxRequired(p_rgh.name());hydrostaticInitialisation(p_rgh, p, rho, U, gh, ghf, pRef, thermo, pimple.dict());Info<< "Creating field dpdt\"<< endl;volScalarField dpdt(IOobject("dpdt", runTime.timeName(), mesh), mesh, dimensionedScalar(p.dimensions()/dimTime, 0));Info<< "Creating field kinetic energy K\"<< endl;volScalarField K("K", 0.5 *magSqr(U));dimensionedScalar initialMass=fvc::domainIntegrate(rho);multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:131
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:318
static tmp< convectionScheme< Type > > New(const fvMesh &mesh, const surfaceScalarField &faceFlux, Istream &schemeData)
Return a pointer to a new convectionScheme created on freestore.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
Abstract base class for multi-variate surface interpolation schemes.
const fvMesh & mesh() const
Return mesh reference.
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.
IOerror FatalIOError
#define InfoInFunction
Report an information message using Foam::Info.