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_angul_C[] = "$Header: /cvsroot/Lorene/C++/Source/Etoile/et_bin_nsbh_angul.C,v 1.3 2005/08/31 09:53:19 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_angul() 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 Cmp yabs (mp) ;
00063 yabs = mp.ya ;
00064
00065
00066 Cmp G_square_cmp (-pow(confpsi(), 4)/nnn()/nnn()*(xabs*xabs+yabs*yabs)) ;
00067 G_square_cmp.std_base_scal() ;
00068 double dG_square = G_square_cmp.dsdx()(0,0,0,0) ;
00069
00070
00071 Cmp G_single_cmp (-2*pow(confpsi(), 4)/nnn()/nnn()*(shift(1)*xabs-shift(0)*yabs)) ;
00072 G_single_cmp.std_base_scal() ;
00073 double dG_single = G_single_cmp.dsdx()(0,0,0,0) ;
00074
00075
00076 Cmp G_const_cmp = 1-pow(confpsi(), 4)/nnn()/nnn()*(shift(0)*shift(0) + shift(1)*shift(1) + shift(2)*shift(2)) ;
00077 G_const_cmp.std_base_scal() ;
00078 double dG_const = G_const_cmp.dsdx()(0,0,0,0) ;
00079
00080
00081 double G_square = G_square_cmp (0,0,0,0) ;
00082 double G_single = G_single_cmp (0,0,0,0) ;
00083 double G_const = G_const_cmp (0,0,0,0) ;
00084
00085
00086 double a_coef = dG_square + 2*G_square*part_dx ;
00087 double b_coef = dG_single + 2*G_single*part_dx ;
00088 double c_coef = dG_const + 2*G_const *part_dx;
00089
00090 double determinant = b_coef*b_coef - 4*a_coef*c_coef ;
00091
00092 bool soluce = true ;
00093 double res ;
00094
00095 if (determinant <0) {
00096 soluce = false ;
00097 }
00098
00099 else {
00100 double sol_un = (-b_coef - sqrt(determinant))/2/a_coef ;
00101 double sol_deux = (-b_coef + sqrt(determinant))/2/a_coef ;
00102
00103 bool signe_un = (sol_un >0) ? true : false ;
00104 bool signe_deux = (sol_deux >0) ? true : false ;
00105
00106
00107 if (signe_un == signe_deux) {
00108 soluce = false ;
00109 }
00110 else {
00111 res = (signe_un) ? sol_un : sol_deux ;
00112 }
00113 }
00114
00115 if (soluce == false)
00116 return 0 ;
00117 else
00118 return res ;
00119 }