OGS
|
Utility class used to assemble global matrices and vectors.
The methods of this class get the global matrices and vectors as input and pass only local data on to the local assemblers.
Definition at line 32 of file VectorMatrixAssembler.h.
#include <VectorMatrixAssembler.h>
Public Member Functions | |
VectorMatrixAssembler (AbstractJacobianAssembler &jacobian_assembler) | |
void | preAssemble (const std::size_t mesh_item_id, LocalAssemblerInterface &local_assembler, const NumLib::LocalToGlobalIndexMap &dof_table, const double t, double const dt, const GlobalVector &x) |
void | assemble (std::size_t const mesh_item_id, LocalAssemblerInterface &local_assembler, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, double const t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalMatrix *M, GlobalMatrix *K, GlobalVector *b) |
void | assembleWithJacobian (std::size_t const mesh_item_id, LocalAssemblerInterface &local_assembler, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalVector *b, GlobalMatrix *Jac) |
Private Attributes | |
std::vector< double > | _local_M_data |
std::vector< double > | _local_K_data |
std::vector< double > | _local_b_data |
std::vector< double > | _local_Jac_data |
AbstractJacobianAssembler & | _jacobian_assembler |
Used to assemble the Jacobian. | |
Assembly::LocalMatrixOutput | _local_output |
|
explicit |
Definition at line 22 of file VectorMatrixAssembler.cpp.
void ProcessLib::VectorMatrixAssembler::assemble | ( | std::size_t const | mesh_item_id, |
LocalAssemblerInterface & | local_assembler, | ||
std::vector< NumLib::LocalToGlobalIndexMap const * > const & | dof_tables, | ||
double const | t, | ||
double const | dt, | ||
std::vector< GlobalVector * > const & | x, | ||
std::vector< GlobalVector * > const & | x_prev, | ||
int const | process_id, | ||
GlobalMatrix * | M, | ||
GlobalMatrix * | K, | ||
GlobalVector * | b ) |
AssemblesM
, K
, and b
.
Definition at line 39 of file VectorMatrixAssembler.cpp.
References _local_b_data, _local_K_data, _local_M_data, _local_output, MathLib::EigenMatrix::add(), MathLib::EigenVector::add(), ProcessLib::LocalAssemblerInterface::assemble(), ProcessLib::LocalAssemblerInterface::assembleForStaggeredScheme(), ProcessLib::getCoupledLocalSolutions(), NumLib::getIndices(), MathLib::toMatrix(), and MathLib::toVector().
Referenced by ProcessLib::AssembledMatrixCache::assemble(), ProcessLib::HeatTransportBHE::HeatTransportBHEProcess::assembleConcreteProcess(), ProcessLib::HT::HTProcess::assembleConcreteProcess(), ProcessLib::HydroMechanics::HydroMechanicsProcess< DisplacementDim >::assembleConcreteProcess(), ProcessLib::LargeDeformation::LargeDeformationProcess< DisplacementDim >::assembleConcreteProcess(), ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< GlobalDim >::assembleConcreteProcess(), ProcessLib::LIE::SmallDeformation::SmallDeformationProcess< DisplacementDim >::assembleConcreteProcess(), ProcessLib::LiquidFlow::LiquidFlowProcess::assembleConcreteProcess(), ProcessLib::PhaseField::PhaseFieldProcess< DisplacementDim >::assembleConcreteProcess(), ProcessLib::RichardsComponentTransport::RichardsComponentTransportProcess::assembleConcreteProcess(), ProcessLib::RichardsFlow::RichardsFlowProcess::assembleConcreteProcess(), ProcessLib::RichardsMechanics::RichardsMechanicsProcess< DisplacementDim >::assembleConcreteProcess(), ProcessLib::SmallDeformation::SmallDeformationProcess< DisplacementDim >::assembleConcreteProcess(), ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalProcess< DisplacementDim >::assembleConcreteProcess(), ProcessLib::SteadyStateDiffusion::SteadyStateDiffusion::assembleConcreteProcess(), ProcessLib::StokesFlow::StokesFlowProcess< GlobalDim >::assembleConcreteProcess(), ProcessLib::TES::TESProcess::assembleConcreteProcess(), ProcessLib::ThermalTwoPhaseFlowWithPP::ThermalTwoPhaseFlowWithPPProcess::assembleConcreteProcess(), ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsProcess< DisplacementDim >::assembleConcreteProcess(), ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldProcess< DisplacementDim >::assembleConcreteProcess(), ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >::assembleConcreteProcess(), ProcessLib::ThermoRichardsFlow::ThermoRichardsFlowProcess::assembleConcreteProcess(), ProcessLib::TwoPhaseFlowWithPP::TwoPhaseFlowWithPPProcess::assembleConcreteProcess(), ProcessLib::TwoPhaseFlowWithPrho::TwoPhaseFlowWithPrhoProcess::assembleConcreteProcess(), and ProcessLib::WellboreSimulator::WellboreSimulatorProcess::assembleConcreteProcess().
void ProcessLib::VectorMatrixAssembler::assembleWithJacobian | ( | std::size_t const | mesh_item_id, |
LocalAssemblerInterface & | local_assembler, | ||
std::vector< NumLib::LocalToGlobalIndexMap const * > const & | dof_tables, | ||
const double | t, | ||
double const | dt, | ||
std::vector< GlobalVector * > const & | x, | ||
std::vector< GlobalVector * > const & | x_prev, | ||
int const | process_id, | ||
GlobalVector * | b, | ||
GlobalMatrix * | Jac ) |
Assembles b
, and the Jacobian Jac
of the residual.
Definition at line 106 of file VectorMatrixAssembler.cpp.
References _jacobian_assembler, _local_b_data, _local_Jac_data, _local_output, MathLib::EigenMatrix::add(), MathLib::EigenVector::add(), ProcessLib::AbstractJacobianAssembler::assembleWithJacobian(), ProcessLib::AbstractJacobianAssembler::assembleWithJacobianForStaggeredScheme(), ProcessLib::getCoupledLocalSolutions(), NumLib::getIndices(), OGS_FATAL, MathLib::toMatrix(), and MathLib::toVector().
Referenced by ProcessLib::ComponentTransport::ComponentTransportProcess::assembleWithJacobianConcreteProcess(), ProcessLib::HeatConduction::HeatConductionProcess::assembleWithJacobianConcreteProcess(), ProcessLib::HeatTransportBHE::HeatTransportBHEProcess::assembleWithJacobianConcreteProcess(), ProcessLib::HT::HTProcess::assembleWithJacobianConcreteProcess(), ProcessLib::LargeDeformation::LargeDeformationProcess< DisplacementDim >::assembleWithJacobianConcreteProcess(), ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< GlobalDim >::assembleWithJacobianConcreteProcess(), ProcessLib::LIE::SmallDeformation::SmallDeformationProcess< DisplacementDim >::assembleWithJacobianConcreteProcess(), ProcessLib::LiquidFlow::LiquidFlowProcess::assembleWithJacobianConcreteProcess(), ProcessLib::PhaseField::PhaseFieldProcess< DisplacementDim >::assembleWithJacobianConcreteProcess(), ProcessLib::RichardsComponentTransport::RichardsComponentTransportProcess::assembleWithJacobianConcreteProcess(), ProcessLib::RichardsFlow::RichardsFlowProcess::assembleWithJacobianConcreteProcess(), ProcessLib::RichardsMechanics::RichardsMechanicsProcess< DisplacementDim >::assembleWithJacobianConcreteProcess(), ProcessLib::SmallDeformation::SmallDeformationProcess< DisplacementDim >::assembleWithJacobianConcreteProcess(), ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalProcess< DisplacementDim >::assembleWithJacobianConcreteProcess(), ProcessLib::SteadyStateDiffusion::SteadyStateDiffusion::assembleWithJacobianConcreteProcess(), ProcessLib::TES::TESProcess::assembleWithJacobianConcreteProcess(), ProcessLib::ThermalTwoPhaseFlowWithPP::ThermalTwoPhaseFlowWithPPProcess::assembleWithJacobianConcreteProcess(), ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsProcess< DisplacementDim >::assembleWithJacobianConcreteProcess(), ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldProcess< DisplacementDim >::assembleWithJacobianConcreteProcess(), ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >::assembleWithJacobianConcreteProcess(), ProcessLib::TwoPhaseFlowWithPP::TwoPhaseFlowWithPPProcess::assembleWithJacobianConcreteProcess(), ProcessLib::TwoPhaseFlowWithPrho::TwoPhaseFlowWithPrhoProcess::assembleWithJacobianConcreteProcess(), and ProcessLib::WellboreSimulator::WellboreSimulatorProcess::assembleWithJacobianConcreteProcess().
void ProcessLib::VectorMatrixAssembler::preAssemble | ( | const std::size_t | mesh_item_id, |
LocalAssemblerInterface & | local_assembler, | ||
const NumLib::LocalToGlobalIndexMap & | dof_table, | ||
const double | t, | ||
double const | dt, | ||
const GlobalVector & | x ) |
Definition at line 28 of file VectorMatrixAssembler.cpp.
References MathLib::EigenVector::get(), NumLib::getIndices(), and ProcessLib::LocalAssemblerInterface::preAssemble().
Referenced by ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalProcess< DisplacementDim >::preAssembleConcreteProcess().
|
private |
Used to assemble the Jacobian.
Definition at line 72 of file VectorMatrixAssembler.h.
Referenced by assembleWithJacobian().
|
private |
Definition at line 68 of file VectorMatrixAssembler.h.
Referenced by assemble(), and assembleWithJacobian().
|
private |
Definition at line 69 of file VectorMatrixAssembler.h.
Referenced by assembleWithJacobian().
|
private |
Definition at line 67 of file VectorMatrixAssembler.h.
Referenced by assemble().
|
private |
Definition at line 66 of file VectorMatrixAssembler.h.
Referenced by assemble().
|
private |
Definition at line 74 of file VectorMatrixAssembler.h.
Referenced by assemble(), and assembleWithJacobian().