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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2017 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "bound.H"
29#include "volFields.H"
30#include "fvc.H"
31
32// * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
33
36{
37 const scalar minVsf = min(vsf).value();
38
39 if (minVsf < lowerBound.value())
40 {
41 Info<< "bounding " << vsf.name()
42 << ", min: " << minVsf
43 << " max: " << max(vsf).value()
44 << " average: " << gAverage(vsf.primitiveField())
45 << endl;
46
48 (
49 max
50 (
51 vsf.primitiveField(),
52 fvc::average(max(vsf, lowerBound))().primitiveField()
53 * pos0(-vsf.primitiveField())
54 ),
55 lowerBound.value()
56 );
57
58 vsf.boundaryFieldRef() = max(vsf.boundaryField(), lowerBound.value());
59 }
60
61 return vsf;
62}
63
64
65// ************************************************************************* //
Y[inertIndex] max(0.0)
Bound the given scalar field if it has gone unbounded.
const Internal::FieldType & primitiveField() const
Return a const-reference to the internal field.
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:65
const Type & value() const
Return const reference to value.
dimensionedScalar pos0(const dimensionedScalar &ds)
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:35
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Type gAverage(const FieldField< Field, Type > &f)