diff --git a/NumLib/DOF/DOFTableUtil.cpp b/NumLib/DOF/DOFTableUtil.cpp index 0fb28cc278b5c0aa9a56eaa6b6e4c597dd14bd78..bebed9bbbf6a7386cfcac4102e80c0190a7efa0a 100644 --- a/NumLib/DOF/DOFTableUtil.cpp +++ b/NumLib/DOF/DOFTableUtil.cpp @@ -12,6 +12,20 @@ namespace NumLib { +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) +{ + MeshLib::Location const l{mesh.getID(), MeshLib::MeshItemType::Node, + node_id}; + + auto const index = dof_table.getLocalIndex( + l, global_component_id, x.getRangeBegin(), x.getRangeEnd()); + + return x.get(index); +} + std::vector<GlobalIndexType> getIndices( std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const& dof_table) @@ -30,4 +44,81 @@ std::vector<GlobalIndexType> getIndices( return indices; } +double norm1(GlobalVector const& x, unsigned const global_component, + LocalToGlobalIndexMap const& dof_table, MeshLib::Mesh const& mesh) +{ + // TODO that also includes ghost nodes. + std::size_t const n_nodes = mesh.getNumberOfNodes(); + + double res = 0.0; + + for (std::size_t node_id = 0; node_id < n_nodes; ++node_id) { + auto const value = + getNodalValue(x, mesh, dof_table, node_id, global_component); + + res += std::abs(value); + } + + // TODO for PETSc some global accumulation is necessary. + return res; +} + +double norm2(GlobalVector const& x, unsigned const global_component, + LocalToGlobalIndexMap const& dof_table, MeshLib::Mesh const& mesh) +{ + // TODO that also includes ghost nodes. + std::size_t const n_nodes = mesh.getNumberOfNodes(); + + double res = 0.0; + + for (std::size_t node_id = 0; node_id < n_nodes; ++node_id) { + auto const value = + getNodalValue(x, mesh, dof_table, node_id, global_component); + + res += value*value; + } + + // TODO for PETSc some global accumulation is necessary. + res = std::sqrt(res); + return res; + +} + +double normMax(GlobalVector const& x, unsigned const global_component, + LocalToGlobalIndexMap const& dof_table, MeshLib::Mesh const& mesh) +{ + // TODO that also includes ghost nodes. + std::size_t const n_nodes = mesh.getNumberOfNodes(); + + double res = 0.0; + + for (std::size_t node_id = 0; node_id < n_nodes; ++node_id) { + auto const value = + getNodalValue(x, mesh, dof_table, node_id, global_component); + auto const abs_value = std::abs(value); + if (abs_value > res) + res = abs_value; + } + + // TODO for PETSc some global accumulation is necessary. + return res; +} + +double norm(GlobalVector const& x, unsigned const global_component, + MathLib::VecNormType norm_type, + LocalToGlobalIndexMap const& dof_table, MeshLib::Mesh const& mesh) +{ + switch(norm_type) + { + case MathLib::VecNormType::NORM1: + return norm1(x, global_component, dof_table, mesh); + case MathLib::VecNormType::NORM2: + return norm2(x, global_component, dof_table, mesh); + case MathLib::VecNormType::INFINITY_N: + return normMax(x, global_component, dof_table, mesh); + default: + OGS_FATAL("An invalid norm type has been passed."); + } +} + } // namespace NumLib diff --git a/NumLib/DOF/DOFTableUtil.h b/NumLib/DOF/DOFTableUtil.h index 546d614b67dd63d2f6a8614b969dfa9c5e627c3e..962079d5e4a63e7ed3dc9043b0c7186839758ff9 100644 --- a/NumLib/DOF/DOFTableUtil.h +++ b/NumLib/DOF/DOFTableUtil.h @@ -10,15 +10,52 @@ #ifndef NUMLIB_DOF_DOFTABLEUTIL_H #define NUMLIB_DOF_DOFTABLEUTIL_H +#include "MathLib/LinAlg/LinAlgEnums.h" #include "NumLib/DOF/LocalToGlobalIndexMap.h" namespace NumLib { +//! Returns the value for the given \c node_id and \c global_component_id from +//! the vector \c x. +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); + //! Returns nodal indices for the item identified by \c mesh_item_id from the //! given \c dof_table. std::vector<GlobalIndexType> getIndices( std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const& dof_table); + +//! Computes the specified norm of the given global component of the given vector x. +//! \remark +//! \c x is typically the solution vector of a monolithically coupled process +//! with several primary variables. +double norm(GlobalVector const& x, unsigned const global_component, + MathLib::VecNormType norm_type, + LocalToGlobalIndexMap const& dof_table, MeshLib::Mesh const& mesh); + +//! Computes the 1-norm of the given global component of the given vector x. +//! \remark +//! \c x is typically the solution vector of a monolithically coupled process +//! with several primary variables. +double norm1(GlobalVector const& x, unsigned const global_component, + LocalToGlobalIndexMap const& dof_table, MeshLib::Mesh const& mesh); + +//! Computes the 2-norm of the given global component of the given vector x. +//! \remark +//! \c x is typically the solution vector of a monolithically coupled process +//! with several primary variables. +double norm2(GlobalVector const& x, unsigned const global_component, + LocalToGlobalIndexMap const& dof_table, MeshLib::Mesh const& mesh); + +//! Computes the maximum norm of the given global component of the given vector x. +//! \remark +//! \c x is typically the solution vector of a monolithically coupled process +//! with several primary variables. +double normMax(GlobalVector const& x, unsigned const global_component, + LocalToGlobalIndexMap const& dof_table, MeshLib::Mesh const& mesh); } // namespace NumLib #endif // NUMLIB_DOF_DOFTABLEUTIL_H