LORENE
et_bfrot_equilibre.C
1/*
2 * Method of class Et_rot_bifluid to compute a static spherical configuration.
3 *
4 * (see file etoile.h for documentation).
5 *
6 */
7
8/*
9 * Copyright (c) 2001 Jerome Novak
10 *
11 * This file is part of LORENE.
12 *
13 * LORENE is free software; you can redistribute it and/or modify
14 * it under the terms of the GNU General Public License as published by
15 * the Free Software Foundation; either version 2 of the License, or
16 * (at your option) any later version.
17 *
18 * LORENE is distributed in the hope that it will be useful,
19 * but WITHOUT ANY WARRANTY; without even the implied warranty of
20 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
21 * GNU General Public License for more details.
22 *
23 * You should have received a copy of the GNU General Public License
24 * along with LORENE; if not, write to the Free Software
25 * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
26 *
27 */
28
29
30
31
32/*
33 * $Id: et_bfrot_equilibre.C,v 1.22 2017/10/06 12:36:34 a_sourie Exp $
34 * $Log: et_bfrot_equilibre.C,v $
35 * Revision 1.22 2017/10/06 12:36:34 a_sourie
36 * Cleaning of tabulated 2-fluid EoS class + superfluid rotating star model.
37 *
38 * Revision 1.21 2016/12/05 16:17:52 j_novak
39 * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
40 *
41 * Revision 1.20 2015/06/10 14:39:17 a_sourie
42 * New class Eos_bf_tabul for tabulated 2-fluid EoSs and associated functions for the computation of rotating stars with such EoSs.
43 *
44 * Revision 1.19 2014/10/13 08:52:54 j_novak
45 * Lorene classes and functions now belong to the namespace Lorene.
46 *
47 * Revision 1.18 2014/10/06 15:13:07 j_novak
48 * Modified #include directives to use c++ syntax.
49 *
50 * Revision 1.17 2006/03/13 10:02:27 j_novak
51 * Added things for triaxial perturbations.
52 *
53 * Revision 1.16 2004/09/01 10:56:05 r_prix
54 * added option of converging baryon-mass to equilibrium_bi()
55 *
56 * Revision 1.15 2004/08/30 09:54:20 r_prix
57 * experimental version of Kepler-limit finder for 2-fluid stars
58 *
59 * Revision 1.14 2004/03/25 10:29:03 j_novak
60 * All LORENE's units are now defined in the namespace Unites (in file unites.h).
61 *
62 * Revision 1.13 2003/12/11 12:43:35 r_prix
63 * activated adaptive grid for 2-fluid star (taken from Etoile_rot)
64 *
65 * Revision 1.12 2003/12/04 14:28:26 r_prix
66 * allow for the case of "slow-rot-style" EOS inversion, in which we need to adapt
67 * the inner domain to n_outer=0 instead of mu_outer=0 ...
68 * (this should only be used for comparison to analytic slow-rot solution!)
69 *
70 * Revision 1.11 2003/11/25 12:49:44 j_novak
71 * Modified headers to compile on IRIX. Changed Mapping to be Map_af (speed
72 * enhancement).
73 *
74 * Revision 1.10 2003/11/20 14:01:25 r_prix
75 * changed member names to better conform to Lorene coding standards:
76 * J_euler -> j_euler, EpS_euler -> enerps_euler, Delta_car -> delta_car
77 *
78 * Revision 1.9 2003/11/19 22:01:57 e_gourgoulhon
79 * -- Relaxation on logn and dzeta performed only if mer >= 10.
80 * -- err_grv2 is now evaluated also in the Newtonian case.
81 *
82 * Revision 1.8 2003/11/18 18:38:11 r_prix
83 * use of new member EpS_euler: matter sources in equilibrium() and global quantities
84 * no longer distinguish Newtonian/relativistic, as all terms should have the right limit...
85 *
86 * Revision 1.7 2003/11/17 13:49:43 r_prix
87 * - moved superluminal check into hydro_euler()
88 * - removed some warnings
89 *
90 * Revision 1.6 2003/11/13 12:14:35 r_prix
91 * *) removed all use of etoile-type specific u_euler, press
92 * and use 3+1 components of Tmunu instead
93 *
94 * Revision 1.5 2002/10/16 14:36:35 j_novak
95 * Reorganization of #include instructions of standard C++, in order to
96 * use experimental version 3 of gcc.
97 *
98 * Revision 1.4 2002/04/05 09:09:36 j_novak
99 * The inversion of the EOS for 2-fluids polytrope has been modified.
100 * Some errors in the determination of the surface were corrected.
101 *
102 * Revision 1.3 2002/01/16 15:03:28 j_novak
103 * *** empty log message ***
104 *
105 * Revision 1.2 2002/01/03 15:30:28 j_novak
106 * Some comments modified.
107 *
108 * Revision 1.1.1.1 2001/11/20 15:19:28 e_gourgoulhon
109 * LORENE
110 *
111 * Revision 1.1 2001/06/22 15:40:06 novak
112 * Initial revision
113 *
114 *
115 * $Header: /cvsroot/Lorene/C++/Source/Etoile/et_bfrot_equilibre.C,v 1.22 2017/10/06 12:36:34 a_sourie Exp $
116 *
117 */
118
119// Headers C
120#include <cmath>
121
122// Headers Lorene
123#include "et_rot_bifluid.h"
124#include "param.h"
125#include "unites.h"
126
127#include "graphique.h"
128#include "utilitaires.h"
129
130namespace Lorene {
131
132//-------------------------------------------------------------------------
133//-------------------------------------------------------------------------
134
135// Axial Equilibrium
136
137//-------------------------------------------------------------------------
138//-------------------------------------------------------------------------
139
141(double ent_c, double ent2_c, double omega0, double omega20,
142 const Tbl& ent_limit, const Tbl& ent2_limit, const Itbl& icontrol,
143 const Tbl& control, Tbl& diff,
144 int mer_mass, double mbar1_wanted, double mbar2_wanted, double aexp_mass)
145{
146
147 // Fundamental constants and units
148 // -------------------------------
149 using namespace Unites ;
150
151 // For the display
152 // ---------------
153 char display_bold[]="x[1m" ; display_bold[0] = 27 ;
154 char display_normal[] = "x[0m" ; display_normal[0] = 27 ;
155
156 // Grid parameters
157 // ---------------
158
159 const Mg3d* mg = mp.get_mg() ;
160 int nz = mg->get_nzone() ; // total number of domains
161
162 // The following is required to initialize mp_prev as a Map_af
163 Map_et& mp_et = dynamic_cast<Map_et&>(mp) ; // reference
164
165 // Index of the point at phi=0, theta=pi/2 at the surface of the star:
166 assert(mg->get_type_t() == SYM) ;
167 int l_b = nzet - 1 ;
168 int i_b = mg->get_nr(l_b) - 1 ;
169 int j_b = mg->get_nt(l_b) - 1 ;
170 int k_b = 0 ;
171
172 // Value of the enthalpies defining the surface of each fluid
173 double ent1_b = ent_limit(nzet-1) ;
174 double ent2_b = ent2_limit(nzet-1) ;
175
176 // This value is chosen so that the grid contain both fluids
177// double ent_b = (ent1_b > ent2_b ? ent1_b : ent2_b) ;
178
179 // Parameters to control the iteration
180 // -----------------------------------
181
182 int mer_max = icontrol(0) ;
183 int mer_rot = icontrol(1) ;
184 int mer_change_omega = icontrol(2) ;
185 int mer_fix_omega = icontrol(3) ;
186 int mermax_poisson = icontrol(4) ;
187 int nzadapt = icontrol(5); // number of domains for adaptive grid
188 int kepler_fluid = icontrol(6); // fluid-index for kepler-limit (0=none,3=both)
189 int kepler_wait_steps = icontrol(7);
190 int mer_triax = icontrol(8) ;
191
192 int niter ;
193
194 // Protections:
195 if (mer_change_omega < mer_rot) {
196 cout << "Et_rot_bifluid::equilibrium: mer_change_omega < mer_rot !" << endl ;
197 cout << " mer_change_omega = " << mer_change_omega << endl ;
198 cout << " mer_rot = " << mer_rot << endl ;
199 abort() ;
200 }
201 if (mer_fix_omega < mer_change_omega) {
202 cout << "Et_rot_bifluid::equilibrium: mer_fix_omega < mer_change_omega !"
203 << endl ;
204 cout << " mer_fix_omega = " << mer_fix_omega << endl ;
205 cout << " mer_change_omega = " << mer_change_omega << endl ;
206 abort() ;
207 }
208
209 double precis = control(0) ;
210 double omega_ini = control(1) ;
211 double omega2_ini = control(2) ;
212 double relax = control(3) ;
213 double relax_prev = double(1) - relax ;
214 double relax_poisson = control(4) ;
215 // some additional stuff for adaptive grid:
216 double thres_adapt = control(5) ;
217 double precis_adapt = control(6) ;
218 double kepler_factor = control(7);
219 if (kepler_factor <= 1.0)
220 {
221 cout << "ERROR: Kepler factor has to be greater than 1!!\n";
222 abort();
223 }
224 double ampli_triax = control(8) ;
225
226 // Error indicators
227 // ----------------
228
229 diff.set_etat_qcq() ;
230 double diff_ent ;
231 double& diff_ent1 = diff.set(0) ;
232 double& diff_ent2 = diff.set(1) ;
233 double& diff_nuf = diff.set(2) ;
234 double& diff_nuq = diff.set(3) ;
235 // double& diff_dzeta = diff.set(4) ;
236 // double& diff_ggg = diff.set(5) ;
237 double& diff_shift_x = diff.set(6) ;
238 double& diff_shift_y = diff.set(7) ;
239 double& vit_triax = diff.set(8) ;
240
241 // Parameters for the function Map_et::adapt
242 // -----------------------------------------
243
244 Param par_adapt ;
245 int nitermax = 100 ;
246 int adapt_flag = 1 ; // 1 = performs the full computation,
247 // 0 = performs only the rescaling by
248 // the factor alpha_r
249 int nz_search = nzet + 1 ; // Number of domains for searching the enthalpy
250 // isosurfaces
251 double alpha_r ;
252 double reg_map = 1. ; // 1 = regular mapping, 0 = contracting mapping
253
254 par_adapt.add_int(nitermax, 0) ; // maximum number of iterations to
255 // locate zeros by the secant method
256 par_adapt.add_int(nzadapt, 1) ; // number of domains where the adjustment
257 // to the isosurfaces of ent is to be
258 // performed
259 par_adapt.add_int(nz_search, 2) ; // number of domains to search for
260 // the enthalpy isosurface
261 par_adapt.add_int(adapt_flag, 3) ; // 1 = performs the full computation,
262 // 0 = performs only the rescaling by
263 // the factor alpha_r
264 par_adapt.add_int(j_b, 4) ; // theta index of the collocation point
265 // (theta_*, phi_*)
266 par_adapt.add_int(k_b, 5) ; // theta index of the collocation point
267 // (theta_*, phi_*)
268
269 par_adapt.add_int_mod(niter, 0) ; // number of iterations actually used in
270 // the secant method
271
272 par_adapt.add_double(precis_adapt, 0) ; // required absolute precision in
273 // the determination of zeros by
274 // the secant method
275 par_adapt.add_double(reg_map, 1) ; // 1. = regular mapping, 0 = contracting mapping
276
277 par_adapt.add_double(alpha_r, 2) ; // factor by which all the radial
278 // distances will be multiplied
279
280 par_adapt.add_tbl(ent_limit, 0) ; // array of values of the field ent
281 // to define the isosurfaces.
282
283
284 // Parameters for the function Map_et::poisson for nuf
285 // ----------------------------------------------------
286
287 double precis_poisson = 1.e-16 ;
288
289 Param par_poisson_nuf ;
290 par_poisson_nuf.add_int(mermax_poisson, 0) ; // maximum number of iterations
291 par_poisson_nuf.add_double(relax_poisson, 0) ; // relaxation parameter
292 par_poisson_nuf.add_double(precis_poisson, 1) ; // required precision
293 par_poisson_nuf.add_int_mod(niter, 0) ; // number of iterations actually used
294 par_poisson_nuf.add_cmp_mod( ssjm1_nuf ) ;
295
296 Param par_poisson_nuq ;
297 par_poisson_nuq.add_int(mermax_poisson, 0) ; // maximum number of iterations
298 par_poisson_nuq.add_double(relax_poisson, 0) ; // relaxation parameter
299 par_poisson_nuq.add_double(precis_poisson, 1) ; // required precision
300 par_poisson_nuq.add_int_mod(niter, 0) ; // number of iterations actually used
301 par_poisson_nuq.add_cmp_mod( ssjm1_nuq ) ;
302
303 Param par_poisson_tggg ;
304 par_poisson_tggg.add_int(mermax_poisson, 0) ; // maximum number of iterations
305 par_poisson_tggg.add_double(relax_poisson, 0) ; // relaxation parameter
306 par_poisson_tggg.add_double(precis_poisson, 1) ; // required precision
307 par_poisson_tggg.add_int_mod(niter, 0) ; // number of iterations actually used
308 par_poisson_tggg.add_cmp_mod( ssjm1_tggg ) ;
309 double lambda_tggg ;
310 par_poisson_tggg.add_double_mod( lambda_tggg ) ;
311
312 Param par_poisson_dzeta ;
313 double lbda_grv2 ;
314 par_poisson_dzeta.add_double_mod( lbda_grv2 ) ;
315
316 // Parameters for the function Tenseur::poisson_vect
317 // -------------------------------------------------
318
319 Param par_poisson_vect ;
320
321 par_poisson_vect.add_int(mermax_poisson, 0) ; // maximum number of iterations
322 par_poisson_vect.add_double(relax_poisson, 0) ; // relaxation parameter
323 par_poisson_vect.add_double(precis_poisson, 1) ; // required precision
324 par_poisson_vect.add_cmp_mod( ssjm1_khi ) ;
325 par_poisson_vect.add_tenseur_mod( ssjm1_wshift ) ;
326 par_poisson_vect.add_int_mod(niter, 0) ;
327
328
329 // Initializations
330 // ---------------
331
332 // Initial angular velocities
333 omega = 0 ;
334 omega2 = 0 ;
335
336 double accrois_omega = (omega0 - omega_ini) /
337 double(mer_fix_omega - mer_change_omega) ;
338 double accrois_omega2 = (omega20 - omega2_ini) /
339 double(mer_fix_omega - mer_change_omega) ;
340
341
342 update_metric() ; // update of the metric coefficients
343
344 equation_of_state() ; // update of the densities, pressure, etc...
345
346 hydro_euler() ; // update of the hydro quantities relative to the
347 // Eulerian observer
348
349 // Quantities at the previous step :
350
351 Map_et mp_prev = mp_et;
352 Tenseur ent_prev = ent ;
353 Tenseur ent2_prev = ent2 ;
354 Tenseur logn_prev = logn ;
355 Tenseur dzeta_prev = dzeta ;
356
357 // Creation of uninitialized tensors:
358 Tenseur source_nuf(mp) ; // source term in the equation for nuf
359 Tenseur source_nuq(mp) ; // source term in the equation for nuq
360 Tenseur source_dzf(mp) ; // matter source term in the eq. for dzeta
361 Tenseur source_dzq(mp) ; // quadratic source term in the eq. for dzeta
362 Tenseur source_tggg(mp) ; // source term in the eq. for tggg
363 Tenseur source_shift(mp, 1, CON, mp.get_bvect_cart()) ;
364 // source term for shift
365 Tenseur mlngamma(mp) ; // centrifugal potential
366 Tenseur mlngamma2(mp) ; // centrifugal potential
367
368 Tenseur *outer_ent_p; // pointer to the enthalpy field of the outer fluid
369
370 // Preparations for the Poisson equations:
371 // --------------------------------------
372 if (nuf.get_etat() == ETATZERO) {
373 nuf.set_etat_qcq() ;
374 nuf.set() = 0 ;
375 }
376
377 if (relativistic) {
378 if (nuq.get_etat() == ETATZERO) {
379 nuq.set_etat_qcq() ;
380 nuq.set() = 0 ;
381 }
382
383 if (tggg.get_etat() == ETATZERO) {
384 tggg.set_etat_qcq() ;
385 tggg.set() = 0 ;
386 }
387
388 if (dzeta.get_etat() == ETATZERO) {
389 dzeta.set_etat_qcq() ;
390 dzeta.set() = 0 ;
391 }
392 }
393
394 ofstream fichconv("convergence.d") ; // Output file for diff_ent
395 fichconv << "# diff_ent GRV2 max_triax vit_triax" << endl ;
396
397 ofstream fichfreq("frequency.d") ; // Output file for omega
398 fichfreq << "# f1 [Hz] f2 [Hz]" << endl ;
399
400 ofstream fichevol("evolution.d") ; // Output file for various quantities
401 fichevol << "# r_pole/r_eq ent_c ent2_c" << endl ;
402
403 diff_ent = 1 ;
404 double err_grv2 = 1 ;
405 double max_triax_prev = 0 ; // Triaxial amplitude at previous step
406
407 //=========================================================================
408 // Start of iteration
409 //=========================================================================
410
411 for(int mer=0 ; (diff_ent > precis) && (mer<mer_max) ; mer++ ) {
412
413 cout << "-----------------------------------------------" << endl ;
414 cout << "step: " << mer << endl ;
415 cout << "diff_ent = " << display_bold << diff_ent << display_normal
416 << endl ;
417 cout << "err_grv2 = " << err_grv2 << endl ;
418 fichconv << mer ;
419 fichfreq << mer ;
420 fichevol << mer ;
421
422 if (mer >= mer_rot) {
423
424 if (mer < mer_change_omega) {
425 omega = omega_ini ;
426 omega2 = omega2_ini ;
427 }
428 else {
429 if (mer <= mer_fix_omega) {
430 omega = omega_ini + accrois_omega *
431 (mer - mer_change_omega) ;
432 omega2 = omega2_ini + accrois_omega2 *
433 (mer - mer_change_omega) ;
434 }
435 }
436
437 }
438
439 //-----------------------------------------------
440 // Sources of the Poisson equations
441 //-----------------------------------------------
442
443 // Source for nu
444 // -------------
445 Tenseur beta = log(bbb) ;
446 beta.set_std_base() ;
447
448 // common source term for relativistic and Newtonian ! (enerps_euler has the right limit)
449 source_nuf = qpig * a_car * enerps_euler;
450
451 if (relativistic)
452 source_nuq = ak_car - flat_scalar_prod(logn.gradient_spher(),logn.gradient_spher() + beta.gradient_spher()) ;
453 else
454 source_nuq = 0 ;
455
456 source_nuf.set_std_base() ;
457 source_nuq.set_std_base() ;
458
459 if (relativistic) {
460 // Source for dzeta
461 // ----------------
462 source_dzf = 2 * qpig * a_car * sphph_euler;
463 source_dzf.set_std_base() ;
464
465 source_dzq = 1.5 * ak_car - flat_scalar_prod(logn.gradient_spher(),logn.gradient_spher() ) ;
466 source_dzq.set_std_base() ;
467
468 // Source for tggg
469 // ---------------
470
471 source_tggg = 2 * qpig * nnn * a_car * bbb * (s_euler - sphph_euler);
472 source_tggg.set_std_base() ;
473
474 (source_tggg.set()).mult_rsint() ;
475
476
477 // Source for shift
478 // ----------------
479
480 // Matter term
481 source_shift = (-4*qpig) * nnn * a_car * j_euler;
482
483 // Quadratic terms:
484 Tenseur vtmp = 3 * beta.gradient_spher() - logn.gradient_spher() ;
485 vtmp.change_triad(mp.get_bvect_cart()) ;
486
487 Tenseur squad = 2 * nnn * flat_scalar_prod(tkij, vtmp) ;
488
489 // The addition of matter terms and quadratic terms is performed
490 // component by component because u_euler is contravariant,
491 // while squad is covariant.
492
493 if (squad.get_etat() == ETATQCQ) {
494 for (int i=0; i<3; i++) {
495 source_shift.set(i) += squad(i) ;
496 }
497 }
498
499 source_shift.set_std_base() ;
500 }
501 //----------------------------------------------
502 // Resolution of the Poisson equation for nuf
503 //----------------------------------------------
504
505 source_nuf().poisson(par_poisson_nuf, nuf.set()) ;
506
507// cout << "Test of the Poisson equation for nuf :" << endl ;
508// Tbl err = source_nuf().test_poisson(nuf(), cout, true) ;
509// diff_nuf = err(0, 0) ;
510 diff_nuf = 0 ;
511
512 //---------------------------------------
513 // Triaxial perturbation of nuf
514 //---------------------------------------
515
516 if (mer == mer_triax) {
517
518 if ( mg->get_np(0) == 1 ) {
519 cout <<
520 "Et_rot_bifluid::equilibrium: np must be stricly greater than 1"
521 << endl << " to set a triaxial perturbation !" << endl ;
522 abort() ;
523 }
524
525 const Coord& phi = mp.phi ;
526 const Coord& sint = mp.sint ;
527 Cmp perturb(mp) ;
528 perturb = 1 + ampli_triax * sint*sint * cos(2*phi) ;
529 nuf.set() = nuf() * perturb ;
530
531 nuf.set_std_base() ; // set the bases for spectral expansions
532 // to be the standard ones for a
533 // scalar field
534
535 }
536
537 // Monitoring of the triaxial perturbation
538 // ---------------------------------------
539
540 Valeur& va_nuf = nuf.set().va ;
541 va_nuf.coef() ; // Computes the spectral coefficients
542 double max_triax = 0 ;
543
544 if ( mg->get_np(0) > 1 ) {
545
546 for (int l=0; l<nz; l++) { // loop on the domains
547 for (int j=0; j<mg->get_nt(l); j++) {
548 for (int i=0; i<mg->get_nr(l); i++) {
549
550 // Coefficient of cos(2 phi) :
551 double xcos2p = (*(va_nuf.c_cf))(l, 2, j, i) ;
552
553 // Coefficient of sin(2 phi) :
554 double xsin2p = (*(va_nuf.c_cf))(l, 3, j, i) ;
555
556 double xx = sqrt( xcos2p*xcos2p + xsin2p*xsin2p ) ;
557
558 max_triax = ( xx > max_triax ) ? xx : max_triax ;
559 }
560 }
561 }
562 }
563 cout << "Triaxial part of nuf : " << max_triax << endl ;
564
565 if (relativistic) {
566
567 //----------------------------------------------
568 // Resolution of the Poisson equation for nuq
569 //----------------------------------------------
570
571 source_nuq().poisson(par_poisson_nuq, nuq.set()) ;
572
573// cout << "Test of the Poisson equation for nuq :" << endl ;
574// err = source_nuq().test_poisson(nuq(), cout, true) ;
575// diff_nuq = err(0, 0) ;
576 diff_nuq = 0 ;
577
578 //---------------------------------------------------------
579 // Resolution of the vector Poisson equation for the shift
580 //---------------------------------------------------------
581
582 if (source_shift.get_etat() != ETATZERO) {
583
584 for (int i=0; i<3; i++) {
585 if(source_shift(i).dz_nonzero()) {
586 assert( source_shift(i).get_dzpuis() == 4 ) ;
587 }
588 else{
589 (source_shift.set(i)).set_dzpuis(4) ;
590 }
591 }
592
593 }
594
595 double lambda_shift = double(1) / double(3) ;
596
597 if ( mg->get_np(0) == 1 ) {
598 lambda_shift = 0 ;
599 }
600
601 source_shift.poisson_vect(lambda_shift, par_poisson_vect,
603
604// cout << "Test of the Poisson equation for shift_x :" << endl ;
605// err = source_shift(0).test_poisson(shift(0), cout, true) ;
606// diff_shift_x = err(0, 0) ;
607 diff_shift_x = 0 ;
608
609// cout << "Test of the Poisson equation for shift_y :" << endl ;
610// err = source_shift(1).test_poisson(shift(1), cout, true) ;
611// diff_shift_y = err(0, 0) ;
612 diff_shift_y = 0 ;
613
614 // Computation of tnphi and nphi from the Cartesian components
615 // of the shift
616 // -----------------------------------------------------------
617
618 fait_nphi() ;
619
620 }
621
622
623 //----------------------------------------
624 // Shall we search for the Kepler limit?
625 //----------------------------------------
626 bool kepler = false;
627 bool too_fast = false;
628
629 if ( (kepler_fluid > 0) && (mer > mer_fix_omega + kepler_wait_steps) )
630 {
631 if (kepler_fluid & 0x01)
632 omega *= kepler_factor;
633 if (kepler_fluid & 0x02)
634 omega2 *= kepler_factor;
635 }
636
637
638 // ============================================================
639 kepler = true;
640 while (kepler)
641 {
642
643 // New computation of delta_car, gam_euler, enerps_euler etc...
644 // ------------------------------------------------------
645 hydro_euler() ;
646
647
648 //------------------------------------------------------
649 // First integral of motion
650 //------------------------------------------------------
651
652 // Centrifugal potential :
653 if (relativistic) {
654 mlngamma = - log( gam_euler ) ;
655 mlngamma2 = - log( gam_euler2) ;
656 }
657 else {
658 mlngamma = - 0.5 * uuu*uuu ;
659 mlngamma2 = -0.5 * uuu2*uuu2 ;
660 }
661
662 // Central values of various potentials :
663 double nuf_c = nuf()(0,0,0,0) ;
664 double nuq_c = nuq()(0,0,0,0) ;
665
666 // Scale factor to ensure that the enthalpy is equal to ent_b at
667 // the equator for the "outer" fluid
668 double alpha_r2 = 0;
669
670 int j=j_b;
671
672 // Boundary values of various potentials :
673 double nuf_b = nuf()(l_b, k_b, j, i_b) ;
674 double nuq_b = nuq()(l_b, k_b, j, i_b) ;
675 double mlngamma_b = mlngamma()(l_b, k_b, j, i_b) ;
676 double mlngamma2_b = mlngamma2()(l_b, k_b, j, i_b) ;
677
678
679 // RP: "hack": adapt the radius correctly if using "slow-rot-style" EOS inversion
680 //
681 if ( eos.identify() == 2 ) // only applies to Eos_bf_poly_newt
682 {
683 const Eos_bf_poly_newt &eos0 = dynamic_cast<const Eos_bf_poly_newt&>(eos);
684 if (eos0.get_typeos() == 5)
685 {
686 double vn_b = uuu()(l_b, k_b, j, i_b);
687 double vp_b = uuu2()(l_b, k_b, j, i_b);
688 double D2_b = (vp_b - vn_b)*(vp_b - vn_b);
689 double kdet = eos0.get_kap3() + eos0.get_beta()*D2_b;
690 double kaps1 = kdet / ( eos0.get_kap2() - kdet );
691 double kaps2 = kdet / ( eos0.get_kap1() - kdet );
692
693 ent1_b = kaps1 * ( ent2_c - ent_c - mlngamma2_b + mlngamma_b );
694 ent2_b = kaps2 * ( ent_c - ent2_c - mlngamma_b + mlngamma2_b );
695
696 cout << "**********************************************************************\n";
697 cout << "DEBUG: Rescaling domain for slow-rot-style EOS inversion \n";
698 cout << "DEBUG: ent1_b = " << ent1_b << "; ent2_b = " << ent2_b << endl;
699 cout << "**********************************************************************\n";
700
701 adapt_flag = 0; // don't do adaptive-grid if using slow-rot-style inversion!
702 }
703 }
704
705 double alpha1_r2 = ( ent_c - ent1_b - mlngamma_b + nuq_c - nuq_b) / ( nuf_b - nuf_c ) ;
706 double alpha2_r2 = ( ent2_c - ent2_b - mlngamma2_b + nuq_c - nuq_b) / ( nuf_b - nuf_c ) ;
707
708 cout << "DEBUG: j= "<< j<<" ; alpha1 = " << alpha1_r2 <<" ; alpha2 = " << alpha2_r2 << endl;
709
710 int outer_fluid = (alpha1_r2 > alpha2_r2) ? 1 : 2; // index of 'outer' fluid (at equator!)
711
712 outer_ent_p = (outer_fluid == 1) ? (&ent) : (&ent2);
713
714 alpha_r2 = (outer_fluid == 1) ? alpha1_r2 : alpha2_r2 ;
715
716 alpha_r = sqrt(alpha_r2);
717
718 cout << "alpha_r = " << alpha_r << endl ;
719
720 // Readjustment of nu :
721 // -------------------
722
723 logn = alpha_r2 * nuf + nuq ;
724 double nu_c = logn()(0,0,0,0) ;
725
726
727 // First integral --> enthalpy in all space
728 //-----------------
729
730 ent = (ent_c + nu_c) - logn - mlngamma ;
731 ent2 = (ent2_c + nu_c) - logn - mlngamma2 ;
732
733
734 // now let's try to figure out if we have overstepped the Kepler-limit
735 // (FIXME) we assume that the enthalpy of the _outer_ fluid being negative
736 // inside the star
737 kepler = false;
738 for (int l=0; l<nzet; l++) {
739 int imax = mg->get_nr(l) - 1 ;
740 if (l == l_b) imax-- ; // The surface point is skipped
741 for (int i=0; i<imax; i++) {
742 if ( (*outer_ent_p)()(l, 0, j_b, i) < 0. ) {
743 kepler = true;
744 cout << "(outer) ent < 0 for l, i : " << l << " " << i
745 << " ent = " << (*outer_ent_p)()(l, 0, j_b, i) << endl ;
746 }
747 }
748 }
749
750 if ( kepler )
751 {
752 cout << "**** KEPLERIAN VELOCITY REACHED ****" << endl ;
753 if (kepler_fluid & 0x01)
754 omega /= kepler_factor ; // Omega is decreased
755 if (kepler_fluid & 0x02)
756 omega2 /= kepler_factor;
757
758 cout << "New rotation frequencies : "
759 << "Omega = " << omega/(2.*M_PI) * f_unit << " Hz; "
760 << "Omega2 = " << omega2/(2.*M_PI) * f_unit << endl ;
761
762 too_fast = true;
763 }
764
765 } /* while kepler */
766
767
768 if ( too_fast )
769 { // fact_omega is decreased for the next step
770 kepler_factor = sqrt( kepler_factor ) ;
771 cout << "**** New fact_omega : " << kepler_factor << endl ;
772 }
773 // ============================================================
774
775
776 // Cusp-check: shall the adaptation still be performed?
777 // ------------------------------------------
778 double dent_eq = (*outer_ent_p)().dsdr()(l_b, k_b, j_b, i_b) ;
779 double dent_pole = (*outer_ent_p)().dsdr()(l_b, k_b, 0, i_b) ;
780 double rap_dent = fabs( dent_eq / dent_pole ) ;
781 cout << "| dH/dr_eq / dH/dr_pole | = " << rap_dent << endl ;
782
783 if ( rap_dent < thres_adapt ) {
784 adapt_flag = 0 ; // No adaptation of the mapping
785 cout << "******* FROZEN MAPPING *********" << endl ;
786 }
787
788 // Rescaling of the grid and adaption to (outer) enthalpy surface
789 //---------------------------------------
790 if (adapt_flag && (nzadapt > 0) )
791 {
792 mp_prev = mp_et ;
793
794 mp.adapt( (*outer_ent_p)(), par_adapt) ;
795
796 mp_prev.homothetie(alpha_r) ;
797
798 mp.reevaluate(&mp_prev, nzet+1, ent.set()) ;
799 mp.reevaluate(&mp_prev, nzet+1, ent2.set()) ;
800 }
801 else
802 mp.homothetie (alpha_r);
803
804
805 //----------------------------------------------------
806 // Equation of state
807 //----------------------------------------------------
808
809 equation_of_state() ; // computes new values for nbar1,2 , ener (e)
810 // and press (p) from the new ent,ent2
811
812 //---------------------------------------------------------
813 // Matter source terms in the gravitational field equations
814 //---------------------------------------------------------
815
816 //## Computation of tnphi and nphi from the Cartesian components
817 // of the shift for the test in hydro_euler():
818
819 fait_nphi() ;
820
821 hydro_euler() ; // computes new values for ener_euler (E),
822 // s_euler (S) and u_euler (U^i)
823
824 if (relativistic) {
825
826 //-------------------------------------------------------
827 // 2-D Poisson equation for tggg
828 //-------------------------------------------------------
829
830 mp.poisson2d(source_tggg(), mp.cmp_zero(), par_poisson_tggg, tggg.set()) ;
831
832 //-------------------------------------------------------
833 // 2-D Poisson equation for dzeta
834 //-------------------------------------------------------
835
836 mp.poisson2d(source_dzf(), source_dzq(), par_poisson_dzeta, dzeta.set()) ;
837
838 err_grv2 = lbda_grv2 - 1;
839 cout << "GRV2: " << err_grv2 << endl ;
840
841 }
842 else {
843 err_grv2 = grv2() ;
844 }
845
846
847 //---------------------------------------
848 // Computation of the metric coefficients (except for N^phi)
849 //---------------------------------------
850
851 // Relaxations on nu and dzeta :
852
853 if (mer >= 10) {
854 logn = relax * logn + relax_prev * logn_prev ;
855
856 dzeta = relax * dzeta + relax_prev * dzeta_prev ;
857 }
858
859 // Update of the metric coefficients N, A, B and computation of K_ij :
860
861 update_metric() ;
862
863 //-----------------------
864 // Informations display
865 //-----------------------
866
867 // partial_display(cout) ;
868 fichfreq << " " << omega / (2*M_PI) * f_unit ;
869 fichfreq << " " << omega2 / (2*M_PI) * f_unit ;
870 fichevol << " " << ray_pole() / ray_eq() ;
871 fichevol << " " << ent_c ;
872 fichevol << " " << ent2_c ;
873
874
875 //-----------------------------------------
876 // Convergence towards given baryon masses (if mer_mass > 0)
877 //-----------------------------------------
878
879
880 cout << "DEBUG MODE : mbar1_wanted : " << mbar1_wanted/msol << endl ;
881 cout << "DEBUG MODE : mbar2_wanted : " << mbar2_wanted/msol << endl ;
882
883 // If we want to impose baryonic masses for both fluids.
884 // Be careful, the code acts on mu_n and mu_p (at the center)
885 // -> beta equilibrium can be not verified
886 // CV towards Mbn and Mbp (without beta equilibrium at the center)
887 if (mbar2_wanted > 0 )
888 {
889 if ((mer_mass>0) && (mer > mer_mass)) {
890
891 double xx, xprog, ax, fact;
892
893 // fluid 1
894 xx = mass_b1() / mbar1_wanted - 1. ;
895 cout << "Discrep. baryon mass1 <-> wanted bar. mass1 : " << xx << endl ;
896
897 xprog = ( mer > 2*mer_mass) ? 1. : double(mer - mer_mass)/double(mer_mass) ;
898 xx *= xprog ;
899 ax = 0.5 * ( 2. + xx ) / (1. + xx ) ;
900 fact = pow(ax, aexp_mass) ;
901 cout << "Fluid1: xprog, xx, ax, fact : " << xprog << " " << xx << " " << ax << " " << fact << endl ;
902 ent_c *= fact ;
903
904 // fluid 2
905 xx = mass_b2() / mbar2_wanted - 1. ;
906 cout << "Discrep. baryon mass2 <-> wanted bar. mass2 : " << xx << endl ;
907
908 xprog = ( mer > 2*mer_mass) ? 1. : double(mer - mer_mass)/double(mer_mass) ;
909 xx *= xprog ;
910 ax = 0.5 * ( 2. + xx ) / (1. + xx ) ;
911 fact = pow(ax, aexp_mass) ;
912 cout << "Fluid2: xprog, xx, ax, fact : " << xprog << " " << xx << " " << ax << " " << fact << endl ;
913 ent2_c *= fact ;
914 cout << "H1c = " << ent_c << " H2c = " << ent2_c << endl ;
915 }
916 }
917
918 // CV towards a given GRAVITATIONAL mass (with beta-eq at the center)
919 else if (mbar2_wanted/msol == -1. ) {
920
921 if ((mer_mass>0) && (mer > mer_mass)) {
922
923 double xx, xprog, ax, fact;
924
925 // total mass
926 xx = mass_g() / mbar1_wanted - 1. ; // mbar1_wanted = "mgrav_wanted"
927 cout << "Discrep. baryon mass <-> wanted bar. mass : " << xx << endl ;
928
929 xprog = ( mer > 2*mer_mass) ? 1. : double(mer - mer_mass)/double(mer_mass) ;
930 xx *= xprog ;
931 ax = 0.5 * ( 2. + xx ) / (1. + xx ) ;
932 fact = pow(ax, aexp_mass) ;
933 cout << "Fluid1: xprog, xx, ax, fact : " << xprog << " " << xx << " " << ax << " " << fact << endl ;
934 ent_c *= fact ;
935
936 double m1 = eos.get_m1() ;
937 double m2 = eos.get_m2() ;
938 cout << "m1 = " << m1 << " m2 = " << m2 << endl;
939 ent2_c = ent_c + log(m1/m2); // to ensure beta_equilibrium at the center
940 cout << "DEBUG MODE : ent_c " << ent_c << endl ;
941 cout << "DEBUG MODE : ent2_c " << ent2_c << endl ;
942 cout << "H1c = " << ent_c << " H2c = " << ent2_c << endl ;
943
944 }
945 }
946 // If we want to impose Mb_tot and beta equilibrium (at the center)
947 // In this case : mbar1_wanted = total baryonic mass wanted
948 // mbar2_wanted should be set to 0 and is not used.
949 else {
950
951 if ((mer_mass>0) && (mer > mer_mass)) {
952
953 double xx, xprog, ax, fact;
954
955 // total mass
956 xx = mass_b() / mbar1_wanted - 1. ; // mbar1_wanted = " mbar_wanted"
957 cout << "Discrep. baryon mass <-> wanted bar. mass : " << xx << endl ;
958
959 xprog = ( mer > 2*mer_mass) ? 1. : double(mer - mer_mass)/double(mer_mass) ;
960 xx *= xprog ;
961 ax = 0.5 * ( 2. + xx ) / (1. + xx ) ;
962 fact = pow(ax, aexp_mass) ;
963 cout << "Fluid1: xprog, xx, ax, fact : " << xprog << " " << xx << " " << ax << " " << fact << endl ;
964 ent_c *= fact ;
965
966 double m1 = eos.get_m1() ;
967 double m2 = eos.get_m2() ;
968 cout << "m1 = " << m1 << " m2 = " << m2 << endl;
969 ent2_c = ent_c + log(m1/m2); // to ensure beta_equilibrium at the center
970 cout << "DEBUG MODE : ent_c " << ent_c << endl ;
971 cout << "DEBUG MODE : ent2_c " << ent2_c << endl ;
972 cout << "H1c = " << ent_c << " H2c = " << ent2_c << endl ;
973
974 }
975
976 } /* if mer > mer_mass */
977
978
979 //-------------------------------------------------------------
980 // Relative change in enthalpies with respect to previous step
981 //-------------------------------------------------------------
982
983 Tbl diff_ent_tbl = diffrel( ent(), ent_prev() ) ;
984 diff_ent1 = diff_ent_tbl(0) ;
985 for (int l=1; l<nzet; l++) {
986 diff_ent1 += diff_ent_tbl(l) ;
987 }
988 diff_ent1 /= nzet ;
989 diff_ent_tbl = diffrel( ent2(), ent2_prev() ) ;
990 diff_ent2 = diff_ent_tbl(0) ;
991 for (int l=1; l<nzet; l++) {
992 diff_ent2 += diff_ent_tbl(l) ;
993 }
994 diff_ent2 /= nzet ;
995 diff_ent = 0.5*(diff_ent1 + diff_ent2) ;
996
997 fichconv << " " << log10( fabs(diff_ent) + 1.e-16 ) ;
998 fichconv << " " << log10( fabs(err_grv2) + 1.e-16 ) ;
999 fichconv << " " << log10( fabs(max_triax) + 1.e-16 ) ;
1000
1001 vit_triax = 0 ;
1002 if ( (mer > mer_triax+1) && (max_triax_prev > 1e-13) ) {
1003 vit_triax = (max_triax - max_triax_prev) / max_triax_prev ;
1004 }
1005
1006 fichconv << " " << vit_triax ;
1007
1008 //------------------------------
1009 // Recycling for the next step
1010 //------------------------------
1011
1012 ent_prev = ent ;
1013 ent2_prev = ent2 ;
1014 logn_prev = logn ;
1015 dzeta_prev = dzeta ;
1016 max_triax_prev = max_triax ;
1017
1018 fichconv << endl ;
1019 fichfreq << endl ;
1020 fichevol << endl ;
1021 fichconv.flush() ;
1022 fichfreq.flush() ;
1023 fichevol.flush() ;
1024
1025 } // End of main loop
1026
1027 //=========================================================================
1028 // End of iteration
1029 //=========================================================================
1030
1031 fichconv.close() ;
1032 fichfreq.close() ;
1033 fichevol.close() ;
1034
1035
1036}
1037
1038}
Component of a tensorial field *** DEPRECATED : use class Scalar instead ***.
Definition cmp.h:446
Active physical coordinates and mapping derivatives.
Definition coord.h:90
Analytic equation of state for two fluids (Newtonian case).
double get_kap1() const
Returns the pressure coefficient [unit: ], where .
double get_kap3() const
Returns the pressure coefficient [unit: ], where .
double get_beta() const
Returns the coefficient [unit: ], where .
double get_kap2() const
Returns the pressure coefficient [unit: ], where .
void equilibrium_bi(double ent_c, double ent_c2, double omega0, double omega20, const Tbl &ent_limit, const Tbl &ent2_limit, const Itbl &icontrol, const Tbl &control, Tbl &diff, int mer_mass, double mbar1_wanted, double mbar2_wanted, double aexp_mass)
Computes an equilibrium configuration.
Tenseur gam_euler2
Lorentz factor between the fluid 2 and Eulerian observers.
Tenseur enerps_euler
the combination : useful because in the Newtonian limit .
Tenseur sphph_euler
The component of the stress tensor .
virtual double mass_b() const
Total Baryon mass.
Tenseur j_euler
Total angular momentum (flat-space!) 3-vector , which is related to of the "3+1" decomposition,...
virtual void equation_of_state()
Computes the proper baryon and energy densities, as well as pressure and the coefficients Knn,...
virtual double mass_g() const
Gravitational mass.
double omega2
Rotation angular velocity for fluid 2 ([f_unit] ).
const Eos_bifluid & eos
Equation of state for two-fluids model.
double mass_b1() const
Baryon mass of fluid 1.
Tenseur uuu2
Norm of the (fluid no.2) 3-velocity with respect to the eulerian observer.
double mass_b2() const
Baryon mass of fluid 2.
virtual double grv2() const
Error on the virial identity GRV2.
Tenseur ent2
Log-enthalpy for the second fluid.
virtual void hydro_euler()
Computes the hydrodynamical quantities relative to the Eulerian observer from those in the fluid fram...
Tenseur ssjm1_wshift
Effective source at the previous step for the resolution of the vector Poisson equation for .
Definition etoile.h:1628
Tenseur uuu
Norm of u_euler.
Definition etoile.h:1521
double omega
Rotation angular velocity ([f_unit] ).
Definition etoile.h:1504
Tenseur & logn
Metric potential = logn_auto.
Definition etoile.h:1524
Tenseur nuq
Part of the Metric potential = logn generated by the quadratic terms.
Definition etoile.h:1534
Tenseur khi_shift
Scalar used in the decomposition of shift , following Shibata's prescription [Prog.
Definition etoile.h:1563
Tenseur tggg
Metric potential .
Definition etoile.h:1540
Tenseur nuf
Part of the Metric potential = logn generated by the matter terms.
Definition etoile.h:1529
Cmp ssjm1_tggg
Effective source at the previous step for the resolution of the Poisson equation for tggg .
Definition etoile.h:1611
Tenseur bbb
Metric factor B.
Definition etoile.h:1507
void update_metric()
Computes metric coefficients from known potentials.
Tenseur ak_car
Scalar .
Definition etoile.h:1589
Tenseur & dzeta
Metric potential = beta_auto.
Definition etoile.h:1537
Cmp ssjm1_nuf
Effective source at the previous step for the resolution of the Poisson equation for nuf by means of ...
Definition etoile.h:1595
void fait_nphi()
Computes tnphi and nphi from the Cartesian components of the shift, stored in shift .
Definition etoile_rot.C:803
Cmp ssjm1_khi
Effective source at the previous step for the resolution of the Poisson equation for the scalar by m...
Definition etoile.h:1619
Tenseur_sym tkij
Tensor related to the extrinsic curvature tensor by .
Definition etoile.h:1570
Cmp ssjm1_nuq
Effective source at the previous step for the resolution of the Poisson equation for nuq by means of ...
Definition etoile.h:1601
Tenseur w_shift
Vector used in the decomposition of shift , following Shibata's prescription [Prog.
Definition etoile.h:1553
int nzet
Number of domains of *mp occupied by the star.
Definition etoile.h:435
double ray_eq() const
Coordinate radius at , [r_unit].
Tenseur nnn
Total lapse function.
Definition etoile.h:512
Tenseur gam_euler
Lorentz factor between the fluid and Eulerian observers.
Definition etoile.h:474
Map & mp
Mapping associated with the star.
Definition etoile.h:432
bool relativistic
Indicator of relativity: true for a relativistic star, false for a Newtonian one.
Definition etoile.h:440
Tenseur shift
Total shift vector.
Definition etoile.h:515
Tenseur s_euler
Trace of the stress tensor in the Eulerian frame.
Definition etoile.h:471
Tenseur ent
Log-enthalpy (relativistic case) or specific enthalpy (Newtonian case).
Definition etoile.h:460
Tenseur a_car
Total conformal factor .
Definition etoile.h:518
double ray_pole() const
Coordinate radius at [r_unit].
Basic integer array class.
Definition itbl.h:122
Radial mapping of rather general form.
Definition map.h:2770
virtual void homothetie(double lambda)
Sets a new radial scale.
Definition map_et.C:928
Multi-domain grid.
Definition grilles.h:279
int get_type_t() const
Returns the type of sampling in the direction: SYM : : symmetry with respect to the equatorial pl...
Definition grilles.h:502
int get_np(int l) const
Returns the number of points in the azimuthal direction ( ) in domain no. l.
Definition grilles.h:479
int get_nt(int l) const
Returns the number of points in the co-latitude direction ( ) in domain no. l.
Definition grilles.h:474
int get_nzone() const
Returns the number of domains.
Definition grilles.h:465
int get_nr(int l) const
Returns the number of points in the radial direction ( ) in domain no. l.
Definition grilles.h:469
Parameter storage.
Definition param.h:125
void add_double(const double &x, int position=0)
Adds the the address of a new double to the list.
Definition param.C:318
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 add_double_mod(double &x, int position=0)
Adds the address of a new modifiable double to the list.
Definition param.C:456
void add_int_mod(int &n, int position=0)
Adds the address of a new modifiable int to the list.
Definition param.C:388
void add_tenseur_mod(Tenseur &ti, int position=0)
Adds the address of a new modifiable Tenseur to the list.
Definition param.C:1145
void add_int(const int &n, int position=0)
Adds the address of a new int to the list.
Definition param.C:249
void add_tbl(const Tbl &ti, int position=0)
Adds the address of a new Tbl to the list.
Definition param.C:525
Basic array class.
Definition tbl.h:161
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition tbl.C:364
double & set(int i)
Read/write of a particular element (index i) (1D case).
Definition tbl.h:281
Tensor handling *** DEPRECATED : use class Tensor instead ***.
Definition tenseur.h:304
Cmp & set()
Read/write for a scalar (see also operator=(const Cmp&) ).
Definition tenseur.C:830
void set_std_base()
Set the standard spectal basis of decomposition for each component.
Definition tenseur.C:1176
const Tenseur & gradient_spher() const
Returns the gradient of *this (Spherical coordinates) (scalar field only).
Definition tenseur.C:1554
void change_triad(const Base_vect &new_triad)
Sets a new vectorial basis (triad) of decomposition and modifies the components accordingly.
Definition tenseur.C:674
void poisson_vect(double lambda, Param &par, Tenseur &shift, Tenseur &vect, Tenseur &scal) const
Solves the vectorial Poisson equation : .
int get_etat() const
Returns the logical state.
Definition tenseur.h:710
Values and coefficients of a (real-value) function.
Definition valeur.h:297
Mtbl_cf * c_cf
Coefficients of the spectral expansion of the function.
Definition valeur.h:312
void coef() const
Computes the coeffcients of *this.
Cmp sqrt(const Cmp &)
Square root.
Definition cmp_math.C:223
Cmp log10(const Cmp &)
Basis 10 logarithm.
Definition cmp_math.C:325
Tbl diffrel(const Cmp &a, const Cmp &b)
Relative difference between two Cmp (norme version).
Definition cmp_math.C:507
Cmp pow(const Cmp &, int)
Power .
Definition cmp_math.C:351
Cmp cos(const Cmp &)
Cosine.
Definition cmp_math.C:97
Cmp log(const Cmp &)
Neperian logarithm.
Definition cmp_math.C:299
Tenseur flat_scalar_prod(const Tenseur &t1, const Tenseur &t2)
Scalar product of two Tenseur when the metric is : performs the contraction of the last index of t1 w...
Lorene prototypes.
Definition app_hor.h:67
virtual void dsdr(const Cmp &ci, Cmp &resu) const =0
Computes of a Cmp .
Coord phi
coordinate centered on the grid
Definition map.h:732
Coord sint
Definition map.h:733
Standard units of space, time and mass.