![]() |
OGS
|
Some general linear algebra functionality.
By using the provided functions linear algebra capabilities can be used for different matrix and vector types in a way that is agnostic towards the specific type used.
For documentation, refer to that of the templated method. All specializations or overload must behave in the same way.
Functions | |
void | setLocalAccessibleVector (PETScVector const &x) |
void | set (PETScVector &x, double const a) |
void | copy (PETScVector const &x, PETScVector &y) |
void | scale (PETScVector &x, double const a) |
void | aypx (PETScVector &y, double const a, PETScVector const &x) |
void | axpy (PETScVector &y, double const a, PETScVector const &x) |
void | axpby (PETScVector &y, double const a, double const b, PETScVector const &x) |
template<> | |
void | componentwiseDivide (PETScVector &w, PETScVector const &x, PETScVector const &y) |
template<> | |
double | norm1 (PETScVector const &x) |
template<> | |
double | norm2 (PETScVector const &x) |
template<> | |
double | normMax (PETScVector const &x) |
void | copy (PETScMatrix const &A, PETScMatrix &B) |
void | scale (PETScMatrix &A, double const a) |
void | aypx (PETScMatrix &Y, double const a, PETScMatrix const &X) |
void | axpy (PETScMatrix &Y, double const a, PETScMatrix const &X) |
void | matMult (PETScMatrix const &A, PETScVector const &x, PETScVector &y) |
void | matMultAdd (PETScMatrix const &A, PETScVector const &v1, PETScVector const &v2, PETScVector &v3) |
void | finalizeAssembly (PETScMatrix &A) |
void | finalizeAssembly (PETScVector &x) |
template<typename MatrixOrVector > | |
void | copy (MatrixOrVector const &x, MatrixOrVector &y) |
Copies x to y . More... | |
template<typename MatrixOrVector > | |
void | scale (MatrixOrVector &x, double const a) |
Scales x by a . More... | |
template<typename MatrixOrVector > | |
void | aypx (MatrixOrVector &y, double const a, MatrixOrVector const &x) |
Computes y = a \cdot y + x . More... | |
template<typename MatrixOrVector > | |
void | axpy (MatrixOrVector &y, double const a, MatrixOrVector const &x) |
Computes y = a \cdot x + y . More... | |
template<typename MatrixOrVector > | |
void | axpby (MatrixOrVector &y, double const a, double const b, MatrixOrVector const &x) |
Computes y = a \cdot x + b \cdot y . More... | |
template<typename MatrixOrVector > | |
void | componentwiseDivide (MatrixOrVector &w, MatrixOrVector const &x, MatrixOrVector const &y) |
Computes w = x/y componentwise. More... | |
template<typename MatrixOrVector > | |
double | norm1 (MatrixOrVector const &x) |
Computes the Manhattan norm of x . More... | |
template<typename MatrixOrVector > | |
double | norm2 (MatrixOrVector const &x) |
Computes the Euclidean norm of x . More... | |
template<typename MatrixOrVector > | |
double | normMax (MatrixOrVector const &x) |
Computes the maximum norm of x . More... | |
template<typename MatrixOrVector > | |
double | norm (MatrixOrVector const &x, MathLib::VecNormType type) |
template<typename Matrix > | |
void | finalizeAssembly (Matrix &) |
template<typename Matrix , typename Vector > | |
void | matMult (Matrix const &A, Vector const &x, Vector &y) |
template<typename Matrix , typename Vector > | |
void | matMultAdd (Matrix const &A, Vector const &v1, Vector const &v2, Vector &v3) |
void MathLib::LinAlg::axpby | ( | MatrixOrVector & | y, |
double const | a, | ||
double const | b, | ||
MatrixOrVector const & | x | ||
) |
void MathLib::LinAlg::axpby | ( | PETScVector & | y, |
double const | a, | ||
double const | b, | ||
PETScVector const & | x | ||
) |
Definition at line 64 of file LinAlg.cpp.
References MathLib::PETScVector::getRawVector().
Referenced by NumLib::TimeDiscretization::getXdot().
void MathLib::LinAlg::axpy | ( | MatrixOrVector & | y, |
double const | a, | ||
MatrixOrVector const & | x | ||
) |
void MathLib::LinAlg::axpy | ( | PETScMatrix & | Y, |
double const | a, | ||
PETScMatrix const & | X | ||
) |
void MathLib::LinAlg::axpy | ( | PETScVector & | y, |
double const | a, | ||
PETScVector const & | x | ||
) |
Definition at line 57 of file LinAlg.cpp.
References MathLib::PETScVector::getRawVector().
Referenced by NumLib::NonlinearSolver< NonlinearSolverTag::Picard >::calculateNonEquilibriumInitialResiduum(), NumLib::TimeDiscretization::computeRelativeChangeFromPreviousTimestep(), NumLib::MatrixTranslatorGeneral< ODESystemTag::FirstOrderImplicitQuasilinear >::computeResidual(), ProcessLib::computeResiduum(), NumLib::NonlinearSolver< NonlinearSolverTag::Newton >::solve(), NumLib::NonlinearSolver< NonlinearSolverTag::Picard >::solve(), and ProcessLib::TimeLoop::solveCoupledEquationSystemsByStaggeredScheme().
void MathLib::LinAlg::aypx | ( | MatrixOrVector & | y, |
double const | a, | ||
MatrixOrVector const & | x | ||
) |
void MathLib::LinAlg::aypx | ( | PETScMatrix & | Y, |
double const | a, | ||
PETScMatrix const & | X | ||
) |
void MathLib::LinAlg::aypx | ( | PETScVector & | y, |
double const | a, | ||
PETScVector const & | x | ||
) |
Definition at line 50 of file LinAlg.cpp.
References MathLib::PETScVector::getRawVector().
Referenced by NumLib::MatrixTranslatorGeneral< ODESystemTag::FirstOrderImplicitQuasilinear >::computeA(), and ProcessLib::ComponentTransport::ComponentTransportProcess::solveReactionEquation().
void MathLib::LinAlg::componentwiseDivide | ( | MatrixOrVector & | w, |
MatrixOrVector const & | x, | ||
MatrixOrVector const & | y | ||
) |
Computes w = x/y componentwise.
void MathLib::LinAlg::componentwiseDivide | ( | PETScVector & | w, |
PETScVector const & | x, | ||
PETScVector const & | y | ||
) |
Definition at line 73 of file LinAlg.cpp.
References MathLib::PETScVector::getRawVector().
Referenced by NumLib::LocalLinearLeastSquaresExtrapolator::extrapolate().
void MathLib::LinAlg::copy | ( | MatrixOrVector const & | x, |
MatrixOrVector & | y | ||
) |
void MathLib::LinAlg::copy | ( | PETScMatrix const & | A, |
PETScMatrix & | B | ||
) |
Definition at line 111 of file LinAlg.cpp.
void MathLib::LinAlg::copy | ( | PETScVector const & | x, |
PETScVector & | y | ||
) |
Definition at line 37 of file LinAlg.cpp.
References MathLib::PETScVector::getRawVector(), and MathLib::PETScVector::shallowCopy().
Referenced by GMSHPrefsDialog::GMSHPrefsDialog(), FileIO::Gocad::GocadSGridReader::GocadSGridReader(), GeoLib::Raster::Raster(), FileIO::Gocad::GocadSGridReader::addFaceSetQuad(), FileIO::Gocad::GocadSGridReader::addGocadPropertiesToMesh(), addIntegrationPointData(), addIntegrationPointMetaData(), MeshLib::MeshLayerMapper::addLayerToMesh(), MeshLib::addLayerToMesh(), MeshLib::addPropertyToMesh(), FileIO::SwmmInterface::addResultsToMesh(), MeshGeoToolsLib::appendLinesAlongPolylines(), ProcessLib::LIE::PostProcessTool::calculateTotalDisplacement(), NumLib::MatrixTranslatorGeneral< ODESystemTag::FirstOrderImplicitQuasilinear >::computeA(), NumLib::MatrixTranslatorGeneral< ODESystemTag::FirstOrderImplicitQuasilinear >::computeJacobian(), NumLib::TimeDiscretization::computeRelativeChangeFromPreviousTimestep(), ProcessLib::TimeLoop::computeTimeStepping(), MeshLib::VtkMeshConverter::convertTypedArray(), ProcessLib::LIE::PostProcessTool::copyPropertyValues(), MeshLib::createQuadraticOrderMesh(), ChemistryLib::createReactionRates(), BaseLib::excludeObjectCopy(), ProcessLib::ComponentTransport::ComponentTransportProcess::extrapolateIntegrationPointValuesToNodes(), flipRaster(), MeshLib::getBaseNodes(), ProcessLib::PhaseFieldIrreversibleDamageOracleBoundaryCondition::getEssentialBCValues(), ProcessLib::LIE::getFractureMatrixDataInMesh(), LayeredMeshGenerator::getMesh(), NumLib::SimpleMatrixVectorProvider::getVector(), NumLib::BackwardEuler::getWeightedOldX(), MaterialLib::Solids::MFront::MFront< DisplacementDim >::integrateStress(), main(), BaseLib::operator<<(), operator<<(), FileIO::Gocad::operator<<(), ChemistryLib::PhreeqcIOData::anonymous_namespace{PhreeqcIO.cpp}::operator<<(), ApplicationUtils::NodeWiseMeshPartitioner::partitionOtherMesh(), ProcessLib::HeatTransportBHE::HeatTransportBHEProcess::postTimestepConcreteProcess(), MeshLib::PropertyVector< T * >::print(), ApplicationUtils::NodeWiseMeshPartitioner::processPartition(), MathLib::TemplatePoint< T, DIM >::read(), FileIO::Gocad::GocadSGridReader::readElementPropertiesBinary(), FileIO::GMSInterface::readGMS3DMMesh(), FileIO::SwmmInterface::readSwmmInputToLineMesh(), FileIO::TetGenInterface::readTetGenMesh(), GeoLib::MinimalBoundingSphere::recurseCalculation(), ProcessLib::TimeLoop::setCoupledSolutions(), ProcessLib::Process::setInitialConditions(), DiagramList::setList(), ProcessLib::Deformation::solidMaterialInternalToSecondaryVariables(), ProcessLib::Deformation::solidMaterialInternalVariablesToIntegrationPointWriter(), NumLib::NonlinearSolver< NonlinearSolverTag::Newton >::solve(), NumLib::NonlinearSolver< NonlinearSolverTag::Picard >::solve(), NumLib::PETScNonlinearSolver::solve(), ProcessLib::TimeLoop::solveCoupledEquationSystemsByStaggeredScheme(), ProcessLib::ComponentTransport::ComponentTransportProcess::solveReactionEquation(), BaseLib::splitString(), ProcessLib::PhaseField::PhaseFieldProcess< DisplacementDim >::updateConstraints(), anonymous_namespace{IdentifySubdomainMesh.cpp}::updateOrCheckExistingSubdomainProperty(), MathLib::TemplatePoint< T, DIM >::write(), and writeDataToMesh().
void MathLib::LinAlg::finalizeAssembly | ( | Matrix & | ) |
void MathLib::LinAlg::finalizeAssembly | ( | PETScMatrix & | A | ) |
Definition at line 162 of file LinAlg.cpp.
References MathLib::PETScMatrix::finalizeAssembly().
Referenced by detail::applyKnownSolutions(), NumLib::TimeDiscretizedODESystem< ODESystemTag::FirstOrderImplicitQuasilinear, NonlinearSolverTag::Newton >::assemble(), NumLib::TimeDiscretizedODESystem< ODESystemTag::FirstOrderImplicitQuasilinear, NonlinearSolverTag::Picard >::assemble(), NumLib::LocalLinearLeastSquaresExtrapolator::calculateResiduals(), NumLib::LocalLinearLeastSquaresExtrapolator::extrapolate(), ProcessLib::ComponentTransport::ComponentTransportProcess::extrapolateIntegrationPointValuesToNodes(), ProcessLib::Process::setInitialConditions(), and ProcessLib::SmallDeformation::writeMaterialForces().
void MathLib::LinAlg::finalizeAssembly | ( | PETScVector & | x | ) |
Definition at line 167 of file LinAlg.cpp.
References MathLib::PETScVector::finalizeAssembly().
void MathLib::LinAlg::matMult | ( | Matrix const & | A, |
Vector const & | x, | ||
Vector & | y | ||
) |
void MathLib::LinAlg::matMult | ( | PETScMatrix const & | A, |
PETScVector const & | x, | ||
PETScVector & | y | ||
) |
Definition at line 141 of file LinAlg.cpp.
References MathLib::PETScMatrix::getRawMatrix(), MathLib::PETScVector::getRawVector(), and MathLib::PETScVector::shallowCopy().
Referenced by NumLib::NonlinearSolver< NonlinearSolverTag::Picard >::calculateNonEquilibriumInitialResiduum(), NumLib::MatrixTranslatorGeneral< ODESystemTag::FirstOrderImplicitQuasilinear >::computeResidual(), ProcessLib::computeResiduum(), NumLib::NonlinearSolver< NonlinearSolverTag::Picard >::solve(), and ProcessLib::ComponentTransport::ComponentTransportProcess::solveReactionEquation().
void MathLib::LinAlg::matMultAdd | ( | Matrix const & | A, |
Vector const & | v1, | ||
Vector const & | v2, | ||
Vector & | v3 | ||
) |
void MathLib::LinAlg::matMultAdd | ( | PETScMatrix const & | A, |
PETScVector const & | v1, | ||
PETScVector const & | v2, | ||
PETScVector & | v3 | ||
) |
Definition at line 151 of file LinAlg.cpp.
References MathLib::PETScMatrix::getRawMatrix(), MathLib::PETScVector::getRawVector(), and MathLib::PETScVector::shallowCopy().
Referenced by NumLib::MatrixTranslatorGeneral< ODESystemTag::FirstOrderImplicitQuasilinear >::computeResidual(), ProcessLib::computeResiduum(), and NumLib::MatrixTranslatorGeneral< ODESystemTag::FirstOrderImplicitQuasilinear >::computeRhs().
double MathLib::LinAlg::norm | ( | MatrixOrVector const & | x, |
MathLib::VecNormType | type | ||
) |
Definition at line 88 of file LinAlg.h.
References MathLib::INFINITY_N, norm1(), MathLib::NORM1, norm2(), MathLib::NORM2, normMax(), and OGS_FATAL.
Referenced by LayeredVolume::addLayerBoundaries(), MathLib::calcProjPntToLineAndDists(), NumLib::ConvergenceCriterionDeltaX::checkDeltaX(), NumLib::ConvergenceCriterionResidual::checkDeltaX(), NumLib::ConvergenceCriterionResidual::checkResidual(), NumLib::TimeDiscretization::computeRelativeChangeFromPreviousTimestep(), MathLib::KelvinVector::Invariants< KelvinVectorSize >::equivalentStress(), GeoLib::Polyline::getDistanceAlongPolyline(), GeoLib::Grid< POINT >::getNearestPoint(), MathLib::KelvinVector::Invariants< KelvinVectorSize >::J2(), MathLib::KelvinVector::Invariants< KelvinVectorSize >::J3(), norm1(), norm2(), normMax(), and GeoLib::PointVec::resetInternalDataStructures().
double MathLib::LinAlg::norm1 | ( | MatrixOrVector const & | x | ) |
Computes the Manhattan norm of x
.
double MathLib::LinAlg::norm1 | ( | PETScVector const & | x | ) |
Definition at line 82 of file LinAlg.cpp.
References MathLib::PETScVector::getRawVector(), and norm().
Referenced by norm().
double MathLib::LinAlg::norm2 | ( | MatrixOrVector const & | x | ) |
Computes the Euclidean norm of x
.
double MathLib::LinAlg::norm2 | ( | PETScVector const & | x | ) |
Definition at line 92 of file LinAlg.cpp.
References MathLib::PETScVector::getRawVector(), and norm().
Referenced by norm().
double MathLib::LinAlg::normMax | ( | MatrixOrVector const & | x | ) |
Computes the maximum norm of x
.
double MathLib::LinAlg::normMax | ( | PETScVector const & | x | ) |
Definition at line 102 of file LinAlg.cpp.
References MathLib::PETScVector::getRawVector(), and norm().
Referenced by norm().
void MathLib::LinAlg::scale | ( | MatrixOrVector & | x, |
double const | a | ||
) |
void MathLib::LinAlg::scale | ( | PETScMatrix & | A, |
double const | a | ||
) |
Definition at line 117 of file LinAlg.cpp.
References MathLib::PETScMatrix::getRawMatrix().
void MathLib::LinAlg::scale | ( | PETScVector & | x, |
double const | a | ||
) |
Definition at line 44 of file LinAlg.cpp.
References MathLib::PETScVector::getRawVector().
Referenced by ProcessLib::computeResiduum(), MaterialPropertyLib::Property::description(), MaterialPropertyLib::ClausiusClapeyron::dMolarMass(), MaterialPropertyLib::RelPermBrooksCorey::dValue(), MaterialPropertyLib::RelPermBrooksCoreyNonwettingPhase::dValue(), MaterialPropertyLib::RelPermLiakopoulos::dValue(), NumLib::BackwardEuler::getWeightedOldX(), VtkVisImageItem::Initialize(), main(), MaterialPropertyLib::ClausiusClapeyron::molarMass(), VtkVisTabWidget::on_scaleZ_textChanged(), ProcessLib::PhaseField::PhaseFieldProcess< DisplacementDim >::postNonLinearSolverConcreteProcess(), VtkVisTabWidget::setActiveItem(), VtkTextureOnSurfaceFilter::SetRaster(), VtkBGImageSource::SetRaster(), MaterialPropertyLib::Property::setScale(), VtkColorByHeightFilter::SetTableRangeScaling(), VtkVisPointSetItem::setTranslation(), MaterialPropertyLib::RelPermBrooksCorey::value(), MaterialPropertyLib::RelPermBrooksCoreyNonwettingPhase::value(), and MaterialPropertyLib::RelPermLiakopoulos::value().
void MathLib::LinAlg::set | ( | PETScVector & | x, |
double const | a | ||
) |
Definition at line 32 of file LinAlg.cpp.
References MathLib::PETScVector::getRawVector().
Referenced by ProcessLib::ComponentTransport::intersection(), and ProcessLib::SmallDeformation::writeMaterialForces().
void MathLib::LinAlg::setLocalAccessibleVector | ( | PETScVector const & | x | ) |
Set local accessible vector in order to get entries. Call this before call operator[] or get(...) of x.
Definition at line 27 of file LinAlg.cpp.
References MathLib::PETScVector::setLocalAccessibleVector().
Referenced by ProcessLib::CoupledSolutionsForStaggeredScheme::CoupledSolutionsForStaggeredScheme(), ProcessLib::Process::assemble(), ProcessLib::Process::assembleWithJacobian(), NumLib::LocalLinearLeastSquaresExtrapolator::calculateResidualElement(), ProcessLib::Process::computeSecondaryVariable(), ProcessLib::Process::postIteration(), ProcessLib::TES::TESProcess::postIterationConcreteProcess(), ProcessLib::Process::postNonLinearSolver(), ProcessLib::Process::postTimestep(), ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< GlobalDim >::postTimestepConcreteProcess(), ProcessLib::Process::preIteration(), ProcessLib::Process::preTimestep(), ProcessLib::Process::setInitialConditions(), NumLib::transformVariableFromGlobalVector(), ProcessLib::PhaseField::PhaseFieldProcess< DisplacementDim >::updateConstraints(), and ProcessLib::SmallDeformation::writeMaterialForces().