![]() |
OGS
|
|
Definition at line 15 of file HeatConductionProcess.h.
#include <HeatConductionProcess.h>
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 ¶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) | |
| 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 ¶meters, 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::Mesh & | getMesh () 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 > |
| 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.
References HeatConductionProcess(), ProcessLib::Process::Process(), ProcessLib::Process::_jacobian_assembler, and ProcessLib::Process::name.
Referenced by HeatConductionProcess(), and AssemblyMixin< HeatConductionProcess >.
|
overrideprivatevirtual |
Implements ProcessLib::Process.
Definition at line 105 of file HeatConductionProcess.cpp.
References ProcessLib::AssemblyMixin< Process >::assemble(), and DBUG().
|
overrideprivatevirtual |
Implements ProcessLib::Process.
Definition at line 116 of file HeatConductionProcess.cpp.
References ProcessLib::AssemblyMixin< Process >::assembleWithJacobian(), and DBUG().
|
overridevirtual |
Reimplemented from ProcessLib::Process.
Definition at line 127 of file HeatConductionProcess.cpp.
References ProcessLib::LocalAssemblerInterface::computeSecondaryVariable(), DBUG(), NumLib::SerialExecutor::executeSelectedMemberOnDereferenced(), ProcessLib::Process::getActiveElementIDs(), and local_assemblers_.
|
overrideprivatevirtual |
Initializes the assembly on submeshes
| meshes | the submeshes on whom the assembly shall proceed. |
meshes must be a must be a non-overlapping cover of the entire simulation domain (bulk mesh)!Reimplemented from ProcessLib::SubmeshAssemblySupport.
Definition at line 92 of file HeatConductionProcess.cpp.
References DBUG(), and ProcessLib::AssemblyMixin< Process >::initializeAssemblyOnSubmeshes().
|
overrideprivatevirtual |
Process specific initialization called by initialize().
Implements ProcessLib::Process.
Definition at line 72 of file HeatConductionProcess.cpp.
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().
|
inlineoverride |
Definition at line 39 of file HeatConductionProcess.h.
References ProcessLib::AssemblyMixinBase::isLinear().
|
overrideprivatevirtual |
Reimplemented from ProcessLib::Process.
Definition at line 153 of file HeatConductionProcess.cpp.
References ProcessLib::AssemblyMixin< Process >::preOutput().
|
overrideprivatevirtual |
Reimplemented from ProcessLib::Process.
Definition at line 144 of file HeatConductionProcess.cpp.
References ProcessLib::AssemblyMixin< Process >::updateActiveElements().
|
inlineoverride |
Definition at line 50 of file HeatConductionProcess.h.
References _ls_compute_only_upon_timestep_change.
|
friend |
Definition at line 15 of file HeatConductionProcess.h.
References HeatConductionProcess(), ProcessLib::Process::Process(), and ProcessLib::Process::name.
|
private |
Definition at line 93 of file HeatConductionProcess.h.
Referenced by shouldLinearSolverComputeOnlyUponTimestepChange().
|
private |
Definition at line 88 of file HeatConductionProcess.h.
Referenced by initializeConcreteProcess().
|
private |
Definition at line 91 of file HeatConductionProcess.h.
Referenced by computeSecondaryVariableConcrete(), and initializeConcreteProcess().