OGS
ProcessLib::HeatConduction::HeatConductionProcess Class Referencefinal

Detailed Description

Definition at line 15 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, 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::vector< std::reference_wrapper< ProcessVariable > > > const & getProcessVariables () 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)
bool requiresNormalization () const override
Public Member Functions inherited from ProcessLib::SubmeshAssemblySupport
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().
std::vector< std::vector< std::string > > initializeAssemblyOnSubmeshes (std::vector< std::reference_wrapper< MeshLib::Mesh > > const &meshes) 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
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 preTimestepConcreteProcess (std::vector< GlobalVector * > const &, const double, const double, const int) 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 Member Functions inherited from ProcessLib::AssemblyMixin< HeatConductionProcess >
void initializeAssemblyOnSubmeshes (std::vector< std::reference_wrapper< MeshLib::Mesh > > const &submeshes, std::vector< std::vector< std::string > > const &residuum_names)
void updateActiveElements ()
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 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)
void preOutput (double const t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id)

Private Attributes

HeatConductionProcessData _process_data
std::vector< std::unique_ptr< HeatConductionLocalAssemblerInterface > > local_assemblers_
bool const _ls_compute_only_upon_timestep_change

Friends

class AssemblyMixin< HeatConductionProcess >

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
void setReleaseNodalForces (GlobalVector const *r_neq, int const process_id) 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
CellAverageData cell_average_data_
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 19 of file HeatConductionProcess.cpp.

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)),
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}
#define OGS_FATAL(...)
Definition Error.h:19
void WARN(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:34
std::string const name
Definition Process.h:361
std::unique_ptr< ProcessLib::AbstractJacobianAssembler > _jacobian_assembler
Definition Process.h:375

References HeatConductionProcess(), ProcessLib::Process::Process(), ProcessLib::Process::_jacobian_assembler, and ProcessLib::Process::name.

Referenced by HeatConductionProcess(), and AssemblyMixin< HeatConductionProcess >.

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 105 of file HeatConductionProcess.cpp.

109{
110 DBUG("Assemble HeatConductionProcess.");
111
112 AssemblyMixin<HeatConductionProcess>::assemble(t, dt, x, x_prev, process_id,
113 M, K, b);
114}
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:22
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)

References ProcessLib::AssemblyMixin< Process >::assemble(), and DBUG().

◆ 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,
GlobalVector & b,
GlobalMatrix & Jac )
overrideprivatevirtual

Implements ProcessLib::Process.

Definition at line 116 of file HeatConductionProcess.cpp.

120{
121 DBUG("AssembleWithJacobian HeatConductionProcess.");
122
124 t, dt, x, x_prev, process_id, b, Jac);
125}
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)

References ProcessLib::AssemblyMixin< Process >::assembleWithJacobian(), and DBUG().

◆ 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 127 of file HeatConductionProcess.cpp.

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}
std::vector< std::unique_ptr< HeatConductionLocalAssemblerInterface > > local_assemblers_
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:160
static void executeSelectedMemberOnDereferenced(Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)

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

◆ initializeAssemblyOnSubmeshes()

std::vector< std::vector< std::string > > ProcessLib::HeatConduction::HeatConductionProcess::initializeAssemblyOnSubmeshes ( std::vector< std::reference_wrapper< MeshLib::Mesh > > const & meshes)
overrideprivatevirtual

Initializes the assembly on submeshes

Parameters
meshesthe submeshes on whom the assembly shall proceed.
Attention
meshes must be a must be a non-overlapping cover of the entire simulation domain (bulk mesh)!
Returns
The names of the residuum vectors that will be assembled for each process: outer vector of size 1 for monolithic schemes and greater for staggered schemes.

Reimplemented from ProcessLib::SubmeshAssemblySupport.

Definition at line 92 of file HeatConductionProcess.cpp.

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}
void initializeAssemblyOnSubmeshes(std::vector< std::reference_wrapper< MeshLib::Mesh > > const &submeshes, std::vector< std::vector< std::string > > const &residuum_names)

References DBUG(), and ProcessLib::AssemblyMixin< Process >::initializeAssemblyOnSubmeshes().

◆ 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 72 of file HeatConductionProcess.cpp.

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}
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:369
NumLib::Extrapolator & getExtrapolator() const
Definition Process.h:201
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)

References _process_data, ProcessLib::Process::_secondary_variables, ProcessLib::createLocalAssemblers(), MeshLib::Mesh::getElements(), ProcessLib::Process::getExtrapolator(), ProcessLib::HeatConduction::HeatConductionLocalAssemblerInterface::getIntPtHeatFlux(), MeshLib::Mesh::isAxiallySymmetric(), local_assemblers_, and ProcessLib::makeExtrapolator().

◆ 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{
161 process_id);
162}
void preOutput(double const t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id)

References ProcessLib::AssemblyMixin< Process >::preOutput().

◆ preTimestepConcreteProcess()

void ProcessLib::HeatConduction::HeatConductionProcess::preTimestepConcreteProcess ( std::vector< GlobalVector * > const & ,
const double ,
const double ,
const int  )
overrideprivatevirtual

◆ shouldLinearSolverComputeOnlyUponTimestepChange()

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

Definition at line 50 of file HeatConductionProcess.h.

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

References _ls_compute_only_upon_timestep_change.

◆ AssemblyMixin< HeatConductionProcess >

Member Data Documentation

◆ _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 88 of file HeatConductionProcess.h.

Referenced by initializeConcreteProcess().

◆ local_assemblers_

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

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