|  | OGS
    | 
Definition at line 22 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 std::vector< std::vector< std::string > > | initializeAssemblyOnSubmeshes (std::vector< std::reference_wrapper< MeshLib::Mesh > > const &meshes) | 
| 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(). | |
| 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 | preOutputConcreteProcess (const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id) override | 
| Private Attributes | |
| HeatConductionProcessData | _process_data | 
| std::vector< std::unique_ptr< HeatConductionLocalAssemblerInterface > > | _local_assemblers | 
| MeshLib::PropertyVector< double > * | _heat_flux = nullptr | 
| AssembledMatrixCache | _asm_mat_cache | 
| bool const | _ls_compute_only_upon_timestep_change | 
| 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 26 of file HeatConductionProcess.cpp.
References _heat_flux, MeshLib::getOrCreateMeshProperty(), MeshLib::Node, OGS_FATAL, and WARN().
| 
 | overrideprivatevirtual | 
Implements ProcessLib::Process.
Definition at line 96 of file HeatConductionProcess.cpp.
References _asm_mat_cache, ProcessLib::Process::_global_assembler, _local_assemblers, ProcessLib::Process::_local_to_global_index_map, ProcessLib::AssembledMatrixCache::assemble(), DBUG(), and ProcessLib::Process::getActiveElementIDs().
Referenced by preOutputConcreteProcess().
| 
 | overrideprivatevirtual | 
Implements ProcessLib::Process.
Definition at line 111 of file HeatConductionProcess.cpp.
References ProcessLib::Process::_global_assembler, _heat_flux, _local_assemblers, ProcessLib::Process::_local_to_global_index_map, ProcessLib::VectorMatrixAssembler::assembleWithJacobian(), DBUG(), NumLib::SerialExecutor::executeSelectedMemberDereferenced(), and ProcessLib::Process::getActiveElementIDs().
| 
 | overridevirtual | 
Reimplemented from ProcessLib::Process.
Definition at line 131 of file HeatConductionProcess.cpp.
References _local_assemblers, ProcessLib::LocalAssemblerInterface::computeSecondaryVariable(), DBUG(), NumLib::SerialExecutor::executeSelectedMemberOnDereferenced(), and ProcessLib::Process::getActiveElementIDs().
| 
 | overrideprivatevirtual | 
Process specific initialization called by initialize().
Implements ProcessLib::Process.
Definition at line 77 of file HeatConductionProcess.cpp.
References _local_assemblers, _process_data, ProcessLib::Process::_secondary_variables, ProcessLib::SecondaryVariableCollection::addSecondaryVariable(), ProcessLib::createLocalAssemblers(), MeshLib::Mesh::getElements(), ProcessLib::Process::getExtrapolator(), ProcessLib::HeatConduction::HeatConductionLocalAssemblerInterface::getIntPtHeatFlux(), MeshLib::Mesh::isAxiallySymmetric(), ProcessLib::makeExtrapolator(), and ProcessLib::HeatConduction::HeatConductionProcessData::mesh_space_dimension.
| 
 | inlineoverride | 
Definition at line 43 of file HeatConductionProcess.h.
References _asm_mat_cache, and ProcessLib::AssembledMatrixCache::isLinear().
| 
 | overrideprivatevirtual | 
Reimplemented from ProcessLib::Process.
Definition at line 148 of file HeatConductionProcess.cpp.
References _heat_flux, ProcessLib::Process::_local_to_global_index_map, assembleConcreteProcess(), ProcessLib::computeResiduum(), BaseLib::RunTime::elapsed(), ProcessLib::Process::getMatrixSpecifications(), INFO(), and BaseLib::RunTime::start().
| 
 | inlineoverride | 
Definition at line 51 of file HeatConductionProcess.h.
References _ls_compute_only_upon_timestep_change.
| 
 | private | 
Definition at line 87 of file HeatConductionProcess.h.
Referenced by assembleConcreteProcess(), and isLinear().
| 
 | private | 
Definition at line 85 of file HeatConductionProcess.h.
Referenced by HeatConductionProcess(), assembleWithJacobianConcreteProcess(), and preOutputConcreteProcess().
| 
 | private | 
Definition at line 83 of file HeatConductionProcess.h.
Referenced by assembleConcreteProcess(), assembleWithJacobianConcreteProcess(), computeSecondaryVariableConcrete(), and initializeConcreteProcess().
| 
 | private | 
Definition at line 89 of file HeatConductionProcess.h.
Referenced by shouldLinearSolverComputeOnlyUponTimestepChange().
| 
 | private | 
Definition at line 80 of file HeatConductionProcess.h.
Referenced by initializeConcreteProcess().