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 }
44}
45
46template <typename IPData, typename FluxVectorType, typename Derived>
47void assembleAdvectionMatrix(IPData const& ip_data_vector,
48 std::vector<FluxVectorType> const& ip_flux_vector,
49 Eigen::MatrixBase<Derived>& laplacian_matrix)
50{
51 for (std::size_t ip = 0; ip < ip_flux_vector.size(); ++ip)
52 {
53 auto const& ip_data = ip_data_vector[ip];
54 auto const w = ip_data.integration_weight;
55 auto const& N = ip_data.N;
56 auto const& dNdx = ip_data.dNdx;
57 laplacian_matrix.noalias() +=
58 N.transpose() * ip_flux_vector[ip].transpose() * dNdx * w;
59 }
60}
61
62template <typename Derived>
63void applyFullUpwind(Eigen::VectorXd const& quasi_nodal_flux,
64 Eigen::MatrixBase<Derived>& laplacian_matrix)
65{
66 Eigen::VectorXd const down_mask =
67 (quasi_nodal_flux.array() < 0).cast<double>();
68 Eigen::VectorXd const down = quasi_nodal_flux.cwiseProduct(down_mask);
69
70 double const q_in = -down.sum();
71 if (q_in < std::numeric_limits<double>::epsilon())
72 {
73 return;
74 }
75
76 Eigen::VectorXd const up_mask =
77 (quasi_nodal_flux.array() >= 0).cast<double>();
78 Eigen::VectorXd const up = quasi_nodal_flux.cwiseProduct(up_mask);
79
80 laplacian_matrix.diagonal().noalias() += up;
81 laplacian_matrix.noalias() += down * up.transpose() / q_in;
82}
83
84template <typename IPData, typename FluxVectorType, typename Derived>
85void applyFullUpwind(IPData const& ip_data_vector,
86 std::vector<FluxVectorType> const& ip_flux_vector,
87 Eigen::MatrixBase<Derived>& laplacian_matrix)
88{
89 Eigen::VectorXd quasi_nodal_flux =
90 Eigen::VectorXd::Zero(laplacian_matrix.rows());
91
92 for (std::size_t ip = 0; ip < ip_flux_vector.size(); ++ip)
93 {
94 auto const& ip_data = ip_data_vector[ip];
95 auto const w = ip_data.integration_weight;
96 auto const& dNdx = ip_data.dNdx;
97
98 quasi_nodal_flux.noalias() -= ip_flux_vector[ip].transpose() * dNdx * w;
99 }
100
101 applyFullUpwind(quasi_nodal_flux, laplacian_matrix);
102}
103} // namespace detail
104
105template <typename MeshElementType,
106 typename IPData,
107 typename FluxVectorType,
108 typename Derived>
110 IPData const& ip_data_vector,
111 NumLib::ShapeMatrixCache const& shape_matrix_cache,
112 std::vector<FluxVectorType> const& ip_flux_vector,
113 double const average_velocity,
114 Eigen::MatrixBase<Derived>& laplacian_matrix)
115{
116 std::visit(
117 [&](auto&& stabilizer)
118 {
119 using Stabilizer = std::decay_t<decltype(stabilizer)>;
120 if constexpr (std::is_same_v<Stabilizer, FullUpwind>)
121 {
122 if (average_velocity > stabilizer.getCutoffVelocity())
123 {
124 detail::applyFullUpwind(ip_data_vector, ip_flux_vector,
125 laplacian_matrix);
126 return;
127 }
128 }
129
131 ip_data_vector, shape_matrix_cache, ip_flux_vector,
132 laplacian_matrix);
133 },
134 stabilizer);
135}
136
137template <typename IPData, typename FluxVectorType, typename Derived>
139 IPData const& ip_data_vector,
140 std::vector<FluxVectorType> const& ip_flux_vector,
141 double const average_velocity,
142 Eigen::MatrixBase<Derived>& laplacian_matrix)
143{
144 std::visit(
145 [&](auto&& stabilizer)
146 {
147 using Stabilizer = std::decay_t<decltype(stabilizer)>;
148 if constexpr (std::is_same_v<Stabilizer, FullUpwind>)
149 {
150 if (average_velocity > stabilizer.getCutoffVelocity())
151 {
152 detail::applyFullUpwind(ip_data_vector, ip_flux_vector,
153 laplacian_matrix);
154 return;
155 }
156 }
157
159 ip_data_vector, ip_flux_vector, laplacian_matrix);
160 },
161 stabilizer);
162}
163} // 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