mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-12-28 03:37:59 +00:00
355 lines
7.6 KiB
C
355 lines
7.6 KiB
C
/*---------------------------------------------------------------------------*\
|
|
========= |
|
|
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
|
\\ / O peration |
|
|
\\ / A nd | Copyright (C) 2011-2013 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/>.
|
|
|
|
Description
|
|
lduMatrix member operations.
|
|
|
|
\*---------------------------------------------------------------------------*/
|
|
|
|
#include "lduMatrix.H"
|
|
|
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
|
|
|
void Foam::lduMatrix::sumDiag()
|
|
{
|
|
const scalarField& Lower = const_cast<const lduMatrix&>(*this).lower();
|
|
const scalarField& Upper = const_cast<const lduMatrix&>(*this).upper();
|
|
scalarField& Diag = diag();
|
|
|
|
const labelUList& l = lduAddr().lowerAddr();
|
|
const labelUList& u = lduAddr().upperAddr();
|
|
|
|
for (register label face=0; face<l.size(); face++)
|
|
{
|
|
Diag[l[face]] += Lower[face];
|
|
Diag[u[face]] += Upper[face];
|
|
}
|
|
}
|
|
|
|
|
|
void Foam::lduMatrix::negSumDiag()
|
|
{
|
|
const scalarField& Lower = const_cast<const lduMatrix&>(*this).lower();
|
|
const scalarField& Upper = const_cast<const lduMatrix&>(*this).upper();
|
|
scalarField& Diag = diag();
|
|
|
|
const labelUList& l = lduAddr().lowerAddr();
|
|
const labelUList& u = lduAddr().upperAddr();
|
|
|
|
for (register label face=0; face<l.size(); face++)
|
|
{
|
|
Diag[l[face]] -= Lower[face];
|
|
Diag[u[face]] -= Upper[face];
|
|
}
|
|
}
|
|
|
|
|
|
void Foam::lduMatrix::sumMagOffDiag
|
|
(
|
|
scalarField& sumOff
|
|
) const
|
|
{
|
|
const scalarField& Lower = const_cast<const lduMatrix&>(*this).lower();
|
|
const scalarField& Upper = const_cast<const lduMatrix&>(*this).upper();
|
|
|
|
const labelUList& l = lduAddr().lowerAddr();
|
|
const labelUList& u = lduAddr().upperAddr();
|
|
|
|
for (register label face = 0; face < l.size(); face++)
|
|
{
|
|
sumOff[u[face]] += mag(Lower[face]);
|
|
sumOff[l[face]] += mag(Upper[face]);
|
|
}
|
|
}
|
|
|
|
|
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
|
|
|
void Foam::lduMatrix::operator=(const lduMatrix& A)
|
|
{
|
|
if (this == &A)
|
|
{
|
|
FatalError
|
|
<< "lduMatrix::operator=(const lduMatrix&) : "
|
|
<< "attempted assignment to self"
|
|
<< abort(FatalError);
|
|
}
|
|
|
|
if (A.lowerPtr_)
|
|
{
|
|
lower() = A.lower();
|
|
}
|
|
else if (lowerPtr_)
|
|
{
|
|
delete lowerPtr_;
|
|
lowerPtr_ = NULL;
|
|
}
|
|
|
|
if (A.upperPtr_)
|
|
{
|
|
upper() = A.upper();
|
|
}
|
|
else if (upperPtr_)
|
|
{
|
|
delete upperPtr_;
|
|
upperPtr_ = NULL;
|
|
}
|
|
|
|
if (A.diagPtr_)
|
|
{
|
|
diag() = A.diag();
|
|
}
|
|
}
|
|
|
|
|
|
void Foam::lduMatrix::negate()
|
|
{
|
|
if (lowerPtr_)
|
|
{
|
|
lowerPtr_->negate();
|
|
}
|
|
|
|
if (upperPtr_)
|
|
{
|
|
upperPtr_->negate();
|
|
}
|
|
|
|
if (diagPtr_)
|
|
{
|
|
diagPtr_->negate();
|
|
}
|
|
}
|
|
|
|
|
|
void Foam::lduMatrix::operator+=(const lduMatrix& A)
|
|
{
|
|
if (A.diagPtr_)
|
|
{
|
|
diag() += A.diag();
|
|
}
|
|
|
|
if (symmetric() && A.symmetric())
|
|
{
|
|
upper() += A.upper();
|
|
}
|
|
else if (symmetric() && A.asymmetric())
|
|
{
|
|
if (upperPtr_)
|
|
{
|
|
lower();
|
|
}
|
|
else
|
|
{
|
|
upper();
|
|
}
|
|
|
|
upper() += A.upper();
|
|
lower() += A.lower();
|
|
}
|
|
else if (asymmetric() && A.symmetric())
|
|
{
|
|
if (A.upperPtr_)
|
|
{
|
|
lower() += A.upper();
|
|
upper() += A.upper();
|
|
}
|
|
else
|
|
{
|
|
lower() += A.lower();
|
|
upper() += A.lower();
|
|
}
|
|
|
|
}
|
|
else if (asymmetric() && A.asymmetric())
|
|
{
|
|
lower() += A.lower();
|
|
upper() += A.upper();
|
|
}
|
|
else if (diagonal())
|
|
{
|
|
if (A.upperPtr_)
|
|
{
|
|
upper() = A.upper();
|
|
}
|
|
|
|
if (A.lowerPtr_)
|
|
{
|
|
lower() = A.lower();
|
|
}
|
|
}
|
|
else if (A.diagonal())
|
|
{
|
|
}
|
|
else
|
|
{
|
|
if (debug > 1)
|
|
{
|
|
WarningIn("lduMatrix::operator+=(const lduMatrix& A)")
|
|
<< "Unknown matrix type combination" << nl
|
|
<< " this :"
|
|
<< " diagonal:" << diagonal()
|
|
<< " symmetric:" << symmetric()
|
|
<< " asymmetric:" << asymmetric() << nl
|
|
<< " A :"
|
|
<< " diagonal:" << A.diagonal()
|
|
<< " symmetric:" << A.symmetric()
|
|
<< " asymmetric:" << A.asymmetric()
|
|
<< endl;
|
|
}
|
|
}
|
|
}
|
|
|
|
|
|
void Foam::lduMatrix::operator-=(const lduMatrix& A)
|
|
{
|
|
if (A.diagPtr_)
|
|
{
|
|
diag() -= A.diag();
|
|
}
|
|
|
|
if (symmetric() && A.symmetric())
|
|
{
|
|
upper() -= A.upper();
|
|
}
|
|
else if (symmetric() && A.asymmetric())
|
|
{
|
|
if (upperPtr_)
|
|
{
|
|
lower();
|
|
}
|
|
else
|
|
{
|
|
upper();
|
|
}
|
|
|
|
upper() -= A.upper();
|
|
lower() -= A.lower();
|
|
}
|
|
else if (asymmetric() && A.symmetric())
|
|
{
|
|
if (A.upperPtr_)
|
|
{
|
|
lower() -= A.upper();
|
|
upper() -= A.upper();
|
|
}
|
|
else
|
|
{
|
|
lower() -= A.lower();
|
|
upper() -= A.lower();
|
|
}
|
|
|
|
}
|
|
else if (asymmetric() && A.asymmetric())
|
|
{
|
|
lower() -= A.lower();
|
|
upper() -= A.upper();
|
|
}
|
|
else if (diagonal())
|
|
{
|
|
if (A.upperPtr_)
|
|
{
|
|
upper() = -A.upper();
|
|
}
|
|
|
|
if (A.lowerPtr_)
|
|
{
|
|
lower() = -A.lower();
|
|
}
|
|
}
|
|
else if (A.diagonal())
|
|
{
|
|
}
|
|
else
|
|
{
|
|
if (debug > 1)
|
|
{
|
|
WarningIn("lduMatrix::operator-=(const lduMatrix& A)")
|
|
<< "Unknown matrix type combination" << nl
|
|
<< " this :"
|
|
<< " diagonal:" << diagonal()
|
|
<< " symmetric:" << symmetric()
|
|
<< " asymmetric:" << asymmetric() << nl
|
|
<< " A :"
|
|
<< " diagonal:" << A.diagonal()
|
|
<< " symmetric:" << A.symmetric()
|
|
<< " asymmetric:" << A.asymmetric()
|
|
<< endl;
|
|
}
|
|
}
|
|
}
|
|
|
|
|
|
void Foam::lduMatrix::operator*=(const scalarField& sf)
|
|
{
|
|
if (diagPtr_)
|
|
{
|
|
*diagPtr_ *= sf;
|
|
}
|
|
|
|
if (upperPtr_)
|
|
{
|
|
scalarField& upper = *upperPtr_;
|
|
|
|
const labelUList& l = lduAddr().lowerAddr();
|
|
|
|
for (register label face=0; face<upper.size(); face++)
|
|
{
|
|
upper[face] *= sf[l[face]];
|
|
}
|
|
}
|
|
|
|
if (lowerPtr_)
|
|
{
|
|
scalarField& lower = *lowerPtr_;
|
|
|
|
const labelUList& u = lduAddr().upperAddr();
|
|
|
|
for (register label face=0; face<lower.size(); face++)
|
|
{
|
|
lower[face] *= sf[u[face]];
|
|
}
|
|
}
|
|
}
|
|
|
|
|
|
void Foam::lduMatrix::operator*=(scalar s)
|
|
{
|
|
if (diagPtr_)
|
|
{
|
|
*diagPtr_ *= s;
|
|
}
|
|
|
|
if (upperPtr_)
|
|
{
|
|
*upperPtr_ *= s;
|
|
}
|
|
|
|
if (lowerPtr_)
|
|
{
|
|
*lowerPtr_ *= s;
|
|
}
|
|
}
|
|
|
|
|
|
// ************************************************************************* //
|