Page Menu
Home
c4science
Search
Configure Global Search
Log In
Files
F102510989
solver_cg_eigen.cc
No One
Temporary
Actions
Download File
Edit File
Delete File
View Transforms
Subscribe
Mute Notifications
Award Token
Subscribers
None
File Metadata
Details
File Info
Storage
Attached
Created
Fri, Feb 21, 11:51
Size
4 KB
Mime Type
text/x-c
Expires
Sun, Feb 23, 11:51 (1 d, 6 h)
Engine
blob
Format
Raw Data
Handle
24338850
Attached To
rMUSPECTRE µSpectre
solver_cg_eigen.cc
View Options
/**
* file solver_cg_eigen.cc
*
* @author Till Junge <till.junge@epfl.ch>
*
* @date 19 Jan 2018
*
* @brief implementation for binding to Eigen's conjugate gradient solver
*
* @section LICENCE
*
* Copyright (C) 2018 Till Junge
*
* µSpectre is free software; you can redistribute it and/or
* modify it under the terms of the GNU General Public License as
* published by the Free Software Foundation, either version 3, or (at
* your option) any later version.
*
* µSpectre is distributed in the hope that it will be useful, but
* WITHOUT ANY WARRANTY; without even the implied warranty of
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
* General Public License for more details.
*
* You should have received a copy of the GNU General Public License
* along with GNU Emacs; see the file COPYING. If not, write to the
* Free Software Foundation, Inc., 59 Temple Place - Suite 330,
* Boston, MA 02111-1307, USA.
*/
#include "solver_cg_eigen.hh"
#include <sstream>
#include <iomanip>
namespace
muSpectre
{
//----------------------------------------------------------------------------//
template
<
class
SolverType
,
Dim_t
DimS
,
Dim_t
DimM
>
SolverEigen
<
SolverType
,
DimS
,
DimM
>::
SolverEigen
(
Sys_t
&
sys
,
Real
tol
,
Uint
maxiter
,
bool
verbose
)
:
Parent
(
sys
,
tol
,
maxiter
,
verbose
),
adaptor
{
sys
.
get_adaptor
()},
solver
{}
{}
//----------------------------------------------------------------------------//
template
<
class
SolverType
,
Dim_t
DimS
,
Dim_t
DimM
>
void
SolverEigen
<
SolverType
,
DimS
,
DimM
>::
initialise
()
{
this
->
solver
.
setTolerance
(
this
->
tol
);
this
->
solver
.
setMaxIterations
(
this
->
maxiter
);
this
->
solver
.
compute
(
this
->
adaptor
);
}
//----------------------------------------------------------------------------//
template
<
class
SolverType
,
Dim_t
DimS
,
Dim_t
DimM
>
typename
SolverEigen
<
SolverType
,
DimS
,
DimM
>::
SolvVectorOut
SolverEigen
<
SolverType
,
DimS
,
DimM
>::
solve
(
const
SolvVectorInC
rhs
,
SolvVectorIn
x_0
)
{
auto
&
this_solver
=
static_cast
<
SolverType
&>
(
*
this
);
SolvVectorOut
retval
=
this
->
solver
.
solveWithGuess
(
rhs
,
x_0
);
this
->
counter
+=
this
->
solver
.
iterations
();
if
(
this
->
solver
.
info
()
!=
Eigen
::
Success
)
{
std
::
stringstream
err
{};
err
<<
this_solver
.
name
()
<<
" has not converged,"
<<
" After "
<<
this
->
solver
.
iterations
()
<<
" steps, the solver "
<<
" FAILED with |r|/|b| = "
<<
std
::
setw
(
15
)
<<
this
->
solver
.
error
()
<<
", cg_tol = "
<<
this
->
tol
<<
std
::
endl
;
throw
ConvergenceError
(
err
.
str
());
}
if
(
this
->
verbose
)
{
std
::
cout
<<
" After "
<<
this
->
solver
.
iterations
()
<<
" "
<<
this_solver
.
name
()
<<
" steps, |r|/|b| = "
<<
std
::
setw
(
15
)
<<
this
->
solver
.
error
()
<<
", cg_tol = "
<<
this
->
tol
<<
std
::
endl
;
}
return
retval
;
}
/* ---------------------------------------------------------------------- */
template
<
class
SolverType
,
Dim_t
DimS
,
Dim_t
DimM
>
typename
SolverEigen
<
SolverType
,
DimS
,
DimM
>::
Tg_req_t
SolverEigen
<
SolverType
,
DimS
,
DimM
>::
get_tangent_req
()
const
{
return
tangent_requirement
;
}
template
class
SolverEigen
<
SolverCGEigen
<
twoD
,
twoD
>
,
twoD
,
twoD
>
;
template
class
SolverEigen
<
SolverCGEigen
<
threeD
,
threeD
>
,
threeD
,
threeD
>
;
template
class
SolverCGEigen
<
twoD
,
twoD
>
;
template
class
SolverCGEigen
<
threeD
,
threeD
>
;
template
class
SolverEigen
<
SolverGMRESEigen
<
twoD
,
twoD
>
,
twoD
,
twoD
>
;
template
class
SolverEigen
<
SolverGMRESEigen
<
threeD
,
threeD
>
,
threeD
,
threeD
>
;
template
class
SolverGMRESEigen
<
twoD
,
twoD
>
;
template
class
SolverGMRESEigen
<
threeD
,
threeD
>
;
template
class
SolverEigen
<
SolverBiCGSTABEigen
<
twoD
,
twoD
>
,
twoD
,
twoD
>
;
template
class
SolverEigen
<
SolverBiCGSTABEigen
<
threeD
,
threeD
>
,
threeD
,
threeD
>
;
template
class
SolverBiCGSTABEigen
<
twoD
,
twoD
>
;
template
class
SolverBiCGSTABEigen
<
threeD
,
threeD
>
;
template
class
SolverEigen
<
SolverDGMRESEigen
<
twoD
,
twoD
>
,
twoD
,
twoD
>
;
template
class
SolverEigen
<
SolverDGMRESEigen
<
threeD
,
threeD
>
,
threeD
,
threeD
>
;
template
class
SolverDGMRESEigen
<
twoD
,
twoD
>
;
template
class
SolverDGMRESEigen
<
threeD
,
threeD
>
;
template
class
SolverEigen
<
SolverMINRESEigen
<
twoD
,
twoD
>
,
twoD
,
twoD
>
;
template
class
SolverEigen
<
SolverMINRESEigen
<
threeD
,
threeD
>
,
threeD
,
threeD
>
;
template
class
SolverMINRESEigen
<
twoD
,
twoD
>
;
template
class
SolverMINRESEigen
<
threeD
,
threeD
>
;
}
// muSpectre
Event Timeline
Log In to Comment