Skip to content
Projects
Groups
Snippets
Help
Loading...
Help
Support
Keyboard shortcuts
?
Submit feedback
Contribute to GitLab
Sign in
Toggle navigation
I
INMOST
Project overview
Project overview
Details
Activity
Releases
Repository
Repository
Files
Commits
Branches
Tags
Contributors
Graph
Compare
Issues
0
Issues
0
List
Boards
Labels
Milestones
Merge Requests
0
Merge Requests
0
CI / CD
CI / CD
Pipelines
Jobs
Schedules
Analytics
Analytics
CI / CD
Repository
Value Stream
Wiki
Wiki
Snippets
Snippets
Members
Members
Collapse sidebar
Close sidebar
Activity
Graph
Create a new issue
Jobs
Commits
Issue Boards
Open sidebar
Kirill Terekhov
INMOST
Commits
3cbe1ef4
Commit
3cbe1ef4
authored
Jan 11, 2015
by
Kirill Terekhov
Browse files
Options
Browse Files
Download
Email Patches
Plain Diff
Fix
Fix introduced bug - parameters are set after inner solvers were initialized.
parent
8fa00d6a
Changes
2
Hide whitespace changes
Inline
Side-by-side
Showing
2 changed files
with
27 additions
and
25 deletions
+27
-25
inmost_solver.h
inmost_solver.h
+6
-6
solver.cpp
solver.cpp
+21
-19
No files found.
inmost_solver.h
View file @
3cbe1ef4
...
...
@@ -441,7 +441,7 @@ namespace INMOST
INMOST_DATA_REAL_TYPE
preconditioner_reuse_tolerance
;
INMOST_DATA_REAL_TYPE
preconditioner_ddpq_tolerance
;
INMOST_DATA_ENUM_TYPE
preconditioner_reorder_nonzero
;
INMOST_DATA_
ENUM
_TYPE
preconditioner_fill_level
;
INMOST_DATA_
REAL
_TYPE
preconditioner_fill_level
;
INMOST_DATA_ENUM_TYPE
preconditioner_rescale_iterations
;
INMOST_DATA_ENUM_TYPE
preconditioner_condition_estimation
;
INMOST_DATA_ENUM_TYPE
preconditioner_adapt_ddpq_tolerance
;
...
...
@@ -473,10 +473,6 @@ namespace INMOST
/// "gmres_substeps" - number of gmres steps performed after each bicgstab step
/// works for:
/// INNER_ILU2, INNER_MLILUC
/// "fill_level" - level of fill for ILU-type preconditioners
/// works for:
/// INNER_ILU2 (if LFILL is defined in solver_ilu2.hpp)
/// Trilinos_Ifpack
/// "reorder_nonzeros" - place sparser rows at the beggining of matrix during reordering
/// works for:
/// INNER_MLILUC
...
...
@@ -504,7 +500,7 @@ namespace INMOST
/// "drop_tolerance" - tolerance for dropping values during incomplete factorization
/// works for:
/// INNER_ILU2, INNER_MLILUC
/// Trilinos_Aztec, Trilinos_
ML, Trilinos_
Ifpack
/// Trilinos_Aztec, Trilinos_Ifpack
/// PETSc
/// "reuse_tolerance" - tolerance for reusing values during incomplete factorization
/// these values are used only during calculation of L and U factors
...
...
@@ -521,6 +517,10 @@ namespace INMOST
/// where on imax, jmax maximum is reached.
/// works for:
/// INNER_MLILUC
/// "fill_level" - level of fill for ILU-type preconditioners
/// works for:
/// INNER_ILU2 (if LFILL is defined in solver_ilu2.hpp)
/// Trilinos, Trilinos_Ifpack
void
SetParameterReal
(
std
::
string
name
,
INMOST_DATA_REAL_TYPE
value
);
/// Get the used defined name of the Solver.
std
::
string
GetName
()
{
return
name
;}
...
...
solver.cpp
View file @
3cbe1ef4
...
...
@@ -1202,12 +1202,10 @@ namespace INMOST
preconditioner_condition_estimation
=
val
;
else
if
(
name
==
"adapt_ddpq_tolerance"
)
preconditioner_adapt_ddpq_tolerance
=
val
;
else
if
(
name
==
"scwartz_overlap"
)
else
if
(
name
==
"sc
h
wartz_overlap"
)
additive_schwartz_overlap
=
val
;
else
if
(
name
==
"gmres_substeps"
)
solver_gmres_substeps
=
val
;
else
if
(
name
==
"fill_level"
)
preconditioner_fill_level
=
val
;
else
if
(
name
==
"reorder_nonzeros"
)
preconditioner_reorder_nonzero
=
val
;
else
if
(
_pack
==
INNER_ILU2
||
_pack
==
INNER_MLILUC
)
...
...
@@ -1241,6 +1239,8 @@ namespace INMOST
preconditioner_reuse_tolerance
=
val
;
else
if
(
name
==
"ddpq_tolerance"
)
preconditioner_ddpq_tolerance
=
val
;
else
if
(
name
==
"fill_level"
)
preconditioner_fill_level
=
val
;
else
if
(
_pack
==
INNER_ILU2
||
_pack
==
INNER_MLILUC
)
{
try
...
...
@@ -1765,10 +1765,10 @@ namespace INMOST
sol
->
ReplaceMAT
(
*
mat
);
if
(
matrix_data
!=
NULL
)
delete
(
Solver
::
Matrix
*
)
matrix_data
;
matrix_data
=
(
void
*
)
mat
;
if
(
!
sol
->
isInitialized
())
{
sol
->
Initialize
();
}
//
if (!sol->isInitialized())
//
{
//
sol->Initialize();
//
}
ok
=
true
;
}
...
...
@@ -1884,16 +1884,16 @@ namespace INMOST
problem
->
SetRHS
(
&
VectorRHS
);
problem
->
SetLHS
(
&
VectorSOL
);
AztecOO
AztecSolver
(
*
problem
);
//
AztecSolver.SetAztecOption(AZ_diagnostics,AZ_none);
//
AztecSolver.SetAztecOption(AZ_output,AZ_none);
AztecSolver
.
SetAztecOption
(
AZ_diagnostics
,
AZ_none
);
AztecSolver
.
SetAztecOption
(
AZ_output
,
AZ_none
);
AztecSolver
.
SetAztecOption
(
AZ_solver
,
AZ_bicgstab
);
AztecSolver
.
SetAztecOption
(
AZ_overlap
,
additive_schwartz_overlap
);
void
*
precinfo
=
NULL
;
if
(
_pack
==
Trilinos_Aztec
)
{
AztecSolver
.
SetAztecParam
(
AZ_
tol
,
preconditioner_drop_tolerance
);
AztecSolver
.
SetAztecParam
(
AZ_
drop
,
preconditioner_fill_level
);
AztecSolver
.
SetAztecParam
(
AZ_
drop
,
preconditioner_drop_tolerance
);
AztecSolver
.
SetAztecParam
(
AZ_
ilut_fill
,
preconditioner_fill_level
);
//AztecSolver.SetAztecOption(AZ_solver,AZ_tfqmr);
//AztecSolver.SetAztecOption(AZ_precond,AZ_Neumann);
//AztecSolver.SetAztecOption(AZ_poly_ord,3);
...
...
@@ -2052,10 +2052,6 @@ namespace INMOST
if
(
_pack
==
INNER_ILU2
)
{
IterativeMethod
*
sol
=
static_cast
<
IterativeMethod
*>
(
solver_data
);
if
(
!
sol
->
isInitialized
())
{
sol
->
Initialize
();
}
sol
->
EnumParameter
(
"maxits"
)
=
maximum_iterations
;
sol
->
EnumParameter
(
"levels"
)
=
solver_gmres_substeps
;
sol
->
RealParameter
(
"rtol"
)
=
relative_tolerance
;
...
...
@@ -2065,6 +2061,11 @@ namespace INMOST
sol
->
RealParameter
(
":tau2"
)
=
preconditioner_reuse_tolerance
;
sol
->
EnumParameter
(
":fill"
)
=
preconditioner_fill_level
;
sol
->
EnumParameter
(
":scale_iters"
)
=
preconditioner_rescale_iterations
;
if
(
!
sol
->
isInitialized
())
{
sol
->
Initialize
();
}
bool
ret
=
sol
->
Solve
(
RHS
,
SOL
);
last_it
=
sol
->
GetIterations
();
last_resid
=
sol
->
GetResidual
();
...
...
@@ -2074,10 +2075,6 @@ namespace INMOST
if
(
_pack
==
INNER_MLILUC
)
{
IterativeMethod
*
sol
=
static_cast
<
IterativeMethod
*>
(
solver_data
);
if
(
!
sol
->
isInitialized
())
{
sol
->
Initialize
();
}
sol
->
EnumParameter
(
"maxits"
)
=
maximum_iterations
;
sol
->
EnumParameter
(
"levels"
)
=
solver_gmres_substeps
;
sol
->
RealParameter
(
"rtol"
)
=
relative_tolerance
;
...
...
@@ -2090,6 +2087,11 @@ namespace INMOST
sol
->
EnumParameter
(
":scale_iters"
)
=
preconditioner_rescale_iterations
;
sol
->
EnumParameter
(
":estimator"
)
=
preconditioner_condition_estimation
;
sol
->
EnumParameter
(
":ddpq_tau_adapt"
)
=
preconditioner_adapt_ddpq_tolerance
;
if
(
!
sol
->
isInitialized
())
{
sol
->
Initialize
();
}
bool
ret
=
sol
->
Solve
(
RHS
,
SOL
);
last_it
=
sol
->
GetIterations
();
last_resid
=
sol
->
GetResidual
();
...
...
Write
Preview
Markdown
is supported
0%
Try again
or
attach a new file
Attach a file
Cancel
You are about to add
0
people
to the discussion. Proceed with caution.
Finish editing this message first!
Cancel
Please
register
or
sign in
to comment