LORENE
et_magnetisation.C
1 /*
2  * Methods for axisymmetric rotating neutron stars with magnetisation in the stress-energy tensor.
3  *
4  * See the file et_rot_mag.h for documentation
5  *
6  */
7 
8 /*
9  * Copyright (c) 2013-2014 Jerome Novak, Debarti Chatterjee, Micaela Oertel
10  *
11  * This file is part of LORENE.
12  *
13  * LORENE is free software; you can redistribute it and/or modify
14  * it under the terms of the GNU General Public License as published by
15  * the Free Software Foundation; either version 2 of the License, or
16  * (at your option) any later version.
17  *
18  * LORENE is distributed in the hope that it will be useful,
19  * but WITHOUT ANY WARRANTY; without even the implied warranty of
20  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
21  * GNU General Public License for more details.
22  *
23  * You should have received a copy of the GNU General Public License
24  * along with LORENE; if not, write to the Free Software
25  * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
26  *
27  */
28 
29 
30 
31 /*
32  * $Id: et_magnetisation.C,v 1.11 2016/12/05 16:17:53 j_novak Exp $
33  * $Log: et_magnetisation.C,v $
34  * Revision 1.11 2016/12/05 16:17:53 j_novak
35  * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
36  *
37  * Revision 1.10 2014/10/21 09:23:53 j_novak
38  * Addition of global functions mass_g(), angu_mom(), grv2/3() and mom_quad().
39  *
40  * Revision 1.9 2014/10/13 08:52:56 j_novak
41  * Lorene classes and functions now belong to the namespace Lorene.
42  *
43  * Revision 1.8 2014/05/28 07:46:06 j_novak
44  * Minor modifications.
45  *
46  * Revision 1.7 2014/05/27 12:32:28 j_novak
47  * Added possibility to converge to a given magnetic moment.
48  *
49  * Revision 1.6 2014/05/14 15:19:05 j_novak
50  * The magnetisation field is now filtered.
51  *
52  * Revision 1.5 2014/04/29 13:46:07 j_novak
53  * Addition of switches 'use_B_in_eos' and 'include_magnetisation' to control the model.
54  *
55  * Revision 1.4 2014/04/28 12:48:13 j_novak
56  * Minor modifications.
57  *
58  * Revision 1.3 2013/12/19 17:05:40 j_novak
59  * Corrected a dzpuis problem.
60  *
61  * Revision 1.2 2013/12/13 16:36:51 j_novak
62  * Addition and computation of magnetisation terms in the Einstein equations.
63  *
64  * Revision 1.1 2013/11/25 13:52:11 j_novak
65  * New class Et_magnetisation to include magnetization terms in the stress energy tensor.
66  *
67  *
68  * $Header: /cvsroot/Lorene/C++/Source/Etoile/et_magnetisation.C,v 1.11 2016/12/05 16:17:53 j_novak Exp $
69  *
70  */
71 
72 // Headers C
73 #include <cmath>
74 
75 // Headers Lorene
76 #include "et_rot_mag.h"
77 #include "utilitaires.h"
78 #include "unites.h"
79 #include "param.h"
80 #include "eos.h"
81 
82  //--------------//
83  // Constructors //
84  //--------------//
85 // Standard constructor
86 // --------------------
87 
88 
89 namespace Lorene {
90 Et_magnetisation::Et_magnetisation(Map& mp_i, int nzet_i, bool relat,
91  const Eos& eos_i, bool magnetisation,
92  bool B_eos)
93  : Et_rot_mag(mp_i, nzet_i, relat, eos_i, 1),
94  use_B_in_eos(B_eos), include_magnetisation(magnetisation),
95  xmag(mp_i),
96  E_I(mp_i),
97  J_I(mp_i, COV, mp_i.get_bvect_spher()),
98  Sij_I(mp_i, COV, mp_i.get_bvect_spher())
99 {
100  const Eos_mag* tmp_p_eos = dynamic_cast<const Eos_mag*>(&eos_i) ;
101  if (tmp_p_eos == 0x0) {
102  cerr << "Et_magnetisation::Et_magnetisation : " << endl ;
103  cerr << "Only magnetised EoS is admitted for this class!" << endl ;
104  cerr << "Aborting ... " << endl ;
105  abort() ;
106  }
107  if ( include_magnetisation && (!use_B_in_eos) ) {
108  cerr << "Et_magnetisation::Et_magnetisation : " << endl ;
109  cerr << "Magnetisation terms can be included only if "
110  << "the magnetic field is used in the EoS!" << endl;
111  cerr << "Aborting ... " << endl ;
112  abort() ;
113  }
114  xmag = 0 ;
115  E_I = 0 ;
116  J_I.set_etat_zero() ;
117  Sij_I.set_etat_zero() ;
118 
119 }
120 
121 
122 
123 Et_magnetisation::Et_magnetisation(Map& mp_i, const Eos& eos_i, FILE* fich)
124  : Et_rot_mag(mp_i, eos_i, fich, 0),
125  use_B_in_eos(true), include_magnetisation(true),
126  xmag(mp_i),
127  E_I(mp_i),
128  J_I(mp_i, COV, mp_i.get_bvect_spher()),
129  Sij_I(mp_i, COV, mp_i.get_bvect_spher())
130 {
131 
132  // Read of the saved fields:
133  // ------------------------
134 
135  fread(&use_B_in_eos, sizeof(bool), 1, fich) ;
136 
137  fread(&include_magnetisation, sizeof(bool), 1, fich) ;
138 
139  Scalar xmag_file(mp, *mp.get_mg(), fich) ;
140  xmag = xmag_file ;
141 
142  Scalar E_I_file(mp, *mp.get_mg(), fich) ;
143  E_I = E_I_file ;
144 
145  Vector J_I_file(mp, mp.get_bvect_spher(), fich) ;
146  J_I = J_I_file ;
147 
148  Sym_tensor Sij_I_file(mp, mp.get_bvect_spher(), fich) ;
149  Sij_I = Sij_I_file ;
150 
151 }
152 
153 
154 // Copy constructor
155 // ----------------
156 
158  : Et_rot_mag(et),
159  use_B_in_eos(et.use_B_in_eos), include_magnetisation(et.include_magnetisation),
160  xmag(et.xmag),
161  E_I(et.E_I),
162  J_I(et.J_I),
163  Sij_I(et.Sij_I)
164 { }
165 
166 
167  //------------//
168  // Destructor //
169  //------------//
170 
172 
173 
174 // Assignment to another Et_magnetisation
175 // --------------------------------
176 
178 
179  // Assignement of quantities common to all the derived classes of Et_rot_mag
183  xmag = et.xmag ;
184  E_I = et.E_I ;
185  J_I = et.J_I ;
186  Sij_I = et.Sij_I ;
187 
188 }
189 
190 
192 
193  Cmp ent_eos = ent() ;
194  const Eos_mag* mageos = dynamic_cast<const Eos_mag*>(&eos) ;
195  assert (mageos != 0x0) ;
196 
197  // Slight rescale of the enthalpy field in case of 2 domains inside the
198  // star
199 
200  double epsilon = 1.e-12 ;
201 
202  const Mg3d* mg = mp.get_mg() ;
203  int nz = mg->get_nzone() ;
204 
205  Mtbl xi(mg) ;
206  xi.set_etat_qcq() ;
207  for (int l=0; l<nz; l++) {
208  xi.t[l]->set_etat_qcq() ;
209  for (int k=0; k<mg->get_np(l); k++) {
210  for (int j=0; j<mg->get_nt(l); j++) {
211  for (int i=0; i<mg->get_nr(l); i++) {
212  xi.set(l,k,j,i) = mg->get_grille3d(l)->x[i] ;
213  }
214  }
215  }
216 
217  }
218 
219  Cmp fact_ent(mp) ;
220  fact_ent.allocate_all() ;
221 
222  fact_ent.set(0) = 1 + epsilon * xi(0) * xi(0) ;
223  fact_ent.set(1) = 1 - 0.25 * epsilon * (xi(1) - 1) * (xi(1) - 1) ;
224 
225  for (int l=nzet; l<nz; l++) {
226  fact_ent.set(l) = 1 ;
227  }
228 
229  if (nzet > 1) {
230 
231  if (nzet == 3) {
232  fact_ent.set(1) = 1 - 0.5 * epsilon * (xi(1) - 0.5) * (xi(1) - 0.5) ;
233  fact_ent.set(2) = 1 - 0.25 * epsilon * (xi(2) - 1) * (xi(2) - 1) ;
234  }
235 
236  if (nzet > 3) {
237  cerr << "Et_magnetisation::equation_of_state: "
238  << "not ready yet for nzet > 3 !" << endl ;
239  }
240 
241  ent_eos = fact_ent * ent_eos ;
242  ent_eos.std_base_scal() ;
243  }
244 
245 
246 
247  // Call to a magnetized EOS
248  // with the norm of the magnetic field passed as an argument to the EoS.
249 
250  double magb0 = 0 ;
251  Cmp norm_b(mp) ;
252  norm_b.set_etat_zero() ;
253  if (use_B_in_eos) {
254  Tenseur Bmag = Magn() ;
255  norm_b = sqrt(a_car() * ( Bmag(0)*Bmag(0) + Bmag(1)*Bmag(1) ))
256  / gam_euler() ; // Use of b^\mu in the fluid frame
257  norm_b.std_base_scal() ;
258  }
259  Param par ;
260  par.add_double_mod(magb0) ;
261 
268  xmag.allocate_all() ;
269 
270  for (int l=0; l< nz; l++) {
271  Tbl* tent = ent_eos.va.c->t[l] ;
272  if ( (tent->get_etat() == ETATZERO) || (l >= nzet) ) {
273  nbar.set().set(l).set_etat_zero() ;
274  ener.set().set(l).set_etat_zero() ;
275  press.set().set(l).set_etat_zero() ;
276  xmag.annule_domain(l) ;
277  }
278  else {
279  nbar.set().set(l).annule_hard() ;
280  ener.set().set(l).annule_hard() ;
281  press.set().set(l).annule_hard() ;
283  for (int k=0; k < mg->get_np(l); k++) {
284  for (int j=0; j < mg->get_nt(l); j++) {
285  for (int i=0; i < mg->get_nr(l); i++) {
286  magb0 = norm_b(l, k, j, i) ;
287  double ent0 = ent_eos(l, k, j, i) ;
288  nbar.set().set(l, k, j, i) = mageos->nbar_ent_p(ent0, &par) ;
289  ener.set().set(l, k, j, i) = mageos->ener_ent_p(ent0, &par) ;
290  press.set().set(l, k, j, i) = mageos->press_ent_p(ent0, &par) ;
291  xmag.set_grid_point(l, k, j, i) = mageos->mag_ent_p(ent0, &par) ;
292  }
293  }
294  }
295  }
296  }
297 
299  xmag.set_etat_zero() ;
300 
301  // Set the bases for spectral expansion
302  nbar.set_std_base() ;
303  ener.set_std_base() ;
304  press.set_std_base() ;
306 
307  Scalar tmp_scal(xmag) ;
308  tmp_scal.exponential_filter_r(0, 0, 1) ;
309  xmag = tmp_scal ;
310 
311  // The derived quantities are obsolete
312  del_deriv() ;
313 
314 }
315 
316 
317 
318  //--------------//
319  // Outputs //
320  //--------------//
321 
322 // Save in a file
323 // --------------
324 void Et_magnetisation::sauve(FILE* fich) const {
325 
326  Et_rot_mag::sauve(fich) ;
327 
328  fwrite(&use_B_in_eos, sizeof(bool), 1, fich) ;
329 
330  fwrite(&include_magnetisation, sizeof(bool), 1, fich) ;
331 
332  xmag.sauve(fich) ;
333  E_I.sauve(fich) ;
334  J_I.sauve(fich) ;
335  Sij_I.sauve(fich) ;
336 
337 }
338 
339 
340 // Printing
341 // --------
342 
343 
344 ostream& Et_magnetisation::operator>>(ostream& ost) const {
345 
346  using namespace Unites_mag ;
347 
349  // int theta_eq = mp.get_mg()->get_nt(nzet-1)-1 ;
350  ost << endl ;
351  ost << "Rotating magnetized neutron star"
352  << endl ;
353  if (use_B_in_eos)
354  ost << "Using magnetic field in the EoS" << endl ;
355  if (include_magnetisation) {
356  ost << "Including magnetisation terms in the equations" << endl ;
357  ost << "Maximal value of the magnetization scalar x : "
358  << max(maxabs(xmag)) << endl ;
359  }
360  return ost ;
361 }
362 
363 //=================================================================================
364 //
365 // Computation of equilibrium configuration
366 //
367 //=================================================================================
368 
369 
370 void Et_magnetisation::equilibrium_mag(double ent_c, double omega0,
371  double fact_omega, int nzadapt, const Tbl& ent_limit,
372  const Itbl& icontrol, const Tbl& control, double mbar_wanted,
373  double magmom_wanted,
374  double aexp_mass, Tbl& diff, double Q0, double a_j0,
375  Cmp (*f_j)(const Cmp&, const double),
376  Cmp (*M_j)(const Cmp& x, const double)) {
377 
378  // Fundamental constants and units
379  // -------------------------------
380  using namespace Unites_mag ;
381 
382  // For the display
383  // ---------------
384  char display_bold[]="x[1m" ; display_bold[0] = 27 ;
385  char display_normal[] = "x[0m" ; display_normal[0] = 27 ;
386 
387  // Grid parameters
388  // ---------------
389 
390  const Mg3d* mg = mp.get_mg() ;
391  int nz = mg->get_nzone() ; // total number of domains
392  int nzm1 = nz - 1 ;
393 
394  // The following is required to initialize mp_prev as a Map_et:
395  Map_et& mp_et = dynamic_cast<Map_et&>(mp) ;
396 
397  // Index of the point at phi=0, theta=pi/2 at the surface of the star:
398  assert(mg->get_type_t() == SYM) ;
399  int l_b = nzet - 1 ;
400  int i_b = mg->get_nr(l_b) - 1 ;
401  int j_b = mg->get_nt(l_b) - 1 ;
402  int k_b = 0 ;
403 
404  // Value of the enthalpy defining the surface of the star
405  double ent_b = ent_limit(nzet-1) ;
406 
407  // Parameters to control the iteration
408  // -----------------------------------
409 
410  int mer_max = icontrol(0) ;
411  int mer_rot = icontrol(1) ;
412  int mer_change_omega = icontrol(2) ;
413  int mer_fix_omega = icontrol(3) ;
414  int mer_mass = icontrol(4) ;
415  int mermax_poisson = icontrol(5) ;
416  int delta_mer_kep = icontrol(6) ;
417  int mer_mag = icontrol(7) ;
418  int mer_change_mag = icontrol(8) ;
419  int mer_fix_mag = icontrol(9) ;
420  int mer_magmom = icontrol(10) ;
421 
422  // Protections:
423  if (mer_change_omega < mer_rot) {
424  cerr << "Et_magnetisation::equilibrium_mag: mer_change_omega < mer_rot !"
425  << endl ;
426  cerr << " mer_change_omega = " << mer_change_omega << endl ;
427  cerr << " mer_rot = " << mer_rot << endl ;
428  abort() ;
429  }
430  if (mer_fix_omega < mer_change_omega) {
431  cerr << "Et_magnetisation::equilibrium_mag: mer_fix_omega < mer_change_omega !"
432  << endl ;
433  cerr << " mer_fix_omega = " << mer_fix_omega << endl ;
434  cerr << " mer_change_omega = " << mer_change_omega << endl ;
435  abort() ;
436  }
437 
438  // In order to converge to a given baryon mass, shall the central
439  // enthalpy be varied or Omega ?
440  bool change_ent = true ;
441  if (mer_mass < 0) {
442  change_ent = false ;
443  mer_mass = abs(mer_mass) ;
444  }
445 
446  double precis = control(0) ;
447  double omega_ini = control(1) ;
448  double relax = control(2) ;
449  double relax_prev = double(1) - relax ;
450  double relax_poisson = control(3) ;
451  double thres_adapt = control(4) ;
452  double precis_adapt = control(5) ;
453  double Q_ini = control(6) ;
454  double a_j_ini = control (7) ;
455 
456  // Error indicators
457  // ----------------
458 
459  diff.set_etat_qcq() ;
460  double& diff_ent = diff.set(0) ;
461 
462  // Parameters for the function Map_et::adapt
463  // -----------------------------------------
464 
465  Param par_adapt ;
466  int nitermax = 100 ;
467  int niter ;
468  int adapt_flag = 1 ; // 1 = performs the full computation,
469  // 0 = performs only the rescaling by
470  // the factor alpha_r
471  int nz_search = nzet + 1 ; // Number of domains for searching the enthalpy
472  // isosurfaces
473  double alpha_r ;
474  double reg_map = 1. ; // 1 = regular mapping, 0 = contracting mapping
475 
476  par_adapt.add_int(nitermax, 0) ; // maximum number of iterations to
477  // locate zeros by the secant method
478  par_adapt.add_int(nzadapt, 1) ; // number of domains where the adjustment
479  // to the isosurfaces of ent is to be
480  // performed
481  par_adapt.add_int(nz_search, 2) ; // number of domains to search for
482  // the enthalpy isosurface
483  par_adapt.add_int(adapt_flag, 3) ; // 1 = performs the full computation,
484  // 0 = performs only the rescaling by
485  // the factor alpha_r
486  par_adapt.add_int(j_b, 4) ; // theta index of the collocation point
487  // (theta_*, phi_*)
488  par_adapt.add_int(k_b, 5) ; // theta index of the collocation point
489  // (theta_*, phi_*)
490 
491  par_adapt.add_int_mod(niter, 0) ; // number of iterations actually used in
492  // the secant method
493 
494  par_adapt.add_double(precis_adapt, 0) ; // required absolute precision in
495  // the determination of zeros by
496  // the secant method
497  par_adapt.add_double(reg_map, 1) ; // 1. = regular mapping,
498  // 0 = contracting mapping
499 
500  par_adapt.add_double(alpha_r, 2) ; // factor by which all the radial
501  // distances will be multiplied
502 
503  par_adapt.add_tbl(ent_limit, 0) ; // array of values of the field ent
504  // to define the isosurfaces.
505 
506 
507  // Parameters for the function Map_et::poisson for nuf
508  // ----------------------------------------------------
509 
510  double precis_poisson = 1.e-16 ;
511 
512  Param par_poisson_nuf ;
513  par_poisson_nuf.add_int(mermax_poisson, 0) ; // maximum number of iterations
514  par_poisson_nuf.add_double(relax_poisson, 0) ; // relaxation parameter
515  par_poisson_nuf.add_double(precis_poisson, 1) ; // required precision
516  par_poisson_nuf.add_int_mod(niter, 0) ; // number of iterations actually used
517  par_poisson_nuf.add_cmp_mod( ssjm1_nuf ) ;
518 
519  Param par_poisson_nuq ;
520  par_poisson_nuq.add_int(mermax_poisson, 0) ; // maximum number of iterations
521  par_poisson_nuq.add_double(relax_poisson, 0) ; // relaxation parameter
522  par_poisson_nuq.add_double(precis_poisson, 1) ; // required precision
523  par_poisson_nuq.add_int_mod(niter, 0) ; // number of iterations actually used
524  par_poisson_nuq.add_cmp_mod( ssjm1_nuq ) ;
525 
526  Param par_poisson_tggg ;
527  par_poisson_tggg.add_int(mermax_poisson, 0) ; // maximum number of iterations
528  par_poisson_tggg.add_double(relax_poisson, 0) ; // relaxation parameter
529  par_poisson_tggg.add_double(precis_poisson, 1) ; // required precision
530  par_poisson_tggg.add_int_mod(niter, 0) ; // number of iterations actually used
531  par_poisson_tggg.add_cmp_mod( ssjm1_tggg ) ;
532  double lambda_tggg ;
533  par_poisson_tggg.add_double_mod( lambda_tggg ) ;
534 
535  Param par_poisson_dzeta ;
536  double lbda_grv2 ;
537  par_poisson_dzeta.add_double_mod( lbda_grv2 ) ;
538 
539  // Parameters for the function Tenseur::poisson_vect
540  // -------------------------------------------------
541 
542  Param par_poisson_vect ;
543 
544  par_poisson_vect.add_int(mermax_poisson, 0) ; // maximum number of iterations
545  par_poisson_vect.add_double(relax_poisson, 0) ; // relaxation parameter
546  par_poisson_vect.add_double(precis_poisson, 1) ; // required precision
547  par_poisson_vect.add_cmp_mod( ssjm1_khi ) ;
548  par_poisson_vect.add_tenseur_mod( ssjm1_wshift ) ;
549  par_poisson_vect.add_int_mod(niter, 0) ;
550 
551  // Parameters for the Maxwell equations
552  // -------------------------------------
553 
554  Param par_poisson_At ; // For scalar At Poisson equation
555  Cmp ssjm1_At(mp) ;
556  ssjm1_At.set_etat_zero() ;
557  par_poisson_At.add_int(mermax_poisson, 0) ; // maximum number of iterations
558  par_poisson_At.add_double(relax_poisson, 0) ; // relaxation parameter
559  par_poisson_At.add_double(precis_poisson, 1) ; // required precision
560  par_poisson_At.add_int_mod(niter, 0) ; // number of iterations actually used
561  par_poisson_At.add_cmp_mod( ssjm1_At ) ;
562 
563  Param par_poisson_Avect ; // For vector Aphi Poisson equation
564 
565  Cmp ssjm1_khi_mag(ssjm1_khi) ;
566  Tenseur ssjm1_w_mag(ssjm1_wshift) ;
567 
568  par_poisson_Avect.add_int(mermax_poisson, 0) ; // maximum number of iterations
569  par_poisson_Avect.add_double(relax_poisson, 0) ; // relaxation parameter
570  par_poisson_Avect.add_double(precis_poisson, 1) ; // required precision
571  par_poisson_Avect.add_cmp_mod( ssjm1_khi_mag ) ;
572  par_poisson_Avect.add_tenseur_mod( ssjm1_w_mag ) ;
573  par_poisson_Avect.add_int_mod(niter, 0) ;
574 
575  // Initializations
576  // ---------------
577 
578  // Initial angular velocity / magnetic quantities
579  omega = 0 ;
580  Q = 0 ;
581  a_j = 0 ;
582 
583  double accrois_omega = (omega0 - omega_ini)
584  / double(mer_fix_omega - mer_change_omega) ;
585  double accrois_Q = (Q0 - Q_ini)
586  / double(mer_fix_mag - mer_change_mag);
587  double accrois_a_j = (a_j0 - a_j_ini)
588  / double(mer_fix_mag - mer_change_mag);
589 
590  update_metric() ; // update of the metric coefficients
591 
592  equation_of_state() ; // update of the density, pressure, etc...
593 
594  hydro_euler() ; // update of the hydro quantities relative to the
595  // Eulerian observer
596 
597  MHD_comput() ; // update of EM contributions to stress-energy tensor
598 
599 
600  // Quantities at the previous step :
601  Map_et mp_prev = mp_et ;
602  Tenseur ent_prev = ent ;
603  Tenseur logn_prev = logn ;
604  Tenseur dzeta_prev = dzeta ;
605 
606  // Creation of uninitialized tensors:
607  Tenseur source_nuf(mp) ; // source term in the equation for nuf
608  Tenseur source_nuq(mp) ; // source term in the equation for nuq
609  Tenseur source_dzf(mp) ; // matter source term in the eq. for dzeta
610  Tenseur source_dzq(mp) ; // quadratic source term in the eq. for dzeta
611  Tenseur source_tggg(mp) ; // source term in the eq. for tggg
612  Tenseur source_shift(mp, 1, CON, mp.get_bvect_cart()) ; // source term for shift
613  Tenseur mlngamma(mp) ; // centrifugal potential
614 
615  // Preparations for the Poisson equations:
616  // --------------------------------------
617  if (nuf.get_etat() == ETATZERO) {
618  nuf.set_etat_qcq() ;
619  nuf.set() = 0 ;
620  }
621 
622  if (relativistic) {
623  if (nuq.get_etat() == ETATZERO) {
624  nuq.set_etat_qcq() ;
625  nuq.set() = 0 ;
626  }
627 
628  if (tggg.get_etat() == ETATZERO) {
629  tggg.set_etat_qcq() ;
630  tggg.set() = 0 ;
631  }
632 
633  if (dzeta.get_etat() == ETATZERO) {
634  dzeta.set_etat_qcq() ;
635  dzeta.set() = 0 ;
636  }
637  }
638 
639  ofstream fichconv("convergence.d") ; // Output file for diff_ent
640  fichconv << "# diff_ent GRV2 " << endl ;
641 
642  ofstream fichfreq("frequency.d") ; // Output file for omega
643  fichfreq << "# f [Hz]" << endl ;
644 
645  ofstream fichevol("evolution.d") ; // Output file for various quantities
646  fichevol <<
647  "# |dH/dr_eq/dH/dr_pole| r_pole/r_eq ent_c"
648  << endl ;
649 
650  diff_ent = 1 ;
651  double err_grv2 = 1 ;
652 
653  //=========================================================================
654  // Start of iteration
655  //=========================================================================
656 
657  for(int mer=0 ; (diff_ent > precis) && (mer<mer_max) ; mer++ ) {
658 
659  cout << "-----------------------------------------------" << endl ;
660  cout << "step: " << mer << endl ;
661  cout << "diff_ent = " << display_bold << diff_ent << display_normal
662  << endl ;
663  cout << "err_grv2 = " << err_grv2 << endl ;
664  fichconv << mer ;
665  fichfreq << mer ;
666  fichevol << mer ;
667 
668  if (mer >= mer_rot) {
669 
670  if (mer < mer_change_omega) {
671  omega = omega_ini ;
672  }
673  else {
674  if (mer <= mer_fix_omega) {
675  omega = omega_ini + accrois_omega *
676  (mer - mer_change_omega) ;
677  }
678  }
679 
680  }
681 
682  if (mer >= mer_mag) {
683  if (mer < mer_change_mag) {
684  Q = Q_ini ;
685  a_j = a_j_ini ;
686  }
687  else {
688  if (mer <= mer_fix_mag) {
689  Q = Q_ini + accrois_Q * (mer - mer_change_mag) ;
690  a_j = a_j_ini + accrois_a_j * (mer - mer_change_mag) ;
691  }
692  }
693  }
694 
695  //-----------------------------------------------
696  // Computation of electromagnetic potentials :
697  // -------------------------------------------
698  magnet_comput(adapt_flag, f_j, par_poisson_At, par_poisson_Avect) ;
699 
700  MHD_comput() ; // computes EM contributions to T_{mu,nu}
701 
702  // S_{rr} + S_{\theta\theta}
703  Tenseur SrrplusStt( Cmp(Sij_I(1, 1) + Sij_I(2, 2)) ) ;
704  SrrplusStt = SrrplusStt / a_car ; // S^r_r + S^\theta_\theta
705 
706  Tenseur Spp (Cmp(Sij_I(3, 3))) ; //S_{\phi\phi}
707  Spp = Spp / b_car ; // S^\phi_\phi
708 
709  Cmp temp(E_I) ;
710  Tenseur E_int(temp) ;
711 
712  //-----------------------------------------------
713  // Sources of the Poisson equations
714  //-----------------------------------------------
715 
716  // Source for nu
717  // -------------
718  Tenseur beta = log(bbb) ;
719  beta.set_std_base() ;
720 
721  if (relativistic) {
722  source_nuf =
723  qpig * a_car *( ener_euler + s_euler + E_int + SrrplusStt + Spp ) ;
724 
725  source_nuq = ak_car
727  + beta.gradient_spher())
728  + qpig * a_car * 2*E_em ;
729  }
730  else {
731  source_nuf = qpig * nbar ;
732 
733  source_nuq = 0 ;
734  }
735  source_nuf.set_std_base() ;
736  source_nuq.set_std_base() ;
737 
738  // Source for dzeta
739  // ----------------
740  source_dzf = 2 * qpig * a_car
741  * (press + (ener_euler+press) * uuu*uuu + Spp ) ;
742  source_dzf.set_std_base() ;
743 
744  source_dzq = 2 * qpig * a_car * E_em + 1.5 * ak_car -
746  source_dzq.set_std_base() ;
747 
748 
749  // Source for tggg
750  // ---------------
751 
752  source_tggg = 2 * qpig * nnn * a_car * bbb
753  * ( 2*press + SrrplusStt ) ;
754  source_tggg.set_std_base() ;
755 
756  (source_tggg.set()).mult_rsint() ;
757 
758 
759  // Source for shift
760  // ----------------
761 
762  // Matter term:
763 
764  Cmp tjpem(J_I(3)) ;
765  tjpem += Jp_em() ;
766  tjpem.div_rsint() ;
767 
768  source_shift = (-4*qpig) * nnn * a_car * (ener_euler + press) * u_euler ;
769 
770  // Quadratic terms:
771  Tenseur vtmp = 6 * beta.gradient_spher() - 2 * logn.gradient_spher() ;
772  Tenseur mtmp(mp, 1, COV, mp.get_bvect_spher()) ;
773  if (tjpem.get_etat() == ETATZERO) mtmp.set_etat_zero() ;
774  else {
775  mtmp.set_etat_qcq() ;
776  mtmp.set(0) = 0 ;
777  mtmp.set(1) = 0 ;
778  mtmp.set(2) = (-4*qpig)*tjpem*nnn()*a_car()/b_car() ;
779  }
780  mtmp.change_triad(mp.get_bvect_cart()) ;
781 
782  vtmp.change_triad(mp.get_bvect_cart()) ;
783 
784  Tenseur squad = nnn * flat_scalar_prod(tkij, vtmp) ;
785 
786  // The addition of matter terms and quadratic terms is performed
787  // component by component because u_euler is contravariant,
788  // while squad is covariant.
789 
790  if (squad.get_etat() == ETATQCQ) {
791  for (int i=0; i<3; i++) {
792  source_shift.set(i) += squad(i) ;
793  }
794  }
795  if (mtmp.get_etat() == ETATQCQ) {
796  if (source_shift.get_etat() == ETATZERO) {
797  source_shift.set_etat_qcq() ;
798  for (int i=0; i<3; i++) {
799  source_shift.set(i) = mtmp(i) ;
800  source_shift.set(i).va.coef_i() ;
801  }
802  }
803  else
804  for (int i=0; i<3; i++)
805  source_shift.set(i) += mtmp(i) ;
806  }
807 
808  source_shift.set_std_base() ;
809 
810  //----------------------------------------------
811  // Resolution of the Poisson equation for nuf
812  //----------------------------------------------
813 
814  source_nuf().poisson(par_poisson_nuf, nuf.set()) ;
815 
816  if (relativistic) {
817 
818  //----------------------------------------------
819  // Resolution of the Poisson equation for nuq
820  //----------------------------------------------
821 
822  source_nuq().poisson(par_poisson_nuq, nuq.set()) ;
823 
824  //---------------------------------------------------------
825  // Resolution of the vector Poisson equation for the shift
826  //---------------------------------------------------------
827 
828 
829  if (source_shift.get_etat() != ETATZERO) {
830 
831  for (int i=0; i<3; i++) {
832  if(source_shift(i).dz_nonzero()) {
833  assert( source_shift(i).get_dzpuis() == 4 ) ;
834  }
835  else{
836  (source_shift.set(i)).set_dzpuis(4) ;
837  }
838  }
839 
840  }
841  //##
842  // source_shift.dec2_dzpuis() ; // dzpuis 4 -> 2
843 
844  double lambda_shift = double(1) / double(3) ;
845 
846  if ( mg->get_np(0) == 1 ) {
847  lambda_shift = 0 ;
848  }
849 
850  source_shift.poisson_vect(lambda_shift, par_poisson_vect,
851  shift, w_shift, khi_shift) ;
852 
853  // Computation of tnphi and nphi from the Cartesian components
854  // of the shift
855  // -----------------------------------------------------------
856 
857  fait_nphi() ;
858 
859  }
860 
861  //-----------------------------------------
862  // Determination of the fluid velociy U
863  //-----------------------------------------
864 
865  if (mer > mer_fix_omega + delta_mer_kep) {
866 
867  omega *= fact_omega ; // Increase of the angular velocity if
868  } // fact_omega != 1
869 
870  bool omega_trop_grand = false ;
871  bool kepler = true ;
872 
873  while ( kepler ) {
874 
875  // Possible decrease of Omega to ensure a velocity < c
876 
877  bool superlum = true ;
878 
879  while ( superlum ) {
880 
881  // New fluid velocity U :
882 
883  Cmp tmp = omega - nphi() ;
884  tmp.annule(nzm1) ;
885  tmp.std_base_scal() ;
886 
887  tmp.mult_rsint() ; // Multiplication by r sin(theta)
888 
889  uuu = bbb() / nnn() * tmp ;
890 
891  if (uuu.get_etat() == ETATQCQ) {
892  // Same basis as (Omega -N^phi) r sin(theta) :
893  ((uuu.set()).va).set_base( (tmp.va).base ) ;
894  }
895 
896  // Is the new velocity larger than c in the equatorial plane ?
897 
898  superlum = false ;
899 
900  for (int l=0; l<nzet; l++) {
901  for (int i=0; i<mg->get_nr(l); i++) {
902 
903  double u1 = uuu()(l, 0, j_b, i) ;
904  if (u1 >= 1.) { // superluminal velocity
905  superlum = true ;
906  cout << "U > c for l, i : " << l << " " << i
907  << " U = " << u1 << endl ;
908  }
909  }
910  }
911  if ( superlum ) {
912  cout << "**** VELOCITY OF LIGHT REACHED ****" << endl ;
913  omega /= fact_omega ; // Decrease of Omega
914  cout << "New rotation frequency : "
915  << omega/(2.*M_PI) * f_unit << " Hz" << endl ;
916  omega_trop_grand = true ;
917  }
918  } // end of while ( superlum )
919 
920 
921  // New computation of U (which this time is not superluminal)
922  // as well as of gam_euler, ener_euler, etc...
923  // -----------------------------------
924 
925  hydro_euler() ;
926 
927  //------------------------------------------------------
928  // First integral of motion
929  //------------------------------------------------------
930 
931  // Centrifugal potential :
932  if (relativistic) {
933  mlngamma = - log( gam_euler ) ;
934  }
935  else {
936  mlngamma = - 0.5 * uuu*uuu ;
937  }
938 
939  Tenseur mag(mp) ;
940  if (is_conduct()) {
941  mag = mu0*M_j(A_phi, a_j) ;}
942  else{
943  mag = mu0*M_j(omega*A_phi-A_t, a_j) ;}
944 
945  // Equatorial values of various potentials :
946  double nuf_b = nuf()(l_b, k_b, j_b, i_b) ;
947  double nuq_b = nuq()(l_b, k_b, j_b, i_b) ;
948  double mlngamma_b = mlngamma()(l_b, k_b, j_b, i_b) ;
949  double mag_b = mag()(l_b, k_b, j_b, i_b) ;
950 
951  // Central values of various potentials :
952  double nuf_c = nuf()(0,0,0,0) ;
953  double nuq_c = nuq()(0,0,0,0) ;
954  double mlngamma_c = 0 ;
955  double mag_c = mag()(0,0,0,0) ;
956 
957  // Scale factor to ensure that the enthalpy is equal to ent_b at
958  // the equator
959  double alpha_r2 = ( ent_c - ent_b + mlngamma_c - mlngamma_b
960  + nuq_c - nuq_b + mag_c - mag_b)
961  / ( nuf_b - nuf_c ) ;
962  alpha_r = sqrt(alpha_r2) ;
963  cout << "alpha_r = " << alpha_r << endl ;
964 
965  // Readjustment of nu :
966  // -------------------
967 
968  logn = alpha_r2 * nuf + nuq ;
969  double nu_c = logn()(0,0,0,0) ;
970 
971  // First integral --> enthalpy in all space
972  //-----------------
973  ent = (ent_c + nu_c + mlngamma_c + mag_c) - logn - mlngamma - mag ;
974 
975  // Test: is the enthalpy negative somewhere in the equatorial plane
976  // inside the star ? If yes, this means that the Keplerian velocity
977  // has been overstep.
978 
979  kepler = false ;
980  for (int l=0; l<nzet; l++) {
981  int imax = mg->get_nr(l) - 1 ;
982  if (l == l_b) imax-- ; // The surface point is skipped
983  for (int i=0; i<imax; i++) {
984  if ( ent()(l, 0, j_b, i) < 0. ) {
985  kepler = true ;
986  cout << "ent < 0 for l, i : " << l << " " << i
987  << " ent = " << ent()(l, 0, j_b, i) << endl ;
988  }
989  }
990  }
991 
992  if ( kepler ) {
993  cout << "**** KEPLERIAN VELOCITY REACHED ****" << endl ;
994  omega /= fact_omega ; // Omega is decreased
995  cout << "New rotation frequency : "
996  << omega/(2.*M_PI) * f_unit << " Hz" << endl ;
997  omega_trop_grand = true ;
998  }
999 
1000  } // End of while ( kepler )
1001 
1002  if ( omega_trop_grand ) { // fact_omega is decreased for the
1003  // next step
1004  fact_omega = sqrt( fact_omega ) ;
1005  cout << "**** New fact_omega : " << fact_omega << endl ;
1006  }
1007 
1008  //----------------------------------------------------
1009  // Adaptation of the mapping to the new enthalpy field
1010  //----------------------------------------------------
1011 
1012  // Shall the adaptation be performed (cusp) ?
1013  // ------------------------------------------
1014 
1015  double dent_eq = ent().dsdr()(l_b, k_b, j_b, i_b) ;
1016  double dent_pole = ent().dsdr()(l_b, k_b, 0, i_b) ;
1017  double rap_dent = fabs( dent_eq / dent_pole ) ;
1018  cout << "| dH/dr_eq / dH/dr_pole | = " << rap_dent << endl ;
1019 
1020  if ( rap_dent < thres_adapt ) {
1021  adapt_flag = 0 ; // No adaptation of the mapping
1022  cout << "******* FROZEN MAPPING *********" << endl ;
1023  }
1024  else{
1025  adapt_flag = 1 ; // The adaptation of the mapping is to be
1026  // performed
1027  }
1028 
1029  mp_prev = mp_et ;
1030 
1031  mp.adapt(ent(), par_adapt) ;
1032 
1033  //----------------------------------------------------
1034  // Computation of the enthalpy at the new grid points
1035  //----------------------------------------------------
1036 
1037  mp_prev.homothetie(alpha_r) ;
1038 
1039  mp.reevaluate(&mp_prev, nzet+1, ent.set()) ;
1040 
1041  //----------------------------------------------------
1042  // Equation of state
1043  //----------------------------------------------------
1044 
1045  equation_of_state() ; // computes new values for nbar (n), ener (e)
1046  // and press (p) from the new ent (H)
1047 
1048  //---------------------------------------------------------
1049  // Matter source terms in the gravitational field equations
1050  //---------------------------------------------------------
1051 
1052  //## Computation of tnphi and nphi from the Cartesian components
1053  // of the shift for the test in hydro_euler():
1054 
1055  fait_nphi() ;
1056 
1057  hydro_euler() ; // computes new values for ener_euler (E),
1058  // s_euler (S) and u_euler (U^i)
1059 
1060  if (relativistic) {
1061 
1062  //-------------------------------------------------------
1063  // 2-D Poisson equation for tggg
1064  //-------------------------------------------------------
1065 
1066  mp.poisson2d(source_tggg(), mp.cmp_zero(), par_poisson_tggg, tggg.set()) ;
1067 
1068  //-------------------------------------------------------
1069  // 2-D Poisson equation for dzeta
1070  //-------------------------------------------------------
1071 
1072  mp.poisson2d(source_dzf(), source_dzq(), par_poisson_dzeta, dzeta.set()) ;
1073 
1074  err_grv2 = lbda_grv2 - 1;
1075  cout << "GRV2: " << err_grv2 << endl ;
1076 
1077  }
1078  else {
1079  err_grv2 = grv2() ;
1080  }
1081 
1082  //---------------------------------------
1083  // Computation of the metric coefficients (except for N^phi)
1084  //---------------------------------------
1085 
1086  // Relaxations on nu and dzeta :
1087 
1088  if (mer >= 10) {
1089  logn = relax * logn + relax_prev * logn_prev ;
1090 
1091  dzeta = relax * dzeta + relax_prev * dzeta_prev ;
1092  }
1093 
1094  // Update of the metric coefficients N, A, B and computation of K_ij :
1095 
1096  update_metric() ;
1097 
1098  //-----------------------
1099  // Informations display
1100  //-----------------------
1101 
1102  // partial_display(cout) ;
1103  fichfreq << " " << omega / (2*M_PI) * f_unit ;
1104  fichevol << " " << rap_dent ;
1105  fichevol << " " << ray_pole() / ray_eq() ;
1106  fichevol << " " << ent_c ;
1107 
1108  //-----------------------------------------
1109  // Convergence towards a given baryon mass
1110  //-----------------------------------------
1111 
1112  if (mer > mer_mass) {
1113 
1114  double xx ;
1115  if (mbar_wanted > 0.) {
1116  xx = mass_b() / mbar_wanted - 1. ;
1117  cout << "Discrep. baryon mass <-> wanted bar. mass : " << xx
1118  << endl ;
1119  }
1120  else{
1121  xx = mass_g() / fabs(mbar_wanted) - 1. ;
1122  cout << "Discrep. grav. mass <-> wanted grav. mass : " << xx
1123  << endl ;
1124  }
1125  double xprog = ( mer > 2*mer_mass) ? 1. :
1126  double(mer-mer_mass)/double(mer_mass) ;
1127  xx *= xprog ;
1128  double ax = .5 * ( 2. + xx ) / (1. + xx ) ;
1129  double fact = pow(ax, aexp_mass) ;
1130  cout << " xprog, xx, ax, fact : " << xprog << " " <<
1131  xx << " " << ax << " " << fact << endl ;
1132 
1133  if ( change_ent ) {
1134  ent_c *= fact ;
1135  }
1136  else {
1137  if (mer%4 == 0) omega *= fact ;
1138  }
1139  }
1140 
1141  //---------------------------------------------
1142  // Convergence towards a given magnetic moment
1143  //---------------------------------------------
1144 
1145  if (mer > mer_magmom) {
1146 
1147  // if(mer > mer_mass) {
1148  // cerr << "et_magnetisation::equilibrium()" << endl ;
1149  // cerr << "Impossible to converge to a given baryon mass" << endl ;
1150  // cerr << "and a given magnetic moment!" << endl ;
1151  // cerr << "Aborting..." << endl ;
1152  // abort() ;
1153  // }
1154  double xx = MagMom() / magmom_wanted - 1. ;
1155  cout << "Discrep. mag. moment <-> wanted mag. moment : " << xx
1156  << endl ;
1157 
1158  double xprog = ( mer > 2*mer_magmom) ? 1. :
1159  double(mer-mer_magmom)/double(mer_magmom) ;
1160  xx *= xprog ;
1161  double ax = .5 * ( 2. + xx ) / (1. + xx ) ;
1162  double fact = pow(ax, aexp_mass) ;
1163  cout << " xprog, xx, ax, fact : " << xprog << " " <<
1164  xx << " " << ax << " " << fact << endl ;
1165 
1166  a_j *= fact ;
1167  }
1168 
1169  //------------------------------------------------------------
1170  // Relative change in enthalpy with respect to previous step
1171  //------------------------------------------------------------
1172 
1173  Tbl diff_ent_tbl = diffrel( ent(), ent_prev() ) ;
1174  diff_ent = diff_ent_tbl(0) ;
1175  for (int l=1; l<nzet; l++) {
1176  diff_ent += diff_ent_tbl(l) ;
1177  }
1178  diff_ent /= nzet ;
1179 
1180  fichconv << " " << log10( fabs(diff_ent) + 1.e-16 ) ;
1181  fichconv << " " << log10( fabs(err_grv2) + 1.e-16 ) ;
1182 
1183  //------------------------------
1184  // Recycling for the next step
1185  //------------------------------
1186 
1187  ent_prev = ent ;
1188  logn_prev = logn ;
1189  dzeta_prev = dzeta ;
1190 
1191  fichconv << endl ;
1192  fichfreq << endl ;
1193  fichevol << endl ;
1194  fichconv.flush() ;
1195  fichfreq.flush() ;
1196  fichevol.flush() ;
1197 
1198  } // End of main loop
1199 
1200  //=========================================================================
1201  // End of iteration
1202  //=========================================================================
1203 
1204  fichconv.close() ;
1205  fichfreq.close() ;
1206  fichevol.close() ;
1207 }
1208 
1209 
1210 
1211 
1212 
1213 
1214 
1215 
1216 
1217 
1218 
1219 
1220 
1221 
1222 }
void annule_domain(int l)
Sets the Tensor to zero in a given domain.
Definition: tensor.C:675
Tbl maxabs(const Tensor &aa, const char *comment=0x0, ostream &ost=cout, bool verb=true)
Maxima in each domain of the absolute values of the tensor components.
const Grille3d * get_grille3d(int l) const
Returns a pointer on the 3D mono-grid for domain no. l.
Definition: grilles.h:517
Tenseur khi_shift
Scalar used in the decomposition of shift , following Shibata&#39;s prescription [Prog.
Definition: etoile.h:1563
Cmp log(const Cmp &)
Neperian logarithm.
Definition: cmp_math.C:299
virtual void sauve(FILE *) const
Save in a binary file.
Definition: tensor_sym.C:375
void add_tenseur_mod(Tenseur &ti, int position=0)
Adds the address of a new modifiable Tenseur to the list.
Definition: param.C:1145
virtual double nbar_ent_p(double ent, const Param *par=0x0) const
Computes the baryon density from the log-enthalpy.
Definition: eos_mag.C:333
Tenseur E_em
electromagnetic energy density in the Eulerian frame
Definition: et_rot_mag.h:161
Cmp ssjm1_khi
Effective source at the previous step for the resolution of the Poisson equation for the scalar by m...
Definition: etoile.h:1619
Component of a tensorial field *** DEPRECATED : use class Scalar instead ***.
Definition: cmp.h:446
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
const Tenseur & gradient_spher() const
Returns the gradient of *this (Spherical coordinates) (scalar field only).
Definition: tenseur.C:1564
const Cmp & cmp_zero() const
Returns the null Cmp defined on *this.
Definition: map.h:825
int get_np(int l) const
Returns the number of points in the azimuthal direction ( ) in domain no. l.
Definition: grilles.h:479
Cmp sqrt(const Cmp &)
Square root.
Definition: cmp_math.C:223
virtual void set_etat_zero()
Sets the logical state to ETATZERO (zero).
Definition: scalar.C:330
void annule(int l)
Sets the Cmp to zero in a given domain.
Definition: cmp.C:351
Et_magnetisation(Map &mp_i, int nzet_i, bool relat, const Eos &eos_i, bool include_mag=true, bool use_B=true)
Standard constructor.
void set_std_base()
Set the standard spectal basis of decomposition for each component.
Definition: tenseur.C:1186
Multi-domain array.
Definition: mtbl.h:118
void fait_nphi()
Computes tnphi and nphi from the Cartesian components of the shift, stored in shift ...
Definition: etoile_rot.C:803
const Base_vect_spher & get_bvect_spher() const
Returns the orthonormal vectorial basis associated with the coordinates of the mapping.
Definition: map.h:801
double Q
In the case of a perfect conductor, the requated baryonic charge.
Definition: et_rot_mag.h:179
Lorene prototypes.
Definition: app_hor.h:67
Tenseur w_shift
Vector used in the decomposition of shift , following Shibata&#39;s prescription [Prog.
Definition: etoile.h:1553
Scalar E_I
Interaction (magnetisation) energy density.
Definition: et_rot_mag.h:624
Equation of state base class.
Definition: eos.h:209
Tenseur ssjm1_wshift
Effective source at the previous step for the resolution of the vector Poisson equation for ...
Definition: etoile.h:1628
const Mg3d * get_mg() const
Gives the Mg3d on which the mapping is defined.
Definition: map.h:783
int get_etat() const
Returns the logical state.
Definition: cmp.h:899
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
Tensor field of valence 0 (or component of a tensorial field).
Definition: scalar.h:393
Tenseur nphi
Metric coefficient .
Definition: etoile.h:1513
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].
Base class for coordinate mappings.
Definition: map.h:688
Tenseur flat_scalar_prod(const Tenseur &t1, const Tenseur &t2)
Scalar product of two Tenseur when the metric is : performs the contraction of the last index of t1 w...
Tenseur b_car
Square of the metric factor B.
Definition: etoile.h:1510
virtual void adapt(const Cmp &ent, const Param &par, int nbr=0)=0
Adaptation of the mapping to a given scalar field.
Tenseur Magn() const
Computes the magnetic field spherical components in Lorene&#39;s units.
int get_type_t() const
Returns the type of sampling in the direction: SYM : : symmetry with respect to the equatorial pl...
Definition: grilles.h:502
Basic integer array class.
Definition: itbl.h:122
virtual void exponential_filter_r(int lzmin, int lzmax, int p, double alpha=-16.)
Applies an exponential filter to the spectral coefficients in the radial direction.
virtual void std_spectral_base()
Sets the spectral bases of the Valeur va to the standard ones for a scalar field. ...
Definition: scalar.C:790
virtual void sauve(FILE *) const
Save in a binary file.
Definition: tensor.C:915
int get_etat() const
Gives the logical state.
Definition: tbl.h:414
virtual void allocate_all()
Sets the logical state to ETATQCQ (ordinary state) and performs the memory allocation of all the elem...
Definition: scalar.C:371
Tenseur press
Fluid pressure.
Definition: etoile.h:464
double * x
Array of values of at the nr collocation points.
Definition: grilles.h:215
Tensor field of valence 1.
Definition: vector.h:188
Tenseur shift
Total shift vector.
Definition: etoile.h:515
Tbl & set_domain(int l)
Read/write of the value in a given domain.
Definition: scalar.h:621
Tbl diffrel(const Cmp &a, const Cmp &b)
Relative difference between two Cmp (norme version).
Definition: cmp_math.C:507
void add_double_mod(double &x, int position=0)
Adds the address of a new modifiable double to the list.
Definition: param.C:456
Cmp A_phi
-component of the electromagnetic potential 1-form divided by .
Definition: et_rot_mag.h:155
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition: tbl.C:364
void update_metric()
Computes metric coefficients from known potentials.
Definition: et_rot_upmetr.C:72
virtual void sauve(FILE *) const
Save in a file.
Definition: scalar.C:692
void set_etat_zero()
Sets the logical state to ETATZERO (zero).
Definition: cmp.C:292
Cmp & set()
Read/write for a scalar (see also operator=(const Cmp&) ).
Definition: tenseur.C:840
Tenseur u_euler
Fluid 3-velocity with respect to the Eulerian observer.
Definition: etoile.h:477
Cmp ssjm1_nuq
Effective source at the previous step for the resolution of the Poisson equation for nuq by means of ...
Definition: etoile.h:1601
void operator=(const Et_rot_mag &)
Assignment to another Et_rot_mag.
Definition: et_rot_mag.C:286
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
virtual double press_ent_p(double ent, const Param *par=0x0) const
Computes the pressure from the log-enthalpy.
Definition: eos_mag.C:422
double a_j
Amplitude of the curent/charge function.
Definition: et_rot_mag.h:180
virtual void poisson2d(const Cmp &source_mat, const Cmp &source_quad, Param &par, Cmp &uu) const =0
Computes the solution of a 2-D Poisson equation.
bool is_conduct() const
Tells if the star is made of conducting or isolating material.
Definition: et_rot_mag.h:241
Vector J_I
Interaction momentum density 3-vector.
Definition: et_rot_mag.h:627
Tenseur nbar
Baryon density in the fluid frame.
Definition: etoile.h:462
Scalar xmag
The magnetisation scalar.
Definition: et_rot_mag.h:622
bool use_B_in_eos
Flag : true if the value of the magnetic field is used in the Eos.
Definition: et_rot_mag.h:617
Class for magnetized (isolator or perfect conductor), rigidly rotating stars.
Definition: et_rot_mag.h:143
double mag_ent_p(double ent, const Param *par=0x0) const
Computes the magnetisation.
Definition: eos_mag.C:463
virtual ostream & operator>>(ostream &) const
Operator >> (virtual function called by the operator <<).
Definition: et_rot_mag.C:339
Tenseur gam_euler
Lorentz factor between the fluid and Eulerian observers.
Definition: etoile.h:474
Mtbl * c
Values of the function at the points of the multi-grid.
Definition: valeur.h:309
bool include_magnetisation
Flag : true if magnetisation terms are included in the equations.
Definition: et_rot_mag.h:620
virtual double mass_b() const
Baryon mass.
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
void mult_rsint()
Multiplication by .
Definition: cmp_r_manip.C:119
Map & mp
Mapping associated with the star.
Definition: etoile.h:432
int get_nzone() const
Returns the number of domains.
Definition: grilles.h:465
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
Cmp ssjm1_nuf
Effective source at the previous step for the resolution of the Poisson equation for nuf by means of ...
Definition: etoile.h:1595
Tenseur ak_car
Scalar .
Definition: etoile.h:1589
virtual void sauve(FILE *) const
Save in a file.
int get_etat() const
Returns the logical state.
Definition: tenseur.h:710
virtual void del_deriv() const
Deletes all the derived quantities.
Definition: et_rot_mag.C:261
void operator=(const Et_magnetisation &)
Assignment to another Et_rot_mag.
Tenseur bbb
Metric factor B.
Definition: etoile.h:1507
virtual void hydro_euler()
Computes the hydrodynamical quantities relative to the Eulerian observer from those in the fluid fram...
Definition: et_rot_hydro.C:86
Cmp pow(const Cmp &, int)
Power .
Definition: cmp_math.C:351
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition: cmp.C:307
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition: mtbl.C:302
Class for a magnetized (tabulated) equation of state.
Definition: eos_mag.h:81
Tenseur uuu
Norm of u_euler.
Definition: etoile.h:1521
Tbl & set(int l)
Read/write of the Tbl in a given domain.
Definition: mtbl.h:221
Tenseur tggg
Metric potential .
Definition: etoile.h:1540
double omega
Rotation angular velocity ([f_unit] )
Definition: etoile.h:1504
virtual ostream & operator>>(ostream &) const
Operator >> (virtual function called by the operator <<).
void set_etat_zero()
Sets the logical state to ETATZERO (zero).
Definition: tbl.C:350
Tenseur nuf
Part of the Metric potential = logn generated by the matter terms.
Definition: etoile.h:1529
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
virtual double grv2() const
Error on the virial identity GRV2.
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
void allocate_all()
Sets the logical state to ETATQCQ (ordinary state) and performs the memory allocation of all the elem...
Definition: cmp.C:326
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].
Tbl & set(int l)
Read/write of the value in a given domain.
Definition: cmp.h:724
Tenseur ener
Total energy density in the fluid frame.
Definition: etoile.h:463
double & set_grid_point(int l, int k, int j, int i)
Setting the value of the field at a given grid point.
Definition: scalar.h:690
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
void equilibrium_mag(double ent_c, double omega0, double fact_omega, int nzadapt, const Tbl &ent_limit, const Itbl &icontrol, const Tbl &control, double mbar_wanted, double magmom_wanted, double aexp_mass, Tbl &diff, double Q0, double a_j0, Cmp(*f_j)(const Cmp &x, const double), Cmp(*M_j)(const Cmp &x, const double))
Computes an equilibrium configuration.
Cmp log10(const Cmp &)
Basis 10 logarithm.
Definition: cmp_math.C:325
Tenseur ent
Log-enthalpy (relativistic case) or specific enthalpy (Newtonian case)
Definition: etoile.h:460
virtual void sauve(FILE *) const
Save in a file.
Definition: et_rot_mag.C:314
virtual double ener_ent_p(double ent, const Param *par=0x0) const
Computes the total energy density from the log-enthalpy.
Definition: eos_mag.C:376
Cmp abs(const Cmp &)
Absolute value.
Definition: cmp_math.C:413
virtual void magnet_comput(const int adapt_flag, Cmp(*f_j)(const Cmp &x, const double), Param &par_poisson_At, Param &par_poisson_Avect)
Computes the electromagnetic quantities solving the Maxwell equations (6) and (7) of [Bocquet...
Sym_tensor Sij_I
Interaction stress 3-tensor.
Definition: et_rot_mag.h:630
virtual void reevaluate(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.
void add_double(const double &x, int position=0)
Adds the the address of a new double to the list.
Definition: param.C:318
virtual double mass_g() const
Gravitational mass.
Cmp A_t
t-component of the elecctromagnetic potential 1-form, divided by .
Definition: et_rot_mag.h:150
Tenseur & logn
Metric potential = logn_auto.
Definition: etoile.h:1524
void set_etat_c_qcq()
Sets the logical state to ETATQCQ (ordinary state) for values in the configuration space (Mtbl c )...
Definition: valeur.C:704
virtual void set_etat_zero()
Sets the logical state of all components to ETATZERO (zero state).
Definition: tensor.C:506
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
virtual void equation_of_state()
Computes the proper baryon and energy density, as well as pressure from the enthalpy.
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
virtual ~Et_magnetisation()
Destructor.
double MagMom() const
Magnetic Momentum in SI units.
Tenseur nuq
Part of the Metric potential = logn generated by the quadratic terms.
Definition: etoile.h:1534
Tbl ** t
Array (size nzone ) of pointers on the Tbl &#39;s.
Definition: mtbl.h:132
Tenseur ener_euler
Total energy density in the Eulerian frame.
Definition: etoile.h:468
void set_etat_zero()
Sets the logical state to ETATZERO (zero state).
Definition: tenseur.C:661
void annule_hard()
Sets the Tbl to zero in a hard way.
Definition: tbl.C:375
Tenseur & dzeta
Metric potential = beta_auto.
Definition: etoile.h:1537
Valeur va
The numerical value of the Cmp.
Definition: cmp.h:464
Cmp ssjm1_tggg
Effective source at the previous step for the resolution of the Poisson equation for tggg ...
Definition: etoile.h:1611
Tenseur_sym tkij
Tensor related to the extrinsic curvature tensor by .
Definition: etoile.h:1570
Class intended to describe valence-2 symmetric tensors.
Definition: sym_tensor.h:226
Standard electro-magnetic units.
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
Tenseur Jp_em
component of the electromagnetic momentum density 3-vector, as measured in the Eulerian frame...
Definition: et_rot_mag.h:167
virtual void MHD_comput()
Computes the electromagnetic part of the stress-energy tensor.