00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023 char et_bin_nsbh_axe_C[] = "$Header: /cvsroot/Lorene/C++/Source/Etoile/et_bin_nsbh_axe.C,v 1.3 2005/08/31 09:53:58 m_saijo Exp $" ;
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046 #include <math.h>
00047
00048
00049 #include "et_bin_nsbh.h"
00050 #include "graphique.h"
00051
00052 double Et_bin_nsbh::compute_axe(double omega) const {
00053
00054
00055 Cmp dx_mu (nnn().dsdx() / nnn());
00056 Cmp dx_loggam (loggam().dsdx()) ;
00057 double part_dx = dx_mu(0,0,0,0) + dx_loggam(0,0,0,0) ;
00058
00059 Cmp xabs (mp) ;
00060 xabs = mp.xa ;
00061
00062
00063 Cmp C_cmp (-pow(confpsi(), 4)/nnn()/nnn()*omega*omega) ;
00064 C_cmp.std_base_scal() ;
00065 double dC = C_cmp.dsdx()(0,0,0,0) ;
00066
00067
00068 Cmp B_cmp (2*pow(confpsi(), 4)/nnn()/nnn()*shift(1)*omega) ;
00069 B_cmp.std_base_scal() ;
00070 double dB = B_cmp.dsdx()(0,0,0,0) ;
00071
00072
00073 Cmp A_cmp = 1-pow(confpsi(), 4)/nnn()/nnn()*(shift(0)*shift(0) + shift(1)*shift(1) + shift(2)*shift(2)) ;
00074 A_cmp.std_base_scal() ;
00075 double dA = A_cmp.dsdx()(0,0,0,0) ;
00076
00077
00078 double A = A_cmp (0,0,0,0) ;
00079 double B = B_cmp (0,0,0,0) ;
00080 double C = C_cmp (0,0,0,0) ;
00081
00082
00083 double a_coef = dC + 2*C*part_dx ;
00084 double b_coef = dB + 2*C + 2*B*part_dx ;
00085 double c_coef = dA + B + 2*A*part_dx;
00086
00087 cout << a_coef << " " << b_coef << " " << c_coef << endl ;
00088
00089 double determinant = b_coef*b_coef - 4*a_coef*c_coef ;
00090
00091 if (determinant <0) {
00092 cout << "No solution for Xabs found in Et_bin_nsbh_compute_axe !" << endl ;
00093 abort() ;
00094 }
00095
00096 double sol_un = (-b_coef - sqrt(determinant))/2/a_coef ;
00097 double sol_deux = (-b_coef + sqrt(determinant))/2/a_coef ;
00098
00099 bool signe_un = (sol_un >0) ? true : false ;
00100 bool signe_deux = (sol_deux >0) ? true : false ;
00101
00102 double res ;
00103
00104 if (signe_un == signe_deux) {
00105 cout << "To many solutions in Et_bin_nsbh_compute_axe !" << endl ;
00106 abort() ;
00107 }
00108 else {
00109 res = (signe_un) ? sol_un : sol_deux ;
00110 }
00111
00112 return res ;
00113 }