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