LORENE
star_bin_equilibrium.C
1 
2 /*
3  * Method of class Star_bin to compute an equilibrium configuration
4  *
5  * (see file star.h for documentation).
6  */
7 /*
8  * Copyright (c) 2004 Francois Limousin
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 version 2
14  * as published by the Free Software Foundation.
15  *
16  * LORENE is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
19  * GNU General Public License for more details.
20  *
21  * You should have received a copy of the GNU General Public License
22  * along with LORENE; if not, write to the Free Software
23  * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
24  *
25  */
26 
27 
28 
29 /*
30  * $Id: star_bin_equilibrium.C,v 1.30 2016/12/05 16:18:14 j_novak Exp $
31  * $Log: star_bin_equilibrium.C,v $
32  * Revision 1.30 2016/12/05 16:18:14 j_novak
33  * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
34  *
35  * Revision 1.29 2014/10/13 08:53:38 j_novak
36  * Lorene classes and functions now belong to the namespace Lorene.
37  *
38  * Revision 1.28 2014/10/06 15:13:16 j_novak
39  * Modified #include directives to use c++ syntax.
40  *
41  * Revision 1.27 2006/08/01 14:26:01 f_limousin
42  * Display
43  *
44  * Revision 1.26 2006/05/31 09:26:04 f_limousin
45  * Modif. of the size of the different domains
46  *
47  * Revision 1.25 2006/04/11 14:24:44 f_limousin
48  * New version of the code : improvement of the computation of some
49  * critical sources, estimation of the dirac gauge, helical symmetry...
50  *
51  * Revision 1.24 2005/11/03 13:27:09 f_limousin
52  * Final version for the letter.
53  *
54  * Revision 1.23 2005/09/14 12:48:02 f_limousin
55  * Comment graphical outputs.
56  *
57  * Revision 1.22 2005/09/14 12:30:52 f_limousin
58  * Saving of fields lnq and logn in class Star.
59  *
60  * Revision 1.21 2005/09/13 19:38:31 f_limousin
61  * Reintroduction of the resolution of the equations in cartesian coordinates.
62  *
63  * Revision 1.20 2005/04/08 12:36:44 f_limousin
64  * Just to avoid warnings...
65  *
66  * Revision 1.19 2005/02/24 16:27:21 f_limousin
67  * Add mermax_poisson and relax_poisson in the parameters of the function.
68  *
69  * Revision 1.18 2005/02/24 16:04:13 f_limousin
70  * Change the name of some variables (for instance dcov_logn --> dlogn).
71  * Improve the resolution of the tensorial poisson equation for hh.
72  *
73  * Revision 1.17 2005/02/18 13:14:18 j_novak
74  * Changing of malloc/free to new/delete + suppression of some unused variables
75  * (trying to avoid compilation warnings).
76  *
77  * Revision 1.16 2005/02/17 17:32:53 f_limousin
78  * Change the name of some quantities to be consistent with other classes
79  * (for instance nnn is changed to nn, shift to beta, beta to lnq...)
80  *
81  * Revision 1.15 2005/02/11 18:13:47 f_limousin
82  * Important modification : all the poisson equations for the metric
83  * quantities are now solved on an affine mapping.
84  *
85  * Revision 1.14 2004/12/17 16:23:19 f_limousin
86  * Modif. comments.
87  *
88  * Revision 1.13 2004/06/22 12:49:12 f_limousin
89  * Change qq, qq_auto and qq_comp to beta, beta_auto and beta_comp.
90  *
91  * Revision 1.12 2004/05/27 12:41:00 p_grandclement
92  * correction of some shadowed variables
93  *
94  * Revision 1.11 2004/05/25 14:18:00 f_limousin
95  * Include filters
96  *
97  * Revision 1.10 2004/05/10 10:26:22 f_limousin
98  * Minor changes to avoid warnings in the compilation of Lorene
99  *
100  * Revision 1.9 2004/04/08 16:32:48 f_limousin
101  * The new variable is ln(Q) instead of Q=psi^2*N. It improves the
102  * convergence of the code.
103  *
104  * Revision 1.8 2004/03/25 10:29:26 j_novak
105  * All LORENE's units are now defined in the namespace Unites (in file unites.h).
106  *
107  * Revision 1.7 2004/03/23 09:56:09 f_limousin
108  * Many minor changes
109  *
110  * Revision 1.6 2004/02/27 21:16:32 e_gourgoulhon
111  * Function contract_desal replaced by function contract with
112  * argument desaliasing set to true.
113  *
114  * Revision 1.5 2004/02/27 09:51:51 f_limousin
115  * Many minor changes.
116  *
117  * Revision 1.4 2004/02/21 17:05:13 e_gourgoulhon
118  * Method Scalar::point renamed Scalar::val_grid_point.
119  * Method Scalar::set_point renamed Scalar::set_grid_point.
120  *
121  * Revision 1.3 2004/01/20 15:17:48 f_limousin
122  * First version
123  *
124  *
125  * $Header: /cvsroot/Lorene/C++/Source/Star/star_bin_equilibrium.C,v 1.30 2016/12/05 16:18:14 j_novak Exp $
126  *
127  */
128 
129 // C headers
130 #include <cmath>
131 
132 // Lorene headers
133 #include "cmp.h"
134 #include "tenseur.h"
135 #include "metrique.h"
136 #include "star.h"
137 #include "param.h"
138 #include "graphique.h"
139 #include "utilitaires.h"
140 #include "tensor.h"
141 #include "nbr_spx.h"
142 #include "unites.h"
143 
144 
145 namespace Lorene {
146 void Star_bin::equilibrium(double ent_c, int mermax, int mermax_potvit,
147  int mermax_poisson, double relax_poisson,
148  double relax_potvit, double thres_adapt,
149  Tbl& diff, double om) {
150 
151  // Fundamental constants and units
152  // -------------------------------
153  using namespace Unites ;
154 
155  // Initializations
156  // ---------------
157 
158  const Mg3d* mg = mp.get_mg() ;
159  int nz = mg->get_nzone() ; // total number of domains
160 
161  // The following is required to initialize mp_prev as a Map_et:
162  Map_et& mp_et = dynamic_cast<Map_et&>(mp) ;
163 
164  // Domain and radial indices of points at the surface of the star:
165  int l_b = nzet - 1 ;
166  int i_b = mg->get_nr(l_b) - 1 ;
167  int k_b ;
168  int j_b ;
169 
170  // Value of the enthalpy defining the surface of the star
171  double ent_b = 0 ;
172 
173  // Error indicators
174  // ----------------
175 
176  double& diff_ent = diff.set(0) ;
177  double& diff_vel_pot = diff.set(1) ;
178  double& diff_logn = diff.set(2) ;
179  double& diff_lnq = diff.set(3) ;
180  double& diff_beta_x = diff.set(4) ;
181  double& diff_beta_y = diff.set(5) ;
182  double& diff_beta_z = diff.set(6) ;
183  double& diff_h11 = diff.set(7) ;
184  double& diff_h21 = diff.set(8) ;
185  double& diff_h31 = diff.set(9) ;
186  double& diff_h22 = diff.set(10) ;
187  double& diff_h32 = diff.set(11) ;
188  double& diff_h33 = diff.set(12) ;
189 
190 
191 
192  // Parameters for te function Map_et::adapt
193  // -----------------------------------------
194 
195  Param par_adapt ;
196  int nitermax = 100 ;
197  int niter ;
198  int adapt_flag = 1 ; // 1 = performs the full computation,
199  // 0 = performs only the rescaling by
200  // the factor alpha_r
201  //## int nz_search = nzet + 1 ; // Number of domains for searching the
202  // enthalpy
203  int nz_search = nzet ; // Number of domains for searching the enthalpy
204  // isosurfaces
205 
206  double precis_secant = 1.e-14 ;
207  double alpha_r ;
208  double reg_map = 1. ; // 1 = regular mapping, 0 = contracting mapping
209 
210  Tbl ent_limit(nz) ;
211 
212 
213  par_adapt.add_int(nitermax, 0) ; // maximum number of iterations to
214  // locate zeros by the secant method
215  par_adapt.add_int(nzet, 1) ; // number of domains where the adjustment
216  // to the isosurfaces of ent is to be
217  // performed
218  par_adapt.add_int(nz_search, 2) ; // number of domains to search
219  // the enthalpy isosurface
220  par_adapt.add_int(adapt_flag, 3) ; // 1 = performs the full computation,
221  // 0 = performs only the rescaling by
222  // the factor alpha_r
223  par_adapt.add_int(j_b, 4) ; // theta index of the collocation point
224  // (theta_*, phi_*)
225  par_adapt.add_int(k_b, 5) ; // theta index of the collocation point
226  // (theta_*, phi_*)
227 
228  par_adapt.add_int_mod(niter, 0) ; // number of iterations actually used in
229  // the secant method
230 
231  par_adapt.add_double(precis_secant, 0) ; // required absolute precision in
232  // the determination of zeros by
233  // the secant method
234  par_adapt.add_double(reg_map, 1) ; // 1. = regular mapping,
235  // 0 = contracting mapping
236 
237  par_adapt.add_double(alpha_r, 2) ; // factor by which all the radial
238  // distances will be multiplied
239 
240  par_adapt.add_tbl(ent_limit, 0) ; // array of values of the field ent
241  // to define the isosurfaces.
242 
243 
244  Cmp ssjm1logn (ssjm1_logn) ;
245  Cmp ssjm1lnq (ssjm1_lnq) ;
246  Cmp ssjm1h11 (ssjm1_h11) ;
247  Cmp ssjm1h21 (ssjm1_h21) ;
248  Cmp ssjm1h31 (ssjm1_h31) ;
249  Cmp ssjm1h22 (ssjm1_h22) ;
250  Cmp ssjm1h32 (ssjm1_h32) ;
251  Cmp ssjm1h33 (ssjm1_h33) ;
252 
253 
254  ssjm1logn.set_etat_qcq() ;
255  ssjm1lnq.set_etat_qcq() ;
256  ssjm1h11.set_etat_qcq() ;
257  ssjm1h21.set_etat_qcq() ;
258  ssjm1h31.set_etat_qcq() ;
259  ssjm1h22.set_etat_qcq() ;
260  ssjm1h32.set_etat_qcq() ;
261  ssjm1h33.set_etat_qcq() ;
262 
263 
264  double precis_poisson = 1.e-16 ;
265 
266  // Parameters for the function Scalar::poisson for logn_auto
267  // ---------------------------------------------------------------
268 
269  Param par_logn ;
270 
271  par_logn.add_int(mermax_poisson, 0) ; // maximum number of iterations
272  par_logn.add_double(relax_poisson, 0) ; // relaxation parameter
273  par_logn.add_double(precis_poisson, 1) ; // required precision
274  par_logn.add_int_mod(niter, 0) ; // number of iterations actually used
275  par_logn.add_cmp_mod( ssjm1logn ) ;
276 
277  // Parameters for the function Scalar::poisson for lnq_auto
278  // ---------------------------------------------------------------
279 
280  Param par_lnq ;
281 
282  par_lnq.add_int(mermax_poisson, 0) ; // maximum number of iterations
283  par_lnq.add_double(relax_poisson, 0) ; // relaxation parameter
284  par_lnq.add_double(precis_poisson, 1) ; // required precision
285  par_lnq.add_int_mod(niter, 0) ; // number of iterations actually used -
286  par_lnq.add_cmp_mod( ssjm1lnq ) ;
287 
288  // Parameters for the function Vector::poisson for beta method 2
289  // ---------------------------------------------------------------
290 
291  Param par_beta2 ;
292 
293  par_beta2.add_int(mermax_poisson, 0) ; // maximum number of iterations
294  par_beta2.add_double(relax_poisson, 0) ; // relaxation parameter
295  par_beta2.add_double(precis_poisson, 1) ; // required precision
296  par_beta2.add_int_mod(niter, 0) ; // number of iterations actually used
297 
298  Cmp ssjm1khi (ssjm1_khi) ;
299  Tenseur ssjm1wbeta(mp, 1, CON, mp.get_bvect_cart()) ;
300  ssjm1wbeta.set_etat_qcq() ;
301  for (int i=0; i<3; i++) {
302  ssjm1wbeta.set(i) = Cmp(ssjm1_wbeta(i+1)) ;
303  }
304 
305  par_beta2.add_cmp_mod(ssjm1khi) ;
306  par_beta2.add_tenseur_mod(ssjm1wbeta) ;
307 
308  // Parameters for the function Scalar::poisson for h11_auto
309  // -------------------------------------------------------------
310 
311  Param par_h11 ;
312 
313  par_h11.add_int(mermax_poisson, 0) ; // maximum number of iterations
314  par_h11.add_double(relax_poisson, 0) ; // relaxation parameter
315  par_h11.add_double(precis_poisson, 1) ; // required precision
316  par_h11.add_int_mod(niter, 0) ; // number of iterations actually used
317  par_h11.add_cmp_mod( ssjm1h11 ) ;
318 
319  // Parameters for the function Scalar::poisson for h21_auto
320  // -------------------------------------------------------------
321 
322  Param par_h21 ;
323 
324  par_h21.add_int(mermax_poisson, 0) ; // maximum number of iterations
325  par_h21.add_double(relax_poisson, 0) ; // relaxation parameter
326  par_h21.add_double(precis_poisson, 1) ; // required precision
327  par_h21.add_int_mod(niter, 0) ; // number of iterations actually used
328  par_h21.add_cmp_mod( ssjm1h21 ) ;
329 
330  // Parameters for the function Scalar::poisson for h31_auto
331  // -------------------------------------------------------------
332 
333  Param par_h31 ;
334 
335  par_h31.add_int(mermax_poisson, 0) ; // maximum number of iterations
336  par_h31.add_double(relax_poisson, 0) ; // relaxation parameter
337  par_h31.add_double(precis_poisson, 1) ; // required precision
338  par_h31.add_int_mod(niter, 0) ; // number of iterations actually used
339  par_h31.add_cmp_mod( ssjm1h31 ) ;
340 
341  // Parameters for the function Scalar::poisson for h22_auto
342  // -------------------------------------------------------------
343 
344  Param par_h22 ;
345 
346  par_h22.add_int(mermax_poisson, 0) ; // maximum number of iterations
347  par_h22.add_double(relax_poisson, 0) ; // relaxation parameter
348  par_h22.add_double(precis_poisson, 1) ; // required precision
349  par_h22.add_int_mod(niter, 0) ; // number of iterations actually used
350  par_h22.add_cmp_mod( ssjm1h22 ) ;
351 
352  // Parameters for the function Scalar::poisson for h32_auto
353  // -------------------------------------------------------------
354 
355  Param par_h32 ;
356 
357  par_h32.add_int(mermax_poisson, 0) ; // maximum number of iterations
358  par_h32.add_double(relax_poisson, 0) ; // relaxation parameter
359  par_h32.add_double(precis_poisson, 1) ; // required precision
360  par_h32.add_int_mod(niter, 0) ; // number of iterations actually used
361  par_h32.add_cmp_mod( ssjm1h32 ) ;
362 
363  // Parameters for the function Scalar::poisson for h33_auto
364  // -------------------------------------------------------------
365 
366  Param par_h33 ;
367 
368  par_h33.add_int(mermax_poisson, 0) ; // maximum number of iterations
369  par_h33.add_double(relax_poisson, 0) ; // relaxation parameter
370  par_h33.add_double(precis_poisson, 1) ; // required precision
371  par_h33.add_int_mod(niter, 0) ; // number of iterations actually used
372  par_h33.add_cmp_mod( ssjm1h33 ) ;
373 
374 
375  // External potential
376  // See Eq (99) from Gourgoulhon et al. (2001)
377  // ------------------
378 
379  cout << "logn_comp" << norme(logn_comp) << endl ;
380  cout << "pot_centri" << norme(pot_centri) << endl ;
381  cout << "loggam" << norme(loggam) << endl ;
382  Scalar pot_ext = logn_comp + pot_centri + loggam ;
383 
384  Scalar ent_jm1 = ent ; // Enthalpy at previous step
385 
386  Scalar source_tot(mp) ; // source term in the equation for hij_auto,
387  // logn_auto and beta_auto
388 
389  Vector source_beta(mp, CON, mp.get_bvect_cart()) ; // source term
390  // in the equation for beta_auto
391 
392 
393 
394  //=========================================================================
395  // Start of iteration
396  //=========================================================================
397 
398  for(int mer=0 ; mer<mermax ; mer++ ) {
399 
400  cout << "-----------------------------------------------" << endl ;
401  cout << "step: " << mer << endl ;
402  cout << "diff_ent = " << diff_ent << endl ;
403 
404  //-----------------------------------------------------
405  // Resolution of the elliptic equation for the velocity
406  // scalar potential
407  //-----------------------------------------------------
408 
409  if (irrotational) {
410  diff_vel_pot = velocity_potential(mermax_potvit, precis_poisson,
411  relax_potvit) ;
412 
413  }
414 
415  diff_vel_pot = 0. ; // to avoid the warning
416 
417  //-----------------------------------------------------
418  // Computation of the new radial scale
419  //--------------------------------------------------
420 
421  // alpha_r (r = alpha_r r') is determined so that the enthalpy
422  // takes the requested value ent_b at the stellar surface
423 
424  // Values at the center of the star:
425  double logn_auto_c = logn_auto.val_grid_point(0, 0, 0, 0) ;
426  double pot_ext_c = pot_ext.val_grid_point(0, 0, 0, 0) ;
427 
428  // Search for the reference point (theta_*, phi_*) [notation of
429  // Bonazzola, Gourgoulhon & Marck PRD 58, 104020 (1998)]
430  // at the surface of the star
431  // ------------------------------------------------------------
432  double alpha_r2 = 0 ;
433  for (int k=0; k<mg->get_np(l_b); k++) {
434  for (int j=0; j<mg->get_nt(l_b); j++) {
435 
436  double pot_ext_b = pot_ext.val_grid_point(l_b, k, j, i_b) ;
437  double logn_auto_b = logn_auto.val_grid_point(l_b, k, j, i_b) ;
438  // See Eq (100) from Gourgoulhon et al. (2001)
439  double alpha_r2_jk = ( ent_c - ent_b + pot_ext_c - pot_ext_b) /
440 
441  ( logn_auto_b - logn_auto_c ) ;
442 
443  if (alpha_r2_jk > alpha_r2) {
444  alpha_r2 = alpha_r2_jk ;
445  k_b = k ;
446  j_b = j ;
447  }
448 
449  }
450  }
451 
452  alpha_r = sqrt(alpha_r2) ;
453 
454  cout << "k_b, j_b, alpha_r: " << k_b << " " << j_b << " "
455  << alpha_r << endl ;
456 
457  // New value of logn_auto
458  // ----------------------
459 
460  logn_auto = alpha_r2 * logn_auto ;
461  logn_auto_c = logn_auto.val_grid_point(0, 0, 0, 0) ;
462 
463  //------------------------------------------------------------
464  // Change the values of the inner points of the second domain
465  // by those of the outer points of the first domain
466  //------------------------------------------------------------
467 
469 
470  //------------------------------------------
471  // First integral --> enthalpy in all space
472  // See Eq (98) from Gourgoulhon et al. (2001)
473  //-------------------------------------------
474 
475  ent = (ent_c + logn_auto_c + pot_ext_c) - logn_auto - pot_ext ;
476  cout.precision(8) ;
477  cout << "pot" << norme(pot_ext) << endl ;
478 
479  (ent.set_spectral_va()).smooth(nzet, ent.set_spectral_va()) ;
480 
481  //----------------------------------------------------
482  // Adaptation of the mapping to the new enthalpy field
483  //----------------------------------------------------
484 
485  // Shall the adaptation be performed (cusp) ?
486  // ------------------------------------------
487 
488  double dent_eq = ent.dsdr().val_point(ray_eq(),M_PI/2.,0.) ;
489  double dent_pole = ent.dsdr().val_point(ray_pole(),0.,0.) ;
490  double rap_dent = fabs( dent_eq / dent_pole ) ;
491  cout << "| dH/dr_eq / dH/dr_pole | = " << rap_dent << endl ;
492 
493  if ( rap_dent < thres_adapt ) {
494  adapt_flag = 0 ; // No adaptation of the mapping
495  cout << "******* FROZEN MAPPING *********" << endl ;
496  }
497  else{
498  adapt_flag = 1 ; // The adaptation of the mapping is to be
499  // performed
500  }
501 
502  ent_limit.set_etat_qcq() ;
503  for (int l=0; l<nzet; l++) { // loop on domains inside the star
504  ent_limit.set(l) = ent.val_grid_point(l, k_b, j_b, i_b) ;
505  }
506  ent_limit.set(nzet-1) = ent_b ;
507 
508  Map_et mp_prev = mp_et ;
509 
510  Cmp ent_cmp(ent) ;
511  mp.adapt(ent_cmp, par_adapt) ;
512  ent = ent_cmp ;
513 
514  // Readjustment of the external boundary of domain l=nzet
515  // to keep a fixed ratio with respect to star's surface
516 
517  if (nz>= 5) {
518  double separation = 2. * fabs(mp.get_ori_x()) ;
519  double ray_eqq = ray_eq() ;
520  double ray_eqq_pi = ray_eq_pi() ;
521  double new_rr_out_2 = (separation - ray_eqq) * 0.95 ;
522  double new_rr_out_3 = (separation + ray_eqq_pi) * 1.05 ;
523 
524  double rr_in_1 = mp.val_r(1,-1., M_PI/2, 0.) ;
525  double rr_out_1 = mp.val_r(1, 1., M_PI/2, 0.) ;
526  double rr_out_2 = mp.val_r(2, 1., M_PI/2, 0.) ;
527  double rr_out_3 = mp.val_r(3, 1., M_PI/2, 0.) ;
528 
529  mp.resize(1, 0.5*(new_rr_out_2 + rr_in_1) / rr_out_1) ;
530  mp.resize(2, new_rr_out_2 / rr_out_2) ;
531  mp.resize(3, new_rr_out_3 / rr_out_3) ;
532 
533  for (int dd=4; dd<=nz-2; dd++) {
534  mp.resize(dd, new_rr_out_3 * pow(4., dd-3) /
535  mp.val_r(dd, 1., M_PI/2, 0.)) ;
536  }
537 
538  }
539  else {
540  cout << "too small number of domains" << endl ;
541  }
542 
543  //----------------------------------------------------
544  // Computation of the enthalpy at the new grid points
545  //----------------------------------------------------
546 
547  mp_prev.homothetie(alpha_r) ;
548 
549  Cmp ent_cmp2 (ent) ;
550  mp.reevaluate_symy(&mp_prev, nzet+1, ent_cmp2) ;
551  ent = ent_cmp2 ;
552 
553  double ent_s_max = -1 ;
554  int k_s_max = -1 ;
555  int j_s_max = -1 ;
556  for (int k=0; k<mg->get_np(l_b); k++) {
557  for (int j=0; j<mg->get_nt(l_b); j++) {
558  double xx = fabs( ent.val_grid_point(l_b, k, j, i_b) ) ;
559  if (xx > ent_s_max) {
560  ent_s_max = xx ;
561  k_s_max = k ;
562  j_s_max = j ;
563  }
564  }
565  }
566  cout << "Max. abs(enthalpy) at the boundary between domains nzet-1"
567  << " and nzet : " << endl ;
568  cout << " " << ent_s_max << " reached for k = " << k_s_max <<
569  " and j = " << j_s_max << endl ;
570 
571  //----------------------------------------------------
572  // Equation of state
573  //----------------------------------------------------
574 
575  equation_of_state() ; // computes new values for nbar (n), ener (e)
576  // and press (p) from the new ent (H)
577 
578  //---------------------------------------------------------
579  // Matter source terms in the gravitational field equations
580  //---------------------------------------------------------
581 
582  hydro_euler() ; // computes new values for ener_euler (E),
583  // s_euler (S) and u_euler (U^i)
584 
585 
586  // -------------------------------
587  // AUXILIARY QUANTITIES
588  // -------------------------------
589 
590  // Derivatives of N and logN
591  //--------------------------
592 
593  const Vector dcov_logn_auto = logn_auto.derive_cov(flat) ;
594 
595  Tensor dcovdcov_logn_auto = (logn_auto.derive_cov(flat))
596  .derive_cov(flat) ;
597  dcovdcov_logn_auto.inc_dzpuis() ;
598 
599  // Derivatives of lnq, phi and Q
600  //-------------------------------
601 
602  const Scalar phi (0.5 * (lnq - logn)) ;
603  const Scalar phi_auto (0.5 * (lnq_auto - logn_auto)) ;
604 
605  const Vector dcov_phi_auto = phi_auto.derive_cov(flat) ;
606 
607  const Vector dcov_lnq = 2*dcov_phi + dcov_logn ;
608  const Vector dcon_lnq = 2*dcon_phi + dcon_logn ;
609  const Vector dcov_lnq_auto = lnq_auto.derive_cov(flat) ;
610  Tensor dcovdcov_lnq_auto = dcov_lnq_auto.derive_cov(flat) ;
611  dcovdcov_lnq_auto.inc_dzpuis() ;
612 
613  Scalar qq = exp(lnq) ;
614  qq.std_spectral_base() ;
615  const Vector& dcov_qq = qq.derive_cov(flat) ;
616 
617  const Scalar& divbeta_auto = beta_auto.divergence(flat) ;
618  const Tensor& dcov_beta_auto = beta_auto.derive_cov(flat) ;
619  Tensor dcovdcov_beta_auto = beta_auto.derive_cov(flat)
620  .derive_cov(flat) ;
621  dcovdcov_beta_auto.inc_dzpuis() ;
622 
623 
624  // Derivatives of hij, gtilde...
625  //------------------------------
626 
627  Scalar psi2 (pow(psi4, 0.5)) ;
628  psi2.std_spectral_base() ;
629 
630  const Tensor& dcov_hij = hij.derive_cov(flat) ;
631  const Tensor& dcon_hij = hij.derive_con(flat) ;
632  const Tensor& dcov_hij_auto = hij_auto.derive_cov(flat) ;
633 
634  const Sym_tensor& gtilde_cov = gtilde.cov() ;
635  const Sym_tensor& gtilde_con = gtilde.con() ;
636  const Tensor& dcov_gtilde = gtilde_cov.derive_cov(flat) ;
637 
638  Connection gamijk (gtilde, flat) ;
639  const Tensor& deltaijk = gamijk.get_delta() ;
640 
641  // H^i and its derivatives ( = O in Dirac gauge)
642  // ---------------------------------------------
643 
644  double lambda_dirac = 0. ;
645 
646  const Vector hdirac = lambda_dirac * hij.divergence(flat) ;
647  const Vector hdirac_auto = lambda_dirac * hij_auto.divergence(flat) ;
648 
649  Tensor dcov_hdirac = lambda_dirac * hdirac.derive_cov(flat) ;
650  dcov_hdirac.inc_dzpuis() ;
651  Tensor dcov_hdirac_auto = lambda_dirac * hdirac_auto.derive_cov(flat) ;
652  dcov_hdirac_auto.inc_dzpuis() ;
653  Tensor dcon_hdirac_auto = lambda_dirac * hdirac_auto.derive_con(flat) ;
654  dcon_hdirac_auto.inc_dzpuis() ;
655 
656 
657  //--------------------------------------------------------
658  // Poisson equation for logn_auto (nu_auto)
659  //--------------------------------------------------------
660 
661  // Source
662  //--------
663 
664  int nr = mp.get_mg()->get_nr(0) ;
665  int nt = mp.get_mg()->get_nt(0) ;
666  int np = mp.get_mg()->get_np(0) ;
667 
668  Scalar source1(mp) ;
669  Scalar source2(mp) ;
670  Scalar source3(mp) ;
671  Scalar source4(mp) ;
672  Scalar source5(mp) ;
673  Scalar source6(mp) ;
674  Scalar source7(mp) ;
675  Scalar source8(mp) ;
676 
677  source1 = qpig * psi4 % (ener_euler + s_euler) ;
678 
679  source2 = psi4 % (kcar_auto + kcar_comp) ;
680 
681  source3 = - contract(dcov_logn_auto, 0, dcon_logn, 0, true)
682  - 2. * contract(contract(gtilde_con, 0, dcov_phi, 0),
683  0, dcov_logn_auto, 0, true) ;
684 
685  source4 = - contract(hij, 0, 1, dcovdcov_logn_auto +
686  dcov_logn_auto*dcov_logn, 0, 1) ;
687 
688  source5 = - contract(hdirac, 0, dcov_logn_auto, 0) ;
689 
690  source_tot = source1 + source2 + source3 + source4 + source5 ;
691 
692 
693  cout << "moyenne de la source 1 pour logn_auto" << endl ;
694  cout << norme(source1/(nr*nt*np)) << endl ;
695  cout << "moyenne de la source 2 pour logn_auto" << endl ;
696  cout << norme(source2/(nr*nt*np)) << endl ;
697  cout << "moyenne de la source 3 pour logn_auto" << endl ;
698  cout << norme(source3/(nr*nt*np)) << endl ;
699  cout << "moyenne de la source 4 pour logn_auto" << endl ;
700  cout << norme(source4/(nr*nt*np)) << endl ;
701  cout << "moyenne de la source 5 pour logn_auto" << endl ;
702  cout << norme(source5/(nr*nt*np)) << endl ;
703  cout << "moyenne de la source pour logn_auto" << endl ;
704  cout << norme(source_tot/(nr*nt*np)) << endl ;
705 
706  // Resolution of the Poisson equation
707  // ----------------------------------
708 
709  source_tot.poisson(par_logn, logn_auto) ;
710  ssjm1_logn = ssjm1logn ;
711 
712  cout << "logn_auto" << endl << norme(logn_auto/(nr*nt*np)) << endl ;
713 
714  // Check: has the Poisson equation been correctly solved ?
715  // -----------------------------------------------------
716 
717  Tbl tdiff_logn = diffrel(logn_auto.laplacian(), source_tot) ;
718  cout <<
719  "Relative error in the resolution of the equation for logn_auto : "
720  << endl ;
721  for (int l=0; l<nz; l++) {
722  cout << tdiff_logn(l) << " " ;
723  }
724  cout << endl ;
725  diff_logn = max(abs(tdiff_logn)) ;
726 
727  //--------------------------------------------------------
728  // Poisson equation for lnq_auto
729  //--------------------------------------------------------
730 
731  // Source
732  //--------
733 
734  source1 = qpig * psi4 % s_euler ;
735 
736  source2 = 0.75 * psi4 % (kcar_auto + kcar_comp) ;
737 
738  source3 = - contract(dcon_lnq, 0, dcov_lnq_auto, 0, true) ;
739 
740  source4 = 2. * contract(contract(gtilde_con, 0, dcov_phi, 0), 0,
741  dcov_phi_auto + dcov_logn_auto, 0, true) ;
742 
743  source5 = 0.0625 * contract(gtilde_con, 0, 1, contract(
744  dcov_hij_auto, 0, 1, dcov_gtilde, 0, 1), 0, 1) ;
745 
746  source6 = - 0.125 * contract(gtilde_con, 0, 1, contract(
747  dcov_hij_auto, 0, 1, dcov_gtilde, 0, 2), 0, 1) ;
748 
749  source7 = - contract(hij, 0, 1, dcovdcov_lnq_auto + dcov_lnq_auto *
750  dcov_lnq, 0, 1) ;
751 
752  source8 = - 0.25 * contract(dcov_hdirac_auto, 0, 1)
753  - contract(hdirac, 0, dcov_lnq_auto, 0) ;
754 
755  source_tot = source1 + source2 + source3 + source4 + source5 +
756  source6 + source7 + source8 ;
757 
758 
759  cout << "moyenne de la source 1 pour lnq_auto" << endl ;
760  cout << norme(source1/(nr*nt*np)) << endl ;
761  cout << "moyenne de la source 2 pour lnq_auto" << endl ;
762  cout << norme(source2/(nr*nt*np)) << endl ;
763  cout << "moyenne de la source 3 pour lnq_auto" << endl ;
764  cout << norme(source3/(nr*nt*np)) << endl ;
765  cout << "moyenne de la source 4 pour lnq_auto" << endl ;
766  cout << norme(source4/(nr*nt*np)) << endl ;
767  cout << "moyenne de la source 5 pour lnq_auto" << endl ;
768  cout << norme(source5/(nr*nt*np)) << endl ;
769  cout << "moyenne de la source 6 pour lnq_auto" << endl ;
770  cout << norme(source6/(nr*nt*np)) << endl ;
771  cout << "moyenne de la source 7 pour lnq_auto" << endl ;
772  cout << norme(source7/(nr*nt*np)) << endl ;
773  cout << "moyenne de la source 8 pour lnq_auto" << endl ;
774  cout << norme(source8/(nr*nt*np)) << endl ;
775  cout << "moyenne de la source pour lnq_auto" << endl ;
776  cout << norme(source_tot/(nr*nt*np)) << endl ;
777 
778 
779  // Resolution of the Poisson equation
780  // ----------------------------------
781 
782  source_tot.poisson(par_lnq, lnq_auto) ;
783  ssjm1_lnq = ssjm1lnq ;
784 
785  cout << "lnq_auto" << endl << norme(lnq_auto/(nr*nt*np)) << endl ;
786 
787  // Check: has the Poisson equation been correctly solved
788  // -----------------------------------------------------
789 
790  Tbl tdiff_lnq = diffrel(lnq_auto.laplacian(), source_tot) ;
791  cout <<
792  "Relative error in the resolution of the equation for lnq : "
793  << endl ;
794  for (int l=0; l<nz; l++) {
795  cout << tdiff_lnq(l) << " " ;
796  }
797  cout << endl ;
798  diff_lnq = max(abs(tdiff_lnq)) ;
799 
800  //--------------------------------------------------------
801  // Vector Poisson equation for beta_auto
802  //--------------------------------------------------------
803 
804  // Source
805  //--------
806 
807  Vector source1_beta(mp, CON, mp.get_bvect_cart()) ;
808  Vector source2_beta(mp, CON, mp.get_bvect_cart()) ;
809  Vector source3_beta(mp, CON, mp.get_bvect_cart()) ;
810  Vector source4_beta(mp, CON, mp.get_bvect_cart()) ;
811  Vector source5_beta(mp, CON, mp.get_bvect_cart()) ;
812  Vector source6_beta(mp, CON, mp.get_bvect_cart()) ;
813  Vector source7_beta(mp, CON, mp.get_bvect_cart()) ;
814 
815  source1_beta = (4.*qpig) * nn % psi4
816  %(ener_euler + press) * u_euler ;
817 
818  source2_beta = 2. * nn * contract(tkij_auto, 1,
819  dcov_logn - 6 * dcov_phi, 0) ;
820 
821  source3_beta = - 2. * nn * contract(tkij_auto, 0, 1, deltaijk,
822  1, 2) ;
823 
824  source4_beta = - contract(hij, 0, 1, dcovdcov_beta_auto, 1, 2) ;
825 
826  source5_beta = - 0.3333333333333333 * contract(hij, 1, contract(
827  dcovdcov_beta_auto, 0, 1), 0) ;
828 
829  source6_beta.set_etat_zero() ; //hdirac_auto.derive_lie(omdsdp) ;
830  source6_beta.inc_dzpuis() ;
831 
832  source7_beta = contract(beta, 0, dcov_hdirac_auto, 1) ;
833  + 2./3. * hdirac * divbeta_auto
834  - contract(hdirac, 0, dcov_beta_auto, 1) ;
835 
836  source_beta = source1_beta + source2_beta + source3_beta
837  + source4_beta + source5_beta + source6_beta + source7_beta ;
838 
839 
840  cout << "moyenne de la source 1 pour beta_auto" << endl ;
841  cout << norme(source1_beta(1)/(nr*nt*np)) << endl ;
842  cout << norme(source1_beta(2)/(nr*nt*np)) << endl ;
843  cout << norme(source1_beta(3)/(nr*nt*np)) << endl ;
844  cout << "moyenne de la source 2 pour beta_auto" << endl ;
845  cout << norme(source2_beta(1)/(nr*nt*np)) << endl ;
846  cout << norme(source2_beta(2)/(nr*nt*np)) << endl ;
847  cout << norme(source2_beta(3)/(nr*nt*np)) << endl ;
848  cout << "moyenne de la source 3 pour beta_auto" << endl ;
849  cout << norme(source3_beta(1)/(nr*nt*np)) << endl ;
850  cout << norme(source3_beta(2)/(nr*nt*np)) << endl ;
851  cout << norme(source3_beta(3)/(nr*nt*np)) << endl ;
852  cout << "moyenne de la source 4 pour beta_auto" << endl ;
853  cout << norme(source4_beta(1)/(nr*nt*np)) << endl ;
854  cout << norme(source4_beta(2)/(nr*nt*np)) << endl ;
855  cout << norme(source4_beta(3)/(nr*nt*np)) << endl ;
856  cout << "moyenne de la source 5 pour beta_auto" << endl ;
857  cout << norme(source5_beta(1)/(nr*nt*np)) << endl ;
858  cout << norme(source5_beta(2)/(nr*nt*np)) << endl ;
859  cout << norme(source5_beta(3)/(nr*nt*np)) << endl ;
860  cout << "moyenne de la source 6 pour beta_auto" << endl ;
861  cout << norme(source6_beta(1)/(nr*nt*np)) << endl ;
862  cout << norme(source6_beta(2)/(nr*nt*np)) << endl ;
863  cout << norme(source6_beta(3)/(nr*nt*np)) << endl ;
864  cout << "moyenne de la source 7 pour beta_auto" << endl ;
865  cout << norme(source7_beta(1)/(nr*nt*np)) << endl ;
866  cout << norme(source7_beta(2)/(nr*nt*np)) << endl ;
867  cout << norme(source7_beta(3)/(nr*nt*np)) << endl ;
868  cout << "moyenne de la source pour beta_auto" << endl ;
869  cout << norme(source_beta(1)/(nr*nt*np)) << endl ;
870  cout << norme(source_beta(2)/(nr*nt*np)) << endl ;
871  cout << norme(source_beta(3)/(nr*nt*np)) << endl ;
872 
873  // Resolution of the Poisson equation
874  // ----------------------------------
875 
876  // Filter for the source of beta vector
877 
878  for (int i=1; i<=3; i++) {
879  if (source_beta(i).get_etat() != ETATZERO)
880  source_beta.set(i).filtre(4) ;
881  }
882 
883  for (int i=1; i<=3; i++) {
884  if(source_beta(i).dz_nonzero()) {
885  assert( source_beta(i).get_dzpuis() == 4 ) ;
886  }
887  else{
888  (source_beta.set(i)).set_dzpuis(4) ;
889  }
890  }
891 
892  double lambda = double(1) / double(3) ;
893 
894  Tenseur source_p(mp, 1, CON, mp.get_bvect_cart() ) ;
895  source_p.set_etat_qcq() ;
896  for (int i=0; i<3; i++) {
897  source_p.set(i) = Cmp(source_beta(i+1)) ;
898  }
899 
900  Tenseur vect_auxi (mp, 1, CON, mp.get_bvect_cart()) ;
901  vect_auxi.set_etat_qcq() ;
902  for (int i=0; i<3 ;i++){
903  vect_auxi.set(i) = 0. ;
904  }
905  Tenseur scal_auxi (mp) ;
906  scal_auxi.set_etat_qcq() ;
907  scal_auxi.set().annule_hard() ;
908  scal_auxi.set_std_base() ;
909 
910  Tenseur resu_p(mp, 1, CON, mp.get_bvect_cart() ) ;
911  resu_p.set_etat_qcq() ;
912  for (int i=0; i<3 ;i++)
913  resu_p.set(i).annule_hard() ;
914  resu_p.set_std_base() ;
915 
916  //source_p.poisson_vect(lambda, par_beta2, resu_p, vect_auxi,
917  // scal_auxi) ;
918 
919  source_p.poisson_vect_oohara(lambda, par_beta2, resu_p, scal_auxi) ;
920 
921  for (int i=1; i<=3; i++)
922  beta_auto.set(i) = resu_p(i-1) ;
923 
924  ssjm1_khi = ssjm1khi ;
925  for (int i=0; i<3; i++){
926  ssjm1_wbeta.set(i+1) = ssjm1wbeta(i) ;
927  }
928 
929  cout << "beta_auto_x" << endl << norme(beta_auto(1)/(nr*nt*np))
930  << endl ;
931  cout << "beta_auto_y" << endl << norme(beta_auto(2)/(nr*nt*np))
932  << endl ;
933  cout << "beta_auto_z" << endl << norme(beta_auto(3)/(nr*nt*np))
934  << endl ;
935 
936 
937  // Check: has the equation for beta_auto been correctly solved ?
938  // --------------------------------------------------------------
939 
940  Vector lap_beta = (beta_auto.derive_con(flat)).divergence(flat)
941  + lambda* beta_auto.divergence(flat).derive_con(flat) ;
942 
943  source_beta.dec_dzpuis() ;
944  Tbl tdiff_beta_x = diffrel(lap_beta(1), source_beta(1)) ;
945  Tbl tdiff_beta_y = diffrel(lap_beta(2), source_beta(2)) ;
946  Tbl tdiff_beta_z = diffrel(lap_beta(3), source_beta(3)) ;
947 
948  cout <<
949  "Relative error in the resolution of the equation for beta_auto : "
950  << endl ;
951  cout << "x component : " ;
952  for (int l=0; l<nz; l++) {
953  cout << tdiff_beta_x(l) << " " ;
954  }
955  cout << endl ;
956  cout << "y component : " ;
957  for (int l=0; l<nz; l++) {
958  cout << tdiff_beta_y(l) << " " ;
959  }
960  cout << endl ;
961  cout << "z component : " ;
962  for (int l=0; l<nz; l++) {
963  cout << tdiff_beta_z(l) << " " ;
964  }
965  cout << endl ;
966 
967  diff_beta_x = max(abs(tdiff_beta_x)) ;
968  diff_beta_y = max(abs(tdiff_beta_y)) ;
969  diff_beta_z = max(abs(tdiff_beta_z)) ;
970 
971 
972  if (!conf_flat) {
973 
974  //--------------------------------------------------------
975  // Poisson equation for hij
976  //--------------------------------------------------------
977 
978 
979  // Declaration of all sources
980  //---------------------------
981 
982  Scalar source_tot_hij(mp) ;
983  Tensor source_Sij(mp, 2, CON, mp.get_bvect_cart()) ;
984  Tensor source_Rij(mp, 2, CON, mp.get_bvect_cart()) ;
985  Tensor tens_temp(mp, 2, CON, mp.get_bvect_cart()) ;
986 
987  Tensor source_1 (mp, 2, CON, mp.get_bvect_cart()) ;
988  Tensor source_2 (mp, 2, CON, mp.get_bvect_cart()) ;
989  Tensor source_3a (mp, 2, CON, mp.get_bvect_cart()) ;
990  Tensor source_3b (mp, 2, CON, mp.get_bvect_cart()) ;
991  Tensor source_4 (mp, 2, CON, mp.get_bvect_cart()) ;
992  Tensor source_5 (mp, 2, CON, mp.get_bvect_cart()) ;
993  Tensor source_6 (mp, 2, CON, mp.get_bvect_cart()) ;
994 
995  // Sources
996  //-----------
997 
998  source_1 = contract(dcon_hij, 1, dcov_lnq_auto, 0) ;
999 
1000  source_2 = - contract(dcon_hij, 2, dcov_lnq_auto, 0)
1001  - 2./3. * contract(hdirac, 0, dcov_lnq_auto, 0) * flat.con() ;
1002 
1003  // Lie derivative of A^{ij}
1004  // --------------------------
1005 
1006  Scalar decouple_logn = (logn_auto - 1.e-8)/(logn - 2.e-8) ;
1007 
1008  // Function exp(-(r-r_0)^2/sigma^2)
1009  // --------------------------------
1010 
1011  double r0 = mp.val_r(nz-2, 1, 0, 0) ;
1012  double sigma = 1.*r0 ;
1013 
1014  Scalar rr (mp) ;
1015  rr = mp.r ;
1016 
1017  Scalar ff (mp) ;
1018  ff = exp( -(rr - r0)*(rr - r0)/sigma/sigma ) ;
1019  for (int ii=0; ii<nz-1; ii++)
1020  ff.set_domain(ii) = 1. ;
1021  ff.set_outer_boundary(nz-1, 0) ;
1022  ff.std_spectral_base() ;
1023 
1024  // ff.annule_domain(nz-1) ;
1025  //des_profile(ff, 0, 20, 0, 0) ;
1026 
1027  // Construction of Omega d/dphi
1028  // ----------------------------
1029 
1030  // Construction of D_k \Phi^i
1031  Itbl type (2) ;
1032  type.set(0) = CON ;
1033  type.set(1) = COV ;
1034 
1035  Tensor dcov_omdsdphi (mp, 2, type, mp.get_bvect_cart()) ;
1036  dcov_omdsdphi.set(1,1) = 0. ;
1037  dcov_omdsdphi.set(2,1) = om * ff ;
1038  dcov_omdsdphi.set(3,1) = 0. ;
1039  dcov_omdsdphi.set(2,2) = 0. ;
1040  dcov_omdsdphi.set(3,2) = 0. ;
1041  dcov_omdsdphi.set(3,3) = 0. ;
1042  dcov_omdsdphi.set(1,2) = -om * ff ;
1043  dcov_omdsdphi.set(1,3) = 0. ;
1044  dcov_omdsdphi.set(2,3) = 0. ;
1045  dcov_omdsdphi.std_spectral_base() ;
1046 
1047  source_3a = contract(tkij_auto, 0, dcov_omdsdphi, 1) ;
1048  source_3a.inc_dzpuis() ;
1049 
1050  // Source 3b
1051  // ------------
1052 
1053  Vector omdsdp (mp, CON, mp.get_bvect_cart()) ;
1054  Scalar yya (mp) ;
1055  yya = mp.ya ;
1056  Scalar xxa (mp) ;
1057  xxa = mp.xa ;
1058  Scalar zza (mp) ;
1059  zza = mp.za ;
1060 
1061  if (fabs(mp.get_rot_phi()) < 1e-10){
1062  omdsdp.set(1) = - om * yya * ff ;
1063  omdsdp.set(2) = om * xxa * ff ;
1064  omdsdp.set(3).annule_hard() ;
1065  }
1066  else{
1067  omdsdp.set(1) = om * yya * ff ;
1068  omdsdp.set(2) = - om * xxa * ff ;
1069  omdsdp.set(3).annule_hard() ;
1070  }
1071 
1072  omdsdp.set(1).set_outer_boundary(nz-1, 0) ;
1073  omdsdp.set(2).set_outer_boundary(nz-1, 0) ;
1074  omdsdp.std_spectral_base() ;
1075 
1076  source_3b = - contract(omdsdp, 0, tkij_auto.derive_cov(flat), 2) ;
1077 
1078  // Source 4
1079  // ---------
1080 
1081  source_4 = - tkij_auto.derive_lie(beta) ;
1082  source_4.inc_dzpuis() ;
1083  source_4 += - 2./3. * beta.divergence(flat) * tkij_auto ;
1084 
1085  source_5 = dcon_hdirac_auto ;
1086 
1087  source_6 = - 2./3. * hdirac_auto.divergence(flat) * flat.con() ;
1088  source_6.inc_dzpuis() ;
1089 
1090  // Source terms for Sij
1091  //---------------------
1092 
1093  source_Sij = 8. * nn / psi4 * phi_auto.derive_con(gtilde) *
1094  contract(gtilde_con, 0, dcov_phi, 0) ;
1095 
1096  source_Sij += 4. / psi4 * phi_auto.derive_con(gtilde) *
1097  nn * contract(gtilde_con, 0, dcov_logn, 0) +
1098  4. / psi4 * nn * contract(gtilde_con, 0, dcov_logn, 0) *
1099  phi_auto.derive_con(gtilde) ;
1100 
1101  source_Sij += - nn / (3.*psi4) * gtilde_con *
1102  ( 0.25 * contract(gtilde_con, 0, 1, contract(dcov_hij_auto, 0, 1,
1103  dcov_gtilde, 0, 1), 0, 1)
1104  - 0.5 * contract(gtilde_con, 0, 1, contract(dcov_hij_auto, 0, 1,
1105  dcov_gtilde, 0, 2), 0, 1)) ;
1106 
1107  source_Sij += - 8.*nn / (3.*psi4) * gtilde_con *
1108  contract(dcov_phi_auto, 0, contract(gtilde_con, 0, dcov_phi, 0), 0) ;
1109 
1110  tens_temp = nn / (3.*psi4) * hdirac.divergence(flat)*hij_auto ;
1111  tens_temp.inc_dzpuis() ;
1112 
1113  source_Sij += tens_temp ;
1114 
1115  source_Sij += - 8./(3.*psi4) * contract(dcov_phi_auto, 0,
1116  nn*contract(gtilde_con, 0, dcov_logn, 0), 0) * gtilde_con ;
1117 
1118  source_Sij += 2.*nn* contract(gtilde_cov, 0, 1, tkij_auto *
1119  (tkij_auto+tkij_comp), 1, 3) ;
1120 
1121  source_Sij += - 2. * qpig * nn * ( psi4 * stress_euler
1122  - 0.33333333333333333 * s_euler * gtilde_con ) ;
1123 
1124  source_Sij += - 1./(psi4*psi2) * contract(gtilde_con, 1,
1125  contract(gtilde_con, 1, qq*dcovdcov_lnq_auto +
1126  qq*dcov_lnq_auto*dcov_lnq, 0), 1) ;
1127 
1128  source_Sij += - 0.5/(psi4*psi2) * contract(contract(hij, 1,
1129  dcov_hij_auto, 2), 1, dcov_qq, 0) -
1130  0.5/(psi4*psi2) * contract(contract(dcov_hij_auto, 2,
1131  hij, 1), 1, dcov_qq, 0) ;
1132 
1133  source_Sij += 0.5/(psi4*psi2) * contract(contract(hij, 0,
1134  dcov_hij_auto, 2), 0, dcov_qq, 0) ;
1135 
1136  source_Sij += 1./(3.*psi4*psi2)*contract(gtilde_con, 0, 1,
1137  qq*dcovdcov_lnq_auto + qq*dcov_lnq_auto*dcov_lnq, 0, 1)
1138  *gtilde_con ;
1139 
1140  source_Sij += 1./(3.*psi4*psi2) * contract(hdirac, 0,
1141  dcov_qq, 0)*hij_auto ;
1142 
1143  // Source terms for Rij
1144  //---------------------
1145 
1146  source_Rij = contract(hij, 0, 1, dcov_hij_auto.derive_cov(flat),
1147  2, 3) ;
1148  source_Rij.inc_dzpuis() ;
1149 
1150 
1151  source_Rij += - contract(hij_auto, 1, dcov_hdirac, 1) -
1152  contract(dcov_hdirac, 1, hij_auto, 1) ;
1153 
1154  source_Rij += contract(hdirac, 0, dcov_hij_auto, 2) ;
1155 
1156  source_Rij += - contract(contract(dcov_hij_auto, 1, dcov_hij, 2),
1157  1, 3) ;
1158 
1159  source_Rij += - contract(gtilde_cov, 0, 1, contract(contract(
1160  gtilde_con, 0, dcov_hij_auto, 2), 0, dcov_hij, 2), 1, 3) ;
1161 
1162  source_Rij += contract(contract(contract(contract(gtilde_cov, 0,
1163  dcov_hij_auto, 1), 2, gtilde_con, 1), 0, dcov_hij, 1), 0, 3) +
1164  contract(contract(contract(contract(gtilde_cov, 0,
1165  dcov_hij_auto, 1), 0, dcov_hij, 1), 0, 3), 0, gtilde_con, 1) ;
1166 
1167  source_Rij += 0.5 * contract(gtilde_con*gtilde_con, 1, 3,
1168  contract(dcov_hij_auto, 0, 1, dcov_gtilde, 0, 1), 0, 1) ;
1169 
1170  source_Rij = source_Rij * 0.5 ;
1171 
1172  for(int i=1; i<=3; i++)
1173  for(int j=1; j<=i; j++) {
1174 
1175  source_tot_hij = source_1(i,j) + source_1(j,i)
1176  + source_2(i,j) + 2.*psi4/nn * (
1177  source_4(i,j) - source_Sij(i,j))
1178  - 2.* source_Rij(i,j) +
1179  source_5(i,j) + source_5(j,i) + source_6(i,j) ;
1180  source_tot_hij.dec_dzpuis() ;
1181 
1182  source3 = 2.*psi4/nn * (source_3a(i,j) + source_3a(j,i) +
1183  source_3b(i,j)) ;
1184 
1185  source_tot_hij = source_tot_hij + source3 ;
1186 
1187  //source_tot_hij.inc_dzpuis() ;
1188 
1189  cout << "source_mat" << endl
1190  << norme((- 2. * qpig * nn * ( psi4 * stress_euler
1191  - 0.33333333333333333 * s_euler * gtilde_con ))
1192  (i,j))/(nr*nt*np) << endl ;
1193  cout << "max source_mat" << endl
1194  << max((- 2. * qpig * nn * ( psi4 * stress_euler
1195  - 0.33333333333333333 * s_euler * gtilde_con ))
1196  (i,j)) << endl ;
1197 
1198  cout << "source1" << endl
1199  << norme(source_1(i,j)/(nr*nt*np)) << endl ;
1200  cout << "max source1" << endl
1201  << max(source_1(i,j)) << endl ;
1202  cout << "source2" << endl
1203  << norme(source_2(i,j)/(nr*nt*np)) << endl ;
1204  cout << "max source2" << endl
1205  << max(source_2(i,j)) << endl ;
1206  cout << "source3a" << endl
1207  << norme(source_3a(i,j)/(nr*nt*np)) << endl ;
1208  cout << "max source3a" << endl
1209  << max(source_3a(i,j)) << endl ;
1210  cout << "source3b" << endl
1211  << norme(source_3b(i,j)/(nr*nt*np)) << endl ;
1212  cout << "max source3b" << endl
1213  << max(source_3b(i,j)) << endl ;
1214  cout << "source4" << endl
1215  << norme(source_4(i,j)/(nr*nt*np)) << endl ;
1216  cout << "max source4" << endl
1217  << max(source_4(i,j)) << endl ;
1218  cout << "source5" << endl
1219  << norme(source_5(i,j)/(nr*nt*np)) << endl ;
1220  cout << "max source5" << endl
1221  << max(source_5(i,j)) << endl ;
1222  cout << "source6" << endl
1223  << norme(source_6(i,j)/(nr*nt*np)) << endl ;
1224  cout << "max source6" << endl
1225  << max(source_6(i,j)) << endl ;
1226  cout << "source_Rij" << endl
1227  << norme(source_Rij(i,j)/(nr*nt*np)) << endl ;
1228  cout << "max source_Rij" << endl
1229  << max(source_Rij(i,j)) << endl ;
1230  cout << "source_Sij" << endl
1231  << norme(source_Sij(i,j)/(nr*nt*np)) << endl ;
1232  cout << "max source_Sij" << endl
1233  << max(source_Sij(i,j)) << endl ;
1234  cout << "source_tot" << endl
1235  << norme(source_tot_hij/(nr*nt*np)) << endl ;
1236  cout << "max source_tot" << endl
1237  << max(source_tot_hij) << endl ;
1238 
1239  // Resolution of the Poisson equations and
1240  // Check: has the Poisson equation been correctly solved ?
1241  // -----------------------------------------------------
1242 
1243  if(i==1 && j==1) {
1244 
1245  source_tot_hij.poisson(par_h11, hij_auto.set(1,1)) ;
1246 
1247  Scalar laplac (hij_auto(1,1).laplacian()) ;
1248  laplac.dec_dzpuis() ;
1249  Tbl tdiff_h11 = diffrel(laplac, source_tot_hij) ;
1250  cout << "Relative error in the resolution of the equation for "
1251  << "h11_auto : " << endl ;
1252  for (int l=0; l<nz; l++) {
1253  cout << tdiff_h11(l) << " " ;
1254  }
1255  cout << endl ;
1256  diff_h11 = max(abs(tdiff_h11)) ;
1257  }
1258 
1259  if(i==2 && j==1) {
1260 
1261  source_tot_hij.poisson(par_h21, hij_auto.set(2,1)) ;
1262 
1263  Scalar laplac (hij_auto(2,1).laplacian()) ;
1264  laplac.dec_dzpuis() ;
1265  Tbl tdiff_h21 = diffrel(laplac, source_tot_hij) ;
1266 
1267  cout <<
1268  "Relative error in the resolution of the equation for "
1269  << "h21_auto : " << endl ;
1270  for (int l=0; l<nz; l++) {
1271  cout << tdiff_h21(l) << " " ;
1272  }
1273  cout << endl ;
1274  diff_h21 = max(abs(tdiff_h21)) ;
1275  }
1276 
1277  if(i==3 && j==1) {
1278 
1279  source_tot_hij.poisson(par_h31, hij_auto.set(3,1)) ;
1280 
1281  Scalar laplac (hij_auto(3,1).laplacian()) ;
1282  laplac.dec_dzpuis() ;
1283  Tbl tdiff_h31 = diffrel(laplac, source_tot_hij) ;
1284 
1285  cout <<
1286  "Relative error in the resolution of the equation for "
1287  << "h31_auto : " << endl ;
1288  for (int l=0; l<nz; l++) {
1289  cout << tdiff_h31(l) << " " ;
1290  }
1291  cout << endl ;
1292  diff_h31 = max(abs(tdiff_h31)) ;
1293  }
1294 
1295  if(i==2 && j==2) {
1296 
1297  source_tot_hij.poisson(par_h22, hij_auto.set(2,2)) ;
1298 
1299  Scalar laplac (hij_auto(2,2).laplacian()) ;
1300  laplac.dec_dzpuis() ;
1301  Tbl tdiff_h22 = diffrel(laplac, source_tot_hij) ;
1302 
1303  cout <<
1304  "Relative error in the resolution of the equation for "
1305  << "h22_auto : " << endl ;
1306  for (int l=0; l<nz; l++) {
1307  cout << tdiff_h22(l) << " " ;
1308  }
1309  cout << endl ;
1310  diff_h22 = max(abs(tdiff_h22)) ;
1311  }
1312 
1313  if(i==3 && j==2) {
1314 
1315  source_tot_hij.poisson(par_h32, hij_auto.set(3,2)) ;
1316 
1317  Scalar laplac (hij_auto(3,2).laplacian()) ;
1318  laplac.dec_dzpuis() ;
1319  Tbl tdiff_h32 = diffrel(laplac, source_tot_hij) ;
1320 
1321  cout <<
1322  "Relative error in the resolution of the equation for "
1323  << "h32_auto : " << endl ;
1324  for (int l=0; l<nz; l++) {
1325  cout << tdiff_h32(l) << " " ;
1326  }
1327  cout << endl ;
1328  diff_h32 = max(abs(tdiff_h32)) ;
1329  }
1330 
1331  if(i==3 && j==3) {
1332 
1333  source_tot_hij.poisson(par_h33, hij_auto.set(3,3)) ;
1334 
1335  Scalar laplac (hij_auto(3,3).laplacian()) ;
1336  laplac.dec_dzpuis() ;
1337  Tbl tdiff_h33 = diffrel(laplac, source_tot_hij) ;
1338 
1339  cout <<
1340  "Relative error in the resolution of the equation for "
1341  << "h33_auto : " << endl ;
1342  for (int l=0; l<nz; l++) {
1343  cout << tdiff_h33(l) << " " ;
1344  }
1345  cout << endl ;
1346  diff_h33 = max(abs(tdiff_h33)) ;
1347  }
1348 
1349  }
1350 
1351  cout << "Tenseur hij auto in cartesian coordinates" << endl ;
1352  for (int i=1; i<=3; i++)
1353  for (int j=1; j<=i; j++) {
1354  cout << " Comp. " << i << " " << j << " : " ;
1355  for (int l=0; l<nz; l++){
1356  cout << norme(hij_auto(i,j)/(nr*nt*np))(l) << " " ;
1357  }
1358  cout << endl ;
1359  }
1360  cout << endl ;
1361 
1362  ssjm1_h11 = ssjm1h11 ;
1363  ssjm1_h21 = ssjm1h21 ;
1364  ssjm1_h31 = ssjm1h31 ;
1365  ssjm1_h22 = ssjm1h22 ;
1366  ssjm1_h32 = ssjm1h32 ;
1367  ssjm1_h33 = ssjm1h33 ;
1368 
1369  }
1370 
1371  // End of relativistic equations
1372 
1373 
1374  //-------------------------------------------------
1375  // Relative change in enthalpy
1376  //-------------------------------------------------
1377 
1378  Tbl diff_ent_tbl = diffrel( ent, ent_jm1 ) ;
1379  diff_ent = diff_ent_tbl(0) ;
1380  for (int l=1; l<nzet; l++) {
1381  diff_ent += diff_ent_tbl(l) ;
1382  }
1383  diff_ent /= nzet ;
1384 
1385 
1386  ent_jm1 = ent ;
1387 
1388 
1389  } // End of main loop
1390 
1391  //=========================================================================
1392  // End of iteration
1393  //=========================================================================
1394 
1395 }
1396 
1397 
1398 
1399 
1400 
1401 
1402 
1403 }
Vector dcov_phi
Covariant derivative of the logarithm of the conformal factor.
Definition: star.h:555
Coord xa
Absolute x coordinate.
Definition: map.h:748
bool irrotational
true for an irrotational star, false for a corotating one
Definition: star.h:491
virtual const Sym_tensor & con() const
Read-only access to the contravariant representation.
Definition: metric.C:293
void add_tenseur_mod(Tenseur &ti, int position=0)
Adds the address of a new modifiable Tenseur to the list.
Definition: param.C:1145
Component of a tensorial field *** DEPRECATED : use class Scalar instead ***.
Definition: cmp.h:446
Cmp exp(const Cmp &)
Exponential.
Definition: cmp_math.C:273
int & set(int i)
Read/write of a particular element (index i ) (1D case)
Definition: itbl.h:247
Radial mapping of rather general form.
Definition: map.h:2783
void add_int(const int &n, int position=0)
Adds the address of a new int to the list.
Definition: param.C:249
Map & mp
Mapping associated with the star.
Definition: star.h:180
Scalar ssjm1_lnq
Effective source at the previous step for the resolution of the Poisson equation for lnq_auto...
Definition: star.h:628
Sym_tensor hij_auto
Deviation of the inverse conformal metric from the inverse flat metric generated principally by the ...
Definition: star.h:588
int get_np(int l) const
Returns the number of points in the azimuthal direction ( ) in domain no. l.
Definition: grilles.h:479
Cmp sqrt(const Cmp &)
Square root.
Definition: cmp_math.C:223
void set_std_base()
Set the standard spectal basis of decomposition for each component.
Definition: tenseur.C:1186
Scalar ssjm1_h31
Effective source at the previous step for the resolution of the Poisson equation for h20_auto...
Definition: star.h:651
Lorene prototypes.
Definition: app_hor.h:67
Scalar ssjm1_h11
Effective source at the previous step for the resolution of the Poisson equation for h00_auto...
Definition: star.h:641
Standard units of space, time and mass.
Sym_tensor stress_euler
Spatial part of the stress-energy tensor with respect to the Eulerian observer.
Definition: star.h:212
const Mg3d * get_mg() const
Gives the Mg3d on which the mapping is defined.
Definition: map.h:783
const Scalar & laplacian(int ced_mult_r=4) const
Returns the Laplacian of *this.
Definition: scalar_deriv.C:436
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
Scalar poisson() const
Solves the scalar Poisson equation with *this as a source.
Definition: scalar_pde.C:139
double get_ori_x() const
Returns the x coordinate of the origin.
Definition: map.h:786
virtual void adapt(const Cmp &ent, const Param &par, int nbr=0)=0
Adaptation of the mapping to a given scalar field.
Scalar ssjm1_h33
Effective source at the previous step for the resolution of the Poisson equation for h22_auto...
Definition: star.h:666
const Tensor_sym & derive_con(const Metric &gam) const
Returns the "contravariant" derivative of this with respect to some metric , by raising the last inde...
Scalar logn_comp
Part of the lapse logarithm (gravitational potential at the Newtonian limit) generated principally by...
Definition: star.h:532
Basic integer array class.
Definition: itbl.h:122
Scalar loggam
Logarithm of the Lorentz factor between the fluid and the co-orbiting observer.
Definition: star.h:512
virtual void std_spectral_base()
Sets the spectral bases of the Valeur va to the standard ones for a scalar field. ...
Definition: scalar.C:790
Scalar kcar_auto
Part of the scalar generated by beta_auto, i.e.
Definition: star.h:612
void annule_hard()
Sets the Cmp to zero in a hard way.
Definition: cmp.C:341
Scalar ent
Log-enthalpy.
Definition: star.h:190
double get_rot_phi() const
Returns the angle between the x –axis and X –axis.
Definition: map.h:793
Scalar s_euler
Trace of the stress scalar in the Eulerian frame.
Definition: star.h:201
Vector beta
Shift vector.
Definition: star.h:228
void set(const Map *mp, Mtbl *(*construct)(const Map *))
Semi-constructor from a mapping and a method.
Definition: coord.C:137
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
Tensor field of valence 1.
Definition: vector.h:188
Scalar logn_auto
Part of the lapse logarithm (gravitational potential at the Newtonian limit) generated principally by...
Definition: star.h:527
const Tensor_sym & derive_cov(const Metric &gam) const
Returns the covariant derivative of this with respect to some metric .
void inc_dzpuis()
dzpuis += 1 ;
Definition: tenseur.C:1133
void smooth(int nzet, Valeur &uuva) const
Changes the function *this as a smooth one when there exists a discontinuity between the nucleus and ...
Definition: valeur_smooth.C:80
const Tensor_sym & get_delta() const
Returns the tensor which defines the connection with respect to the flat one: is the difference bet...
Definition: connection.h:271
Tbl diffrel(const Cmp &a, const Cmp &b)
Relative difference between two Cmp (norme version).
Definition: cmp_math.C:507
Scalar ssjm1_h21
Effective source at the previous step for the resolution of the Poisson equation for h10_auto...
Definition: star.h:646
double val_grid_point(int l, int k, int j, int i) const
Returns the value of the field at a specified grid point.
Definition: scalar.h:643
virtual void dec_dzpuis(int dec=1)
Decreases by dec units the value of dzpuis and changes accordingly the values in the compactified ext...
Definition: tensor.C:817
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition: tbl.C:364
int nzet
Number of domains of *mp occupied by the star.
Definition: star.h:183
virtual const Sym_tensor & con() const
Read-only access to the contravariant representation.
Definition: metric_flat.C:156
virtual double val_r(int l, double xi, double theta, double pphi) const =0
Returns the value of the radial coordinate r for a given in a given domain.
Cmp & set()
Read/write for a scalar (see also operator=(const Cmp&) ).
Definition: tenseur.C:840
const Vector & divergence(const Metric &) const
Returns the divergence of this with respect to a Metric .
Definition: sym_tensor.C:352
Tbl norme(const Cmp &)
Sums of the absolute values of all the values of the Cmp in each domain.
Definition: cmp_math.C:484
Sym_tensor hij
Total deviation of the inverse conformal metric from the inverse flat metric.
Definition: star.h:581
Scalar pot_centri
Centrifugal potential.
Definition: star.h:521
Scalar press
Fluid pressure.
Definition: star.h:194
Class Connection.
Definition: connection.h:113
Metric_flat flat
Flat metric defined on the mapping (Spherical components with respect to the mapping of the star) ...
Definition: star.h:562
Scalar lnq_auto
Scalar field generated principally by the star.
Definition: star.h:543
void equilibrium(double ent_c, int mermax, int mermax_potvit, int mermax_poisson, double relax_poisson, double relax_potvit, double thres_adapt, Tbl &diff, double om)
Computes an equilibrium configuration.
Parameter storage.
Definition: param.h:125
void add_tbl(const Tbl &ti, int position=0)
Adds the address of a new Tbl to the list.
Definition: param.C:525
int get_nzone() const
Returns the number of domains.
Definition: grilles.h:465
virtual void homothetie(double lambda)
Sets a new radial scale.
Definition: map_et.C:928
Tbl max(const Cmp &)
Maximum values of a Cmp in each domain.
Definition: cmp_math.C:438
const Scalar & divergence(const Metric &) const
The divergence of this with respect to a Metric .
Definition: vector.C:387
Vector u_euler
Fluid 3-velocity with respect to the Eulerian observer.
Definition: star.h:207
Sym_tensor tkij_comp
Part of the extrinsic curvature tensor generated by beta_comp.
Definition: star.h:606
double ray_eq_pi() const
Coordinate radius at , [r_unit].
Definition: star_global.C:189
Cmp pow(const Cmp &, int)
Power .
Definition: cmp_math.C:351
Tensor handling.
Definition: tensor.h:294
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition: cmp.C:307
Vector dcov_logn
Covariant derivative of the total logarithm of the lapse.
Definition: star.h:535
Tenseur contract(const Tenseur &, int id1, int id2)
Self contraction of two indices of a Tenseur .
const Tensor & derive_cov(const Metric &gam) const
Returns the covariant derivative of this with respect to some metric .
Definition: tensor.C:1011
Scalar logn
Logarithm of the lapse N .
Definition: star.h:222
virtual void resize(int l, double lambda)=0
Rescales the outer boundary of one domain.
double val_point(double r, double theta, double phi) const
Computes the value of the field at an arbitrary point , by means of the spectral expansion.
Definition: scalar.C:896
Metric gtilde
Conformal metric .
Definition: star.h:565
virtual void inc_dzpuis(int inc=1)
Increases by inc units the value of dzpuis and changes accordingly the values in the compactified ext...
Definition: tensor.C:825
virtual void hydro_euler()
Computes the hydrodynamical quantities relative to the Eulerian observer from those in the fluid fram...
Vector dcon_logn
Contravariant derivative of the total logarithm of the lapse.
Definition: star.h:538
double ray_pole() const
Coordinate radius at [r_unit].
Definition: star_global.C:281
int get_nr(int l) const
Returns the number of points in the radial direction ( ) in domain no. l.
Definition: grilles.h:469
virtual const Sym_tensor & cov() const
Read-only access to the covariant representation.
Definition: metric.C:283
Coord ya
Absolute y coordinate.
Definition: map.h:749
const Tensor & derive_con(const Metric &gam) const
Returns the "contravariant" derivative of this with respect to some metric , by raising the last inde...
Definition: tensor.C:1023
Multi-domain grid.
Definition: grilles.h:279
double velocity_potential(int mermax, double precis, double relax)
Computes the non-translational part of the velocity scalar potential by solving the continuity equat...
Tbl & set(int l)
Read/write of the value in a given domain.
Definition: cmp.h:724
virtual void reevaluate_symy(const Map *mp_prev, int nzet, Cmp &uu) const =0
Recomputes the values of a Cmp at the collocation points after a change in the mapping.
const Base_vect_cart & get_bvect_cart() const
Returns the Cartesian basis associated with the coordinates (x,y,z) of the mapping, i.e.
Definition: map.h:809
Scalar & set(const Itbl &ind)
Returns the value of a component (read/write version).
Definition: tensor.C:663
Scalar nn
Lapse function N .
Definition: star.h:225
Scalar ssjm1_h22
Effective source at the previous step for the resolution of the Poisson equation for h11_auto...
Definition: star.h:656
Cmp abs(const Cmp &)
Absolute value.
Definition: cmp_math.C:413
Sym_tensor tkij_auto
Part of the extrinsic curvature tensor generated by beta_auto.
Definition: star.h:600
Coord za
Absolute z coordinate.
Definition: map.h:750
void add_double(const double &x, int position=0)
Adds the the address of a new double to the list.
Definition: param.C:318
bool conf_flat
true if the 3-metric is conformally flat, false for a more general metric.
Definition: star.h:681
double ray_eq() const
Coordinate radius at , [r_unit].
Definition: star_global.C:111
Scalar ssjm1_h32
Effective source at the previous step for the resolution of the Poisson equation for h21_auto...
Definition: star.h:661
const Scalar & dsdr() const
Returns of *this .
Definition: scalar_deriv.C:113
Vector beta_auto
Part of the shift vector generated principally by the star (Spherical components with respect to the ...
Definition: star.h:570
Scalar kcar_comp
Part of the scalar generated by beta_auto and beta_comp, i.e.
Definition: star.h:618
Basic array class.
Definition: tbl.h:164
int get_nt(int l) const
Returns the number of points in the co-latitude direction ( ) in domain no. l.
Definition: grilles.h:474
Valeur & set_spectral_va()
Returns va (read/write version)
Definition: scalar.h:610
Scalar ssjm1_khi
Effective source at the previous step for the resolution of the Poisson equation for khi...
Definition: star.h:634
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
void add_cmp_mod(Cmp &ti, int position=0)
Adds the address of a new modifiable Cmp to the list.
Definition: param.C:1007
void equation_of_state()
Computes the proper baryon and energy density, as well as pressure from the enthalpy.
Definition: star.C:465
const Vector & derive_cov(const Metric &gam) const
Returns the gradient (1-form = covariant vector) of *this
Definition: scalar_deriv.C:390
Scalar ener_euler
Total energy density in the Eulerian frame.
Definition: star.h:198
virtual void dec_dzpuis(int dec=1)
Decreases by dec units the value of dzpuis and changes accordingly the values of the Scalar in the co...
Sym_tensor derive_lie(const Vector &v) const
Computes the Lie derivative of this with respect to some vector field v.
Definition: sym_tensor.C:363
Scalar psi4
Conformal factor .
Definition: star.h:552
Class intended to describe valence-2 symmetric tensors.
Definition: sym_tensor.h:226
Tensor handling *** DEPRECATED : use class Tensor instead ***.
Definition: tenseur.h:304
void add_int_mod(int &n, int position=0)
Adds the address of a new modifiable int to the list.
Definition: param.C:388
Vector dcon_phi
Contravariant derivative of the logarithm of the conformal factor.
Definition: star.h:557
Coord r
r coordinate centered on the grid
Definition: map.h:736
Scalar ssjm1_logn
Effective source at the previous step for the resolution of the Poisson equation for logn_auto...
Definition: star.h:623