Skip to content
Snippets Groups Projects
Commit 20b821be authored by Christoph Grüninger's avatar Christoph Grüninger Committed by Christoph Grüninger
Browse files

[cleanup] Replace SLList by std::forward_list

parent 2e71688e
No related branches found
No related tags found
1 merge request!8[cleanup] Replace SLList by std::forward_list
......@@ -9,8 +9,9 @@
#include <dune/grid/uggrid.hh>
/** \todo Remove the following include once getAllSubfaces... is gone */
#include <dune/common/sllist.hh>
/** \todo Remove the following two includes once getAllSubfaces... is gone */
#include <forward_list>
#include <iterator>
#include <dune/common/stdstreams.hh>
#include <dune/grid/common/mcmgmapper.hh>
......@@ -420,7 +421,7 @@ Dune::UGGrid<dim>::getChildrenOfSubface(const T & e,
typedef std::pair<typename UG_NS<dim>::Element*,int> ListEntryType;
SLList<ListEntryType> list;
std::forward_list<ListEntryType> list;
// //////////////////////////////////////////////////////////////////////
// Change the input face number from Dune numbering to UG numbering
......@@ -451,7 +452,7 @@ Dune::UGGrid<dim>::getChildrenOfSubface(const T & e,
DUNE_THROW(GridError, "Get_Sons_of_ElementSide returned with error value " << rv);
for (int i=0; i<Sons_of_Side; i++)
list.push_back(ListEntryType(SonList[i],SonSides[i]));
list.push_front(ListEntryType(SonList[i],SonSides[i]));
}
......@@ -459,11 +460,10 @@ Dune::UGGrid<dim>::getChildrenOfSubface(const T & e,
// Traverse and collect all children of the side
// //////////////////////////////////////////////////
typename SLList<ListEntryType>::iterator f = list.begin();
for (; f!=list.end(); ++f) {
typename UG_NS<dim>::Element* theElement = f->first;
int side = f->second;
for (const auto& f : list)
{
typename UG_NS<dim>::Element* theElement = f.first;
int side = f.second;
int Sons_of_Side = 0;
typename UG_NS<dim>::Element* SonList[UG_NS<dim>::MAX_SONS];
......@@ -480,7 +480,7 @@ Dune::UGGrid<dim>::getChildrenOfSubface(const T & e,
true);
for (int i=0; i<Sons_of_Side; i++)
list.push_back(ListEntryType(SonList[i],SonSides[i]));
list.push_front(ListEntryType(SonList[i],SonSides[i]));
}
......@@ -492,23 +492,22 @@ Dune::UGGrid<dim>::getChildrenOfSubface(const T & e,
// Use reserve / push_back since EntityPointer is not default constructable
childElements.clear();
childElements.reserve( list.size() );
childElementSides.resize(list.size());
childElements.reserve(std::distance(list.begin(), list.end()));
childElementSides.resize(childElements.size());
int i=0;
for (f = list.begin(); f!=list.end(); ++f, ++i)
for (const auto &f : list)
{
// Set element
typedef typename Traits::template Codim< 0 >::EntityPointer EntityPointer;
childElements.push_back( EntityPointer( UGGridEntityPointer< 0, const UGGrid< dim > >( f->first, this ) ) );
childElements.push_back( EntityPointer( UGGridEntityPointer< 0, const UGGrid< dim > >( f.first, this ) ) );
int side = f->second;
int side = f.second;
// Dune numbers the faces of several elements differently than UG.
// The following switch does the transformation
childElementSides[i] = UGGridRenumberer<dim>::facesUGtoDUNE(side, childElements[i]->type());
++i;
}
}
......@@ -523,7 +522,7 @@ void Dune::UGGrid<dim>::getChildrenOfSubface(const typename Traits::template Cod
typedef std::pair<typename UG_NS<dim>::Element*,int> ListEntryType;
SLList<ListEntryType> list;
std::forward_list<ListEntryType> list;
// //////////////////////////////////////////////////////////////////////
// Change the input face number from Dune numbering to UG numbering
......@@ -554,7 +553,7 @@ void Dune::UGGrid<dim>::getChildrenOfSubface(const typename Traits::template Cod
DUNE_THROW(GridError, "Get_Sons_of_ElementSide returned with error value " << rv);
for (int i=0; i<Sons_of_Side; i++)
list.push_back(ListEntryType(SonList[i],SonSides[i]));
list.push_front(ListEntryType(SonList[i],SonSides[i]));
}
......@@ -562,11 +561,10 @@ void Dune::UGGrid<dim>::getChildrenOfSubface(const typename Traits::template Cod
// Traverse and collect all children of the side
// //////////////////////////////////////////////////
typename SLList<ListEntryType>::iterator f = list.begin();
for (; f!=list.end(); ++f) {
typename UG_NS<dim>::Element* theElement = f->first;
int side = f->second;
for (const auto& f : list)
{
typename UG_NS<dim>::Element* theElement = f.first;
int side = f.second;
UG::INT Sons_of_Side = 0;
typename UG_NS<dim>::Element* SonList[UG_NS<dim>::MAX_SONS];
......@@ -583,7 +581,7 @@ void Dune::UGGrid<dim>::getChildrenOfSubface(const typename Traits::template Cod
true);
for (int i=0; i<Sons_of_Side; i++)
list.push_back(ListEntryType(SonList[i],SonSides[i]));
list.push_front(ListEntryType(SonList[i],SonSides[i]));
}
......@@ -595,23 +593,22 @@ void Dune::UGGrid<dim>::getChildrenOfSubface(const typename Traits::template Cod
// Use reserve / push_back since EntityPointer is not default constructable
childElements.clear();
childElements.reserve( list.size() );
childElementSides.resize(list.size());
childElements.reserve(std::distance(list.begin(), list.end()));
childElementSides.resize(childElements.size());
int i=0;
for (f = list.begin(); f!=list.end(); ++f, ++i)
for (const auto f : list)
{
// Set element
typedef typename Traits::template Codim< 0 >::Entity Entity;
childElements.push_back( Entity( UGGridEntity< 0, dim, const UGGrid< dim > >( f->first, this ) ) );
childElements.push_back( Entity( UGGridEntity< 0, dim, const UGGrid< dim > >( f.first, this ) ) );
int side = f->second;
int side = f.second;
// Dune numbers the faces of several elements differently than UG.
// The following switch does the transformation
childElementSides[i] = UGGridRenumberer<dim>::facesUGtoDUNE(side, childElements[i].type());
++i;
}
}
......
......@@ -6,7 +6,7 @@
#include <dune/grid/uggrid/uggridintersections.hh>
#include <set>
#include <forward_list>
template<class GridImp>
const typename Dune::UGGridLevelIntersection<GridImp>::WorldVector&
......@@ -696,7 +696,7 @@ void Dune::UGGridLeafIntersection<GridImp>::constructLeafSubfaces() {
else {
SLList<Face> list;
std::forward_list<Face> list;
int levelNeighborSide = numberInNeighbor(center_, levelNeighbor);
UG::INT Sons_of_Side = 0;
......@@ -716,7 +716,7 @@ void Dune::UGGridLeafIntersection<GridImp>::constructLeafSubfaces() {
DUNE_THROW(GridError, "Get_Sons_of_ElementSide returned with error value " << rv);
for (int i=0; i<Sons_of_Side; i++)
list.push_back(Face(SonList[i],SonSides[i]));
list.push_front(Face(SonList[i],SonSides[i]));
// //////////////////////////////////////////////////
// Get_Sons_of_ElementSide only computes direct sons. Therefore in order to get all
......@@ -724,10 +724,9 @@ void Dune::UGGridLeafIntersection<GridImp>::constructLeafSubfaces() {
// We do a breadth-first search.
// //////////////////////////////////////////////////
typename SLList<Face>::iterator f = list.begin();
for (; f!=list.end(); ++f) {
const typename UG_NS<dim>::Element* theElement = f->first;
for (const auto& f : list)
{
const typename UG_NS<dim>::Element* theElement = f.first;
UG::INT Sons_of_Side = 0;
typename UG_NS<dim>::Element* SonList[UG_NS<dim>::MAX_SONS];
......@@ -736,7 +735,7 @@ void Dune::UGGridLeafIntersection<GridImp>::constructLeafSubfaces() {
if (!UG_NS<dim>::isLeaf(theElement)) {
Get_Sons_of_ElementSide(theElement,
f->second, // Input element side number
f.second, // Input element side number
&Sons_of_Side, // Number of topological sons of the element side
SonList, // Output elements
SonSides, // Output element side numbers
......@@ -745,7 +744,7 @@ void Dune::UGGridLeafIntersection<GridImp>::constructLeafSubfaces() {
true);
for (int i=0; i<Sons_of_Side; i++)
list.push_back(Face(SonList[i],SonSides[i]));
list.push_front(Face(SonList[i],SonSides[i]));
}
......@@ -756,11 +755,11 @@ void Dune::UGGridLeafIntersection<GridImp>::constructLeafSubfaces() {
// //////////////////////////////
leafSubFaces_.resize(0);
for (f = list.begin(); f!=list.end(); ++f) {
for (const auto& f : list)
{
// Set element
if (UG_NS<dim>::isLeaf(f->first))
leafSubFaces_.push_back(*f);
if (UG_NS<dim>::isLeaf(f.first))
leafSubFaces_.push_back(f);
}
}
......
......@@ -5,8 +5,6 @@
#include <memory>
#include <dune/common/sllist.hh>
#include <dune/grid/uggrid/uggridrenumberer.hh>
/** \file
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment