Skip to content
GitLab
Explore
Sign in
Primary navigation
Search or go to…
Project
D
dune-istl
Manage
Activity
Members
Labels
Plan
Issues
Issue boards
Milestones
Wiki
Code
Merge requests
Repository
Branches
Commits
Tags
Repository graph
Compare revisions
Build
Pipelines
Jobs
Pipeline schedules
Artifacts
Deploy
Releases
Model registry
Operate
Environments
Monitor
Incidents
Service Desk
Analyze
Value stream analytics
Contributor analytics
CI/CD analytics
Repository analytics
Model experiments
Help
Help
Support
GitLab documentation
Compare GitLab plans
Community forum
Contribute to GitLab
Provide feedback
Keyboard shortcuts
?
Snippets
Groups
Projects
Show more breadcrumbs
Core Modules
dune-istl
Commits
e4b1b35e
Commit
e4b1b35e
authored
15 years ago
by
Markus Blatt
Browse files
Options
Downloads
Patches
Plain Diff
BugFix: Random access containers should have random access iterators
as expected! [[Imported from SVN: r1033]]
parent
8f4d867c
Branches
Branches containing commit
Tags
Tags containing commit
No related merge requests found
Changes
1
Hide whitespace changes
Inline
Side-by-side
Showing
1 changed file
istl/basearray.hh
+13
-3
13 additions, 3 deletions
istl/basearray.hh
with
13 additions
and
3 deletions
istl/basearray.hh
+
13
−
3
View file @
e4b1b35e
...
...
@@ -76,14 +76,14 @@ namespace Dune {
/** \brief Iterator implementation class */
template
<
class
T
>
class
RealIterator
:
public
Bidirectional
IteratorFacade
<
RealIterator
<
T
>
,
T
>
:
public
RandomAccess
IteratorFacade
<
RealIterator
<
T
>
,
T
>
{
public:
//! \brief The unqualified value type
typedef
typename
remove_const
<
T
>::
type
ValueType
;
friend
class
Bidirectional
IteratorFacade
<
RealIterator
<
const
ValueType
>
,
const
ValueType
>
;
friend
class
Bidirectional
IteratorFacade
<
RealIterator
<
ValueType
>
,
ValueType
>
;
friend
class
RandomAccess
IteratorFacade
<
RealIterator
<
const
ValueType
>
,
const
ValueType
>
;
friend
class
RandomAccess
IteratorFacade
<
RealIterator
<
ValueType
>
,
ValueType
>
;
friend
class
RealIterator
<
const
ValueType
>
;
friend
class
RealIterator
<
ValueType
>
;
...
...
@@ -119,6 +119,11 @@ namespace Dune {
return
i
==
other
.
i
;
}
std
::
ptrdiff_t
distanceTo
(
const
RealIterator
&
o
)
const
{
return
o
.
i
-
i
;
}
private
:
//! prefix increment
void
increment
()
...
...
@@ -138,6 +143,11 @@ namespace Dune {
return
*
i
;
}
void
advance
(
std
::
ptrdiff_t
d
)
{
i
+=
d
;
}
const
B
*
p
;
B
*
i
;
};
...
...
This diff is collapsed.
Click to expand it.
Preview
0%
Loading
Try again
or
attach a new file
.
Cancel
You are about to add
0
people
to the discussion. Proceed with caution.
Finish editing this message first!
Save comment
Cancel
Please
register
or
sign in
to comment