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  Copyright (C) 2023 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 "bound.H"
30 #include "volFields.H"
31 #include "fvc.H"
32 
33 // * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
34 
36 Foam::bound(volScalarField& vsf, const dimensionedScalar& lowerBound)
37 {
38  const scalar minVsf = min(vsf).value();
39 
40  if (minVsf < lowerBound.value())
41  {
42  Info<< "bounding " << vsf.name()
43  << ", min: " << minVsf
44  << " max: " << max(vsf).value()
45  << " average: " << gAverage(vsf.primitiveField())
46  << endl;
47 
48  vsf.primitiveFieldRef() = max
49  (
50  max
51  (
52  vsf.primitiveField(),
53  fvc::average(max(vsf, lowerBound))().primitiveField()
54  * pos0(-vsf.primitiveField())
55  ),
56  lowerBound.value()
57  );
58 
59  vsf.boundaryFieldRef() = max(vsf.boundaryField(), lowerBound.value());
60 
61  // Give coupled bcs chance to update since cell values changed
63  }
64 
65  return vsf;
66 }
67 
68 
69 // ************************************************************************* //
const Type & value() const noexcept
Return const reference to value.
const Internal::FieldType & primitiveField() const noexcept
Return a const-reference to the internal field values.
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:195
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
An abstract base class for patches that couple regions of the computational domain e...
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
tmp< GeometricField< Type, faPatchField, areaMesh > > average(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Area-weighted average a edgeField creating a areaField.
Definition: facAverage.C:40
void evaluateCoupled()
Evaluate boundary conditions on a subset of coupled patches.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
dimensionedScalar pos0(const dimensionedScalar &ds)
Bound the given scalar field if it has gone unbounded.
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:29
Type gAverage(const FieldField< Field, Type > &f)
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
messageStream Info
Information stream (stdout output on master, null elsewhere)
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.