OGS
DOFTableUtil.h
Go to the documentation of this file.
1 
11 #pragma once
12 
13 #include "MathLib/LinAlg/LinAlg.h"
16 
17 namespace NumLib
18 {
21 double getNonGhostNodalValue(GlobalVector const& x, MeshLib::Mesh const& mesh,
22  NumLib::LocalToGlobalIndexMap const& dof_table,
23  std::size_t const node_id,
24  std::size_t const global_component_id);
25 
28 double getNodalValue(GlobalVector const& x, MeshLib::Mesh const& mesh,
29  NumLib::LocalToGlobalIndexMap const& dof_table,
30  std::size_t const node_id,
31  std::size_t const global_component_id);
32 
35 std::vector<GlobalIndexType> getIndices(
36  std::size_t const mesh_item_id,
37  NumLib::LocalToGlobalIndexMap const& dof_table);
38 
42  std::size_t const id,
43  NumLib::LocalToGlobalIndexMap const& dof_table,
44  std::vector<GlobalIndexType>& indices);
45 
49 double norm(GlobalVector const& x, unsigned const global_component,
50  MathLib::VecNormType norm_type,
51  LocalToGlobalIndexMap const& dof_table, MeshLib::Mesh const& mesh);
52 
58 template <typename Functor>
60  GlobalVector const& input_vector, int const variable_id,
61  NumLib::LocalToGlobalIndexMap const& local_to_global_index_map,
62  MeshLib::PropertyVector<double>& output_vector, Functor mapFunction)
63 {
65 
66  std::fill(begin(output_vector), end(output_vector),
67  std::numeric_limits<double>::quiet_NaN());
68 
69  int const n_components =
70  local_to_global_index_map.getNumberOfVariableComponents(variable_id);
71  for (int component = 0; component < n_components; ++component)
72  {
73  auto const& mesh_subset =
74  local_to_global_index_map.getMeshSubset(variable_id, component);
75  auto const mesh_id = mesh_subset.getMeshID();
76  for (auto const& node : mesh_subset.getNodes())
77  {
78  auto const node_id = node->getID();
80  node_id);
81  output_vector.getComponent(node_id, component) = mapFunction(
82  input_vector[local_to_global_index_map.getGlobalIndex(
83  l, variable_id, component)]);
84  }
85  }
86 }
87 } // namespace NumLib
Global vector based on Eigen vector.
Definition: EigenVector.h:26
std::size_t getMeshID() const
return this mesh ID
Definition: MeshSubset.h:71
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
void setLocalAccessibleVector(PETScVector const &x)
Definition: LinAlg.cpp:27
VecNormType
Norm type. Not declared as class type in order to use the members as integers.
Definition: LinAlgEnums.h:22
double getNodalValue(GlobalVector const &x, MeshLib::Mesh const &mesh, NumLib::LocalToGlobalIndexMap const &dof_table, std::size_t const node_id, std::size_t const global_component_id)
double norm(GlobalVector const &x, unsigned const global_component, MathLib::VecNormType norm_type, LocalToGlobalIndexMap const &dof_table, MeshLib::Mesh const &mesh)
double getNonGhostNodalValue(GlobalVector const &x, MeshLib::Mesh const &mesh, NumLib::LocalToGlobalIndexMap const &dof_table, std::size_t const node_id, 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)
void transformVariableFromGlobalVector(GlobalVector const &input_vector, int const variable_id, NumLib::LocalToGlobalIndexMap const &local_to_global_index_map, MeshLib::PropertyVector< double > &output_vector, Functor mapFunction)
Definition: DOFTableUtil.h:59