bound.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) 2011-2016 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 "bound.H"
27 #include "volFields.H"
28 #include "fvc.H"
29 
30 // * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
31 
34 {
35  const scalar minVsf = min(vsf).value();
36 
37  if (minVsf < lowerBound.value())
38  {
39  Info<< "bounding " << vsf.name()
40  << ", min: " << minVsf
41  << " max: " << max(vsf).value()
42  << " average: " << gAverage(vsf.primitiveField())
43  << endl;
44 
45  vsf.primitiveFieldRef() = max
46  (
47  max
48  (
49  vsf.primitiveField(),
50  fvc::average(max(vsf, lowerBound))().primitiveField()
51  * pos(-vsf.primitiveField())
52  ),
53  lowerBound.value()
54  );
55 
56  vsf.boundaryFieldRef() = max(vsf.boundaryField(), lowerBound.value());
57  }
58 
59  return vsf;
60 }
61 
62 
63 // ************************************************************************* //
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
const Internal::FieldType & primitiveField() const
Return a const-reference to the internal field.
const Type & value() const
Return const reference to value.
dimensionedScalar pos(const dimensionedScalar &ds)
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Internal::FieldType & primitiveFieldRef()
Return a reference to the internal field.
Bound the given scalar field if it has gone unbounded.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:33
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
Type gAverage(const FieldField< Field, Type > &f)
messageStream Info
tmp< GeometricField< Type, fvPatchField, volMesh > > average(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Area-weighted average a surfaceField creating a volField.
Definition: fvcAverage.C:46
const word & name() const
Return name.
Definition: IOobject.H:260