OGS
DOFTableUtil.cpp
Go to the documentation of this file.
1
11#include "DOFTableUtil.h"
12
13#include <algorithm>
14#include <cassert>
15#include <functional>
16
18namespace NumLib
19{
20namespace
21{
22template <typename CalculateNorm>
23double norm(GlobalVector const& x, unsigned const global_component,
24 LocalToGlobalIndexMap const& dof_table, MeshLib::Mesh const& mesh,
25 CalculateNorm calculate_norm)
26{
27#ifdef USE_PETSC
28 x.setLocalAccessibleVector();
29#endif
30
31 double res = 0.0;
32 auto const& ms = dof_table.getMeshSubset(global_component);
33
34 assert(ms.getMeshID() == mesh.getID());
35 for (MeshLib::Node const* node : ms.getNodes())
36 {
37 auto const value = getNonGhostNodalValue(
38 x, mesh, dof_table, node->getID(), global_component);
39 res = calculate_norm(res, value);
40 }
41 return res;
42}
43
44double norm1(GlobalVector const& x, unsigned const global_component,
45 LocalToGlobalIndexMap const& dof_table, MeshLib::Mesh const& mesh)
46{
47 double res =
48 norm(x, global_component, dof_table, mesh,
49 [](double res, double value) { return res + std::abs(value); });
50
51#ifdef USE_PETSC
52 double global_result = 0.0;
53 MPI_Allreduce(&res, &global_result, 1, MPI_DOUBLE, MPI_SUM,
54 PETSC_COMM_WORLD);
55 res = global_result;
56#endif
57 return res;
58}
59
60double norm2(GlobalVector const& x, unsigned const global_component,
61 LocalToGlobalIndexMap const& dof_table, MeshLib::Mesh const& mesh)
62{
63 double res =
64 norm(x, global_component, dof_table, mesh,
65 [](double res, double value) { return res + value * value; });
66
67#ifdef USE_PETSC
68 double global_result = 0.0;
69 MPI_Allreduce(&res, &global_result, 1, MPI_DOUBLE, MPI_SUM,
70 PETSC_COMM_WORLD);
71 res = global_result;
72#endif
73 return std::sqrt(res);
74}
75
76double normInfinity(GlobalVector const& x, unsigned const global_component,
77 LocalToGlobalIndexMap const& dof_table,
78 MeshLib::Mesh const& mesh)
79{
80 double res = norm(x, global_component, dof_table, mesh,
81 [](double res, double value)
82 { return std::max(res, std::abs(value)); });
83
84#ifdef USE_PETSC
85 double global_result = 0.0;
86 MPI_Allreduce(&res, &global_result, 1, MPI_DOUBLE, MPI_MAX,
87 PETSC_COMM_WORLD);
88 res = global_result;
89#endif
90 return res;
91}
92} // anonymous namespace
93
95 NumLib::LocalToGlobalIndexMap const& dof_table,
96 std::size_t const node_id,
97 std::size_t const global_component_id)
98{
100 node_id};
101
102 auto const index = dof_table.getGlobalIndex(l, global_component_id);
103 assert(index != NumLib::MeshComponentMap::nop);
104
105 if (index < 0)
106 { // ghost node value
107 return 0.0;
108 }
109
110 return x.get(index);
111}
112
113double getNodalValue(GlobalVector const& x, MeshLib::Mesh const& mesh,
114 NumLib::LocalToGlobalIndexMap const& dof_table,
115 std::size_t const node_id,
116 std::size_t const global_component_id)
117{
119 node_id};
120
121 auto const index = dof_table.getGlobalIndex(l, global_component_id);
122 assert(index != NumLib::MeshComponentMap::nop);
123
124 return x.get(index);
125}
126
127std::vector<GlobalIndexType> getIndices(
128 std::size_t const mesh_item_id,
129 NumLib::LocalToGlobalIndexMap const& dof_table)
130{
131 assert(dof_table.size() > mesh_item_id);
132 std::vector<GlobalIndexType> indices;
133
134 // Local matrices and vectors will always be ordered by component
135 // no matter what the order of the global matrix is.
136 for (int c = 0; c < dof_table.getNumberOfGlobalComponents(); ++c)
137 {
138 auto const& idcs = dof_table(mesh_item_id, c).rows;
139 indices.reserve(indices.size() + idcs.size());
140 indices.insert(indices.end(), idcs.begin(), idcs.end());
141 }
142
143 return indices;
144}
145
147 std::size_t const id,
148 NumLib::LocalToGlobalIndexMap const& dof_table,
149 std::vector<GlobalIndexType>& indices)
150{
151 assert(dof_table.size() > id);
152 indices.clear();
153
154 // Local matrices and vectors will always be ordered by component,
155 // no matter what the order of the global matrix is.
156 for (int c = 0; c < dof_table.getNumberOfGlobalComponents(); ++c)
157 {
158 auto const& idcs = dof_table(id, c).rows;
159 indices.reserve(indices.size() + idcs.size());
160 indices.insert(indices.end(), idcs.begin(), idcs.end());
161 }
162
164}
165
166double norm(GlobalVector const& x, unsigned const global_component,
167 MathLib::VecNormType norm_type,
168 LocalToGlobalIndexMap const& dof_table, MeshLib::Mesh const& mesh)
169{
170 switch (norm_type)
171 {
173 return norm1(x, global_component, dof_table, mesh);
175 return norm2(x, global_component, dof_table, mesh);
177 return normInfinity(x, global_component, dof_table, mesh);
178 default:
179 OGS_FATAL("An invalid norm type has been passed.");
180 }
181}
182
183Eigen::VectorXd getLocalX(
184 std::size_t const mesh_item_id,
185 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_tables,
186 std::vector<GlobalVector*> const& x)
187{
188 Eigen::VectorXd local_x_vec;
189
190 auto const n_processes = x.size();
191 for (std::size_t process_id = 0; process_id < n_processes; ++process_id)
192 {
193 auto const indices =
194 NumLib::getIndices(mesh_item_id, *dof_tables[process_id]);
195 assert(!indices.empty());
196 auto const last = local_x_vec.size();
197 local_x_vec.conservativeResize(last + indices.size());
198 auto const local_solution = x[process_id]->get(indices);
199 assert(indices.size() == local_solution.size());
200 local_x_vec.tail(local_solution.size()).noalias() =
201 MathLib::toVector(local_solution);
202 }
203 return local_x_vec;
204}
205
206} // namespace NumLib
#define OGS_FATAL(...)
Definition Error.h:26
Global vector based on Eigen vector.
Definition EigenVector.h:25
double get(IndexType rowId) const
get entry
Definition EigenVector.h:58
std::size_t getID() const
Get id of the mesh.
Definition Mesh.h:121
MathLib::RowColumnIndices< GlobalIndexType > RowColumnIndices
GlobalIndexType getGlobalIndex(MeshLib::Location const &l, int const variable_id, int const component_id) const
MeshLib::MeshSubset const & getMeshSubset(int const variable_id, int const component_id) const
static constexpr NUMLIB_EXPORT GlobalIndexType const nop
Eigen::Map< const Vector > toVector(std::vector< double > const &data, Eigen::VectorXd::Index size)
Creates an Eigen mapped vector from the given data vector.
double normInfinity(GlobalVector const &x, unsigned const global_component, LocalToGlobalIndexMap const &dof_table, MeshLib::Mesh const &mesh)
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)
Eigen::VectorXd getLocalX(std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x)