OGS
ProcessLib::HeatConduction::HeatConductionProcess Class Referencefinal

Detailed Description

Definition at line 22 of file HeatConductionProcess.h.

#include <HeatConductionProcess.h>

Inheritance diagram for ProcessLib::HeatConduction::HeatConductionProcess:
[legend]
Collaboration diagram for ProcessLib::HeatConduction::HeatConductionProcess:
[legend]

Public Member Functions

 HeatConductionProcess (std::string name, MeshLib::Mesh &mesh, std::unique_ptr< ProcessLib::AbstractJacobianAssembler > &&jacobian_assembler, std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const &parameters, unsigned const integration_order, std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > &&process_variables, HeatConductionProcessData &&process_data, SecondaryVariableCollection &&secondary_variables, bool const is_linear, bool const ls_compute_only_upon_timestep_change)
 
void computeSecondaryVariableConcrete (double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_prev, int const process_id) override
 
bool shouldLinearSolverComputeOnlyUponTimestepChange () const override
 
ODESystem interface
bool isLinear () const override
 
- Public Member Functions inherited from ProcessLib::Process
 Process (std::string name_, MeshLib::Mesh &mesh, std::unique_ptr< AbstractJacobianAssembler > &&jacobian_assembler, std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const &parameters, unsigned const integration_order, std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > &&process_variables, SecondaryVariableCollection &&secondary_variables, const bool use_monolithic_scheme=true)
 
void preTimestep (std::vector< GlobalVector * > const &x, const double t, const double delta_t, const int process_id)
 Preprocessing before starting assembly for new timestep.
 
void postTimestep (std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, const double t, const double delta_t, int const process_id)
 Postprocessing after a complete timestep.
 
void postNonLinearSolver (std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, const double t, double const dt, int const process_id)
 
void preIteration (const unsigned iter, GlobalVector const &x) final
 
void computeSecondaryVariable (double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_prev, int const process_id)
 compute secondary variables for the coupled equations or for output.
 
NumLib::IterationResult postIteration (GlobalVector const &x) final
 
void initialize (std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
 
void setInitialConditions (std::vector< GlobalVector * > &process_solutions, std::vector< GlobalVector * > const &process_solutions_prev, double const t, int const process_id)
 
MathLib::MatrixSpecifications getMatrixSpecifications (const int process_id) const override
 
void updateDeactivatedSubdomains (double const time, const int process_id)
 
virtual bool isMonolithicSchemeUsed () const
 
virtual void extrapolateIntegrationPointValuesToNodes (const double, std::vector< GlobalVector * > const &, std::vector< GlobalVector * > &)
 
void preAssemble (const double t, double const dt, GlobalVector const &x) final
 
void assemble (const double 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) final
 
void assembleWithJacobian (const double 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, GlobalMatrix &Jac) final
 
void preOutput (const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id)
 
std::vector< NumLib::IndexValueVector< GlobalIndexType > > const * getKnownSolutions (double const t, GlobalVector const &x, int const process_id) const final
 
virtual NumLib::LocalToGlobalIndexMap const & getDOFTable (const int) const
 
MeshLib::MeshgetMesh () const
 
std::vector< std::reference_wrapper< ProcessVariable > > const & getProcessVariables (const int process_id) const
 
std::vector< std::size_t > const & getActiveElementIDs () const
 
SecondaryVariableCollection const & getSecondaryVariables () const
 
std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > const & getIntegrationPointWriters () const
 
virtual Eigen::Vector3d getFlux (std::size_t, MathLib::Point3d const &, double const, std::vector< GlobalVector * > const &) const
 
virtual void solveReactionEquation (std::vector< GlobalVector * > &, std::vector< GlobalVector * > const &, double const, double const, NumLib::EquationSystem &, int const)
 
- Public Member Functions inherited from ProcessLib::SubmeshAssemblySupport
virtual std::vector< std::string > initializeAssemblyOnSubmeshes (std::vector< std::reference_wrapper< MeshLib::Mesh > > const &meshes)
 
virtual ~SubmeshAssemblySupport ()=default
 

Private Member Functions

void initializeConcreteProcess (NumLib::LocalToGlobalIndexMap const &dof_table, MeshLib::Mesh const &mesh, unsigned const integration_order) override
 Process specific initialization called by initialize().
 
void assembleConcreteProcess (const double t, double const, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b) override
 
void assembleWithJacobianConcreteProcess (const double t, double const, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b, GlobalMatrix &Jac) override
 
void preOutputConcreteProcess (const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id) override
 

Private Attributes

HeatConductionProcessData _process_data
 
std::vector< std::unique_ptr< HeatConductionLocalAssemblerInterface > > _local_assemblers
 
MeshLib::PropertyVector< double > * _heat_flux = nullptr
 
AssembledMatrixCache _asm_mat_cache
 
bool const _ls_compute_only_upon_timestep_change
 

Additional Inherited Members

- Public Attributes inherited from ProcessLib::Process
std::string const name
 
- Static Public Attributes inherited from ProcessLib::Process
static PROCESSLIB_EXPORT const std::string constant_one_parameter_name = "constant_one"
 
- Protected Member Functions inherited from ProcessLib::Process
std::vector< NumLib::LocalToGlobalIndexMap const * > getDOFTables (int const number_of_processes) const
 
NumLib::ExtrapolatorgetExtrapolator () const
 
NumLib::LocalToGlobalIndexMap const & getSingleComponentDOFTable () const
 
void initializeProcessBoundaryConditionsAndSourceTerms (const NumLib::LocalToGlobalIndexMap &dof_table, const int process_id, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
 
virtual void constructDofTable ()
 
void constructMonolithicProcessDofTable ()
 
void constructDofTableOfSpecifiedProcessStaggeredScheme (const int specified_process_id)
 
virtual std::tuple< NumLib::LocalToGlobalIndexMap *, bool > getDOFTableForExtrapolatorData () const
 
std::vector< GlobalIndexTypegetIndicesOfResiduumWithoutInitialCompensation () const override
 
- Protected Attributes inherited from ProcessLib::Process
MeshLib::Mesh_mesh
 
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_all_nodes
 
std::unique_ptr< NumLib::LocalToGlobalIndexMap_local_to_global_index_map
 
SecondaryVariableCollection _secondary_variables
 
std::unique_ptr< ProcessLib::AbstractJacobianAssembler_jacobian_assembler
 
VectorMatrixAssembler _global_assembler
 
const bool _use_monolithic_scheme
 
unsigned const _integration_order
 
std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > _integration_point_writer
 
GlobalSparsityPattern _sparsity_pattern
 
std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > _process_variables
 
std::vector< BoundaryConditionCollection_boundary_conditions
 

Constructor & Destructor Documentation

◆ HeatConductionProcess()

ProcessLib::HeatConduction::HeatConductionProcess::HeatConductionProcess ( std::string name,
MeshLib::Mesh & mesh,
std::unique_ptr< ProcessLib::AbstractJacobianAssembler > && jacobian_assembler,
std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const & parameters,
unsigned const integration_order,
std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > && process_variables,
HeatConductionProcessData && process_data,
SecondaryVariableCollection && secondary_variables,
bool const is_linear,
bool const ls_compute_only_upon_timestep_change )

Definition at line 26 of file HeatConductionProcess.cpp.

38 : Process(std::move(name), mesh, std::move(jacobian_assembler), parameters,
39 integration_order, std::move(process_variables),
40 std::move(secondary_variables)),
41 _process_data(std::move(process_data)),
42 _asm_mat_cache{is_linear, true /*use_monolithic_scheme*/},
44 ls_compute_only_upon_timestep_change}
45{
46 if (ls_compute_only_upon_timestep_change)
47 {
48 // TODO move this feature to some common location for all processes.
49 if (!is_linear)
50 {
52 "Using the linear solver compute() method only upon timestep "
53 "change only makes sense for linear model equations.");
54 }
55
56 WARN(
57 "You specified that the HeatConduction linear solver will do "
58 "the compute() step only upon timestep change. This is an expert "
59 "option. It is your responsibility to ensure that "
60 "the conditions for the correct use of this feature are met! "
61 "Otherwise OGS might compute garbage without being recognized. "
62 "There is no safety net!");
63
64 WARN(
65 "You specified that the HeatConduction linear solver will do "
66 "the compute() step only upon timestep change. This option will "
67 "only be used by the Picard non-linear solver. The Newton-Raphson "
68 "non-linear solver will silently ignore this setting, i.e., it "
69 "won't do any harm, there, but you won't observe the speedup you "
70 "probably expect.");
71 }
72
73 _heat_flux = MeshLib::getOrCreateMeshProperty<double>(
74 mesh, "HeatFlowRate", MeshLib::MeshItemType::Node, 1);
75}
#define OGS_FATAL(...)
Definition Error.h:26
void WARN(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:40
std::string const name
Definition Process.h:354
Process(std::string name_, MeshLib::Mesh &mesh, std::unique_ptr< AbstractJacobianAssembler > &&jacobian_assembler, std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const &parameters, unsigned const integration_order, std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > &&process_variables, SecondaryVariableCollection &&secondary_variables, const bool use_monolithic_scheme=true)
Definition Process.cpp:44

References _heat_flux, MeshLib::Node, OGS_FATAL, and WARN().

Member Function Documentation

◆ assembleConcreteProcess()

void ProcessLib::HeatConduction::HeatConductionProcess::assembleConcreteProcess ( const double 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 )
overrideprivatevirtual

Implements ProcessLib::Process.

Definition at line 96 of file HeatConductionProcess.cpp.

100{
101 DBUG("Assemble HeatConductionProcess.");
102
103 ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
104
105 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_table = {
107
108 _asm_mat_cache.assemble(t, dt, x, x_prev, process_id, M, K, b, dof_table,
111}
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:30
std::vector< std::unique_ptr< HeatConductionLocalAssemblerInterface > > _local_assemblers
std::vector< std::size_t > const & getActiveElementIDs() const
std::vector< std::reference_wrapper< ProcessVariable > > const & getProcessVariables(const int process_id) const
Definition Process.h:155
VectorMatrixAssembler _global_assembler
Definition Process.h:367
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map
Definition Process.h:360
void assemble(const double 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, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, VectorMatrixAssembler &global_assembler, VectorOfLocalAssemblers const &local_assemblers, std::vector< std::size_t > const &active_element_ids)

References _asm_mat_cache, ProcessLib::Process::_global_assembler, _local_assemblers, ProcessLib::Process::_local_to_global_index_map, ProcessLib::AssembledMatrixCache::assemble(), DBUG(), ProcessLib::ProcessVariable::getActiveElementIDs(), and ProcessLib::Process::getProcessVariables().

Referenced by preOutputConcreteProcess().

◆ assembleWithJacobianConcreteProcess()

void ProcessLib::HeatConduction::HeatConductionProcess::assembleWithJacobianConcreteProcess ( const double 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,
GlobalMatrix & Jac )
overrideprivatevirtual

Implements ProcessLib::Process.

Definition at line 113 of file HeatConductionProcess.cpp.

117{
118 DBUG("AssembleWithJacobian HeatConductionProcess.");
119
120 ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
121
122 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_table = {
124 // Call global assembler for each local assembly item.
127 _local_assemblers, pv.getActiveElementIDs(), dof_table, t, dt, x,
128 x_prev, process_id, M, K, b, Jac);
129
130 transformVariableFromGlobalVector(b, 0 /*variable id*/,
132 std::negate<double>());
133}
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, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b, GlobalMatrix &Jac)
void transformVariableFromGlobalVector(GlobalVector const &input_vector, int const variable_id, NumLib::LocalToGlobalIndexMap const &local_to_global_index_map, MeshLib::PropertyVector< double > &output_vector, Functor map_function)
static void executeSelectedMemberDereferenced(Object &object, Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)

References ProcessLib::Process::_global_assembler, _heat_flux, _local_assemblers, ProcessLib::Process::_local_to_global_index_map, ProcessLib::VectorMatrixAssembler::assembleWithJacobian(), DBUG(), NumLib::SerialExecutor::executeSelectedMemberDereferenced(), ProcessLib::ProcessVariable::getActiveElementIDs(), and ProcessLib::Process::getProcessVariables().

◆ computeSecondaryVariableConcrete()

void ProcessLib::HeatConduction::HeatConductionProcess::computeSecondaryVariableConcrete ( double const t,
double const dt,
std::vector< GlobalVector * > const & x,
GlobalVector const & x_prev,
int const process_id )
overridevirtual

Reimplemented from ProcessLib::Process.

Definition at line 135 of file HeatConductionProcess.cpp.

138{
139 DBUG("Compute heat flux for HeatConductionProcess.");
140
141 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
142 dof_tables.reserve(x.size());
143 std::generate_n(std::back_inserter(dof_tables), x.size(),
144 [&]() { return _local_to_global_index_map.get(); });
145
146 ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
149 _local_assemblers, pv.getActiveElementIDs(), dof_tables, t, dt, x,
150 x_prev, process_id);
151}
virtual void computeSecondaryVariable(std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_prev, int const process_id)
static void executeSelectedMemberOnDereferenced(Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)

References _local_assemblers, ProcessLib::LocalAssemblerInterface::computeSecondaryVariable(), DBUG(), NumLib::SerialExecutor::executeSelectedMemberOnDereferenced(), ProcessLib::ProcessVariable::getActiveElementIDs(), and ProcessLib::Process::getProcessVariables().

◆ initializeConcreteProcess()

void ProcessLib::HeatConduction::HeatConductionProcess::initializeConcreteProcess ( NumLib::LocalToGlobalIndexMap const & dof_table,
MeshLib::Mesh const & mesh,
unsigned const integration_order )
overrideprivatevirtual

Process specific initialization called by initialize().

Implements ProcessLib::Process.

Definition at line 77 of file HeatConductionProcess.cpp.

81{
82 int const mesh_space_dimension = _process_data.mesh_space_dimension;
83
84 ProcessLib::createLocalAssemblers<LocalAssemblerData>(
85 mesh_space_dimension, mesh.getElements(), dof_table, _local_assemblers,
86 NumLib::IntegrationOrder{integration_order}, mesh.isAxiallySymmetric(),
88
90 "heat_flux",
92 mesh_space_dimension, getExtrapolator(), _local_assemblers,
94}
virtual std::vector< double > const & getIntPtHeatFlux(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
SecondaryVariableCollection _secondary_variables
Definition Process.h:362
NumLib::Extrapolator & getExtrapolator() const
Definition Process.h:199
void addSecondaryVariable(std::string const &internal_name, SecondaryVariableFunctions &&fcts)
SecondaryVariableFunctions makeExtrapolator(const unsigned num_components, NumLib::Extrapolator &extrapolator, LocalAssemblerCollection const &local_assemblers, typename NumLib::ExtrapolatableLocalAssemblerCollection< LocalAssemblerCollection >::IntegrationPointValuesMethod integration_point_values_method)

References _local_assemblers, _process_data, ProcessLib::Process::_secondary_variables, ProcessLib::SecondaryVariableCollection::addSecondaryVariable(), MeshLib::Mesh::getElements(), ProcessLib::Process::getExtrapolator(), ProcessLib::HeatConduction::HeatConductionLocalAssemblerInterface::getIntPtHeatFlux(), MeshLib::Mesh::isAxiallySymmetric(), ProcessLib::makeExtrapolator(), and ProcessLib::HeatConduction::HeatConductionProcessData::mesh_space_dimension.

◆ isLinear()

bool ProcessLib::HeatConduction::HeatConductionProcess::isLinear ( ) const
inlineoverride

◆ preOutputConcreteProcess()

void ProcessLib::HeatConduction::HeatConductionProcess::preOutputConcreteProcess ( const double t,
double const dt,
std::vector< GlobalVector * > const & x,
std::vector< GlobalVector * > const & x_prev,
int const process_id )
overrideprivatevirtual

Reimplemented from ProcessLib::Process.

Definition at line 153 of file HeatConductionProcess.cpp.

159{
160 auto const matrix_specification = getMatrixSpecifications(process_id);
161
163 matrix_specification);
165 matrix_specification);
167 matrix_specification);
168
169 M->setZero();
170 K->setZero();
171 b->setZero();
172
173 assembleConcreteProcess(t, dt, x, x_prev, process_id, *M, *K, *b);
174
175 BaseLib::RunTime time_residuum;
176 time_residuum.start();
177
178 auto const residuum = computeResiduum(dt, *x[0], *x_prev[0], *M, *K, *b);
179
181 *_heat_flux, std::negate<double>());
182
183 INFO("[time] Computing residuum flow rates took {:g} s",
184 time_residuum.elapsed());
185}
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:35
Count the running time.
Definition RunTime.h:29
double elapsed() const
Get the elapsed time in seconds.
Definition RunTime.h:42
void start()
Start the timer.
Definition RunTime.h:32
void assembleConcreteProcess(const double t, double const, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b) override
MathLib::MatrixSpecifications getMatrixSpecifications(const int process_id) const override
Definition Process.cpp:210
GlobalVector computeResiduum(double const dt, GlobalVector const &x, GlobalVector const &x_prev, GlobalMatrix const &M, GlobalMatrix const &K, GlobalVector const &b)

References _heat_flux, ProcessLib::Process::_local_to_global_index_map, assembleConcreteProcess(), ProcessLib::computeResiduum(), BaseLib::RunTime::elapsed(), ProcessLib::Process::getMatrixSpecifications(), INFO(), and BaseLib::RunTime::start().

◆ shouldLinearSolverComputeOnlyUponTimestepChange()

bool ProcessLib::HeatConduction::HeatConductionProcess::shouldLinearSolverComputeOnlyUponTimestepChange ( ) const
inlineoverride

Definition at line 51 of file HeatConductionProcess.h.

52 {
53 // TODO move this feature to some common location for all processes.
55 }

References _ls_compute_only_upon_timestep_change.

Member Data Documentation

◆ _asm_mat_cache

AssembledMatrixCache ProcessLib::HeatConduction::HeatConductionProcess::_asm_mat_cache
private

Definition at line 88 of file HeatConductionProcess.h.

Referenced by assembleConcreteProcess(), and isLinear().

◆ _heat_flux

MeshLib::PropertyVector<double>* ProcessLib::HeatConduction::HeatConductionProcess::_heat_flux = nullptr
private

◆ _local_assemblers

std::vector<std::unique_ptr<HeatConductionLocalAssemblerInterface> > ProcessLib::HeatConduction::HeatConductionProcess::_local_assemblers
private

◆ _ls_compute_only_upon_timestep_change

bool const ProcessLib::HeatConduction::HeatConductionProcess::_ls_compute_only_upon_timestep_change
private

◆ _process_data

HeatConductionProcessData ProcessLib::HeatConduction::HeatConductionProcess::_process_data
private

Definition at line 81 of file HeatConductionProcess.h.

Referenced by initializeConcreteProcess().


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