OGS
ProcessLib::HeatConduction::HeatConductionProcess Class Referencefinal

Detailed Description

Definition at line 21 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)
 
- 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. More...
 
void postTimestep (std::vector< GlobalVector * > const &x, const double t, const double delta_t, int const process_id)
 Postprocessing after a complete timestep. More...
 
void postNonLinearSolver (GlobalVector const &x, GlobalVector const &xdot, 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_dot, int const process_id)
 compute secondary variables for the coupled equations or for output. More...
 
NumLib::IterationResult postIteration (GlobalVector const &x) final
 
void initialize ()
 
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 setCoupledSolutionsForStaggeredScheme (CoupledSolutionsForStaggeredScheme *const coupled_solutions)
 
void updateDeactivatedSubdomains (double const time, const int process_id)
 
bool isMonolithicSchemeUsed () const
 
virtual void setCoupledTermForTheStaggeredSchemeToLocalAssemblers (int 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 &xdot, 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 &xdot, const double dxdot_dx, const double dx_dx, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b, GlobalMatrix &Jac) final
 
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
 
SecondaryVariableCollection const & getSecondaryVariables () const
 
std::vector< std::unique_ptr< IntegrationPointWriter > > const * getIntegrationPointWriter (MeshLib::Mesh const &mesh) 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)
 

ODESystem interface

HeatConductionProcessData _process_data
 
std::vector< std::unique_ptr< HeatConductionLocalAssemblerInterface > > _local_assemblers
 
std::unique_ptr< GlobalVector_x_previous_timestep = nullptr
 Solution of the previous time step. More...
 
MeshLib::PropertyVector< double > * _heat_flux = nullptr
 
bool isLinear () const override
 
void computeSecondaryVariableConcrete (double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_dot, int const process_id) override
 
void initializeConcreteProcess (NumLib::LocalToGlobalIndexMap const &dof_table, MeshLib::Mesh const &mesh, unsigned const integration_order) override
 Process specific initialization called by initialize(). More...
 
void assembleConcreteProcess (const double t, double const, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &xdot, 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 &xdot, const double dxdot_dx, const double dx_dx, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b, GlobalMatrix &Jac) override
 

Additional Inherited Members

- Public Types inherited from ProcessLib::Process
using NonlinearSolver = NumLib::NonlinearSolverBase
 
using TimeDiscretization = NumLib::TimeDiscretization
 
- Public Attributes inherited from ProcessLib::Process
std::string const name
 
- Protected Member Functions inherited from ProcessLib::Process
NumLib::ExtrapolatorgetExtrapolator () const
 
NumLib::LocalToGlobalIndexMap const & getSingleComponentDOFTable () const
 
void initializeProcessBoundaryConditionsAndSourceTerms (const NumLib::LocalToGlobalIndexMap &dof_table, const int process_id)
 
virtual void constructDofTable ()
 
void constructMonolithicProcessDofTable ()
 
void constructDofTableOfSpecifiedProcessStaggeredScheme (const int specified_prosess_id)
 
virtual std::tuple< NumLib::LocalToGlobalIndexMap *, bool > getDOFTableForExtrapolatorData () const
 
- 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
 
VectorMatrixAssembler _global_assembler
 
const bool _use_monolithic_scheme
 
CoupledSolutionsForStaggeredScheme_coupled_solutions
 
unsigned const _integration_order
 
std::vector< std::unique_ptr< 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 
)

Definition at line 25 of file HeatConductionProcess.cpp.

35  : Process(std::move(name), mesh, std::move(jacobian_assembler), parameters,
36  integration_order, std::move(process_variables),
37  std::move(secondary_variables)),
38  _process_data(std::move(process_data))
39 {
40  _heat_flux = MeshLib::getOrCreateMeshProperty<double>(
41  mesh, "HeatFlux", MeshLib::MeshItemType::Node, 1);
42 }
MeshLib::PropertyVector< double > * _heat_flux
std::string const name
Definition: Process.h:323
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:22

References _heat_flux, and MeshLib::Node.

Member Function Documentation

◆ assembleConcreteProcess()

void ProcessLib::HeatConduction::HeatConductionProcess::assembleConcreteProcess ( const double  t,
double const  dt,
std::vector< GlobalVector * > const &  x,
std::vector< GlobalVector * > const &  xdot,
int const  process_id,
GlobalMatrix M,
GlobalMatrix K,
GlobalVector b 
)
overrideprivatevirtual

Implements ProcessLib::Process.

Definition at line 80 of file HeatConductionProcess.cpp.

84 {
85  DBUG("Assemble HeatConductionProcess.");
86 
87  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
88 
89  std::vector<std::reference_wrapper<NumLib::LocalToGlobalIndexMap>>
90  dof_table = {std::ref(*_local_to_global_index_map)};
91  // Call global assembler for each local assembly item.
94  pv.getActiveElementIDs(), dof_table, t, dt, x, xdot, process_id, M, K,
95  b);
96 
100 
101  auto const residuum = computeResiduum(*x[0], *xdot[0], M, K, b);
102  transformVariableFromGlobalVector(residuum, 0 /*variable id*/,
104  std::negate<double>());
105 }
void DBUG(char const *fmt, Args const &... args)
Definition: Logging.h:27
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:145
VectorMatrixAssembler _global_assembler
Definition: Process.h:333
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map
Definition: Process.h:329
void assemble(std::size_t const mesh_item_id, LocalAssemblerInterface &local_assembler, std::vector< std::reference_wrapper< NumLib::LocalToGlobalIndexMap >> const &dof_tables, double const t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &xdot, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b)
void finalizeVectorAssembly(VEC_T &)
General function to finalize the vector assembly.
bool finalizeMatrixAssembly(MAT_T &)
void transformVariableFromGlobalVector(GlobalVector const &input_vector, int const variable_id, NumLib::LocalToGlobalIndexMap const &local_to_global_index_map, MeshLib::PropertyVector< double > &output_vector, Functor mapFunction)
Definition: DOFTableUtil.h:59
GlobalVector computeResiduum(GlobalVector const &x, GlobalVector const &xdot, GlobalMatrix const &M, GlobalMatrix const &K, GlobalVector const &b)
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::assemble(), ProcessLib::computeResiduum(), DBUG(), NumLib::SerialExecutor::executeSelectedMemberDereferenced(), MathLib::finalizeMatrixAssembly(), MathLib::finalizeVectorAssembly(), ProcessLib::ProcessVariable::getActiveElementIDs(), ProcessLib::Process::getProcessVariables(), and NumLib::transformVariableFromGlobalVector().

◆ assembleWithJacobianConcreteProcess()

void ProcessLib::HeatConduction::HeatConductionProcess::assembleWithJacobianConcreteProcess ( const double  t,
double const  dt,
std::vector< GlobalVector * > const &  x,
std::vector< GlobalVector * > const &  xdot,
const double  dxdot_dx,
const double  dx_dx,
int const  process_id,
GlobalMatrix M,
GlobalMatrix K,
GlobalVector b,
GlobalMatrix Jac 
)
overrideprivatevirtual

Implements ProcessLib::Process.

Definition at line 107 of file HeatConductionProcess.cpp.

112 {
113  DBUG("AssembleWithJacobian HeatConductionProcess.");
114 
115  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
116 
117  std::vector<std::reference_wrapper<NumLib::LocalToGlobalIndexMap>>
118  dof_table = {std::ref(*_local_to_global_index_map)};
119  // Call global assembler for each local assembly item.
122  _local_assemblers, pv.getActiveElementIDs(), dof_table, t, dt, x, xdot,
123  dxdot_dx, dx_dx, process_id, M, K, b, Jac);
124 
125  transformVariableFromGlobalVector(b, 0 /*variable id*/,
127  std::negate<double>());
128 }
void assembleWithJacobian(std::size_t const mesh_item_id, LocalAssemblerInterface &local_assembler, std::vector< std::reference_wrapper< NumLib::LocalToGlobalIndexMap >> const &dof_tables, const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &xdot, const double dxdot_dx, const double dx_dx, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b, GlobalMatrix &Jac)

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(), ProcessLib::Process::getProcessVariables(), and NumLib::transformVariableFromGlobalVector().

◆ computeSecondaryVariableConcrete()

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

Reimplemented from ProcessLib::Process.

Definition at line 130 of file HeatConductionProcess.cpp.

133 {
134  DBUG("Compute heat flux for HeatConductionProcess.");
135 
136  std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
137  dof_tables.reserve(x.size());
138  std::generate_n(std::back_inserter(dof_tables), x.size(),
139  [&]() { return _local_to_global_index_map.get(); });
140 
141  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
144  _local_assemblers, pv.getActiveElementIDs(), dof_tables, t, dt, x,
145  x_dot, process_id);
146 }
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_dot, 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 44 of file HeatConductionProcess.cpp.

48 {
49  const int process_id = 0;
50  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
51  ProcessLib::createLocalAssemblers<LocalAssemblerData>(
52  mesh.getDimension(), mesh.getElements(), dof_table,
54  mesh.isAxiallySymmetric(), integration_order, _process_data);
55 
57  "heat_flux_x",
61 
62  if (mesh.getDimension() > 1)
63  {
65  "heat_flux_y",
69  }
70  if (mesh.getDimension() > 2)
71  {
73  "heat_flux_z",
77  }
78 }
virtual std::vector< double > const & getIntPtHeatFluxY(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtHeatFluxX(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtHeatFluxZ(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
unsigned getShapeFunctionOrder() const
NumLib::Extrapolator & getExtrapolator() const
Definition: Process.h:185
SecondaryVariableCollection _secondary_variables
Definition: Process.h:331
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::getDimension(), MeshLib::Mesh::getElements(), ProcessLib::Process::getExtrapolator(), ProcessLib::HeatConduction::HeatConductionLocalAssemblerInterface::getIntPtHeatFluxX(), ProcessLib::HeatConduction::HeatConductionLocalAssemblerInterface::getIntPtHeatFluxY(), ProcessLib::HeatConduction::HeatConductionLocalAssemblerInterface::getIntPtHeatFluxZ(), ProcessLib::Process::getProcessVariables(), ProcessLib::ProcessVariable::getShapeFunctionOrder(), MeshLib::Mesh::isAxiallySymmetric(), and ProcessLib::makeExtrapolator().

◆ isLinear()

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

Definition at line 41 of file HeatConductionProcess.h.

41 { return true; }

Member Data Documentation

◆ _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

◆ _process_data

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

Definition at line 67 of file HeatConductionProcess.h.

Referenced by initializeConcreteProcess().

◆ _x_previous_timestep

std::unique_ptr<GlobalVector> ProcessLib::HeatConduction::HeatConductionProcess::_x_previous_timestep = nullptr
private

Solution of the previous time step.

Definition at line 73 of file HeatConductionProcess.h.


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