LORENE
et_rot_mag_global.C
1 /*
2  * Methods for computing global quantities within the class Etoile_rot
3  *
4  * (see file etoile.h for documentation)
5  */
6 
7 /*
8  * Copyright (c) 2000-2001 Eric Gourgoulhon
9  * Copyright (c) 2002 Emmanuel Marcq
10  * Copyright (c) 2002 Jerome Novak
11  *
12  * This file is part of LORENE.
13  *
14  * LORENE is free software; you can redistribute it and/or modify
15  * it under the terms of the GNU General Public License as published by
16  * the Free Software Foundation; either version 2 of the License, or
17  * (at your option) any later version.
18  *
19  * LORENE is distributed in the hope that it will be useful,
20  * but WITHOUT ANY WARRANTY; without even the implied warranty of
21  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
22  * GNU General Public License for more details.
23  *
24  * You should have received a copy of the GNU General Public License
25  * along with LORENE; if not, write to the Free Software
26  * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
27  *
28  */
29 
30 
31 char et_rot_mag_global_C[] = "$Header: /cvsroot/Lorene/C++/Source/Etoile/et_rot_mag_global.C,v 1.22 2015/06/12 12:38:25 j_novak Exp $" ;
32 
33 /*
34  * $Id: et_rot_mag_global.C,v 1.22 2015/06/12 12:38:25 j_novak Exp $
35  * $Log: et_rot_mag_global.C,v $
36  * Revision 1.22 2015/06/12 12:38:25 j_novak
37  * Implementation of the corrected formula for the quadrupole momentum.
38  *
39  * Revision 1.21 2014/10/13 08:52:58 j_novak
40  * Lorene classes and functions now belong to the namespace Lorene.
41  *
42  * Revision 1.20 2014/05/13 10:06:13 j_novak
43  * Change of magnetic units, to make the Lorene unit system coherent. Magnetic field is now expressed in Lorene units. Improvement on the comments on units.
44  *
45  * Revision 1.19 2012/08/12 17:48:35 p_cerda
46  * Magnetstar: New classes for magnetstar. Allowing for non-equatorial symmetry in Etoile et al. Adding B_phi in Et_rot_mag.
47  *
48  * Revision 1.18 2006/01/31 15:54:57 j_novak
49  * Corrected a missing '-' sign for the theta component of the magnetic field in
50  * Et_rot_mag::Magn(). This had no influence in the calculations, only in the
51  * display of B values.
52  *
53  * Revision 1.17 2004/03/25 10:29:06 j_novak
54  * All LORENE's units are now defined in the namespace Unites (in file unites.h).
55  *
56  * Revision 1.16 2003/10/27 10:52:19 e_gourgoulhon
57  * Suppressed the global #include "unites.h"
58  * and made it local to each function.
59  *
60  * Revision 1.15 2002/10/17 11:30:54 j_novak
61  * Corrected mistake in angu_mom()
62  *
63  * Revision 1.14 2002/06/03 13:00:45 e_marcq
64  *
65  * conduc parameter read in parmag.d
66  *
67  * Revision 1.12 2002/05/22 12:20:17 j_novak
68  * *** empty log message ***
69  *
70  * Revision 1.11 2002/05/20 15:44:55 e_marcq
71  *
72  * Dimension errors corrected, parmag.d input file created and read
73  *
74  * Revision 1.10 2002/05/20 10:31:59 j_novak
75  * *** empty log message ***
76  *
77  * Revision 1.9 2002/05/20 08:27:59 j_novak
78  * *** empty log message ***
79  *
80  * Revision 1.8 2002/05/17 15:08:01 e_marcq
81  *
82  * Rotation progressive plug-in, units corrected, Q and a_j new member data
83  *
84  * Revision 1.7 2002/05/16 13:27:11 j_novak
85  * *** empty log message ***
86  *
87  * Revision 1.6 2002/05/16 10:02:09 j_novak
88  * Errors in stress energy tensor corrected
89  *
90  * Revision 1.5 2002/05/15 09:53:59 j_novak
91  * First operational version
92  *
93  * Revision 1.4 2002/05/14 13:45:30 e_marcq
94  *
95  * Correction de la formule du rapport gyromagnetique
96  *
97  * Revision 1.1 2002/05/10 09:26:52 j_novak
98  * Added new class Et_rot_mag for magnetized rotating neutron stars (under development)
99  *
100  *
101  * $Header: /cvsroot/Lorene/C++/Source/Etoile/et_rot_mag_global.C,v 1.22 2015/06/12 12:38:25 j_novak Exp $
102  *
103  */
104 
105 // Headers C
106 #include <cstdlib>
107 #include <cmath>
108 
109 // Headers Lorene
110 #include "et_rot_mag.h"
111 #include "unites.h"
112 
113 // Definition des fonctions membres differentes ou nouvelles
114 
115 namespace Lorene {
117  // Calcule les grandeurs du tenseur impulsion-energie EM a partir des champs
118 
119  using namespace Unites_mag ;
120 
121  Tenseur ATTENS(A_t) ;
122 
123  Tenseur APTENS(A_phi) ;
124 
126  APTENS.gradient_spher())() );
128  ATTENS.gradient_spher())() );
130  ATTENS.gradient_spher())() );
131 
132  if (ApAp.get_etat() != ETATZERO) {
133  ApAp.set().div_rsint() ;
134  ApAp.set().div_rsint() ;
135  }
136  if (ApAt.get_etat() != ETATZERO)
137  ApAt.set().div_rsint() ;
138 
139  E_em = 0.5*mu0 * ( 1/(a_car*nnn*nnn) * (AtAt + 2*tnphi*ApAt)
140  + ( (tnphi*tnphi/(a_car*nnn*nnn)) + 1/(a_car*b_car) )*ApAp );
141  Jp_em = -mu0 * (ApAt + tnphi*ApAp) /(a_car*nnn) ;
142  if (Jp_em.get_etat() != ETATZERO) Jp_em.set().mult_rsint() ;
143  Srr_em = 0 ;
144  // Stt_em = -Srr_em
145  Spp_em = E_em ;
146 }
147 
149 
150  using namespace Unites_mag ;
151 
152  Cmp E_r(mp); Cmp E_t(mp);
153  E_r = 1/(sqrt(a_car())*nnn())*(A_t.dsdr()+nphi()*A_phi.dsdr()) ;
154  E_t = 1/(sqrt(a_car())*nnn())*(A_t.srdsdt()+nphi()*A_phi.srdsdt()) ;
155  E_r.va.set_base((A_t.dsdr()).va.base) ;
156  E_t.va.set_base((A_t.srdsdt()).va.base) ;
157  Tenseur Elect(mp, 1, CON, mp.get_bvect_spher()) ;
158  Elect.set_etat_qcq() ;
159  Elect.set(0) = E_r ;
160  Elect.set(1) = E_t ;
161  Elect.set(2) = 0. ;
162 
163  return Elect*mu0 ;
164 
165 }
166 
168 
169  using namespace Unites_mag ;
170 
171  Cmp B_r(mp); Cmp B_t(mp);
172  B_r = 1/(sqrt(a_car())*bbb())*A_phi.srdsdt();
173  B_r.va.set_base((A_phi.srdsdt()).va.base) ;
174  B_r.div_rsint();
175  B_t = 1/(sqrt(a_car())*bbb())*A_phi.dsdr();
176  B_t.va.set_base((A_phi.dsdr()).va.base) ;
177  B_t.div_rsint();
178 
179  Tenseur Bmag(mp, 1, CON, mp.get_bvect_spher()) ;
180  Bmag.set_etat_qcq() ;
181  Bmag.set(0) = B_r ;
182  Bmag.set(1) = -B_t ;
183  Bmag.set(2) = B_phi ;
184 
185  return Bmag*mu0 ;
186 
187 }
188 
189 double Et_rot_mag::MagMom() const {
190 
191  using namespace Unites_mag ;
192 
193  int Z = mp.get_mg()->get_nzone();
194  double mm ;
195 
196  if(A_phi.get_etat()==ETATZERO) {
197 
198  mm = 0 ;
199  }else{
200 
201  Valeur** asymp = A_phi.asymptot(1) ;
202  mm = 4*M_PI*(*asymp[1])(Z-1,0,mp.get_mg()->get_nt(Z-1)-1,0) ;
203 
204  delete asymp[0] ;
205  delete asymp[1] ;
206 
207  delete [] asymp ;
208  }
209 
210  return mm*j_unit*pow(r_unit,4) ;
211 
212 }
213 
214 double Et_rot_mag::Q_comput() const {
215 
216  using namespace Unites_mag ;
217 
218  int Z = mp.get_mg()->get_nzone();
219 
220  if(A_t.get_etat()==ETATZERO) {
221  return 0 ;
222  }else{
223  Valeur** asymp = A_t.asymptot(1) ;
224 
225  double Q_c = -4*M_PI*(*asymp[1])(Z-1,0,0,0) ;
226  delete asymp[0] ;
227  delete asymp[1] ;
228 
229  delete [] asymp ;
230 
231  return Q_c *(j_unit/v_unit*pow(r_unit,3)) ;}
232  }
233 
234 
235 double Et_rot_mag::Q_int() const {
236 
237  using namespace Unites_mag ;
238 
239  double Qi = 0. ;
240 
241  if (relativistic) {
242 
243  Cmp dens = a_car() * bbb() * nnn() * j_t ;
244 
245  dens.std_base_scal() ;
246 
247  Qi = dens.integrale() ;
248 
249 
250  }
251  else{ // Newtonian case
252  assert(nbar.get_etat() == ETATQCQ) ;
253 
254  Qi = ( j_t.integrale() ) ;
255 
256  }
257 
258 
259 
260  return Qi * (j_unit/v_unit*pow(r_unit,3)) ;
261 
262 }
263 
264 
265 double Et_rot_mag::GyroMag() const {
266 
267  using namespace Unites_mag ;
268 
269  return 2*MagMom()*mass_g()/(Q_comput()*angu_mom()*v_unit*r_unit);
270 
271 }
272  //----------------------------//
273  // Gravitational mass //
274  //----------------------------//
275 
276 double Et_rot_mag::mass_g() const {
277 
278  if (p_mass_g == 0x0) { // a new computation is required
279 
280  if (relativistic) {
281 
282  Tenseur source = nnn * (ener_euler + E_em + s_euler + Spp_em) +
283  nphi * Jp_em + 2 * bbb * (ener_euler + press) * tnphi * uuu ;
284 
285  source = a_car * bbb * source ;
286 
287  source.set_std_base() ;
288 
289  p_mass_g = new double( source().integrale() ) ;
290 
291 
292  }
293  else{ // Newtonian case
294  p_mass_g = new double( mass_b() ) ; // in the Newtonian case
295  // M_g = M_b
296  }
297  }
298 
299  return *p_mass_g ;
300 
301 }
302 
303  //----------------------------//
304  // Angular momentum //
305  //----------------------------//
306 
307 double Et_rot_mag::angu_mom() const {
308 
309  if (p_angu_mom == 0x0) { // a new computation is required
310 
311  Cmp dens = uuu() ;
312 
313  dens.mult_r() ; // Multiplication by
314  dens.va = (dens.va).mult_st() ; // r sin(theta)
315 
316  if (relativistic) {
317  dens = a_car() * (b_car() * (ener_euler() + press())
318  * dens + bbb() * Jp_em()) ;
319  }
320  else { // Newtonian case
321  dens = nbar() * dens ;
322  }
323 
324  dens.std_base_scal() ;
325 
326  p_angu_mom = new double( dens.integrale() ) ;
327 
328  }
329 
330  return *p_angu_mom ;
331 
332 }
333 
334 
335  //----------------------------//
336  // T/W //
337  //----------------------------//
338 
339 // Redefini en virtual dans le .h : A CHANGER
340 
341 double Et_rot_mag::tsw() const {
342 
343  if (p_tsw == 0x0) { // a new computation is required
344 
345  double tcin = 0.5 * omega * angu_mom() ;
346 
347  if (relativistic) {
348 
349  Cmp dens = a_car() * bbb() * gam_euler() * ener() ;
350  dens.std_base_scal() ;
351  double mass_p = dens.integrale() ;
352 
353  p_tsw = new double( tcin / ( mass_p + tcin - mass_g() ) ) ;
354 
355  }
356  else { // Newtonian case
357  Cmp dens = 0.5 * nbar() * logn() ;
358  dens.std_base_scal() ;
359  double wgrav = dens.integrale() ;
360  p_tsw = new double( tcin / fabs(wgrav) ) ;
361  }
362 
363 
364  }
365 
366  return *p_tsw ;
367 
368 }
369 
370 
371  //----------------------------//
372  // GRV2 //
373  //----------------------------//
374 
375 double Et_rot_mag::grv2() const {
376 
377  if (p_grv2 == 0x0) { // a new computation is required
378 
379  // To get qpig:
380  using namespace Unites ;
381 
382  Tenseur sou_m = 2 * qpig * a_car * (press + (ener_euler+press)
383  * uuu*uuu ) ;
384 
385  Tenseur sou_q = 2 * qpig * a_car * Spp_em + 1.5 * ak_car
387 
388  p_grv2 = new double( double(1) - lambda_grv2(sou_m(), sou_q()) ) ;
389 
390  }
391 
392  return *p_grv2 ;
393 
394 }
395 
396 
397  //----------------------------//
398  // GRV3 //
399  //----------------------------//
400 
401 double Et_rot_mag::grv3(ostream* ost) const {
402 
403  if (p_grv3 == 0x0) { // a new computation is required
404 
405  // To get qpig:
406  using namespace Unites ;
407 
408  Tenseur source(mp) ;
409 
410  // Gravitational term [cf. Eq. (43) of Gourgoulhon & Bonazzola
411  // ------------------ Class. Quantum Grav. 11, 443 (1994)]
412 
413  if (relativistic) {
414  Tenseur alpha = dzeta - logn ;
415  Tenseur beta = log( bbb ) ;
416  beta.set_std_base() ;
417 
418  source = 0.75 * ak_car
419  - flat_scalar_prod(logn.gradient_spher(),
420  logn.gradient_spher() )
421  + 0.5 * flat_scalar_prod(alpha.gradient_spher(),
422  beta.gradient_spher() ) ;
423 
424  Cmp aa = alpha() - 0.5 * beta() ;
425  Cmp daadt = aa.srdsdt() ; // 1/r d/dth
426 
427  // What follows is valid only for a mapping of class Map_radial :
428  const Map_radial* mpr = dynamic_cast<const Map_radial*>(&mp) ;
429  if (mpr == 0x0) {
430  cout << "Etoile_rot::grv3: the mapping does not belong"
431  << " to the class Map_radial !" << endl ;
432  abort() ;
433  }
434 
435  // Computation of 1/tan(theta) * 1/r daa/dtheta
436  if (daadt.get_etat() == ETATQCQ) {
437  Valeur& vdaadt = daadt.va ;
438  vdaadt = vdaadt.ssint() ; // division by sin(theta)
439  vdaadt = vdaadt.mult_ct() ; // multiplication by cos(theta)
440  }
441 
442  Cmp temp = aa.dsdr() + daadt ;
443  temp = ( bbb() - a_car()/bbb() ) * temp ;
444  temp.std_base_scal() ;
445 
446  // Division by r
447  Valeur& vtemp = temp.va ;
448  vtemp = vtemp.sx() ; // division by xi in the nucleus
449  // Id in the shells
450  // division by xi-1 in the ZEC
451  vtemp = (mpr->xsr) * vtemp ; // multiplication by xi/r in the nucleus
452  // by 1/r in the shells
453  // by r(xi-1) in the ZEC
454 
455  // In the ZEC, a multiplication by r has been performed instead
456  // of the division:
457  temp.set_dzpuis( temp.get_dzpuis() + 2 ) ;
458 
459  source = bbb() * source() + 0.5 * temp ;
460 
461  }
462  else{
463  source = - 0.5 * flat_scalar_prod(logn.gradient_spher(),
464  logn.gradient_spher() ) ;
465  }
466 
467  source.set_std_base() ;
468 
469  double int_grav = source().integrale() ;
470 
471  // Matter term
472  // -----------
473 
474  if (relativistic) {
475  source = qpig * a_car * bbb * ( s_euler + Spp_em ) ;
476  }
477  else{
478  source = qpig * ( 3 * press + nbar * uuu * uuu ) ;
479  }
480 
481  source.set_std_base() ;
482 
483  double int_mat = source().integrale() ;
484 
485  // Virial error
486  // ------------
487  if (ost != 0x0) {
488  *ost << "Etoile_rot::grv3 : gravitational term : " << int_grav
489  << endl ;
490  *ost << "Etoile_rot::grv3 : matter term : " << int_mat
491  << endl ;
492  }
493 
494  p_grv3 = new double( (int_grav + int_mat) / int_mat ) ;
495 
496  }
497 
498  return *p_grv3 ;
499 
500 }
501 
502  //----------------------------//
503  // Quadrupole moment //
504  //----------------------------//
505 
506  double Et_rot_mag::mom_quad_old() const {
507 
508  if (p_mom_quad_old == 0x0) { // a new computation is required
509 
510  // To get qpig:
511  using namespace Unites ;
512 
513  // Source for of the Poisson equation for nu
514  // -----------------------------------------
515 
516  Tenseur source(mp) ;
517 
518  if (relativistic) {
519  Tenseur beta = log(bbb) ;
520  beta.set_std_base() ;
521  source = qpig * a_car *( ener_euler + s_euler + Spp_em )
523  logn.gradient_spher() + beta.gradient_spher()) ;
524  }
525  else {
526  source = qpig * nbar ;
527  }
528  source.set_std_base() ;
529 
530  // Multiplication by -r^2 P_2(cos(theta))
531  // [cf Eq.(7) of Salgado et al. Astron. Astrophys. 291, 155 (1994) ]
532  // ------------------------------------------------------------------
533 
534  // Multiplication by r^2 :
535  // ----------------------
536  Cmp& csource = source.set() ;
537  csource.mult_r() ;
538  csource.mult_r() ;
539  if (csource.check_dzpuis(2)) {
540  csource.inc2_dzpuis() ;
541  }
542 
543  // Muliplication by cos^2(theta) :
544  // -----------------------------
545  Cmp temp = csource ;
546 
547  // What follows is valid only for a mapping of class Map_radial :
548  assert( dynamic_cast<const Map_radial*>(&mp) != 0x0 ) ;
549 
550  if (temp.get_etat() == ETATQCQ) {
551  Valeur& vtemp = temp.va ;
552  vtemp = vtemp.mult_ct() ; // multiplication by cos(theta)
553  vtemp = vtemp.mult_ct() ; // multiplication by cos(theta)
554  }
555 
556  // Muliplication by -P_2(cos(theta)) :
557  // ----------------------------------
558  source = 0.5 * source() - 1.5 * temp ;
559 
560  // Final result
561  // ------------
562 
563  p_mom_quad_old = new double( source().integrale() / qpig ) ;
564 
565  }
566 
567  return *p_mom_quad_old ;
568 
569  }
570 
571 }
const Cmp & dsdr() const
Returns of *this .
Definition: cmp_deriv.C:84
double * p_mom_quad_old
Part of the quadrupole moment.
Definition: etoile.h:1642
Cmp log(const Cmp &)
Neperian logarithm.
Definition: cmp_math.C:296
virtual double grv3(ostream *ost=0x0) const
Error on the virial identity GRV3.
Tenseur E_em
electromagnetic energy density in the Eulerian frame
Definition: et_rot_mag.h:161
Component of a tensorial field *** DEPRECATED : use class Scalar instead ***.
Definition: cmp.h:443
Tenseur tnphi
Component of the shift vector.
Definition: etoile.h:1515
double GyroMag() const
Gyromagnetic ratio .
const Tenseur & gradient_spher() const
Returns the gradient of *this (Spherical coordinates) (scalar field only).
Definition: tenseur.C:1548
virtual double angu_mom() const
Angular momentum.
Cmp sqrt(const Cmp &)
Square root.
Definition: cmp_math.C:220
void set_std_base()
Set the standard spectal basis of decomposition for each component.
Definition: tenseur.C:1170
const Base_vect_spher & get_bvect_spher() const
Returns the orthonormal vectorial basis associated with the coordinates of the mapping.
Definition: map.h:783
Lorene prototypes.
Definition: app_hor.h:64
Standard units of space, time and mass.
const Mg3d * get_mg() const
Gives the Mg3d on which the mapping is defined.
Definition: map.h:765
int get_etat() const
Returns the logical state.
Definition: cmp.h:896
Tenseur nnn
Total lapse function.
Definition: etoile.h:509
Tenseur nphi
Metric coefficient .
Definition: etoile.h:1510
Tenseur s_euler
Trace of the stress tensor in the Eulerian frame.
Definition: etoile.h:468
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...
const Cmp & srdsdt() const
Returns of *this .
Definition: cmp_deriv.C:105
Tenseur b_car
Square of the metric factor B.
Definition: etoile.h:1507
virtual double mass_g() const
Gravitational mass.
Tenseur Magn() const
Computes the magnetic field spherical components in Lorene&#39;s units.
const Valeur & sx() const
Returns (r -sampling = RARE ) \ Id (r sampling = FIN ) \ (r -sampling = UNSURR ) ...
Definition: valeur_sx.C:110
Values and coefficients of a (real-value) function.
Definition: valeur.h:287
Tenseur press
Fluid pressure.
Definition: etoile.h:461
void set_base(const Base_val &)
Sets the bases for spectral expansions (member base )
Definition: valeur.C:810
Cmp A_phi
-component of the electromagnetic potential 1-form divided by .
Definition: et_rot_mag.h:155
virtual double grv2() const
Error on the virial identity GRV2.
Tenseur Srr_em
rr component of the electromagnetic stress 3-tensor, as measured in the Eulerian frame. (not used and set to 0, should be supressed)
Definition: et_rot_mag.h:170
void mult_r()
Multiplication by r everywhere.
Definition: cmp_r_manip.C:91
Cmp & set()
Read/write for a scalar (see also operator=(const Cmp&) ).
Definition: tenseur.C:824
double * p_grv3
Error on the virial identity GRV3.
Definition: etoile.h:1634
Tenseur flat_scalar_prod_desal(const Tenseur &t1, const Tenseur &t2)
Same as flat_scalar_prod but with desaliasing.
const Valeur & ssint() const
Returns of *this.
Definition: valeur_ssint.C:112
Tenseur nbar
Baryon density in the fluid frame.
Definition: etoile.h:459
Tenseur gam_euler
Lorentz factor between the fluid and Eulerian observers.
Definition: etoile.h:471
double * p_tsw
Ratio T/W.
Definition: etoile.h:1632
virtual double mass_b() const
Baryon mass.
Base class for pure radial mappings.
Definition: map.h:1536
Map & mp
Mapping associated with the star.
Definition: etoile.h:429
int get_nzone() const
Returns the number of domains.
Definition: grilles.h:448
double integrale() const
Computes the integral over all space of *this .
Definition: cmp_integ.C:55
Tenseur ak_car
Scalar .
Definition: etoile.h:1586
int get_etat() const
Returns the logical state.
Definition: tenseur.h:704
static double lambda_grv2(const Cmp &sou_m, const Cmp &sou_q)
Computes the coefficient which ensures that the GRV2 virial identity is satisfied.
Coord xsr
in the nucleus; \ 1/R in the non-compactified shells; \ in the compactified outer domain...
Definition: map.h:1549
Tenseur bbb
Metric factor B.
Definition: etoile.h:1504
Cmp pow(const Cmp &, int)
Power .
Definition: cmp_math.C:348
void inc2_dzpuis()
Increases by 2 the value of dzpuis and changes accordingly the values of the Cmp in the external comp...
Definition: cmp_r_manip.C:192
Tenseur uuu
Norm of u_euler.
Definition: etoile.h:1518
double * p_mass_g
Gravitational mass.
Definition: etoile.h:548
double omega
Rotation angular velocity ([f_unit] )
Definition: etoile.h:1501
void std_base_scal()
Sets the spectral bases of the Valeur va to the standard ones for a scalar.
Definition: cmp.C:644
Tenseur a_car
Total conformal factor .
Definition: etoile.h:515
Cmp j_t
t-component of the current 4-vector
Definition: et_rot_mag.h:158
bool relativistic
Indicator of relativity: true for a relativistic star, false for a Newtonian one. ...
Definition: etoile.h:437
Tbl & set(int l)
Read/write of the value in a given domain.
Definition: cmp.h:721
Tenseur ener
Total energy density in the fluid frame.
Definition: etoile.h:460
virtual double tsw() const
Ratio T/W.
int get_dzpuis() const
Returns dzpuis.
Definition: cmp.h:900
double Q_int() const
Computed charge from the integration of charge density over the star (i.e.
Tenseur Elec() const
Computes the electric field spherical components in Lorene&#39;s units.
bool check_dzpuis(int dzi) const
Returns false if the last domain is compactified and *this is not zero in this domain and dzpuis is n...
Definition: cmp.C:715
Valeur ** asymptot(int n, const int flag=0) const
Asymptotic expansion at r = infinity.
Definition: cmp_asymptot.C:71
const Valeur & mult_ct() const
Returns applied to *this.
Cmp A_t
t-component of the elecctromagnetic potential 1-form, divided by .
Definition: et_rot_mag.h:150
virtual void MHD_comput()
Computes the electromagnetic part of the stress-energy tensor.
Tenseur & logn
Metric potential = logn_auto.
Definition: etoile.h:1521
void set_dzpuis(int)
Set a value to dzpuis.
Definition: cmp.C:654
int get_nt(int l) const
Returns the number of points in the co-latitude direction ( ) in domain no. l.
Definition: grilles.h:457
double * p_grv2
Error on the virial identity GRV2.
Definition: etoile.h:1633
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition: tenseur.C:636
double Q_comput() const
Computed charge deduced from the asymptotic behaviour of At [SI units].
double MagMom() const
Magnetic Momentum in SI units.
Tenseur ener_euler
Total energy density in the Eulerian frame.
Definition: etoile.h:465
Tenseur & dzeta
Metric potential = beta_auto.
Definition: etoile.h:1534
Cmp B_phi
-component of the magnetic field
Definition: et_rot_mag.h:157
Valeur va
The numerical value of the Cmp.
Definition: cmp.h:461
virtual double mom_quad_old() const
Part of the quadrupole moment.
Standard electro-magnetic units.
Tensor handling *** DEPRECATED : use class Tensor instead ***.
Definition: tenseur.h:298
Tenseur Jp_em
component of the electromagnetic momentum density 3-vector, as measured in the Eulerian frame...
Definition: et_rot_mag.h:167
Tenseur Spp_em
component of the electromagnetic stress 3-tensor, as measured in the Eulerian frame.
Definition: et_rot_mag.h:173
double * p_angu_mom
Angular momentum.
Definition: etoile.h:1631