OGS
ThermoRichardsFlowFEM.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) OpenGeoSys Community (opengeosys.org)
2// SPDX-License-Identifier: BSD-3-Clause
3
4#pragma once
5
6#include <memory>
7#include <vector>
8
21
22namespace ProcessLib
23{
24namespace ThermoRichardsFlow
25{
26namespace MPL = MaterialPropertyLib;
27
28template <typename ShapeFunction, int GlobalDim>
30{
31public:
32 // Note: temperature variable uses the same shape functions as that are used
33 // by pressure variable.
35
38
40
42 delete;
44 delete;
45
47 MeshLib::Element const& e,
48 std::size_t const /*local_matrix_size*/,
49 NumLib::GenericIntegrationMethod const& integration_method,
50 bool const is_axially_symmetric,
51 ThermoRichardsFlowProcessData& process_data);
52
55 std::string_view const name,
56 double const* values,
57 int const integration_order) override;
58
59 void setInitialConditionsConcrete(Eigen::VectorXd const local_x,
60 double const t,
61 int const process_id) override;
62
63 void assembleWithJacobian(double const t, double const dt,
64 std::vector<double> const& local_x,
65 std::vector<double> const& local_x_prev,
66 std::vector<double>& local_rhs_data,
67 std::vector<double>& local_Jac_data) override;
68
69 void assemble(double const t, double const dt,
70 std::vector<double> const& local_x,
71 std::vector<double> const& local_x_prev,
72 std::vector<double>& local_M_data,
73 std::vector<double>& local_K_data,
74 std::vector<double>& local_rhs_data) override;
75
76 void initializeConcrete() override
77 {
78 unsigned const n_integration_points =
79 _integration_method.getNumberOfPoints();
80
81 for (unsigned ip = 0; ip < n_integration_points; ip++)
82 {
83 auto& ip_data = _ip_data[ip];
84 ip_data.pushBackState();
85 }
86 }
87
88 void postTimestepConcrete(Eigen::VectorXd const& /*local_x*/,
89 Eigen::VectorXd const& /*local_x_prev*/,
90 double const /*t*/, double const /*dt*/,
91 int const /*process_id*/) override
92 {
93 unsigned const n_integration_points =
94 _integration_method.getNumberOfPoints();
95
96 for (unsigned ip = 0; ip < n_integration_points; ip++)
97 {
98 _ip_data[ip].pushBackState();
99 }
100 }
101
103 double const t, double const dt, Eigen::VectorXd const& local_x,
104 Eigen::VectorXd const& local_x_prev) override;
105
106 Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
107 const unsigned integration_point) const override
108 {
109 auto const& N = _ip_data[integration_point].N;
110
111 // assumes N is stored contiguously in memory
112 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
113 }
114
115 std::vector<double> const& getIntPtDarcyVelocity(
116 const double t,
117 std::vector<GlobalVector*> const& x,
118 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
119 std::vector<double>& cache) const override;
120
121 std::vector<double> getSaturation() const override;
122 std::vector<double> const& getIntPtSaturation(
123 const double t,
124 std::vector<GlobalVector*> const& x,
125 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
126 std::vector<double>& cache) const override;
127
128 std::vector<double> getPorosity() const override;
129 std::vector<double> const& getIntPtPorosity(
130 const double t,
131 std::vector<GlobalVector*> const& x,
132 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
133 std::vector<double>& cache) const override;
134
135 std::vector<double> const& getIntPtDryDensitySolid(
136 const double t,
137 std::vector<GlobalVector*> const& x,
138 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
139 std::vector<double>& cache) const override;
140
141private:
142 unsigned getNumberOfIntegrationPoints() const override;
143
144private:
146
147 std::vector<IpData, Eigen::aligned_allocator<IpData>> _ip_data;
148
152
153 static const int temperature_index = 0;
154 static const int temperature_size = ShapeFunction::NPOINTS;
156 static const int pressure_size = ShapeFunction::NPOINTS;
157};
158
159} // namespace ThermoRichardsFlow
160} // namespace ProcessLib
161
EigenFixedShapeMatrixPolicy< ShapeFunction, GlobalDim > ShapeMatrixPolicyType
void postTimestepConcrete(Eigen::VectorXd const &, Eigen::VectorXd const &, double const, double const, int const) override
void assembleWithJacobian(double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &local_x_prev, std::vector< double > &local_rhs_data, std::vector< double > &local_Jac_data) override
std::vector< double > const & getIntPtSaturation(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
std::size_t setIPDataInitialConditions(std::string_view const name, double const *values, int const integration_order) override
ThermoRichardsFlowLocalAssembler(ThermoRichardsFlowLocalAssembler &&)=delete
std::vector< double > const & getIntPtDryDensitySolid(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
std::vector< double > const & getIntPtDarcyVelocity(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
typename ShapeMatricesType::GlobalDimVectorType GlobalDimVectorType
void assemble(double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &local_x_prev, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_rhs_data) override
ThermoRichardsFlowLocalAssembler(ThermoRichardsFlowLocalAssembler const &)=delete
std::vector< double > const & getIntPtPorosity(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
void setInitialConditionsConcrete(Eigen::VectorXd const local_x, double const t, int const process_id) override
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix(const unsigned integration_point) const override
Provides the shape matrix at the given integration point.
typename ShapeMatricesType::GlobalDimMatrixType GlobalDimMatrixType
ShapeMatrixPolicyType< ShapeFunction, GlobalDim > ShapeMatricesType
void computeSecondaryVariableConcrete(double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev) override
std::vector< IpData, Eigen::aligned_allocator< IpData > > _ip_data
MatrixType< GlobalDim, GlobalDim > GlobalDimMatrixType
VectorType< GlobalDim > GlobalDimVectorType