Loading [MathJax]/extensions/tex2jax.js
OGS
DOFTableUtil.h
Go to the documentation of this file.
1
11#pragma once
12
13#include <Eigen/Core>
14#include <range/v3/algorithm/fill.hpp>
15
19
20namespace NumLib
21{
25 MeshLib::Location const& location,
26 NumLib::LocalToGlobalIndexMap const& dof_table,
27 std::size_t const global_component_id);
28
31double getNodalValue(GlobalVector const& x, MeshLib::Location const& location,
32 NumLib::LocalToGlobalIndexMap const& dof_table,
33 std::size_t const global_component_id);
34
37std::vector<GlobalIndexType> getIndices(
38 std::size_t const mesh_item_id,
39 NumLib::LocalToGlobalIndexMap const& dof_table);
40
43std::vector<GlobalIndexType> getIndices(
44 std::size_t const variable_id,
45 std::size_t const component_id,
46 std::size_t const bulk_element_id,
47 NumLib::LocalToGlobalIndexMap const& dof_table);
48
52 std::size_t const id,
53 NumLib::LocalToGlobalIndexMap const& dof_table,
54 std::vector<GlobalIndexType>& indices);
55
59double norm(GlobalVector const& x, unsigned const global_component,
60 MathLib::VecNormType norm_type,
61 LocalToGlobalIndexMap const& dof_table);
62
68template <typename Functor>
70 GlobalVector const& input_vector, int const variable_id,
71 NumLib::LocalToGlobalIndexMap const& local_to_global_index_map,
72 MeshLib::PropertyVector<double>& output_vector, Functor map_function)
73{
75
76 // We fill the output with zeroes, because filling with NaN
77 // "breaks" visualization, e.g. of heat or mass flow rate
78 // with Taylor-Hood elements. I.e., all lower order
79 // properties would have NaN on the higher order nodes.
80 ranges::fill(output_vector, 0);
81
82 int const n_components =
83 local_to_global_index_map.getNumberOfVariableComponents(variable_id);
84 for (int component = 0; component < n_components; ++component)
85 {
86 auto const& mesh_subset =
87 local_to_global_index_map.getMeshSubset(variable_id, component);
88 for (auto const& l : MeshLib::views::meshLocations(
89 mesh_subset, MeshLib::MeshItemType::Node))
90 {
91 auto const node_id = l.item_id;
92 auto const input_index = local_to_global_index_map.getGlobalIndex(
93 l, variable_id, component);
94 double const value = input_vector[input_index];
95
96 output_vector.getComponent(node_id, component) =
97 map_function(value);
98 }
99 }
100}
101
104template <typename Functor>
106 GlobalVector const& input_vector_on_bulk_mesh, int const variable_id,
107 NumLib::LocalToGlobalIndexMap const& local_to_global_index_map,
108 MeshLib::PropertyVector<double>& output_vector_on_submesh,
109 std::span<std::size_t const> const&
110 map_submesh_node_id_to_bulk_mesh_node_id,
111 Functor map_function)
112{
113 MathLib::LinAlg::setLocalAccessibleVector(input_vector_on_bulk_mesh);
114
115 // We fill the output with zeroes, because filling with NaN
116 // "breaks" visualization, e.g. of heat or mass flow rate
117 // with Taylor-Hood elements. I.e., all lower order
118 // properties would have NaN on the higher order nodes.
119 ranges::fill(output_vector_on_submesh, 0);
120
121 int const n_components =
122 local_to_global_index_map.getNumberOfVariableComponents(variable_id);
123 for (int component = 0; component < n_components; ++component)
124 {
125 auto const& mesh_subset =
126 local_to_global_index_map.getMeshSubset(variable_id, component);
127 auto const mesh_id = mesh_subset.getMeshID();
128
129 for (std::size_t submesh_node_id = 0;
130 submesh_node_id < map_submesh_node_id_to_bulk_mesh_node_id.size();
131 ++submesh_node_id)
132 {
133 std::size_t const bulk_node_id =
134 map_submesh_node_id_to_bulk_mesh_node_id[submesh_node_id];
136 bulk_node_id);
137 auto const input_index = local_to_global_index_map.getGlobalIndex(
138 l, variable_id, component);
139
140 // On higher-order meshes some d.o.f.s might not be defined on all
141 // nodes. We silently ignore the d.o.f.s we cannot find.
142 [[unlikely]] if (input_index == NumLib::MeshComponentMap::nop)
143 {
144 continue;
145 }
146
147 double const value = input_vector_on_bulk_mesh[input_index];
148
149 output_vector_on_submesh.getComponent(submesh_node_id, component) =
150 map_function(value);
151 }
152 }
153}
154
155Eigen::VectorXd getLocalX(
156 std::size_t const mesh_item_id,
157 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_tables,
158 std::vector<GlobalVector*> const& x);
159} // namespace NumLib
Global vector based on Eigen vector.
Definition EigenVector.h:26
std::size_t getMeshID() const
return this mesh ID
Definition MeshSubset.h:76
PROP_VAL_TYPE & getComponent(std::size_t tuple_index, int component)
Returns the value for the given component stored in the given tuple.
MathLib::RowColumnIndices< GlobalIndexType > RowColumnIndices
GlobalIndexType getGlobalIndex(MeshLib::Location const &l, int const variable_id, int const component_id) const
int getNumberOfVariableComponents(int variable_id) const
MeshLib::MeshSubset const & getMeshSubset(int const variable_id, int const component_id) const
static constexpr NUMLIB_EXPORT GlobalIndexType const nop
void setLocalAccessibleVector(PETScVector const &x)
Definition LinAlg.cpp:27
auto meshLocations(Mesh const &mesh, MeshItemType const item_type)
Definition Mesh.h:238
void transformVariableFromGlobalVector(GlobalVector const &input_vector, int const variable_id, NumLib::LocalToGlobalIndexMap const &local_to_global_index_map, MeshLib::PropertyVector< double > &output_vector, Functor map_function)
double getNodalValue(GlobalVector const &x, MeshLib::Location const &location, NumLib::LocalToGlobalIndexMap const &dof_table, std::size_t const global_component_id)
std::vector< GlobalIndexType > getIndices(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
NumLib::LocalToGlobalIndexMap::RowColumnIndices getRowColumnIndices(std::size_t const id, NumLib::LocalToGlobalIndexMap const &dof_table, std::vector< GlobalIndexType > &indices)
Eigen::VectorXd getLocalX(std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x)
double getNonGhostNodalValue(GlobalVector const &x, MeshLib::Location const &location, NumLib::LocalToGlobalIndexMap const &dof_table, std::size_t const global_component_id)
double norm(GlobalVector const &x, unsigned const global_component, MathLib::VecNormType norm_type, LocalToGlobalIndexMap const &dof_table)