OGS
HeatConductionProcess.cpp
Go to the documentation of this file.
1
12
13#include <cassert>
14
21
22namespace ProcessLib
23{
24namespace HeatConduction
25{
27 std::string name,
28 MeshLib::Mesh& mesh,
29 std::unique_ptr<ProcessLib::AbstractJacobianAssembler>&& jacobian_assembler,
30 std::vector<std::unique_ptr<ParameterLib::ParameterBase>> const& parameters,
31 unsigned const integration_order,
32 std::vector<std::vector<std::reference_wrapper<ProcessVariable>>>&&
33 process_variables,
34 HeatConductionProcessData&& process_data,
35 SecondaryVariableCollection&& secondary_variables,
36 bool const is_linear,
37 bool const ls_compute_only_upon_timestep_change)
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*/},
43 _ls_compute_only_upon_timestep_change{
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
74 mesh, "HeatFlowRate", MeshLib::MeshItemType::Node, 1);
75}
76
78 NumLib::LocalToGlobalIndexMap const& dof_table,
79 MeshLib::Mesh const& mesh,
80 unsigned const integration_order)
81{
82 int const mesh_space_dimension = _process_data.mesh_space_dimension;
83
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}
95
97 const double t, double const dt, std::vector<GlobalVector*> const& x,
98 std::vector<GlobalVector*> const& x_prev, int const process_id,
100{
101 DBUG("Assemble HeatConductionProcess.");
102
103 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_table = {
105
106 _asm_mat_cache.assemble(t, dt, x, x_prev, process_id, &M, &K, &b, dof_table,
109}
110
112 const double t, double const dt, std::vector<GlobalVector*> const& x,
113 std::vector<GlobalVector*> const& x_prev, int const process_id,
114 GlobalVector& b, GlobalMatrix& Jac)
115{
116 DBUG("AssembleWithJacobian HeatConductionProcess.");
117
118 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_table = {
120 // Call global assembler for each local assembly item.
123 _local_assemblers, getActiveElementIDs(), dof_table, t, dt, x, x_prev,
124 process_id, &b, &Jac);
125
126 transformVariableFromGlobalVector(b, 0 /*variable id*/,
128 std::negate<double>());
129}
130
132 double const t, double const dt, std::vector<GlobalVector*> const& x,
133 GlobalVector const& x_prev, int const process_id)
134{
135 DBUG("Compute heat flux for HeatConductionProcess.");
136
137 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
138 dof_tables.reserve(x.size());
139 std::generate_n(std::back_inserter(dof_tables), x.size(),
140 [&]() { return _local_to_global_index_map.get(); });
141
144 _local_assemblers, getActiveElementIDs(), dof_tables, t, dt, x, x_prev,
145 process_id);
146}
147
149 const double t,
150 double const dt,
151 std::vector<GlobalVector*> const& x,
152 std::vector<GlobalVector*> const& x_prev,
153 int const process_id)
154{
155 auto const matrix_specification = getMatrixSpecifications(process_id);
156
158 matrix_specification);
160 matrix_specification);
162 matrix_specification);
163
164 M->setZero();
165 K->setZero();
166 b->setZero();
167
168 assembleConcreteProcess(t, dt, x, x_prev, process_id, *M, *K, *b);
169
170 BaseLib::RunTime time_residuum;
171 time_residuum.start();
172
173 auto const residuum = computeResiduum(dt, *x[0], *x_prev[0], *M, *K, *b);
174
175 transformVariableFromGlobalVector(residuum, 0, *_local_to_global_index_map,
176 *_heat_flux, std::negate<double>());
177
178 INFO("[time] Computing residuum flow rates took {:g} s",
179 time_residuum.elapsed());
180}
181} // namespace HeatConduction
182} // namespace ProcessLib
#define OGS_FATAL(...)
Definition Error.h:26
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:35
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:30
void WARN(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:40
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
Global vector based on Eigen vector.
Definition EigenVector.h:25
bool isAxiallySymmetric() const
Definition Mesh.h:137
std::vector< Element * > const & getElements() const
Get the element-vector for the mesh.
Definition Mesh.h:109
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
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)
std::vector< std::unique_ptr< HeatConductionLocalAssemblerInterface > > _local_assemblers
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::vector< std::size_t > const & getActiveElementIDs() const
Definition Process.h:167
SecondaryVariableCollection _secondary_variables
Definition Process.h:370
VectorMatrixAssembler _global_assembler
Definition Process.h:377
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map
Definition Process.h:368
MathLib::MatrixSpecifications getMatrixSpecifications(const int process_id) const override
Definition Process.cpp:211
NumLib::Extrapolator & getExtrapolator() const
Definition Process.h:208
Handles configuration of several secondary variables from the project file.
void addSecondaryVariable(std::string const &internal_name, SecondaryVariableFunctions &&fcts)
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)
PropertyVector< T > * getOrCreateMeshProperty(Mesh &mesh, std::string const &property_name, MeshItemType const item_type, int const number_of_components)
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)
GlobalVector computeResiduum(double const dt, GlobalVector const &x, GlobalVector const &x_prev, GlobalMatrix const &M, GlobalMatrix const &K, GlobalVector const &b)
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)
static void executeSelectedMemberDereferenced(Object &object, Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)
void assemble(const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalMatrix *const M, GlobalMatrix *const K, GlobalVector *const 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)