Loading [MathJax]/extensions/MathMenu.js
OGS
DOFTableUtil.h
Go to the documentation of this file.
1
11#pragma once
12
13#include <Eigen/Core>
14
18
19namespace NumLib
20{
24 MeshLib::Location const& location,
25 NumLib::LocalToGlobalIndexMap const& dof_table,
26 std::size_t const global_component_id);
27
30double getNodalValue(GlobalVector const& x, MeshLib::Location const& location,
31 NumLib::LocalToGlobalIndexMap const& dof_table,
32 std::size_t const global_component_id);
33
36std::vector<GlobalIndexType> getIndices(
37 std::size_t const mesh_item_id,
38 NumLib::LocalToGlobalIndexMap const& dof_table);
39
43 std::size_t const id,
44 NumLib::LocalToGlobalIndexMap const& dof_table,
45 std::vector<GlobalIndexType>& indices);
46
50double norm(GlobalVector const& x, unsigned const global_component,
51 MathLib::VecNormType norm_type,
52 LocalToGlobalIndexMap const& dof_table);
53
59template <typename Functor>
61 GlobalVector const& input_vector, int const variable_id,
62 NumLib::LocalToGlobalIndexMap const& local_to_global_index_map,
63 MeshLib::PropertyVector<double>& output_vector, Functor map_function)
64{
66
67 // We fill the output with zeroes, because filling with NaN
68 // "breaks" visualization, e.g. of heat or mass flow rate
69 // with Taylor-Hood elements. I.e., all lower order
70 // properties would have NaN on the higher order nodes.
71 std::fill(begin(output_vector), end(output_vector), 0);
72
73 int const n_components =
74 local_to_global_index_map.getNumberOfVariableComponents(variable_id);
75 for (int component = 0; component < n_components; ++component)
76 {
77 auto const& mesh_subset =
78 local_to_global_index_map.getMeshSubset(variable_id, component);
79 for (auto const& l : MeshLib::views::meshLocations(
80 mesh_subset, MeshLib::MeshItemType::Node))
81 {
82 auto const node_id = l.item_id;
83 auto const input_index = local_to_global_index_map.getGlobalIndex(
84 l, variable_id, component);
85 double const value = input_vector[input_index];
86
87 output_vector.getComponent(node_id, component) =
88 map_function(value);
89 }
90 }
91}
92
95template <typename Functor>
97 GlobalVector const& input_vector_on_bulk_mesh, int const variable_id,
98 NumLib::LocalToGlobalIndexMap const& local_to_global_index_map,
99 MeshLib::PropertyVector<double>& output_vector_on_submesh,
100 std::vector<std::size_t> const& map_submesh_node_id_to_bulk_mesh_node_id,
101 Functor map_function)
102{
103 MathLib::LinAlg::setLocalAccessibleVector(input_vector_on_bulk_mesh);
104
105 // We fill the output with zeroes, because filling with NaN
106 // "breaks" visualization, e.g. of heat or mass flow rate
107 // with Taylor-Hood elements. I.e., all lower order
108 // properties would have NaN on the higher order nodes.
109 std::fill(begin(output_vector_on_submesh), end(output_vector_on_submesh),
110 0);
111
112 int const n_components =
113 local_to_global_index_map.getNumberOfVariableComponents(variable_id);
114 for (int component = 0; component < n_components; ++component)
115 {
116 auto const& mesh_subset =
117 local_to_global_index_map.getMeshSubset(variable_id, component);
118 auto const mesh_id = mesh_subset.getMeshID();
119
120 for (std::size_t submesh_node_id = 0;
121 submesh_node_id < map_submesh_node_id_to_bulk_mesh_node_id.size();
122 ++submesh_node_id)
123 {
124 std::size_t const bulk_node_id =
125 map_submesh_node_id_to_bulk_mesh_node_id[submesh_node_id];
127 bulk_node_id);
128 auto const input_index = local_to_global_index_map.getGlobalIndex(
129 l, variable_id, component);
130
131 // On higher-order meshes some d.o.f.s might not be defined on all
132 // nodes. We silently ignore the d.o.f.s we cannot find.
133 [[unlikely]] if (input_index == NumLib::MeshComponentMap::nop)
134 {
135 continue;
136 }
137
138 double const value = input_vector_on_bulk_mesh[input_index];
139
140 output_vector_on_submesh.getComponent(submesh_node_id, component) =
141 map_function(value);
142 }
143 }
144}
145
146Eigen::VectorXd getLocalX(
147 std::size_t const mesh_item_id,
148 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_tables,
149 std::vector<GlobalVector*> const& x);
150} // namespace NumLib
Global vector based on Eigen vector.
Definition EigenVector.h:25
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)