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
264d819b
Commit
264d819b
authored
18 years ago
by
Oliver Sander
Committed by
sander@PCPOOL.MI.FU-BERLIN.DE
18 years ago
Browse files
Options
Downloads
Patches
Plain Diff
moved interpolation of derivatives to the quaternion class
[[Imported from SVN: r1116]]
parent
2663701b
No related branches found
Branches containing commit
No related tags found
No related merge requests found
Changes
2
Hide whitespace changes
Inline
Side-by-side
Showing
2 changed files
src/quaternion.hh
+44
-0
44 additions, 0 deletions
src/quaternion.hh
src/rodassembler.cc
+6
-12
6 additions, 12 deletions
src/rodassembler.cc
with
50 additions
and
12 deletions
src/quaternion.hh
+
44
−
0
View file @
264d819b
...
...
@@ -121,9 +121,20 @@ public:
return
result
;
}
/** \brief Invert the quaternion */
void
invert
()
{
(
*
this
)[
0
]
*=
-
1
;
(
*
this
)[
1
]
*=
-
1
;
(
*
this
)[
2
]
*=
-
1
;
(
*
this
)
/=
this
->
two_norm2
();
}
/** \brief Interpolate between two rotations */
static
Quaternion
<
T
>
interpolate
(
const
Quaternion
<
T
>&
a
,
const
Quaternion
<
T
>&
b
,
double
omega
)
{
#if 0 // Interpolation in H plus normalization
Quaternion<T> result;
for (int i=0; i<4; i++)
...
...
@@ -131,6 +142,39 @@ public:
result.normalize();
return result;
#else
// Interpolation on the geodesic in SO(3) from a to b
Quaternion
<
T
>
diff
=
a
;
diff
.
invert
();
diff
.
mult
(
b
);
T
dist
=
2
*
std
::
acos
(
diff
[
3
]);
T
invSinc
=
(
dist
<
1e-4
)
?
1
/
(
0.5
+
(
dist
*
dist
/
48
))
:
dist
/
std
::
sin
(
dist
/
2
);
// Compute difference on T_a SO(3)
Dune
::
FieldVector
<
double
,
3
>
v
;
v
[
0
]
=
diff
[
0
]
*
invSinc
;
v
[
1
]
=
diff
[
1
]
*
invSinc
;
v
[
2
]
=
diff
[
2
]
*
invSinc
;
v
*=
omega
;
return
a
.
mult
(
exp
(
v
[
0
],
v
[
1
],
v
[
2
]));
#endif
}
/** \brief Interpolate between two rotations */
static
Quaternion
<
T
>
interpolateDerivative
(
const
Quaternion
<
T
>&
a
,
const
Quaternion
<
T
>&
b
,
double
omega
,
double
intervallLength
)
{
Quaternion
<
T
>
result
;
for
(
int
i
=
0
;
i
<
4
;
i
++
)
result
[
i
]
=
a
[
i
]
/
(
-
intervallLength
)
+
b
[
i
]
/
intervallLength
;
return
result
;
}
...
...
This diff is collapsed.
Click to expand it.
src/rodassembler.cc
+
6
−
12
View file @
264d819b
...
...
@@ -303,9 +303,8 @@ getLocalMatrix( EntityPointer &entity,
}
// Get the derivative of the rotation at the quadrature point by interpolating in $H$
Quaternion
<
double
>
hatq_s
;
for
(
int
i
=
0
;
i
<
4
;
i
++
)
hatq_s
[
i
]
=
localSolution
[
0
].
q
[
i
]
*
shapeGrad
[
0
]
+
localSolution
[
1
].
q
[
i
]
*
shapeGrad
[
1
];
Quaternion
<
double
>
hatq_s
=
Quaternion
<
double
>::
interpolateDerivative
(
localSolution
[
0
].
q
,
localSolution
[
1
].
q
,
quadPos
,
1
/
shapeGrad
[
1
]);
// The current strain
FieldVector
<
double
,
blocksize
>
strain
=
getStrain
(
globalSolution
,
entity
,
quadPos
);
...
...
@@ -402,7 +401,6 @@ getLocalMatrix( EntityPointer &entity,
// \partial W^2 \partial v^i_j \partial v^k_l
// All other derivatives are zero
//double sum = duLocal_dvij[k][l][m] * (duCan_dvij[i][j] * hatq.director(m) + darbouxCan*dd_dvj[m][j]*shapeFunction[i]);
double
sum
=
du_dvij
[
k
][
l
][
m
]
*
du_dvij
[
i
][
j
][
m
];
sum
+=
(
strain
[
m
+
3
]
-
referenceStrain
[
m
+
3
])
*
du_dvij_dvkl
[
i
][
j
][
k
][
l
][
m
];
...
...
@@ -501,9 +499,8 @@ assembleGradient(const std::vector<Configuration>& sol,
Quaternion
<
double
>
hatq
=
Quaternion
<
double
>::
interpolate
(
localSolution
[
0
].
q
,
localSolution
[
1
].
q
,
quadPos
[
0
]);
// Get the derivative of the rotation at the quadrature point by interpolating in $H$
Quaternion
<
double
>
hatq_s
;
for
(
int
i
=
0
;
i
<
4
;
i
++
)
hatq_s
[
i
]
=
localSolution
[
0
].
q
[
i
]
*
shapeGrad
[
0
]
+
localSolution
[
1
].
q
[
i
]
*
shapeGrad
[
1
];
Quaternion
<
double
>
hatq_s
=
Quaternion
<
double
>::
interpolateDerivative
(
localSolution
[
0
].
q
,
localSolution
[
1
].
q
,
quadPos
,
1
/
shapeGrad
[
1
]);
// The current strain
FieldVector
<
double
,
blocksize
>
strain
=
getStrain
(
sol
,
it
,
quadPos
);
...
...
@@ -811,11 +808,8 @@ Dune::FieldVector<double, 6> Dune::RodAssembler<GridType>::getStrain(const std::
Quaternion
<
double
>
q
=
Quaternion
<
double
>::
interpolate
(
localSolution
[
0
].
q
,
localSolution
[
1
].
q
,
pos
);
// Get the derivative of the rotation at the quadrature point by interpolating in $H$
Quaternion
<
double
>
q_s
;
q_s
[
0
]
=
localSolution
[
0
].
q
[
0
]
*
shapeGrad
[
0
][
0
]
+
localSolution
[
1
].
q
[
0
]
*
shapeGrad
[
1
][
0
];
q_s
[
1
]
=
localSolution
[
0
].
q
[
1
]
*
shapeGrad
[
0
][
0
]
+
localSolution
[
1
].
q
[
1
]
*
shapeGrad
[
1
][
0
];
q_s
[
2
]
=
localSolution
[
0
].
q
[
2
]
*
shapeGrad
[
0
][
0
]
+
localSolution
[
1
].
q
[
2
]
*
shapeGrad
[
1
][
0
];
q_s
[
3
]
=
localSolution
[
0
].
q
[
3
]
*
shapeGrad
[
0
][
0
]
+
localSolution
[
1
].
q
[
3
]
*
shapeGrad
[
1
][
0
];
Quaternion
<
double
>
q_s
=
Quaternion
<
double
>::
interpolateDerivative
(
localSolution
[
0
].
q
,
localSolution
[
1
].
q
,
pos
,
1
/
shapeGrad
[
1
]);
// /////////////////////////////////////////////
// Sum it all up
...
...
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