LORENE
boson_star_equil.C
1 /*
2  * Method Boson_star::equilibrium
3  *
4  * (see file boson_star.h for documentation).
5  */
6 
7 /*
8  * Copyright (c) 2012 Claire Some, Eric Gourgoulhon
9  *
10  * This file is part of LORENE.
11  *
12  * LORENE is free software; you can redistribute it and/or modify
13  * it under the terms of the GNU General Public License as published by
14  * the Free Software Foundation; either version 2 of the License, or
15  * (at your option) any later version.
16  *
17  * LORENE is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
20  * GNU General Public License for more details.
21  *
22  * You should have received a copy of the GNU General Public License
23  * along with LORENE; if not, write to the Free Software
24  * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
25  *
26  */
27 
28 
29 char boson_star_equil_C[] = "$Header: /cvsroot/Lorene/C++/Source/Compobj/boson_star_equil.C,v 1.6 2014/10/13 08:52:49 j_novak Exp $" ;
30 
31 /*
32  * $Id: boson_star_equil.C,v 1.6 2014/10/13 08:52:49 j_novak Exp $
33  * $Log: boson_star_equil.C,v $
34  * Revision 1.6 2014/10/13 08:52:49 j_novak
35  * Lorene classes and functions now belong to the namespace Lorene.
36  *
37  * Revision 1.5 2014/10/06 15:13:04 j_novak
38  * Modified #include directives to use c++ syntax.
39  *
40  * Revision 1.4 2013/04/03 12:10:13 e_gourgoulhon
41  * Added member kk to Compobj; suppressed tkij
42  *
43  * Revision 1.3 2012/12/03 15:27:30 c_some
44  * Small changes
45  *
46  * Revision 1.2 2012/11/23 15:43:05 c_some
47  * Small changes
48  *
49  * Revision 1.1 2012/11/22 16:04:12 c_some
50  * New class Boson_star
51  *
52  *
53  * $Header: /cvsroot/Lorene/C++/Source/Compobj/boson_star_equil.C,v 1.6 2014/10/13 08:52:49 j_novak Exp $
54  *
55  */
56 
57 // Headers C
58 #include <cmath>
59 
60 // Headers Lorene
61 #include "boson_star.h"
62 #include "param.h"
63 #include "tenseur.h"
64 
65 #include "graphique.h"
66 #include "utilitaires.h"
67 #include "unites.h"
68 
69 namespace Lorene {
70 void Boson_star::equilibrium(double, double,
71  int nzadapt, const Tbl& phi_limit, const Itbl& icontrol,
72  const Tbl& control, Tbl& diff, Param*) {
73 
74  // Fundamental constants and units
75  // -------------------------------
76 
77  using namespace Unites ;
78 
79  // For the display
80  // ---------------
81  char display_bold[]="x[1m" ; display_bold[0] = 27 ;
82  char display_normal[] = "x[0m" ; display_normal[0] = 27 ;
83 
84  // Grid parameters
85  // ---------------
86 
87  const Mg3d* mg = mp.get_mg() ;
88  int nz = mg->get_nzone() ; // total number of domains
89 
90  // The following is required to initialize mp_prev as a Map_et:
91  Map_et& mp_et = dynamic_cast<Map_et&>(mp) ;
92 
93  // Index of the point at phi=0, theta=pi/2 at the surface of the star:
94 
95  int nzet = nzadapt ; //## to be checked
96 
97  assert(mg->get_type_t() == SYM) ;
98  int l_b = nzet - 1 ;
99  int j_b = mg->get_nt(l_b) - 1 ;
100  int k_b = 0 ;
101 
102  // Value of the enthalpy defining the surface of the star
103  // double ent_b = phi_limit(nzet-1) ;
104 
105  // Parameters to control the iteration
106  // -----------------------------------
107 
108  int mer_max = icontrol(0) ;
109 // int mer_rot = icontrol(1) ;
110 // int mer_change_omega = icontrol(2) ;
111 // int mer_fix_omega = icontrol(3) ;
112 //## int mer_mass = icontrol(4) ;
113  int mermax_poisson = icontrol(5) ;
114  int mer_triax = icontrol(6) ;
115 //## int delta_mer_kep = icontrol(7) ;
116 
117 
118  double precis = control(0) ;
119 // double omega_ini = control(1) ;
120  double relax = control(2) ;
121  double relax_prev = double(1) - relax ;
122  double relax_poisson = control(3) ;
123 // double thres_adapt = control(4) ;
124  double ampli_triax = control(5) ;
125  double precis_adapt = control(6) ;
126 
127 
128  // Error indicators
129  // ----------------
130 
131  diff.set_etat_qcq() ;
132  double& diff_phi = diff.set(0) ;
133  double& diff_nuf = diff.set(1) ;
134  double& diff_nuq = diff.set(2) ;
135 // double& diff_dzeta = diff.set(3) ;
136 // double& diff_ggg = diff.set(4) ;
137  double& diff_shift_x = diff.set(5) ;
138  double& diff_shift_y = diff.set(6) ;
139  double& vit_triax = diff.set(7) ;
140 
141  // Parameters for the function Map_et::adapt
142  // -----------------------------------------
143 
144  Param par_adapt ;
145  int nitermax = 100 ;
146  int niter ;
147  int adapt_flag = 1 ; // 1 = performs the full computation,
148  // 0 = performs only the rescaling by
149  // the factor alpha_r
150  int nz_search = nzet + 1 ; // Number of domains for searching the enthalpy
151  // isosurfaces
152  double alpha_r ;
153  double reg_map = 1. ; // 1 = regular mapping, 0 = contracting mapping
154 
155  par_adapt.add_int(nitermax, 0) ; // maximum number of iterations to
156  // locate zeros by the secant method
157  par_adapt.add_int(nzadapt, 1) ; // number of domains where the adjustment
158  // to the isosurfaces of ent is to be
159  // performed
160  par_adapt.add_int(nz_search, 2) ; // number of domains to search for
161  // the enthalpy isosurface
162  par_adapt.add_int(adapt_flag, 3) ; // 1 = performs the full computation,
163  // 0 = performs only the rescaling by
164  // the factor alpha_r
165  par_adapt.add_int(j_b, 4) ; // theta index of the collocation point
166  // (theta_*, phi_*)
167  par_adapt.add_int(k_b, 5) ; // theta index of the collocation point
168  // (theta_*, phi_*)
169 
170  par_adapt.add_int_mod(niter, 0) ; // number of iterations actually used in
171  // the secant method
172 
173  par_adapt.add_double(precis_adapt, 0) ; // required absolute precision in
174  // the determination of zeros by
175  // the secant method
176  par_adapt.add_double(reg_map, 1) ; // 1. = regular mapping, 0 = contracting mapping
177 
178  par_adapt.add_double(alpha_r, 2) ; // factor by which all the radial
179  // distances will be multiplied
180 
181  par_adapt.add_tbl(phi_limit, 0) ; // array of values of the field Phi
182  // to define the isosurfaces.
183 
184  // Parameters for the function Map_et::poisson for nuf
185  // ----------------------------------------------------
186 
187  double precis_poisson = 1.e-16 ;
188 
189  // Preparations
190  // ------------
191 
192  // Cartesian components of the shift vector are required
194 
195 
196  Cmp cssjm1_nuf(ssjm1_nuf) ;
197  Cmp cssjm1_nuq(ssjm1_nuq) ;
198  Cmp cssjm1_tggg(ssjm1_tggg) ;
199  Cmp cssjm1_khi(ssjm1_khi) ;
200  Tenseur cssjm1_wshift(mp, 1, CON, mp.get_bvect_cart() ) ;
201  cssjm1_wshift.set_etat_qcq() ;
202  for (int i=1; i<=3; i++) {
203  cssjm1_wshift.set(i-1) = ssjm1_wshift(i) ;
204  }
205 
206  Tenseur cshift(mp, 1, CON, mp.get_bvect_cart() ) ;
207  cshift.set_etat_qcq() ;
208  for (int i=1; i<=3; i++) {
209  cshift.set(i-1) = -beta(i) ;
210  }
211 
212  Tenseur cw_shift(mp, 1, CON, mp.get_bvect_cart() ) ;
213  cw_shift.set_etat_qcq() ;
214  for (int i=1; i<=3; i++) {
215  cw_shift.set(i-1) = w_shift(i) ;
216  }
217 
218  Tenseur ckhi_shift(mp) ;
219  ckhi_shift.set_etat_qcq() ;
220  ckhi_shift.set() = khi_shift ;
221 
222  Param par_poisson_nuf ;
223  par_poisson_nuf.add_int(mermax_poisson, 0) ; // maximum number of iterations
224  par_poisson_nuf.add_double(relax_poisson, 0) ; // relaxation parameter
225  par_poisson_nuf.add_double(precis_poisson, 1) ; // required precision
226  par_poisson_nuf.add_int_mod(niter, 0) ; // number of iterations actually used
227  par_poisson_nuf.add_cmp_mod( cssjm1_nuf ) ;
228 
229  Param par_poisson_nuq ;
230  par_poisson_nuq.add_int(mermax_poisson, 0) ; // maximum number of iterations
231  par_poisson_nuq.add_double(relax_poisson, 0) ; // relaxation parameter
232  par_poisson_nuq.add_double(precis_poisson, 1) ; // required precision
233  par_poisson_nuq.add_int_mod(niter, 0) ; // number of iterations actually used
234  par_poisson_nuq.add_cmp_mod( cssjm1_nuq ) ;
235 
236  Param par_poisson_tggg ;
237  par_poisson_tggg.add_int(mermax_poisson, 0) ; // maximum number of iterations
238  par_poisson_tggg.add_double(relax_poisson, 0) ; // relaxation parameter
239  par_poisson_tggg.add_double(precis_poisson, 1) ; // required precision
240  par_poisson_tggg.add_int_mod(niter, 0) ; // number of iterations actually used
241  par_poisson_tggg.add_cmp_mod( cssjm1_tggg ) ;
242  double lambda_tggg ;
243  par_poisson_tggg.add_double_mod( lambda_tggg ) ;
244 
245  Param par_poisson_dzeta ;
246  double lbda_grv2 ;
247  par_poisson_dzeta.add_double_mod( lbda_grv2 ) ;
248 
249  // Parameters for the function Scalar::poisson_vect
250  // -------------------------------------------------
251 
252  Param par_poisson_vect ;
253 
254  par_poisson_vect.add_int(mermax_poisson, 0) ; // maximum number of iterations
255  par_poisson_vect.add_double(relax_poisson, 0) ; // relaxation parameter
256  par_poisson_vect.add_double(precis_poisson, 1) ; // required precision
257  par_poisson_vect.add_cmp_mod( cssjm1_khi ) ;
258  par_poisson_vect.add_tenseur_mod( cssjm1_wshift ) ;
259  par_poisson_vect.add_int_mod(niter, 0) ;
260 
261 
262  // Initializations
263  // ---------------
264 
265  //## Spherical components of the shift vector are restored
267  update_metric() ;
268  //## Back to Cartesian components
270 
271 
272  // Quantities at the previous step :
273  Map_et mp_prev = mp_et ;
274  Scalar rphi_prev = rphi ;
275  Scalar logn_prev = logn ;
276  Scalar dzeta_prev = dzeta ;
277 
278  // Creation of uninitialized tensors:
279  Scalar source_nuf(mp) ; // source term in the equation for nuf
280  Scalar source_nuq(mp) ; // source term in the equation for nuq
281  Scalar source_dzf(mp) ; // matter source term in the eq. for dzeta
282  Scalar source_dzq(mp) ; // quadratic source term in the eq. for dzeta
283  Scalar source_tggg(mp) ; // source term in the eq. for tggg
284  Vector source_shift(mp, CON, mp.get_bvect_cart()) ;
285  // source term for shift
286 
287 
288  ofstream fichconv("convergence.d") ; // Output file for diff_phi
289  fichconv << "# diff_phi GRV2 max_triax vit_triax" << endl ;
290 
291 
292  ofstream fichevol("evolution.d") ; // Output file for various quantities
293  fichevol <<
294  "# |dH/dr_eq/dH/dr_pole| r_pole/r_eq rphi_c"
295  << endl ;
296 
297  diff_phi = 1 ;
298  double err_grv2 = 1 ;
299  double max_triax_prev = 0 ; // Triaxial amplitude at previous step
300 
301  //=========================================================================
302  // Start of iteration
303  //=========================================================================
304 
305  for(int mer=0 ; (diff_phi > precis) && (mer<mer_max) ; mer++ ) {
306 
307  cout << "-----------------------------------------------" << endl ;
308  cout << "step: " << mer << endl ;
309  cout << "diff_phi = " << display_bold << diff_phi << display_normal
310  << endl ;
311  cout << "err_grv2 = " << err_grv2 << endl ;
312  fichconv << mer ;
313  fichevol << mer ;
314 
315 
316  //-----------------------------------------------
317  // Sources of the Poisson equations
318  //-----------------------------------------------
319 
320  // Source for nu
321  // -------------
322  Scalar bet = log(bbb) ;
323  bet.std_spectral_base() ;
324 
325  Vector d_logn = logn.derive_cov( mp.flat_met_spher() ) ;
326  Vector d_bet = bet.derive_cov( mp.flat_met_spher() ) ;
327 
328  Scalar s_euler = stress_euler.trace(gamma) ;
329  source_nuf = qpig * a_car *( ener_euler + s_euler ) ;
330 
331  source_nuq = ak_car - d_logn(1)*(d_logn(1)+d_bet(1))
332  - d_logn(2)*(d_logn(2)+d_bet(2))
333  - d_logn(3)*(d_logn(3)+d_bet(3)) ;
334 
335  source_nuf.std_spectral_base() ;
336  source_nuq.std_spectral_base() ;
337 
338  // Source for dzeta
339  // ----------------
340  source_dzf = 2 * qpig * a_car * b_car * stress_euler(3,3) ;
341  source_dzf.std_spectral_base() ;
342 
343  source_dzq = 1.5 * ak_car
344  - d_logn(1)*d_logn(1) - d_logn(2)*d_logn(2) - d_logn(3)*d_logn(3) ;
345  source_dzq.std_spectral_base() ;
346 
347  // Source for tggg
348  // ---------------
349 
350  source_tggg = 2 * qpig * nn * a_car * bbb * (s_euler - b_car * stress_euler(3,3)) ;
351  source_tggg.std_spectral_base() ;
352 
353  source_tggg.mult_rsint() ;
354 
355 
356  // Source for shift
357  // ----------------
358 
359  // Matter term:
360  Vector mom_euler_cart = mom_euler ;
361  mom_euler_cart.change_triad(mp.get_bvect_cart()) ;
362  source_shift = (-4*qpig) * nn * a_car * mom_euler_cart ;
363 
364  // Quadratic terms:
365  Vector vtmp = 6 * bet.derive_con( mp.flat_met_spher() )
366  - 2 * logn.derive_con( mp.flat_met_spher() ) ;
367 
368  Vector squad = nn * contract(kk, 1, vtmp, 0) / b_car ;
369  squad.change_triad(mp.get_bvect_cart()) ;
370 
371  source_shift = source_shift + squad.up(0, mp.flat_met_cart() ) ;
372 
373  //----------------------------------------------
374  // Resolution of the Poisson equation for nuf
375  //----------------------------------------------
376 
377  source_nuf.poisson(par_poisson_nuf, nuf) ;
378 
379  cout << "Test of the Poisson equation for nuf :" << endl ;
380  Tbl err = source_nuf.test_poisson(nuf, cout, true) ;
381  diff_nuf = err(0, 0) ;
382 
383  //---------------------------------------
384  // Triaxial perturbation of nuf
385  //---------------------------------------
386 
387  if (mer == mer_triax) {
388 
389  if ( mg->get_np(0) == 1 ) {
390  cout <<
391  "Boson_star::equilibrium: np must be stricly greater than 1"
392  << endl << " to set a triaxial perturbation !" << endl ;
393  abort() ;
394  }
395 
396  const Coord& phi = mp.phi ;
397  const Coord& sint = mp.sint ;
398  Scalar perturb(mp) ;
399  perturb = 1 + ampli_triax * sint*sint * cos(2*phi) ;
400  nuf = nuf * perturb ;
401 
402  nuf.std_spectral_base() ; // set the bases for spectral expansions
403  // to be the standard ones for a
404  // scalar field
405 
406  }
407 
408  // Monitoring of the triaxial perturbation
409  // ---------------------------------------
410 
411  const Valeur& va_nuf = nuf.get_spectral_va() ;
412  va_nuf.coef() ; // Computes the spectral coefficients
413  double max_triax = 0 ;
414 
415  if ( mg->get_np(0) > 1 ) {
416 
417  for (int l=0; l<nz; l++) { // loop on the domains
418  for (int j=0; j<mg->get_nt(l); j++) {
419  for (int i=0; i<mg->get_nr(l); i++) {
420 
421  // Coefficient of cos(2 phi) :
422  double xcos2p = (*(va_nuf.c_cf))(l, 2, j, i) ;
423 
424  // Coefficient of sin(2 phi) :
425  double xsin2p = (*(va_nuf.c_cf))(l, 3, j, i) ;
426 
427  double xx = sqrt( xcos2p*xcos2p + xsin2p*xsin2p ) ;
428 
429  max_triax = ( xx > max_triax ) ? xx : max_triax ;
430  }
431  }
432  }
433 
434  }
435 
436  cout << "Triaxial part of nuf : " << max_triax << endl ;
437 
438  //----------------------------------------------
439  // Resolution of the Poisson equation for nuq
440  //----------------------------------------------
441 
442  source_nuq.poisson(par_poisson_nuq, nuq) ;
443 
444  cout << "Test of the Poisson equation for nuq :" << endl ;
445  err = source_nuq.test_poisson(nuq, cout, true) ;
446  diff_nuq = err(0, 0) ;
447 
448  //---------------------------------------------------------
449  // Resolution of the vector Poisson equation for the shift
450  //---------------------------------------------------------
451 
452 
453  for (int i=1; i<=3; i++) {
454  if(source_shift(i).get_etat() != ETATZERO) {
455  if(source_shift(i).dz_nonzero()) {
456  assert( source_shift(i).get_dzpuis() == 4 ) ;
457  }
458  else{
459  (source_shift.set(i)).set_dzpuis(4) ;
460  }
461  }
462  }
463 
464  double lambda_shift = double(1) / double(3) ;
465 
466  if ( mg->get_np(0) == 1 ) {
467  lambda_shift = 0 ;
468  }
469 
470  Tenseur csource_shift(mp, 1, CON, mp.get_bvect_cart() ) ;
471  csource_shift.set_etat_qcq() ;
472  for (int i=1; i<=3; i++) {
473  csource_shift.set(i-1) = source_shift(i) ;
474  }
475  csource_shift.set(2).set_etat_zero() ; //## bizarre...
476 
477  csource_shift.poisson_vect(lambda_shift, par_poisson_vect,
478  cshift, cw_shift, ckhi_shift) ;
479 
480  for (int i=1; i<=3; i++) {
481  beta.set(i) = - cshift(i-1) ;
482  beta.set(i).set_dzpuis(0) ; //## bizarre...
483  w_shift.set(i) = cw_shift(i-1) ;
484  }
485  khi_shift = ckhi_shift() ;
486 
487  cout << "Test of the Poisson equation for shift_x :" << endl ;
488  err = source_shift(1).test_poisson(-beta(1), cout, true) ;
489  diff_shift_x = err(0, 0) ;
490 
491  cout << "Test of the Poisson equation for shift_y :" << endl ;
492  err = source_shift(2).test_poisson(-beta(2), cout, true) ;
493  diff_shift_y = err(0, 0) ;
494 
495  // Computation of tnphi and nphi from the Cartesian components
496  // of the shift
497  // -----------------------------------------------------------
498 
499  fait_nphi() ;
500 
501 
502 
503  //----------------------------------------------------
504  // Adaptation of the mapping to the new enthalpy field
505  //----------------------------------------------------
506 
507  // Shall the adaptation be performed ?
508  // ---------------------------------
509 
510  adapt_flag = 0 ; // No adaptation of the mapping
511 
512  mp_prev = mp_et ;
513 
514 
515  //---------------------------------------------------------
516  // Matter source terms in the gravitational field equations
517  //---------------------------------------------------------
518 
519  //## Computation of tnphi and nphi from the Cartesian components
520  // of the shift for the test in hydro_euler():
521 
522  fait_nphi() ;
523 
524  // Update of the energy-momentum tensor
525  update_ener_mom() ;
526 
527 
528  //-------------------------------------------------------
529  // 2-D Poisson equation for tggg
530  //-------------------------------------------------------
531 
532  Cmp csource_tggg(source_tggg) ;
533  Cmp ctggg(tggg) ;
534  mp.poisson2d(csource_tggg, mp.cmp_zero(), par_poisson_tggg,
535  ctggg) ;
536  tggg = ctggg ;
537 
538 
539  //-------------------------------------------------------
540  // 2-D Poisson equation for dzeta
541  //-------------------------------------------------------
542 
543  Cmp csource_dzf(source_dzf) ;
544  Cmp csource_dzq(source_dzq) ;
545  Cmp cdzeta(dzeta) ;
546  mp.poisson2d(csource_dzf, csource_dzq, par_poisson_dzeta,
547  cdzeta) ;
548  dzeta = cdzeta ;
549 
550  err_grv2 = lbda_grv2 - 1;
551  cout << "GRV2: " << err_grv2 << endl ;
552 
553 
554  //---------------------------------------
555  // Computation of the metric coefficients (except for N^phi)
556  //---------------------------------------
557 
558  // Relaxations on nu and dzeta :
559 
560  if (mer >= 10) {
561  logn = relax * logn + relax_prev * logn_prev ;
562 
563  dzeta = relax * dzeta + relax_prev * dzeta_prev ;
564  }
565 
566  // Update of the metric coefficients N, A, B and computation of K_ij :
567 
568  //## Spherical components of the shift vector are restored
570  update_metric() ;
571  //## Back to Cartesian components
573 
574  //-----------------------
575  // Informations display
576  //-----------------------
577 
578  cout << *this << endl ;
579 
580 
581  //------------------------------------------------------------
582  // Relative change in Phi with respect to previous step
583  //------------------------------------------------------------
584 
585  Tbl diff_phi_tbl = diffrel( rphi, rphi_prev ) ;
586  diff_phi = diff_phi_tbl(0) ;
587  for (int l=1; l<nzet; l++) {
588  diff_phi += diff_phi_tbl(l) ;
589  }
590  diff_phi /= nzet ;
591 
592  fichconv << " " << log10( fabs(diff_phi) + 1.e-16 ) ;
593  fichconv << " " << log10( fabs(err_grv2) + 1.e-16 ) ;
594  fichconv << " " << log10( fabs(max_triax) + 1.e-16 ) ;
595 
596  vit_triax = 0 ;
597  if ( (mer > mer_triax+1) && (max_triax_prev > 1e-13) ) {
598  vit_triax = (max_triax - max_triax_prev) / max_triax_prev ;
599  }
600 
601  fichconv << " " << vit_triax ;
602 
603  //------------------------------
604  // Recycling for the next step
605  //------------------------------
606 
607  rphi_prev = rphi ;
608  logn_prev = logn ;
609  dzeta_prev = dzeta ;
610  max_triax_prev = max_triax ;
611 
612  fichconv << endl ;
613  fichevol << endl ;
614  fichconv.flush() ;
615  fichevol.flush() ;
616 
617  } // End of main loop
618 
619  //=========================================================================
620  // End of iteration
621  //=========================================================================
622 
623  ssjm1_nuf = cssjm1_nuf ;
624  ssjm1_nuq = cssjm1_nuq ;
625  ssjm1_tggg = cssjm1_tggg ;
626  ssjm1_khi = cssjm1_khi ;
627  for (int i=1; i<=3; i++) {
628  ssjm1_wshift.set(i) = cssjm1_wshift(i-1) ;
629  }
630 
631  // Spherical components of the shift vector are restored
633 
634  fichconv.close() ;
635  fichevol.close() ;
636 
637 }
638 }
Lorene::Star_QI::ssjm1_nuq
Scalar ssjm1_nuq
Effective source at the previous step for the resolution of the Poisson equation for nuq by means of ...
Definition: compobj.h:551
Lorene::Map_et
Radial mapping of rather general form.
Definition: map.h:2752
Lorene::Tbl::set
double & set(int i)
Read/write of a particular element (index i) (1D case)
Definition: tbl.h:281
Lorene::Mg3d::get_np
int get_np(int l) const
Returns the number of points in the azimuthal direction ( ) in domain no. l.
Definition: grilles.h:462
Lorene::Tenseur
Tensor handling *** DEPRECATED : use class Tensor instead ***.
Definition: tenseur.h:301
Lorene::Tensor::up
Tensor up(int ind, const Metric &gam) const
Computes a new tensor by raising an index of *this.
Definition: tensor_calculus.C:225
Lorene::Star_QI::dzeta
Scalar dzeta
Metric potential .
Definition: compobj.h:513
Lorene::Valeur
Values and coefficients of a (real-value) function.
Definition: valeur.h:287
Lorene::Compobj_QI::bbb
Scalar bbb
Metric factor B.
Definition: compobj.h:290
Lorene::Compobj::mp
Map & mp
Mapping describing the coordinate system (r,theta,phi)
Definition: compobj.h:132
Lorene::Boson_star::equilibrium
virtual void equilibrium(double rphi_c, double iphi_c, int nzadapt, const Tbl &ent_limit, const Itbl &icontrol, const Tbl &control, Tbl &diff, Param *=0x0)
Solves the equation satisfied by the scalar field.
Definition: boson_star_equil.C:70
Lorene::Mg3d
Multi-domain grid.
Definition: grilles.h:273
Lorene::Tenseur::set
Cmp & set()
Read/write for a scalar (see also operator=(const Cmp&) ).
Definition: tenseur.C:824
Lorene::log
Cmp log(const Cmp &)
Neperian logarithm.
Definition: cmp_math.C:296
Lorene::Valeur::c_cf
Mtbl_cf * c_cf
Coefficients of the spectral expansion of the function.
Definition: valeur.h:302
Lorene::Scalar::mult_rsint
void mult_rsint()
Multiplication by everywhere; dzpuis is not changed.
Definition: scalar_r_manip.C:281
Lorene::Param::add_double
void add_double(const double &x, int position=0)
Adds the the address of a new double to the list.
Definition: param.C:315
Lorene::Map::cmp_zero
const Cmp & cmp_zero() const
Returns the null Cmp defined on *this.
Definition: map.h:807
Lorene::diffrel
Tbl diffrel(const Cmp &a, const Cmp &b)
Relative difference between two Cmp (norme version).
Definition: cmp_math.C:504
Lorene::Mg3d::get_nt
int get_nt(int l) const
Returns the number of points in the co-latitude direction ( ) in domain no. l.
Definition: grilles.h:457
Lorene
Lorene prototypes.
Definition: app_hor.h:64
Lorene::Boson_star::rphi
Scalar rphi
Real part of the scalar field Phi.
Definition: boson_star.h:72
Lorene::Map::get_bvect_spher
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::Map::get_mg
const Mg3d * get_mg() const
Gives the Mg3d on which the mapping is defined.
Definition: map.h:765
Lorene::Map::flat_met_spher
const Metric_flat & flat_met_spher() const
Returns the flat metric associated with the spherical coordinates and with components expressed in th...
Definition: map.C:321
Lorene::Tenseur::poisson_vect
void poisson_vect(double lambda, Param &par, Tenseur &shift, Tenseur &vect, Tenseur &scal) const
Solves the vectorial Poisson equation : .
Definition: tenseur_pde.C:118
Lorene::Compobj::stress_euler
Sym_tensor stress_euler
Stress tensor with respect to the Eulerian observer.
Definition: compobj.h:150
Lorene::Scalar::derive_cov
const Vector & derive_cov(const Metric &gam) const
Returns the gradient (1-form = covariant vector) of *this
Definition: scalar_deriv.C:390
Lorene::Map::sint
Coord sint
Definition: map.h:721
Lorene::Star_QI::khi_shift
Scalar khi_shift
Scalar used in the decomposition of shift , following Shibata's prescription [Prog.
Definition: compobj.h:539
Lorene::Cmp::set_etat_zero
void set_etat_zero()
Sets the logical state to ETATZERO (zero).
Definition: cmp.C:289
Lorene::Itbl
Basic integer array class.
Definition: itbl.h:122
Lorene::Compobj::beta
Vector beta
Shift vector .
Definition: compobj.h:138
Lorene::Star_QI::ssjm1_wshift
Vector ssjm1_wshift
Effective source at the previous step for the resolution of the vector Poisson equation for .
Definition: compobj.h:578
Lorene::Scalar
Tensor field of valence 0 (or component of a tensorial field).
Definition: scalar.h:387
Lorene::Param::add_double_mod
void add_double_mod(double &x, int position=0)
Adds the address of a new modifiable double to the list.
Definition: param.C:453
Lorene::Vector::change_triad
virtual void change_triad(const Base_vect &)
Sets a new vectorial basis (triad) of decomposition and modifies the components accordingly.
Definition: vector_change_triad.C:75
Lorene::Valeur::coef
void coef() const
Computes the coeffcients of *this.
Definition: valeur_coef.C:148
Lorene::Scalar::get_spectral_va
const Valeur & get_spectral_va() const
Returns va (read only version)
Definition: scalar.h:601
Lorene::Star_QI::w_shift
Vector w_shift
Vector used in the decomposition of shift , following Shibata's prescription [Prog.
Definition: compobj.h:529
Unites
Standard units of space, time and mass.
Lorene::Tbl::set_etat_qcq
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition: tbl.C:361
Lorene::Compobj::kk
Sym_tensor kk
Extrinsic curvature tensor
Definition: compobj.h:153
Lorene::Mg3d::get_type_t
int get_type_t() const
Returns the type of sampling in the direction: SYM : : symmetry with respect to the equatorial pl...
Definition: grilles.h:485
Lorene::Star_QI::nuq
Scalar nuq
Part of the Metric potential = logn generated by the quadratic terms.
Definition: compobj.h:510
Lorene::Coord
Active physical coordinates and mapping derivatives.
Definition: coord.h:90
Lorene::Star_QI::logn
Scalar logn
Logarithm of the lapse N .
Definition: compobj.h:495
Lorene::Star_QI::update_metric
void update_metric()
Computes metric coefficients from known potentials.
Definition: star_QI.C:410
Lorene::Vector::set
Scalar & set(int)
Read/write access to a component.
Definition: vector.C:296
Lorene::Compobj::ener_euler
Scalar ener_euler
Total energy density E in the Eulerian frame.
Definition: compobj.h:144
Lorene::Map::poisson2d
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.
Lorene::Star_QI::fait_nphi
void fait_nphi()
Computes tnphi and nphi from the Cartesian components of the shift, stored in shift .
Definition: star_QI.C:386
Lorene::Cmp
Component of a tensorial field *** DEPRECATED : use class Scalar instead ***.
Definition: cmp.h:446
Lorene::Star_QI::ssjm1_nuf
Scalar ssjm1_nuf
Effective source at the previous step for the resolution of the Poisson equation for nuf by means of ...
Definition: compobj.h:545
Lorene::Tbl
Basic array class.
Definition: tbl.h:161
Lorene::Map::flat_met_cart
const Metric_flat & flat_met_cart() const
Returns the flat metric associated with the Cartesian coordinates and with components expressed in th...
Definition: map.C:331
Lorene::Mg3d::get_nzone
int get_nzone() const
Returns the number of domains.
Definition: grilles.h:448
Lorene::Mg3d::get_nr
int get_nr(int l) const
Returns the number of points in the radial direction ( ) in domain no. l.
Definition: grilles.h:452
Lorene::Param::add_int_mod
void add_int_mod(int &n, int position=0)
Adds the address of a new modifiable int to the list.
Definition: param.C:385
Lorene::log10
Cmp log10(const Cmp &)
Basis 10 logarithm.
Definition: cmp_math.C:322
Lorene::Map::get_bvect_cart
const Base_vect_cart & get_bvect_cart() const
Returns the Cartesian basis associated with the coordinates (x,y,z) of the mapping,...
Definition: map.h:791
Lorene::Map::phi
Coord phi
coordinate centered on the grid
Definition: map.h:720
Lorene::sqrt
Cmp sqrt(const Cmp &)
Square root.
Definition: cmp_math.C:220
Lorene::Scalar::set_dzpuis
void set_dzpuis(int)
Modifies the dzpuis flag.
Definition: scalar.C:808
Lorene::Vector
Tensor field of valence 1.
Definition: vector.h:188
Lorene::Star_QI::ssjm1_khi
Scalar ssjm1_khi
Effective source at the previous step for the resolution of the Poisson equation for the scalar by m...
Definition: compobj.h:569
Lorene::Tensor::trace
Tensor trace(int ind1, int ind2) const
Trace on two different type indices.
Definition: tensor_calculus.C:94
Lorene::Tenseur::set_etat_qcq
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition: tenseur.C:636
Lorene::Scalar::derive_con
const Vector & derive_con(const Metric &gam) const
Returns the "contravariant" derivative of *this with respect to some metric , by raising the index of...
Definition: scalar_deriv.C:402
Lorene::Scalar::test_poisson
Tbl test_poisson(const Scalar &uu, ostream &ostr, bool detail=false) const
Checks if a Poisson equation with *this as a source has been correctly solved.
Definition: scalar_test_poisson.C:60
Lorene::Star_QI::ssjm1_tggg
Scalar ssjm1_tggg
Effective source at the previous step for the resolution of the Poisson equation for tggg .
Definition: compobj.h:561
Lorene::cos
Cmp cos(const Cmp &)
Cosine.
Definition: cmp_math.C:94
Lorene::Param::add_tbl
void add_tbl(const Tbl &ti, int position=0)
Adds the address of a new Tbl to the list.
Definition: param.C:522
Lorene::Scalar::std_spectral_base
virtual void std_spectral_base()
Sets the spectral bases of the Valeur va to the standard ones for a scalar field.
Definition: scalar.C:784
Lorene::Compobj_QI::ak_car
Scalar ak_car
Scalar .
Definition: compobj.h:315
Lorene::Star_QI::nuf
Scalar nuf
Part of the Metric potential = logn generated by the matter terms.
Definition: compobj.h:505
Lorene::Param
Parameter storage.
Definition: param.h:125
Lorene::Compobj_QI::b_car
Scalar b_car
Square of the metric factor B.
Definition: compobj.h:293
Lorene::Star_QI::tggg
Scalar tggg
Metric potential .
Definition: compobj.h:516
Lorene::Scalar::poisson
Scalar poisson() const
Solves the scalar Poisson equation with *this as a source.
Definition: scalar_pde.C:136
Lorene::Param::add_cmp_mod
void add_cmp_mod(Cmp &ti, int position=0)
Adds the address of a new modifiable Cmp to the list.
Definition: param.C:1004
Lorene::Boson_star::update_ener_mom
void update_ener_mom()
Computes the 3+1 components of the energy-momentum tensor (E, P_i and S_{ij}) from the values of the ...
Definition: boson_star.C:211
Lorene::Compobj::mom_euler
Vector mom_euler
Total 3-momentum density in the Eulerian frame.
Definition: compobj.h:147
Lorene::Compobj_QI::a_car
Scalar a_car
Square of the metric factor A.
Definition: compobj.h:287
Lorene::Param::add_int
void add_int(const int &n, int position=0)
Adds the address of a new int to the list.
Definition: param.C:246
Lorene::Compobj::nn
Scalar nn
Lapse function N .
Definition: compobj.h:135
Lorene::Compobj::gamma
Metric gamma
3-metric
Definition: compobj.h:141
Lorene::contract
Tenseur contract(const Tenseur &, int id1, int id2)
Self contraction of two indices of a Tenseur .
Definition: tenseur_operateur.C:279
Lorene::Param::add_tenseur_mod
void add_tenseur_mod(Tenseur &ti, int position=0)
Adds the address of a new modifiable Tenseur to the list.
Definition: param.C:1142