source: trunk/source/global/HEPNumerics/include/G4JTPolynomialSolver.hh@ 1036

Last change on this file since 1036 was 850, checked in by garnier, 17 years ago

geant4.8.2 beta

File size: 4.6 KB
Line 
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26//
27// $Id: G4JTPolynomialSolver.hh,v 1.6 2006/06/29 18:59:49 gunter Exp $
28// GEANT4 tag $Name: HEAD $
29//
30// Class description:
31//
32// G4JTPolynomialSolver implements the Jenkins-Traub algorithm
33// for real polynomial root finding.
34// The solver returns -1, if the leading coefficient is zero,
35// the number of roots found, otherwise.
36//
37// ----------------------------- INPUT --------------------------------
38//
39// op - double precision vector of coefficients in order of
40// decreasing powers
41// degree - integer degree of polynomial
42//
43// ----------------------------- OUTPUT -------------------------------
44//
45// zeror,zeroi - double precision vectors of the
46// real and imaginary parts of the zeros
47//
48// ---------------------------- EXAMPLE -------------------------------
49//
50// G4JTPolynomialSolver trapEq ;
51// G4double coef[8] ;
52// G4double zr[7] , zi[7] ;
53// G4int num = trapEq.FindRoots(coef,7,zr,zi);
54
55// ---------------------------- HISTORY -------------------------------
56//
57// Translated from original TOMS493 Fortran77 routine (ANSI C, by C.Bond).
58// Translated to C++ and adapted to use STL vectors,
59// by Oliver Link (Oliver.Link@cern.ch)
60//
61// --------------------------------------------------------------------
62
63#ifndef G4JTPOLYNOMIALSOLVER_HH
64#define G4JTPOLYNOMIALSOLVER_HH
65
66#include <cmath>
67#include <vector>
68
69#include "globals.hh"
70
71class G4JTPolynomialSolver
72{
73
74 public:
75
76 G4JTPolynomialSolver();
77 ~G4JTPolynomialSolver();
78
79 G4int FindRoots(G4double *op, G4int degree,
80 G4double *zeror, G4double *zeroi);
81
82 private:
83
84 std::vector<G4double> p;
85 std::vector<G4double> qp;
86 std::vector<G4double> k;
87 std::vector<G4double> qk;
88 std::vector<G4double> svk;
89
90 G4double sr;
91 G4double si;
92 G4double u,v;
93 G4double a,b,c,d;
94 G4double a1,a2,a3,a6,a7;
95 G4double e,f,g,h;
96 G4double szr,szi;
97 G4double lzr,lzi;
98 G4int n,nmi;
99
100 /* The following statements set machine constants */
101
102 static const G4double base;
103 static const G4double eta;
104 static const G4double infin;
105 static const G4double smalno;
106 static const G4double are;
107 static const G4double mre;
108 static const G4double lo;
109
110 void Quadratic(G4double a,G4double b1,G4double c,
111 G4double *sr,G4double *si, G4double *lr,G4double *li);
112 void ComputeFixedShiftPolynomial(G4int l2, G4int *nz);
113 void QuadraticPolynomialIteration(G4double *uu,G4double *vv,G4int *nz);
114 void RealPolynomialIteration(G4double *sss, G4int *nz, G4int *iflag);
115 void ComputeScalarFactors(G4int *type);
116 void ComputeNextPolynomial(G4int *type);
117 void ComputeNewEstimate(G4int type,G4double *uu,G4double *vv);
118 void QuadraticSyntheticDivision(G4int n, G4double *u, G4double *v,
119 std::vector<G4double> &p,
120 std::vector<G4double> &q,
121 G4double *a, G4double *b);
122};
123
124#endif
Note: See TracBrowser for help on using the repository browser.