Skip to content
GitLab
Explore
Sign in
Primary navigation
Search or go to…
Project
D
dune-subgrid
Manage
Activity
Members
Labels
Plan
Issues
Issue boards
Milestones
Wiki
Requirements
Code
Merge requests
Repository
Branches
Commits
Tags
Repository graph
Compare revisions
Snippets
Locked files
Build
Pipelines
Jobs
Pipeline schedules
Test cases
Artifacts
Deploy
Releases
Container registry
Model registry
Operate
Environments
Monitor
Incidents
Analyze
Value stream analytics
Contributor analytics
CI/CD analytics
Repository analytics
Code review analytics
Issue analytics
Insights
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
neubertw
dune-subgrid
Commits
8ce8b3c0
Commit
8ce8b3c0
authored
7 years ago
by
Jonathan Youett
Browse files
Options
Downloads
Patches
Plain Diff
Cleanup
parent
17321126
No related branches found
No related tags found
No related merge requests found
Changes
2
Show whitespace changes
Inline
Side-by-side
Showing
2 changed files
dune/subgrid/subgrid/subgridgeometry.hh
+10
-27
10 additions, 27 deletions
dune/subgrid/subgrid/subgridgeometry.hh
dune/subgrid/subgrid/subgridintersection.hh
+50
-49
50 additions, 49 deletions
dune/subgrid/subgrid/subgridintersection.hh
with
60 additions
and
76 deletions
dune/subgrid/subgrid/subgridgeometry.hh
+
10
−
27
View file @
8ce8b3c0
...
...
@@ -28,12 +28,10 @@ template<int mydim, int coorddim, class GridImp> class SubGridLocalGeometry;
because the Dune grid interface specifies the number and types of template
parameters of that class.
*/
template
<
int
mydim
,
int
coorddim
,
class
GridImp
>
//, bool isLocal>
template
<
int
mydim
,
int
coorddim
,
class
GridImp
>
class
SubGridGeometry
{
typedef
typename
GridImp
::
ctype
ctype
;
using
ctype
=
typename
GridImp
::
ctype
;
public:
...
...
@@ -41,12 +39,6 @@ class SubGridGeometry
enum
{
CodimInHostGrid
=
GridImp
::
HostGridType
::
dimension
-
mydim
};
enum
{
DimensionWorld
=
GridImp
::
HostGridType
::
dimensionworld
};
// Select the type we are wrapping depending on the isLocal flag
// - isLocal == true: LocalGeometry
// - isLocal == false: Geometry
// typedef typename std::conditional<isLocal,
//typename GridImp::HostGridType::Traits::template Codim<CodimInHostGrid>::LocalGeometry,
// typename GridImp::HostGridType::Traits::template Codim<CodimInHostGrid>::Geometry>::type HostGridGeometry;
using
HostGridGeometry
=
typename
GridImp
::
HostGridType
::
Traits
::
template
Codim
<
CodimInHostGrid
>
::
Geometry
;
//! The type of the Jacobian matrix of the mapping from the reference element to this element
...
...
@@ -126,24 +118,21 @@ class SubGridGeometry
return
hostGeometry_
.
center
();
}
//! The Jacobian matrix of the mapping from the reference element to this element
const
JacobianTransposed
jacobianTransposed
(
const
FieldVector
<
ctype
,
mydim
>&
local
)
const
{
return
hostGeometry_
.
jacobianTransposed
(
local
);
}
//! The inverse of the Jacobian matrix of the mapping from the reference element to this element
const
JacobianInverseTransposed
jacobianInverseTransposed
(
const
FieldVector
<
ctype
,
mydim
>&
local
)
const
{
return
hostGeometry_
.
jacobianInverseTransposed
(
local
);
}
private
:
const
HostGridGeometry
hostGeometry_
;
};
/** \brief
Encapsulate a
host
grid
LocalG
eometry */
/** \brief
Local SubGrid geometries. They can be either a wrapped
hostgrid
geometry or a multilinearg
eometry
.
*/
template
<
int
mydim
,
int
coorddim
,
class
GridImp
>
class
SubGridLocalGeometry
{
...
...
@@ -160,7 +149,8 @@ class SubGridLocalGeometry
// TODO does this work for arbitratry co-dimensions?
using
HostGridLocalGeometry
=
typename
GridImp
::
HostGridType
::
Traits
::
template
Codim
<
CodimInHostGrid
>
::
LocalGeometry
;
using
MultiLinGeometry
=
MultiLinearGeometry
<
ctype
,
mydim
,
coorddim
>
;
// variant can be either of the template types
// use a variant to store either a HostGridLocalGeometry or a MultiLinearGeometry
using
GeometryContainer
=
Std
::
variant
<
HostGridLocalGeometry
,
MultiLinGeometry
>
;
public:
...
...
@@ -178,7 +168,7 @@ class SubGridLocalGeometry
using
JacobianInverseTransposed
=
typename
std
::
common_type
<
Dune
::
FieldMatrix
<
ctype
,
coorddim
,
mydim
>
,
typename
HostGridLocalGeometry
::
JacobianInverseTransposed
>::
type
;
/** Construct
or for
a given host grid geometry*/
/**
\brief
Construct
from
a given host grid geometry*/
SubGridLocalGeometry
(
HostGridLocalGeometry
hostLocalGeometry
)
:
localGeometry_
(
std
::
move
(
hostLocalGeometry
))
{}
...
...
@@ -194,7 +184,6 @@ class SubGridLocalGeometry
return
Std
::
visit
([
&
](
auto
&&
geom
)
{
return
geom
.
type
();},
localGeometry_
);
}
/** \brief Return true if the geometry mapping is affine and false otherwise */
bool
affine
()
const
{
return
Std
::
visit
([
&
](
auto
&&
geom
)
{
return
geom
.
affine
();},
localGeometry_
);
...
...
@@ -205,27 +194,23 @@ class SubGridLocalGeometry
return
Std
::
visit
([
&
](
auto
&&
geom
)
{
return
geom
.
corners
();},
localGeometry_
);
}
//! access to coordinates of corners. Index is the number of the corner
FieldVector
<
ctype
,
coorddim
>
corner
(
int
i
)
const
{
return
Std
::
visit
([
&
](
auto
&&
geom
)
{
return
geom
.
corner
(
i
);},
localGeometry_
);
}
/** \brief Maps a local coordinate within reference element to
* global coordinate in element */
FieldVector
<
ctype
,
coorddim
>
global
(
const
FieldVector
<
ctype
,
mydim
>&
local
)
const
{
return
Std
::
visit
([
&
](
auto
&&
geom
)
{
return
geom
.
global
(
local
);},
localGeometry_
);
}
/** \brief Maps a global coordinate within the element to a
* local coordinate in its reference element */
FieldVector
<
ctype
,
mydim
>
local
(
const
FieldVector
<
ctype
,
coorddim
>&
global
)
const
{
return
Std
::
visit
([
&
](
auto
&&
geom
)
{
return
geom
.
local
(
global
);},
localGeometry_
);
}
/**
*/
ctype
integrationElement
(
const
FieldVector
<
ctype
,
mydim
>&
local
)
const
{
...
...
@@ -250,19 +235,17 @@ class SubGridLocalGeometry
return
Std
::
visit
([
&
](
auto
&&
geom
)
{
return
geom
.
center
();},
localGeometry_
);
}
//! The Jacobian matrix of the mapping from the reference element to this element
const
JacobianTransposed
jacobianTransposed
(
const
FieldVector
<
ctype
,
mydim
>&
local
)
const
{
return
Std
::
visit
([
&
](
auto
&&
geom
)
{
return
geom
.
jacobianTransposed
(
local
);},
localGeometry_
);
}
//! The inverse of the Jacobian matrix of the mapping from the reference element to this element
const
JacobianInverseTransposed
jacobianInverseTransposed
(
const
FieldVector
<
ctype
,
mydim
>&
local
)
const
{
return
Std
::
visit
([
&
](
auto
&&
geom
)
{
return
geom
.
jacobianInverseTransposed
(
local
);},
localGeometry_
);
}
private
:
const
GeometryContainer
localGeometry_
;
};
...
...
This diff is collapsed.
Click to expand it.
dune/subgrid/subgrid/subgridintersection.hh
+
50
−
49
View file @
8ce8b3c0
...
...
@@ -92,7 +92,6 @@ class SubGridLeafIntersection :
using
ctype
=
typename
GridImp
::
ctype
;
using
Base
::
insideIntersect_
;
using
ctype
=
typename
GridImp
::
ctype
;
public
:
typedef
typename
GridImp
::
template
Codim
<
0
>
::
Entity
Entity
;
...
...
@@ -233,38 +232,39 @@ public:
return
boundary
()
or
((
outside_
.
level
()
==
inside_
.
level
())
and
insideIntersect_
.
conforming
());
}
/
/! intersection of codimension 1 of this neighbor with element where
//! iteration started
.
//! Here returned element is in LOCAL coordinates of the element
//! where iteration started.
/
** \brief Return the local geometry, mapping local coordinates of the intersection to
* local coordinates of the inside element
.
* If the intersection is not part of the host grid, construct one using MultiLinearGeometry.
*/
LocalGeometry
geometryInInside
()
const
{
if
(
boundary
())
return
LocalGeometry
(
LocalGeometryImpl
(
insideIntersect_
.
geometryInInside
()));
else
{
// Case 0: inside_'s the smaller element, or both elements are on the same level.
if
(
inside_
.
level
()
>=
outside_
.
level
())
{
if
(
inside_
.
level
()
>=
outside_
.
level
())
return
LocalGeometry
(
LocalGeometryImpl
(
insideIntersect_
.
geometryInInside
()));
}
// Case 1: outsideIntersect_ is the intersection we want, embedded in inside_'s geometry
auto
outsideGeo
=
outsideIntersect_
.
geometry
();
auto
outsideIsGeo
=
outsideIntersect_
.
geometry
();
auto
insideGeo
=
inside_
.
geometry
();
using
Coords
=
std
::
decay_t
<
decltype
(
outsideGeo
.
corner
(
0
))
>
;
auto
corners
=
std
::
vector
<
Coords
>
(
outsideGeo
.
corners
());
using
Coords
=
std
::
decay_t
<
decltype
(
insideGeo
.
local
(
outside
Is
Geo
.
corner
(
0
))
)
>
;
auto
corners
=
std
::
vector
<
Coords
>
(
outside
Is
Geo
.
corners
());
// Map the coordinates of the corners from the outside's local coordinates into inside's
// local coordinates
for
(
size_t
i
=
0
;
i
<
corners
.
size
();
i
++
)
corners
[
i
]
=
inside
_
.
geometry
()
.
local
(
outsideGeo
.
corner
(
i
));
corners
[
i
]
=
inside
Geo
.
local
(
outside
Is
Geo
.
corner
(
i
));
// create a new geometry with the coordinates we just computed
auto
geo
=
Dune
::
MultiLinearGeometry
<
ctype
,
dim
-
1
,
dim
>
(
outsideGeo
.
type
(),
corners
);
auto
geo
=
Dune
::
MultiLinearGeometry
<
ctype
,
dim
-
1
,
dim
>
(
outside
Is
Geo
.
type
(),
corners
);
return
LocalGeometry
(
LocalGeometryImpl
(
std
::
move
(
geo
)));
}
}
//! intersection of codimension 1 of this neighbor with element where iteration started.
//! Here returned element is in LOCAL coordinates of neighbor
/** \brief Return the local geometry, mapping local coordinates of the intersection to
* local coordinates of the outside element.
*/
LocalGeometry
geometryInOutside
()
const
{
if
(
boundary
())
DUNE_THROW
(
Dune
::
Exception
,
"No outside geometry available for boundary intersections"
);
...
...
@@ -272,32 +272,33 @@ public:
// Case 0: outside_'s the smaller element, or both elements are on the same level.
if
(
inside_
.
level
()
<
outside_
.
level
())
{
auto
outsideGeo
=
outsideIntersect_
.
geometry
();
auto
outsideIsGeo
=
outsideIntersect_
.
geometry
();
auto
outsideGeo
=
outside_
.
geometry
();
using
Coords
=
std
::
decay_t
<
decltype
(
outside
_
.
geometry
()
.
local
(
outsideGeo
.
corner
(
0
)))
>
;
auto
corners
=
std
::
vector
<
Coords
>
(
outsideGeo
.
corners
());
using
Coords
=
std
::
decay_t
<
decltype
(
outside
Geo
.
local
(
outside
Is
Geo
.
corner
(
0
)))
>
;
auto
corners
=
std
::
vector
<
Coords
>
(
outside
Is
Geo
.
corners
());
for
(
size_t
i
=
0
;
i
<
corners
.
size
();
i
++
)
corners
[
i
]
=
outside
_
.
geometry
()
.
local
(
outsideGeo
.
corner
(
i
));
corners
[
i
]
=
outside
Geo
.
local
(
outside
Is
Geo
.
corner
(
i
));
// create a new geometry with the coordinates we just computed
auto
geo
=
Dune
::
MultiLinearGeometry
<
ctype
,
dim
-
1
,
dim
>
(
outsideGeo
.
type
(),
corners
);
auto
geo
=
Dune
::
MultiLinearGeometry
<
ctype
,
dim
-
1
,
dim
>
(
outside
Is
Geo
.
type
(),
corners
);
return
LocalGeometry
(
LocalGeometryImpl
(
geo
));
}
// Case 1: outsideIntersect_ is the intersection we want, embedded in inside_'s geometry
auto
insideIsGeo
=
insideIntersect_
.
geometry
();
auto
outsideGeo
=
outside_
.
geometry
();
auto
insideGeo
=
insideIntersect_
.
geometry
();
using
Coords
=
std
::
decay_t
<
decltype
(
outside_
.
geometry
().
local
(
insideGeo
.
corner
(
0
)))
>
;
auto
corners
=
std
::
vector
<
Coords
>
(
insideGeo
.
corners
());
using
Coords
=
std
::
decay_t
<
decltype
(
outsideGeo
.
local
(
insideIsGeo
.
corner
(
0
)))
>
;
auto
corners
=
std
::
vector
<
Coords
>
(
insideIsGeo
.
corners
());
// Map the coordinates of the corners from the insides's local coordinates into outside's
// local coordinates
for
(
size_t
i
=
0
;
i
<
corners
.
size
();
i
++
)
corners
[
i
]
=
outside
_
.
geometry
()
.
local
(
insideGeo
.
corner
(
i
));
corners
[
i
]
=
outside
Geo
.
local
(
inside
Is
Geo
.
corner
(
i
));
// create a new geometry with the coordinates we just computed
auto
geo
=
Dune
::
MultiLinearGeometry
<
ctype
,
dim
-
1
,
dim
>
(
insideGeo
.
type
(),
corners
);
auto
geo
=
Dune
::
MultiLinearGeometry
<
ctype
,
dim
-
1
,
dim
>
(
inside
Is
Geo
.
type
(),
corners
);
return
LocalGeometry
(
LocalGeometryImpl
(
std
::
move
(
geo
)));
}
...
...
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