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 25 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 15 of file VectorMatrixAssembler.cpp.

17 : _jacobian_assembler(jacobian_assembler)
18{
19}
AbstractJacobianAssembler & _jacobian_assembler
Used to assemble the Jacobian.

References _jacobian_assembler.

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 32 of file VectorMatrixAssembler.cpp.

38{
39 std::vector<std::vector<GlobalIndexType>> indices_of_processes;
40 indices_of_processes.reserve(dof_tables.size());
41 transform(cbegin(dof_tables), cend(dof_tables),
42 back_inserter(indices_of_processes),
43 [&](auto const dof_table)
44 { return NumLib::getIndices(mesh_item_id, *dof_table); });
45
46 auto const& indices = indices_of_processes[process_id];
47 _local_M_data.clear();
48 _local_K_data.clear();
49 _local_b_data.clear();
50
51 std::size_t const number_of_processes = x.size();
52 // Monolithic scheme
53 if (number_of_processes == 1)
54 {
55 auto const local_x = x[process_id]->get(indices);
56 auto const local_x_prev = x_prev[process_id]->get(indices);
57 local_assembler.assemble(t, dt, local_x, local_x_prev, _local_M_data,
59 }
60 else // Staggered scheme
61 {
62 auto local_coupled_xs =
63 getCoupledLocalSolutions(x, indices_of_processes);
64 auto const local_x = MathLib::toVector(local_coupled_xs);
65
66 auto local_coupled_x_prevs =
67 getCoupledLocalSolutions(x_prev, indices_of_processes);
68 auto const local_x_prev = MathLib::toVector(local_coupled_x_prevs);
69
70 local_assembler.assembleForStaggeredScheme(
71 t, dt, local_x, local_x_prev, process_id, _local_M_data,
73 }
74
75 auto const num_r_c = indices.size();
76 auto const r_c_indices =
78
79 if (M && !_local_M_data.empty())
80 {
81 auto const local_M = MathLib::toMatrix(_local_M_data, num_r_c, num_r_c);
82 M->add(r_c_indices, local_M);
83 }
84 if (K && !_local_K_data.empty())
85 {
86 auto const local_K = MathLib::toMatrix(_local_K_data, num_r_c, num_r_c);
87 K->add(r_c_indices, local_K);
88 }
89 if (b && !_local_b_data.empty())
90 {
91 assert(_local_b_data.size() == num_r_c);
92 b->add(indices, _local_b_data);
93 }
94
95 _local_output(t, process_id, mesh_item_id, _local_M_data, _local_K_data,
97}
int add(IndexType row, IndexType col, double val)
Definition EigenMatrix.h:80
void add(IndexType rowId, double v)
add entry
Definition EigenVector.h:70
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::HeatTransportBHE::HeatTransportBHEProcess::assembleConcreteProcess(), ProcessLib::HT::HTProcess::assembleConcreteProcess(), ProcessLib::LargeDeformation::LargeDeformationProcess< DisplacementDim >::assembleConcreteProcess(), ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< DisplacementDim >::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::SteadyStateDiffusion::SteadyStateDiffusion::assembleConcreteProcess(), ProcessLib::ThermalTwoPhaseFlowWithPP::ThermalTwoPhaseFlowWithPPProcess::assembleConcreteProcess(), ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >::assembleConcreteProcess(), ProcessLib::TwoPhaseFlowWithPP::TwoPhaseFlowWithPPProcess::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 99 of file VectorMatrixAssembler.cpp.

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

References _jacobian_assembler, _local_b_data, _local_Jac_data, _local_output, MathLib::EigenMatrix::add(), MathLib::EigenVector::add(), ProcessLib::getCoupledLocalSolutions(), NumLib::getIndices(), OGS_FATAL, MathLib::toMatrix(), and MathLib::toVector().

Referenced by ProcessLib::HMPhaseField::HMPhaseFieldProcess< DisplacementDim >::assembleWithJacobianConcreteProcess(), ProcessLib::HT::HTProcess::assembleWithJacobianConcreteProcess(), ProcessLib::LargeDeformation::LargeDeformationProcess< DisplacementDim >::assembleWithJacobianConcreteProcess(), ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< DisplacementDim >::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::SteadyStateDiffusion::SteadyStateDiffusion::assembleWithJacobianConcreteProcess(), ProcessLib::ThermalTwoPhaseFlowWithPP::ThermalTwoPhaseFlowWithPPProcess::assembleWithJacobianConcreteProcess(), ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >::assembleWithJacobianConcreteProcess(), ProcessLib::TwoPhaseFlowWithPP::TwoPhaseFlowWithPPProcess::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 21 of file VectorMatrixAssembler.cpp.

25{
26 auto const indices = NumLib::getIndices(mesh_item_id, dof_table);
27 auto const local_x = x.get(indices);
28
29 local_assembler.preAssemble(t, dt, local_x);
30}
double get(IndexType rowId) const
get entry
Definition EigenVector.h:52

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

Member Data Documentation

◆ _jacobian_assembler

AbstractJacobianAssembler& ProcessLib::VectorMatrixAssembler::_jacobian_assembler
private

Used to assemble the Jacobian.

Definition at line 65 of file VectorMatrixAssembler.h.

Referenced by VectorMatrixAssembler(), and assembleWithJacobian().

◆ _local_b_data

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

Definition at line 61 of file VectorMatrixAssembler.h.

Referenced by assemble(), and assembleWithJacobian().

◆ _local_Jac_data

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

Definition at line 62 of file VectorMatrixAssembler.h.

Referenced by assembleWithJacobian().

◆ _local_K_data

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

Definition at line 60 of file VectorMatrixAssembler.h.

Referenced by assemble().

◆ _local_M_data

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

Definition at line 59 of file VectorMatrixAssembler.h.

Referenced by assemble().

◆ _local_output

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

Definition at line 67 of file VectorMatrixAssembler.h.

Referenced by assemble(), and assembleWithJacobian().


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