Skip to content
GitLab
Menu
Projects
Groups
Snippets
Loading...
Help
Help
Support
Community forum
Keyboard shortcuts
?
Submit feedback
Contribute to GitLab
Sign in
Toggle navigation
Menu
Open sidebar
Kirill Terekhov
INMOST
Commits
0af669bb
Commit
0af669bb
authored
Dec 10, 2016
by
Dmitry Bagaev
Browse files
FCBIILU2 interface improvments
1)SetParameter for fcbiilu2 2) deleted void*
parent
7034199e
Changes
4
Expand all
Hide whitespace changes
Inline
Side-by-side
Source/Solver/solver_fcbiilu2/SolverFCBIILU2.cpp
View file @
0af669bb
#include "SolverFCBIILU2.h"
#include "solver_fcbiilu2.h"
namespace
INMOST
{
SolverFCBIILU2
::
SolverFCBIILU2
()
{
solver_data
=
NULL
;
matrix_data
=
NULL
;
}
SolverInterface
*
SolverFCBIILU2
::
Copy
(
const
SolverInterface
*
other
)
{
...
...
@@ -47,7 +49,12 @@ namespace INMOST {
void
SolverFCBIILU2
::
Setup
(
int
*
argc
,
char
***
argv
,
SolverParameters
&
p
)
{
SolverInitDataFcbiilu2
(
&
solver_data
,
communicator
,
p
.
solverPrefix
.
c_str
());
SolverInitializeFcbiilu2
(
argc
,
argv
,
p
.
internalFile
.
c_str
());
solver_data
->
kovl
=
0
;
// number of overlap layers: kovl=0,1,2,...
solver_data
->
tau
=
3e-3
;
// the ILU2 precision (for the submatrix factorization); tau=3e-3
solver_data
->
eps
=
1e-5
;
// the residual precision: ||r|| < eps * ||b||; eps=1e-6
solver_data
->
nit
=
999
;
// number of iterations permitted; nit=999
solver_data
->
msglev
=
2
;
// messages level; msglev=0 for silent; msglev=1 to output solution statistics
SolverInitializeFcbiilu2
(
solver_data
,
argc
,
argv
,
p
.
internalFile
.
c_str
());
}
void
SolverFCBIILU2
::
SetMatrix
(
Sparse
::
Matrix
&
A
,
bool
ModifiedPattern
,
bool
OldPreconditioner
)
{
...
...
@@ -117,11 +124,11 @@ namespace INMOST {
INMOST_DATA_ENUM_TYPE
vbeg
,
vend
;
RHS
.
GetInterval
(
vbeg
,
vend
);
v
oid
*
rhs_data
=
NULL
;
v
ector
*
rhs_data
=
NULL
;
VectorInitDataFcbiilu2
(
&
rhs_data
,
RHS
.
GetCommunicator
(),
RHS
.
GetName
().
c_str
());
VectorPreallocateFcbiilu2
(
rhs_data
,
local_size
);
v
oid
*
solution_data
=
NULL
;
v
ector
*
solution_data
=
NULL
;
VectorInitDataFcbiilu2
(
&
solution_data
,
SOL
.
GetCommunicator
(),
SOL
.
GetName
().
c_str
());
VectorPreallocateFcbiilu2
(
solution_data
,
local_size
);
VectorFillFcbiilu2
(
rhs_data
,
&
RHS
[
vbeg
]);
...
...
@@ -153,16 +160,22 @@ namespace INMOST {
}
void
SolverFCBIILU2
::
SetParameter
(
std
::
string
name
,
std
::
string
value
)
{
std
::
cout
<<
"SolverFCBIILU2::SetParameter unsupported operation"
<<
std
::
endl
;
const
char
*
val
=
value
.
c_str
();
if
(
name
==
"kovl"
)
solver_data
->
kovl
=
atoi
(
val
);
else
if
(
name
==
"tau"
)
solver_data
->
tau
=
atof
(
val
);
else
if
(
name
==
"eps"
)
solver_data
->
eps
=
atof
(
val
);
else
if
(
name
==
"nit"
)
solver_data
->
nit
=
atoi
(
val
);
else
if
(
name
==
"msglev"
)
solver_data
->
msglev
=
atoi
(
val
);
else
std
::
cout
<<
"Parameter "
<<
name
<<
" is unknown"
<<
std
::
endl
;
//throw INMOST::SolverUnsupportedOperation;
}
const
INMOST_DATA_ENUM_TYPE
SolverFCBIILU2
::
Iterations
()
const
{
return
static_cast
<
INMOST_DATA_ENUM_TYPE
>
(
SolverIterationNumberFcbiilu2
(
solver_data
)
);
return
static_cast
<
INMOST_DATA_ENUM_TYPE
>
(
solver_data
->
ITER
);
}
const
INMOST_DATA_REAL_TYPE
SolverFCBIILU2
::
Residual
()
const
{
return
SolverResidualNormFcbiilu2
(
solver_data
)
;
return
solver_data
->
RESID
;
}
const
std
::
string
SolverFCBIILU2
::
ReturnReason
()
const
{
...
...
Source/Solver/solver_fcbiilu2/SolverFCBIILU2.h
View file @
0af669bb
...
...
@@ -8,8 +8,8 @@ namespace INMOST {
class
SolverFCBIILU2
:
public
SolverInterface
{
private:
void
*
solver_data
;
void
*
matrix_data
;
bcg
*
solver_data
;
matrix
*
matrix_data
;
INMOST_DATA_ENUM_TYPE
local_size
,
global_size
;
public:
SolverFCBIILU2
();
...
...
Source/Solver/solver_fcbiilu2/solver_fcbiilu2.cpp
View file @
0af669bb
This diff is collapsed.
Click to expand it.
Source/Solver/solver_fcbiilu2/solver_fcbiilu2.h
View file @
0af669bb
...
...
@@ -2,39 +2,71 @@
#define SOLVER_FCBIILU2_H_INCLUDED
#include "inmost_solver.h"
/* BiCGStab solver structure */
typedef
struct
{
int
n
;
// local number of unknowns at the local processor
int
nproc
;
// number of processors
int
*
ibl
;
// block splitting: ibl[0]=0; ibl[nproc]=nglob
int
*
ia
;
// row pointers: ia[0]=0; ia[nloc]=nzloc
int
*
ja
;
// column numbers (NOTE: starting from 0 or 1); ja[nzloc]
double
*
a
;
// matrix A coefficients; a[nzloc]
int
len_r8
;
// size of the working memory, set len_r8=nbl for the first call
double
*
W
;
// poiter to the working memory W[len_r8]
int
kovl
;
// number of overlap layers: kovl=0,1,2,...
double
tau
;
// the ILU2 precision (for the submatrix factorization); tau=3e-3
double
eps
;
// the residual precision: ||r|| < eps * ||b||; eps=1e-6
int
nit
;
// number of iterations permitted; nit=999
int
msglev
;
// messages level; msglev=0 for silent; msglev=1 to output solution statistics
int
ierr
;
// error flag on return; ierr=0 for success
int
istat
[
16
];
// integer statistics array on return
double
dstat
[
16
];
// double statistics array on return
double
RESID
;
// residual norm
int
ITER
;
// number of BiCGStab iterations performed
}
bcg
;
//#define USE_SOLVER_FCBIILU2 //just for test; see flag HAVE_SOLVER_FCBIILU2
//#if defined(USE_SOLVER_FCBIILU2)
void
MatrixCopyDataFcbiilu2
(
void
**
ppA
,
void
*
pB
);
void
MatrixAssignDataFcbiilu2
(
void
*
pA
,
void
*
pB
);
void
MatrixInitDataFcbiilu2
(
void
**
ppA
,
INMOST_MPI_Comm
comm
,
const
char
*
name
);
void
MatrixDestroyDataFcbiilu2
(
void
**
pA
);
void
MatrixFillFcbiilu2
(
void
*
pA
,
int
size
,
int
nproc
,
int
*
ibl
,
int
*
ia
,
int
*
ja
,
double
*
values
);
void
MatrixFillValuesFcbiilu2
(
void
*
pA
,
double
*
values
);
void
MatrixFinalizeFcbiilu2
(
void
*
data
);
void
VectorInitDataFcbiilu2
(
void
**
ppA
,
INMOST_MPI_Comm
comm
,
const
char
*
name
);
void
VectorCopyDataFcbiilu2
(
void
**
ppA
,
void
*
pB
);
void
VectorAssignDataFcbiilu2
(
void
*
pA
,
void
*
pB
);
void
VectorPreallocateFcbiilu2
(
void
*
pA
,
int
size
);
void
VectorFillFcbiilu2
(
void
*
pA
,
double
*
values
);
void
VectorLoadFcbiilu2
(
void
*
pA
,
double
*
values
);
void
VectorFinalizeFcbiilu2
(
void
*
data
);
void
VectorDestroyDataFcbiilu2
(
void
**
ppA
);
void
SolverInitializeFcbiilu2
(
int
*
argc
,
char
***
argv
,
const
char
*
file_options
);
typedef
struct
{
int
n
;
// local number of unknowns at the local processor
int
nproc
;
// number of processors
int
*
ibl
;
// block splitting: ibl[0]=0; ibl[nproc]=nglob
int
*
ia
;
int
*
ja
;
double
*
A
;
}
matrix
;
typedef
struct
{
int
n
;
// local number of unknowns at the local processor
double
*
v
;
}
vector
;
void
MatrixCopyDataFcbiilu2
(
matrix
**
pA
,
matrix
*
B
);
void
MatrixAssignDataFcbiilu2
(
matrix
*
A
,
matrix
*
B
);
void
MatrixInitDataFcbiilu2
(
matrix
**
ppA
,
INMOST_MPI_Comm
comm
,
const
char
*
name
);
void
MatrixDestroyDataFcbiilu2
(
matrix
**
pA
);
void
MatrixFillFcbiilu2
(
matrix
*
pA
,
int
size
,
int
nproc
,
int
*
ibl
,
int
*
ia
,
int
*
ja
,
double
*
values
);
void
MatrixFillValuesFcbiilu2
(
matrix
*
pA
,
double
*
values
);
void
MatrixFinalizeFcbiilu2
(
matrix
*
data
);
void
VectorInitDataFcbiilu2
(
vector
**
ppA
,
INMOST_MPI_Comm
comm
,
const
char
*
name
);
void
VectorCopyDataFcbiilu2
(
vector
**
ppA
,
vector
*
pB
);
void
VectorAssignDataFcbiilu2
(
vector
*
pA
,
vector
*
pB
);
void
VectorPreallocateFcbiilu2
(
vector
*
pA
,
int
size
);
void
VectorFillFcbiilu2
(
vector
*
pA
,
double
*
values
);
void
VectorLoadFcbiilu2
(
vector
*
pA
,
double
*
values
);
void
VectorFinalizeFcbiilu2
(
vector
*
data
);
void
VectorDestroyDataFcbiilu2
(
vector
**
ppA
);
void
SolverInitializeFcbiilu2
(
bcg
*
data
,
int
*
argc
,
char
***
argv
,
const
char
*
file_options
);
bool
SolverIsFinalizedFcbiilu2
();
void
SolverFinalizeFcbiilu2
();
void
SolverDestroyDataFcbiilu2
(
void
**
data
);
void
SolverInitDataFcbiilu2
(
void
**
data
,
INMOST_MPI_Comm
comm
,
const
char
*
name
);
void
SolverCopyDataFcbiilu2
(
void
**
data
,
void
*
other_data
,
INMOST_MPI_Comm
comm
);
void
SolverAssignDataFcbiilu2
(
void
*
data
,
void
*
other_data
);
void
SolverSetMatrixFcbiilu2
(
void
*
data
,
void
*
matrix_data
,
bool
same_pattern
,
bool
reuse_preconditioner
);
bool
SolverSolveFcbiilu2
(
void
*
data
,
void
*
rhs_data
,
void
*
sol_data
);
int
SolverIterationNumberFcbiilu2
(
void
*
data
);
double
SolverResidualNormFcbiilu2
(
void
*
data
);
void
SolverAddOtherStatFcbiilu2
(
void
*
data
,
unsigned
int
*
pivmod
,
double
*
prdens
,
double
*
t_prec
,
double
*
t_iter
);
//#endif //USE_SOLVER_FCBIILU2
void
SolverDestroyDataFcbiilu2
(
bcg
**
data
);
void
SolverInitDataFcbiilu2
(
bcg
**
data
,
INMOST_MPI_Comm
comm
,
const
char
*
name
);
void
SolverCopyDataFcbiilu2
(
bcg
**
data
,
bcg
*
other_data
,
INMOST_MPI_Comm
comm
);
void
SolverAssignDataFcbiilu2
(
bcg
*
data
,
bcg
*
other_data
);
void
SolverSetMatrixFcbiilu2
(
bcg
*
data
,
matrix
*
matrix_data
,
bool
same_pattern
,
bool
reuse_preconditioner
);
bool
SolverSolveFcbiilu2
(
bcg
*
data
,
vector
*
rhs_data
,
vector
*
sol_data
);
#endif //SOLVER_FCBIILU2_H_INCLUDED
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