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
411d0fae
Commit
411d0fae
authored
11 years ago
by
Oliver Sander
Committed by
sander
11 years ago
Browse files
Options
Downloads
Patches
Plain Diff
Allow to specify a minimal number of trust-region iterations
[[Imported from SVN: r9450]]
parent
46fe369d
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
dune/gfe/targetspacertrsolver.cc
+5
-1
5 additions, 1 deletion
dune/gfe/targetspacertrsolver.cc
dune/gfe/targetspacertrsolver.hh
+9
-0
9 additions, 0 deletions
dune/gfe/targetspacertrsolver.hh
with
14 additions
and
1 deletion
dune/gfe/targetspacertrsolver.cc
+
5
−
1
View file @
411d0fae
...
...
@@ -24,6 +24,7 @@ setup(const AverageDistanceAssembler<TargetSpace>* assembler,
innerIterations_
=
innerIterations
;
innerTolerance_
=
innerTolerance
;
this
->
verbosity_
=
NumProc
::
QUIET
;
minNumberOfIterations_
=
4
;
// ////////////////////////////////
// Create a projected gauss-seidel solver
...
...
@@ -48,6 +49,8 @@ setup(const AverageDistanceAssembler<TargetSpace>* assembler,
template
<
class
TargetSpace
>
void
TargetSpaceRiemannianTRSolver
<
TargetSpace
>::
solve
()
{
assert
(
minNumberOfIterations_
>
0
);
MaxNormTrustRegion
<
blocksize
,
field_type
>
trustRegion
(
1
,
// we have only one block
initialTrustRegionRadius_
);
...
...
@@ -93,7 +96,7 @@ void TargetSpaceRiemannianTRSolver<TargetSpace>::solve()
if
(
this
->
verbosity_
==
NumProc
::
FULL
)
std
::
cout
<<
"Infinity norm of the correction: "
<<
corr
.
infinity_norm
()
<<
std
::
endl
;
if
(
corr
.
infinity_norm
()
<
this
->
tolerance_
)
{
if
(
corr
.
infinity_norm
()
<
this
->
tolerance_
and
i
>=
minNumberOfIterations_
-
1
)
{
if
(
this
->
verbosity_
==
NumProc
::
FULL
)
std
::
cout
<<
"CORRECTION IS SMALL ENOUGH"
<<
std
::
endl
;
...
...
@@ -136,6 +139,7 @@ void TargetSpaceRiemannianTRSolver<TargetSpace>::solve()
}
if
(
energy
>=
oldEnergy
&&
i
>
minNumberOfIterations_
-
1
&&
(
std
::
abs
(
oldEnergy
-
energy
)
/
energy
<
1e-9
||
modelDecrease
/
energy
<
1e-9
))
{
if
(
this
->
verbosity_
==
NumProc
::
FULL
)
std
::
cout
<<
"Suspecting rounding problems"
<<
std
::
endl
;
...
...
This diff is collapsed.
Click to expand it.
dune/gfe/targetspacertrsolver.hh
+
9
−
0
View file @
411d0fae
...
...
@@ -81,6 +81,15 @@ protected:
/** \brief Norm for the quadratic inner problems */
std
::
auto_ptr
<
EnergyNorm
<
MatrixType
,
CorrectionType
>
>
energyNorm_
;
/** \brief Specify a minimal number of iterations the trust-region solver has to do
*
* This is needed when working with automatic differentiation. While a very low
* number of iterations may be enough to precisely compute the value of a
* geodesic finite element function, a higher number may be needed to make an AD
* system compute a derivative with sufficient precision.
*/
size_t
minNumberOfIterations_
;
};
#include
"targetspacertrsolver.cc"
...
...
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