blendingFactorTemplates.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) 2013-2022 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 "gaussConvectionScheme.H"
27 #include "blendedSchemeBase.H"
28 #include "fvcCellReduce.H"
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 template<class Type>
33 bool Foam::functionObjects::blendingFactor::calcBF()
34 {
35  typedef GeometricField<Type, fvPatchField, volMesh> FieldType;
36 
37  if (!foundObject<FieldType>(fieldName_))
38  {
39  return false;
40  }
41 
42  const FieldType& field = lookupObject<FieldType>(fieldName_);
43 
44  const word divScheme("div(" + phiName_ + ',' + fieldName_ + ')');
45  ITstream& its = mesh_.schemes().div(divScheme);
46 
47  const surfaceScalarField& phi = lookupObject<surfaceScalarField>(phiName_);
48 
49  tmp<fv::convectionScheme<Type>> cs =
51 
52  const fv::gaussConvectionScheme<Type>& gcs =
53  refCast<const fv::gaussConvectionScheme<Type>>(cs());
54 
55  const surfaceInterpolationScheme<Type>& interpScheme =
56  gcs.interpScheme();
57 
58  if (!isA<blendedSchemeBase<Type>>(interpScheme))
59  {
61  << interpScheme.typeName << " is not a blended scheme"
62  << exit(FatalError);
63  }
64 
65  // Retrieve the face-based blending factor
66  const blendedSchemeBase<Type>& blendedScheme =
67  refCast<const blendedSchemeBase<Type>>(interpScheme);
68  tmp<surfaceScalarField> factorf(blendedScheme.blendingFactor(field));
69 
70  // Convert into vol field whose values represent the local face maxima
71  return store
72  (
74  fvc::cellReduce(factorf, maxEqOp<scalar>())
75  );
76 }
77 
78 
79 // ************************************************************************* //
bool isA(const Type &t)
Check if a dynamic_cast to typeid is possible.
Definition: typeInfo.H:134
word resultName_
Name of result field.
const fvSchemes & schemes() const
Return the fvSchemes.
Definition: fvMesh.C:1672
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
ITstream & div(const word &name) const
Definition: fvSchemes.C:411
Construct a volume field from a surface field using a combine operator.
const word fieldName_
Name of field to process.
bool store(const tmp< ObjectType > &tfield)
Store the given field in the objectRegistry.
static tmp< convectionScheme< Type > > New(const fvMesh &mesh, const surfaceScalarField &faceFlux, Istream &schemeData)
Return a pointer to a new convectionScheme created on freestore.
tmp< GeometricField< Type, fvPatchField, volMesh > > cellReduce(const GeometricField< Type, fvsPatchField, surfaceMesh > &, const CombineOp &cop, const Type &nullValue=pTraits< Type >::zero)
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
const fvMesh & mesh_
Reference to the fvMesh.