The Gaudi Framework  v28r3 (cc1cf868)
EqSolver.cpp
Go to the documentation of this file.
1 // Include files
2 
3 #include <stdlib.h>
4 #include <stdio.h>
5 // from Gaudi
7 //from GSL
8 #include "gsl/gsl_vector.h"
9 #include "gsl/gsl_multimin.h"
10 #include "gsl/gsl_math.h"
11 #include "gsl/gsl_matrix.h"
12 #include "gsl/gsl_linalg.h"
13 #include "gsl/gsl_blas.h"
14 #include "gsl/gsl_errno.h"
15 
16 // local
17 #include "EqSolver.h"
18 
26 // Declaration of the Tool Factory
28 
29 
30 // ============================================================================
31 // ============================================================================
32 EqSolver::EqSolverMisc::EqSolverMisc
33 ( const EqSolver::Equations& funcs ,
34  EqSolver::Arg& arg )
35  : m_argum ( arg )
36  , m_eqs ( &funcs )
37 {
38  const size_t N = funcs.size () ;
39  for( size_t i = 0 ; i < N ; ++i ) {
40  Equations last;
41  for( size_t j = 0 ; j < N ; ++j ) {
42  Genfun::GENFUNCTION fij = funcs[i]->partial(j);
43  last.push_back( fij.clone() ) ;
44  }
45  m_jac.push_back( last );
46  }
47 }
48 // ============================================================================
49 
50 // ============================================================================
52 {
53  while( !m_jac.empty() )
54  {
55  Equations& last = m_jac.back() ;
56  while( !last.empty() ) {
57  delete last.back() ;
58  last.pop_back () ;
59  }
60  m_jac.pop_back();
61  }
62 }
63 // ============================================================================
64 
65 namespace
66 {
67  using namespace Genfun;
68 
71  int fun_gsl ( const gsl_vector* v ,
72  void* params ,
73  gsl_vector* f )
74  {
75  EqSolver::EqSolverMisc* local =
76  static_cast <EqSolver::EqSolverMisc*> (params);
77  const EqSolver::Equations& eqs = *(local->equations()) ;
78  Argument& arg = local->argument() ;
79 
80  for ( unsigned int i = 0; i < v->size; ++i) {
81  arg[i] = gsl_vector_get (v, i);
82  }
83 
84  for ( unsigned int i = 0; i < v->size; ++i) {
85  gsl_vector_set(f, i, (*eqs[i])(arg));
86  }
87  return GSL_SUCCESS;
88  }
89 
90 
93  int dfun_gsl ( const gsl_vector* v ,
94  void* params ,
95  gsl_matrix* df )
96  {
97  EqSolver::EqSolverMisc* local =
98  static_cast <EqSolver::EqSolverMisc*> (params);
99  const EqSolver::Jacobi& jac = local -> jacobi() ;
100  Argument& arg = local -> argument() ;
101 
102  for ( unsigned int i = 0; i < v->size; ++i) {
103  arg[i] = gsl_vector_get (v, i);
104  }
105 
106 
107  for( unsigned int i = 0 ; i < v->size ; ++i ) {
108  for (unsigned int j = 0; j < v->size; ++j) {
109  Genfun::GENFUNCTION f = *(jac[i][j]) ;
110  gsl_matrix_set ( df , i , j , f ( arg ) ) ;
111  }
112  }
113  return GSL_SUCCESS;
114  }
115 
119  int fdfun_gsl ( const gsl_vector* v ,
120  void* params ,
121  gsl_vector* f ,
122  gsl_matrix* df )
123  {
124  EqSolver::EqSolverMisc* local =
125  static_cast <EqSolver::EqSolverMisc*> (params);
126  const EqSolver::Equations& eqs = *(local->equations()) ;
127  const EqSolver::Jacobi& jac = local->jacobi() ;
128  Argument& arg = local->argument() ;
129 
130  for ( unsigned int i = 0; i < v->size; ++i) {
131  arg[i] = gsl_vector_get (v, i);
132  }
133 
134  for( unsigned int i = 0 ; i < v->size ; ++i ) {
135  gsl_vector_set(f, i, (*eqs[i])(arg ) ) ;
136  for (unsigned int j = 0; j < v->size; ++j) {
137  Genfun::GENFUNCTION f1 = *(jac[i][j]) ;
138  gsl_matrix_set ( df , i , j , f1(arg) ) ;
139  }
140  }
141  return GSL_SUCCESS;
142  }
143 }
144 
145 //=============================================================================
150 StatusCode EqSolver::solver (const Equations& funcs ,
151  Arg& arg ) const
152 //=============================================================================
153 {
154  using namespace Genfun;
155 
156  gsl_vector_view vect = gsl_vector_view_array ( &arg[0] ,
157  arg.dimension() );
158  EqSolverMisc local (funcs, arg);
159 
160  const gsl_multiroot_fdfsolver_type *T = m_type;
161  gsl_multiroot_fdfsolver *s = nullptr;
162  int status = 0 ;
163  size_t iter = 0 ;
164 
165  gsl_multiroot_function_fdf function;
166 
167  function.f = &fun_gsl;
168  function.df = &dfun_gsl;
169  function.fdf = &fdfun_gsl;
170  function.n = vect.vector.size;
171  function.params = &local;
172 
173  s = gsl_multiroot_fdfsolver_alloc(T, vect.vector.size);
174  gsl_multiroot_fdfsolver_set (s, &function, &vect.vector);
175 
176  for (iter = 0; iter < m_max_iter; ++iter) {
177  status = gsl_multiroot_fdfsolver_iterate (s);
178  if (status) {
179  return Error
180  ("Error from gsl_gsl_multiroot_fdfsolver_iterate '"
181  + std::string(gsl_strerror(status)) + "'") ;
182  }
183 
184  status = gsl_multiroot_test_residual (s->f,
186 
187  if ( status != GSL_CONTINUE ) { break; }
188  }
189 
190  for (unsigned int i = 0; i < vect.vector.size; ++i) {
191  gsl_vector_set (&vect.vector, i, gsl_vector_get (s->x, i));
192  }
193 
194  if (status == GSL_SUCCESS) {
195  debug()
196  << "We stopped in the method on the " << iter
197  << " iteration (we have maximum " << m_max_iter
198  << " iterations)" << endmsg;
199  } else if (status == GSL_CONTINUE && iter <= m_max_iter ) {
200  return Error ( "Method finished with '"
201  + std::string(gsl_strerror(status))
202  + "' error" );
203  } else {
204  return Error ( "Method finished with '" +
205  std::string(gsl_strerror(status))
206  + "' error" );
207  }
208 
209  gsl_multiroot_fdfsolver_free (s);
210 
211  if (status) {
212  return Error ( "Method finished with '"
213  + std::string(gsl_strerror(status))
214  + "' error" );
215  }
216 
217  return GSL_SUCCESS;
218 }
219 
220 //=============================================================================
222 {
224  if( sc.isFailure() ) {
225  return Error ("Could not initiliaze base class GaudiTool", sc);
226  }
227 
228  /* The algorithm for multiional root-finding
229  (solving nonlinear systems with n equations in n unknowns)*/
230 
231  if( "fdfsolver_hybridsj" == m_algType ) {
232  m_type = gsl_multiroot_fdfsolver_hybridsj ;
233  debug()
234  << "Root finding algorithm to be used: "
235  << "'gsl_multiroot_fdfsolver_hybridsj'"
236  << endmsg;
237  } else if ( "fdfsolver_hybridj" == m_algType ) {
238  m_type = gsl_multiroot_fdfsolver_hybridj ;
239  debug()
240  << "Root finding algorithm to be used: "
241  << "'gsl_multiroot_fdfsolver_hybridj'"
242  << endmsg;
243  } else if ( "fdfsolver_newton" == m_algType ) {
244  m_type = gsl_multiroot_fdfsolver_newton ;
245  debug()
246  << "Root findind algorithm to be used: "
247  << "'gsl_multiroot_fdfsolver_newton'"
248  << endmsg;
249  } else if ( "fdfsolver_gnewton" == m_algType ) {
250  m_type = gsl_multiroot_fdfsolver_gnewton ;
251  debug()
252  << "Root findind algorithm to be used: "
253  << "'gsl_multiroot_fdfsolver_gnewton'"
254  << endmsg;
255  } else {
256  return Error(" Unknown algorithm type '"
257  + std::string(m_algType) + "'");
258  }
259 
260  return StatusCode::SUCCESS;
261 }
262 //=============================================================================
T empty(T...args)
const Jacobi & jacobi() const
Definition: EqSolver.h:53
bool isFailure() const
Test for a status code of FAILURE.
Definition: StatusCode.h:84
#define DECLARE_COMPONENT(type)
Definition: PluginService.h:36
STL class.
Gaudi::Property< std::string > m_algType
Definition: EqSolver.h:80
int N
Definition: IOTest.py:90
The simplest concrete implementation of IEqSolver interface.
Definition: EqSolver.h:23
Gaudi::Property< double > m_max_iter
Definition: EqSolver.h:82
StatusCode initialize() override
Overriding initialize.
Definition: EqSolver.cpp:221
const Arg & argument() const
Definition: EqSolver.h:50
This class is used for returning status codes from appropriate routines.
Definition: StatusCode.h:26
const Equations * equations() const
Definition: EqSolver.h:52
T pop_back(T...args)
StatusCode solver(const Equations &funcs, Arg &arg) const override
Solving nonlinear system with N equations in N unknowns of the function "GenFunc".
Definition: EqSolver.cpp:150
MsgStream & debug() const
shortcut for the method msgStream(MSG::DEBUG)
StatusCode Error(const std::string &msg, const StatusCode st=StatusCode::FAILURE, const size_t mx=10) const
Print the error message and return with the given StatusCode.
Definition: GaudiTool.h:671
T back(T...args)
string s
Definition: gaudirun.py:245
const gsl_multiroot_fdfsolver_type * m_type
Definition: EqSolver.h:85
CLHEP.
Definition: IEqSolver.h:13
MsgStream & endmsg(MsgStream &s)
MsgStream Modifier: endmsg. Calls the output method of the MsgStream.
Definition: MsgStream.h:244
Gaudi::Property< double > m_norm_residual
Definition: EqSolver.h:83
StatusCode initialize() override
standard initialization method
Definition: GaudiTool.cpp:157