Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4JTPolynomialSolver.hh
Go to the documentation of this file.
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$
28//
29// Class description:
30//
31// G4JTPolynomialSolver implements the Jenkins-Traub algorithm
32// for real polynomial root finding.
33// The solver returns -1, if the leading coefficient is zero,
34// the number of roots found, otherwise.
35//
36// ----------------------------- INPUT --------------------------------
37//
38// op - double precision vector of coefficients in order of
39// decreasing powers
40// degree - integer degree of polynomial
41//
42// ----------------------------- OUTPUT -------------------------------
43//
44// zeror,zeroi - double precision vectors of the
45// real and imaginary parts of the zeros
46//
47// ---------------------------- EXAMPLE -------------------------------
48//
49// G4JTPolynomialSolver trapEq ;
50// G4double coef[8] ;
51// G4double zr[7] , zi[7] ;
52// G4int num = trapEq.FindRoots(coef,7,zr,zi);
53
54// ---------------------------- HISTORY -------------------------------
55//
56// Translated from original TOMS493 Fortran77 routine (ANSI C, by C.Bond).
57// Translated to C++ and adapted to use STL vectors,
58// by Oliver Link ([email protected])
59//
60// --------------------------------------------------------------------
61
62#ifndef G4JTPOLYNOMIALSOLVER_HH
63#define G4JTPOLYNOMIALSOLVER_HH
64
65#include <cmath>
66#include <vector>
67
68#include "globals.hh"
69
71{
72
73 public:
74
77
78 G4int FindRoots(G4double *op, G4int degree,
79 G4double *zeror, G4double *zeroi);
80
81 private:
82
83 std::vector<G4double> p;
84 std::vector<G4double> qp;
85 std::vector<G4double> k;
86 std::vector<G4double> qk;
87 std::vector<G4double> svk;
88
89 G4double sr;
90 G4double si;
91 G4double u,v;
92 G4double a,b,c,d;
93 G4double a1,a2,a3,a6,a7;
94 G4double e,f,g,h;
95 G4double szr,szi;
96 G4double lzr,lzi;
97 G4int n,nmi;
98
99 /* The following statements set machine constants */
100
101 static const G4double base;
102 static const G4double eta;
103 static const G4double infin;
104 static const G4double smalno;
105 static const G4double are;
106 static const G4double mre;
107 static const G4double lo;
108
109 void Quadratic(G4double a,G4double b1,G4double c,
110 G4double *sr,G4double *si, G4double *lr,G4double *li);
111 void ComputeFixedShiftPolynomial(G4int l2, G4int *nz);
112 void QuadraticPolynomialIteration(G4double *uu,G4double *vv,G4int *nz);
113 void RealPolynomialIteration(G4double *sss, G4int *nz, G4int *iflag);
114 void ComputeScalarFactors(G4int *type);
115 void ComputeNextPolynomial(G4int *type);
116 void ComputeNewEstimate(G4int type,G4double *uu,G4double *vv);
117 void QuadraticSyntheticDivision(G4int n, G4double *u, G4double *v,
118 std::vector<G4double> &p,
119 std::vector<G4double> &q,
120 G4double *a, G4double *b);
121};
122
123#endif
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
G4int FindRoots(G4double *op, G4int degree, G4double *zeror, G4double *zeroi)