blendingFactorTemplates.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) 2013-2016 OpenFOAM Foundation
9 Copyright (C) 2016-2022 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
31#include "blendedSchemeBase.H"
32#include "fvcCellReduce.H"
33#include "DEShybrid.H"
34
35// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36
37template<class Type>
38void Foam::functionObjects::blendingFactor::calcBlendingFactor
39(
40 const GeometricField<Type, fvPatchField, volMesh>& field,
41 const typename fv::convectionScheme<Type>& cs
42)
43{
44 if (!isA<fv::gaussConvectionScheme<Type>>(cs))
45 {
47 << "Scheme for field " << field.name() << " is not a "
49 << " scheme. Not calculating " << resultName_ << endl;
50
51 return;
52 }
53
54 const fv::gaussConvectionScheme<Type>& gcs =
55 refCast<const fv::gaussConvectionScheme<Type>>(cs);
56
57
58 const surfaceInterpolationScheme<Type>& interpScheme = gcs.interpScheme();
59
60 if (!isA<blendedSchemeBase<Type>>(interpScheme))
61 {
63 << interpScheme.type() << " is not a blended scheme"
64 << ". Not calculating " << resultName_ << endl;
65
66 return;
67 }
68
69 // Retrieve the face-based blending factor
70 const blendedSchemeBase<Type>& blendedScheme =
71 refCast<const blendedSchemeBase<Type>>(interpScheme);
72 const surfaceScalarField factorf(blendedScheme.blendingFactor(field));
73
74 // Convert into vol field whose values represent the local face minima
75 // Note:
76 // - not using the store(...) mechanism due to need to correct BCs
77 auto& indicator = lookupObjectRef<volScalarField>(resultName_);
78
79 if (isA<DEShybrid<Type>>(blendedScheme))
80 {
81 indicator = fvc::cellReduce(factorf, minEqOp<scalar>(), GREAT);
82 }
83 else
84 {
85 indicator = 1 - fvc::cellReduce(factorf, minEqOp<scalar>(), GREAT);
86 }
87
88 indicator.correctBoundaryConditions();
89}
90
91
92template<class Type>
93bool Foam::functionObjects::blendingFactor::calcScheme()
94{
95 typedef GeometricField<Type, fvPatchField, volMesh> FieldType;
96
97 if (!foundObject<FieldType>(fieldName_, false))
98 {
99 return false;
100 }
101
102 const FieldType& field = lookupObject<FieldType>(fieldName_);
103
104 const word divScheme("div(" + phiName_ + ',' + fieldName_ + ')');
105 ITstream& its = mesh_.divScheme(divScheme);
106
107 const surfaceScalarField& phi = lookupObject<surfaceScalarField>(phiName_);
108
109 tmp<fv::convectionScheme<Type>> tcs =
111
112 if (isA<fv::boundedConvectionScheme<Type>>(tcs()))
113 {
114 const fv::boundedConvectionScheme<Type>& bcs =
115 refCast<const fv::boundedConvectionScheme<Type>>(tcs());
116
117 calcBlendingFactor(field, bcs.scheme());
118 }
119 else
120 {
121 const fv::gaussConvectionScheme<Type>& gcs =
122 refCast<const fv::gaussConvectionScheme<Type>>(tcs());
123
124 calcBlendingFactor(field, gcs);
125 }
126
127 return true;
128}
129
130
131// ************************************************************************* //
surfaceScalarField & phi
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
word resultName_
Name of result field.
rDeltaTY field()
Construct a volume field from a surface field using a combine operator.
#define WarningInFunction
Report a warning using Foam::Warning.
tmp< GeometricField< Type, fvPatchField, volMesh > > cellReduce(const GeometricField< Type, fvsPatchField, surfaceMesh > &, const CombineOp &cop, const Type &nullValue=pTraits< Type >::zero)
const TargetType * isA(const Type &t)
Check if dynamic_cast to TargetType is possible.
Definition: typeInfo.H:197
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
static const char *const typeName
The type name used in ensight case files.