initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / finiteVolume / cfdTools / general / bound / bound.C
blobb6290337546252e9c547907cde00c7989a5b5e05
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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 the
13     Free Software Foundation; either version 2 of the License, or (at your
14     option) any later version.
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.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM; if not, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 \*---------------------------------------------------------------------------*/
27 #include "bound.H"
28 #include "volFields.H"
29 #include "fvc.H"
31 // * * * * * * * * * * * * * * * Global Functions  * * * * * * * * * * * * * //
33 void Foam::bound(volScalarField& vsf, const dimensionedScalar& vsf0)
35     scalar minVsf = min(vsf).value();
37     if (minVsf < vsf0.value())
38     {
39         Info<< "bounding " << vsf.name()
40             << ", min: " << minVsf
41             << " max: " << max(vsf).value()
42             << " average: " << gAverage(vsf.internalField())
43             << endl;
45         vsf.internalField() = max
46         (
47             max
48             (
49                 vsf.internalField(),
50                 fvc::average(max(vsf, vsf0))().internalField()
51                 *pos(-vsf.internalField())
52             ),
53             vsf0.value()
54         );
56         vsf.boundaryField() = max(vsf.boundaryField(), vsf0.value());
57     }
61 // ************************************************************************* //