mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
ENH: Added position output to fieldMinMax function object
This commit is contained in:
@ -156,7 +156,9 @@ void Foam::fieldMinMax::writeFileHeader()
|
||||
if (fieldMinMaxFilePtr_.valid())
|
||||
{
|
||||
fieldMinMaxFilePtr_()
|
||||
<< "# Time" << tab << "field" << tab << "min" << tab << "max"
|
||||
<< "# Time" << token::TAB << "field" << token::TAB
|
||||
<< "min" << token::TAB << "position(min)" << token::TAB
|
||||
<< "max" << token::TAB << "position(max)" << token::TAB
|
||||
<< endl;
|
||||
}
|
||||
}
|
||||
@ -184,50 +186,27 @@ void Foam::fieldMinMax::write()
|
||||
makeFile();
|
||||
}
|
||||
|
||||
if (log_)
|
||||
{
|
||||
Info<< type() << " output:" << nl;
|
||||
}
|
||||
|
||||
forAll(fieldSet_, fieldI)
|
||||
{
|
||||
calcMinMaxFields<scalar>(fieldSet_[fieldI]);
|
||||
calcMinMaxFields<vector>(fieldSet_[fieldI]);
|
||||
calcMinMaxFields<sphericalTensor>(fieldSet_[fieldI]);
|
||||
calcMinMaxFields<symmTensor>(fieldSet_[fieldI]);
|
||||
calcMinMaxFields<tensor>(fieldSet_[fieldI]);
|
||||
calcMinMaxFields<scalar>(fieldSet_[fieldI], mdCmpt);
|
||||
calcMinMaxFields<vector>(fieldSet_[fieldI], mode_);
|
||||
calcMinMaxFields<sphericalTensor>(fieldSet_[fieldI], mode_);
|
||||
calcMinMaxFields<symmTensor>(fieldSet_[fieldI], mode_);
|
||||
calcMinMaxFields<tensor>(fieldSet_[fieldI], mode_);
|
||||
}
|
||||
}
|
||||
}
|
||||
|
||||
|
||||
template<>
|
||||
void Foam::fieldMinMax::calcMinMaxFields<Foam::scalar>
|
||||
(
|
||||
const word& fieldName
|
||||
)
|
||||
{
|
||||
if (obr_.foundObject<volScalarField>(fieldName))
|
||||
{
|
||||
const volScalarField& field =
|
||||
obr_.lookupObject<volScalarField>(fieldName);
|
||||
const scalar minValue = min(field).value();
|
||||
const scalar maxValue = max(field).value();
|
||||
|
||||
if (Pstream::master())
|
||||
if (log_)
|
||||
{
|
||||
if (write_)
|
||||
{
|
||||
fieldMinMaxFilePtr_()
|
||||
<< obr_.time().value() << tab
|
||||
<< fieldName << tab << minValue << tab << maxValue << endl;
|
||||
}
|
||||
|
||||
if (log_)
|
||||
{
|
||||
Info<< "fieldMinMax output:" << nl
|
||||
<< " min(" << fieldName << ") = " << minValue << nl
|
||||
<< " max(" << fieldName << ") = " << maxValue << nl
|
||||
<< endl;
|
||||
}
|
||||
Info<< endl;
|
||||
}
|
||||
}
|
||||
}
|
||||
|
||||
|
||||
|
||||
// ************************************************************************* //
|
||||
|
||||
@ -163,7 +163,11 @@ public:
|
||||
|
||||
//- Calculate the field min/max
|
||||
template<class Type>
|
||||
void calcMinMaxFields(const word& fieldName);
|
||||
void calcMinMaxFields
|
||||
(
|
||||
const word& fieldName,
|
||||
const modeType& mode
|
||||
);
|
||||
|
||||
//- Write the fieldMinMax
|
||||
virtual void write();
|
||||
@ -178,11 +182,6 @@ public:
|
||||
};
|
||||
|
||||
|
||||
// Template specialisation for scalar fields
|
||||
template<>
|
||||
void fieldMinMax::calcMinMaxFields<scalar>(const word& fieldName);
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
} // End namespace Foam
|
||||
|
||||
@ -27,81 +27,140 @@ License
|
||||
#include "volFields.H"
|
||||
#include "dictionary.H"
|
||||
#include "Time.H"
|
||||
#include "ListOps.H"
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
template<class Type>
|
||||
void Foam::fieldMinMax::calcMinMaxFields(const word& fieldName)
|
||||
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);
|
||||
switch (mode_)
|
||||
switch (mode)
|
||||
{
|
||||
case mdMag:
|
||||
{
|
||||
const scalar minValue = min(mag(field)).value();
|
||||
const scalar maxValue = max(mag(field)).value();
|
||||
const scalarField magField(mag(field));
|
||||
|
||||
labelList minIs(Pstream::nProcs());
|
||||
scalarList minVs(Pstream::nProcs());
|
||||
List<vector> minCs(Pstream::nProcs());
|
||||
minIs[procI] = findMin(magField);
|
||||
minVs[procI] = magField[minIs[procI]];
|
||||
minCs[procI] = field.mesh().C()[minIs[procI]];
|
||||
|
||||
Pstream::gatherList(minIs);
|
||||
Pstream::gatherList(minVs);
|
||||
Pstream::gatherList(minCs);
|
||||
|
||||
labelList maxIs(Pstream::nProcs());
|
||||
scalarList maxVs(Pstream::nProcs());
|
||||
List<vector> maxCs(Pstream::nProcs());
|
||||
maxIs[procI] = findMax(magField);
|
||||
maxVs[procI] = magField[maxIs[procI]];
|
||||
maxCs[procI] = field.mesh().C()[maxIs[procI]];
|
||||
|
||||
Pstream::gatherList(maxIs);
|
||||
Pstream::gatherList(maxVs);
|
||||
Pstream::gatherList(maxCs);
|
||||
|
||||
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];
|
||||
|
||||
if (Pstream::master())
|
||||
{
|
||||
if (write_)
|
||||
{
|
||||
fieldMinMaxFilePtr_()
|
||||
<< obr_.time().value() << tab
|
||||
<< fieldName << tab << minValue << tab << maxValue
|
||||
<< endl;
|
||||
<< obr_.time().value() << token::TAB
|
||||
<< fieldName << token::TAB
|
||||
<< minValue << token::TAB << minC << token::TAB
|
||||
<< maxValue << token::TAB << maxC << endl;
|
||||
}
|
||||
|
||||
if (log_)
|
||||
{
|
||||
Info<< "fieldMinMax output:" << nl
|
||||
<< " min(mag(" << fieldName << ")) = "
|
||||
<< minValue << nl
|
||||
Info<< " min(mag(" << fieldName << ")) = "
|
||||
<< minValue << " at position " << minC << nl
|
||||
<< " max(mag(" << fieldName << ")) = "
|
||||
<< maxValue << nl
|
||||
<< endl;
|
||||
<< maxValue << " at position " << maxC << nl;
|
||||
}
|
||||
}
|
||||
break;
|
||||
}
|
||||
case mdCmpt:
|
||||
{
|
||||
const Type minValue = min(field).value();
|
||||
const Type maxValue = max(field).value();
|
||||
List<Type> minVs(Pstream::nProcs());
|
||||
labelList minIs(Pstream::nProcs());
|
||||
List<vector> minCs(Pstream::nProcs());
|
||||
minIs[procI] = findMin(field);
|
||||
minVs[procI] = field[minIs[procI]];
|
||||
minCs[procI] = field.mesh().C()[minIs[procI]];
|
||||
|
||||
Pstream::gatherList(minIs);
|
||||
Pstream::gatherList(minVs);
|
||||
Pstream::gatherList(minCs);
|
||||
|
||||
List<Type> maxVs(Pstream::nProcs());
|
||||
labelList maxIs(Pstream::nProcs());
|
||||
List<vector> maxCs(Pstream::nProcs());
|
||||
maxIs[procI] = findMax(field);
|
||||
maxVs[procI] = field[maxIs[procI]];
|
||||
maxCs[procI] = field.mesh().C()[maxIs[procI]];
|
||||
|
||||
Pstream::gatherList(maxIs);
|
||||
Pstream::gatherList(maxVs);
|
||||
Pstream::gatherList(maxCs);
|
||||
|
||||
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];
|
||||
|
||||
|
||||
if (Pstream::master())
|
||||
{
|
||||
if (write_)
|
||||
{
|
||||
fieldMinMaxFilePtr_()
|
||||
<< obr_.time().value() << tab
|
||||
<< fieldName << tab << minValue << tab << maxValue
|
||||
<< endl;
|
||||
<< obr_.time().value() << token::TAB
|
||||
<< fieldName << token::TAB
|
||||
<< minValue << token::TAB << minC << token::TAB
|
||||
<< maxValue << token::TAB << maxC << endl;
|
||||
}
|
||||
|
||||
if (log_)
|
||||
{
|
||||
Info<< "fieldMinMax output:" << nl
|
||||
<< " cmptMin(" << fieldName << ") = "
|
||||
<< minValue << nl
|
||||
<< " cmptMax(" << fieldName << ") = "
|
||||
<< maxValue << nl
|
||||
<< endl;
|
||||
Info<< " min(" << fieldName << ") = "
|
||||
<< minValue << " at position " << minC << nl
|
||||
<< " max(" << fieldName << ") = "
|
||||
<< maxValue << " at position " << maxC << nl;
|
||||
}
|
||||
}
|
||||
break;
|
||||
}
|
||||
default:
|
||||
{
|
||||
FatalErrorIn
|
||||
(
|
||||
"Foam::fieldMinMax::calcMinMaxFields"
|
||||
"(const word& fieldName)"
|
||||
)<< "Unknown min/max mode: " << modeTypeNames_[mode_]
|
||||
<< exit(FatalError);
|
||||
FatalErrorIn("Foam::fieldMinMax::calcMinMaxFields(const word&)")
|
||||
<< "Unknown min/max mode: " << modeTypeNames_[mode_]
|
||||
<< exit(FatalError);
|
||||
}
|
||||
}
|
||||
}
|
||||
|
||||
Reference in New Issue
Block a user