OGS
|
The implementation uses a monolithic approach, i.e., both processes are assembled within one global system of equations.
The advective term of the heat conduction equation is given by the confined groundwater flow process, i.e., the heat conduction depends on darcy velocity of the groundwater flow process. On the other hand the temperature dependencies of the viscosity and density in the groundwater flow couples the H process to the T process.
Definition at line 52 of file HTProcess.h.
#include <HTProcess.h>
Public Member Functions | |
HTProcess (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, HTProcessData &&process_data, SecondaryVariableCollection &&secondary_variables, bool const use_monolithic_scheme, std::unique_ptr< ProcessLib::SurfaceFluxData > &&surfaceflux) | |
Eigen::Vector3d | getFlux (std::size_t element_id, MathLib::Point3d const &p, double const t, std::vector< GlobalVector * > const &x) const override |
void | postTimestepConcreteProcess (std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, const double t, const double delta_t, int const process_id) 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 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 dt, 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 dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalVector &b, GlobalMatrix &Jac) override |
std::tuple< NumLib::LocalToGlobalIndexMap *, bool > | getDOFTableForExtrapolatorData () const override |
Private Attributes | |
HTProcessData | _process_data |
std::vector< std::unique_ptr< HTLocalAssemblerInterface > > | _local_assemblers |
std::unique_ptr< ProcessLib::SurfaceFluxData > | _surfaceflux |
ProcessLib::HT::HTProcess::HTProcess | ( | 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, | ||
HTProcessData && | process_data, | ||
SecondaryVariableCollection && | secondary_variables, | ||
bool const | use_monolithic_scheme, | ||
std::unique_ptr< ProcessLib::SurfaceFluxData > && | surfaceflux ) |
Definition at line 27 of file HTProcess.cpp.
|
overrideprivatevirtual |
Implements ProcessLib::Process.
Definition at line 76 of file HTProcess.cpp.
References ProcessLib::Process::_global_assembler, _local_assemblers, ProcessLib::Process::_local_to_global_index_map, _process_data, ProcessLib::Process::_use_monolithic_scheme, ProcessLib::VectorMatrixAssembler::assemble(), DBUG(), NumLib::SerialExecutor::executeSelectedMemberDereferenced(), ProcessLib::Process::getActiveElementIDs(), and ProcessLib::HT::HTProcessData::heat_transport_process_id.
|
overrideprivatevirtual |
Implements ProcessLib::Process.
Definition at line 112 of file HTProcess.cpp.
References ProcessLib::Process::_global_assembler, _local_assemblers, ProcessLib::Process::_local_to_global_index_map, ProcessLib::Process::_use_monolithic_scheme, ProcessLib::VectorMatrixAssembler::assembleWithJacobian(), DBUG(), NumLib::SerialExecutor::executeSelectedMemberDereferenced(), and ProcessLib::Process::getActiveElementIDs().
|
overrideprivatevirtual |
Get the address of a LocalToGlobalIndexMap, and the status of its memory. If the LocalToGlobalIndexMap is created as new in this function, the function also returns a true boolean value to let Extrapolator manage the memory by the address returned by this function.
Reimplemented from ProcessLib::Process.
Definition at line 138 of file HTProcess.cpp.
References ProcessLib::Process::_local_to_global_index_map, ProcessLib::Process::_mesh_subset_all_nodes, ProcessLib::Process::_use_monolithic_scheme, and NumLib::BY_LOCATION.
|
overridevirtual |
Reimplemented from ProcessLib::Process.
Definition at line 161 of file HTProcess.cpp.
References _local_assemblers, ProcessLib::Process::_local_to_global_index_map, ProcessLib::getCoupledLocalSolutions(), and NumLib::getRowColumnIndices().
|
overrideprivatevirtual |
Process specific initialization called by initialize().
Implements ProcessLib::Process.
Definition at line 47 of file HTProcess.cpp.
References _local_assemblers, _process_data, ProcessLib::Process::_secondary_variables, ProcessLib::Process::_use_monolithic_scheme, ProcessLib::SecondaryVariableCollection::addSecondaryVariable(), ProcessLib::createLocalAssemblers(), MeshLib::Mesh::getElements(), ProcessLib::Process::getExtrapolator(), ProcessLib::HT::HTLocalAssemblerInterface::getIntPtDarcyVelocity(), MeshLib::Mesh::isAxiallySymmetric(), ProcessLib::makeExtrapolator(), and ProcessLib::HT::HTProcessData::mesh_space_dimension.
|
inlineoverride |
Definition at line 72 of file HTProcess.h.
|
overridevirtual |
Reimplemented from ProcessLib::Process.
Definition at line 179 of file HTProcess.cpp.
References ProcessLib::Process::_integration_order, ProcessLib::Process::_mesh, _process_data, _surfaceflux, ProcessLib::Process::_use_monolithic_scheme, DBUG(), ProcessLib::Process::getActiveElementIDs(), ProcessLib::HT::HTProcessData::hydraulic_process_id, and OGS_FATAL.
|
private |
Definition at line 111 of file HTProcess.h.
Referenced by assembleConcreteProcess(), assembleWithJacobianConcreteProcess(), getFlux(), and initializeConcreteProcess().
|
private |
Definition at line 109 of file HTProcess.h.
Referenced by assembleConcreteProcess(), initializeConcreteProcess(), and postTimestepConcreteProcess().
|
private |
Definition at line 113 of file HTProcess.h.
Referenced by postTimestepConcreteProcess().