Loading...
Searching...
No Matches
fusedGaussConvectionScheme.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-2016 OpenFOAM Foundation
9 Copyright (C) 2024 M. Janssens
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
29#include "fvcSurfaceOps.H"
31#include "fvcSurfaceIntegrate.H"
32#include "fvMatrices.H"
33
34// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35
36namespace Foam
37{
38
39// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
40
41namespace fv
42{
43
44// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45
46template<class Type>
50 return tinterpScheme_();
51}
52
53
54template<class Type>
57(
58 const surfaceScalarField&,
60) const
62 return tinterpScheme_().interpolate(vf);
63}
64
65
66template<class Type>
69(
70 const surfaceScalarField& faceFlux,
72) const
74 return faceFlux*interpolate(faceFlux, vf);
75}
76
77
78template<class Type>
81(
82 const surfaceScalarField& faceFlux,
84) const
85{
86 DebugPout<< "fusedGaussConvectionScheme<Type>::fvmDiv on " << vf.name()
87 << " with flux " << faceFlux.name() << endl;
88
89 tmp<surfaceScalarField> tweights = tinterpScheme_().weights(vf);
90 const surfaceScalarField& weights = tweights();
91
93 (
95 (
96 vf,
97 faceFlux.dimensions()*vf.dimensions()
98 )
99 );
100 fvMatrix<Type>& fvm = tfvm.ref();
101
102 //fvm.lower() = -weights.primitiveField()*faceFlux.primitiveField();
104 (
105 fvm.lower(),
106 weights.primitiveField(),
107 faceFlux.primitiveField()
108 );
109
110 //fvm.upper() = fvm.lower() + faceFlux.primitiveField();
111 add(fvm.upper(), fvm.lower(), faceFlux.primitiveField());
112
113 fvm.negSumDiag();
114
115 forAll(vf.boundaryField(), patchi)
116 {
117 const fvPatchField<Type>& psf = vf.boundaryField()[patchi];
118 const fvsPatchScalarField& patchFlux = faceFlux.boundaryField()[patchi];
119 const fvsPatchScalarField& pw = weights.boundaryField()[patchi];
120
121 auto& intCoeffs = fvm.internalCoeffs()[patchi];
122 auto& bouCoeffs = fvm.boundaryCoeffs()[patchi];
123
124 //fvm.internalCoeffs()[patchi] = patchFlux*psf.valueInternalCoeffs(pw);
125 multiply(intCoeffs, patchFlux, psf.valueInternalCoeffs(pw)());
126
127 //fvm.boundaryCoeffs()[patchi] = -patchFlux*psf.valueBoundaryCoeffs(pw);
128 multiply(bouCoeffs, patchFlux, psf.valueBoundaryCoeffs(pw)());
129 bouCoeffs.negate();
130 }
131
132 if (tinterpScheme_().corrected())
133 {
134 fvm += fvc::surfaceIntegrate(faceFlux*tinterpScheme_().correction(vf));
135 }
137 return tfvm;
138}
139
140
141template<class Type>
144(
145 const surfaceScalarField& faceFlux,
147) const
148{
149 DebugPout<< "fusedGaussConvectionScheme<Type>::fvcDiv on " << vf.name()
150 << " with flux " << faceFlux.name() << endl;
151
153
154 const fvMesh& mesh = vf.mesh();
155
156 tmp<FieldType> tConvection
157 (
158 new FieldType
159 (
161 (
162 "convection(" + faceFlux.name() + ',' + vf.name() + ')',
163 vf.instance(),
164 mesh,
167 ),
168 mesh,
170 (
171 faceFlux.dimensions()
172 *vf.dimensions()
173 /dimVol,
174 Zero
175 ),
177 )
178 );
179
180 if (this->tinterpScheme_().corrected())
181 {
182 const auto tfaceCorr(this->tinterpScheme_().correction(vf));
183 auto& faceCorr = tfaceCorr();
184
185 const auto interpolator = [&]
186 (
187 const vector& area,
188 const scalar lambda,
189
190 const Type& ownVal,
191 const Type& neiVal,
192
193 const scalar& faceVal,
194 const Type& correction
195
196 ) -> Type
197 {
198 return faceVal*((lambda*(ownVal - neiVal) + neiVal) + correction);
199 };
200
202 (
203 // interpolation factors for volume field
204 this->tinterpScheme_().weights(vf),
205
206 // volume field(s)
207 vf,
208
209 // surface field(s)
210 faceFlux,
211 faceCorr,
212
213 // operation
214 interpolator,
215
216 tConvection.ref(),
217 false
218 );
219 }
220 else
221 {
222 const auto interpolator = [&]
223 (
224 const vector& area,
225 const scalar lambda,
226 const Type& ownVal,
227 const Type& neiVal,
228
229 const scalar& faceVal
230 ) -> Type
231 {
232 return faceVal*(lambda*(ownVal - neiVal) + neiVal);
233 };
234
236 (
237 tinterpScheme_().weights(vf),
238
239 vf,
240
241 faceFlux,
242
243 interpolator,
244 tConvection.ref(),
245 false
246 );
247 }
248
249 tConvection.ref().primitiveFieldRef() /= mesh.Vsc();
250
251 tConvection.ref().correctBoundaryConditions();
252
253 return tConvection;
254}
255
256
257// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
258
259} // End namespace fv
260
261// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
262
263} // End namespace Foam
264
265// ************************************************************************* //
const Mesh & mesh() const noexcept
Return const reference to mesh.
const dimensionSet & dimensions() const noexcept
Return dimensions.
Generic GeometricField class.
const Internal::FieldType & primitiveField() const noexcept
Return a const-reference to the internal field values.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
const fileName & instance() const noexcept
Read access to instance path component.
Definition IOobjectI.H:289
Generic dimensioned Type class.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition fvMatrix.H:118
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
static const word & extrapolatedCalculatedType() noexcept
The type name for extrapolatedCalculated patch fields combines zero-gradient and calculated.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
virtual tmp< Field< Type > > valueBoundaryCoeffs(const tmp< scalarField > &) const
Return the matrix source coefficients corresponding to the evaluation of the value of this patchField...
virtual tmp< Field< Type > > valueInternalCoeffs(const tmp< scalarField > &) const
Return the matrix diagonal coefficients corresponding to the evaluation of the value of this patchFie...
const fvMesh & mesh() const
Return mesh reference.
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > flux(const surfaceScalarField &, const GeometricField< Type, fvPatchField, volMesh > &) const
const surfaceInterpolationScheme< Type > & interpScheme() const
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const surfaceScalarField &, const GeometricField< Type, fvPatchField, volMesh > &) const
tmp< fvMatrix< Type > > fvmDiv(const surfaceScalarField &, const GeometricField< Type, fvPatchField, volMesh > &) const
tmp< GeometricField< Type, fvPatchField, volMesh > > fvcDiv(const surfaceScalarField &, const GeometricField< Type, fvPatchField, volMesh > &) const
Abstract base class for surface interpolation schemes.
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
dynamicFvMesh & mesh
A special matrix type and solver, designed for finite volume solutions of scalar equations.
Surface integrate surfaceField creating a volField. Surface sum a surfaceField creating a volField.
Surface integrate surfaceField creating a volField. Surface sum a surfaceField creating a volField.
#define DebugPout
Report an information message using Foam::Pout.
Namespace for finite-volume.
tmp< GeometricField< Type, fvPatchField, volMesh > > surfaceSum(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
void surfaceIntegrate(Field< Type > &ivf, const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Namespace of functions to calculate implicit derivatives returning a matrix.
Namespace for OpenFOAM.
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
void add(DimensionedField< scalar, GeoMesh > &result, const dimensioned< scalar > &dt1, const DimensionedField< scalar, GeoMesh > &f2)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
void multiplySubtract(Field< typename outerProduct< Type1, Type2 >::type > &result, const UList< Type1 > &f1, const UList< Type2 > &f2)
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition curveTools.C:75
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix by subtracting the matrix multiplied by the current fi...
Vector< scalar > vector
Definition vector.H:57
void multiply(DimensionedField< Type, GeoMesh > &result, const DimensionedField< Type, GeoMesh > &f1, const DimensionedField< scalar, GeoMesh > &f2)
fvsPatchField< scalar > fvsPatchScalarField
const dimensionSet dimVol(dimVolume)
Older spelling for dimVolume.
dimensionedScalar lambda("lambda", dimTime/sqr(dimLength), laminarTransport)
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299