Files
OpenFOAM-5.x/src/postProcessing/functionObjects/field/fieldMinMax/fieldMinMaxTemplates.C
Henry Weller 75ea76187b GeometricField::GeometricBoundaryField -> GeometricField::Boundary
When the GeometricBoundaryField template class was originally written it
was a separate class in the Foam namespace rather than a sub-class of
GeometricField as it is now.  Without loss of clarity and simplifying
code which access the boundary field of GeometricFields it is better
that GeometricBoundaryField be renamed Boundary for consistency with the
new naming convention for the type of the dimensioned internal field:
Internal, see commit a25a449c9e

This is a very simple text substitution change which can be applied to
any code which compiles with the OpenFOAM-dev libraries.
2016-04-28 07:22:02 +01:00

280 lines
8.5 KiB
C

/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-2016 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
\*---------------------------------------------------------------------------*/
#include "fieldMinMax.H"
#include "volFields.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
template<class Type>
void Foam::fieldMinMax::output
(
const word& fieldName,
const word& outputName,
const vector& minC,
const vector& maxC,
const label minProci,
const label maxProci,
const Type& minValue,
const Type& maxValue
)
{
OFstream& file = this->file();
if (location_)
{
writeTime(file());
writeTabbed(file, fieldName);
file<< token::TAB << minValue
<< token::TAB << minC;
if (Pstream::parRun())
{
file<< token::TAB << minProci;
}
file<< token::TAB << maxValue
<< token::TAB << maxC;
if (Pstream::parRun())
{
file<< token::TAB << maxProci;
}
file<< endl;
if (log_) Info
<< " min(" << outputName << ") = " << minValue
<< " at location " << minC;
if (Pstream::parRun())
{
if (log_) Info<< " on processor " << minProci;
}
if (log_) Info
<< nl << " max(" << outputName << ") = " << maxValue
<< " at location " << maxC;
if (Pstream::parRun())
{
if (log_) Info<< " on processor " << maxProci;
}
}
else
{
file<< token::TAB << minValue << token::TAB << maxValue;
if (log_) Info
<< " min/max(" << outputName << ") = "
<< minValue << ' ' << maxValue;
}
if (log_) Info<< endl;
}
template<class Type>
void Foam::fieldMinMax::calcMinMaxFields
(
const word& fieldName,
const modeType& mode
)
{
typedef GeometricField<Type, fvPatchField, volMesh> fieldType;
if (obr_.foundObject<fieldType>(fieldName))
{
const label proci = Pstream::myProcNo();
const fieldType& field = obr_.lookupObject<fieldType>(fieldName);
const fvMesh& mesh = field.mesh();
const volVectorField::Boundary& CfBoundary =
mesh.C().boundaryField();
switch (mode)
{
case mdMag:
{
const volScalarField magField(mag(field));
const volScalarField::Boundary& magFieldBoundary =
magField.boundaryField();
scalarList minVs(Pstream::nProcs());
List<vector> minCs(Pstream::nProcs());
label minProci = findMin(magField);
minVs[proci] = magField[minProci];
minCs[proci] = field.mesh().C()[minProci];
labelList maxIs(Pstream::nProcs());
scalarList maxVs(Pstream::nProcs());
List<vector> maxCs(Pstream::nProcs());
label maxProci = findMax(magField);
maxVs[proci] = magField[maxProci];
maxCs[proci] = field.mesh().C()[maxProci];
forAll(magFieldBoundary, patchi)
{
const scalarField& mfp = magFieldBoundary[patchi];
if (mfp.size())
{
const vectorField& Cfp = CfBoundary[patchi];
label minPI = findMin(mfp);
if (mfp[minPI] < minVs[proci])
{
minVs[proci] = mfp[minPI];
minCs[proci] = Cfp[minPI];
}
label maxPI = findMax(mfp);
if (mfp[maxPI] > maxVs[proci])
{
maxVs[proci] = mfp[maxPI];
maxCs[proci] = Cfp[maxPI];
}
}
}
Pstream::gatherList(minVs);
Pstream::gatherList(minCs);
Pstream::gatherList(maxVs);
Pstream::gatherList(maxCs);
if (Pstream::master())
{
label minI = findMin(minVs);
scalar minValue = minVs[minI];
const vector& minC = minCs[minI];
label maxI = findMax(maxVs);
scalar maxValue = maxVs[maxI];
const vector& maxC = maxCs[maxI];
output
(
fieldName,
word("mag(" + fieldName + ")"),
minC,
maxC,
minI,
maxI,
minValue,
maxValue
);
}
break;
}
case mdCmpt:
{
const typename fieldType::Boundary&
fieldBoundary = field.boundaryField();
List<Type> minVs(Pstream::nProcs());
List<vector> minCs(Pstream::nProcs());
label minProci = findMin(field);
minVs[proci] = field[minProci];
minCs[proci] = field.mesh().C()[minProci];
Pstream::gatherList(minVs);
Pstream::gatherList(minCs);
List<Type> maxVs(Pstream::nProcs());
List<vector> maxCs(Pstream::nProcs());
label maxProci = findMax(field);
maxVs[proci] = field[maxProci];
maxCs[proci] = field.mesh().C()[maxProci];
forAll(fieldBoundary, patchi)
{
const Field<Type>& fp = fieldBoundary[patchi];
if (fp.size())
{
const vectorField& Cfp = CfBoundary[patchi];
label minPI = findMin(fp);
if (fp[minPI] < minVs[proci])
{
minVs[proci] = fp[minPI];
minCs[proci] = Cfp[minPI];
}
label maxPI = findMax(fp);
if (fp[maxPI] > maxVs[proci])
{
maxVs[proci] = fp[maxPI];
maxCs[proci] = Cfp[maxPI];
}
}
}
Pstream::gatherList(minVs);
Pstream::gatherList(minCs);
Pstream::gatherList(maxVs);
Pstream::gatherList(maxCs);
if (Pstream::master())
{
label minI = findMin(minVs);
Type minValue = minVs[minI];
const vector& minC = minCs[minI];
label maxI = findMax(maxVs);
Type maxValue = maxVs[maxI];
const vector& maxC = maxCs[maxI];
output
(
fieldName,
fieldName,
minC,
maxC,
minI,
maxI,
minValue,
maxValue
);
}
break;
}
default:
{
FatalErrorInFunction
<< "Unknown min/max mode: " << modeTypeNames_[mode_]
<< exit(FatalError);
}
}
}
}
// ************************************************************************* //