Skip to content
GitLab
Explore
Sign in
Primary navigation
Search or go to…
Project
D
dune-gfe
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
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
Sander, Oliver
dune-gfe
Commits
a06fe3af
Commit
a06fe3af
authored
14 years ago
by
Oliver Sander
Committed by
sander@PCPOOL.MI.FU-BERLIN.DE
14 years ago
Browse files
Options
Downloads
Patches
Plain Diff
add method to get an orthonormal frame at each point (for dim==2 only)
[[Imported from SVN: r5964]]
parent
0b5a0097
No related branches found
Branches containing commit
No related tags found
No related merge requests found
Changes
1
Hide whitespace changes
Inline
Side-by-side
Showing
1 changed file
dune/gfe/unitvector.hh
+32
-0
32 additions, 0 deletions
dune/gfe/unitvector.hh
with
32 additions
and
0 deletions
dune/gfe/unitvector.hh
+
32
−
0
View file @
a06fe3af
...
...
@@ -39,7 +39,11 @@ public:
/** \brief The type used for coordinates */
typedef
double
ctype
;
/** \brief Global coordinates wrt an isometric embedding function are available */
static
const
bool
isometricallyEmbedded
=
true
;
typedef
Dune
::
FieldVector
<
double
,
dim
-
1
>
TangentVector
;
typedef
Dune
::
FieldVector
<
double
,
dim
>
EmbeddedTangentVector
;
UnitVector
<
dim
>&
operator
=
(
const
Dune
::
FieldVector
<
double
,
dim
>&
vector
)
...
...
@@ -49,6 +53,17 @@ public:
return
*
this
;
}
/** \brief The exponential map */
static
UnitVector
exp
(
const
UnitVector
&
p
,
const
TangentVector
&
v
)
{
Dune
::
FieldMatrix
<
double
,
dim
-
1
,
dim
>
frame
=
p
.
orthonormalFrame
();
EmbeddedTangentVector
ev
;
frame
.
mtv
(
v
,
ev
);
return
exp
(
p
,
ev
);
}
/** \brief The exponential map */
static
UnitVector
exp
(
const
UnitVector
&
p
,
const
EmbeddedTangentVector
&
v
)
{
...
...
@@ -146,6 +161,23 @@ public:
return
data_
;
}
/** \brief Compute an orthonormal basis of the tangent space of S^n.
This basis is of course not globally continuous.
*/
Dune
::
FieldMatrix
<
double
,
dim
-
1
,
dim
>
orthonormalFrame
()
const
{
Dune
::
FieldMatrix
<
double
,
dim
-
1
,
dim
>
result
;
if
(
dim
==
2
)
{
result
[
0
][
0
]
=
-
data_
[
1
];
result
[
0
][
1
]
=
data_
[
0
];
}
else
DUNE_THROW
(
Dune
::
NotImplemented
,
"orthonormalFrame for dim!=2!"
);
return
result
;
}
/** \brief Write LocalKey object to output stream */
friend
std
::
ostream
&
operator
<<
(
std
::
ostream
&
s
,
const
UnitVector
&
unitVector
)
{
...
...
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