SolverPerformance: Complete the integration of the templated SolverPerformance<Type>

Now solvers return solver performance information for all components
with backward compatibility provided by the "max" function which created
the scalar solverPerformance from the maximum component residuals from
the SolverPerformance<Type>.

The residuals functionObject has been upgraded to support
SolverPerformance<Type> so that now the initial residuals for all
(valid) components are tabulated, e.g. for the cavity tutorial case the
residuals for p, Ux and Uy are listed vs time.

Currently the residualControl option of pimpleControl and simpleControl
is supported in backward compatibility mode (only the maximum component
residual is considered) but in the future this will be upgraded to
support convergence control for the components individually.

This development started from patches provided by Bruno Santos, See
http://www.openfoam.org/mantisbt/view.php?id=1824
This commit is contained in:
Henry Weller
2015-11-10 08:50:11 +00:00
parent 5841ee3481
commit 1944b09bb5
31 changed files with 428 additions and 142 deletions

View File

@ -93,7 +93,7 @@ int main(int argc, char *argv[])
//DEqn.setComponentReference(1, 0, vector::X, 0);
//DEqn.setComponentReference(1, 0, vector::Z, 0);
initialResidual = DEqn.solve().initialResidual();
initialResidual = DEqn.solve().max().initialResidual();
if (!compactNormalStress)
{

View File

@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-2013 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2011-2015 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
@ -83,13 +83,35 @@ int main(int argc, char *argv[])
zeroGradientFvPatchSymmTensorField::typeName
);
solve
SolverPerformance<symmTensor> sP =
(
fvm::ddt(st)
+ fvm::div(phi, st)
- fvm::laplacian(dimensionedScalar("D", sqr(dimLength)/dimTime, 1), st)
solve
(
fvm::ddt(st)
+ fvm::div(phi, st)
- fvm::laplacian
(
dimensionedScalar("D", sqr(dimLength)/dimTime, 1),
st
)
==
dimensioned<symmTensor>
(
"source",
dimless/dimTime,
symmTensor(0, 2, 0, 1, 1.5, 0)
)
)
);
Info<< nl
<< "Detailed SolverPerformance<symmTensor>: " << nl
<< " " << sP << endl;
Info<< nl
<< "solverPerformanceDict: "
<< mesh.solverPerformanceDict() << endl;
return 0;
}

View File

@ -231,7 +231,7 @@ inline const word& token::wordToken() const
}
else
{
parseError("word");
parseError(word::typeName);
return word::null;
}
}
@ -254,7 +254,7 @@ inline const string& token::stringToken() const
}
else
{
parseError("string");
parseError(string::typeName);
return string::null;
}
}
@ -272,7 +272,7 @@ inline label token::labelToken() const
}
else
{
parseError("label");
parseError(pTraits<label>::typeName);
return 0;
}
}
@ -332,7 +332,7 @@ inline scalar token::scalarToken() const
}
else
{
parseError("scalar");
parseError(pTraits<scalar>::typeName);
return 0.0;
}
}

View File

@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-2013 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2011-2015 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
@ -126,6 +126,36 @@ void Foam::SolverPerformance<Type>::print
}
template<class Type>
void Foam::SolverPerformance<Type>::replace
(
const Foam::label cmpt,
const Foam::SolverPerformance<typename pTraits<Type>::cmptType>& sp
)
{
initialResidual_.replace(cmpt, sp.initialResidual());
finalResidual_.replace(cmpt, sp.finalResidual());
singular_[cmpt] = sp.singular();
}
template<class Type>
Foam::SolverPerformance<typename Foam::pTraits<Type>::cmptType>
Foam::SolverPerformance<Type>::max()
{
return SolverPerformance<typename pTraits<Type>::cmptType>
(
solverName_,
fieldName_,
cmptMax(initialResidual_),
cmptMax(finalResidual_),
noIterations_,
converged_,
singular()
);
}
template<class Type>
bool Foam::SolverPerformance<Type>::operator!=
(

View File

@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2012-2013 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2012-2015 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
@ -78,13 +78,15 @@ Ostream& operator<<
template<class Type>
class SolverPerformance
{
word solverName_;
word fieldName_;
Type initialResidual_;
Type finalResidual_;
label noIterations_;
bool converged_;
FixedList<bool, pTraits<Type>::nComponents> singular_;
// Private data
word solverName_;
word fieldName_;
Type initialResidual_;
Type finalResidual_;
label noIterations_;
bool converged_;
FixedList<bool, pTraits<Type>::nComponents> singular_;
public:
@ -220,6 +222,17 @@ public:
//- Print summary of solver performance to the given stream
void print(Ostream& os) const;
//- Replace component based on the minimal SolverPerformance
void replace
(
const label cmpt,
const SolverPerformance<typename pTraits<Type>::cmptType>& sp
);
//- Return the summary maximum of SolverPerformance<Type>
// Effectively it will mostly return solverPerformanceScalar
SolverPerformance<typename pTraits<Type>::cmptType> max();
// Member Operators
@ -229,7 +242,7 @@ public:
// Friend functions
//- Return the element-wise maximum of two SolverPerformance<Type>s
friend SolverPerformance<Type> max <Type>
friend SolverPerformance<Type> Foam::max <Type>
(
const SolverPerformance<Type>&,
const SolverPerformance<Type>&

View File

@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2012 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2012-2015 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
@ -36,4 +36,12 @@ namespace Foam
};
template<>
Foam::SolverPerformance<Foam::scalar>
Foam::SolverPerformance<Foam::scalar>::max()
{
return *this;
}
// ************************************************************************* //

View File

@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2012 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2012-2015 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
@ -42,8 +42,13 @@ SourceFiles
namespace Foam
{
typedef SolverPerformance<scalar> solverPerformance;
// Specialization of the max function for scalar object
template<>
SolverPerformance<scalar> SolverPerformance<scalar>::max();
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif

View File

@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-2012 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2011-2015 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
@ -25,7 +25,6 @@ License
#include "data.H"
#include "Time.H"
#include "solverPerformance.H"
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
@ -60,41 +59,4 @@ const Foam::dictionary& Foam::data::solverPerformanceDict() const
}
void Foam::data::setSolverPerformance
(
const word& name,
const solverPerformance& sp
) const
{
dictionary& dict = const_cast<dictionary&>(solverPerformanceDict());
List<solverPerformance> perfs;
if (prevTimeIndex_ != this->time().timeIndex())
{
// reset solver performance between iterations
prevTimeIndex_ = this->time().timeIndex();
dict.clear();
}
else
{
dict.readIfPresent(name, perfs);
}
// append to list
perfs.setSize(perfs.size()+1, sp);
dict.set(name, perfs);
}
void Foam::data::setSolverPerformance
(
const solverPerformance& sp
) const
{
setSolverPerformance(sp.fieldName(), sp);
}
// ************************************************************************* //

View File

@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-2012 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2011-2015 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
@ -91,16 +91,18 @@ public:
const dictionary& solverPerformanceDict() const;
//- Add/set the solverPerformance entry for the named field
template<class Type>
void setSolverPerformance
(
const word& name,
const solverPerformance&
const SolverPerformance<Type>&
) const;
//- Add/set the solverPerformance entry, using its fieldName
template<class Type>
void setSolverPerformance
(
const solverPerformance&
const SolverPerformance<Type>&
) const;
};
@ -111,6 +113,12 @@ public:
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#ifdef NoRepository
# include "dataTemplates.C"
#endif
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //

View File

@ -0,0 +1,71 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2015 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 "data.H"
#include "Time.H"
#include "solverPerformance.H"
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
template<class Type>
void Foam::data::setSolverPerformance
(
const word& name,
const SolverPerformance<Type>& sp
) const
{
dictionary& dict = const_cast<dictionary&>(solverPerformanceDict());
List<SolverPerformance<Type> > perfs;
if (prevTimeIndex_ != this->time().timeIndex())
{
// Reset solver performance between iterations
prevTimeIndex_ = this->time().timeIndex();
dict.clear();
}
else
{
dict.readIfPresent(name, perfs);
}
// Append to list
perfs.setSize(perfs.size()+1, sp);
dict.set(name, perfs);
}
template<class Type>
void Foam::data::setSolverPerformance
(
const SolverPerformance<Type>& sp
) const
{
setSolverPerformance(sp.fieldName(), sp);
}
// ************************************************************************* //

View File

@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-2013 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2011-2015 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
@ -38,7 +38,7 @@ const Scalar pTraits<Scalar>::max = ScalarVGREAT;
const Scalar pTraits<Scalar>::rootMin = -ScalarROOTVGREAT;
const Scalar pTraits<Scalar>::rootMax = ScalarROOTVGREAT;
const char* pTraits<Scalar>::componentNames[] = { "x" };
const char* pTraits<Scalar>::componentNames[] = { "" };
pTraits<Scalar>::pTraits(const Scalar& p)
:

View File

@ -50,6 +50,9 @@ public:
//- Component type
typedef Scalar cmptType;
//- Equivalent type of labels used for valid component indexing
typedef label labelType;
// Member constants
enum

View File

@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2011-2015 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
@ -31,7 +31,7 @@ const char* const Foam::pTraits<bool>::typeName = "bool";
const bool Foam::pTraits<bool>::zero = false;
const bool Foam::pTraits<bool>::one = true;
const char* Foam::pTraits<bool>::componentNames[] = { "x" };
const char* Foam::pTraits<bool>::componentNames[] = { "" };
Foam::pTraits<bool>::pTraits(const bool& p)
:

View File

@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2014 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2014-2015 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
@ -34,7 +34,7 @@ const int32_t Foam::pTraits<int32_t>::max = INT32_MAX;
const int32_t Foam::pTraits<int32_t>::rootMin = pTraits<int32_t>::min;
const int32_t Foam::pTraits<int32_t>::rootMax = pTraits<int32_t>::max;
const char* Foam::pTraits<int32_t>::componentNames[] = { "x" };
const char* Foam::pTraits<int32_t>::componentNames[] = { "" };
Foam::pTraits<int32_t>::pTraits(const int32_t& p)
:

View File

@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2014 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2014-2015 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
@ -34,7 +34,7 @@ const int64_t Foam::pTraits<int64_t>::max = INT64_MAX;
const int64_t Foam::pTraits<int64_t>::rootMin = pTraits<int64_t>::min;
const int64_t Foam::pTraits<int64_t>::rootMax = pTraits<int64_t>::max;
const char* Foam::pTraits<int64_t>::componentNames[] = { "x" };
const char* Foam::pTraits<int64_t>::componentNames[] = { "" };
Foam::pTraits<int64_t>::pTraits(const int64_t& p)
:

View File

@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2014 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2014-2015 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
@ -34,7 +34,7 @@ const uint32_t Foam::pTraits<uint32_t>::max = INT32_MAX;
const uint32_t Foam::pTraits<uint32_t>::rootMin = pTraits<uint32_t>::min;
const uint32_t Foam::pTraits<uint32_t>::rootMax = pTraits<uint32_t>::max;
const char* Foam::pTraits<uint32_t>::componentNames[] = { "x" };
const char* Foam::pTraits<uint32_t>::componentNames[] = { "" };
Foam::pTraits<uint32_t>::pTraits(const uint32_t& p)
:

View File

@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2014 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2014-2015 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
@ -34,7 +34,7 @@ const uint64_t Foam::pTraits<uint64_t>::max = INT64_MAX;
const uint64_t Foam::pTraits<uint64_t>::rootMin = pTraits<uint64_t>::min;
const uint64_t Foam::pTraits<uint64_t>::rootMax = pTraits<uint64_t>::max;
const char* Foam::pTraits<uint64_t>::componentNames[] = { "x" };
const char* Foam::pTraits<uint64_t>::componentNames[] = { "" };
Foam::pTraits<uint64_t>::pTraits(const uint64_t& p)
:

View File

@ -70,15 +70,15 @@ bool Foam::pimpleControl::criteriaSatisfied()
const label fieldI = applyToField(variableName);
if (fieldI != -1)
{
const List<solverPerformance> sp(iter().stream());
const scalar residual = sp.last().initialResidual();
scalar residual = 0;
const scalar firstResidual =
maxResidual(variableName, iter().stream(), residual);
checked = true;
if (storeIni)
{
residualControl_[fieldI].initialResidual =
sp.first().initialResidual();
residualControl_[fieldI].initialResidual = firstResidual;
}
const bool absCheck = residual < residualControl_[fieldI].absTol;

View File

@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-2012 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2011-2015 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
@ -59,8 +59,9 @@ bool Foam::simpleControl::criteriaSatisfied()
const label fieldI = applyToField(variableName);
if (fieldI != -1)
{
const List<solverPerformance> sp(iter().stream());
const scalar residual = sp.first().initialResidual();
scalar lastResidual = 0;
const scalar residual =
maxResidual(variableName, iter().stream(), lastResidual);
checked = true;

View File

@ -172,6 +172,45 @@ void Foam::solutionControl::storePrevIterFields() const
}
template<class Type>
void Foam::solutionControl::maxTypeResidual
(
const word& fieldName,
ITstream& data,
scalar& firstRes,
scalar& lastRes
) const
{
typedef GeometricField<Type, fvPatchField, volMesh> fieldType;
if (mesh_.foundObject<fieldType>(fieldName))
{
const List<SolverPerformance<Type> > sp(data);
firstRes = cmptMax(sp.first().initialResidual());
lastRes = cmptMax(sp.last().initialResidual());
}
}
Foam::scalar Foam::solutionControl::maxResidual
(
const word& fieldName,
ITstream& data,
scalar& lastRes
) const
{
scalar firstRes = 0;
maxTypeResidual<scalar>(fieldName, data, firstRes, lastRes);
maxTypeResidual<vector>(fieldName, data, firstRes, lastRes);
maxTypeResidual<sphericalTensor>(fieldName, data, firstRes, lastRes);
maxTypeResidual<symmTensor>(fieldName, data, firstRes, lastRes);
maxTypeResidual<tensor>(fieldName, data, firstRes, lastRes);
return firstRes;
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::solutionControl::solutionControl(fvMesh& mesh, const word& algorithmName)

View File

@ -122,6 +122,25 @@ protected:
template<class Type>
void storePrevIter() const;
template<class Type>
void maxTypeResidual
(
const word& fieldName,
ITstream& data,
scalar& firstRes,
scalar& lastRes
) const;
scalar maxResidual
(
const word& fieldName,
ITstream& data,
scalar& lastRes
) const;
private:
//- Disallow default bitwise copy construct
solutionControl(const solutionControl&);

View File

@ -808,11 +808,7 @@ Foam::fvMatrix<Type>::H() const
typename Type::labelType validComponents
(
pow
(
psi_.mesh().solutionD(),
pTraits<typename powProduct<Vector<label>, Type::rank>::type>::zero
)
psi_.mesh().template validComponents<Type>()
);
for (direction cmpt=0; cmpt<Type::nComponents; cmpt++)
@ -1350,7 +1346,7 @@ void Foam::checkMethod
template<class Type>
Foam::solverPerformance Foam::solve
Foam::SolverPerformance<Type> Foam::solve
(
fvMatrix<Type>& fvm,
const dictionary& solverControls
@ -1360,13 +1356,13 @@ Foam::solverPerformance Foam::solve
}
template<class Type>
Foam::solverPerformance Foam::solve
Foam::SolverPerformance<Type> Foam::solve
(
const tmp<fvMatrix<Type> >& tfvm,
const dictionary& solverControls
)
{
solverPerformance solverPerf =
SolverPerformance<Type> solverPerf =
const_cast<fvMatrix<Type>&>(tfvm()).solve(solverControls);
tfvm.clear();
@ -1376,15 +1372,15 @@ Foam::solverPerformance Foam::solve
template<class Type>
Foam::solverPerformance Foam::solve(fvMatrix<Type>& fvm)
Foam::SolverPerformance<Type> Foam::solve(fvMatrix<Type>& fvm)
{
return fvm.solve();
}
template<class Type>
Foam::solverPerformance Foam::solve(const tmp<fvMatrix<Type> >& tfvm)
Foam::SolverPerformance<Type> Foam::solve(const tmp<fvMatrix<Type> >& tfvm)
{
solverPerformance solverPerf =
SolverPerformance<Type> solverPerf =
const_cast<fvMatrix<Type>&>(tfvm()).solve();
tfvm.clear();

View File

@ -239,11 +239,11 @@ public:
//- Solve returning the solution statistics.
// Use the given solver controls
solverPerformance solve(const dictionary&);
SolverPerformance<Type> solve(const dictionary&);
//- Solve returning the solution statistics.
// Solver controls read from fvSolution
solverPerformance solve();
SolverPerformance<Type> solve();
};
@ -389,19 +389,19 @@ public:
//- Solve segregated or coupled returning the solution statistics.
// Use the given solver controls
solverPerformance solve(const dictionary&);
SolverPerformance<Type> solve(const dictionary&);
//- Solve segregated returning the solution statistics.
// Use the given solver controls
solverPerformance solveSegregated(const dictionary&);
SolverPerformance<Type> solveSegregated(const dictionary&);
//- Solve coupled returning the solution statistics.
// Use the given solver controls
solverPerformance solveCoupled(const dictionary&);
SolverPerformance<Type> solveCoupled(const dictionary&);
//- Solve returning the solution statistics.
// Solver controls read from fvSolution
solverPerformance solve();
SolverPerformance<Type> solve();
//- Return the matrix residual
tmp<Field<Type> > residual() const;
@ -549,14 +549,14 @@ void checkMethod
//- Solve returning the solution statistics given convergence tolerance
// Use the given solver controls
template<class Type>
solverPerformance solve(fvMatrix<Type>&, const dictionary&);
SolverPerformance<Type> solve(fvMatrix<Type>&, const dictionary&);
//- Solve returning the solution statistics given convergence tolerance,
// deleting temporary matrix after solution.
// Use the given solver controls
template<class Type>
solverPerformance solve
SolverPerformance<Type> solve
(
const tmp<fvMatrix<Type> >&,
const dictionary&
@ -566,14 +566,14 @@ solverPerformance solve
//- Solve returning the solution statistics given convergence tolerance
// Solver controls read fvSolution
template<class Type>
solverPerformance solve(fvMatrix<Type>&);
SolverPerformance<Type> solve(fvMatrix<Type>&);
//- Solve returning the solution statistics given convergence tolerance,
// deleting temporary matrix after solution.
// Solver controls read fvSolution
template<class Type>
solverPerformance solve(const tmp<fvMatrix<Type> >&);
SolverPerformance<Type> solve(const tmp<fvMatrix<Type> >&);
//- Return the correction form of the given matrix

View File

@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-2014 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2011-2015 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
@ -53,7 +53,7 @@ void Foam::fvMatrix<Type>::setComponentReference
template<class Type>
Foam::solverPerformance Foam::fvMatrix<Type>::solve
Foam::SolverPerformance<Type> Foam::fvMatrix<Type>::solve
(
const dictionary& solverControls
)
@ -71,7 +71,7 @@ Foam::solverPerformance Foam::fvMatrix<Type>::solve
{
if (maxIter == 0)
{
return solverPerformance();
return SolverPerformance<Type>();
}
}
@ -95,13 +95,13 @@ Foam::solverPerformance Foam::fvMatrix<Type>::solve
<< "; currently supported solver types are segregated and coupled"
<< exit(FatalIOError);
return solverPerformance();
return SolverPerformance<Type>();
}
}
template<class Type>
Foam::solverPerformance Foam::fvMatrix<Type>::solveSegregated
Foam::SolverPerformance<Type> Foam::fvMatrix<Type>::solveSegregated
(
const dictionary& solverControls
)
@ -118,7 +118,7 @@ Foam::solverPerformance Foam::fvMatrix<Type>::solveSegregated
GeometricField<Type, fvPatchField, volMesh>& psi =
const_cast<GeometricField<Type, fvPatchField, volMesh>&>(psi_);
solverPerformance solverPerfVec
SolverPerformance<Type> solverPerfVec
(
"fvMatrix<Type>::solveSegregated",
psi.name()
@ -135,11 +135,7 @@ Foam::solverPerformance Foam::fvMatrix<Type>::solveSegregated
typename Type::labelType validComponents
(
pow
(
psi.mesh().solutionD(),
pTraits<typename powProduct<Vector<label>, Type::rank>::type>::zero
)
psi.mesh().template validComponents<Type>()
);
for (direction cmpt=0; cmpt<Type::nComponents; cmpt++)
@ -200,13 +196,12 @@ Foam::solverPerformance Foam::fvMatrix<Type>::solveSegregated
solverControls
)->solve(psiCmpt, sourceCmpt, cmpt);
if (solverPerformance::debug)
if (SolverPerformance<Type>::debug)
{
solverPerf.print(Info.masterStream(this->mesh().comm()));
}
solverPerfVec = max(solverPerfVec, solverPerf);
solverPerfVec.solverName() = solverPerf.solverName();
solverPerfVec.replace(cmpt, solverPerf);
psi.internalField().replace(cmpt, psiCmpt);
diag() = saveDiag;
@ -221,7 +216,7 @@ Foam::solverPerformance Foam::fvMatrix<Type>::solveSegregated
template<class Type>
Foam::solverPerformance Foam::fvMatrix<Type>::solveCoupled
Foam::SolverPerformance<Type> Foam::fvMatrix<Type>::solveCoupled
(
const dictionary& solverControls
)
@ -274,9 +269,9 @@ Foam::solverPerformance Foam::fvMatrix<Type>::solveCoupled
psi.correctBoundaryConditions();
// psi.mesh().setSolverPerformance(psi.name(), solverPerf);
psi.mesh().setSolverPerformance(psi.name(), solverPerf);
return solverPerformance();
return solverPerf;
}
@ -299,7 +294,7 @@ Foam::fvMatrix<Type>::solver()
template<class Type>
Foam::solverPerformance Foam::fvMatrix<Type>::fvSolver::solve()
Foam::SolverPerformance<Type> Foam::fvMatrix<Type>::fvSolver::solve()
{
return solve
(
@ -316,7 +311,7 @@ Foam::solverPerformance Foam::fvMatrix<Type>::fvSolver::solve()
template<class Type>
Foam::solverPerformance Foam::fvMatrix<Type>::solve()
Foam::SolverPerformance<Type> Foam::fvMatrix<Type>::solve()
{
return solve
(

View File

@ -884,6 +884,14 @@ bool Foam::fvMesh::write() const
}
template<>
typename Foam::pTraits<Foam::sphericalTensor>::labelType
Foam::fvMesh::validComponents<Foam::sphericalTensor>() const
{
return Foam::pTraits<Foam::sphericalTensor>::labelType(1);
}
// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
bool Foam::fvMesh::operator!=(const fvMesh& bm) const

View File

@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-2014 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2011-2015 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
@ -324,6 +324,12 @@ public:
//- Return face deltas as surfaceVectorField
tmp<surfaceVectorField> delta() const;
//- Return a labelType of valid component indicators
// 1 : valid (solved)
// -1 : invalid (not solved)
template<class Type>
typename pTraits<Type>::labelType validComponents() const;
// Edit
@ -371,6 +377,11 @@ public:
};
template<>
typename pTraits<sphericalTensor>::labelType
fvMesh::validComponents<sphericalTensor>() const;
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace Foam
@ -378,6 +389,7 @@ public:
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#ifdef NoRepository
# include "fvMeshTemplates.C"
# include "fvPatchFvMeshTemplates.C"
#endif

View File

@ -0,0 +1,45 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2015 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 "fvMesh.H"
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
template<class Type>
typename Foam::pTraits<Type>::labelType Foam::fvMesh::validComponents() const
{
return pow
(
this->solutionD(),
pTraits
<
typename powProduct<Vector<label>,
pTraits<Type>::rank>::type
>::zero
);
}
// ************************************************************************* //

View File

@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2011-2015 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License

View File

@ -35,6 +35,7 @@ namespace Foam
defineTypeNameAndDebug(residuals, 0);
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::residuals::residuals
@ -98,7 +99,13 @@ void Foam::residuals::writeFileHeader(const label i)
forAll(fieldSet_, fieldI)
{
writeTabbed(file(), fieldSet_[fieldI]);
const word& fieldName = fieldSet_[fieldI];
writeFileHeader<scalar>(fieldName);
writeFileHeader<vector>(fieldName);
writeFileHeader<sphericalTensor>(fieldName);
writeFileHeader<symmTensor>(fieldName);
writeFileHeader<tensor>(fieldName);
}
file() << endl;
@ -107,21 +114,15 @@ void Foam::residuals::writeFileHeader(const label i)
void Foam::residuals::execute()
{
// Do nothing - only valid on write
}
{}
void Foam::residuals::end()
{
// Do nothing - only valid on write
}
{}
void Foam::residuals::timeSet()
{
// Do nothing - only valid on write
}
{}
void Foam::residuals::write()
@ -150,4 +151,5 @@ void Foam::residuals::write()
}
}
// ************************************************************************* //

View File

@ -115,6 +115,10 @@ protected:
//- Disallow default bitwise assignment
void operator=(const residuals&);
//- Output field header information
template<class Type>
void writeFileHeader(const word& fieldName);
//- Output file header information
virtual void writeFileHeader(const label i);

View File

@ -32,10 +32,37 @@ License
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
template<class Type>
void Foam::residuals::writeResidual
(
const word& fieldName
)
void Foam::residuals::writeFileHeader(const word& fieldName)
{
typedef GeometricField<Type, fvPatchField, volMesh> fieldType;
if (obr_.foundObject<fieldType>(fieldName))
{
const fieldType& field = obr_.lookupObject<fieldType>(fieldName);
const fvMesh& mesh = field.mesh();
typename pTraits<Type>::labelType validComponents
(
mesh.validComponents<Type>()
);
for (direction cmpt=0; cmpt<pTraits<Type>::nComponents; cmpt++)
{
if (component(validComponents, cmpt) != -1)
{
writeTabbed
(
file(),
fieldName + word(pTraits<Type>::componentNames[cmpt])
);
}
}
}
}
template<class Type>
void Foam::residuals::writeResidual(const word& fieldName)
{
typedef GeometricField<Type, fvPatchField, volMesh> fieldType;
@ -47,9 +74,25 @@ void Foam::residuals::writeResidual
if (solverDict.found(fieldName))
{
const List<solverPerformance> sp(solverDict.lookup(fieldName));
const scalar residual = sp.first().initialResidual();
file() << token::TAB << residual;
const List<SolverPerformance<Type> > sp
(
solverDict.lookup(fieldName)
);
const Type& residual = sp.first().initialResidual();
typename pTraits<Type>::labelType validComponents
(
mesh.validComponents<Type>()
);
for (direction cmpt=0; cmpt<pTraits<Type>::nComponents; cmpt++)
{
if (component(validComponents, cmpt) != -1)
{
file() << token::TAB << component(residual, cmpt);
}
}
}
}
}