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 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 "gaussConvectionScheme.H"
31 #include "blendedSchemeBase.H"
32 #include "fvcCellReduce.H"
33 
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 
36 template<class Type>
37 void Foam::functionObjects::blendingFactor::calcBlendingFactor
38 (
39  const GeometricField<Type, fvPatchField, volMesh>& field,
40  const typename fv::convectionScheme<Type>& cs
41 )
42 {
43  if (!isA<fv::gaussConvectionScheme<Type>>(cs))
44  {
46  << "Scheme for field " << field.name() << " is not a "
47  << fv::gaussConvectionScheme<Type>::typeName
48  << " scheme. Not calculating " << resultName_ << endl;
49 
50  return;
51  }
52 
53  const fv::gaussConvectionScheme<Type>& gcs =
54  refCast<const fv::gaussConvectionScheme<Type>>(cs);
55 
56 
57  const surfaceInterpolationScheme<Type>& interpScheme = gcs.interpScheme();
58 
59  if (!isA<blendedSchemeBase<Type>>(interpScheme))
60  {
62  << interpScheme.type() << " is not a blended scheme"
63  << ". Not calculating " << resultName_ << endl;
64 
65  return;
66  }
67 
68  // Retrieve the face-based blending factor
69  const blendedSchemeBase<Type>& blendedScheme =
70  refCast<const blendedSchemeBase<Type>>(interpScheme);
71  const surfaceScalarField factorf(blendedScheme.blendingFactor(field));
72 
73  // Convert into vol field whose values represent the local face minima
74  // Note:
75  // - factor applied to 1st scheme, and (1-factor) to 2nd scheme
76  // - not using the store(...) mechanism due to need to correct BCs
77  volScalarField& indicator =
78  lookupObjectRef<volScalarField>(resultName_);
79 
80  indicator = 1 - fvc::cellReduce(factorf, minEqOp<scalar>(), GREAT);
81  indicator.correctBoundaryConditions();
82 }
83 
84 
85 template<class Type>
86 bool Foam::functionObjects::blendingFactor::calcScheme()
87 {
88  typedef GeometricField<Type, fvPatchField, volMesh> FieldType;
89 
90  if (!foundObject<FieldType>(fieldName_, false))
91  {
92  return false;
93  }
94 
95  const FieldType& field = lookupObject<FieldType>(fieldName_);
96 
97  const word divScheme("div(" + phiName_ + ',' + fieldName_ + ')');
98  ITstream& its = mesh_.divScheme(divScheme);
99 
100  const surfaceScalarField& phi = lookupObject<surfaceScalarField>(phiName_);
101 
102  tmp<fv::convectionScheme<Type>> tcs =
104 
105  if (isA<fv::boundedConvectionScheme<Type>>(tcs()))
106  {
107  const fv::boundedConvectionScheme<Type>& bcs =
108  refCast<const fv::boundedConvectionScheme<Type>>(tcs());
109 
110  calcBlendingFactor(field, bcs.scheme());
111  }
112  else
113  {
114  const fv::gaussConvectionScheme<Type>& gcs =
115  refCast<const fv::gaussConvectionScheme<Type>>(tcs());
116 
117  calcBlendingFactor(field, gcs);
118  }
119 
120  return true;
121 }
122 
123 
124 // ************************************************************************* //
blendedSchemeBase.H
Foam::fv::convectionScheme::New
static tmp< convectionScheme< Type > > New(const fvMesh &mesh, const surfaceScalarField &faceFlux, Istream &schemeData)
Return a pointer to a new convectionScheme created on freestore.
Definition: convectionScheme.C:61
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
field
rDeltaTY field()
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::fvc::cellReduce
tmp< GeometricField< Type, fvPatchField, volMesh > > cellReduce(const GeometricField< Type, fvsPatchField, surfaceMesh > &, const CombineOp &cop, const Type &nullValue=pTraits< Type >::zero)
gaussConvectionScheme.H
Foam::isA
const TargetType * isA(const Type &t)
Check if dynamic_cast to TargetType is possible.
Definition: typeInfo.H:197
boundedConvectionScheme.H
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:54
Foam::functionObjects::fieldExpression::resultName_
word resultName_
Name of result field.
Definition: fieldExpression.H:132
fvcCellReduce.H
Construct a volume field from a surface field using a combine operator.
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:328