LORENE
et_magnetisation_comp.C
1 /*
2  * Computational functions for magnetized rotating equilibrium
3  *
4  * (see file et_rot_mag.h for documentation)
5  *
6  */
7 
8 /*
9  * Copyright (c) 2013 Debarati Chatterjee, Jerome Novak
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_comp.C,v 1.15 2016/12/05 16:17:53 j_novak Exp $
33  * $Log: et_magnetisation_comp.C,v $
34  * Revision 1.15 2016/12/05 16:17:53 j_novak
35  * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
36  *
37  * Revision 1.14 2016/11/01 09:12:59 j_novak
38  * Correction of a missing '-' in mom_quad_old().
39  *
40  * Revision 1.13 2015/06/12 12:38:25 j_novak
41  * Implementation of the corrected formula for the quadrupole momentum.
42  *
43  * Revision 1.12 2014/10/21 09:23:54 j_novak
44  * Addition of global functions mass_g(), angu_mom(), grv2/3() and mom_quad().
45  *
46  * Revision 1.11 2014/10/13 08:52:57 j_novak
47  * Lorene classes and functions now belong to the namespace Lorene.
48  *
49  * Revision 1.10 2014/07/04 12:08:02 j_novak
50  * Added some filtering.
51  *
52  * Revision 1.9 2014/05/14 15:19:05 j_novak
53  * The magnetisation field is now filtered.
54  *
55  * Revision 1.8 2014/05/13 15:37:12 j_novak
56  * Updated to new magnetic units.
57  *
58  * Revision 1.7 2014/05/01 13:07:16 j_novak
59  * Fixed two bugs: in the computation of F31,F32 and the triad of U_up.
60  *
61  * Revision 1.6 2014/04/29 13:46:07 j_novak
62  * Addition of switches 'use_B_in_eos' and 'include_magnetisation' to control the model.
63  *
64  * Revision 1.5 2014/04/28 14:53:29 j_novak
65  * Minor modif.
66  *
67  * Revision 1.4 2014/04/28 12:48:13 j_novak
68  * Minor modifications.
69  *
70  * Revision 1.2 2013/12/19 17:05:40 j_novak
71  * Corrected a dzpuis problem.
72  *
73  * Revision 1.1 2013/12/13 16:36:51 j_novak
74  * Addition and computation of magnetisation terms in the Einstein equations.
75  *
76  *
77  *
78  * $Header: /cvsroot/Lorene/C++/Source/Etoile/et_magnetisation_comp.C,v 1.15 2016/12/05 16:17:53 j_novak Exp $
79  *
80  */
81 
82 // Headers C
83 #include <cstdlib>
84 #include <cmath>
85 
86 // Headers Lorene
87 #include "et_rot_mag.h"
88 #include "metric.h"
89 #include "utilitaires.h"
90 #include "param.h"
91 #include "proto_f77.h"
92 #include "unites.h"
93 
94 namespace Lorene {
95 
96  using namespace Unites_mag ;
97 
98 // Algo du papier de 1995
99 
100 void Et_magnetisation::magnet_comput(const int adapt_flag,
101  Cmp (*f_j)(const Cmp&, const double),
102  Param& par_poisson_At,
103  Param& par_poisson_Avect){
104  double relax_mag = 0.5 ;
105 
106  int Z = mp.get_mg()->get_nzone();
107 
108  bool adapt(adapt_flag) ;
109  /****************************************************************
110  * Assertion that all zones have same number of points in theta
111  ****************************************************************/
112  int nt = mp.get_mg()->get_nt(nzet-1) ;
113  for (int l=0; l<Z; l++) assert(mp.get_mg()->get_nt(l) == nt) ;
114 
115  Tbl Rsurf(nt) ;
116  Rsurf.set_etat_qcq() ;
117  mp.r.fait() ;
118  mp.tet.fait() ;
119  Mtbl* theta = mp.tet.c ;
120  const Map_radial* mpr = dynamic_cast<const Map_radial*>(&mp) ;
121  assert (mpr != 0x0) ;
122  for (int j=0; j<nt; j++)
123  Rsurf.set(j) = mpr->val_r_jk(l_surf()(0,j), xi_surf()(0,j), j, 0) ;
124 
125 
126  // Calcul de A_0t dans l'etoile (conducteur parfait)
127 
128  Cmp A_0t(- omega * A_phi) ;
129  A_0t.annule(nzet,Z-1) ;
130 
131  Tenseur ATTENS(A_t) ;
132  Tenseur APTENS(A_phi) ;
133  Tenseur BMN(-logn) ;
134  BMN = BMN + log(bbb) ;
135  BMN.set_std_base() ;
136 
137 
139  nphi.gradient_spher())());
141  nphi.gradient_spher())()) ;
143  BMN.gradient_spher())()
144  + 2*nphi()*flat_scalar_prod_desal(APTENS.gradient_spher(),
145  BMN.gradient_spher())()) ;
146 
147  Cmp ATANT(A_phi.srdsdt()); // Constrction par copie pour mapping
148 
149  ATANT.va = ATANT.va.mult_ct().ssint() ;
150 
151  Cmp ttnphi(tnphi()) ;
152  ttnphi.mult_rsint() ;
153  Cmp BLAH(- b_car()/(nnn()*nnn())*ttnphi*grad1) ;
154  BLAH -= (1+b_car()/(nnn()*nnn())*tnphi()*tnphi())*grad2 ;
155  Cmp nphisr(nphi()) ;
156  nphisr.div_r() ;
157  Cmp npgrada(2*nphisr*(A_phi.dsdr()+ATANT )) ;
158  npgrada.inc2_dzpuis() ;
159  BLAH -= grad3 + npgrada ;
160  Cmp gtt(-nnn()*nnn()+b_car()*tnphi()*tnphi()) ;
161  Cmp gtphi( - b_car()*ttnphi) ;
162 
163  // Computation of j_t thanks to Maxwell-Gauss
164  // modified to include Magnetisation
165  // components of F
166  Cmp F01 = 1/(a_car()*nnn()*nnn())*A_0t.dsdr()
167  + 1/(a_car()*nnn()*nnn())*nphi()*A_phi.dsdr() ;
168 
169  Cmp F02 = 1/(a_car()*nnn()*nnn())*A_0t.srdsdt()
170  + 1/(a_car()*nnn()*nnn())*nphi()*A_phi.srdsdt() ;
171 
172  Cmp tmp = A_phi.dsdr() / (bbb() * bbb() * a_car() );
173  tmp.div_rsint() ;
174  tmp.div_rsint() ;
175  Cmp F31 = 1/(a_car()*nnn()*nnn())*nphi()*nphi()*A_phi.dsdr()
176  + 1/(a_car()*nnn()*nnn())*nphi()*A_0t.dsdr()
177  + tmp ;
178 
179  tmp = A_phi.srdsdt() / (bbb() * bbb() * a_car() );
180  tmp.div_rsint() ;
181  tmp.div_rsint() ;
182  Cmp F32 = 1/(a_car()*nnn()*nnn())*nphi()*nphi()*A_phi.srdsdt()
183  + 1/(a_car()*nnn()*nnn())*nphi()*A_0t.srdsdt()
184  + tmp ;
185 
186  Cmp x = get_magnetisation();
187  Cmp one_minus_x = 1 - x ;
188  one_minus_x.std_base_scal() ;
189 
190  tmp = ((BLAH - A_0t.laplacien())*one_minus_x/a_car()
191  - gtphi*j_phi
192  - gtt*(F01*x.dsdr()+F02*x.srdsdt())
193  - gtphi*(F31*x.dsdr()+F32*x.srdsdt()) ) / gtt ;
194 
195  tmp.annule(nzet, Z-1) ;
196  if (adapt) {
197  j_t = tmp ;
198  }
199  else {
200  j_t.allocate_all() ;
201  for (int j=0; j<nt; j++)
202  for (int l=0; l<nzet; l++)
203  for (int i=0; i<mp.get_mg()->get_nr(l); i++)
204  j_t.set(l,0,j,i) = ( (*mp.r.c)(l,0,j,i) > Rsurf(j) ?
205  0. : tmp(l,0,j,i) ) ;
206  j_t.annule(nzet,Z-1) ;
207  }
208  j_t.std_base_scal() ;
209 
210  // Calcul du courant j_phi
211  j_phi = omega * j_t + (ener() + press())*f_j(A_phi, a_j) ;
212  j_phi.std_base_scal() ;
213 
214  // Resolution de Maxwell Ampere (-> A_phi)
215  // Calcul des termes sources avec A-t du pas precedent.
216 
218  BMN.gradient_spher())());
219 
220  Tenseur source_tAphi(mp, 1, CON, mp.get_bvect_spher()) ;
221 
222  source_tAphi.set_etat_qcq() ;
223  Cmp tjphi(j_phi) ;
224  tjphi.mult_rsint() ;
225  Cmp tgrad1(grad1) ;
226  tgrad1.mult_rsint() ;
227  Cmp d_grad4(grad4) ;
228  d_grad4.div_rsint() ;
229  source_tAphi.set(0)=0 ;
230  source_tAphi.set(1)=0 ;
231 
232 // modified to include Magnetisation
233  Cmp phifac = (F31-nphi()*F01)*x.dsdr()
234  + (F32-nphi()*F02)*x.srdsdt() ;
235  phifac.mult_rsint();
236  source_tAphi.set(2)= -b_car()*a_car()/one_minus_x
237  *(tjphi-tnphi()*j_t + phifac)
238  + b_car()/(nnn()*nnn())*(tgrad1+tnphi()*grad2)
239  + d_grad4 ;
240 
241  source_tAphi.change_triad(mp.get_bvect_cart());
242 
243  // Filtering
244  for (int i=0; i<3; i++) {
245  Scalar tmp_filter = source_tAphi(i) ;
246  tmp_filter.exponential_filter_r(0, 2, 1) ;
247  tmp_filter.exponential_filter_ylm(0, 2, 1) ;
248  source_tAphi.set(i) = tmp_filter ;
249  }
250 
251  Tenseur WORK_VECT(mp, 1, CON, mp.get_bvect_cart()) ;
252  WORK_VECT.set_etat_qcq() ;
253  for (int i=0; i<3; i++) {
254  WORK_VECT.set(i) = 0 ;
255  }
256  Tenseur WORK_SCAL(mp) ;
257  WORK_SCAL.set_etat_qcq() ;
258  WORK_SCAL.set() = 0 ;
259 
260  double lambda_mag = 0. ; // No 3D version !
261 
262  Tenseur AVECT(source_tAphi) ;
263  if (source_tAphi.get_etat() != ETATZERO) {
264 
265  for (int i=0; i<3; i++) {
266  if(source_tAphi(i).dz_nonzero()) {
267  assert( source_tAphi(i).get_dzpuis() == 4 ) ;
268  }
269  else{
270  (source_tAphi.set(i)).set_dzpuis(4) ;
271  }
272  }
273 
274  }
275 
276  source_tAphi.poisson_vect(lambda_mag, par_poisson_Avect, AVECT, WORK_VECT,
277  WORK_SCAL) ;
278  AVECT.change_triad(mp.get_bvect_spher());
279  Cmp A_phi_n(AVECT(2));
280  A_phi_n.mult_rsint() ;
281 
282  // Solution to Maxwell-Ampere : A_1
283  // modified to include Magnetisation
284  Cmp source_A_1t(-a_car()*( j_t*gtt + j_phi*gtphi
285  + gtt*(F01*x.dsdr()+F02*x.srdsdt())
286  + gtphi*(F31*x.dsdr()+F32*x.srdsdt()) )/one_minus_x
287  + BLAH);
288  Scalar tmp_filter = source_A_1t ;
289  tmp_filter.exponential_filter_r(0, 2, 1) ;
290  tmp_filter.exponential_filter_ylm(0, 2, 1) ;
291  source_A_1t = tmp_filter ;
292 
293  Cmp A_1t(mp);
294  A_1t = 0 ;
295  source_A_1t.poisson(par_poisson_At, A_1t) ;
296 
297  int L = mp.get_mg()->get_nt(0);
298 
299  Tbl MAT(L,L) ;
300  Tbl MAT_PHI(L,L);
301  Tbl VEC(L) ;
302 
303  MAT.set_etat_qcq() ;
304  VEC.set_etat_qcq() ;
305  MAT_PHI.set_etat_qcq() ;
306 
307  Tbl leg(L,2*L) ;
308  leg.set_etat_qcq() ;
309 
310  Cmp psi(mp);
311  Cmp psi2(mp);
312  psi.allocate_all() ;
313  psi2.allocate_all() ;
314 
315  for (int p=0; p<mp.get_mg()->get_np(0); p++) {
316  // leg[k,l] : legendre_l(cos(theta_k))
317  // Construction par recurrence de degre 2
318  for(int k=0;k<L;k++){
319  for(int l=0;l<2*L;l++){
320 
321  if(l==0) leg.set(k,l)=1. ;
322  if(l==1) leg.set(k,l)=cos((*theta)(l_surf()(p,k),p,k,0)) ;
323  if(l>=2) leg.set(k,l) = double(2*l-1)/double(l)
324  * cos((*theta)(l_surf()(p,k),p,k,0))
325  * leg(k,l-1)-double(l-1)/double(l)*leg(k,l-2) ;
326  }
327  }
328 
329  for(int k=0;k<L;k++){
330 
331  // Valeurs a la surface trouvees via va.val_point_jk(l,xisurf,k,p)
332 
333  VEC.set(k) = A_0t.va.val_point_jk(l_surf()(p,k), xi_surf()(p,k), k, p)
334  -A_1t.va.val_point_jk(l_surf()(p,k), xi_surf()(p,k), k, p);
335 
336  for(int l=0;l<L;l++) MAT.set(l,k) = leg(k,2*l)/pow(Rsurf(k),2*l+1);
337 
338  }
339  // appel fortran :
340 
341  int* IPIV=new int[L] ;
342  int INFO ;
343 
344  Tbl MAT_SAVE(MAT) ;
345  Tbl VEC2(L) ;
346  VEC2.set_etat_qcq() ;
347  int un = 1 ;
348 
349  F77_dgesv(&L, &un, MAT.t, &L, IPIV, VEC.t, &L, &INFO) ;
350 
351  // coeffs a_l dans VEC
352 
353  for(int k=0;k<L;k++) {VEC2.set(k)=1. ; }
354 
355  F77_dgesv(&L, &un, MAT_SAVE.t, &L, IPIV, VEC2.t, &L, &INFO) ;
356 
357  delete [] IPIV ;
358 
359  for(int nz=0;nz < Z; nz++){
360  for(int i=0;i< mp.get_mg()->get_nr(nz);i++){
361  for(int k=0;k<L;k++){
362  psi.set(nz,p,k,i) = 0. ;
363  psi2.set(nz,p,k,i) = 0. ;
364  for(int l=0;l<L;l++){
365  psi.set(nz,p,k,i) += VEC(l)*leg(k,2*l) /
366  pow((*mp.r.c)(nz,p,k,i),2*l+1);
367  psi2.set(nz,p,k,i) += VEC2(l)*leg(k,2*l)/
368  pow((*mp.r.c)(nz, p, k,i),2*l+1);
369  }
370  }
371  }
372  }
373  }
374  psi.std_base_scal() ;
375  psi2.std_base_scal() ;
376 
377  assert(psi.get_dzpuis() == 0) ;
378  int dif = A_1t.get_dzpuis() ;
379  if (dif > 0) {
380  for (int d=0; d<dif; d++) A_1t.dec_dzpuis() ;
381  }
382 
383  if (adapt) {
384  Cmp A_t_ext(A_1t + psi) ;
385  A_t_ext.annule(0,nzet-1) ;
386  A_0t += A_t_ext ;
387  }
388  else {
389  tmp = A_0t ;
390  A_0t.allocate_all() ;
391  for (int j=0; j<nt; j++)
392  for (int l=0; l<Z; l++)
393  for (int i=0; i<mp.get_mg()->get_nr(l); i++)
394  A_0t.set(l,0,j,i) = ( (*mp.r.c)(l,0,j,i) > Rsurf(j) ?
395  A_1t(l,0,j,i) + psi(l,0,j,i) : tmp(l,0,j,i) ) ;
396  }
397  A_0t.std_base_scal() ;
398 
399  tmp_filter = A_0t ;
400  tmp_filter.exponential_filter_r(0, 2, 1) ;
401  tmp_filter.exponential_filter_ylm(0, 2, 1) ;
402  A_0t = tmp_filter ;
403 
404  Valeur** asymp = A_0t.asymptot(1) ;
405 
406  double Q_0 = -4*M_PI*(*asymp[1])(Z-1,0,0,0) ; // utilise A_0t plutot que E
407  delete asymp[0] ;
408  delete asymp[1] ;
409 
410  delete [] asymp ;
411 
412  asymp = psi2.asymptot(1) ;
413 
414  double Q_2 = -4*M_PI*(*asymp[1])(Z-1,0,0,0) ; // A_2t = psi2 a l'infini
415  delete asymp[0] ;
416  delete asymp[1] ;
417 
418  delete [] asymp ;
419 
420  // solution definitive de A_t:
421 
422  double C = (Q-Q_0)/Q_2 ;
423 
424  assert(psi2.get_dzpuis() == 0) ;
425  dif = A_0t.get_dzpuis() ;
426  if (dif > 0) {
427  for (int d=0; d<dif; d++) A_0t.dec_dzpuis() ;
428  }
429  Cmp A_t_n(mp) ;
430  if (adapt) {
431  A_t_n = A_0t + C ;
432  Cmp A_t_ext(A_0t + C*psi2) ;
433  A_t_ext.annule(0,nzet-1) ;
434  A_t_n.annule(nzet,Z-1) ;
435  A_t_n += A_t_ext ;
436  }
437  else {
438  A_t_n.allocate_all() ;
439  for (int j=0; j<nt; j++)
440  for (int l=0; l<Z; l++)
441  for (int i=0; i<mp.get_mg()->get_nr(l); i++)
442  A_t_n.set(l,0,j,i) = ( (*mp.r.c)(l,0,j,i) > Rsurf(j) ?
443  A_0t(l,0,j,i) + C*psi2(l,0,j,i) :
444  A_0t(l,0,j,i) + C ) ;
445  }
446  A_t_n.std_base_scal() ;
447  tmp_filter = A_t_n ;
448  tmp_filter.exponential_filter_r(0, 2, 1) ;
449  tmp_filter.exponential_filter_ylm(0, 2, 1) ;
450  A_t_n = tmp_filter ;
451 
452  asymp = A_t_n.asymptot(1) ;
453 
454  delete asymp[0] ;
455  delete asymp[1] ;
456 
457  delete [] asymp ;
458  A_t = relax_mag*A_t_n + (1.-relax_mag)*A_t ;
459  A_phi = relax_mag*A_phi_n + (1. - relax_mag)*A_phi ;
460 
461 }
462 
463 
465  // Computes the E-M terms of the stress-energy tensor...
466 
467  Tenseur ATTENS(A_t) ;
468 
469  Tenseur APTENS(A_phi) ;
470 
472  APTENS.gradient_spher())() );
474  ATTENS.gradient_spher())() );
476  ATTENS.gradient_spher())() );
477 
478  if (ApAp.get_etat() != ETATZERO) {
479  ApAp.set().div_rsint() ;
480  ApAp.set().div_rsint() ;
481  }
482  if (ApAt.get_etat() != ETATZERO)
483  ApAt.set().div_rsint() ;
484 
485  E_em = 0.5*mu0 * ( 1/(a_car*nnn*nnn) * (AtAt + 2*tnphi*ApAt)
486  + ( (tnphi*tnphi/(a_car*nnn*nnn)) + 1/(a_car*b_car) )*ApAp );
487  Jp_em = -mu0 * (ApAt + tnphi*ApAp) /(a_car*nnn) ;
488  if (Jp_em.get_etat() != ETATZERO) Jp_em.set().mult_rsint() ;
489  Srr_em = 0 ;
490  // Stt_em = -Srr_em
491  Spp_em = E_em ;
492 
493  // ... and those corresponding to the magnetization.
494  Tenseur Efield = Elec() ;
495  Tenseur Bfield = Magn() ;
496 
497  Scalar EiEi ( flat_scalar_prod(Efield, Efield)() ) ;
498  Scalar BiBi ( flat_scalar_prod(Bfield, Bfield)() ) ;
499 
500  Vector U_up(mp, CON, mp.get_bvect_cart()) ;
501  for (int i=1; i<=3; i++)
502  U_up.set(i) = u_euler(i-1) ;
503  U_up.change_triad(mp.get_bvect_spher()) ;
504 
505  Sym_tensor gamij(mp, COV, mp.get_bvect_spher()) ;
506  for (int i=1; i<=3; i++)
507  for (int j=1; j<i; j++) {
508  gamij.set(i,j) = 0 ;
509  }
510  gamij.set(1,1) = a_car() ;
511  gamij.set(2,2) = a_car() ;
512  gamij.set(3,3) = b_car() ;
513  Metric met(gamij) ;
514  Vector Ui = U_up.down(0, met) ;
515 
516  Scalar fac = sqrt(a_car()) ;
517  Vector B_up(mp, CON, mp.get_bvect_spher()) ;
518  B_up.set(1) = Scalar(Bfield(0)) / fac ;
519  B_up.set(2) = Scalar(Bfield(1)) / fac ;
520  B_up.set(3) = 0 ;
521  Vector Bi = B_up.down(0, met) ;
522 
523  fac = Scalar(gam_euler()*gam_euler()) ;
524 
525  E_I = get_magnetisation() * EiEi / mu0 ;
526 
527  J_I = get_magnetisation() * BiBi * Ui / mu0 ;
528  Sij_I = get_magnetisation()
529  * ( (BiBi / fac) * gamij + BiBi*Ui*Ui - Bi*Bi / fac ) / mu0 ;
530 
531  for (int i=1; i<=3; i++)
532  for (int j=i; j<=3; j++)
533  Sij_I.set(i,j).set_dzpuis(0) ;
534 
535 }
536 
537  //----------------------------//
538  // Gravitational mass //
539  //----------------------------//
540 
541 double Et_magnetisation::mass_g() const {
542 
543  if (p_mass_g == 0x0) { // a new computation is required
544 
545  if (relativistic) {
546 
547  // Magnetisation: S_{rr} + S_{\theta\theta}
548  Tenseur SrrplusStt( Cmp(Sij_I(1, 1) + Sij_I(2, 2)) ) ;
549  SrrplusStt = SrrplusStt / a_car ; // S^r_r + S^\theta_\theta
550 
551  Tenseur Spp (Cmp(Sij_I(3, 3))) ; // Magnetisation: S_{\phi\phi}
552  Spp = Spp / b_car ; // S^\phi_\phi
553 
554  Cmp temp(E_I) ;
555  Tenseur E_i (temp) ;
556  Tenseur J_i (Cmp(J_I(3))) ;
557 
558  Tenseur source = nnn * (ener_euler + E_em + E_i
559  + s_euler + Spp_em + SrrplusStt + Spp) +
560  nphi * (Jp_em + J_i)
561  + 2 * bbb * (ener_euler + press) * tnphi * uuu ;
562 
563  source = a_car * bbb * source ;
564 
565  source.set_std_base() ;
566 
567  p_mass_g = new double( source().integrale() ) ;
568 
569 
570  }
571  else{ // Newtonian case
572  p_mass_g = new double( mass_b() ) ; // in the Newtonian case
573  // M_g = M_b
574  }
575  }
576 
577  return *p_mass_g ;
578 
579 }
580 
581  //----------------------------//
582  // Angular momentum //
583  //----------------------------//
584 
586 
587  if (p_angu_mom == 0x0) { // a new computation is required
588 
589  Cmp dens = uuu() ;
590 
591  dens.mult_r() ; // Multiplication by
592  dens.va = (dens.va).mult_st() ; // r sin(theta)
593 
594  if (relativistic) {
595  dens = a_car() * (b_car() * (ener_euler() + press())
596  * dens + bbb() * (Jp_em() + Cmp(J_I(3)) ) ) ;
597  }
598  else { // Newtonian case
599  dens = nbar() * dens ;
600  }
601 
602  dens.std_base_scal() ;
603 
604  p_angu_mom = new double( dens.integrale() ) ;
605 
606  }
607 
608  return *p_angu_mom ;
609 
610 }
611 
612  //----------------------------//
613  // GRV2 //
614  //----------------------------//
615 
616 double Et_magnetisation::grv2() const {
617 
618  if (p_grv2 == 0x0) { // a new computation is required
619 
620  // To get qpig:
621  using namespace Unites ;
622 
623  Tenseur Spp (Cmp(Sij_I(3, 3))) ; //S_{\phi\phi}
624  Spp = Spp / b_car ; // S^\phi_\phi
625 
626  Tenseur sou_m = 2 * qpig * a_car * (press + (ener_euler+press)
627  * uuu*uuu + Spp) ;
628 
629  Tenseur sou_q = 2 * qpig * a_car * Spp_em + 1.5 * ak_car
630  - flat_scalar_prod(logn.gradient_spher(), logn.gradient_spher() ) ;
631 
632  p_grv2 = new double( double(1) - lambda_grv2(sou_m(), sou_q()) ) ;
633 
634  }
635 
636  return *p_grv2 ;
637 
638 }
639 
640 
641  //----------------------------//
642  // GRV3 //
643  //----------------------------//
644 
645 double Et_magnetisation::grv3(ostream* ost) const {
646 
647  if (p_grv3 == 0x0) { // a new computation is required
648 
649  // To get qpig:
650  using namespace Unites ;
651 
652  Tenseur source(mp) ;
653 
654  // Gravitational term [cf. Eq. (43) of Gourgoulhon & Bonazzola
655  // ------------------ Class. Quantum Grav. 11, 443 (1994)]
656 
657  if (relativistic) {
658  Tenseur alpha = dzeta - logn ;
659  Tenseur beta = log( bbb ) ;
660  beta.set_std_base() ;
661 
662  source = 0.75 * ak_car
663  - flat_scalar_prod(logn.gradient_spher(),
664  logn.gradient_spher() )
665  + 0.5 * flat_scalar_prod(alpha.gradient_spher(),
666  beta.gradient_spher() ) ;
667 
668  Cmp aa = alpha() - 0.5 * beta() ;
669  Cmp daadt = aa.srdsdt() ; // 1/r d/dth
670 
671  // What follows is valid only for a mapping of class Map_radial :
672  const Map_radial* mpr = dynamic_cast<const Map_radial*>(&mp) ;
673  if (mpr == 0x0) {
674  cout << "Etoile_rot::grv3: the mapping does not belong"
675  << " to the class Map_radial !" << endl ;
676  abort() ;
677  }
678 
679  // Computation of 1/tan(theta) * 1/r daa/dtheta
680  if (daadt.get_etat() == ETATQCQ) {
681  Valeur& vdaadt = daadt.va ;
682  vdaadt = vdaadt.ssint() ; // division by sin(theta)
683  vdaadt = vdaadt.mult_ct() ; // multiplication by cos(theta)
684  }
685 
686  Cmp temp = aa.dsdr() + daadt ;
687  temp = ( bbb() - a_car()/bbb() ) * temp ;
688  temp.std_base_scal() ;
689 
690  // Division by r
691  Valeur& vtemp = temp.va ;
692  vtemp = vtemp.sx() ; // division by xi in the nucleus
693  // Id in the shells
694  // division by xi-1 in the ZEC
695  vtemp = (mpr->xsr) * vtemp ; // multiplication by xi/r in the nucleus
696  // by 1/r in the shells
697  // by r(xi-1) in the ZEC
698 
699  // In the ZEC, a multiplication by r has been performed instead
700  // of the division:
701  temp.set_dzpuis( temp.get_dzpuis() + 2 ) ;
702 
703  source = bbb() * source() + 0.5 * temp ;
704 
705  }
706  else{
707  source = - 0.5 * flat_scalar_prod(logn.gradient_spher(),
708  logn.gradient_spher() ) ;
709  }
710 
711  source.set_std_base() ;
712 
713  double int_grav = source().integrale() ;
714 
715  // Matter term
716  // -----------
717 
718  if (relativistic) {
719 
720  // S_{rr} + S_{\theta\theta}
721  Tenseur SrrplusStt( Cmp(Sij_I(1, 1) + Sij_I(2, 2)) ) ;
722  SrrplusStt = SrrplusStt / a_car ; // S^r_r + S^\theta_\theta
723 
724  Tenseur Spp (Cmp(Sij_I(3, 3))) ; //S_{\phi\phi}
725  Spp = Spp / b_car ; // S^\phi_\phi
726 
727  source = qpig * a_car * bbb * ( s_euler + Spp_em + SrrplusStt + Spp ) ;
728  }
729  else{
730  source = qpig * ( 3 * press + nbar * uuu * uuu ) ;
731  }
732 
733  source.set_std_base() ;
734 
735  double int_mat = source().integrale() ;
736 
737  // Virial error
738  // ------------
739  if (ost != 0x0) {
740  *ost << "Et_magnetisation::grv3 : gravitational term : " << int_grav
741  << endl ;
742  *ost << "Et_magnetisation::grv3 : matter term : " << int_mat
743  << endl ;
744  }
745 
746  p_grv3 = new double( (int_grav + int_mat) / int_mat ) ;
747 
748  }
749 
750  return *p_grv3 ;
751 
752 }
753 
754  //----------------------------//
755  // Quadrupole moment //
756  //----------------------------//
757 
759 
760  if (p_mom_quad_old == 0x0) { // a new computation is required
761 
762  // To get qpig:
763  using namespace Unites ;
764 
765  // Source for of the Poisson equation for nu
766  // -----------------------------------------
767 
768  Tenseur source(mp) ;
769 
770  if (relativistic) {
771  // S_{rr} + S_{\theta\theta}
772  Tenseur SrrplusStt( Cmp(Sij_I(1, 1) + Sij_I(2, 2)) ) ;
773  SrrplusStt = SrrplusStt / a_car ; // S^r_r + S^\theta_\theta
774 
775  Tenseur Spp (Cmp(Sij_I(3, 3))) ; //S_{\phi\phi}
776  Spp = Spp / b_car ; // S^\phi_\phi
777 
778  Cmp temp(E_I) ;
779  Tenseur E_i(temp) ;
780 
781  Tenseur beta = log(bbb) ;
782  beta.set_std_base() ;
783  source = qpig * a_car *( ener_euler + E_em + E_i
784  + s_euler + Spp_em + SrrplusStt + Spp)
785  + ak_car - flat_scalar_prod(logn.gradient_spher(),
786  logn.gradient_spher() + beta.gradient_spher()) ;
787  }
788  else {
789  source = qpig * nbar ;
790  }
791  source.set_std_base() ;
792 
793  // Multiplication by -r^2 P_2(cos(theta))
794  // [cf Eq.(7) of Salgado et al. Astron. Astrophys. 291, 155 (1994) ]
795  // ------------------------------------------------------------------
796 
797  // Multiplication by r^2 :
798  // ----------------------
799  Cmp& csource = source.set() ;
800  csource.mult_r() ;
801  csource.mult_r() ;
802  if (csource.check_dzpuis(2)) {
803  csource.inc2_dzpuis() ;
804  }
805 
806  // Muliplication by cos^2(theta) :
807  // -----------------------------
808  Cmp temp = csource ;
809 
810  // What follows is valid only for a mapping of class Map_radial :
811  assert( dynamic_cast<const Map_radial*>(&mp) != 0x0 ) ;
812 
813  if (temp.get_etat() == ETATQCQ) {
814  Valeur& vtemp = temp.va ;
815  vtemp = vtemp.mult_ct() ; // multiplication by cos(theta)
816  vtemp = vtemp.mult_ct() ; // multiplication by cos(theta)
817  }
818 
819  // Muliplication by -P_2(cos(theta)) :
820  // ----------------------------------
821  source = 0.5 * source() - 1.5 * temp ;
822 
823  // Final result
824  // ------------
825  p_mom_quad_old = new double( - source().integrale() / qpig ) ;
826  }
827  return *p_mom_quad_old ;
828  }
829 
830 
832 
833  using namespace Unites ;
834 
835  if (p_mom_quad_Bo == 0x0) { // a new computation is required
836 
837  // S_{rr} + S_{\theta\theta} = A^2*(S^r_r + S^\theta_\theta)
838  Tenseur SrrplusStt( Cmp(Sij_I(1, 1) + Sij_I(2, 2)) ) ;
839 
840  Cmp dens = a_car() * press() ;
841  dens = bbb() * nnn() * (SrrplusStt() + 2*dens) ;
842  dens.mult_rsint() ;
843  dens.std_base_scal() ;
844 
845  p_mom_quad_Bo = new double( - 16. * dens.integrale() / qpig ) ;
846  }
847  return *p_mom_quad_Bo ;
848  }
849 
850 
851 
852 }
const Cmp & dsdr() const
Returns of *this .
Definition: cmp_deriv.C:87
Metric for tensor calculation.
Definition: metric.h:90
virtual double mom_quad_Bo() const
Part of the quadrupole moment.
Cmp log(const Cmp &)
Neperian logarithm.
Definition: cmp_math.C:299
Component of a tensorial field *** DEPRECATED : use class Scalar instead ***.
Definition: cmp.h:446
const Tenseur & gradient_spher() const
Returns the gradient of *this (Spherical coordinates) (scalar field only).
Definition: tenseur.C:1564
Cmp sqrt(const Cmp &)
Square root.
Definition: cmp_math.C:223
void dec_dzpuis()
Decreases by 1 the value of dzpuis and changes accordingly the values of the Cmp in the external comp...
Definition: cmp_r_manip.C:157
void annule(int l)
Sets the Cmp to zero in a given domain.
Definition: cmp.C:351
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
Lorene prototypes.
Definition: app_hor.h:67
Standard units of space, time and mass.
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
Tensor field of valence 0 (or component of a tensorial field).
Definition: scalar.h:393
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:108
const Valeur & sx() const
Returns (r -sampling = RARE ) \ Id (r sampling = FIN ) \ (r -sampling = UNSURR ) ...
Definition: valeur_sx.C:113
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.
Values and coefficients of a (real-value) function.
Definition: valeur.h:297
Tensor field of valence 1.
Definition: vector.h:188
Cmp cos(const Cmp &)
Cosine.
Definition: cmp_math.C:97
void div_r()
Division by r everywhere.
Definition: cmp_r_manip.C:81
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition: tbl.C:364
virtual double angu_mom() const
Angular momentum.
void mult_r()
Multiplication by r everywhere.
Definition: cmp_r_manip.C:94
Cmp & set()
Read/write for a scalar (see also operator=(const Cmp&) ).
Definition: tenseur.C:840
void change_triad(const Base_vect &new_triad)
Sets a new vectorial basis (triad) of decomposition and modifies the components accordingly.
Definition: tenseur.C:684
Tenseur flat_scalar_prod_desal(const Tenseur &t1, const Tenseur &t2)
Same as flat_scalar_prod but with desaliasing.
const Valeur & ssint() const
Returns of *this.
Definition: valeur_ssint.C:115
double * t
The array of double.
Definition: tbl.h:176
Parameter storage.
Definition: param.h:125
Base class for pure radial mappings.
Definition: map.h:1557
void mult_rsint()
Multiplication by .
Definition: cmp_r_manip.C:119
virtual double val_r_jk(int l, double xi, int j, int k) const =0
Returns the value of the radial coordinate r for a given and a given collocation point in in a give...
double integrale() const
Computes the integral over all space of *this .
Definition: cmp_integ.C:58
double val_point_jk(int l, double x, int j, int k) const
Computes the value of the field represented by *this at an arbitrary point in , but collocation point...
Definition: valeur.C:903
Coord xsr
in the nucleus; \ 1/R in the non-compactified shells; \ in the compactified outer domain...
Definition: map.h:1570
Cmp pow(const Cmp &, int)
Power .
Definition: cmp_math.C:351
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:195
void std_base_scal()
Sets the spectral bases of the Valeur va to the standard ones for a scalar.
Definition: cmp.C:647
virtual double mom_quad_old() const
Part of the quadrupole moment.
Cmp poisson() const
Solves the scalar Poisson equation with *this as a source.
Definition: cmp_pde.C:97
virtual double grv2() const
Error on the virial identity GRV2.
void allocate_all()
Sets the logical state to ETATQCQ (ordinary state) and performs the memory allocation of all the elem...
Definition: cmp.C:326
Tbl & set(int l)
Read/write of the value in a given domain.
Definition: cmp.h:724
int get_dzpuis() const
Returns dzpuis.
Definition: cmp.h:903
Scalar & set(const Itbl &ind)
Returns the value of a component (read/write version).
Definition: tensor.C:663
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...
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:718
virtual double mass_g() const
Gravitational mass.
Valeur ** asymptot(int n, const int flag=0) const
Asymptotic expansion at r = infinity.
Definition: cmp_asymptot.C:74
const Valeur & mult_ct() const
Returns applied to *this.
virtual double grv3(ostream *ost=0x0) const
Error on the virial identity GRV3.
void set_dzpuis(int)
Set a value to dzpuis.
Definition: cmp.C:657
Basic array class.
Definition: tbl.h:164
Scalar & set(int)
Read/write access to a component.
Definition: vector.C:302
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition: tenseur.C:652
virtual void exponential_filter_ylm(int lzmin, int lzmax, int p, double alpha=-16.)
Applies an exponential filter to the spectral coefficients in the angular directions.
const Cmp & laplacien(int zec_mult_r=4) const
Returns the Laplacian of *this.
Definition: cmp_deriv.C:245
Valeur va
The numerical value of the Cmp.
Definition: cmp.h:464
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
virtual void MHD_comput()
Computes the electromagnetic part of the stress-energy tensor.
Tensor down(int ind, const Metric &gam) const
Computes a new tensor by lowering an index of *this.