LORENE
et_bin_equil_regu.C
1/*
2 * Method of class Etoile_bin to compute an equilibrium configuration
3 * by regularizing source.
4 *
5 * (see file etoile.h for documentation).
6 *
7 */
8
9/*
10 * Copyright (c) 2000-2001 Eric Gourgoulhon
11 * Copyright (c) 2000-2001 Keisuke Taniguchi
12 *
13 * This file is part of LORENE.
14 *
15 * LORENE is free software; you can redistribute it and/or modify
16 * it under the terms of the GNU General Public License as published by
17 * the Free Software Foundation; either version 2 of the License, or
18 * (at your option) any later version.
19 *
20 * LORENE is distributed in the hope that it will be useful,
21 * but WITHOUT ANY WARRANTY; without even the implied warranty of
22 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
23 * GNU General Public License for more details.
24 *
25 * You should have received a copy of the GNU General Public License
26 * along with LORENE; if not, write to the Free Software
27 * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
28 *
29 */
30
31
32char et_bin_equil_regu_C[] = "$Header: /cvsroot/Lorene/C++/Source/Etoile/et_bin_equil_regu.C,v 1.8 2014/10/13 08:52:55 j_novak Exp $" ;
33
34/*
35 * $Id: et_bin_equil_regu.C,v 1.8 2014/10/13 08:52:55 j_novak Exp $
36 * $Log: et_bin_equil_regu.C,v $
37 * Revision 1.8 2014/10/13 08:52:55 j_novak
38 * Lorene classes and functions now belong to the namespace Lorene.
39 *
40 * Revision 1.7 2014/10/06 15:13:08 j_novak
41 * Modified #include directives to use c++ syntax.
42 *
43 * Revision 1.6 2009/06/15 09:26:57 k_taniguchi
44 * Improved the rescaling of the domains.
45 *
46 * Revision 1.5 2004/03/25 10:29:03 j_novak
47 * All LORENE's units are now defined in the namespace Unites (in file unites.h).
48 *
49 * Revision 1.4 2003/09/01 06:48:08 k_taniguchi
50 * Change of the domain which should be resized.
51 *
52 * Revision 1.3 2003/08/31 05:35:38 k_taniguchi
53 * Addition of the specification of the domain
54 * which is resized.
55 *
56 * Revision 1.2 2002/12/11 12:51:26 k_taniguchi
57 * Change the multiplication "*" to "%"
58 * and flat_scalar_prod to flat_scalar_prod_desal.
59 *
60 * Revision 1.1.1.1 2001/11/20 15:19:28 e_gourgoulhon
61 * LORENE
62 *
63 * Revision 2.17 2001/08/07 09:49:00 keisuke
64 * Change of the method to set the longest radius of a star
65 * on the first domain.
66 * Addition of a new argument in Etoile_bin::equil_regular : Tbl fact.
67 *
68 * Revision 2.16 2001/06/22 08:54:53 keisuke
69 * Set the inner values of the second domain of ent
70 * by using the outer ones of the first domain.
71 *
72 * Revision 2.15 2001/05/17 12:22:26 keisuke
73 * Change of the method to calculate chi from setting position in map
74 * to val_point.
75 *
76 * Revision 2.14 2001/02/07 09:47:28 eric
77 * unsgam1 est desormais donne par Eos::der_nbar_ent (cas newtonien)
78 * ou Eos::der_ener_ent (cas relativiste).
79 *
80 * Revision 2.13 2001/01/16 17:02:32 keisuke
81 * *** empty log message ***
82 *
83 * Revision 2.12 2001/01/16 16:58:08 keisuke
84 * Change the method to set the values on the surface.
85 *
86 * Revision 2.11 2001/01/10 16:45:34 keisuke
87 * Set the inner values of the second domain of logn_auto
88 * by using the outer ones of the first domain.
89 *
90 * Revision 2.10 2000/12/20 10:33:14 eric
91 * Changement important : nz_search = nzet ---> nz_search = nzet + 1
92 *
93 * Revision 2.9 2000/10/25 14:01:03 keisuke
94 * Modif de Map_et::adapt: on y rentre desormais avec nz_search
95 * (dans le cas present nz_search = nzet).
96 *
97 * Revision 2.8 2000/10/06 15:29:01 keisuke
98 * Change poisson_vect into poisson_vect_regu.
99 *
100 * Revision 2.7 2000/09/25 15:01:10 keisuke
101 * Suppress "int" from the declaration of k_div.
102 *
103 * Revision 2.6 2000/09/22 15:51:39 keisuke
104 * d_logn_auto est desormais calcule en dehors (dans update_metric).
105 *
106 * Revision 2.5 2000/09/13 09:50:33 keisuke
107 * Minor change on change_triad.
108 *
109 * Revision 2.4 2000/09/08 15:57:31 keisuke
110 * Change the basis of d_logn_auto_div from the spherical coordinate
111 * to the Cartesian one with respect to ref_triad.
112 *
113 * Revision 2.3 2000/09/07 15:47:19 keisuke
114 * Minor change.
115 *
116 * Revision 2.2 2000/09/07 15:43:41 keisuke
117 * Add a new argument in poisson_regular and suppress logn_auto_total.
118 *
119 * Revision 2.1 2000/08/29 14:01:43 keisuke
120 * Modify the arguments of poisson_regular.
121 *
122 * Revision 2.0 2000/08/29 11:39:02 eric
123 * Version provisoire.
124 *
125 *
126 * $Header: /cvsroot/Lorene/C++/Source/Etoile/et_bin_equil_regu.C,v 1.8 2014/10/13 08:52:55 j_novak Exp $
127 *
128 */
129
130// Headers C
131#include <cmath>
132
133// Headers Lorene
134#include "etoile.h"
135#include "param.h"
136#include "eos.h"
137#include "utilitaires.h"
138#include "unites.h"
139
140namespace Lorene {
141
142//********************************************************************
143
144void Etoile_bin::equil_regular(double ent_c, int mermax, int mermax_poisson,
145 double relax_poisson, int mermax_potvit,
146 double relax_potvit, double thres_adapt,
147 const Tbl& fact_resize, Tbl& diff) {
148
149 // Fundamental constants and units
150 // -------------------------------
151 using namespace Unites ;
152
153 // Initializations
154 // ---------------
155
156 k_div = 2 ; // Regularity parameter for poisson_regular
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_beta = diff.set(3) ;
180 double& diff_shift_x = diff.set(4) ;
181 double& diff_shift_y = diff.set(5) ;
182 double& diff_shift_z = diff.set(6) ;
183
184 // Parameters for the function Map_et::adapt
185 // -----------------------------------------
186
187 Param par_adapt ;
188 int nitermax = 100 ;
189 int niter ;
190 int adapt_flag = 1 ; // 1 = performs the full computation,
191 // 0 = performs only the rescaling by
192 // the factor alpha_r
193 //## int nz_search = nzet + 1 ; // Number of domains for searching the enthalpy
194 int nz_search = nzet ; // Number of domains for searching the enthalpy
195 // isosurfaces
196
197 double precis_secant = 1.e-14 ;
198 double alpha_r ;
199 double reg_map = 1. ; // 1 = regular mapping, 0 = contracting mapping
200
201 Tbl ent_limit(nz) ;
202
203 par_adapt.add_int(nitermax, 0) ; // maximum number of iterations to
204 // locate zeros by the secant method
205 par_adapt.add_int(nzet, 1) ; // number of domains where the adjustment
206 // to the isosurfaces of ent is to be
207 // performed
208 par_adapt.add_int(nz_search, 2) ; // number of domains to search for
209 // the enthalpy isosurface
210 par_adapt.add_int(adapt_flag, 3) ; // 1 = performs the full computation,
211 // 0 = performs only the rescaling by
212 // the factor alpha_r
213 par_adapt.add_int(j_b, 4) ; // theta index of the collocation point
214 // (theta_*, phi_*)
215 par_adapt.add_int(k_b, 5) ; // theta index of the collocation point
216 // (theta_*, phi_*)
217 par_adapt.add_int_mod(niter, 0) ; // number of iterations actually
218 // used in the secant method
219 par_adapt.add_double(precis_secant, 0) ; // required absolute precision in
220 // the determination of zeros by
221 // the secant method
222 par_adapt.add_double(reg_map, 1) ; // 1. = regular mapping,
223 // 0 = contracting mapping
224 par_adapt.add_double(alpha_r, 2) ; // factor by which all the radial
225 // distances will be multiplied
226 par_adapt.add_tbl(ent_limit, 0) ; // array of values of the field ent
227 // to define the isosurfaces.
228
229 // Parameters for the function Map_et::poisson_regular for logn_auto
230 // -----------------------------------------------------------------
231
232 double precis_poisson = 1.e-16 ;
233
234 Param par_poisson1 ;
235
236 par_poisson1.add_int(mermax_poisson, 0) ; // maximum number of iterations
237 par_poisson1.add_double(relax_poisson, 0) ; // relaxation parameter
238 par_poisson1.add_double(precis_poisson, 1) ; // required precision
239 par_poisson1.add_int_mod(niter, 0) ; // number of iterations actually
240 // used
241 par_poisson1.add_cmp_mod( ssjm1_logn ) ;
242
243 // Parameters for the function Map_et::poisson for beta_auto
244 // ---------------------------------------------------------
245
246 Param par_poisson2 ;
247
248 par_poisson2.add_int(mermax_poisson, 0) ; // maximum number of iterations
249 par_poisson2.add_double(relax_poisson, 0) ; // relaxation parameter
250 par_poisson2.add_double(precis_poisson, 1) ; // required precision
251 par_poisson2.add_int_mod(niter, 0) ; // number of iterations actually
252 // used
253 par_poisson2.add_cmp_mod( ssjm1_beta ) ;
254
255 // Parameters for the function Tenseur::poisson_vect_regu
256 // ------------------------------------------------------
257
258 Param par_poisson_vect ;
259
260 par_poisson_vect.add_int(mermax_poisson, 0) ; // maximum number of
261 // iterations
262 par_poisson_vect.add_double(relax_poisson, 0) ; // relaxation parameter
263 par_poisson_vect.add_double(precis_poisson, 1) ; // required precision
264 par_poisson_vect.add_cmp_mod( ssjm1_khi ) ;
265 par_poisson_vect.add_tenseur_mod( ssjm1_wshift ) ;
266 par_poisson_vect.add_int_mod(niter, 0) ;
267
268
269 // External potential
270 // ------------------
271
272 Tenseur pot_ext = logn_comp + pot_centri + loggam ;
273//##
274// des_coupe_z(pot_ext(), 0., 1, "pot_ext", &(ent()) ) ;
275//##
276
277 Tenseur ent_jm1 = ent ; // Enthalpy at previous step
278
279 Tenseur source(mp) ; // source term in the equation for logn_auto
280 // and beta_auto
281
282 Tenseur source_shift(mp, 1, CON, ref_triad) ; // source term in the
283 // equation for shift_auto
284
285 Cmp source_regu(mp) ;
286 Cmp source_div(mp) ;
287
288 //=========================================================================
289 // Start of iteration
290 //=========================================================================
291
292 for(int mer=0 ; mer<mermax ; mer++ ) {
293
294 cout << "-----------------------------------------------" << endl ;
295 cout << "step: " << mer << endl ;
296 cout << "diff_ent = " << diff_ent << endl ;
297
298 //-----------------------------------------------------
299 // Resolution of the elliptic equation for the velocity
300 // scalar potential
301 //-----------------------------------------------------
302
303 if (irrotational) {
304
305 diff_vel_pot = velocity_potential(mermax_potvit, precis_poisson,
306 relax_potvit) ;
307
308 }
309
310 //-----------------------------------------------------
311 // Computation of the new radial scale
312 //-----------------------------------------------------
313
314 // alpha_r (r = alpha_r r') is determined so that the enthalpy
315 // takes the requested value ent_b at the stellar surface
316
317 // Values at the center of the star:
318 double logn_auto_c = logn_auto()(0, 0, 0, 0) ;
319 double pot_ext_c = pot_ext()(0, 0, 0, 0) ;
320
321 // Search for the reference point (theta_*, phi_*) [notation of
322 // Bonazzola, Gourgoulhon & Marck PRD 58, 104020 (1998)]
323 // at the surface of the star
324 // ------------------------------------------------------------
325 double alpha_r2 = 0 ;
326 for (int k=0; k<mg->get_np(l_b); k++) {
327 for (int j=0; j<mg->get_nt(l_b); j++) {
328
329 double pot_ext_b = pot_ext()(l_b, k, j, i_b) ;
330 double logn_auto_b = logn_auto()(l_b, k, j, i_b) ;
331
332 double alpha_r2_jk = ( ent_c - ent_b + pot_ext_c - pot_ext_b) /
333 ( logn_auto_b - logn_auto_c ) ;
334
335// cout << "k, j, alpha_r2_jk : " << k << " " << j << " "
336// << alpha_r2_jk << endl ;
337
338 if (alpha_r2_jk > alpha_r2) {
339 alpha_r2 = alpha_r2_jk ;
340 k_b = k ;
341 j_b = j ;
342 }
343
344 }
345 }
346
347 alpha_r = sqrt(alpha_r2) ;
348
349 cout << "k_b, j_b, alpha_r: " << k_b << " " << j_b << " "
350 << alpha_r << endl ;
351
352 // New value of logn_auto
353 // ----------------------
354
355 logn_auto = alpha_r2 * logn_auto ;
356 logn_auto_regu = alpha_r2 * logn_auto_regu ;
357 logn_auto_c = logn_auto()(0, 0, 0, 0) ;
358
359
360 //------------------------------------------------------------
361 // Change the values of the inner points of the second domain
362 // by those of the outer points of the first domain
363 //------------------------------------------------------------
364
365 (logn_auto().va).smooth(nzet, (logn_auto.set()).va) ;
366
367
368 //--------------------
369 // First integral --> enthalpy in all space
370 //--------------------
371
372 ent = (ent_c + logn_auto_c + pot_ext_c) - logn_auto - pot_ext ;
373
374 (ent().va).smooth(nzet, (ent.set()).va) ;
375
376 //----------------------------------------------------
377 // Adaptation of the mapping to the new enthalpy field
378 //----------------------------------------------------
379
380 // Shall the adaptation be performed (cusp) ?
381 // ------------------------------------------
382
383 double dent_eq = ent().dsdr().val_point(ray_eq(),M_PI/2.,0.) ;
384 double dent_pole = ent().dsdr().val_point(ray_pole(),0.,0.) ;
385 double rap_dent = fabs( dent_eq / dent_pole ) ;
386 cout << "| dH/dr_eq / dH/dr_pole | = " << rap_dent << endl ;
387
388 if ( rap_dent < thres_adapt ) {
389 adapt_flag = 0 ; // No adaptation of the mapping
390 cout << "******* FROZEN MAPPING *********" << endl ;
391 }
392 else{
393 adapt_flag = 1 ; // The adaptation of the mapping is to be
394 // performed
395 }
396
397
398 ent_limit.set_etat_qcq() ;
399 for (int l=0; l<nzet; l++) { // loop on domains inside the star
400 ent_limit.set(l) = ent()(l, k_b, j_b, i_b) ;
401 }
402 ent_limit.set(nzet-1) = ent_b ;
403
404 Map_et mp_prev = mp_et ;
405
406//##
407// des_coupe_z(ent(), 0., 1, "ent before adapt", &(ent()) ) ;
408//##
409
410 mp.adapt(ent(), par_adapt) ;
411
412 // Readjustment of the external boundary of domain l=nzet
413 // to keep a fixed ratio with respect to star's surface
414
415 double rr_in_1 = mp.val_r(nzet, -1., M_PI/2., 0.) ;
416
417 // Resizes the outer boundary of the shell including the comp. NS
418 double rr_out_nm2 = mp.val_r(nz-2, 1., M_PI/2., 0.) ;
419
420 mp.resize(nz-2, rr_in_1/rr_out_nm2 * fact_resize(1)) ;
421
422 // Resizes the inner boundary of the shell including the comp. NS
423 double rr_out_nm3 = mp.val_r(nz-3, 1., M_PI/2., 0.) ;
424
425 mp.resize(nz-3, rr_in_1/rr_out_nm3 * fact_resize(0)) ;
426
427 if (nz > nzet+3) {
428
429 // Resize of the domain from 1(nzet) to N-4
430 double rr_out_nm3_new = mp.val_r(nz-3, 1., M_PI/2., 0.) ;
431
432 for (int i=nzet-1; i<nz-4; i++) {
433
434 double rr_out_i = mp.val_r(i, 1., M_PI/2., 0.) ;
435
436 double rr_mid = rr_out_i
437 + (rr_out_nm3_new - rr_out_i) / double(nz - 3 - i) ;
438
439 double rr_2timesi = 2. * rr_out_i ;
440
441 if (rr_2timesi < rr_mid) {
442
443 double rr_out_ip1 = mp.val_r(i+1, 1., M_PI/2., 0.) ;
444
445 mp.resize(i+1, rr_2timesi / rr_out_ip1) ;
446
447 }
448 else {
449
450 double rr_out_ip1 = mp.val_r(i+1, 1., M_PI/2., 0.) ;
451
452 mp.resize(i+1, rr_mid / rr_out_ip1) ;
453
454 } // End of else
455
456 } // End of i loop
457
458 } // End of (nz > nzet+3) loop
459
460//##
461// des_coupe_z(ent(), 0., 1, "ent after adapt", &(ent()) ) ;
462//##
463 //----------------------------------------------------
464 // Computation of the enthalpy at the new grid points
465 //----------------------------------------------------
466
467 mp_prev.homothetie(alpha_r) ;
468
469 mp.reevaluate_symy(&mp_prev, nzet+1, ent.set()) ;
470
471// des_coupe_z(ent(), 0., 1, "ent after reevaluate", &(ent()) ) ;
472
473 double ent_s_max = -1 ;
474 int k_s_max = -1 ;
475 int j_s_max = -1 ;
476 for (int k=0; k<mg->get_np(l_b); k++) {
477 for (int j=0; j<mg->get_nt(l_b); j++) {
478 double xx = fabs( ent()(l_b, k, j, i_b) ) ;
479 if (xx > ent_s_max) {
480 ent_s_max = xx ;
481 k_s_max = k ;
482 j_s_max = j ;
483 }
484 }
485 }
486 cout << "Max. abs(enthalpy) at the boundary between domains nzet-1"
487 << " and nzet : " << endl ;
488 cout << " " << ent_s_max << " reached for k = " << k_s_max <<
489 " and j = " << j_s_max << endl ;
490
491 //----------------------------------------------------
492 // Equation of state
493 //----------------------------------------------------
494
495 equation_of_state() ; // computes new values for nbar (n), ener (e)
496 // and press (p) from the new ent (H)
497
498 //---------------------------------------------------------
499 // Matter source terms in the gravitational field equations
500 //---------------------------------------------------------
501
502 hydro_euler() ; // computes new values for ener_euler (E),
503 // s_euler (S) and u_euler (U^i)
504
505 //--------------------------------------------------------
506 // Poisson equation for logn_auto (nu_auto)
507 //--------------------------------------------------------
508
509 // Source
510 // ------
511
512 double unsgam1 ; // effective power of H in the source
513 // close to the surface
514
515 if (relativistic) {
516 source = qpig * a_car % (ener_euler + s_euler)
520
521 // 1/(gam-1) = dln(e)/dln(H) close to the surface :
522 unsgam1 = eos.der_ener_ent_p(ent_b + 1e-10*(ent_c-ent_b)) ;
523
524 }
525 else {
526 source = qpig * nbar ;
527
528 // 1/(gam-1) = dln(n)/dln(H) close to the surface :
529 unsgam1 = eos.der_nbar_ent_p(ent_b + 1e-10*(ent_c-ent_b)) ;
530 }
531
532 source.set_std_base() ;
533
534 // Resolution of the Poisson equation
535 // ----------------------------------
536
540
541 source_regu.std_base_scal() ;
542 source_div.std_base_scal() ;
543
544 source().poisson_regular(k_div, nzet, unsgam1, par_poisson1,
548 source_regu, source_div) ;
549
550 // Check: has the Poisson equation been correctly solved ?
551 // -----------------------------------------------------
552
553 Tbl tdiff_logn = diffrel(logn_auto().laplacien(), source()) ;
554 cout <<
555 "Relative error in the resolution of the equation for logn_auto : "
556 << endl ;
557 for (int l=0; l<nz; l++) {
558 cout << tdiff_logn(l) << " " ;
559 }
560 cout << endl ;
561 diff_logn = max(abs(tdiff_logn)) ;
562
563
564 if (relativistic) {
565
566 //--------------------------------------------------------
567 // Poisson equation for beta_auto
568 //--------------------------------------------------------
569
570 // Source
571 // ------
572
573 source = qpig * a_car % s_euler
574 + .75 * ( akcar_auto + akcar_comp )
579
580 source.set_std_base() ;
581
582 // Resolution of the Poisson equation
583 // ----------------------------------
584
585 source().poisson(par_poisson2, beta_auto.set()) ;
586
587
588 // Check: has the Poisson equation been correctly solved ?
589 // -----------------------------------------------------
590
591 Tbl tdiff_beta = diffrel(beta_auto().laplacien(), source()) ;
592 cout <<
593 "Relative error in the resolution of the equation for beta_auto : "
594 << endl ;
595 for (int l=0; l<nz; l++) {
596 cout << tdiff_beta(l) << " " ;
597 }
598 cout << endl ;
599 diff_beta = max(abs(tdiff_beta)) ;
600
601 //--------------------------------------------------------
602 // Vector Poisson equation for shift_auto
603 //--------------------------------------------------------
604
605 // Source
606 // ------
607
608 Tenseur vtmp = 6. * ( d_beta_auto + d_beta_comp )
609 -8. * ( d_logn_auto + d_logn_comp ) ;
610
611 source_shift = (-4.*qpig) * nnn % a_car % (ener_euler + press)
612 % u_euler
614
615 source_shift.set_std_base() ;
616
617 // Resolution of the Poisson equation
618 // ----------------------------------
619
620 // Filter for the source of shift vector
621
622 for (int i=0; i<3; i++) {
623
624 if (source_shift(i).get_etat() != ETATZERO)
625 source_shift.set(i).filtre(4) ;
626
627 }
628
629 // For Tenseur::poisson_vect_regu,
630 // the triad must be the mapping triad,
631 // not the reference one:
632
633 source_shift.change_triad( mp.get_bvect_cart() ) ;
634
635 for (int i=0; i<3; i++) {
636 if(source_shift(i).dz_nonzero()) {
637 assert( source_shift(i).get_dzpuis() == 4 ) ;
638 }
639 else{
640 (source_shift.set(i)).set_dzpuis(4) ;
641 }
642 }
643
644 //##
645 // source_shift.dec2_dzpuis() ; // dzpuis 4 -> 2
646
647 double lambda_shift = double(1) / double(3) ;
648
649 source_shift.poisson_vect_regu(k_div, nzet, unsgam1,
650 lambda_shift, par_poisson_vect,
652
653
654 // Check: has the equation for shift_auto been correctly solved ?
655 // --------------------------------------------------------------
656
657 // Divergence of shift_auto :
658 Tenseur divn = contract(shift_auto.gradient(), 0, 1) ;
659 divn.dec2_dzpuis() ; // dzpuis 2 -> 0
660
661 // Grad(div) :
662 Tenseur graddivn = divn.gradient() ;
663 graddivn.inc2_dzpuis() ; // dzpuis 2 -> 4
664
665 // Full operator :
666 Tenseur lap_shift(mp, 1, CON, mp.get_bvect_cart() ) ;
667 lap_shift.set_etat_qcq() ;
668 for (int i=0; i<3; i++) {
669 lap_shift.set(i) = shift_auto(i).laplacien()
670 + lambda_shift * graddivn(i) ;
671 }
672
673 Tbl tdiff_shift_x = diffrel(lap_shift(0), source_shift(0)) ;
674 Tbl tdiff_shift_y = diffrel(lap_shift(1), source_shift(1)) ;
675 Tbl tdiff_shift_z = diffrel(lap_shift(2), source_shift(2)) ;
676
677 cout <<
678 "Relative error in the resolution of the equation "
679 "for shift_auto : "
680 << endl ;
681 cout << "x component : " ;
682 for (int l=0; l<nz; l++) {
683 cout << tdiff_shift_x(l) << " " ;
684 }
685 cout << endl ;
686 cout << "y component : " ;
687 for (int l=0; l<nz; l++) {
688 cout << tdiff_shift_y(l) << " " ;
689 }
690 cout << endl ;
691 cout << "z component : " ;
692 for (int l=0; l<nz; l++) {
693 cout << tdiff_shift_z(l) << " " ;
694 }
695 cout << endl ;
696
697 diff_shift_x = max(abs(tdiff_shift_x)) ;
698 diff_shift_y = max(abs(tdiff_shift_y)) ;
699 diff_shift_z = max(abs(tdiff_shift_z)) ;
700
701 // Final result
702 // ------------
703 // The output of Tenseur::poisson_vect is on the mapping triad,
704 // it should therefore be transformed to components on the
705 // reference triad :
706
708
709
710 } // End of relativistic equations
711
712
713 //-------------------------------------------------
714 // Relative change in enthalpy
715 //-------------------------------------------------
716
717 Tbl diff_ent_tbl = diffrel( ent(), ent_jm1() ) ;
718 diff_ent = diff_ent_tbl(0) ;
719 for (int l=1; l<nzet; l++) {
720 diff_ent += diff_ent_tbl(l) ;
721 }
722 diff_ent /= nzet ;
723
724 ent_jm1 = ent ;
725
726
727 } // End of main loop
728
729 //=========================================================================
730 // End of iteration
731 //=========================================================================
732
733
734}
735
736}
Component of a tensorial field *** DEPRECATED : use class Scalar instead ***.
Definition cmp.h:446
void std_base_scal()
Sets the spectral bases of the Valeur va to the standard ones for a scalar.
Definition cmp.C:644
void filtre(int n)
Sets the n lasts coefficients in r to 0 in the external domain.
Definition cmp_manip.C:74
virtual double der_ener_ent_p(double ent, const Param *par=0x0) const =0
Computes the logarithmic derivative from the log-enthalpy with extra parameters (virtual function im...
virtual double der_nbar_ent_p(double ent, const Param *par=0x0) const =0
Computes the logarithmic derivative from the log-enthalpy and extra parameters (virtual function imp...
Tenseur d_beta_auto
Gradient of beta_auto (Cartesian components with respect to ref_triad )
Definition etoile.h:879
Tenseur d_logn_auto
Gradient of logn_auto (Cartesian components with respect to ref_triad )
Definition etoile.h:859
Tenseur akcar_comp
Part of the scalar generated by shift_auto and shift_comp , i.e.
Definition etoile.h:944
Cmp ssjm1_khi
Effective source at the previous step for the resolution of the Poisson equation for the scalar by m...
Definition etoile.h:973
double velocity_potential(int mermax, double precis, double relax)
Computes the non-translational part of the velocity scalar potential by solving the continuity equat...
const Base_vect & ref_triad
Reference triad ("absolute frame"), with respect to which the components of all the member Tenseur 's...
Definition etoile.h:828
bool irrotational
true for an irrotational star, false for a corotating one
Definition etoile.h:822
Tenseur w_shift
Vector used in the decomposition of shift_auto , following Shibata's prescription [Prog.
Definition etoile.h:908
Tenseur logn_comp
Part of the lapse logarithm (gravitational potential at the Newtonian limit) generated principaly by ...
Definition etoile.h:854
Tenseur d_logn_comp
Gradient of logn_comp (Cartesian components with respect to ref_triad )
Definition etoile.h:869
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 by means...
Definition etoile.h:983
Cmp ssjm1_logn
Effective source at the previous step for the resolution of the Poisson equation for logn_auto by mea...
Definition etoile.h:959
Tenseur pot_centri
Centrifugal potential.
Definition etoile.h:953
void equil_regular(double ent_c, int mermax, int mermax_poisson, double relax_poisson, int mermax_potvit, double relax_potvit, double thres_adapt, const Tbl &fact, Tbl &diff)
Computes an equilibrium configuration by regularizing the diverging source.
Tenseur loggam
Logarithm of the Lorentz factor between the fluid and the co-orbiting observer.
Definition etoile.h:849
Tenseur shift_auto
Part of the shift vector generated principaly by the star.
Definition etoile.h:889
Tenseur_sym tkij_auto
Part of the extrinsic curvature tensor generated by shift_auto .
Definition etoile.h:925
Tenseur akcar_auto
Part of the scalar generated by shift_auto , i.e.
Definition etoile.h:938
Cmp ssjm1_beta
Effective source at the previous step for the resolution of the Poisson equation for beta_auto by mea...
Definition etoile.h:965
Tenseur d_beta_comp
Gradient of beta_comp (Cartesian components with respect to ref_triad )
Definition etoile.h:884
Tenseur khi_shift
Scalar used in the decomposition of shift_auto , following Shibata's prescription [Prog.
Definition etoile.h:918
Tenseur logn_auto_div
Divergent part (if k_div!=0 ) of the logarithm of the part of the lapse N generated principaly by t...
Definition etoile.h:497
int nzet
Number of domains of *mp occupied by the star.
Definition etoile.h:432
Tenseur logn_auto_regu
Regular part of the logarithm of the part of the lapse N generated principaly by the star.
Definition etoile.h:491
int k_div
Index of regularity of the gravitational potential logn_auto .
Definition etoile.h:449
double ray_eq() const
Coordinate radius at , [r_unit].
Tenseur nnn
Total lapse function.
Definition etoile.h:509
Tenseur logn_auto
Total of the logarithm of the part of the lapse N generated principaly by the star.
Definition etoile.h:484
Tenseur nbar
Baryon density in the fluid frame.
Definition etoile.h:459
const Eos & eos
Equation of state of the stellar matter.
Definition etoile.h:451
virtual void equation_of_state()
Computes the proper baryon and energy density, as well as pressure from the enthalpy.
Definition etoile.C:566
Tenseur u_euler
Fluid 3-velocity with respect to the Eulerian observer.
Definition etoile.h:474
Map & mp
Mapping associated with the star.
Definition etoile.h:429
Tenseur d_logn_auto_div
Gradient of logn_auto_div (if k_div!=0 )
Definition etoile.h:501
Tenseur press
Fluid pressure.
Definition etoile.h:461
bool relativistic
Indicator of relativity: true for a relativistic star, false for a Newtonian one.
Definition etoile.h:437
Tenseur ener_euler
Total energy density in the Eulerian frame.
Definition etoile.h:465
Tenseur s_euler
Trace of the stress tensor in the Eulerian frame.
Definition etoile.h:468
Tenseur ent
Log-enthalpy (relativistic case) or specific enthalpy (Newtonian case)
Definition etoile.h:457
Tenseur beta_auto
Logarithm of the part of the product AN generated principaly by by the star.
Definition etoile.h:506
Tenseur a_car
Total conformal factor .
Definition etoile.h:515
double ray_pole() const
Coordinate radius at [r_unit].
Radial mapping of rather general form.
Definition map.h:2752
virtual void homothetie(double lambda)
Sets a new radial scale.
Definition map_et.C:905
const Base_vect_cart & get_bvect_cart() const
Returns the Cartesian basis associated with the coordinates (x,y,z) of the mapping,...
Definition map.h:791
virtual void resize(int l, double lambda)=0
Rescales the outer boundary of one domain.
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.
virtual void adapt(const Cmp &ent, const Param &par, int nbr=0)=0
Adaptation of the mapping to a given scalar field.
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.
const Mg3d * get_mg() const
Gives the Mg3d on which the mapping is defined.
Definition map.h:765
Multi-domain grid.
Definition grilles.h:273
int get_np(int l) const
Returns the number of points in the azimuthal direction ( ) in domain no. l.
Definition grilles.h:462
int get_nt(int l) const
Returns the number of points in the co-latitude direction ( ) in domain no. l.
Definition grilles.h:457
int get_nzone() const
Returns the number of domains.
Definition grilles.h:448
int get_nr(int l) const
Returns the number of points in the radial direction ( ) in domain no. l.
Definition grilles.h:452
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:315
void add_cmp_mod(Cmp &ti, int position=0)
Adds the address of a new modifiable Cmp to the list.
Definition param.C:1004
void add_int_mod(int &n, int position=0)
Adds the address of a new modifiable int to the list.
Definition param.C:385
void add_tenseur_mod(Tenseur &ti, int position=0)
Adds the address of a new modifiable Tenseur to the list.
Definition param.C:1142
void add_int(const int &n, int position=0)
Adds the address of a new int to the list.
Definition param.C:246
void add_tbl(const Tbl &ti, int position=0)
Adds the address of a new Tbl to the list.
Definition param.C:522
Basic array class.
Definition tbl.h:161
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition tbl.C:361
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:301
Cmp & set()
Read/write for a scalar (see also operator=(const Cmp&) ).
Definition tenseur.C:824
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition tenseur.C:636
void dec2_dzpuis()
dzpuis -= 2 ;
Definition tenseur.C:1130
const Tenseur & gradient() const
Returns the gradient of *this (Cartesian coordinates)
Definition tenseur.C:1542
void set_std_base()
Set the standard spectal basis of decomposition for each component.
Definition tenseur.C:1170
void inc2_dzpuis()
dzpuis += 2 ;
Definition tenseur.C:1143
void change_triad(const Base_vect &new_triad)
Sets a new vectorial basis (triad) of decomposition and modifies the components accordingly.
Definition tenseur.C:668
void poisson_vect_regu(int k_div, int nzet, double unsgam1, double lambda, Param &par, Tenseur &shift, Tenseur &vect, Tenseur &scal) const
Solves the vectorial Poisson equation : .
Cmp sqrt(const Cmp &)
Square root.
Definition cmp_math.C:220
Tbl diffrel(const Cmp &a, const Cmp &b)
Relative difference between two Cmp (norme version).
Definition cmp_math.C:504
Tbl max(const Cmp &)
Maximum values of a Cmp in each domain.
Definition cmp_math.C:435
Cmp abs(const Cmp &)
Absolute value.
Definition cmp_math.C:410
Tenseur contract(const Tenseur &, int id1, int id2)
Self contraction of two indices of a Tenseur .
Tenseur flat_scalar_prod_desal(const Tenseur &t1, const Tenseur &t2)
Same as flat_scalar_prod but with desaliasing.
Lorene prototypes.
Definition app_hor.h:64
Standard units of space, time and mass.