File: boson_star_equil.C

package info (click to toggle)
lorene 0.0.0~cvs20161116%2Bdfsg-2
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid, trixie
  • size: 26,472 kB
  • sloc: cpp: 212,946; fortran: 21,645; makefile: 1,750; sh: 4
file content (638 lines) | stat: -rw-r--r-- 20,315 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
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
/*
 * Method Boson_star::equilibrium
 *
 *    (see file boson_star.h for documentation).
 */

/*
 *   Copyright (c) 2012 Claire Some, Eric Gourgoulhon
 *
 *   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 boson_star_equil_C[] = "$Header: /cvsroot/Lorene/C++/Source/Compobj/boson_star_equil.C,v 1.6 2014/10/13 08:52:49 j_novak Exp $" ;

/*
 * $Id: boson_star_equil.C,v 1.6 2014/10/13 08:52:49 j_novak Exp $
 * $Log: boson_star_equil.C,v $
 * Revision 1.6  2014/10/13 08:52:49  j_novak
 * Lorene classes and functions now belong to the namespace Lorene.
 *
 * Revision 1.5  2014/10/06 15:13:04  j_novak
 * Modified #include directives to use c++ syntax.
 *
 * Revision 1.4  2013/04/03 12:10:13  e_gourgoulhon
 * Added member kk to Compobj; suppressed tkij
 *
 * Revision 1.3  2012/12/03 15:27:30  c_some
 * Small changes
 *
 * Revision 1.2  2012/11/23 15:43:05  c_some
 * Small changes
 *
 * Revision 1.1  2012/11/22 16:04:12  c_some
 * New class Boson_star
 *
 *
 * $Header: /cvsroot/Lorene/C++/Source/Compobj/boson_star_equil.C,v 1.6 2014/10/13 08:52:49 j_novak Exp $
 *
 */

// Headers C
#include <cmath>

// Headers Lorene
#include "boson_star.h"
#include "param.h"
#include "tenseur.h"

#include "graphique.h"
#include "utilitaires.h"
#include "unites.h"

namespace Lorene {
void Boson_star::equilibrium(double, double,
			     int nzadapt, const Tbl& phi_limit, const Itbl& icontrol,
			     const Tbl& control, Tbl& diff, Param*) {
			     
    // Fundamental constants and units
    // -------------------------------

    using namespace Unites ;	    
    
    // For the display 
    // ---------------
    char display_bold[]="x[1m" ; display_bold[0] = 27 ;
    char display_normal[] = "x[0m" ; display_normal[0] = 27 ;

    // Grid parameters
    // ---------------
    
    const Mg3d* mg = mp.get_mg() ; 
    int nz = mg->get_nzone() ;	    // total number of domains
    
    // The following is required to initialize mp_prev as a Map_et:
    Map_et& mp_et = dynamic_cast<Map_et&>(mp) ; 
        
    // Index of the point at phi=0, theta=pi/2 at the surface of the star:

    int nzet = nzadapt ; //## to be checked 
    
    assert(mg->get_type_t() == SYM) ; 
    int l_b = nzet - 1 ; 
    int j_b = mg->get_nt(l_b) - 1 ; 
    int k_b = 0 ; 
    
    // Value of the enthalpy defining the surface of the star
    // double ent_b = phi_limit(nzet-1) ;
    
    // Parameters to control the iteration
    // -----------------------------------
    
    int mer_max = icontrol(0) ; 
//    int mer_rot = icontrol(1) ;
//    int mer_change_omega = icontrol(2) ; 
//    int mer_fix_omega = icontrol(3) ; 
//##    int mer_mass = icontrol(4) ; 
    int mermax_poisson = icontrol(5) ; 
    int mer_triax = icontrol(6) ; 
//##    int delta_mer_kep = icontrol(7) ; 


    double precis = control(0) ; 
//    double omega_ini = control(1) ; 
    double relax = control(2) ;
    double relax_prev = double(1) - relax ;  
    double relax_poisson = control(3) ; 
//    double thres_adapt = control(4) ; 
    double ampli_triax = control(5) ; 
    double precis_adapt = control(6) ; 


    // Error indicators
    // ----------------
    
    diff.set_etat_qcq() ; 
    double& diff_phi = diff.set(0) ; 
    double& diff_nuf = diff.set(1) ; 
    double& diff_nuq = diff.set(2) ; 
//    double& diff_dzeta = diff.set(3) ; 
//    double& diff_ggg = diff.set(4) ; 
    double& diff_shift_x = diff.set(5) ; 
    double& diff_shift_y = diff.set(6) ; 
    double& vit_triax = diff.set(7) ; 
    
    // Parameters for the function Map_et::adapt
    // -----------------------------------------
    
    Param par_adapt ; 
    int nitermax = 100 ;  
    int niter ; 
    int adapt_flag = 1 ;    //  1 = performs the full computation, 
			    //  0 = performs only the rescaling by 
			    //      the factor alpha_r
    int nz_search = nzet + 1 ;  // Number of domains for searching the enthalpy
				//  isosurfaces
    double alpha_r ; 
    double reg_map = 1. ; // 1 = regular mapping, 0 = contracting mapping

    par_adapt.add_int(nitermax, 0) ; // maximum number of iterations to
				     // locate zeros by the secant method
    par_adapt.add_int(nzadapt, 1) ; // number of domains where the adjustment 
				    // to the isosurfaces of ent is to be 
				    // performed
    par_adapt.add_int(nz_search, 2) ;	// number of domains to search for
					// the enthalpy isosurface
    par_adapt.add_int(adapt_flag, 3) ; //  1 = performs the full computation, 
				       //  0 = performs only the rescaling by 
				       //      the factor alpha_r
    par_adapt.add_int(j_b, 4) ; //  theta index of the collocation point 
			        //  (theta_*, phi_*)
    par_adapt.add_int(k_b, 5) ; //  theta index of the collocation point 
			        //  (theta_*, phi_*)

    par_adapt.add_int_mod(niter, 0) ;  //  number of iterations actually used in 
				       //  the secant method
    
    par_adapt.add_double(precis_adapt, 0) ; // required absolute precision in 
					     // the determination of zeros by 
					     // the secant method
    par_adapt.add_double(reg_map, 1)	;  // 1. = regular mapping, 0 = contracting mapping
    
    par_adapt.add_double(alpha_r, 2) ;	// factor by which all the radial 
					   // distances will be multiplied 
    	   
    par_adapt.add_tbl(phi_limit, 0) ;	// array of values of the field Phi
				        // to define the isosurfaces. 			   

    // Parameters for the function Map_et::poisson for nuf
    // ----------------------------------------------------

    double precis_poisson = 1.e-16 ;     

	// Preparations 
	// ------------
	
	// Cartesian components of the shift vector are required
	beta.change_triad( mp.get_bvect_cart() ) ; 


    Cmp cssjm1_nuf(ssjm1_nuf) ; 
    Cmp cssjm1_nuq(ssjm1_nuq) ; 
    Cmp cssjm1_tggg(ssjm1_tggg) ; 
    Cmp cssjm1_khi(ssjm1_khi) ; 
    Tenseur cssjm1_wshift(mp, 1, CON, mp.get_bvect_cart() ) ;
    cssjm1_wshift.set_etat_qcq() ;
    for (int i=1; i<=3; i++) {
	cssjm1_wshift.set(i-1) = ssjm1_wshift(i) ; 
    }

    Tenseur cshift(mp, 1, CON, mp.get_bvect_cart() ) ;
    cshift.set_etat_qcq() ;
    for (int i=1; i<=3; i++) {
      cshift.set(i-1) = -beta(i) ; 
    }

    Tenseur cw_shift(mp, 1, CON, mp.get_bvect_cart() ) ;
    cw_shift.set_etat_qcq() ;
    for (int i=1; i<=3; i++) {
      cw_shift.set(i-1) = w_shift(i) ; 
    }

    Tenseur ckhi_shift(mp) ;
    ckhi_shift.set_etat_qcq() ;
    ckhi_shift.set() = khi_shift ; 

    Param par_poisson_nuf ; 
    par_poisson_nuf.add_int(mermax_poisson,  0) ;  // maximum number of iterations
    par_poisson_nuf.add_double(relax_poisson,  0) ; // relaxation parameter
    par_poisson_nuf.add_double(precis_poisson, 1) ; // required precision
    par_poisson_nuf.add_int_mod(niter, 0) ;  //  number of iterations actually used 
    par_poisson_nuf.add_cmp_mod( cssjm1_nuf ) ; 
					   
    Param par_poisson_nuq ; 
    par_poisson_nuq.add_int(mermax_poisson,  0) ;  // maximum number of iterations
    par_poisson_nuq.add_double(relax_poisson,  0) ; // relaxation parameter
    par_poisson_nuq.add_double(precis_poisson, 1) ; // required precision
    par_poisson_nuq.add_int_mod(niter, 0) ;  //  number of iterations actually used 
    par_poisson_nuq.add_cmp_mod( cssjm1_nuq ) ; 
					   
    Param par_poisson_tggg ; 
    par_poisson_tggg.add_int(mermax_poisson,  0) ;  // maximum number of iterations
    par_poisson_tggg.add_double(relax_poisson,  0) ; // relaxation parameter
    par_poisson_tggg.add_double(precis_poisson, 1) ; // required precision
    par_poisson_tggg.add_int_mod(niter, 0) ;  //  number of iterations actually used 
    par_poisson_tggg.add_cmp_mod( cssjm1_tggg ) ; 
    double lambda_tggg ;
    par_poisson_tggg.add_double_mod( lambda_tggg ) ; 
    
    Param par_poisson_dzeta ; 
    double lbda_grv2 ;
    par_poisson_dzeta.add_double_mod( lbda_grv2 ) ; 
 					   
    // Parameters for the function Scalar::poisson_vect
    // -------------------------------------------------

    Param par_poisson_vect ; 

    par_poisson_vect.add_int(mermax_poisson,  0) ;  // maximum number of iterations
    par_poisson_vect.add_double(relax_poisson,  0) ; // relaxation parameter
    par_poisson_vect.add_double(precis_poisson, 1) ; // required precision
    par_poisson_vect.add_cmp_mod( cssjm1_khi ) ; 
    par_poisson_vect.add_tenseur_mod( cssjm1_wshift ) ; 
    par_poisson_vect.add_int_mod(niter, 0) ;   

 					   
    // Initializations
    // ---------------

	//## Spherical components of the shift vector are restored
	beta.change_triad( mp.get_bvect_spher() ) ; 
	update_metric() ; 
	//## Back to Cartesian components 
	beta.change_triad( mp.get_bvect_cart() ) ; 
	

    // Quantities at the previous step : 	
    Map_et mp_prev = mp_et ; 
    Scalar rphi_prev = rphi ;	    
    Scalar logn_prev = logn ;	    
    Scalar dzeta_prev = dzeta ;	    
    
    // Creation of uninitialized tensors:
    Scalar source_nuf(mp) ;    // source term in the equation for nuf
    Scalar source_nuq(mp) ;    // source term in the equation for nuq
    Scalar source_dzf(mp) ;	// matter source term in the eq. for dzeta
    Scalar source_dzq(mp) ;	// quadratic source term in the eq. for dzeta
    Scalar source_tggg(mp) ;	// source term in the eq. for tggg
    Vector source_shift(mp, CON, mp.get_bvect_cart()) ;  
						    // source term for shift
    
		    
    ofstream fichconv("convergence.d") ;    // Output file for diff_phi
    fichconv << "#     diff_phi     GRV2       max_triax      vit_triax" << endl ; 
    
    
    ofstream fichevol("evolution.d") ;    // Output file for various quantities
    fichevol << 
    "#       |dH/dr_eq/dH/dr_pole|      r_pole/r_eq	rphi_c" 
    << endl ; 
    
    diff_phi = 1 ; 
    double err_grv2 = 1 ; 
    double max_triax_prev = 0 ;	    // Triaxial amplitude at previous step

    //=========================================================================
    // 			Start of iteration
    //=========================================================================

    for(int mer=0 ; (diff_phi > precis) && (mer<mer_max) ; mer++ ) {

 	cout << "-----------------------------------------------" << endl ;
	cout << "step: " << mer << endl ;
	cout << "diff_phi = " << display_bold << diff_phi << display_normal
	     << endl ;    
	cout << "err_grv2 = " << err_grv2 << endl ;    
	fichconv << mer ;
	fichevol << mer ;


	//-----------------------------------------------
	//  Sources of the Poisson equations
	//-----------------------------------------------
	
	// Source for nu
	// -------------
	Scalar bet = log(bbb) ; 
	bet.std_spectral_base() ; 

	Vector d_logn = logn.derive_cov( mp.flat_met_spher() ) ; 
	Vector d_bet = bet.derive_cov( mp.flat_met_spher() ) ; 

	Scalar s_euler = stress_euler.trace(gamma) ; 
	source_nuf =  qpig * a_car *( ener_euler + s_euler ) ; 

	source_nuq = ak_car - d_logn(1)*(d_logn(1)+d_bet(1))
		- d_logn(2)*(d_logn(2)+d_bet(2))
		- d_logn(3)*(d_logn(3)+d_bet(3)) ; 

	source_nuf.std_spectral_base() ; 	
	source_nuq.std_spectral_base() ; 	

	// Source for dzeta
	// ----------------
	source_dzf = 2 * qpig * a_car * b_car * stress_euler(3,3) ;
	source_dzf.std_spectral_base() ; 
  
	source_dzq = 1.5 * ak_car 
	    - d_logn(1)*d_logn(1) - d_logn(2)*d_logn(2) - d_logn(3)*d_logn(3) ;	    
	source_dzq.std_spectral_base() ; 	
	
	// Source for tggg
	// ---------------
	
	source_tggg = 2 * qpig * nn * a_car * bbb * (s_euler - b_car * stress_euler(3,3)) ;
	source_tggg.std_spectral_base() ; 
	
	source_tggg.mult_rsint() ; 
	

	// Source for shift
	// ----------------
	    
	// Matter term: 
	Vector mom_euler_cart = mom_euler ;
	mom_euler_cart.change_triad(mp.get_bvect_cart()) ;
	source_shift = (-4*qpig) * nn * a_car  * mom_euler_cart ;

	// Quadratic terms:
	Vector vtmp =  6 * bet.derive_con( mp.flat_met_spher() ) 
	    - 2 * logn.derive_con( mp.flat_met_spher() ) ;

	Vector squad  = nn * contract(kk, 1, vtmp, 0) / b_car ; 
    squad.change_triad(mp.get_bvect_cart()) ;
    
	source_shift = source_shift + squad.up(0, mp.flat_met_cart() ) ; 

	//----------------------------------------------
	// Resolution of the Poisson equation for nuf 
	//----------------------------------------------

	source_nuf.poisson(par_poisson_nuf, nuf) ; 
	
	cout << "Test of the Poisson equation for nuf :" << endl ; 
	Tbl err = source_nuf.test_poisson(nuf, cout, true) ; 
	diff_nuf = err(0, 0) ; 

	//---------------------------------------
	// Triaxial perturbation of nuf
	//---------------------------------------
	
	if (mer == mer_triax) {
	
	    if ( mg->get_np(0) == 1 ) {
		cout << 
		"Boson_star::equilibrium: np must be stricly greater than 1"
		<< endl << " to set a triaxial perturbation !" << endl ; 
		abort() ; 
	    }
	
	    const Coord& phi = mp.phi ; 
	    const Coord& sint = mp.sint ; 
	    Scalar perturb(mp) ; 
	    perturb = 1 + ampli_triax * sint*sint * cos(2*phi) ;
	    nuf = nuf * perturb ;  
	    
	    nuf.std_spectral_base() ;    // set the bases for spectral expansions
				    //  to be the standard ones for a 
				    //  scalar field

	}
	
	// Monitoring of the triaxial perturbation
	// ---------------------------------------
	
	const Valeur& va_nuf = nuf.get_spectral_va() ; 
	va_nuf.coef() ;		// Computes the spectral coefficients
	double max_triax = 0 ; 

	if ( mg->get_np(0) > 1 ) {

	    for (int l=0; l<nz; l++) {	// loop on the domains
		for (int j=0; j<mg->get_nt(l); j++) {
		    for (int i=0; i<mg->get_nr(l); i++) {
		
			// Coefficient of cos(2 phi) : 
			double xcos2p = (*(va_nuf.c_cf))(l, 2, j, i) ; 

			// Coefficient of sin(2 phi) : 
			double xsin2p = (*(va_nuf.c_cf))(l, 3, j, i) ; 
		    
			double xx = sqrt( xcos2p*xcos2p + xsin2p*xsin2p ) ; 
		    
			max_triax = ( xx > max_triax ) ? xx : max_triax ; 		    
		    }
		} 
	    }
	    
	}
	
	cout << "Triaxial part of nuf : " << max_triax << endl ; 

	    //----------------------------------------------
	    // Resolution of the Poisson equation for nuq 
	    //----------------------------------------------

	    source_nuq.poisson(par_poisson_nuq, nuq) ; 
	    
	    cout << "Test of the Poisson equation for nuq :" << endl ; 
	    err = source_nuq.test_poisson(nuq, cout, true) ;
	    diff_nuq = err(0, 0) ; 
	
	    //---------------------------------------------------------
	    // Resolution of the vector Poisson equation for the shift
	    //---------------------------------------------------------


	    for (int i=1; i<=3; i++) {
	      if(source_shift(i).get_etat() != ETATZERO) {
		    if(source_shift(i).dz_nonzero()) {
			assert( source_shift(i).get_dzpuis() == 4 ) ; 
		    }
		    else{
			(source_shift.set(i)).set_dzpuis(4) ; 
		    }
		}
	    }

	    double lambda_shift = double(1) / double(3) ; 

	    if ( mg->get_np(0) == 1 ) {
		lambda_shift = 0 ; 
	    }
	
	    Tenseur csource_shift(mp, 1, CON, mp.get_bvect_cart() ) ;
	    csource_shift.set_etat_qcq() ;
	    for (int i=1; i<=3; i++) {
		csource_shift.set(i-1) = source_shift(i) ; 
	    }
	    csource_shift.set(2).set_etat_zero() ;  //## bizarre...

	    csource_shift.poisson_vect(lambda_shift, par_poisson_vect, 
	    			      cshift, cw_shift, ckhi_shift) ;      	    

	    for (int i=1; i<=3; i++) {
		beta.set(i) = - cshift(i-1) ;
		beta.set(i).set_dzpuis(0) ;     //## bizarre...
		w_shift.set(i) = cw_shift(i-1) ; 
	    }
	    khi_shift = ckhi_shift() ; 

	    cout << "Test of the Poisson equation for shift_x :" << endl ; 
	    err = source_shift(1).test_poisson(-beta(1), cout, true) ;
	    diff_shift_x = err(0, 0) ; 
	
	    cout << "Test of the Poisson equation for shift_y :" << endl ; 
	    err = source_shift(2).test_poisson(-beta(2), cout, true) ;
	    diff_shift_y = err(0, 0) ; 
	    
	    // Computation of tnphi and nphi from the Cartesian components
	    //  of the shift
	    // -----------------------------------------------------------
	    
	    fait_nphi() ; 
	


	//----------------------------------------------------
	// Adaptation of the mapping to the new enthalpy field
	//----------------------------------------------------
    
	// Shall the adaptation be performed ?
	// ---------------------------------
	
	adapt_flag = 0 ;	// No adaptation of the mapping 

	mp_prev = mp_et ; 

	
	//---------------------------------------------------------
	// Matter source terms in the gravitational field equations	
	//---------------------------------------------------------

	//## Computation of tnphi and nphi from the Cartesian components
	//  of the shift for the test in hydro_euler():
	    
    fait_nphi() ; 

	// Update of the energy-momentum tensor
	update_ener_mom() ; 
	

	    //-------------------------------------------------------
	    //	2-D Poisson equation for tggg
	    //-------------------------------------------------------

	    Cmp csource_tggg(source_tggg) ; 
	    Cmp ctggg(tggg) ; 
	    mp.poisson2d(csource_tggg, mp.cmp_zero(), par_poisson_tggg,
			 ctggg) ; 
	    tggg = ctggg ; 

	    
	    //-------------------------------------------------------
	    //	2-D Poisson equation for dzeta
	    //-------------------------------------------------------

	    Cmp csource_dzf(source_dzf) ; 
	    Cmp csource_dzq(source_dzq) ; 
	    Cmp cdzeta(dzeta) ; 
	    mp.poisson2d(csource_dzf, csource_dzq, par_poisson_dzeta,
			 cdzeta) ;
	    dzeta = cdzeta ; 
	    
	    err_grv2 = lbda_grv2 - 1; 
	    cout << "GRV2: " << err_grv2 << endl ; 
	    

	//---------------------------------------
	// Computation of the metric coefficients (except for N^phi)
	//---------------------------------------

	// Relaxations on nu and dzeta :  

	if (mer >= 10) {
		logn = relax * logn + relax_prev * logn_prev ;

		dzeta = relax * dzeta + relax_prev * dzeta_prev ; 
	}

	// Update of the metric coefficients N, A, B and computation of K_ij :

	//## Spherical components of the shift vector are restored
	beta.change_triad( mp.get_bvect_spher() ) ; 
	update_metric() ; 
	//## Back to Cartesian components 
	beta.change_triad( mp.get_bvect_cart() ) ; 
	
	//-----------------------
	//  Informations display
	//-----------------------

	cout << *this << endl ; 


	//------------------------------------------------------------
	//  Relative change in Phi with respect to previous step 
	//------------------------------------------------------------

	Tbl diff_phi_tbl = diffrel( rphi, rphi_prev ) ; 
	diff_phi = diff_phi_tbl(0) ; 
	for (int l=1; l<nzet; l++) {
	    diff_phi += diff_phi_tbl(l) ; 
	}
	diff_phi /= nzet ; 
	
	fichconv << "  " << log10( fabs(diff_phi) + 1.e-16 ) ;
	fichconv << "  " << log10( fabs(err_grv2) + 1.e-16 ) ;
	fichconv << "  " << log10( fabs(max_triax) + 1.e-16 ) ;

	vit_triax = 0 ;
	if ( (mer > mer_triax+1) && (max_triax_prev > 1e-13) ) {
		vit_triax = (max_triax - max_triax_prev) / max_triax_prev ; 
	}

	fichconv << "  " << vit_triax ;
	
	//------------------------------
	//  Recycling for the next step
	//------------------------------
	
	rphi_prev = rphi ; 
	logn_prev = logn ; 
	dzeta_prev = dzeta ; 
	max_triax_prev = max_triax ; 
	
	fichconv << endl ;
	fichevol << endl ;
	fichconv.flush() ; 
	fichevol.flush() ; 
	
    } // End of main loop
    
    //=========================================================================
    // 			End of iteration
    //=========================================================================

    ssjm1_nuf = cssjm1_nuf ; 
    ssjm1_nuq = cssjm1_nuq ; 
    ssjm1_tggg = cssjm1_tggg ; 
    ssjm1_khi = cssjm1_khi ; 
    for (int i=1; i<=3; i++) {
		ssjm1_wshift.set(i) = cssjm1_wshift(i-1) ; 
    }

	// Spherical components of the shift vector are restored
	beta.change_triad( mp.get_bvect_spher() ) ; 

    fichconv.close() ; 
    fichevol.close() ; 
    
}
}