Loading...
Searching...
No Matches
limitedSurfaceInterpolationScheme.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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2019-2021 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
30#include "volFields.H"
32#include "coupledFvPatchField.H"
33
34// * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
35
36template<class Type>
39(
40 const fvMesh& mesh,
41 Istream& schemeData
42)
43{
44 if (surfaceInterpolation::debug)
45 {
47 << "Constructing limitedSurfaceInterpolationScheme<Type>" << endl;
48 }
49
50 if (schemeData.eof())
51 {
52 FatalIOErrorInFunction(schemeData)
53 << "Discretisation scheme not specified"
54 << endl << endl
55 << "Valid schemes are :" << endl
56 << MeshConstructorTablePtr_->sortedToc()
58 }
59
60 const word schemeName(schemeData);
61
62 auto* ctorPtr = MeshConstructorTable(schemeName);
63
64 if (!ctorPtr)
65 {
67 (
68 schemeData,
69 "discretisation",
70 schemeName,
71 *MeshConstructorTablePtr_
72 ) << exit(FatalIOError);
73 }
75 return ctorPtr(mesh, schemeData);
76}
77
78
79template<class Type>
82(
83 const fvMesh& mesh,
84 const surfaceScalarField& faceFlux,
85 Istream& schemeData
86)
87{
88 if (surfaceInterpolation::debug)
89 {
91 << "Constructing limitedSurfaceInterpolationScheme<Type>"
92 << endl;
93 }
94
95 if (schemeData.eof())
96 {
97 FatalIOErrorInFunction(schemeData)
98 << "Discretisation scheme not specified"
99 << endl << endl
100 << "Valid schemes are :" << endl
101 << MeshConstructorTablePtr_->sortedToc()
102 << exit(FatalIOError);
103 }
104
105 const word schemeName(schemeData);
106
107 auto* ctorPtr = MeshFluxConstructorTable(schemeName);
108
109 if (!ctorPtr)
110 {
112 (
113 schemeData,
114 "discretisation",
115 schemeName,
116 *MeshFluxConstructorTablePtr_
117 ) << exit(FatalIOError);
118 }
119
120 return ctorPtr(mesh, faceFlux, schemeData);
121}
122
123
124// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
125
126template<class Type>
130
131
132// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
133
134template<class Type>
137(
139 const surfaceScalarField& CDweights,
141) const
142{
143 // Note that here the weights field is initialised as the limiter
144 // from which the weight is calculated using the limiter value
145 surfaceScalarField& Weights = tLimiter.ref();
146
147 scalarField& pWeights = Weights.primitiveFieldRef();
148
149 forAll(pWeights, face)
150 {
151 pWeights[face] =
152 pWeights[face]*CDweights[face]
153 + (1.0 - pWeights[face])*pos0(faceFlux_[face]);
154 }
155
156 surfaceScalarField::Boundary& bWeights =
157 Weights.boundaryFieldRef();
158
159 forAll(bWeights, patchi)
160 {
161 scalarField& pWeights = bWeights[patchi];
162
163 const scalarField& pCDweights = CDweights.boundaryField()[patchi];
164 const scalarField& pFaceFlux = faceFlux_.boundaryField()[patchi];
165
166 forAll(pWeights, face)
167 {
168 pWeights[face] =
169 pWeights[face]*pCDweights[face]
170 + (1.0 - pWeights[face])*pos0(pFaceFlux[face]);
171 }
173
174 return tLimiter;
175}
176
177template<class Type>
180(
181 const GeometricField<Type, fvPatchField, volMesh>& phi
182) const
183{
184 return this->weights
185 (
186 phi,
195(
197) const
198{
199 return faceFlux_*this->interpolate(phi);
200}
201
202
203// ************************************************************************* //
Generic GeometricField class.
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
GeometricBoundaryField< scalar, fvsPatchField, surfaceMesh > Boundary
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
bool eof() const noexcept
True if end of input seen.
Definition IOstream.H:289
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition Istream.H:60
A face is a list of labels corresponding to mesh vertices.
Definition face.H:71
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
static tmp< limitedSurfaceInterpolationScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return new tmp interpolation scheme.
virtual tmp< surfaceScalarField > limiter(const GeometricField< Type, fvPatchField, volMesh > &) const =0
Return the interpolation weighting factors.
tmp< surfaceScalarField > weights(const GeometricField< Type, fvPatchField, volMesh > &, const surfaceScalarField &CDweights, tmp< surfaceScalarField > tLimiter) const
Return the interpolation weighting factors for the given field,.
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > flux(const GeometricField< Type, fvPatchField, volMesh > &) const
Return the interpolation weighting factors.
const fvMesh & mesh() const
Return mesh reference.
virtual const surfaceScalarField & weights() const
Return reference to linear difference weighting factors.
A class for managing temporary objects.
Definition tmp.H:75
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
A class for handling words, derived from Foam::string.
Definition word.H:66
mesh interpolate(rAU)
dynamicFvMesh & mesh
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition error.H:637
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition error.H:629
#define InfoInFunction
Report an information message using Foam::Info.
dimensionedScalar pos0(const dimensionedScalar &ds)
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Foam::surfaceFields.