mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
ENH: lduMatrix: added I/O
This commit is contained in:
@ -2,7 +2,7 @@
|
|||||||
========= |
|
========= |
|
||||||
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||||
\\ / O peration |
|
\\ / O peration |
|
||||||
\\ / A nd | Copyright (C) 2011-2012 OpenFOAM Foundation
|
\\ / A nd | Copyright (C) 2011-2013 OpenFOAM Foundation
|
||||||
\\/ M anipulation |
|
\\/ M anipulation |
|
||||||
-------------------------------------------------------------------------------
|
-------------------------------------------------------------------------------
|
||||||
License
|
License
|
||||||
@ -25,6 +25,7 @@ License
|
|||||||
|
|
||||||
#include "lduMatrix.H"
|
#include "lduMatrix.H"
|
||||||
#include "IOstreams.H"
|
#include "IOstreams.H"
|
||||||
|
#include "Switch.H"
|
||||||
|
|
||||||
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
|
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
|
||||||
|
|
||||||
@ -116,17 +117,27 @@ Foam::lduMatrix::lduMatrix(lduMatrix& A, bool reUse)
|
|||||||
}
|
}
|
||||||
|
|
||||||
|
|
||||||
Foam::lduMatrix::lduMatrix
|
Foam::lduMatrix::lduMatrix(const lduMesh& mesh, Istream& is)
|
||||||
(
|
|
||||||
const lduMesh& mesh,
|
|
||||||
Istream& is
|
|
||||||
)
|
|
||||||
:
|
:
|
||||||
lduMesh_(mesh),
|
lduMesh_(mesh)
|
||||||
lowerPtr_(new scalarField(is)),
|
{
|
||||||
diagPtr_(new scalarField(is)),
|
Switch hasLow(is);
|
||||||
upperPtr_(new scalarField(is))
|
Switch hasDiag(is);
|
||||||
{}
|
Switch hasUp(is);
|
||||||
|
|
||||||
|
if (hasLow)
|
||||||
|
{
|
||||||
|
lowerPtr_ = new scalarField(is);
|
||||||
|
}
|
||||||
|
if (hasDiag)
|
||||||
|
{
|
||||||
|
diagPtr_ = new scalarField(is);
|
||||||
|
}
|
||||||
|
if (hasUp)
|
||||||
|
{
|
||||||
|
upperPtr_ = new scalarField(is);
|
||||||
|
}
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
Foam::lduMatrix::~lduMatrix()
|
Foam::lduMatrix::~lduMatrix()
|
||||||
@ -252,25 +263,86 @@ const Foam::scalarField& Foam::lduMatrix::upper() const
|
|||||||
|
|
||||||
Foam::Ostream& Foam::operator<<(Ostream& os, const lduMatrix& ldum)
|
Foam::Ostream& Foam::operator<<(Ostream& os, const lduMatrix& ldum)
|
||||||
{
|
{
|
||||||
if (ldum.lowerPtr_)
|
Switch hasLow = ldum.hasLower();
|
||||||
|
Switch hasDiag = ldum.hasDiag();
|
||||||
|
Switch hasUp = ldum.hasUpper();
|
||||||
|
|
||||||
|
os << hasLow << token::SPACE << hasDiag << token::SPACE
|
||||||
|
<< hasUp << token::SPACE;
|
||||||
|
|
||||||
|
if (hasLow)
|
||||||
{
|
{
|
||||||
os << "Lower triangle = "
|
os << ldum.lower();
|
||||||
<< *ldum.lowerPtr_
|
|
||||||
<< endl << endl;
|
|
||||||
}
|
}
|
||||||
|
|
||||||
if (ldum.diagPtr_)
|
if (hasDiag)
|
||||||
{
|
{
|
||||||
os << "diagonal = "
|
os << ldum.diag();
|
||||||
<< *ldum.diagPtr_
|
|
||||||
<< endl << endl;
|
|
||||||
}
|
}
|
||||||
|
|
||||||
if (ldum.upperPtr_)
|
if (hasUp)
|
||||||
{
|
{
|
||||||
os << "Upper triangle = "
|
os << ldum.upper();
|
||||||
<< *ldum.upperPtr_
|
}
|
||||||
<< endl << endl;
|
|
||||||
|
os.check("Ostream& operator<<(Ostream&, const lduMatrix&");
|
||||||
|
|
||||||
|
return os;
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
Foam::Ostream& Foam::operator<<(Ostream& os, const InfoProxy<lduMatrix>& ip)
|
||||||
|
{
|
||||||
|
const lduMatrix& ldum = ip.t_;
|
||||||
|
|
||||||
|
Switch hasLow = ldum.hasLower();
|
||||||
|
Switch hasDiag = ldum.hasDiag();
|
||||||
|
Switch hasUp = ldum.hasUpper();
|
||||||
|
|
||||||
|
os << "Lower:" << hasLow
|
||||||
|
<< " Diag:" << hasDiag
|
||||||
|
<< " Upper:" << hasUp << endl;
|
||||||
|
|
||||||
|
if (hasLow)
|
||||||
|
{
|
||||||
|
os << "lower:" << ldum.lower().size() << endl;
|
||||||
|
}
|
||||||
|
if (hasDiag)
|
||||||
|
{
|
||||||
|
os << "diag :" << ldum.diag().size() << endl;
|
||||||
|
}
|
||||||
|
if (hasUp)
|
||||||
|
{
|
||||||
|
os << "upper:" << ldum.upper().size() << endl;
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
if (hasLow)
|
||||||
|
{
|
||||||
|
os << "lower contents:" << endl;
|
||||||
|
forAll(ldum.lower(), i)
|
||||||
|
{
|
||||||
|
os << "i:" << i << "\t" << ldum.lower()[i] << endl;
|
||||||
|
}
|
||||||
|
os << endl;
|
||||||
|
}
|
||||||
|
if (hasDiag)
|
||||||
|
{
|
||||||
|
os << "diag contents:" << endl;
|
||||||
|
forAll(ldum.diag(), i)
|
||||||
|
{
|
||||||
|
os << "i:" << i << "\t" << ldum.diag()[i] << endl;
|
||||||
|
}
|
||||||
|
os << endl;
|
||||||
|
}
|
||||||
|
if (hasUp)
|
||||||
|
{
|
||||||
|
os << "upper contents:" << endl;
|
||||||
|
forAll(ldum.upper(), i)
|
||||||
|
{
|
||||||
|
os << "i:" << i << "\t" << ldum.upper()[i] << endl;
|
||||||
|
}
|
||||||
|
os << endl;
|
||||||
}
|
}
|
||||||
|
|
||||||
os.check("Ostream& operator<<(Ostream&, const lduMatrix&");
|
os.check("Ostream& operator<<(Ostream&, const lduMatrix&");
|
||||||
|
|||||||
@ -2,7 +2,7 @@
|
|||||||
========= |
|
========= |
|
||||||
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||||
\\ / O peration |
|
\\ / O peration |
|
||||||
\\ / A nd | Copyright (C) 2011-2012 OpenFOAM Foundation
|
\\ / A nd | Copyright (C) 2011-2013 OpenFOAM Foundation
|
||||||
\\/ M anipulation |
|
\\/ M anipulation |
|
||||||
-------------------------------------------------------------------------------
|
-------------------------------------------------------------------------------
|
||||||
License
|
License
|
||||||
@ -58,6 +58,7 @@ SourceFiles
|
|||||||
#include "autoPtr.H"
|
#include "autoPtr.H"
|
||||||
#include "runTimeSelectionTables.H"
|
#include "runTimeSelectionTables.H"
|
||||||
#include "solverPerformance.H"
|
#include "solverPerformance.H"
|
||||||
|
#include "InfoProxy.H"
|
||||||
|
|
||||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
@ -691,6 +692,16 @@ public:
|
|||||||
tmp<Field<Type> > faceH(const tmp<Field<Type> >&) const;
|
tmp<Field<Type> > faceH(const tmp<Field<Type> >&) const;
|
||||||
|
|
||||||
|
|
||||||
|
// Info
|
||||||
|
|
||||||
|
//- Return info proxy.
|
||||||
|
// Used to print matrix information to a stream
|
||||||
|
InfoProxy<lduMatrix> info() const
|
||||||
|
{
|
||||||
|
return *this;
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
// Member operators
|
// Member operators
|
||||||
|
|
||||||
void operator=(const lduMatrix&);
|
void operator=(const lduMatrix&);
|
||||||
@ -707,6 +718,7 @@ public:
|
|||||||
// Ostream operator
|
// Ostream operator
|
||||||
|
|
||||||
friend Ostream& operator<<(Ostream&, const lduMatrix&);
|
friend Ostream& operator<<(Ostream&, const lduMatrix&);
|
||||||
|
friend Ostream& operator<<(Ostream&, const InfoProxy<lduMatrix>&);
|
||||||
};
|
};
|
||||||
|
|
||||||
|
|
||||||
|
|||||||
Reference in New Issue
Block a user