Loading [MathJax]/extensions/MathZoom.js
OGS
AdvectionMatrixAssembler.h
Go to the documentation of this file.
1
10#pragma once
11
12#include <Eigen/Core>
13#include <limits>
14#include <memory>
15#include <vector>
16
19
20namespace NumLib
21{
22namespace detail
23{
24template <typename MeshElementType,
25 typename IPData,
26 typename FluxVectorType,
27 typename Derived>
28void assembleAdvectionMatrix(IPData const& ip_data_vector,
29 NumLib::ShapeMatrixCache const& shape_matrix_cache,
30 std::vector<FluxVectorType> const& ip_flux_vector,
31 Eigen::MatrixBase<Derived>& laplacian_matrix)
32{
33 auto const& Ns = shape_matrix_cache.NsHigherOrder<MeshElementType>();
34
35 for (std::size_t ip = 0; ip < ip_flux_vector.size(); ++ip)
36 {
37 auto const& ip_data = ip_data_vector[ip];
38 auto const w = ip_data.integration_weight;
39 auto const& dNdx = ip_data.dNdx;
40 auto const& N = Ns[ip];
41 laplacian_matrix.noalias() +=
42 (N.transpose() * ip_flux_vector[ip].transpose() * dNdx * w)
43#if defined(__GNUC__) && (__GNUC__ == 14) && (__GNUC_MINOR__ == 2) && \
44 (__GNUC_PATCHLEVEL__ == 1)
45 .eval()
46#endif
47 ;
48 }
49}
50
51template <typename IPData, typename FluxVectorType, typename Derived>
52void assembleAdvectionMatrix(IPData const& ip_data_vector,
53 std::vector<FluxVectorType> const& ip_flux_vector,
54 Eigen::MatrixBase<Derived>& laplacian_matrix)
55{
56 for (std::size_t ip = 0; ip < ip_flux_vector.size(); ++ip)
57 {
58 auto const& ip_data = ip_data_vector[ip];
59 auto const w = ip_data.integration_weight;
60 auto const& N = ip_data.N;
61 auto const& dNdx = ip_data.dNdx;
62 laplacian_matrix.noalias() +=
63 N.transpose() * ip_flux_vector[ip].transpose() * dNdx * w;
64 }
65}
66
67template <typename Derived>
68void applyFullUpwind(Eigen::VectorXd const& quasi_nodal_flux,
69 Eigen::MatrixBase<Derived>& laplacian_matrix)
70{
71 Eigen::VectorXd const down_mask =
72 (quasi_nodal_flux.array() < 0).cast<double>();
73 Eigen::VectorXd const down = quasi_nodal_flux.cwiseProduct(down_mask);
74
75 double const q_in = -down.sum();
76 if (q_in < std::numeric_limits<double>::epsilon())
77 {
78 return;
79 }
80
81 Eigen::VectorXd const up_mask =
82 (quasi_nodal_flux.array() >= 0).cast<double>();
83 Eigen::VectorXd const up = quasi_nodal_flux.cwiseProduct(up_mask);
84
85 laplacian_matrix.diagonal().noalias() += up;
86 laplacian_matrix.noalias() += down * up.transpose() / q_in;
87}
88
89template <typename IPData, typename FluxVectorType, typename Derived>
90void applyFullUpwind(IPData const& ip_data_vector,
91 std::vector<FluxVectorType> const& ip_flux_vector,
92 Eigen::MatrixBase<Derived>& laplacian_matrix)
93{
94 Eigen::VectorXd quasi_nodal_flux =
95 Eigen::VectorXd::Zero(laplacian_matrix.rows());
96
97 for (std::size_t ip = 0; ip < ip_flux_vector.size(); ++ip)
98 {
99 auto const& ip_data = ip_data_vector[ip];
100 auto const w = ip_data.integration_weight;
101 auto const& dNdx = ip_data.dNdx;
102
103 quasi_nodal_flux.noalias() -= ip_flux_vector[ip].transpose() * dNdx * w;
104 }
105
106 applyFullUpwind(quasi_nodal_flux, laplacian_matrix);
107}
108} // namespace detail
109
110template <typename MeshElementType,
111 typename IPData,
112 typename FluxVectorType,
113 typename Derived>
115 IPData const& ip_data_vector,
116 NumLib::ShapeMatrixCache const& shape_matrix_cache,
117 std::vector<FluxVectorType> const& ip_flux_vector,
118 double const average_velocity,
119 Eigen::MatrixBase<Derived>& laplacian_matrix)
120{
121 std::visit(
122 [&](auto&& stabilizer)
123 {
124 using Stabilizer = std::decay_t<decltype(stabilizer)>;
125 if constexpr (std::is_same_v<Stabilizer, FullUpwind>)
126 {
127 if (average_velocity > stabilizer.getCutoffVelocity())
128 {
129 detail::applyFullUpwind(ip_data_vector, ip_flux_vector,
130 laplacian_matrix);
131 return;
132 }
133 }
134
136 ip_data_vector, shape_matrix_cache, ip_flux_vector,
137 laplacian_matrix);
138 },
139 stabilizer);
140}
141
142template <typename IPData, typename FluxVectorType, typename Derived>
144 IPData const& ip_data_vector,
145 std::vector<FluxVectorType> const& ip_flux_vector,
146 double const average_velocity,
147 Eigen::MatrixBase<Derived>& laplacian_matrix)
148{
149 std::visit(
150 [&](auto&& stabilizer)
151 {
152 using Stabilizer = std::decay_t<decltype(stabilizer)>;
153 if constexpr (std::is_same_v<Stabilizer, FullUpwind>)
154 {
155 if (average_velocity > stabilizer.getCutoffVelocity())
156 {
157 detail::applyFullUpwind(ip_data_vector, ip_flux_vector,
158 laplacian_matrix);
159 return;
160 }
161 }
162
164 ip_data_vector, ip_flux_vector, laplacian_matrix);
165 },
166 stabilizer);
167}
168} // namespace NumLib
auto const & NsHigherOrder() const
void assembleAdvectionMatrix(IPData const &ip_data_vector, NumLib::ShapeMatrixCache const &shape_matrix_cache, std::vector< FluxVectorType > const &ip_flux_vector, Eigen::MatrixBase< Derived > &laplacian_matrix)
void applyFullUpwind(Eigen::VectorXd const &quasi_nodal_flux, Eigen::MatrixBase< Derived > &laplacian_matrix)
std::variant< NoStabilization, IsotropicDiffusionStabilization, FullUpwind, FluxCorrectedTransport > NumericalStabilization