version 3.11-dev
Loading...
Searching...
No Matches
localdof.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightText: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_CVFE_LOCAL_DOF_HH
13#define DUMUX_CVFE_LOCAL_DOF_HH
14
15#include <ranges>
16
18
19namespace Dumux::CVFE {
20
25template<class LocalIndex, class GridIndex>
27{
28public:
29 using LocalIndexType = LocalIndex;
30 using GridIndexType = GridIndex;
31
32 LocalDof(LocalIndex index, GridIndex dofIndex, GridIndex eIdx)
33 : index_(index), dofIndex_(dofIndex), eIdx_(eIdx) {}
34 LocalIndex index() const { return index_; }
35 GridIndex dofIndex() const { return dofIndex_; }
36 GridIndex elementIndex() const { return eIdx_; }
37private:
38 LocalIndex index_;
39 GridIndex dofIndex_;
40 GridIndex eIdx_;
41};
42
43
44} // end namespace Dumux::CVFE
45
46namespace Dumux {
47
49template<class FVElementGeometry>
50inline auto localDofs(const FVElementGeometry& fvGeometry)
51{
54
55 auto scvRange = scvs(fvGeometry);
56 return std::ranges::transform_view(
57 std::ranges::subrange(scvRange.begin(), scvRange.end()),
58 [](const auto& scv) {
59 return CVFE::LocalDof<LocalIndexType, GridIndexType>{
60 static_cast<LocalIndexType>(scv.localDofIndex()),
61 static_cast<GridIndexType>(scv.dofIndex()),
62 static_cast<GridIndexType>(scv.elementIndex())
63 };
64 }
65 );
66}
67
69template<class FVElementGeometry>
70inline auto cvLocalDofs(const FVElementGeometry& fvGeometry)
71{
72 // As default all dofs are cv dofs
73 return localDofs(fvGeometry);
74}
75
80template<class FVElementGeometry, class LocalDof>
81inline std::ranges::range auto
82scvs(const FVElementGeometry& fvGeometry, const LocalDof& localDof)
83{
84 assert(fvGeometry.numScv() > localDof.index());
85 return std::views::single(1) | std::views::transform(
86 [&](const auto i) -> const typename FVElementGeometry::SubControlVolume& { return fvGeometry.scv(localDof.index()); }
87 );
88}
89
90} // end namespace Dumux
91#endif
GridIndex GridIndexType
Definition localdof.hh:30
GridIndex dofIndex() const
Definition localdof.hh:35
LocalIndex index() const
Definition localdof.hh:34
GridIndex elementIndex() const
Definition localdof.hh:36
LocalDof(LocalIndex index, GridIndex dofIndex, GridIndex eIdx)
Definition localdof.hh:32
LocalIndex LocalIndexType
Definition localdof.hh:29
Defines the index types used for grid and local indices.
Definition interpolate.hh:27
Definition adapt.hh:17
auto cvLocalDofs(const FVElementGeometry &fvGeometry)
range over control-volume local dofs
Definition localdof.hh:70
std::ranges::range auto scvs(const FVElementGeometry &fvGeometry, const LocalDof &localDof)
Definition localdof.hh:82
auto localDofs(const FVElementGeometry &fvGeometry)
range over local dofs
Definition localdof.hh:50
typename GridView::IndexSet::IndexType GridIndex
Definition indextraits.hh:27
unsigned int LocalIndex
Definition indextraits.hh:28