laplacianScheme.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 \*---------------------------------------------------------------------------*/
25 
26 #include "fv.H"
27 #include "HashTable.H"
28 #include "linear.H"
29 #include "fvMatrix.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 namespace fv
39 {
40 
41 // * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
42 
43 template<class Type, class GType>
44 tmp<laplacianScheme<Type, GType>> laplacianScheme<Type, GType>::New
45 (
46  const fvMesh& mesh,
47  Istream& schemeData
48 )
49 {
50  if (fv::debug)
51  {
52  InfoInFunction << "Constructing laplacianScheme<Type, GType>" << endl;
53  }
54 
55  if (schemeData.eof())
56  {
58  (
59  schemeData
60  ) << "Laplacian scheme not specified" << endl << endl
61  << "Valid laplacian schemes are :" << endl
62  << IstreamConstructorTablePtr_->sortedToc()
63  << exit(FatalIOError);
64  }
65 
66  const word schemeName(schemeData);
67 
68  typename IstreamConstructorTable::iterator cstrIter =
69  IstreamConstructorTablePtr_->find(schemeName);
70 
71  if (cstrIter == IstreamConstructorTablePtr_->end())
72  {
74  (
75  schemeData
76  ) << "Unknown laplacian scheme " << schemeName << nl << nl
77  << "Valid laplacian schemes are :" << endl
78  << IstreamConstructorTablePtr_->sortedToc()
79  << exit(FatalIOError);
80  }
81 
82  return cstrIter()(mesh, schemeData);
83 }
84 
85 
86 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
87 
88 template<class Type, class GType>
90 {}
91 
92 
93 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
94 
95 template<class Type, class GType>
98 (
101 )
102 {
103  return fvmLaplacian(tinterpGammaScheme_().interpolate(gamma)(), vf);
104 }
105 
106 
107 template<class Type, class GType>
110 (
113 )
114 {
115  return fvcLaplacian(tinterpGammaScheme_().interpolate(gamma)(), vf);
116 }
117 
118 
119 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
120 
121 } // End namespace fv
122 
123 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
124 
125 } // End namespace Foam
126 
127 // ************************************************************************* //
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition: curveTools.C:75
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
virtual ~laplacianScheme()
Destructor.
Generic GeometricField class.
fvMesh & mesh
A class for handling words, derived from string.
Definition: word.H:59
labelList fv(nPoints)
bool eof() const
Return true if end of input seen.
Definition: IOstream.H:336
static const char nl
Definition: Ostream.H:260
virtual tmp< fvMatrix< Type > > fvmLaplacian(const GeometricField< GType, fvsPatchField, surfaceMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)=0
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:318
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
A class for managing temporary objects.
Definition: PtrList.H:53
static tmp< laplacianScheme< Type, GType > > New(const fvMesh &mesh, Istream &schemeData)
Return a pointer to a new laplacianScheme created on freestore.
virtual tmp< GeometricField< Type, fvPatchField, volMesh > > fvcLaplacian(const GeometricField< Type, fvPatchField, volMesh > &)=0
Namespace for OpenFOAM.
IOerror FatalIOError
#define InfoInFunction
Report an information message using Foam::Info.