22template <
typename CalculateNorm>
25 CalculateNorm calculate_norm)
28 x.setLocalAccessibleVector();
34 assert(ms.getMeshID() == mesh.
getID());
38 x, mesh, dof_table, node->getID(), global_component);
39 res = calculate_norm(res, value);
44double norm1(
GlobalVector const& x,
unsigned const global_component,
48 norm(x, global_component, dof_table, mesh,
49 [](
double res,
double value) {
return res + std::abs(value); });
52 double global_result = 0.0;
53 MPI_Allreduce(&res, &global_result, 1, MPI_DOUBLE, MPI_SUM,
60double norm2(
GlobalVector const& x,
unsigned const global_component,
64 norm(x, global_component, dof_table, mesh,
65 [](
double res,
double value) {
return res + value * value; });
68 double global_result = 0.0;
69 MPI_Allreduce(&res, &global_result, 1, MPI_DOUBLE, MPI_SUM,
73 return std::sqrt(res);
80 double res =
norm(x, global_component, dof_table, mesh,
81 [](
double res,
double value)
82 {
return std::max(res, std::abs(value)); });
85 double global_result = 0.0;
86 MPI_Allreduce(&res, &global_result, 1, MPI_DOUBLE, MPI_MAX,
96 std::size_t
const node_id,
97 std::size_t
const global_component_id)
102 auto const index = dof_table.
getGlobalIndex(l, global_component_id);
115 std::size_t
const node_id,
116 std::size_t
const global_component_id)
121 auto const index = dof_table.
getGlobalIndex(l, global_component_id);
128 std::size_t
const mesh_item_id,
131 assert(dof_table.
size() > mesh_item_id);
132 std::vector<GlobalIndexType> indices;
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());
147 std::size_t
const id,
149 std::vector<GlobalIndexType>& indices)
151 assert(dof_table.
size() >
id);
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());
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);
179 OGS_FATAL(
"An invalid norm type has been passed.");
184 std::size_t
const mesh_item_id,
185 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& dof_tables,
186 std::vector<GlobalVector*>
const& x)
188 Eigen::VectorXd local_x_vec;
190 auto const n_processes = x.size();
191 for (std::size_t process_id = 0; process_id < n_processes; ++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() =
Global vector based on Eigen vector.
double get(IndexType rowId) const
get entry
std::size_t getID() const
Get id of the mesh.
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
int getNumberOfGlobalComponents() 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)