OGS
ProcessLib::VectorMatrixAssembler Class Referencefinal

Detailed Description

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>

Collaboration diagram for ProcessLib::VectorMatrixAssembler:
[legend]

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
 

Constructor & Destructor Documentation

◆ VectorMatrixAssembler()

ProcessLib::VectorMatrixAssembler::VectorMatrixAssembler ( AbstractJacobianAssembler & jacobian_assembler)
explicit

Definition at line 22 of file VectorMatrixAssembler.cpp.

24 : _jacobian_assembler(jacobian_assembler)
25{
26}
AbstractJacobianAssembler & _jacobian_assembler
Used to assemble the Jacobian.

Member Function Documentation

◆ assemble()

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.

Remarks
Jacobian is not assembled here, see assembleWithJacobian().

Definition at line 39 of file VectorMatrixAssembler.cpp.

45{
46 std::vector<std::vector<GlobalIndexType>> indices_of_processes;
47 indices_of_processes.reserve(dof_tables.size());
48 transform(cbegin(dof_tables), cend(dof_tables),
49 back_inserter(indices_of_processes),
50 [&](auto const dof_table)
51 { return NumLib::getIndices(mesh_item_id, *dof_table); });
52
53 auto const& indices = indices_of_processes[process_id];
54 _local_M_data.clear();
55 _local_K_data.clear();
56 _local_b_data.clear();
57
58 std::size_t const number_of_processes = x.size();
59 // Monolithic scheme
60 if (number_of_processes == 1)
61 {
62 auto const local_x = x[process_id]->get(indices);
63 auto const local_x_prev = x_prev[process_id]->get(indices);
64 local_assembler.assemble(t, dt, local_x, local_x_prev, _local_M_data,
66 }
67 else // Staggered scheme
68 {
69 auto local_coupled_xs =
70 getCoupledLocalSolutions(x, indices_of_processes);
71 auto const local_x = MathLib::toVector(local_coupled_xs);
72
73 auto local_coupled_x_prevs =
74 getCoupledLocalSolutions(x_prev, indices_of_processes);
75 auto const local_x_prev = MathLib::toVector(local_coupled_x_prevs);
76
77 local_assembler.assembleForStaggeredScheme(
78 t, dt, local_x, local_x_prev, process_id, _local_M_data,
80 }
81
82 auto const num_r_c = indices.size();
83 auto const r_c_indices =
85
86 if (M && !_local_M_data.empty())
87 {
88 auto const local_M = MathLib::toMatrix(_local_M_data, num_r_c, num_r_c);
89 M->add(r_c_indices, local_M);
90 }
91 if (K && !_local_K_data.empty())
92 {
93 auto const local_K = MathLib::toMatrix(_local_K_data, num_r_c, num_r_c);
94 K->add(r_c_indices, local_K);
95 }
96 if (b && !_local_b_data.empty())
97 {
98 assert(_local_b_data.size() == num_r_c);
99 b->add(indices, _local_b_data);
100 }
101
102 _local_output(t, process_id, mesh_item_id, _local_M_data, _local_K_data,
104}
int add(IndexType row, IndexType col, double val)
Definition EigenMatrix.h:87
void add(IndexType rowId, double v)
add entry
Definition EigenVector.h:76
MathLib::RowColumnIndices< GlobalIndexType > RowColumnIndices
Assembly::LocalMatrixOutput _local_output
Eigen::Map< const Vector > toVector(std::vector< double > const &data, Eigen::VectorXd::Index size)
Creates an Eigen mapped vector from the given data vector.
Eigen::Map< const Matrix > toMatrix(std::vector< double > const &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)
std::vector< GlobalIndexType > getIndices(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
std::vector< double > getCoupledLocalSolutions(std::vector< GlobalVector * > const &global_solutions, std::vector< std::vector< GlobalIndexType > > const &indices)

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().

◆ assembleWithJacobian()

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.

Note
The Jacobian must be assembled.

Definition at line 106 of file VectorMatrixAssembler.cpp.

112{
113 std::vector<std::vector<GlobalIndexType>> indices_of_processes;
114 indices_of_processes.reserve(dof_tables.size());
115 transform(cbegin(dof_tables), cend(dof_tables),
116 back_inserter(indices_of_processes),
117 [&](auto const dof_table)
118 { return NumLib::getIndices(mesh_item_id, *dof_table); });
119
120 auto const& indices = indices_of_processes[process_id];
121
122 _local_b_data.clear();
123 _local_Jac_data.clear();
124
125 std::size_t const number_of_processes = x.size();
126 // Monolithic scheme
127 if (number_of_processes == 1)
128 {
129 auto const local_x = x[process_id]->get(indices);
130 auto const local_x_prev = x_prev[process_id]->get(indices);
132 local_assembler, t, dt, local_x, local_x_prev, _local_b_data,
134 }
135 else // Staggered scheme
136 {
137 auto local_coupled_xs =
138 getCoupledLocalSolutions(x, indices_of_processes);
139 auto const local_x = MathLib::toVector(local_coupled_xs);
140
141 auto local_coupled_x_prevs =
142 getCoupledLocalSolutions(x_prev, indices_of_processes);
143 auto const local_x_prev = MathLib::toVector(local_coupled_x_prevs);
144
146 local_assembler, t, dt, local_x, local_x_prev, process_id,
148 }
149
150 auto const num_r_c = indices.size();
151 auto const r_c_indices =
153
154 if (b && !_local_b_data.empty())
155 {
156 assert(_local_b_data.size() == num_r_c);
157 b->add(indices, _local_b_data);
158 }
159 if (Jac && !_local_Jac_data.empty())
160 {
161 auto const local_Jac =
162 MathLib::toMatrix(_local_Jac_data, num_r_c, num_r_c);
163 Jac->add(r_c_indices, local_Jac);
164 }
165 else
166 {
167 OGS_FATAL(
168 "No Jacobian has been assembled! This might be due to programming "
169 "errors in the local assembler of the current process.");
170 }
171
172 _local_output(t, process_id, mesh_item_id, _local_b_data, _local_Jac_data);
173}
#define OGS_FATAL(...)
Definition Error.h:26
virtual void assembleWithJacobian(LocalAssemblerInterface &local_assembler, double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &local_x_prev, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data)=0
virtual void assembleWithJacobianForStaggeredScheme(LocalAssemblerInterface &, double const, double const, Eigen::VectorXd const &, Eigen::VectorXd const &, int const, std::vector< double > &, std::vector< double > &)

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().

◆ preAssemble()

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.

32{
33 auto const indices = NumLib::getIndices(mesh_item_id, dof_table);
34 auto const local_x = x.get(indices);
35
36 local_assembler.preAssemble(t, dt, local_x);
37}
double get(IndexType rowId) const
get entry
Definition EigenVector.h:58

References MathLib::EigenVector::get(), NumLib::getIndices(), and ProcessLib::LocalAssemblerInterface::preAssemble().

Referenced by ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalProcess< DisplacementDim >::preAssembleConcreteProcess().

Member Data Documentation

◆ _jacobian_assembler

AbstractJacobianAssembler& ProcessLib::VectorMatrixAssembler::_jacobian_assembler
private

Used to assemble the Jacobian.

Definition at line 72 of file VectorMatrixAssembler.h.

Referenced by assembleWithJacobian().

◆ _local_b_data

std::vector<double> ProcessLib::VectorMatrixAssembler::_local_b_data
private

Definition at line 68 of file VectorMatrixAssembler.h.

Referenced by assemble(), and assembleWithJacobian().

◆ _local_Jac_data

std::vector<double> ProcessLib::VectorMatrixAssembler::_local_Jac_data
private

Definition at line 69 of file VectorMatrixAssembler.h.

Referenced by assembleWithJacobian().

◆ _local_K_data

std::vector<double> ProcessLib::VectorMatrixAssembler::_local_K_data
private

Definition at line 67 of file VectorMatrixAssembler.h.

Referenced by assemble().

◆ _local_M_data

std::vector<double> ProcessLib::VectorMatrixAssembler::_local_M_data
private

Definition at line 66 of file VectorMatrixAssembler.h.

Referenced by assemble().

◆ _local_output

Assembly::LocalMatrixOutput ProcessLib::VectorMatrixAssembler::_local_output
private

Definition at line 74 of file VectorMatrixAssembler.h.

Referenced by assemble(), and assembleWithJacobian().


The documentation for this class was generated from the following files: