File: bhole_coal.C

package info (click to toggle)
lorene 0.0.0~cvs20161116%2Bdfsg-1
  • links: PTS, VCS
  • area: main
  • in suites: bookworm, bullseye, buster, stretch
  • size: 26,444 kB
  • ctags: 13,953
  • sloc: cpp: 212,946; fortran: 21,645; makefile: 1,750; sh: 4
file content (359 lines) | stat: -rw-r--r-- 10,970 bytes parent folder | download | duplicates (2)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
/*
 *   Copyright (c) 2000-2001 Philippe Grandclement
 *
 *   This file is part of LORENE.
 *
 *   LORENE is free software; you can redistribute it and/or modify
 *   it under the terms of the GNU General Public License as published by
 *   the Free Software Foundation; either version 2 of the License, or
 *   (at your option) any later version.
 *
 *   LORENE is distributed in the hope that it will be useful,
 *   but WITHOUT ANY WARRANTY; without even the implied warranty of
 *   MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
 *   GNU General Public License for more details.
 *
 *   You should have received a copy of the GNU General Public License
 *   along with LORENE; if not, write to the Free Software
 *   Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
 *
 */


char bhole_coal_C[] = "$Header: /cvsroot/Lorene/C++/Source/Bhole_binaire/bhole_coal.C,v 1.8 2014/10/13 08:52:40 j_novak Exp $" ;

/*
 * $Id: bhole_coal.C,v 1.8 2014/10/13 08:52:40 j_novak Exp $
 * $Log: bhole_coal.C,v $
 * Revision 1.8  2014/10/13 08:52:40  j_novak
 * Lorene classes and functions now belong to the namespace Lorene.
 *
 * Revision 1.7  2014/10/06 15:12:58  j_novak
 * Modified #include directives to use c++ syntax.
 *
 * Revision 1.6  2005/08/31 09:48:00  m_saijo
 * Delete one <math.h>
 *
 * Revision 1.5  2005/08/31 09:06:18  p_grandclement
 * add math.h in bhole_coal.C
 *
 * Revision 1.4  2005/08/29 15:10:14  p_grandclement
 * Addition of things needed :
 *   1) For BBH with different masses
 *   2) Provisory files for the mixted binaries (Bh and NS) : THIS IS NOT
 *   WORKING YET !!!
 *
 * Revision 1.3  2003/11/13 13:43:54  p_grandclement
 * Addition of things needed for Bhole::update_metric (const Etoile_bin&, double, double)
 *
 * Revision 1.2  2002/10/16 14:36:32  j_novak
 * Reorganization of #include instructions of standard C++, in order to
 * use experimental version 3 of gcc.
 *
 * Revision 1.1.1.1  2001/11/20 15:19:28  e_gourgoulhon
 * LORENE
 *
 * Revision 2.15  2001/05/07  09:12:17  phil
 * *** empty log message ***
 *
 * Revision 2.14  2001/04/26  12:23:06  phil
 * *** empty log message ***
 *
 * Revision 2.13  2001/04/26  12:04:17  phil
 * *** empty log message ***
 *
 * Revision 2.12  2001/03/22  10:49:42  phil
 * *** empty log message ***
 *
 * Revision 2.11  2001/02/28  13:23:54  phil
 * vire kk_auto
 *
 * Revision 2.10  2001/01/29  14:31:04  phil
 * ajout tuype rotation
 *
 * Revision 2.9  2001/01/22  09:29:34  phil
 * vire convergence vers bare masse
 *
 * Revision 2.8  2001/01/10  09:31:52  phil
 * ajoute fait_kk_auto
 *
 * Revision 2.7  2000/12/20  15:02:57  phil
 * *** empty log message ***
 *
 * Revision 2.6  2000/12/20  09:09:48  phil
 * ajout set_statiques
 *
 * Revision 2.5  2000/12/18  17:43:06  phil
 * ajout sortie pour le rayon
 *
 * Revision 2.4  2000/12/18  16:38:39  phil
 * ajout convergence vers une masse donneee
 *
 * Revision 2.3  2000/12/14  10:45:38  phil
 * ATTENTION : PASSAGE DE PHI A PSI
 *
 * Revision 2.2  2000/12/04  14:29:17  phil
 * changement nom omega pour eviter interference avec membre prive
 *
 * Revision 2.1  2000/11/17  10:07:14  phil
 * corrections diverses
 *
 * Revision 2.0  2000/11/17  10:04:08  phil
 * *** empty log message ***
 *
 *
 * $Header: /cvsroot/Lorene/C++/Source/Bhole_binaire/bhole_coal.C,v 1.8 2014/10/13 08:52:40 j_novak Exp $
 *
 */

//standard
#include <cmath>
#include <cstdlib>

// Lorene
#include "tenseur.h"
#include "bhole.h"

namespace Lorene {
void Bhole_binaire::set_statiques (double precis, double relax) {
    
    int nz_un = hole1.mp.get_mg()->get_nzone() ;
    int nz_deux = hole2.mp.get_mg()->get_nzone() ;
    
    set_omega(0) ;
    init_bhole_binaire() ;
    
    int indic = 1 ;
    int conte = 0 ;
 
    cout << "TROUS STATIQUES : " << endl ;
    while (indic == 1) {
	Cmp lapse_un_old (hole1.get_n_auto()()) ;
	Cmp lapse_deux_old (hole2.get_n_auto()()) ;
	
	solve_psi (precis, relax) ;
	solve_lapse (precis, relax) ;
	
	double erreur = 0 ;
	Tbl diff_un (diffrelmax (lapse_un_old, hole1.get_n_auto()())) ;
	for (int i=1 ; i<nz_un ; i++)
	    if (diff_un(i) > erreur)
		erreur = diff_un(i) ;
	
	Tbl diff_deux (diffrelmax (lapse_deux_old, hole2.get_n_auto()())) ;
	for (int i=1 ; i<nz_deux ; i++)
	    if (diff_deux(i) > erreur)
		erreur = diff_deux(i) ;
	
		
	cout << "PAS TOTAL : " << conte << " DIFFERENCE : " << erreur << endl ;
	
	if (erreur < precis)
	    indic = -1 ;
	conte ++ ;
    }
}

void Bhole_binaire::coal (double precis, double relax, int nbre_ome, double seuil_search, double m1, double m2, const int sortie) {

    assert (omega == 0) ;
    int nz1 = hole1.mp.get_mg()->get_nzone() ;
    int nz2 = hole1.mp.get_mg()->get_nzone() ;
    
    // Distance initiale
    double distance = hole1.mp.get_ori_x()-hole2.mp.get_ori_x() ;
    set_pos_axe (distance*(hole2.rayon/(hole1.rayon+hole2.rayon))) ;
    double scale_linear = (hole1.rayon + hole2.rayon)/2.*distance ;
    
    // Omega initial :
    double angulaire = sqrt((hole1.rayon+hole2.rayon)/distance/distance/distance) ;
    
    int indic = 1 ;
    int conte = 0 ;
    
    char name_iteration[40] ;
    char name_correction[40] ;
    char name_viriel[40] ;
    char name_ome [40] ;
    char name_linear[40] ;
    char name_axe[40] ;
    char name_error_m1[40] ;
    char name_error_m2[40] ;
    char name_r1[40] ;
    char name_r2[40] ;
    
    sprintf(name_iteration, "ite.dat") ;
    sprintf(name_correction, "cor.dat") ;
    sprintf(name_viriel, "vir.dat") ;
    sprintf(name_ome, "ome.dat") ;
    sprintf(name_linear, "linear.dat") ;
    sprintf(name_axe, "axe.dat") ;
    sprintf(name_error_m1, "error_m1.dat") ;
    sprintf(name_error_m2, "error_m2.dat") ; 
    sprintf(name_r1, "r1.dat") ;
    sprintf(name_r2, "r2.dat") ; 
    
    ofstream fiche_iteration(name_iteration) ;
    fiche_iteration.precision(8) ; 

    ofstream fiche_correction(name_correction) ;
    fiche_correction.precision(8) ; 
    
    ofstream fiche_viriel(name_viriel) ;
    fiche_viriel.precision(8) ; 
    
    ofstream fiche_ome(name_ome) ;
    fiche_ome.precision(8) ; 
   
    ofstream fiche_linear(name_linear) ;
    fiche_linear.precision(8) ; 
     
    ofstream fiche_axe(name_axe) ;
    fiche_axe.precision(8) ; 
    
    ofstream fiche_error_m1 (name_error_m1) ;
    fiche_error_m1.precision(8) ;
      
    ofstream fiche_error_m2 (name_error_m2) ;
    fiche_error_m2.precision(8) ;
  
    ofstream fiche_r1 (name_r1) ;
    fiche_r1.precision(8) ;
      
    ofstream fiche_r2 (name_r2) ;
    fiche_r2.precision(8) ;
    
    // LA BOUCLE EN AUGMENTANT OMEGA A LA MAIN PROGRESSIVEMENT : 
    cout << "OMEGA AUGMENTE A LA MAIN." << endl ;
    double homme = 0 ;
    for (int pas = 0 ; pas <nbre_ome ; pas ++) {
	
	homme += angulaire/nbre_ome ;
	set_omega (homme) ;
	
	Cmp shift_un_old (hole1.get_shift_auto()(0)) ;
	Cmp shift_deux_old (hole2.get_shift_auto()(0)) ;
	
	solve_shift (precis, relax) ;
	fait_tkij() ;
	
	solve_psi (precis, relax) ;
	solve_lapse (precis, relax) ;
	
	double erreur = 0 ;
	Tbl diff_un (diffrelmax (shift_un_old, hole1.get_shift_auto()(0))) ;
	for (int i=1 ; i<nz1 ; i++)
	    if (diff_un(i) > erreur)
		erreur = diff_un(i) ;
	
	Tbl diff_deux (diffrelmax (shift_deux_old, hole2.get_shift_auto()(0))) ;
	for (int i=1 ; i<nz2 ; i++)
	    if (diff_deux(i) > erreur)
		erreur = diff_deux(i) ;
	
	double error_viriel = viriel() ;
	double error_linear = linear_momentum_systeme_inf()(1)/scale_linear ;
	double error_m1 = 1.-sqrt(hole1.area()/16./M_PI)/m1 ;
	double error_m2 = 1.-sqrt(hole2.area()/16./M_PI)/m2 ;
	double r1 = hole1.mp.val_r(0, 1, 0, 0) ;
	double r2 = hole2.mp.val_r(0, 1, 0, 0) ;
		
	if (sortie != 0) {
	    fiche_iteration << conte << " " << erreur << endl ;
	    fiche_correction << conte << " " << hole1.get_regul() << " " << hole2.get_regul() << endl ;
	    fiche_viriel << conte << " " << error_viriel << endl ;
	    fiche_ome << conte << " " << homme << endl ;
	    fiche_linear << conte << " " << error_linear << endl ;
	    fiche_axe << conte << " " << pos_axe << endl ;
	    fiche_error_m1 << conte << " " << error_m1 << endl ;
	    fiche_error_m2 << conte << " " << error_m2 << endl ;
	    fiche_r1 << conte << " " << r1 << endl ;
	    fiche_r2 << conte << " " << r2 << endl ;
	    }
	    
	cout << "PAS TOTAL : " << conte << " DIFFERENCE : " << erreur << endl ;
	conte ++ ;
    }
    
    // BOUCLE AVEC BLOQUE :
    cout << "OMEGA VARIABLE" << endl ;
    indic = 1 ;     
    bool scale = false ;
    
    while (indic == 1) {
	
	Cmp shift_un_old (hole1.get_shift_auto()(0)) ;
	Cmp shift_deux_old (hole2.get_shift_auto()(0)) ;
	
	solve_shift (precis, relax) ;
	fait_tkij() ;
	
	solve_psi (precis, relax) ;
	solve_lapse (precis, relax) ;

	double erreur = 0 ;
	Tbl diff_un (diffrelmax (shift_un_old, hole1.get_shift_auto()(0))) ;
	for (int i=1 ; i<nz1 ; i++)
	    if (diff_un(i) > erreur)
		erreur = diff_un(i) ;
	
	Tbl diff_deux (diffrelmax (shift_deux_old, hole2.get_shift_auto()(0))) ;
	for (int i=1 ; i<nz2 ; i++)
	    if (diff_deux(i) > erreur)
		erreur = diff_deux(i) ;
	
        double error_viriel = viriel() ;
	double error_linear = linear_momentum_systeme_inf()(1)/scale_linear ;
	double error_m1 = 1.-sqrt(hole1.area()/16./M_PI)/m1 ;
	double error_m2 = 1.-sqrt(hole2.area()/16./M_PI)/m2 ;
	double r1 = hole1.mp.val_r(0, 1, 0, 0) ;
	double r2 = hole2.mp.val_r(0, 1, 0, 0) ;
	
	if (sortie != 0) {
	    fiche_iteration << conte << " " << erreur << endl ;
	    fiche_correction << conte << " " << hole1.regul << " " << hole2.regul << endl ;
	    fiche_viriel << conte << " " << error_viriel << endl ;
	    fiche_ome << conte << " " << omega << endl ;
	    fiche_linear << conte << " " << error_linear << endl ;
	    fiche_axe << conte << " " << pos_axe << endl ; 
	    fiche_error_m1 << conte << " " << error_m1 << endl ;
	    fiche_error_m2 << conte << " " << error_m2 << endl ; 
	    fiche_r1 << conte << " " << r1 << endl ;
	    fiche_r2 << conte << " " << r2 << endl ;
	    }
	    
	// On modifie omega, position de l'axe et les masses !
	if (erreur <= seuil_search)
	    scale = true ;
	if (scale) {
	    double scaling_ome = pow((2-error_viriel)/(2-2*error_viriel), 1.) ;
	    set_omega (omega*scaling_ome) ; 
	    
	    double scaling_axe = pow((2-error_linear)/(2-2*error_linear), 0.1) ;
	    set_pos_axe (pos_axe*scaling_axe) ;
	    
	    double scaling_r1 = pow((2-error_m1)/(2-2*error_m1), 0.1) ;
	    hole1.mp.homothetie_interne(scaling_r1) ;
	  
	    double scaling_r2 = pow((2-error_m2)/(2-2*error_m2), 0.1) ;
	    hole2.mp.homothetie_interne(scaling_r2) ;
	}
	
	cout << "PAS TOTAL : " << conte << " DIFFERENCE : " << erreur << endl ;
	if (erreur < precis)
	    indic = -1 ;
	conte ++ ;
    }
    
    fiche_iteration.close() ;
    fiche_correction.close() ;
    fiche_viriel.close() ;
    fiche_ome.close() ;
    fiche_linear.close() ;
    fiche_axe.close() ;
    fiche_error_m1.close() ;
    fiche_error_m2.close() ;
    fiche_r1.close() ;
    fiche_r2.close() ;
}
}