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 | Copyright (C) 2013-2015 OpenFOAM Foundation
6  \\/ M anipulation | Copyright (C) 2015 OpenCFD Ltd.
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>
34 {
36 
37  if (!obr_.foundObject<fieldType>(fieldName_))
38  {
39  return;
40  }
41 
42  const fvMesh& mesh = refCast<const fvMesh>(obr_);
43 
44  const fieldType& field = mesh.lookupObject<fieldType>(fieldName_);
45 
46  const word divScheme("div(" + phiName_ + ',' + fieldName_ + ')');
47  ITstream& its = mesh.divScheme(divScheme);
48 
49  const surfaceScalarField& phi =
50  mesh.lookupObject<surfaceScalarField>(phiName_);
51 
54 
56  refCast<const fv::gaussConvectionScheme<Type> >(cs());
57 
58  const surfaceInterpolationScheme<Type>& interpScheme =
59  gcs.interpScheme();
60 
61  if (!isA<blendedSchemeBase<Type> >(interpScheme))
62  {
64  << interpScheme.typeName << " is not a blended scheme"
65  << exit(FatalError);
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: factor applied to 1st scheme, and (1-factor) to 2nd scheme
75  volScalarField& indicator =
76  const_cast<volScalarField&>
77  (
79  );
80  indicator = 1 - fvc::cellReduce(factorf, minEqOp<scalar>(), GREAT);
81  indicator.correctBoundaryConditions();
82 
83  // Generate scheme statistics
84  label nCellsScheme1 = 0;
85  label nCellsScheme2 = 0;
86  label nCellsBlended = 0;
87  forAll(indicator, cellI)
88  {
89  scalar i = indicator[cellI];
90 
91  if (i < tolerance_)
92  {
93  nCellsScheme1++;
94  }
95  else if (i > (1 - tolerance_))
96  {
97  nCellsScheme2++;
98  }
99  else
100  {
101  nCellsBlended++;
102  }
103  }
104 
105  reduce(nCellsScheme1, sumOp<label>());
106  reduce(nCellsScheme2, sumOp<label>());
107  reduce(nCellsBlended, sumOp<label>());
108 
109  if (log_) Info
110  << type() << " " << name_ << " output:" << nl
111  << " scheme 1 cells : " << nCellsScheme1 << nl
112  << " scheme 2 cells : " << nCellsScheme2 << nl
113  << " blended cells : " << nCellsBlended << nl
114  << endl;
115 
116  file()
117  << obr_.time().time().value()
118  << token::TAB << nCellsScheme1
119  << token::TAB << nCellsScheme2
120  << token::TAB << nCellsBlended
121  << endl;
122 }
123 
124 
125 // ************************************************************************* //
Foam::token::TAB
@ TAB
Definition: token.H:96
Foam::fvc::cellReduce
tmp< GeometricField< Type, fvPatchField, volMesh > > cellReduce(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf, const CombineOp &cop, const Type &nullValue)
Definition: fvcCellReduce.C:46
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
blendedSchemeBase.H
Foam::fv::gaussConvectionScheme::interpScheme
const surfaceInterpolationScheme< Type > & interpScheme() const
Definition: gaussConvectionScheme.C:44
Foam::blendingFactor::log_
Switch log_
Switch to send output to Info as well as to file.
Definition: blendingFactor.H:161
Foam::blendingFactor::tolerance_
scalar tolerance_
Tolerance used when calculating the number of blended cells.
Definition: blendingFactor.H:158
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::objectRegistry::time
const Time & time() const
Return time.
Definition: objectRegistry.H:117
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:58
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:261
Foam::fv::gaussConvectionScheme
Basic second-order convection using face-gradients and Gauss' theorem.
Definition: gaussConvectionScheme.H:58
Foam::blendingFactor::resultName_
word resultName_
Result field name.
Definition: blendingFactor.H:155
Foam::functionObjectFile::file
OFstream & file()
Return access to the file (if only 1)
Definition: functionObjectFile.C:189
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:43
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::ITstream
Input token stream.
Definition: ITstream.H:49
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::minEqOp
Definition: ops.H:78
Foam::blendingFactor::name_
const word name_
Name.
Definition: blendingFactor.H:143
Foam::FatalError
error FatalError
Foam::blendingFactor::obr_
const objectRegistry & obr_
Reference to the database.
Definition: blendingFactor.H:146
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
gaussConvectionScheme.H
Foam::blendingFactor::fieldName_
word fieldName_
Field name.
Definition: blendingFactor.H:152
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::GeometricField::correctBoundaryConditions
void correctBoundaryConditions()
Correct boundary field.
Definition: GeometricField.C:885
Foam::objectRegistry::foundObject
bool foundObject(const word &name) const
Is the named Type found?
Definition: objectRegistryTemplates.C:142
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::sumOp
Definition: ops.H:162
Foam::isA
bool isA(const Type &t)
Check if a dynamic_cast to typeid is possible.
Definition: typeInfo.H:134
Foam::blendedSchemeBase::blendingFactor
virtual tmp< surfaceScalarField > blendingFactor(const GeometricField< Type, fvPatchField, volMesh > &vf) const =0
Return the face-based blending factor.
Foam::surfaceInterpolationScheme
Abstract base class for surface interpolation schemes.
Definition: surfaceInterpolationScheme.H:55
Foam::blendingFactor::phiName_
word phiName_
Name of flux field, default is "phi".
Definition: blendingFactor.H:149
Foam::objectRegistry::lookupObject
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
Definition: objectRegistryTemplates.C:165
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
Foam::blendingFactor::calc
void calc()
Calculate the blending factor.
Definition: blendingFactorTemplates.C:33
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
fvcCellReduce.H
Construct a volume field from a surface field using a combine operator.
Foam::blendedSchemeBase
Base class for blended schemes to provide access to the blending factor surface field.
Definition: blendedSchemeBase.H:52