bound.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) 2011-2023 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 "fvcAverage.H"
28 
29 // * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
30 
32 {
33  const scalar minVsf = Foam::min(vsf).value();
34 
35  if (minVsf < min.value())
36  {
37  scalarField& isf = vsf.primitiveFieldRef();
38 
39  Info<< "bounding " << vsf.name()
40  << ", min: " << minVsf
41  << " max: " << max(vsf).value()
42  << " average: " << gAverage(isf)
43  << endl;
44 
45  isf = max
46  (
47  max
48  (
49  isf,
50  fvc::average(max(vsf, min))().primitiveField()
51  *pos0(-isf)
52  ),
53  min.value()
54  );
55 
57  forAll(bsf, patchi)
58  {
59  bsf[patchi] == max
60  (
61  max
62  (
63  bsf[patchi],
64  bsf[patchi].patchInternalField()
65  *pos0(-bsf[patchi])
66  ),
67  min.value()
68  );
69  }
70 
71  return true;
72  }
73  else
74  {
75  return false;
76  }
77 }
78 
79 
80 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Bound the given scalar field where it is below the specified minimum.
Generic GeometricBoundaryField class.
Generic GeometricField class.
Internal::FieldType & primitiveFieldRef()
Return a reference to the primitive field.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
const word & name() const
Return name.
Definition: IOobject.H:310
Area-weighted average a surfaceField creating a volField.
label patchi
tmp< VolField< Type > > average(const SurfaceField< Type > &ssf)
Area-weighted average a surfaceField creating a volField.
Definition: fvcAverage.C:46
dimensionedScalar pos0(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
messageStream Info
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
bool bound(volScalarField &, const dimensionedScalar &min)
Bound the given scalar field where it is below the specified min value.
Definition: bound.C:31
Type gAverage(const FieldField< Field, Type > &f)