00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023 char bhole_coal_C[] = "$Header: /cvsroot/Lorene/C++/Source/Bhole_binaire/bhole_coal.C,v 1.6 2005/08/31 09:48:00 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
00047
00048
00049
00050
00051
00052
00053
00054
00055
00056
00057
00058
00059
00060
00061
00062
00063
00064
00065
00066
00067
00068
00069
00070
00071
00072
00073
00074
00075
00076
00077
00078
00079
00080
00081
00082
00083
00084
00085
00086
00087
00088
00089
00090
00091
00092
00093
00094
00095
00096
00097
00098
00099
00100
00101
00102
00103
00104 #include <math.h>
00105 #include <stdlib.h>
00106
00107
00108 #include "tenseur.h"
00109 #include "bhole.h"
00110
00111 void Bhole_binaire::set_statiques (double precis, double relax) {
00112
00113 int nz_un = hole1.mp.get_mg()->get_nzone() ;
00114 int nz_deux = hole2.mp.get_mg()->get_nzone() ;
00115
00116 set_omega(0) ;
00117 init_bhole_binaire() ;
00118
00119 int indic = 1 ;
00120 int conte = 0 ;
00121
00122 cout << "TROUS STATIQUES : " << endl ;
00123 while (indic == 1) {
00124 Cmp lapse_un_old (hole1.get_n_auto()()) ;
00125 Cmp lapse_deux_old (hole2.get_n_auto()()) ;
00126
00127 solve_psi (precis, relax) ;
00128 solve_lapse (precis, relax) ;
00129
00130 double erreur = 0 ;
00131 Tbl diff_un (diffrelmax (lapse_un_old, hole1.get_n_auto()())) ;
00132 for (int i=1 ; i<nz_un ; i++)
00133 if (diff_un(i) > erreur)
00134 erreur = diff_un(i) ;
00135
00136 Tbl diff_deux (diffrelmax (lapse_deux_old, hole2.get_n_auto()())) ;
00137 for (int i=1 ; i<nz_deux ; i++)
00138 if (diff_deux(i) > erreur)
00139 erreur = diff_deux(i) ;
00140
00141
00142 cout << "PAS TOTAL : " << conte << " DIFFERENCE : " << erreur << endl ;
00143
00144 if (erreur < precis)
00145 indic = -1 ;
00146 conte ++ ;
00147 }
00148 }
00149
00150 void Bhole_binaire::coal (double precis, double relax, int nbre_ome, double seuil_search, double m1, double m2, const int sortie) {
00151
00152 assert (omega == 0) ;
00153 int nz1 = hole1.mp.get_mg()->get_nzone() ;
00154 int nz2 = hole1.mp.get_mg()->get_nzone() ;
00155
00156
00157 double distance = hole1.mp.get_ori_x()-hole2.mp.get_ori_x() ;
00158 set_pos_axe (distance*(hole2.rayon/(hole1.rayon+hole2.rayon))) ;
00159 double scale_linear = (hole1.rayon + hole2.rayon)/2.*distance ;
00160
00161
00162 double angulaire = sqrt((hole1.rayon+hole2.rayon)/distance/distance/distance) ;
00163
00164 int indic = 1 ;
00165 int conte = 0 ;
00166
00167 char name_iteration[40] ;
00168 char name_correction[40] ;
00169 char name_viriel[40] ;
00170 char name_ome [40] ;
00171 char name_linear[40] ;
00172 char name_axe[40] ;
00173 char name_error_m1[40] ;
00174 char name_error_m2[40] ;
00175 char name_r1[40] ;
00176 char name_r2[40] ;
00177
00178 sprintf(name_iteration, "ite.dat") ;
00179 sprintf(name_correction, "cor.dat") ;
00180 sprintf(name_viriel, "vir.dat") ;
00181 sprintf(name_ome, "ome.dat") ;
00182 sprintf(name_linear, "linear.dat") ;
00183 sprintf(name_axe, "axe.dat") ;
00184 sprintf(name_error_m1, "error_m1.dat") ;
00185 sprintf(name_error_m2, "error_m2.dat") ;
00186 sprintf(name_r1, "r1.dat") ;
00187 sprintf(name_r2, "r2.dat") ;
00188
00189 ofstream fiche_iteration(name_iteration) ;
00190 fiche_iteration.precision(8) ;
00191
00192 ofstream fiche_correction(name_correction) ;
00193 fiche_correction.precision(8) ;
00194
00195 ofstream fiche_viriel(name_viriel) ;
00196 fiche_viriel.precision(8) ;
00197
00198 ofstream fiche_ome(name_ome) ;
00199 fiche_ome.precision(8) ;
00200
00201 ofstream fiche_linear(name_linear) ;
00202 fiche_linear.precision(8) ;
00203
00204 ofstream fiche_axe(name_axe) ;
00205 fiche_axe.precision(8) ;
00206
00207 ofstream fiche_error_m1 (name_error_m1) ;
00208 fiche_error_m1.precision(8) ;
00209
00210 ofstream fiche_error_m2 (name_error_m2) ;
00211 fiche_error_m2.precision(8) ;
00212
00213 ofstream fiche_r1 (name_r1) ;
00214 fiche_r1.precision(8) ;
00215
00216 ofstream fiche_r2 (name_r2) ;
00217 fiche_r2.precision(8) ;
00218
00219
00220 cout << "OMEGA AUGMENTE A LA MAIN." << endl ;
00221 double homme = 0 ;
00222 for (int pas = 0 ; pas <nbre_ome ; pas ++) {
00223
00224 homme += angulaire/nbre_ome ;
00225 set_omega (homme) ;
00226
00227 Cmp shift_un_old (hole1.get_shift_auto()(0)) ;
00228 Cmp shift_deux_old (hole2.get_shift_auto()(0)) ;
00229
00230 solve_shift (precis, relax) ;
00231 fait_tkij() ;
00232
00233 solve_psi (precis, relax) ;
00234 solve_lapse (precis, relax) ;
00235
00236 double erreur = 0 ;
00237 Tbl diff_un (diffrelmax (shift_un_old, hole1.get_shift_auto()(0))) ;
00238 for (int i=1 ; i<nz1 ; i++)
00239 if (diff_un(i) > erreur)
00240 erreur = diff_un(i) ;
00241
00242 Tbl diff_deux (diffrelmax (shift_deux_old, hole2.get_shift_auto()(0))) ;
00243 for (int i=1 ; i<nz2 ; i++)
00244 if (diff_deux(i) > erreur)
00245 erreur = diff_deux(i) ;
00246
00247 double error_viriel = viriel() ;
00248 double error_linear = linear_momentum_systeme_inf()(1)/scale_linear ;
00249 double error_m1 = 1.-sqrt(hole1.area()/16./M_PI)/m1 ;
00250 double error_m2 = 1.-sqrt(hole2.area()/16./M_PI)/m2 ;
00251 double r1 = hole1.mp.val_r(0, 1, 0, 0) ;
00252 double r2 = hole2.mp.val_r(0, 1, 0, 0) ;
00253
00254 if (sortie != 0) {
00255 fiche_iteration << conte << " " << erreur << endl ;
00256 fiche_correction << conte << " " << hole1.get_regul() << " " << hole2.get_regul() << endl ;
00257 fiche_viriel << conte << " " << error_viriel << endl ;
00258 fiche_ome << conte << " " << homme << endl ;
00259 fiche_linear << conte << " " << error_linear << endl ;
00260 fiche_axe << conte << " " << pos_axe << endl ;
00261 fiche_error_m1 << conte << " " << error_m1 << endl ;
00262 fiche_error_m2 << conte << " " << error_m2 << endl ;
00263 fiche_r1 << conte << " " << r1 << endl ;
00264 fiche_r2 << conte << " " << r2 << endl ;
00265 }
00266
00267 cout << "PAS TOTAL : " << conte << " DIFFERENCE : " << erreur << endl ;
00268 conte ++ ;
00269 }
00270
00271
00272 cout << "OMEGA VARIABLE" << endl ;
00273 indic = 1 ;
00274 bool scale = false ;
00275
00276 while (indic == 1) {
00277
00278 Cmp shift_un_old (hole1.get_shift_auto()(0)) ;
00279 Cmp shift_deux_old (hole2.get_shift_auto()(0)) ;
00280
00281 solve_shift (precis, relax) ;
00282 fait_tkij() ;
00283
00284 solve_psi (precis, relax) ;
00285 solve_lapse (precis, relax) ;
00286
00287 double erreur = 0 ;
00288 Tbl diff_un (diffrelmax (shift_un_old, hole1.get_shift_auto()(0))) ;
00289 for (int i=1 ; i<nz1 ; i++)
00290 if (diff_un(i) > erreur)
00291 erreur = diff_un(i) ;
00292
00293 Tbl diff_deux (diffrelmax (shift_deux_old, hole2.get_shift_auto()(0))) ;
00294 for (int i=1 ; i<nz2 ; i++)
00295 if (diff_deux(i) > erreur)
00296 erreur = diff_deux(i) ;
00297
00298 double error_viriel = viriel() ;
00299 double error_linear = linear_momentum_systeme_inf()(1)/scale_linear ;
00300 double error_m1 = 1.-sqrt(hole1.area()/16./M_PI)/m1 ;
00301 double error_m2 = 1.-sqrt(hole2.area()/16./M_PI)/m2 ;
00302 double r1 = hole1.mp.val_r(0, 1, 0, 0) ;
00303 double r2 = hole2.mp.val_r(0, 1, 0, 0) ;
00304
00305 if (sortie != 0) {
00306 fiche_iteration << conte << " " << erreur << endl ;
00307 fiche_correction << conte << " " << hole1.regul << " " << hole2.regul << endl ;
00308 fiche_viriel << conte << " " << error_viriel << endl ;
00309 fiche_ome << conte << " " << omega << endl ;
00310 fiche_linear << conte << " " << error_linear << endl ;
00311 fiche_axe << conte << " " << pos_axe << endl ;
00312 fiche_error_m1 << conte << " " << error_m1 << endl ;
00313 fiche_error_m2 << conte << " " << error_m2 << endl ;
00314 fiche_r1 << conte << " " << r1 << endl ;
00315 fiche_r2 << conte << " " << r2 << endl ;
00316 }
00317
00318
00319 if (erreur <= seuil_search)
00320 scale = true ;
00321 if (scale) {
00322 double scaling_ome = pow((2-error_viriel)/(2-2*error_viriel), 1.) ;
00323 set_omega (omega*scaling_ome) ;
00324
00325 double scaling_axe = pow((2-error_linear)/(2-2*error_linear), 0.1) ;
00326 set_pos_axe (pos_axe*scaling_axe) ;
00327
00328 double scaling_r1 = pow((2-error_m1)/(2-2*error_m1), 0.1) ;
00329 hole1.mp.homothetie_interne(scaling_r1) ;
00330
00331 double scaling_r2 = pow((2-error_m2)/(2-2*error_m2), 0.1) ;
00332 hole2.mp.homothetie_interne(scaling_r2) ;
00333 }
00334
00335 cout << "PAS TOTAL : " << conte << " DIFFERENCE : " << erreur << endl ;
00336 if (erreur < precis)
00337 indic = -1 ;
00338 conte ++ ;
00339 }
00340
00341 fiche_iteration.close() ;
00342 fiche_correction.close() ;
00343 fiche_viriel.close() ;
00344 fiche_ome.close() ;
00345 fiche_linear.close() ;
00346 fiche_axe.close() ;
00347 fiche_error_m1.close() ;
00348 fiche_error_m2.close() ;
00349 fiche_r1.close() ;
00350 fiche_r2.close() ;
00351 }