OGS
HeatConductionProcess.cpp
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) OpenGeoSys Community (opengeosys.org)
2// SPDX-License-Identifier: BSD-3-Clause
3
5
6#include <cassert>
7
14
15namespace ProcessLib
16{
17namespace HeatConduction
18{
20 std::string name,
21 MeshLib::Mesh& mesh,
22 std::unique_ptr<ProcessLib::AbstractJacobianAssembler>&& jacobian_assembler,
23 std::vector<std::unique_ptr<ParameterLib::ParameterBase>> const& parameters,
24 unsigned const integration_order,
25 std::vector<std::vector<std::reference_wrapper<ProcessVariable>>>&&
26 process_variables,
27 HeatConductionProcessData&& process_data,
28 SecondaryVariableCollection&& secondary_variables,
29 bool const is_linear,
30 bool const ls_compute_only_upon_timestep_change)
31 : Process(std::move(name), mesh, std::move(jacobian_assembler), parameters,
32 integration_order, std::move(process_variables),
33 std::move(secondary_variables)),
35 true /*use_monolithic_scheme*/},
36 _process_data(std::move(process_data)),
37 _ls_compute_only_upon_timestep_change{
38 ls_compute_only_upon_timestep_change}
39{
40 // For numerical Jacobian assembler
41 this->_jacobian_assembler->setNonDeformationComponentIDs(
42 {0} /* only one variable: temperature */);
43
44 if (ls_compute_only_upon_timestep_change)
45 {
46 // TODO move this feature to some common location for all processes.
47 if (!is_linear)
48 {
50 "Using the linear solver compute() method only upon timestep "
51 "change only makes sense for linear model equations.");
52 }
53
54 WARN(
55 "You specified that the HeatConduction linear solver will do "
56 "the compute() step only upon timestep change. This is an expert "
57 "option. It is your responsibility to ensure that "
58 "the conditions for the correct use of this feature are met! "
59 "Otherwise OGS might compute garbage without being recognized. "
60 "There is no safety net!");
61
62 WARN(
63 "You specified that the HeatConduction linear solver will do "
64 "the compute() step only upon timestep change. This option will "
65 "only be used by the Picard non-linear solver. The Newton-Raphson "
66 "non-linear solver will silently ignore this setting, i.e., it "
67 "won't do any harm, there, but you won't observe the speedup you "
68 "probably expect.");
69 }
70}
71
73 NumLib::LocalToGlobalIndexMap const& dof_table,
74 MeshLib::Mesh const& mesh,
75 unsigned const integration_order)
76{
77 int const mesh_space_dimension = _process_data.mesh_space_dimension;
78
80 mesh_space_dimension, mesh.getElements(), dof_table, local_assemblers_,
81 NumLib::IntegrationOrder{integration_order}, mesh.isAxiallySymmetric(),
83
84 _secondary_variables.addSecondaryVariable(
85 "heat_flux",
87 mesh_space_dimension, getExtrapolator(), local_assemblers_,
89}
90
91std::vector<std::vector<std::string>>
93 std::vector<std::reference_wrapper<MeshLib::Mesh>> const& meshes)
94{
95 DBUG("HeatConductionProcess initializeSubmeshOutput().");
96
97 std::vector<std::vector<std::string>> residuum_names{{"HeatFlowRate"}};
98
100 meshes, residuum_names);
101
102 return residuum_names;
103}
104
106 const double t, double const dt, std::vector<GlobalVector*> const& x,
107 std::vector<GlobalVector*> const& x_prev, int const process_id,
109{
110 DBUG("Assemble HeatConductionProcess.");
111
112 AssemblyMixin<HeatConductionProcess>::assemble(t, dt, x, x_prev, process_id,
113 M, K, b);
114}
115
117 const double t, double const dt, std::vector<GlobalVector*> const& x,
118 std::vector<GlobalVector*> const& x_prev, int const process_id,
119 GlobalVector& b, GlobalMatrix& Jac)
120{
121 DBUG("AssembleWithJacobian HeatConductionProcess.");
122
124 t, dt, x, x_prev, process_id, b, Jac);
125}
126
128 double const t, double const dt, std::vector<GlobalVector*> const& x,
129 GlobalVector const& x_prev, int const process_id)
130{
131 DBUG("Compute heat flux for HeatConductionProcess.");
132
133 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
134 dof_tables.reserve(x.size());
135 std::generate_n(std::back_inserter(dof_tables), x.size(),
136 [&]() { return _local_to_global_index_map.get(); });
137
140 local_assemblers_, getActiveElementIDs(), dof_tables, t, dt, x, x_prev,
141 process_id);
142}
143
145 std::vector<GlobalVector*> const& /*x*/,
146 const double /*t*/,
147 const double /*dt*/,
148 const int /*process_id*/)
149{
151}
152
154 const double t,
155 double const dt,
156 std::vector<GlobalVector*> const& x,
157 std::vector<GlobalVector*> const& x_prev,
158 int const process_id)
159{
161 process_id);
162}
163} // namespace HeatConduction
164} // namespace ProcessLib
#define OGS_FATAL(...)
Definition Error.h:19
MathLib::EigenMatrix GlobalMatrix
MathLib::EigenVector GlobalVector
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:22
void WARN(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:34
bool isAxiallySymmetric() const
Definition Mesh.h:128
std::vector< Element * > const & getElements() const
Get the element-vector for the mesh.
Definition Mesh.h:100
void assemble(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, std::vector< std::size_t > const *const sorted_element_subset=nullptr, bool const copy_residua_to_mesh=false)
void preOutput(double const t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id)
void initializeAssemblyOnSubmeshes(std::vector< std::reference_wrapper< MeshLib::Mesh > > const &submeshes, std::vector< std::vector< std::string > > const &residuum_names)
void assembleWithJacobian(double const t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalVector &b, GlobalMatrix &Jac, std::vector< std::size_t > const *const sorted_element_subset=nullptr, bool const copy_residua_to_mesh=false)
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
void computeSecondaryVariableConcrete(double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_prev, int const process_id) override
std::vector< std::unique_ptr< HeatConductionLocalAssemblerInterface > > local_assemblers_
std::vector< std::vector< std::string > > initializeAssemblyOnSubmeshes(std::vector< std::reference_wrapper< MeshLib::Mesh > > const &meshes) override
void preTimestepConcreteProcess(std::vector< GlobalVector * > const &, const double, const double, const int) override
void assembleWithJacobianConcreteProcess(const double t, double const, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalVector &b, GlobalMatrix &Jac) override
void initializeConcreteProcess(NumLib::LocalToGlobalIndexMap const &dof_table, MeshLib::Mesh const &mesh, unsigned const integration_order) override
Process specific initialization called by initialize().
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 preOutputConcreteProcess(const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id) override
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
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)
std::string const name
Definition Process.h:361
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:37
std::vector< std::size_t > const & getActiveElementIDs() const
Definition Process.h:160
SecondaryVariableCollection _secondary_variables
Definition Process.h:369
std::unique_ptr< ProcessLib::AbstractJacobianAssembler > _jacobian_assembler
Definition Process.h:375
NumLib::Extrapolator & getExtrapolator() const
Definition Process.h:201
Handles configuration of several secondary variables from the project file.
void createLocalAssemblers(std::vector< MeshLib::Element * > const &mesh_elements, NumLib::LocalToGlobalIndexMap const &dof_table, std::vector< std::unique_ptr< LocalAssemblerInterface > > &local_assemblers, ProviderOrOrder const &provider_or_order, ExtraCtorArgs &&... extra_ctor_args)
SecondaryVariableFunctions makeExtrapolator(const unsigned num_components, NumLib::Extrapolator &extrapolator, LocalAssemblerCollection const &local_assemblers, typename NumLib::ExtrapolatableLocalAssemblerCollection< LocalAssemblerCollection >::IntegrationPointValuesMethod integration_point_values_method)
static void executeSelectedMemberOnDereferenced(Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)