Skip to content
Snippets Groups Projects
Commit cb64779e authored by Jonathan Youett's avatar Jonathan Youett
Browse files

Clean-up using range-based for-loops

parent 45220813
No related branches found
No related tags found
No related merge requests found
...@@ -28,10 +28,10 @@ class DOFConstraints { ...@@ -28,10 +28,10 @@ class DOFConstraints {
{ {
public: public:
LinearFactor() : index(0), factor(1.0) LinearFactor() : index(0), factor(1.0)
{}; {}
LinearFactor(int index, double factor) : index(index), factor(factor) LinearFactor(int index, double factor) : index(index), factor(factor)
{}; {}
int index; int index;
double factor; double factor;
...@@ -43,7 +43,7 @@ class DOFConstraints { ...@@ -43,7 +43,7 @@ class DOFConstraints {
//! create hanging node set //! create hanging node set
DOFConstraints() DOFConstraints()
{}; {}
/* /*
...@@ -52,7 +52,7 @@ class DOFConstraints { ...@@ -52,7 +52,7 @@ class DOFConstraints {
DOFConstraints(int size) DOFConstraints(int size)
{ {
resize(size); resize(size);
}; }
/* /*
...@@ -87,29 +87,22 @@ class DOFConstraints { ...@@ -87,29 +87,22 @@ class DOFConstraints {
void setupFromBasis(const B& basis) void setupFromBasis(const B& basis)
{ {
typedef typename B::GridView::template Codim<0>::Entity Element; typedef typename B::GridView::template Codim<0>::Entity Element;
typedef typename B::GridView::template Codim<0>::Iterator ElementIterator;
typedef typename B::GridView::IntersectionIterator NeighborIterator;
typedef typename B::LocalFiniteElement FE; typedef typename B::LocalFiniteElement FE;
typedef typename FE::Traits::LocalInterpolationType LocalInterpolationType;
typedef typename FE::Traits::LocalBasisType LocalBasisType;
typedef std::vector< std::map<int, double> > InterpolationMap; typedef std::vector< std::map<int, double> > InterpolationMap;
const typename B::GridView gridview = basis.getGridView(); const typename B::GridView& gridview = basis.getGridView();
resize(basis.size()); resize(basis.size());
InterpolationMap interpolationMap(basis.size()); InterpolationMap interpolationMap(basis.size());
ElementIterator it = gridview.template begin<0>(); for(const auto& element : elements(gridview))
ElementIterator end = gridview.template end<0>();
for(; it != end; ++it)
{ {
const auto& element = *it;
int level = element.level(); int level = element.level();
const LocalBasisType& localBasis = basis.getLocalFiniteElement(element).localBasis(); const auto& localBasis = basis.getLocalFiniteElement(element).localBasis();
const LocalInterpolationType& localInterpolation = basis.getLocalFiniteElement(element).localInterpolation(); const auto& localInterpolation = basis.getLocalFiniteElement(element).localInterpolation();
// check if all dofs on element are already constrained // check if all dofs on element are already constrained
bool allDofsConstrained = true; bool allDofsConstrained = true;
...@@ -119,19 +112,17 @@ class DOFConstraints { ...@@ -119,19 +112,17 @@ class DOFConstraints {
if (allDofsConstrained) if (allDofsConstrained)
continue; continue;
NeighborIterator nIt = gridview.ibegin(element); for(const auto& intersection : intersections(gridview, element))
NeighborIterator nEnd = gridview.iend(element);
for(; nIt != nEnd; ++nIt)
{ {
if (not(nIt->neighbor())) if (not(intersection.neighbor()))
continue; continue;
const Element& neighbor = nIt->outside(); const Element& neighbor = intersection.outside();
if (not(neighbor.level()<level)) if (not(neighbor.level()<level))
continue; continue;
const LocalBasisType& nLocalBasis = basis.getLocalFiniteElement(neighbor).localBasis(); const auto& nLocalBasis = basis.getLocalFiniteElement(neighbor).localBasis();
typedef typename Dune::LocalFiniteElementFunctionBase<FE>::type FunctionBaseClass; typedef typename Dune::LocalFiniteElementFunctionBase<FE>::type FunctionBaseClass;
...@@ -180,10 +171,8 @@ class DOFConstraints { ...@@ -180,10 +171,8 @@ class DOFConstraints {
// copy linear combinations // copy linear combinations
for(size_t i=0; i<interpolationMap.size(); ++i) for(size_t i=0; i<interpolationMap.size(); ++i)
{ {
InterpolationMap::value_type::const_iterator it = interpolationMap[i].begin(); for(const auto& it : interpolationMap[i])
InterpolationMap::value_type::const_iterator end = interpolationMap[i].end(); interpolation_[i].push_back(LinearFactor(it.first,it.second));
for(; it != end; ++it)
interpolation_[i].push_back(LinearFactor((*it).first,(*it).second));
} }
resolveDependencies(); resolveDependencies();
...@@ -250,14 +239,10 @@ class DOFConstraints { ...@@ -250,14 +239,10 @@ class DOFConstraints {
void report(int dofIndex) const void report(int dofIndex) const
{ {
typedef LinearCombination::const_iterator ConstIt;
std::cout << "x_" << dofIndex << " = "; std::cout << "x_" << dofIndex << " = ";
ConstIt it = interpolation_[dofIndex].begin(); for(const auto& it : interpolation_[dofIndex])
ConstIt end = interpolation_[dofIndex].end(); std::cout << it.factor << " * x_" << it.index << " + " ;
for(; it!=end; ++it)
std::cout << it->factor << " * x_" << it->index << " + " ;
std::cout << std::endl; std::cout << std::endl;
} }
...@@ -267,7 +252,6 @@ class DOFConstraints { ...@@ -267,7 +252,6 @@ class DOFConstraints {
*/ */
bool check() const bool check() const
{ {
typedef LinearCombination::const_iterator ConstIt;
bool r = true; bool r = true;
int size = isConstrained_.size(); int size = isConstrained_.size();
...@@ -278,14 +262,12 @@ class DOFConstraints { ...@@ -278,14 +262,12 @@ class DOFConstraints {
double sum = 0.0; double sum = 0.0;
bool errorForDof = false; bool errorForDof = false;
ConstIt it = interpolation_[i].begin(); for(const auto& it : interpolation_[i])
ConstIt end = interpolation_[i].end();
for(; it!=end; ++it)
{ {
sum += it->factor; sum += it.factor;
if (isConstrained_[it->index][0]) if (isConstrained_[it.index][0])
{ {
std::cout << "Error: interpolation for constrained DOF " << i << " depends on constrained DOF " << it->index << std::endl; std::cout << "Error: interpolation for constrained DOF " << i << " depends on constrained DOF " << it.index << std::endl;
errorForDof = true; errorForDof = true;
} }
} }
...@@ -329,12 +311,9 @@ class DOFConstraints { ...@@ -329,12 +311,9 @@ class DOFConstraints {
{ {
if (isConstrained_[i][0] and not(dependenciesResolved[i][0])) if (isConstrained_[i][0] and not(dependenciesResolved[i][0]))
{ {
typedef LinearCombination::const_iterator ConstIt;
typedef LinearCombination::iterator It;
LinearCombination inserted; LinearCombination inserted;
It it = interpolation_[i].begin(); auto it = interpolation_[i].begin();
while (it != interpolation_[i].end()) while (it != interpolation_[i].end())
{ {
// only factors which are hanging itself have to be inserted recursively // only factors which are hanging itself have to be inserted recursively
...@@ -347,10 +326,8 @@ class DOFConstraints { ...@@ -347,10 +326,8 @@ class DOFConstraints {
// now insert the factors of the constrained DOFs // now insert the factors of the constrained DOFs
// these dofs of the inserted are all nonhanging itself, since they have // these dofs of the inserted are all nonhanging itself, since they have
// already been resolved by the recursive resolveDOF above // already been resolved by the recursive resolveDOF above
ConstIt depIt = interpolation_[it->index].begin(); for (const auto& depIt : interpolation_[it->index])
ConstIt depEnd = interpolation_[it->index].end(); inserted.push_back(LinearFactor(depIt.index, it->factor * depIt.factor));
for (; depIt!=depEnd; ++depIt)
inserted.push_back(LinearFactor(depIt->index, it->factor * depIt->factor));
it = interpolation_[i].erase(it); it = interpolation_[i].erase(it);
} }
...@@ -367,5 +344,4 @@ class DOFConstraints { ...@@ -367,5 +344,4 @@ class DOFConstraints {
Interpolation interpolation_; Interpolation interpolation_;
}; };
#endif #endif
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment