LORENE
et_bin_equil_regu.C
1 /*
2  * Method of class Etoile_bin to compute an equilibrium configuration
3  * by regularizing source.
4  *
5  * (see file etoile.h for documentation).
6  *
7  */
8 
9 /*
10  * Copyright (c) 2000-2001 Eric Gourgoulhon
11  * Copyright (c) 2000-2001 Keisuke Taniguchi
12  *
13  * This file is part of LORENE.
14  *
15  * LORENE is free software; you can redistribute it and/or modify
16  * it under the terms of the GNU General Public License as published by
17  * the Free Software Foundation; either version 2 of the License, or
18  * (at your option) any later version.
19  *
20  * LORENE is distributed in the hope that it will be useful,
21  * but WITHOUT ANY WARRANTY; without even the implied warranty of
22  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
23  * GNU General Public License for more details.
24  *
25  * You should have received a copy of the GNU General Public License
26  * along with LORENE; if not, write to the Free Software
27  * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
28  *
29  */
30 
31 
32 
33 
34 /*
35  * $Id: et_bin_equil_regu.C,v 1.9 2016/12/05 16:17:52 j_novak Exp $
36  * $Log: et_bin_equil_regu.C,v $
37  * Revision 1.9 2016/12/05 16:17:52 j_novak
38  * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
39  *
40  * Revision 1.8 2014/10/13 08:52:55 j_novak
41  * Lorene classes and functions now belong to the namespace Lorene.
42  *
43  * Revision 1.7 2014/10/06 15:13:08 j_novak
44  * Modified #include directives to use c++ syntax.
45  *
46  * Revision 1.6 2009/06/15 09:26:57 k_taniguchi
47  * Improved the rescaling of the domains.
48  *
49  * Revision 1.5 2004/03/25 10:29:03 j_novak
50  * All LORENE's units are now defined in the namespace Unites (in file unites.h).
51  *
52  * Revision 1.4 2003/09/01 06:48:08 k_taniguchi
53  * Change of the domain which should be resized.
54  *
55  * Revision 1.3 2003/08/31 05:35:38 k_taniguchi
56  * Addition of the specification of the domain
57  * which is resized.
58  *
59  * Revision 1.2 2002/12/11 12:51:26 k_taniguchi
60  * Change the multiplication "*" to "%"
61  * and flat_scalar_prod to flat_scalar_prod_desal.
62  *
63  * Revision 1.1.1.1 2001/11/20 15:19:28 e_gourgoulhon
64  * LORENE
65  *
66  * Revision 2.17 2001/08/07 09:49:00 keisuke
67  * Change of the method to set the longest radius of a star
68  * on the first domain.
69  * Addition of a new argument in Etoile_bin::equil_regular : Tbl fact.
70  *
71  * Revision 2.16 2001/06/22 08:54:53 keisuke
72  * Set the inner values of the second domain of ent
73  * by using the outer ones of the first domain.
74  *
75  * Revision 2.15 2001/05/17 12:22:26 keisuke
76  * Change of the method to calculate chi from setting position in map
77  * to val_point.
78  *
79  * Revision 2.14 2001/02/07 09:47:28 eric
80  * unsgam1 est desormais donne par Eos::der_nbar_ent (cas newtonien)
81  * ou Eos::der_ener_ent (cas relativiste).
82  *
83  * Revision 2.13 2001/01/16 17:02:32 keisuke
84  * *** empty log message ***
85  *
86  * Revision 2.12 2001/01/16 16:58:08 keisuke
87  * Change the method to set the values on the surface.
88  *
89  * Revision 2.11 2001/01/10 16:45:34 keisuke
90  * Set the inner values of the second domain of logn_auto
91  * by using the outer ones of the first domain.
92  *
93  * Revision 2.10 2000/12/20 10:33:14 eric
94  * Changement important : nz_search = nzet ---> nz_search = nzet + 1
95  *
96  * Revision 2.9 2000/10/25 14:01:03 keisuke
97  * Modif de Map_et::adapt: on y rentre desormais avec nz_search
98  * (dans le cas present nz_search = nzet).
99  *
100  * Revision 2.8 2000/10/06 15:29:01 keisuke
101  * Change poisson_vect into poisson_vect_regu.
102  *
103  * Revision 2.7 2000/09/25 15:01:10 keisuke
104  * Suppress "int" from the declaration of k_div.
105  *
106  * Revision 2.6 2000/09/22 15:51:39 keisuke
107  * d_logn_auto est desormais calcule en dehors (dans update_metric).
108  *
109  * Revision 2.5 2000/09/13 09:50:33 keisuke
110  * Minor change on change_triad.
111  *
112  * Revision 2.4 2000/09/08 15:57:31 keisuke
113  * Change the basis of d_logn_auto_div from the spherical coordinate
114  * to the Cartesian one with respect to ref_triad.
115  *
116  * Revision 2.3 2000/09/07 15:47:19 keisuke
117  * Minor change.
118  *
119  * Revision 2.2 2000/09/07 15:43:41 keisuke
120  * Add a new argument in poisson_regular and suppress logn_auto_total.
121  *
122  * Revision 2.1 2000/08/29 14:01:43 keisuke
123  * Modify the arguments of poisson_regular.
124  *
125  * Revision 2.0 2000/08/29 11:39:02 eric
126  * Version provisoire.
127  *
128  *
129  * $Header: /cvsroot/Lorene/C++/Source/Etoile/et_bin_equil_regu.C,v 1.9 2016/12/05 16:17:52 j_novak Exp $
130  *
131  */
132 
133 // Headers C
134 #include <cmath>
135 
136 // Headers Lorene
137 #include "etoile.h"
138 #include "param.h"
139 #include "eos.h"
140 #include "utilitaires.h"
141 #include "unites.h"
142 
143 namespace Lorene {
144 
145 //********************************************************************
146 
147 void Etoile_bin::equil_regular(double ent_c, int mermax, int mermax_poisson,
148  double relax_poisson, int mermax_potvit,
149  double relax_potvit, double thres_adapt,
150  const Tbl& fact_resize, Tbl& diff) {
151 
152  // Fundamental constants and units
153  // -------------------------------
154  using namespace Unites ;
155 
156  // Initializations
157  // ---------------
158 
159  k_div = 2 ; // Regularity parameter for poisson_regular
160 
161  const Mg3d* mg = mp.get_mg() ;
162  int nz = mg->get_nzone() ; // total number of domains
163 
164  // The following is required to initialize mp_prev as a Map_et:
165  Map_et& mp_et = dynamic_cast<Map_et&>(mp) ;
166 
167  // Domain and radial indices of points at the surface of the star:
168  int l_b = nzet - 1 ;
169  int i_b = mg->get_nr(l_b) - 1 ;
170  int k_b ;
171  int j_b ;
172 
173  // Value of the enthalpy defining the surface of the star
174  double ent_b = 0 ;
175 
176  // Error indicators
177  // ----------------
178 
179  double& diff_ent = diff.set(0) ;
180  double& diff_vel_pot = diff.set(1) ;
181  double& diff_logn = diff.set(2) ;
182  double& diff_beta = diff.set(3) ;
183  double& diff_shift_x = diff.set(4) ;
184  double& diff_shift_y = diff.set(5) ;
185  double& diff_shift_z = diff.set(6) ;
186 
187  // Parameters for the function Map_et::adapt
188  // -----------------------------------------
189 
190  Param par_adapt ;
191  int nitermax = 100 ;
192  int niter ;
193  int adapt_flag = 1 ; // 1 = performs the full computation,
194  // 0 = performs only the rescaling by
195  // the factor alpha_r
196  //## int nz_search = nzet + 1 ; // Number of domains for searching the enthalpy
197  int nz_search = nzet ; // Number of domains for searching the enthalpy
198  // isosurfaces
199 
200  double precis_secant = 1.e-14 ;
201  double alpha_r ;
202  double reg_map = 1. ; // 1 = regular mapping, 0 = contracting mapping
203 
204  Tbl ent_limit(nz) ;
205 
206  par_adapt.add_int(nitermax, 0) ; // maximum number of iterations to
207  // locate zeros by the secant method
208  par_adapt.add_int(nzet, 1) ; // number of domains where the adjustment
209  // to the isosurfaces of ent is to be
210  // performed
211  par_adapt.add_int(nz_search, 2) ; // number of domains to search for
212  // the enthalpy isosurface
213  par_adapt.add_int(adapt_flag, 3) ; // 1 = performs the full computation,
214  // 0 = performs only the rescaling by
215  // the factor alpha_r
216  par_adapt.add_int(j_b, 4) ; // theta index of the collocation point
217  // (theta_*, phi_*)
218  par_adapt.add_int(k_b, 5) ; // theta index of the collocation point
219  // (theta_*, phi_*)
220  par_adapt.add_int_mod(niter, 0) ; // number of iterations actually
221  // used in the secant method
222  par_adapt.add_double(precis_secant, 0) ; // required absolute precision in
223  // the determination of zeros by
224  // the secant method
225  par_adapt.add_double(reg_map, 1) ; // 1. = regular mapping,
226  // 0 = contracting mapping
227  par_adapt.add_double(alpha_r, 2) ; // factor by which all the radial
228  // distances will be multiplied
229  par_adapt.add_tbl(ent_limit, 0) ; // array of values of the field ent
230  // to define the isosurfaces.
231 
232  // Parameters for the function Map_et::poisson_regular for logn_auto
233  // -----------------------------------------------------------------
234 
235  double precis_poisson = 1.e-16 ;
236 
237  Param par_poisson1 ;
238 
239  par_poisson1.add_int(mermax_poisson, 0) ; // maximum number of iterations
240  par_poisson1.add_double(relax_poisson, 0) ; // relaxation parameter
241  par_poisson1.add_double(precis_poisson, 1) ; // required precision
242  par_poisson1.add_int_mod(niter, 0) ; // number of iterations actually
243  // used
244  par_poisson1.add_cmp_mod( ssjm1_logn ) ;
245 
246  // Parameters for the function Map_et::poisson for beta_auto
247  // ---------------------------------------------------------
248 
249  Param par_poisson2 ;
250 
251  par_poisson2.add_int(mermax_poisson, 0) ; // maximum number of iterations
252  par_poisson2.add_double(relax_poisson, 0) ; // relaxation parameter
253  par_poisson2.add_double(precis_poisson, 1) ; // required precision
254  par_poisson2.add_int_mod(niter, 0) ; // number of iterations actually
255  // used
256  par_poisson2.add_cmp_mod( ssjm1_beta ) ;
257 
258  // Parameters for the function Tenseur::poisson_vect_regu
259  // ------------------------------------------------------
260 
261  Param par_poisson_vect ;
262 
263  par_poisson_vect.add_int(mermax_poisson, 0) ; // maximum number of
264  // iterations
265  par_poisson_vect.add_double(relax_poisson, 0) ; // relaxation parameter
266  par_poisson_vect.add_double(precis_poisson, 1) ; // required precision
267  par_poisson_vect.add_cmp_mod( ssjm1_khi ) ;
268  par_poisson_vect.add_tenseur_mod( ssjm1_wshift ) ;
269  par_poisson_vect.add_int_mod(niter, 0) ;
270 
271 
272  // External potential
273  // ------------------
274 
275  Tenseur pot_ext = logn_comp + pot_centri + loggam ;
276 //##
277 // des_coupe_z(pot_ext(), 0., 1, "pot_ext", &(ent()) ) ;
278 //##
279 
280  Tenseur ent_jm1 = ent ; // Enthalpy at previous step
281 
282  Tenseur source(mp) ; // source term in the equation for logn_auto
283  // and beta_auto
284 
285  Tenseur source_shift(mp, 1, CON, ref_triad) ; // source term in the
286  // equation for shift_auto
287 
288  Cmp source_regu(mp) ;
289  Cmp source_div(mp) ;
290 
291  //=========================================================================
292  // Start of iteration
293  //=========================================================================
294 
295  for(int mer=0 ; mer<mermax ; mer++ ) {
296 
297  cout << "-----------------------------------------------" << endl ;
298  cout << "step: " << mer << endl ;
299  cout << "diff_ent = " << diff_ent << endl ;
300 
301  //-----------------------------------------------------
302  // Resolution of the elliptic equation for the velocity
303  // scalar potential
304  //-----------------------------------------------------
305 
306  if (irrotational) {
307 
308  diff_vel_pot = velocity_potential(mermax_potvit, precis_poisson,
309  relax_potvit) ;
310 
311  }
312 
313  //-----------------------------------------------------
314  // Computation of the new radial scale
315  //-----------------------------------------------------
316 
317  // alpha_r (r = alpha_r r') is determined so that the enthalpy
318  // takes the requested value ent_b at the stellar surface
319 
320  // Values at the center of the star:
321  double logn_auto_c = logn_auto()(0, 0, 0, 0) ;
322  double pot_ext_c = pot_ext()(0, 0, 0, 0) ;
323 
324  // Search for the reference point (theta_*, phi_*) [notation of
325  // Bonazzola, Gourgoulhon & Marck PRD 58, 104020 (1998)]
326  // at the surface of the star
327  // ------------------------------------------------------------
328  double alpha_r2 = 0 ;
329  for (int k=0; k<mg->get_np(l_b); k++) {
330  for (int j=0; j<mg->get_nt(l_b); j++) {
331 
332  double pot_ext_b = pot_ext()(l_b, k, j, i_b) ;
333  double logn_auto_b = logn_auto()(l_b, k, j, i_b) ;
334 
335  double alpha_r2_jk = ( ent_c - ent_b + pot_ext_c - pot_ext_b) /
336  ( logn_auto_b - logn_auto_c ) ;
337 
338 // cout << "k, j, alpha_r2_jk : " << k << " " << j << " "
339 // << alpha_r2_jk << endl ;
340 
341  if (alpha_r2_jk > alpha_r2) {
342  alpha_r2 = alpha_r2_jk ;
343  k_b = k ;
344  j_b = j ;
345  }
346 
347  }
348  }
349 
350  alpha_r = sqrt(alpha_r2) ;
351 
352  cout << "k_b, j_b, alpha_r: " << k_b << " " << j_b << " "
353  << alpha_r << endl ;
354 
355  // New value of logn_auto
356  // ----------------------
357 
358  logn_auto = alpha_r2 * logn_auto ;
359  logn_auto_regu = alpha_r2 * logn_auto_regu ;
360  logn_auto_c = logn_auto()(0, 0, 0, 0) ;
361 
362 
363  //------------------------------------------------------------
364  // Change the values of the inner points of the second domain
365  // by those of the outer points of the first domain
366  //------------------------------------------------------------
367 
368  (logn_auto().va).smooth(nzet, (logn_auto.set()).va) ;
369 
370 
371  //--------------------
372  // First integral --> enthalpy in all space
373  //--------------------
374 
375  ent = (ent_c + logn_auto_c + pot_ext_c) - logn_auto - pot_ext ;
376 
377  (ent().va).smooth(nzet, (ent.set()).va) ;
378 
379  //----------------------------------------------------
380  // Adaptation of the mapping to the new enthalpy field
381  //----------------------------------------------------
382 
383  // Shall the adaptation be performed (cusp) ?
384  // ------------------------------------------
385 
386  double dent_eq = ent().dsdr().val_point(ray_eq(),M_PI/2.,0.) ;
387  double dent_pole = ent().dsdr().val_point(ray_pole(),0.,0.) ;
388  double rap_dent = fabs( dent_eq / dent_pole ) ;
389  cout << "| dH/dr_eq / dH/dr_pole | = " << rap_dent << endl ;
390 
391  if ( rap_dent < thres_adapt ) {
392  adapt_flag = 0 ; // No adaptation of the mapping
393  cout << "******* FROZEN MAPPING *********" << endl ;
394  }
395  else{
396  adapt_flag = 1 ; // The adaptation of the mapping is to be
397  // performed
398  }
399 
400 
401  ent_limit.set_etat_qcq() ;
402  for (int l=0; l<nzet; l++) { // loop on domains inside the star
403  ent_limit.set(l) = ent()(l, k_b, j_b, i_b) ;
404  }
405  ent_limit.set(nzet-1) = ent_b ;
406 
407  Map_et mp_prev = mp_et ;
408 
409 //##
410 // des_coupe_z(ent(), 0., 1, "ent before adapt", &(ent()) ) ;
411 //##
412 
413  mp.adapt(ent(), par_adapt) ;
414 
415  // Readjustment of the external boundary of domain l=nzet
416  // to keep a fixed ratio with respect to star's surface
417 
418  double rr_in_1 = mp.val_r(nzet, -1., M_PI/2., 0.) ;
419 
420  // Resizes the outer boundary of the shell including the comp. NS
421  double rr_out_nm2 = mp.val_r(nz-2, 1., M_PI/2., 0.) ;
422 
423  mp.resize(nz-2, rr_in_1/rr_out_nm2 * fact_resize(1)) ;
424 
425  // Resizes the inner boundary of the shell including the comp. NS
426  double rr_out_nm3 = mp.val_r(nz-3, 1., M_PI/2., 0.) ;
427 
428  mp.resize(nz-3, rr_in_1/rr_out_nm3 * fact_resize(0)) ;
429 
430  if (nz > nzet+3) {
431 
432  // Resize of the domain from 1(nzet) to N-4
433  double rr_out_nm3_new = mp.val_r(nz-3, 1., M_PI/2., 0.) ;
434 
435  for (int i=nzet-1; i<nz-4; i++) {
436 
437  double rr_out_i = mp.val_r(i, 1., M_PI/2., 0.) ;
438 
439  double rr_mid = rr_out_i
440  + (rr_out_nm3_new - rr_out_i) / double(nz - 3 - i) ;
441 
442  double rr_2timesi = 2. * rr_out_i ;
443 
444  if (rr_2timesi < rr_mid) {
445 
446  double rr_out_ip1 = mp.val_r(i+1, 1., M_PI/2., 0.) ;
447 
448  mp.resize(i+1, rr_2timesi / rr_out_ip1) ;
449 
450  }
451  else {
452 
453  double rr_out_ip1 = mp.val_r(i+1, 1., M_PI/2., 0.) ;
454 
455  mp.resize(i+1, rr_mid / rr_out_ip1) ;
456 
457  } // End of else
458 
459  } // End of i loop
460 
461  } // End of (nz > nzet+3) loop
462 
463 //##
464 // des_coupe_z(ent(), 0., 1, "ent after adapt", &(ent()) ) ;
465 //##
466  //----------------------------------------------------
467  // Computation of the enthalpy at the new grid points
468  //----------------------------------------------------
469 
470  mp_prev.homothetie(alpha_r) ;
471 
472  mp.reevaluate_symy(&mp_prev, nzet+1, ent.set()) ;
473 
474 // des_coupe_z(ent(), 0., 1, "ent after reevaluate", &(ent()) ) ;
475 
476  double ent_s_max = -1 ;
477  int k_s_max = -1 ;
478  int j_s_max = -1 ;
479  for (int k=0; k<mg->get_np(l_b); k++) {
480  for (int j=0; j<mg->get_nt(l_b); j++) {
481  double xx = fabs( ent()(l_b, k, j, i_b) ) ;
482  if (xx > ent_s_max) {
483  ent_s_max = xx ;
484  k_s_max = k ;
485  j_s_max = j ;
486  }
487  }
488  }
489  cout << "Max. abs(enthalpy) at the boundary between domains nzet-1"
490  << " and nzet : " << endl ;
491  cout << " " << ent_s_max << " reached for k = " << k_s_max <<
492  " and j = " << j_s_max << endl ;
493 
494  //----------------------------------------------------
495  // Equation of state
496  //----------------------------------------------------
497 
498  equation_of_state() ; // computes new values for nbar (n), ener (e)
499  // and press (p) from the new ent (H)
500 
501  //---------------------------------------------------------
502  // Matter source terms in the gravitational field equations
503  //---------------------------------------------------------
504 
505  hydro_euler() ; // computes new values for ener_euler (E),
506  // s_euler (S) and u_euler (U^i)
507 
508  //--------------------------------------------------------
509  // Poisson equation for logn_auto (nu_auto)
510  //--------------------------------------------------------
511 
512  // Source
513  // ------
514 
515  double unsgam1 ; // effective power of H in the source
516  // close to the surface
517 
518  if (relativistic) {
519  source = qpig * a_car % (ener_euler + s_euler)
523 
524  // 1/(gam-1) = dln(e)/dln(H) close to the surface :
525  unsgam1 = eos.der_ener_ent_p(ent_b + 1e-10*(ent_c-ent_b)) ;
526 
527  }
528  else {
529  source = qpig * nbar ;
530 
531  // 1/(gam-1) = dln(n)/dln(H) close to the surface :
532  unsgam1 = eos.der_nbar_ent_p(ent_b + 1e-10*(ent_c-ent_b)) ;
533  }
534 
535  source.set_std_base() ;
536 
537  // Resolution of the Poisson equation
538  // ----------------------------------
539 
543 
544  source_regu.std_base_scal() ;
545  source_div.std_base_scal() ;
546 
547  source().poisson_regular(k_div, nzet, unsgam1, par_poisson1,
549  logn_auto_div.set(),
551  source_regu, source_div) ;
552 
553  // Check: has the Poisson equation been correctly solved ?
554  // -----------------------------------------------------
555 
556  Tbl tdiff_logn = diffrel(logn_auto().laplacien(), source()) ;
557  cout <<
558  "Relative error in the resolution of the equation for logn_auto : "
559  << endl ;
560  for (int l=0; l<nz; l++) {
561  cout << tdiff_logn(l) << " " ;
562  }
563  cout << endl ;
564  diff_logn = max(abs(tdiff_logn)) ;
565 
566 
567  if (relativistic) {
568 
569  //--------------------------------------------------------
570  // Poisson equation for beta_auto
571  //--------------------------------------------------------
572 
573  // Source
574  // ------
575 
576  source = qpig * a_car % s_euler
577  + .75 * ( akcar_auto + akcar_comp )
582 
583  source.set_std_base() ;
584 
585  // Resolution of the Poisson equation
586  // ----------------------------------
587 
588  source().poisson(par_poisson2, beta_auto.set()) ;
589 
590 
591  // Check: has the Poisson equation been correctly solved ?
592  // -----------------------------------------------------
593 
594  Tbl tdiff_beta = diffrel(beta_auto().laplacien(), source()) ;
595  cout <<
596  "Relative error in the resolution of the equation for beta_auto : "
597  << endl ;
598  for (int l=0; l<nz; l++) {
599  cout << tdiff_beta(l) << " " ;
600  }
601  cout << endl ;
602  diff_beta = max(abs(tdiff_beta)) ;
603 
604  //--------------------------------------------------------
605  // Vector Poisson equation for shift_auto
606  //--------------------------------------------------------
607 
608  // Source
609  // ------
610 
611  Tenseur vtmp = 6. * ( d_beta_auto + d_beta_comp )
612  -8. * ( d_logn_auto + d_logn_comp ) ;
613 
614  source_shift = (-4.*qpig) * nnn % a_car % (ener_euler + press)
615  % u_euler
617 
618  source_shift.set_std_base() ;
619 
620  // Resolution of the Poisson equation
621  // ----------------------------------
622 
623  // Filter for the source of shift vector
624 
625  for (int i=0; i<3; i++) {
626 
627  if (source_shift(i).get_etat() != ETATZERO)
628  source_shift.set(i).filtre(4) ;
629 
630  }
631 
632  // For Tenseur::poisson_vect_regu,
633  // the triad must be the mapping triad,
634  // not the reference one:
635 
636  source_shift.change_triad( mp.get_bvect_cart() ) ;
637 
638  for (int i=0; i<3; i++) {
639  if(source_shift(i).dz_nonzero()) {
640  assert( source_shift(i).get_dzpuis() == 4 ) ;
641  }
642  else{
643  (source_shift.set(i)).set_dzpuis(4) ;
644  }
645  }
646 
647  //##
648  // source_shift.dec2_dzpuis() ; // dzpuis 4 -> 2
649 
650  double lambda_shift = double(1) / double(3) ;
651 
652  source_shift.poisson_vect_regu(k_div, nzet, unsgam1,
653  lambda_shift, par_poisson_vect,
655 
656 
657  // Check: has the equation for shift_auto been correctly solved ?
658  // --------------------------------------------------------------
659 
660  // Divergence of shift_auto :
661  Tenseur divn = contract(shift_auto.gradient(), 0, 1) ;
662  divn.dec2_dzpuis() ; // dzpuis 2 -> 0
663 
664  // Grad(div) :
665  Tenseur graddivn = divn.gradient() ;
666  graddivn.inc2_dzpuis() ; // dzpuis 2 -> 4
667 
668  // Full operator :
669  Tenseur lap_shift(mp, 1, CON, mp.get_bvect_cart() ) ;
670  lap_shift.set_etat_qcq() ;
671  for (int i=0; i<3; i++) {
672  lap_shift.set(i) = shift_auto(i).laplacien()
673  + lambda_shift * graddivn(i) ;
674  }
675 
676  Tbl tdiff_shift_x = diffrel(lap_shift(0), source_shift(0)) ;
677  Tbl tdiff_shift_y = diffrel(lap_shift(1), source_shift(1)) ;
678  Tbl tdiff_shift_z = diffrel(lap_shift(2), source_shift(2)) ;
679 
680  cout <<
681  "Relative error in the resolution of the equation "
682  "for shift_auto : "
683  << endl ;
684  cout << "x component : " ;
685  for (int l=0; l<nz; l++) {
686  cout << tdiff_shift_x(l) << " " ;
687  }
688  cout << endl ;
689  cout << "y component : " ;
690  for (int l=0; l<nz; l++) {
691  cout << tdiff_shift_y(l) << " " ;
692  }
693  cout << endl ;
694  cout << "z component : " ;
695  for (int l=0; l<nz; l++) {
696  cout << tdiff_shift_z(l) << " " ;
697  }
698  cout << endl ;
699 
700  diff_shift_x = max(abs(tdiff_shift_x)) ;
701  diff_shift_y = max(abs(tdiff_shift_y)) ;
702  diff_shift_z = max(abs(tdiff_shift_z)) ;
703 
704  // Final result
705  // ------------
706  // The output of Tenseur::poisson_vect is on the mapping triad,
707  // it should therefore be transformed to components on the
708  // reference triad :
709 
711 
712 
713  } // End of relativistic equations
714 
715 
716  //-------------------------------------------------
717  // Relative change in enthalpy
718  //-------------------------------------------------
719 
720  Tbl diff_ent_tbl = diffrel( ent(), ent_jm1() ) ;
721  diff_ent = diff_ent_tbl(0) ;
722  for (int l=1; l<nzet; l++) {
723  diff_ent += diff_ent_tbl(l) ;
724  }
725  diff_ent /= nzet ;
726 
727  ent_jm1 = ent ;
728 
729 
730  } // End of main loop
731 
732  //=========================================================================
733  // End of iteration
734  //=========================================================================
735 
736 
737 }
738 
739 }
int k_div
Index of regularity of the gravitational potential logn_auto .
Definition: etoile.h:452
void add_tenseur_mod(Tenseur &ti, int position=0)
Adds the address of a new modifiable Tenseur to the list.
Definition: param.C:1145
const Base_vect & ref_triad
Reference triad ("absolute frame"), with respect to which the components of all the member Tenseur &#39;s...
Definition: etoile.h:831
Component of a tensorial field *** DEPRECATED : use class Scalar instead ***.
Definition: cmp.h:446
void dec2_dzpuis()
dzpuis -= 2 ;
Definition: tenseur.C:1146
Radial mapping of rather general form.
Definition: map.h:2783
void add_int(const int &n, int position=0)
Adds the address of a new int to the list.
Definition: param.C:249
int get_np(int l) const
Returns the number of points in the azimuthal direction ( ) in domain no. l.
Definition: grilles.h:479
virtual double der_nbar_ent_p(double ent, const Param *par=0x0) const =0
Computes the logarithmic derivative from the log-enthalpy and extra parameters (virtual function imp...
Cmp sqrt(const Cmp &)
Square root.
Definition: cmp_math.C:223
void set_std_base()
Set the standard spectal basis of decomposition for each component.
Definition: tenseur.C:1186
Tenseur pot_centri
Centrifugal potential.
Definition: etoile.h:956
Cmp ssjm1_logn
Effective source at the previous step for the resolution of the Poisson equation for logn_auto by mea...
Definition: etoile.h:962
Tenseur logn_auto_regu
Regular part of the logarithm of the part of the lapse N generated principaly by the star...
Definition: etoile.h:494
Lorene prototypes.
Definition: app_hor.h:67
Standard units of space, time and mass.
const Mg3d * get_mg() const
Gives the Mg3d on which the mapping is defined.
Definition: map.h:783
double & set(int i)
Read/write of a particular element (index i) (1D case)
Definition: tbl.h:301
Tenseur nnn
Total lapse function.
Definition: etoile.h:512
Tenseur s_euler
Trace of the stress tensor in the Eulerian frame.
Definition: etoile.h:471
double ray_eq() const
Coordinate radius at , [r_unit].
virtual void adapt(const Cmp &ent, const Param &par, int nbr=0)=0
Adaptation of the mapping to a given scalar field.
Tenseur press
Fluid pressure.
Definition: etoile.h:464
Tenseur d_beta_auto
Gradient of beta_auto (Cartesian components with respect to ref_triad )
Definition: etoile.h:882
Tbl diffrel(const Cmp &a, const Cmp &b)
Relative difference between two Cmp (norme version).
Definition: cmp_math.C:507
void inc2_dzpuis()
dzpuis += 2 ;
Definition: tenseur.C:1159
Tenseur shift_auto
Part of the shift vector generated principaly by the star.
Definition: etoile.h:892
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition: tbl.C:364
Tenseur logn_auto_div
Divergent part (if k_div!=0 ) of the logarithm of the part of the lapse N generated principaly by t...
Definition: etoile.h:500
bool irrotational
true for an irrotational star, false for a corotating one
Definition: etoile.h:825
virtual double val_r(int l, double xi, double theta, double pphi) const =0
Returns the value of the radial coordinate r for a given in a given domain.
Cmp & set()
Read/write for a scalar (see also operator=(const Cmp&) ).
Definition: tenseur.C:840
Tenseur w_shift
Vector used in the decomposition of shift_auto , following Shibata&#39;s prescription [Prog...
Definition: etoile.h:911
Tenseur u_euler
Fluid 3-velocity with respect to the Eulerian observer.
Definition: etoile.h:477
Tenseur d_beta_comp
Gradient of beta_comp (Cartesian components with respect to ref_triad )
Definition: etoile.h:887
void change_triad(const Base_vect &new_triad)
Sets a new vectorial basis (triad) of decomposition and modifies the components accordingly.
Definition: tenseur.C:684
Tenseur flat_scalar_prod_desal(const Tenseur &t1, const Tenseur &t2)
Same as flat_scalar_prod but with desaliasing.
Tenseur logn_comp
Part of the lapse logarithm (gravitational potential at the Newtonian limit) generated principaly by ...
Definition: etoile.h:857
Tenseur d_logn_auto
Gradient of logn_auto (Cartesian components with respect to ref_triad )
Definition: etoile.h:862
Tenseur nbar
Baryon density in the fluid frame.
Definition: etoile.h:462
void filtre(int n)
Sets the n lasts coefficients in r to 0 in the external domain.
Definition: cmp_manip.C:77
Tenseur akcar_comp
Part of the scalar generated by shift_auto and shift_comp , i.e.
Definition: etoile.h:947
Parameter storage.
Definition: param.h:125
void add_tbl(const Tbl &ti, int position=0)
Adds the address of a new Tbl to the list.
Definition: param.C:525
Map & mp
Mapping associated with the star.
Definition: etoile.h:432
int get_nzone() const
Returns the number of domains.
Definition: grilles.h:465
virtual void equation_of_state()
Computes the proper baryon and energy density, as well as pressure from the enthalpy.
Definition: etoile.C:569
const Eos & eos
Equation of state of the stellar matter.
Definition: etoile.h:454
virtual void homothetie(double lambda)
Sets a new radial scale.
Definition: map_et.C:928
Tbl max(const Cmp &)
Maximum values of a Cmp in each domain.
Definition: cmp_math.C:438
Tenseur khi_shift
Scalar used in the decomposition of shift_auto , following Shibata&#39;s prescription [Prog...
Definition: etoile.h:921
Tenseur contract(const Tenseur &, int id1, int id2)
Self contraction of two indices of a Tenseur .
Tenseur_sym tkij_auto
Part of the extrinsic curvature tensor generated by shift_auto .
Definition: etoile.h:928
virtual void resize(int l, double lambda)=0
Rescales the outer boundary of one domain.
void std_base_scal()
Sets the spectral bases of the Valeur va to the standard ones for a scalar.
Definition: cmp.C:647
int nzet
Number of domains of *mp occupied by the star.
Definition: etoile.h:435
int get_nr(int l) const
Returns the number of points in the radial direction ( ) in domain no. l.
Definition: grilles.h:469
Tenseur a_car
Total conformal factor .
Definition: etoile.h:518
bool relativistic
Indicator of relativity: true for a relativistic star, false for a Newtonian one. ...
Definition: etoile.h:440
Multi-domain grid.
Definition: grilles.h:279
double ray_pole() const
Coordinate radius at [r_unit].
Cmp ssjm1_beta
Effective source at the previous step for the resolution of the Poisson equation for beta_auto by mea...
Definition: etoile.h:968
virtual void reevaluate_symy(const Map *mp_prev, int nzet, Cmp &uu) const =0
Recomputes the values of a Cmp at the collocation points after a change in the mapping.
const Base_vect_cart & get_bvect_cart() const
Returns the Cartesian basis associated with the coordinates (x,y,z) of the mapping, i.e.
Definition: map.h:809
Tenseur logn_auto
Total of the logarithm of the part of the lapse N generated principaly by the star.
Definition: etoile.h:487
void equil_regular(double ent_c, int mermax, int mermax_poisson, double relax_poisson, int mermax_potvit, double relax_potvit, double thres_adapt, const Tbl &fact, Tbl &diff)
Computes an equilibrium configuration by regularizing the diverging source.
Tenseur loggam
Logarithm of the Lorentz factor between the fluid and the co-orbiting observer.
Definition: etoile.h:852
Tenseur ent
Log-enthalpy (relativistic case) or specific enthalpy (Newtonian case)
Definition: etoile.h:460
Cmp abs(const Cmp &)
Absolute value.
Definition: cmp_math.C:413
void add_double(const double &x, int position=0)
Adds the the address of a new double to the list.
Definition: param.C:318
double velocity_potential(int mermax, double precis, double relax)
Computes the non-translational part of the velocity scalar potential by solving the continuity equat...
Basic array class.
Definition: tbl.h:164
int get_nt(int l) const
Returns the number of points in the co-latitude direction ( ) in domain no. l.
Definition: grilles.h:474
Tenseur beta_auto
Logarithm of the part of the product AN generated principaly by by the star.
Definition: etoile.h:509
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition: tenseur.C:652
void add_cmp_mod(Cmp &ti, int position=0)
Adds the address of a new modifiable Cmp to the list.
Definition: param.C:1007
Tenseur akcar_auto
Part of the scalar generated by shift_auto , i.e.
Definition: etoile.h:941
Tenseur ener_euler
Total energy density in the Eulerian frame.
Definition: etoile.h:468
void poisson_vect_regu(int k_div, int nzet, double unsgam1, double lambda, Param &par, Tenseur &shift, Tenseur &vect, Tenseur &scal) const
Solves the vectorial Poisson equation : .
Tenseur d_logn_comp
Gradient of logn_comp (Cartesian components with respect to ref_triad )
Definition: etoile.h:872
Tenseur d_logn_auto_div
Gradient of logn_auto_div (if k_div!=0 )
Definition: etoile.h:504
virtual double der_ener_ent_p(double ent, const Param *par=0x0) const =0
Computes the logarithmic derivative from the log-enthalpy with extra parameters (virtual function im...
Tensor handling *** DEPRECATED : use class Tensor instead ***.
Definition: tenseur.h:304
void add_int_mod(int &n, int position=0)
Adds the address of a new modifiable int to the list.
Definition: param.C:388
const Tenseur & gradient() const
Returns the gradient of *this (Cartesian coordinates)
Definition: tenseur.C:1558
Cmp ssjm1_khi
Effective source at the previous step for the resolution of the Poisson equation for the scalar by m...
Definition: etoile.h:976
virtual void hydro_euler()
Computes the hydrodynamical quantities relative to the Eulerian observer from those in the fluid fram...
Definition: et_bin_hydro.C:109
Tenseur ssjm1_wshift
Effective source at the previous step for the resolution of the vector Poisson equation for by means...
Definition: etoile.h:986