24namespace HeatConduction
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>>>&&
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 },
43 _ls_compute_only_upon_timestep_change{
44 ls_compute_only_upon_timestep_change}
46 if (ls_compute_only_upon_timestep_change)
52 "Using the linear solver compute() method only upon timestep "
53 "change only makes sense for linear model equations.");
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!");
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 "
80 unsigned const integration_order)
97 const double t,
double const dt, std::vector<GlobalVector*>
const& x,
98 std::vector<GlobalVector*>
const& x_prev,
int const process_id,
101 DBUG(
"Assemble HeatConductionProcess.");
103 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_table = {
112 const double t,
double const dt, std::vector<GlobalVector*>
const& x,
113 std::vector<GlobalVector*>
const& x_prev,
int const process_id,
116 DBUG(
"AssembleWithJacobian HeatConductionProcess.");
118 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_table = {
124 process_id, &b, &Jac);
126 transformVariableFromGlobalVector(b, 0 ,
128 std::negate<double>());
132 double const t,
double const dt, std::vector<GlobalVector*>
const& x,
135 DBUG(
"Compute heat flux for HeatConductionProcess.");
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(); });
151 std::vector<GlobalVector*>
const& x,
152 std::vector<GlobalVector*>
const& x_prev,
153 int const process_id)
158 matrix_specification);
160 matrix_specification);
162 matrix_specification);
171 time_residuum.
start();
173 auto const residuum =
computeResiduum(dt, *x[0], *x_prev[0], *M, *K, *b);
178 INFO(
"[time] Computing residuum flow rates took {:g} s",
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
void WARN(fmt::format_string< Args... > fmt, Args &&... args)
double elapsed() const
Get the elapsed time in seconds.
void start()
Start the timer.
Global vector based on Eigen vector.
bool isAxiallySymmetric() const
std::vector< Element * > const & getElements() const
Get the element-vector for the mesh.
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 ¶meters, 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)
HeatConductionProcessData _process_data
std::vector< std::unique_ptr< HeatConductionLocalAssemblerInterface > > _local_assemblers
AssembledMatrixCache _asm_mat_cache
void preOutputConcreteProcess(const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id) override
MeshLib::PropertyVector< double > * _heat_flux
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
SecondaryVariableCollection _secondary_variables
VectorMatrixAssembler _global_assembler
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map
MathLib::MatrixSpecifications getMatrixSpecifications(const int process_id) const override
NumLib::Extrapolator & getExtrapolator() const
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)
int const mesh_space_dimension