LORENE
et_bin_vel_pot.C
1/*
2 * Method of class Etoile_bin to compute the velocity scalar potential $\psi$
3 * by solving the continuity equation.
4 *
5 * (see file etoile.h for documentation).
6 *
7 */
8
9/*
10 * Copyright (c) 2000-2001 Eric Gourgoulhon
11 *
12 * This file is part of LORENE.
13 *
14 * LORENE is free software; you can redistribute it and/or modify
15 * it under the terms of the GNU General Public License as published by
16 * the Free Software Foundation; either version 2 of the License, or
17 * (at your option) any later version.
18 *
19 * LORENE is distributed in the hope that it will be useful,
20 * but WITHOUT ANY WARRANTY; without even the implied warranty of
21 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
22 * GNU General Public License for more details.
23 *
24 * You should have received a copy of the GNU General Public License
25 * along with LORENE; if not, write to the Free Software
26 * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
27 *
28 */
29
30
31
32
33/*
34 * $Id: et_bin_vel_pot.C,v 1.16 2016/12/05 16:17:53 j_novak Exp $
35 * $Log: et_bin_vel_pot.C,v $
36 * Revision 1.16 2016/12/05 16:17:53 j_novak
37 * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
38 *
39 * Revision 1.15 2014/10/13 08:52:56 j_novak
40 * Lorene classes and functions now belong to the namespace Lorene.
41 *
42 * Revision 1.14 2007/10/18 14:26:43 e_gourgoulhon
43 * Changed the call to Eos::der_nbar_ent in order to allow for MEos type
44 * of equation of state.
45 *
46 * Revision 1.13 2007/10/16 21:56:26 e_gourgoulhon
47 * Can deal with more than one domain into the star,
48 * thanks to the new function Map_radial::poisson_compact.
49 *
50 * Revision 1.12 2005/10/18 13:12:33 p_grandclement
51 * update of the mixted binary codes
52 *
53 * Revision 1.11 2004/05/25 15:38:38 f_limousin
54 * Minor modifs.
55 *
56 * Revision 1.10 2004/05/10 10:17:27 f_limousin
57 * Add a new member ssjm1_psi of class Etoile for the resolution of the
58 * oisson_interne equation
59 *
60 * Revision 1.9 2004/04/19 11:26:17 f_limousin
61 * Add a new function Etoile_bin::velocity_potential( , , , ) for the
62 * case of strange stars
63 *
64 * Revision 1.8 2004/04/08 17:02:00 f_limousin
65 * Modif to avoid an error in the compilation
66 *
67 * Revision 1.7 2004/04/08 16:52:58 f_limousin
68 * Minor change
69 *
70 * Revision 1.6 2004/04/08 16:36:36 f_limousin
71 * Implement the resolution of the continuity equation for strange
72 * stars.
73 *
74 * Revision 1.5 2003/10/24 11:43:57 e_gourgoulhon
75 * beta is now computed as ln(AN) in the case beta_auto
76 * is undefined (for instance, if the companion is black hole).
77 *
78 * Revision 1.4 2003/01/17 13:38:56 f_limousin
79 * Add comments
80 *
81 * Revision 1.3 2003/01/13 15:31:50 e_gourgoulhon
82 * Suppressed the desaliasing
83 * (did not worked due to missing basis in ylm).
84 *
85 * Revision 1.2 2002/12/10 14:44:21 k_taniguchi
86 * Change the multiplication "*" to "%"
87 * and flat_scalar_prod to flat_scalar_prod_desal.
88 *
89 * Revision 1.1.1.1 2001/11/20 15:19:28 e_gourgoulhon
90 * LORENE
91 *
92 * Revision 2.9 2001/02/23 15:18:59 eric
93 * Modification du calcul de zeta_h pour eviter division par zero
94 * dans les domaines externes a l'etoile.
95 *
96 * Revision 2.8 2001/02/07 09:47:42 eric
97 * zeta_h est desormais donne par Eos::der_nbar_ent.
98 *
99 * Revision 2.7 2000/12/22 13:10:03 eric
100 * Prolongement C^1 de dpsi en dehors de l'etoile.
101 *
102 * Revision 2.6 2000/03/22 12:56:44 eric
103 * Nouveau prototype d'Etoile_bin::velocity_potential : l'erreur est
104 * retournee en double.
105 *
106 * Revision 2.5 2000/02/25 17:35:29 eric
107 * Annulation de la source dans les zones externes avant l'appel a
108 * poisson_compact.
109 *
110 * Revision 2.4 2000/02/22 11:42:55 eric
111 * Test resolution de l'equation.
112 *
113 * Revision 2.3 2000/02/22 10:42:25 eric
114 * Correction erreur dans les termes sources: multiplication par unsurc2 de
115 * termes relativistes.
116 *
117 * Revision 2.2 2000/02/21 15:05:50 eric
118 * Traitement du cas psi0 = 0 .
119 *
120 * Revision 2.1 2000/02/21 13:59:39 eric
121 * Remplacement du membre psi par psi0.
122 * Modif calcul de d_psi a la fin.
123 *
124 * Revision 2.0 2000/02/17 18:50:44 eric
125 * *** empty log message ***
126 *
127 *
128 * $Header: /cvsroot/Lorene/C++/Source/Etoile/et_bin_vel_pot.C,v 1.16 2016/12/05 16:17:53 j_novak Exp $
129 *
130 */
131
132// Headers Lorene
133#include "scalar.h"
134#include "metrique.h"
135#include "etoile.h"
136#include "eos.h"
137#include "param.h"
138#include "et_bin_nsbh.h"
139#include "utilitaires.h"
140
141// Local prototype
142namespace Lorene {
143Cmp raccord_c1(const Cmp& uu, int l1) ;
144
145double Etoile_bin::velocity_potential(int mermax, double precis, double relax) {
146
147 // Which star is that ?
148 const Et_bin_nsbh* pnsbh = dynamic_cast<const Et_bin_nsbh*>(this) ;
149
150 if (eos.identify() == 5 || eos.identify() == 4 ||
151 eos.identify() == 3) {
152
153 // Routine used for binary strange stars.
154
155 int nzm1 = mp.get_mg()->get_nzone() - 1 ;
156
157 //----------------------------------
158 // Specific relativistic enthalpy ---> hhh
159 //----------------------------------
160
161 Tenseur hhh = exp(unsurc2 * ent) ; // = 1 at the Newtonian limit
162 hhh.set_std_base() ;
163
164 //----------------------------------------------
165 // Computation of W^i = - A^2 h Gamma_n B^i/N
166 // See Eq (62) from Gourgoulhon et al. (2001)
167 //----------------------------------------------
168
169 Tenseur www = - a_car * hhh * gam_euler * bsn ;
170
171 www.change_triad( mp.get_bvect_cart() ) ; // components on the mapping
172 // Cartesian basis
173
174 //-------------------------------------------------
175 // Constant value of W^i at the center of the star
176 //-------------------------------------------------
177
178 Tenseur v_orb(mp, 1, CON, mp.get_bvect_cart()) ;
179
180 v_orb.set_etat_qcq() ;
181 for (int i=0; i<3; i++) {
182 v_orb.set(i) = www(i)(0, 0, 0, 0) ;
183 }
184
185 v_orb.annule(nzm1, nzm1) ; // set to zero in the ZEC
186
187
188 v_orb.set_triad( *(www.get_triad()) ) ;
189 v_orb.set_std_base() ;
190
191
192 //-------------------------------------------------
193 // Source and coefficients a,b for poisson_compact (idenpendent from psi0)
194 //-------------------------------------------------
195
196 Cmp dndh_log = eos.der_nbar_ent(ent(), nzet) ;
197
198 // In order to avoid any division by zero in the computation of zeta_h
199 // the value of dndh_log is set to 1 in the external domains:
200 for (int l=nzet; l <= nzm1; l++) {
201 dndh_log.set(l) = 1 ;
202 }
203
204 double erreur ;
205
206 Tenseur zeta_h( ent() / dndh_log ) ;
207 zeta_h.set_std_base() ;
208
209 Scalar zeta_h_scalar (zeta_h()) ;
210 zeta_h_scalar.set_outer_boundary(0, (ent() / dndh_log)(0,0,0,0)) ;
211 for (int l=1; l<=nzm1; l++)
212 zeta_h_scalar.set_domain(l) = 1 ;
213
214 Cmp zeta_h_cmp (zeta_h_scalar) ;
215 zeta_h.set() = zeta_h_cmp ;
216 zeta_h.set_std_base() ;
217
218
219
220 Tenseur beta(mp) ;
221
222 if (pnsbh!=0x0) {
223 beta = log( sqrt(a_car) * nnn ) ;
224 beta.set_std_base() ;
225 }
226 else {
227 beta = beta_auto + beta_comp ;
228 }
229
230 Tenseur tmp_zeta = 1 - unsurc2 * zeta_h ;
231 tmp_zeta.set_std_base() ;
232
233 Tenseur bb = tmp_zeta * ent.gradient_spher()
234 + unsurc2 * zeta_h * beta.gradient_spher() ;
235
236 Tenseur entmb = ent - beta ;
237
238 Tenseur grad_ent (ent.gradient()) ;
239 grad_ent.change_triad(mp.get_bvect_spher()) ;
240
241 // Source for the poisson equation
242 // See Eq (63) from Gourgoulhon et al. (2001)
243 Tenseur source = flat_scalar_prod( www - v_orb, ent.gradient() )
244 + unsurc2 * zeta_h * (
245 flat_scalar_prod( v_orb, entmb.gradient() )
246 + flat_scalar_prod( www, gam_euler.gradient() )
247 / gam_euler ) ;
248
249 for (int l=1; l<=nzm1; l++)
250 source.set().annule(l) ;
251
252 source = (source - flat_scalar_prod(bb, psi0.gradient_spher()))
253 / zeta_h ;
254 source.annule(nzet, nzm1) ;
255
256 Param par ;
257 int niter ;
258 par.add_int(mermax) ;
259 par.add_double(precis, 0) ;
260 par.add_double(relax, 1) ;
261 par.add_int_mod(niter) ;
262
263 par.add_cmp_mod(ssjm1_psi, 0) ;
264
265 if (psi0.get_etat() == ETATZERO) {
266 psi0.set_etat_qcq() ;
267 psi0.set() = 0 ;
268 }
269
270 int nr = mp.get_mg()->get_nr(0);
271 int nt = mp.get_mg()->get_nt(0);
272 int np = mp.get_mg()->get_np(0);
273
274 cout << "nr = " << nr << " nt = " << nt << " np = " << np << endl ;
275
276 cout << "psi0" << endl << norme(psi0()/(nr*nt*np)) << endl ;
277 cout << "d(psi)/dr" << endl << norme(psi0.set().dsdr()/(nr*nt*np)) << endl ;
278
279 Valeur lim(mp.get_mg()->get_angu()) ;
280 lim.annule_hard() ;
281
282 Tenseur normal (mp, 1, CON, mp.get_bvect_cart()) ;
283 Tenseur normal2 (mp, 1, COV, mp.get_bvect_cart()) ;
284 normal.set_etat_qcq() ;
285 normal2.set_etat_qcq() ;
286
287 const Coord& rr0 = mp.r ;
288 Tenseur rr(mp) ;
289 rr.set_etat_qcq() ;
290 rr.set() = rr0 ;
291 rr.set_std_base() ;
292
293 Tenseur_sym plat(mp, 2, COV, mp.get_bvect_cart() ) ;
294 plat.set_etat_qcq() ;
295 for (int i=0; i<3; i++) {
296 for (int j=0; j<i; j++) {
297 plat.set(i,j) = 0 ;
298 }
299 plat.set(i,i) = 1 ;
300 }
301 plat.set_std_base() ;
302
303 Metrique flat(plat, true) ;
304 Tenseur dcov_r = rr.derive_cov(flat) ;
305
306
307 for (int i=0; i<3; i++) {
308 normal.set(i) = dcov_r(i) ;
309 normal2.set(i) = dcov_r(i) ;
310 }
311
312 normal.change_triad(mp.get_bvect_spher()) ;
313 normal2.change_triad(mp.get_bvect_spher()) ;
314
315
316
317 Tenseur bsn0 (bsn) ;
318 bsn0.change_triad(mp.get_bvect_cart()) ;
319 Tenseur aa (mp, 1, CON, mp.get_bvect_cart()) ;
320 aa = - v_orb - a_car * gam_euler * hhh * bsn0 ;
321 aa.change_triad(mp.get_bvect_spher()) ;
322
323
324 Tenseur dcov_psi = psi0.derive_cov(flat) ;
325 dcov_psi.change_triad(mp.get_bvect_spher()) ;
326
327 Cmp limite (mp) ;
328 limite = ( - dcov_psi(1) * normal(1) - dcov_psi(2) * normal(2)
329 + contract(aa, 0, normal2, 0)())
330 /normal(0) ;
331
332 for (int j=0; j<nt; j++)
333 for (int k=0; k<np; k++)
334 lim.set(0, k, j, 0) = limite(0, k, j, nr-1) ;
335
336// cout << "lim" << endl << lim << endl ;
337
338 lim.std_base_scal() ;
339 Cmp resu (psi0()) ;
340 source().poisson_neumann_interne(lim, par, resu) ;
341 psi0 = resu ;
342
343/*
344 resu.va.ylm() ;
345 Scalar psi00(resu) ;
346 psi00.spectral_display("psi00") ;
347
348 cout << "value of d(psi)/dr at the surface after poisson" << endl ;
349 for (int j=0; j<nt; j++)
350 for (int k=0; k<np; k++)
351 cout << "j = " << j << " ; k = " << k << " : " <<
352 psi0.set().dsdr()(0, k, j, nr-1) << endl ;
353*/
354 for (int l=1; l<=nzm1; l++)
355 psi0.set().annule(l) ;
356
357
358 //---------------------------------------------------
359 // Check of the solution
360 //---------------------------------------------------
361
362 Cmp laplacien_psi0 = psi0().laplacien() ;
363
364 erreur = diffrel(laplacien_psi0, source())(0) ;
365
366 cout << "Check of the resolution of the continuity equation for strange stars: "
367 << endl ;
368 cout << "norme(source) : " << norme(source())(0) << endl
369 << "Error in the solution : " << erreur << endl ;
370
371 //--------------------------------
372 // Computation of grad(psi)
373 //--------------------------------
374
375 // The computation is done component by component because psi0.gradient()
376 // is a covariant vector, whereas v_orb is a contravariant one.
377
378 d_psi.set_etat_qcq() ;
379
380 for (int i=0; i<3; i++) {
381 d_psi.set(i) = (psi0.gradient())(i) + v_orb(i) ;
382 }
383
384 d_psi.set_triad( *(v_orb.get_triad()) ) ;
385
386 // C^1 continuation of d_psi outside the star
387 // (to ensure a smooth enthalpy field accross the stellar surface)
388 // ----------------------------------------------------------------
389
390 d_psi.annule(nzet, nzm1) ;
391 for (int i=0; i<3; i++) {
392 d_psi.set(i) = raccord_c1(d_psi(i), nzet) ;
393 }
394
395 assert( d_psi.get_triad() == &(mp.get_bvect_cart()) ) ;
396
397 d_psi.change_triad(ref_triad) ;
398
399 return erreur ;
400
401
402 } // End of strange stars case
403
404//=============================================================================
405
406 else {
407
408 int nzm1 = mp.get_mg()->get_nzone() - 1 ;
409
410 //----------------------------------
411 // Specific relativistic enthalpy ---> hhh
412 //----------------------------------
413
414 Tenseur hhh = exp(unsurc2 * ent) ; // = 1 at the Newtonian limit
415 hhh.set_std_base() ;
416
417 //----------------------------------------------
418 // Computation of W^i = - A^2 h Gamma_n B^i/N
419 // See Eq (62) from Gourgoulhon et al. (2001)
420 //----------------------------------------------
421
422 Tenseur www = - a_car * hhh * gam_euler * bsn ;
423
424 www.change_triad( mp.get_bvect_cart() ) ; // components on the mapping
425 // Cartesian basis
426
427 //-------------------------------------------------
428 // Constant value of W^i at the center of the star
429 //-------------------------------------------------
430
431 Tenseur v_orb(mp, 1, CON, mp.get_bvect_cart()) ;
432
433 v_orb.set_etat_qcq() ;
434 for (int i=0; i<3; i++) {
435 v_orb.set(i) = www(i)(0, 0, 0, 0) ;
436 }
437
438 v_orb.set_triad( *(www.get_triad()) ) ;
439 v_orb.set_std_base() ;
440
441
442 //-------------------------------------------------
443 // Source and coefficients a,b for poisson_compact (idenpendent from psi0)
444 //-------------------------------------------------
445
446 Cmp dndh_log(mp) ;
447 dndh_log = 0 ;
448
449 for (int l=0; l<nzet; l++) {
450
451 Param par ; // Paramater for multi-domain equation of state
452 par.add_int(l) ;
453
454 dndh_log = dndh_log + eos.der_nbar_ent(ent(), 1, l, &par) ;
455
456 }
457
458 // Cmp dndh_log = eos.der_nbar_ent(ent(), nzet) ;
459
460 // In order to avoid any division by zero in the computation of zeta_h
461 // the value of dndh_log is set to 1 in the external domains:
462 for (int l=nzet; l <= nzm1; l++) {
463 dndh_log.set(l) = 1 ;
464 }
465
466 Tenseur zeta_h( ent() / dndh_log ) ;
467 zeta_h.set_std_base() ;
468
469 Tenseur beta(mp) ;
470
471 if (pnsbh!=0x0) {
472 beta = log( sqrt(a_car) * nnn ) ;
473 beta.set_std_base() ;
474 }
475 else {
476 beta = beta_auto + beta_comp ;
477 }
478
479 Tenseur tmp_zeta = 1 - unsurc2 * zeta_h ;
480 tmp_zeta.set_std_base() ;
481
482 Tenseur bb = tmp_zeta * ent.gradient_spher()
483 + unsurc2 * zeta_h * beta.gradient_spher() ;
484
485 Tenseur entmb = ent - beta ;
486
487 // See Eq (63) from Gourgoulhon et al. (2001)
488 Tenseur source = flat_scalar_prod( www - v_orb, ent.gradient() )
489 + unsurc2 * zeta_h * (
490 flat_scalar_prod( v_orb, entmb.gradient() )
491 + flat_scalar_prod( www, gam_euler.gradient() )
492 / gam_euler ) ;
493
494
495 source.annule(nzet, nzm1) ;
496
497 //---------------------------------------------------
498 // Resolution by means of Map_radial::poisson_compact
499 //---------------------------------------------------
500
501 Param par ;
502 int niter ;
503 par.add_int(mermax) ;
504 par.add_double(precis, 0) ;
505 par.add_double(relax, 1) ;
506 par.add_int_mod(niter) ;
507
508
509 if (psi0.get_etat() == ETATZERO) {
510 psi0.set_etat_qcq() ;
511 psi0.set() = 0 ;
512 }
513
514 source.set().va.ylm() ;
515
516 mp.poisson_compact(nzet, source(), zeta_h(), bb, par, psi0.set() ) ;
517
518 //---------------------------------------------------
519 // Check of the solution
520 //---------------------------------------------------
521
522 Tenseur bb_dpsi0 = flat_scalar_prod( bb, psi0.gradient_spher() ) ;
523
524 Cmp oper = zeta_h() * psi0().laplacien() + bb_dpsi0() ;
525
526 source.set().va.ylm_i() ;
527
528 cout << "Check of the resolution of the continuity equation : " << endl ;
529 Tbl terr = diffrel(oper, source()) ;
530 double erreur = 0 ;
531 for (int l=0; l<nzet; l++) {
532 double err = terr(l) ;
533 cout << " domain " << l << " : norme(source) : " << norme(source())(l)
534 << " relative error : " << err << endl ;
535 if (err > erreur) erreur = err ;
536 }
537 // arrete() ;
538
539 //--------------------------------
540 // Computation of grad(psi)
541 //--------------------------------
542
543 // The computation is done component by component because psi0.gradient()
544 // is a covariant vector, whereas v_orb is a contravariant one.
545
546 d_psi.set_etat_qcq() ;
547
548 for (int i=0; i<3; i++) {
549 d_psi.set(i) = (psi0.gradient())(i) + v_orb(i) ;
550 }
551
552 d_psi.set_triad( *(v_orb.get_triad()) ) ;
553
554 // C^1 continuation of d_psi outside the star
555 // (to ensure a smooth enthalpy field accross the stellar surface)
556 // ----------------------------------------------------------------
557
558 d_psi.annule(nzet, nzm1) ;
559 for (int i=0; i<3; i++) {
560 d_psi.set(i) = raccord_c1(d_psi(i), nzet) ;
561 }
562
563 assert( d_psi.get_triad() == &(mp.get_bvect_cart()) ) ;
564
565 d_psi.change_triad(ref_triad) ;
566
567 return erreur ;
568
569
570 }
571}
572}
Component of a tensorial field *** DEPRECATED : use class Scalar instead ***.
Definition cmp.h:446
Valeur va
The numerical value of the Cmp.
Definition cmp.h:464
void annule(int l)
Sets the Cmp to zero in a given domain.
Definition cmp.C:351
Tbl & set(int l)
Read/write of the value in a given domain.
Definition cmp.h:724
Active physical coordinates and mapping derivatives.
Definition coord.h:90
Class for a star in a NS-BH binary system.
Definition et_bin_nsbh.h:79
Tenseur psi0
Scalar potential of the non-translational part of fluid 4-velocity (in the irrotational case).
Definition etoile.h:836
Tenseur bsn
3-vector shift, divided by N , of the rotating coordinates, .
Definition etoile.h:953
Tenseur d_psi
Gradient of (in the irrotational case) (Cartesian components with respect to ref_triad ).
Definition etoile.h:841
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:831
Tenseur beta_comp
Part of the logarithm of AN generated principaly by the companion star.
Definition etoile.h:877
Cmp ssjm1_psi
Effective source at the previous step for the resolution of the Poisson equation for the scalar by m...
Definition etoile.h:992
int nzet
Number of domains of *mp occupied by the star.
Definition etoile.h:435
Tenseur nnn
Total lapse function.
Definition etoile.h:512
const Eos & eos
Equation of state of the stellar matter.
Definition etoile.h:454
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
Tenseur ent
Log-enthalpy (relativistic case) or specific enthalpy (Newtonian case).
Definition etoile.h:460
Tenseur beta_auto
Logarithm of the part of the product AN generated principaly by by the star.
Definition etoile.h:509
Tenseur a_car
Total conformal factor .
Definition etoile.h:518
double unsurc2
: unsurc2=1 for a relativistic star, 0 for a Newtonian one.
Definition etoile.h:445
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_int_mod(int &n, int position=0)
Adds the address of a new modifiable int to the list.
Definition param.C:388
void add_int(const int &n, int position=0)
Adds the address of a new int to the list.
Definition param.C:249
Tensor field of valence 0 (or component of a tensorial field).
Definition scalar.h:393
Tbl & set_domain(int l)
Read/write of the value in a given domain.
Definition scalar.h:621
void set_outer_boundary(int l, double x)
Sets the value of the Scalar at the outer boundary of a given domain.
Basic array class.
Definition tbl.h:161
Class intended to describe tensors with a symmetry on the two last indices *** DEPRECATED : use class...
Definition tenseur.h:1256
Tensor handling *** DEPRECATED : use class Tensor instead ***.
Definition tenseur.h:304
const Base_vect * get_triad() const
Returns the vectorial basis (triad) on which the components are defined.
Definition tenseur.h:707
Cmp & set()
Read/write for a scalar (see also operator=(const Cmp&) ).
Definition tenseur.C:830
const Tenseur & derive_cov(const Metrique &met) const
Returns the covariant derivative of *this , with respect to met .
Definition tenseur.C:1560
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition tenseur.C:642
void annule(int l)
Sets the Tenseur to zero in a given domain.
Definition tenseur.C:906
const Tenseur & gradient() const
Returns the gradient of *this (Cartesian coordinates).
Definition tenseur.C:1548
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 set_triad(const Base_vect &new_triad)
Assigns a new vectorial basis (triad) of decomposition.
Definition tenseur.C:680
Values and coefficients of a (real-value) function.
Definition valeur.h:297
void ylm()
Computes the coefficients of *this.
Definition valeur_ylm.C:141
Tbl & set(int l)
Read/write of the value in a given domain (configuration space).
Definition valeur.h:373
void ylm_i()
Inverse of ylm().
void std_base_scal()
Sets the bases for spectral expansions (member base ) to the standard ones for a scalar.
Definition valeur.C:827
void annule_hard()
Sets the Valeur to zero in a hard way.
Definition valeur.C:726
Cmp sqrt(const Cmp &)
Square root.
Definition cmp_math.C:223
Cmp exp(const Cmp &)
Exponential.
Definition cmp_math.C:273
Tbl diffrel(const Cmp &a, const Cmp &b)
Relative difference between two Cmp (norme version).
Definition cmp_math.C:507
Tbl norme(const Cmp &)
Sums of the absolute values of all the values of the Cmp in each domain.
Definition cmp_math.C:484
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...
Tenseur contract(const Tenseur &, int id1, int id2)
Self contraction of two indices of a Tenseur .
Lorene prototypes.
Definition app_hor.h:67