LORENE
bhole.C
1/*
2 * Methods of class Bhole
3 *
4 * (see file bhole.h for documentation)
5 *
6 */
7
8/*
9 * Copyright (c) 2000-2001 Philippe Grandclement
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
30char bhole_C[] = "$Header: /cvsroot/Lorene/C++/Source/Bhole/bhole.C,v 1.15 2014/10/13 08:52:39 j_novak Exp $" ;
31
32/*
33 * $Id: bhole.C,v 1.15 2014/10/13 08:52:39 j_novak Exp $
34 * $Log: bhole.C,v $
35 * Revision 1.15 2014/10/13 08:52:39 j_novak
36 * Lorene classes and functions now belong to the namespace Lorene.
37 *
38 * Revision 1.14 2014/10/06 15:12:57 j_novak
39 * Modified #include directives to use c++ syntax.
40 *
41 * Revision 1.13 2007/04/26 13:16:21 f_limousin
42 * Correction of an error in the computation of grad_n_tot and grad_psi_tot
43 *
44 * Revision 1.12 2007/04/24 20:14:04 f_limousin
45 * Implementation of Dirichlet and Neumann BC for the lapse
46 *
47 * Revision 1.11 2006/09/25 10:01:48 p_grandclement
48 * Addition of N-dimensional Tbl
49 *
50 * Revision 1.10 2006/06/01 12:47:51 p_grandclement
51 * update of the Bin_ns_bh project
52 *
53 * Revision 1.9 2006/04/27 09:12:31 p_grandclement
54 * First try at irrotational black holes
55 *
56 * Revision 1.8 2005/08/29 15:10:13 p_grandclement
57 * Addition of things needed :
58 * 1) For BBH with different masses
59 * 2) Provisory files for the mixted binaries (Bh and NS) : THIS IS NOT
60 * WORKING YET !!!
61 *
62 * Revision 1.7 2003/12/16 05:27:07 k_taniguchi
63 * Change some arguments.
64 *
65 * Revision 1.6 2003/11/25 07:10:05 k_taniguchi
66 * Change some arguments from the class Etoile_bin to Et_bin_nsbh.
67 *
68 * Revision 1.5 2003/02/13 16:40:25 p_grandclement
69 * Addition of various things for the Bin_ns_bh project, non of them being
70 * completely tested
71 *
72 * Revision 1.4 2003/01/31 16:57:12 p_grandclement
73 * addition of the member Cmp decouple used to compute the K_ij auto, once
74 * the K_ij total is known
75 *
76 * Revision 1.3 2002/10/16 14:36:31 j_novak
77 * Reorganization of #include instructions of standard C++, in order to
78 * use experimental version 3 of gcc.
79 *
80 * Revision 1.2 2001/12/04 21:27:52 e_gourgoulhon
81 *
82 * All writing/reading to a binary file are now performed according to
83 * the big endian convention, whatever the system is big endian or
84 * small endian, thanks to the functions fwrite_be and fread_be
85 *
86 * Revision 1.1.1.1 2001/11/20 15:19:27 e_gourgoulhon
87 * LORENE
88 *
89 * Revision 2.23 2001/06/28 15:08:27 eric
90 * Ajout de la fonction area().
91 *
92 * Revision 2.22 2001/05/16 14:22:14 phil
93 * modif init_bhole_seul
94 *
95 * Revision 2.21 2001/05/16 11:33:08 phil
96 * ajout de init_bhole_seul
97 *
98 * Revision 2.20 2001/05/07 12:45:29 phil
99 * *** empty log message ***
100 *
101 * Revision 2.19 2001/05/07 09:30:34 phil
102 * *** empty log message ***
103 *
104 * Revision 2.18 2001/05/07 09:11:36 phil
105 * *** empty log message ***
106 *
107 * Revision 2.17 2001/04/26 12:04:27 phil
108 * *** empty log message ***
109 *
110 * Revision 2.16 2001/04/02 12:17:21 phil
111 * *** empty log message ***
112 *
113 * Revision 2.15 2001/02/28 13:22:43 phil
114 * vire kk_auto
115 *
116 * Revision 2.14 2001/01/29 14:30:37 phil
117 * ajout type_rotation
118 *
119 * Revision 2.13 2001/01/10 09:31:31 phil
120 * vire fait_kk_auto
121 *
122 * Revision 2.12 2001/01/09 13:38:26 phil
123 * ajout memebre kk_auto
124 *
125 * Revision 2.11 2000/12/21 10:48:59 phil
126 * retour a version 2.9
127 *
128 * Revision 2.9 2000/12/18 16:40:45 phil
129 * *** empty log message ***
130 *
131 * Revision 2.8 2000/12/14 10:44:48 phil
132 * ATTENTION : PASSAGE DE PHI A PSI
133 *
134 * Revision 2.7 2000/12/04 14:29:47 phil
135 * ajout de grad_n_tot
136 *
137 * Revision 2.6 2000/12/01 16:39:25 phil
138 * correction impoetation grad_phi_comp
139 *
140 * Revision 2.5 2000/12/01 16:21:25 phil
141 * modification init_bhole
142 *
143 * Revision 2.4 2000/12/01 16:17:42 phil
144 * correction base import de grad_phi_comp
145 *
146 * Revision 2.3 2000/12/01 16:14:33 phil
147 * *** empty log message ***
148 *
149 * Revision 2.2 2000/12/01 16:12:55 phil
150 * *** empty log message ***
151 *
152 * Revision 2.1 2000/10/23 09:22:07 phil
153 * rearrangement dans constructeurs.
154 *
155 * Revision 2.0 2000/10/20 09:18:47 phil
156 * *** empty log message ***
157 *
158 *
159 * $Header: /cvsroot/Lorene/C++/Source/Bhole/bhole.C,v 1.15 2014/10/13 08:52:39 j_novak Exp $
160 *
161 */
162
163//standard
164#include <cstdlib>
165#include <cmath>
166
167// Lorene
168#include "tenseur.h"
169#include "bhole.h"
170#include "proto.h"
171#include "utilitaires.h"
172#include "et_bin_nsbh.h"
173#include "graphique.h"
174
175// Constructeur standard
176namespace Lorene {
177Bhole::Bhole (Map_af& mpi) : mp(mpi),
178 rayon ((mp.get_alpha())[0]),
179 omega(0), omega_local(0), rot_state(COROT), boost (new double[3]), regul(0),
180 n_auto(mpi), n_comp(mpi), n_tot(mpi),
181 psi_auto(mpi), psi_comp(mpi), psi_tot(mpi),
182 grad_n_tot (mpi, 1, COV, mpi.get_bvect_cart()),
183 grad_psi_tot (mpi, 1, COV, mpi.get_bvect_cart()),
184 shift_auto(mpi, 1, CON, mpi.get_bvect_cart()),
185 taij_auto(mpi, 2, CON, mpi.get_bvect_cart()),
186 taij_comp(mpi, 2, CON, mpi.get_bvect_cart()),
187 taij_tot (mpi, 2, CON, mpi.get_bvect_cart()),
188 tkij_auto (mpi, 2, CON, mpi.get_bvect_cart()),
189 tkij_tot (mpi, 2, CON, mpi.get_bvect_cart()),
190 decouple (mpi) {
191 for (int i=0 ; i<3 ; i++)
192 boost[i] = 0 ;
193}
194
195
196// Constructeur par recopie
197Bhole::Bhole (const Bhole& source) :
198 mp(source.mp), rayon(source.rayon), omega (source.omega), omega_local(source.omega_local), rot_state(source.rot_state),
199 boost (new double [3]), regul(source.regul), n_auto(source.n_auto),
200 n_comp(source.n_comp), n_tot(source.n_tot), psi_auto(source.psi_auto),
201 psi_comp(source.psi_comp), psi_tot(source.psi_tot),
202 grad_n_tot (source.grad_n_tot), grad_psi_tot (source.grad_psi_tot),
203 shift_auto(source.shift_auto),
204 taij_auto (source.taij_auto),
205 taij_comp(source.taij_comp), taij_tot(source.taij_tot),
206 tkij_auto(source.tkij_auto), tkij_tot(source.tkij_tot), decouple(source.decouple) {
207
208 for (int i=0 ; i<3 ; i++)
209 boost[i] = source.boost[i] ;
210}
211
212// Destructeur
214 delete [] boost ;
215}
216
217// Affectation
218void Bhole::operator= (const Bhole& source) {
219
220 assert (&mp == &source.mp) ;
221
222 rayon = source.rayon ;
223 omega = source.omega ;
224 omega_local = source.omega_local ;
225 rot_state = source.rot_state ;
226 for (int i=0 ; i<3 ; i++)
227 boost[i] = source.boost[i] ;
228 regul = regul ;
229
230 n_auto = source.n_auto ;
231 n_comp = source.n_comp ;
232 n_tot = source.n_tot ;
233
234 psi_auto = source.psi_auto ;
235 psi_comp = source.psi_comp ;
236 psi_tot = source.psi_tot ;
237
238 grad_n_tot = source.grad_n_tot ;
239 grad_psi_tot = source.grad_psi_tot ;
240
241 shift_auto = source.shift_auto ;
242
243 taij_auto = source.taij_auto ;
244 taij_comp = source.taij_comp ;
245 taij_tot = source.taij_tot ;
246
247 tkij_auto = source.tkij_auto ;
248 tkij_tot = source.tkij_tot ;
249 decouple = source.decouple ;
250}
251
252// Importe le lapse du compagnon (Bhole case)
253
254void Bhole::fait_n_comp (const Bhole& comp) {
255
257 n_comp.set().import_symy(comp.n_auto()) ;
259
260 n_tot = n_comp + n_auto ;
262
263 Tenseur auxi (comp.n_auto.gradient()) ;
264 auxi.dec2_dzpuis() ;
265
266 Tenseur grad_comp (mp, 1, COV, *auxi.get_triad()) ;
267 grad_comp.set_etat_qcq() ;
268 grad_comp.set(0).import_symy(auxi(0)) ;
269 grad_comp.set(1).import_asymy(auxi(1)) ;
270 grad_comp.set(2).import_symy(auxi(2)) ;
271 grad_comp.set_std_base() ;
272 grad_comp.inc2_dzpuis() ;
273 grad_comp.change_triad(mp.get_bvect_cart()) ;
274
275 grad_n_tot = n_auto.gradient() + grad_comp ;
276}
277
278// Importe le facteur conforme du compagnon (Bhole case)
279
280void Bhole::fait_psi_comp (const Bhole& comp) {
281
283 psi_comp.set().import_symy(comp.psi_auto()) ;
285
288
289
290 Tenseur auxi (comp.psi_auto.gradient()) ;
291 auxi.dec2_dzpuis() ;
292
293 Tenseur grad_comp (mp, 1, COV, *auxi.get_triad()) ;
294 grad_comp.set_etat_qcq() ;
295 grad_comp.set(0).import_symy(auxi(0)) ;
296 grad_comp.set(1).import_asymy(auxi(1)) ;
297 grad_comp.set(2).import_symy(auxi(2)) ;
298 grad_comp.set_std_base() ;
299 grad_comp.inc2_dzpuis() ;
300 grad_comp.change_triad(mp.get_bvect_cart()) ;
301
302 grad_psi_tot = psi_auto.gradient() + grad_comp ;
303}
304
305
306// Importe le lapse du compagnon (NS case)
307void Bhole::fait_n_comp (const Et_bin_nsbh& comp) {
308
310 if (regul == 0)
311 n_comp.set().import(comp.get_n_auto()()) ;
312 else
313 n_comp.set().import_symy(comp.get_n_auto()()) ;
314
316
317 n_tot = n_comp + n_auto ;
319
320 Tenseur grad_comp (mp, 1, COV, mp.get_bvect_cart()) ;
321 Tenseur auxi (comp.get_d_n_auto()) ;
322 auxi.dec2_dzpuis() ;
323 grad_comp.set_etat_qcq() ;
324 if (regul == 0){
325 grad_comp.set(0).import(auxi(0)) ;
326 grad_comp.set(1).import(auxi(1)) ;
327 grad_comp.set(2).import(auxi(2)) ;
328 }
329 else{
330 grad_comp.set(0).import_symy(auxi(0)) ;
331 grad_comp.set(1).import_asymy(auxi(1)) ;
332 grad_comp.set(2).import_symy(auxi(2)) ;
333 }
334 grad_comp.set_std_base() ;
335 grad_comp.inc2_dzpuis() ;
336 grad_comp.set_triad( *(auxi.get_triad()) ) ;
337 grad_comp.change_triad(mp.get_bvect_cart()) ;
338
339 grad_n_tot = n_auto.gradient() + grad_comp ;
340}
341
342// Importe le facteur conforme du compagnon (Bhole case)
343
345
347 if (regul == 0)
348 psi_comp.set().import(comp.get_confpsi_auto()()) ;
349 else
352
355
356 Tenseur grad_comp (mp, 1, COV, mp.get_bvect_cart()) ;
357 Tenseur auxi (comp.get_d_confpsi_auto()) ;
358 auxi.dec2_dzpuis() ;
359 grad_comp.set_etat_qcq() ;
360 if (regul == 0){
361 grad_comp.set(0).import(auxi(0)) ;
362 grad_comp.set(1).import(auxi(1)) ;
363 grad_comp.set(2).import(auxi(2)) ;
364 }
365 else{
366 grad_comp.set(0).import_symy(auxi(0)) ;
367 grad_comp.set(1).import_asymy(auxi(1)) ;
368 grad_comp.set(2).import_symy(auxi(2)) ;
369 }
370 grad_comp.set_std_base() ;
371 grad_comp.inc2_dzpuis() ;
372 grad_comp.set_triad( *(auxi.get_triad()) ) ;
373 grad_comp.change_triad(mp.get_bvect_cart()) ;
374
375 grad_psi_tot = psi_auto.gradient() + grad_comp ;
376}
377
378// Calcul Taij auto (nul sur H que si la regularisation a ete faite)
380
381 if (shift_auto.get_etat() == ETATZERO)
383 else {
384 Tenseur grad (shift_auto.gradient()) ;
385 Tenseur trace (mp) ;
386 trace = grad(0, 0)+grad(1, 1)+grad(2, 2) ;
387
390 for (int i=0 ; i<3 ; i++) {
391 for (int j=i+1 ; j<3 ; j++)
392 taij_auto.set(i, j) = grad(i, j)+grad(j, i) ;
393 taij_auto.set(i, i) = 2*grad(i, i) -2./3.*trace() ;
394 }
395
396 for (int i=0 ; i<3 ; i++)
397 for (int j=0 ; j<3 ; j++)
398 taij_auto.set(i, j).raccord(1) ;
399 }
400}
401
402// Regularise le shift, untilisant le compagnon.
403void Bhole::regularise_shift (const Tenseur& shift_comp) {
404 regul = regle (shift_auto, shift_comp, omega, omega_local) ;
405}
406
407//Initialise a Schwartz
408// Compagnon a 0
410
411 Cmp auxi(mp) ;
412
413 auxi = 1./2.-rayon/mp.r ;
414 auxi.annule(0);
415 auxi.set_dzpuis(0) ;
416 n_auto = auxi;
418 n_auto.set().raccord(1) ;
419 n_comp =0.5;
422
423 auxi = 0.5+rayon/mp.r ;
424 auxi.annule(0);
425 auxi.set_dzpuis(0) ;
426 psi_auto = auxi;
428 psi_auto.set().raccord(1) ;
429 psi_comp = 0.5;
432
435
443}
444
446
447 Cmp auxi(mp) ;
448
449 auxi = (1-rayon/mp.r)/(1+rayon/mp.r) ;
450 auxi.annule(0);
451 auxi.set_val_inf(1) ;
452 auxi.set_dzpuis(0) ;
453 n_auto = auxi;
455 n_auto.set().raccord(1) ;
458
459 auxi = 1+rayon/mp.r ;
460 auxi.annule(0);
461 auxi.set_val_inf(1) ;
462 auxi.set_dzpuis(0) ;
463 psi_auto = auxi;
465 psi_auto.set().raccord(1) ;
468
471
479}
480
481// Sauvegarde dans fichier
482void Bhole::sauve (FILE* fich) const {
483
484 fwrite_be (&omega, sizeof(double), 1, fich) ;
485 fwrite_be (&omega_local, sizeof(double), 1, fich) ;
486 fwrite_be (&rot_state, sizeof(int), 1, fich) ;
487 fwrite_be (boost, sizeof(double), 3, fich) ;
488 fwrite_be (&regul, sizeof(double), 1, fich) ;
489
490 n_auto.sauve(fich) ;
491 psi_auto.sauve(fich) ;
492 shift_auto.sauve(fich) ;
493}
494
495//Constructeur par fichier
496Bhole::Bhole(Map_af& mpi, FILE* fich, bool old) :
497 mp(mpi), rayon ((mp.get_alpha())[0]),
498 boost (new double[3]), n_auto(mpi), n_comp(mpi), n_tot(mpi),
499 psi_auto(mpi), psi_comp(mpi), psi_tot(mpi),
500 grad_n_tot (mpi, 1, COV, mpi.get_bvect_cart()),
501 grad_psi_tot (mpi, 1, COV, mpi.get_bvect_cart()),
502 shift_auto(mpi, 1, CON, mpi.get_bvect_cart()),
503 taij_auto(mpi, 2, CON, mpi.get_bvect_cart()),
504 taij_comp(mpi, 2, CON, mpi.get_bvect_cart()),
505 taij_tot (mpi, 2, CON, mpi.get_bvect_cart()),
506 tkij_auto (mpi, 2, CON, mpi.get_bvect_cart()) ,
507 tkij_tot (mpi, 2, CON, mpi.get_bvect_cart()) ,
508 decouple (mpi) {
509
510 fread_be(&omega, sizeof(double), 1, fich) ;
511
512 if (!old) {
513 fread_be(&omega_local, sizeof(double), 1, fich) ;
514 fread_be(&rot_state, sizeof(int), 1, fich) ;
515 }
516 fread_be(boost, sizeof(double), 3, fich) ;
517 fread_be(&regul, sizeof(double), 1, fich) ;
518
519 Tenseur n_auto_file (mp, fich) ;
520 n_auto = n_auto_file ;
521
522 Tenseur psi_auto_file (mp, fich) ;
523 psi_auto = psi_auto_file ;
524
525 Tenseur shift_auto_file (mp, mp.get_bvect_cart(), fich) ;
526 shift_auto = shift_auto_file ;
527
530
534
538}
539
540
541 //---------------------------------//
542 // Computation of throat area //
543 //---------------------------------//
544
545double Bhole::area() const {
546
547 Cmp integrant( pow(psi_tot(), 4) ) ; // Psi^4
548
549 integrant.std_base_scal() ;
550 integrant.raccord(1) ;
551
552 return mp.integrale_surface(integrant, rayon) ;
553
554}
555
556 //---------------------------------//
557 // Moment local //
558 //---------------------------------//
559
560double Bhole::local_momentum() const {
561
562 Cmp xx (mp) ;
563 xx = mp.x ;
564 xx.std_base_scal() ;
565
566 Cmp yy (mp) ;
567 yy = mp.y ;
568 yy.std_base_scal() ;
569
570 Tenseur vecteur (mp, 1, CON, mp.get_bvect_cart()) ;
571 vecteur.set_etat_qcq() ;
572 for (int i=0 ; i<3 ; i++)
573 vecteur.set(i) = (-yy*tkij_tot(0, i)+xx*tkij_tot(1, i)) ;
574 vecteur.set_std_base() ;
575 vecteur.annule(mp.get_mg()->get_nzone()-1) ;
576 vecteur.change_triad (mp.get_bvect_spher()) ;
577
578 Cmp integrant (pow(psi_tot(), 6)*vecteur(0)) ;
579 integrant.std_base_scal() ;
580 double moment = mp.integrale_surface(integrant, rayon)/8/M_PI ;
581
582 return moment ;
583}
584
585}
Black hole.
Definition bhole.h:268
Tenseur psi_auto
Part of generated by the hole.
Definition bhole.h:290
Tenseur shift_auto
Part of generated by the hole.
Definition bhole.h:297
double omega_local
local angular velocity
Definition bhole.h:276
Tenseur psi_tot
Total .
Definition bhole.h:292
Tenseur psi_comp
Part of generated by the companion hole.
Definition bhole.h:291
Tenseur_sym tkij_auto
Auto .
Definition bhole.h:307
double omega
Angular velocity in LORENE's units.
Definition bhole.h:275
Tenseur grad_psi_tot
Total gradient of .
Definition bhole.h:295
double regul
Intensity of the correction on the shift vector.
Definition bhole.h:284
Tenseur_sym taij_auto
Part of generated by the hole.
Definition bhole.h:299
int rot_state
State of rotation.
Definition bhole.h:277
Tenseur grad_n_tot
Total gradient of N .
Definition bhole.h:294
Tenseur n_auto
Part of N generated by the hole.
Definition bhole.h:286
~Bhole()
Destructor.
Definition bhole.C:213
Tenseur_sym taij_comp
Part of generated by the companion hole.
Definition bhole.h:300
void init_bhole()
Sets the values of the fields to :
Definition bhole.C:409
double * boost
Lapse on the horizon.
Definition bhole.h:283
double rayon
Radius of the horizon in LORENE's units.
Definition bhole.h:274
void init_bhole_seul()
Initiates for a single the black hole.
Definition bhole.C:445
double area() const
Computes the area of the throat.
Definition bhole.C:545
Tenseur n_comp
Part of N generated by the companion hole.
Definition bhole.h:287
Map_af & mp
Affine mapping.
Definition bhole.h:273
Tenseur_sym taij_tot
Total , which must be zero on the horizon of the regularisation on the shift has been done.
Definition bhole.h:305
Tenseur n_tot
Total N .
Definition bhole.h:288
void sauve(FILE *fich) const
Write on a file.
Definition bhole.C:482
void fait_psi_comp(const Bhole &comp)
Imports the part of due to the companion hole comp .
Definition bhole.C:280
Bhole(Map_af &mapping)
Standard constructor.
Definition bhole.C:177
Tenseur_sym tkij_tot
Total .
Definition bhole.h:308
void fait_taij_auto()
Calculates the part of generated by shift_auto .
Definition bhole.C:379
void fait_n_comp(const Bhole &comp)
Imports the part of N due to the companion hole comp .
Definition bhole.C:254
void operator=(const Bhole &)
Affectation.
Definition bhole.C:218
Cmp decouple
Function used to construct the part of generated by the hole from the total .
Definition bhole.h:318
void regularise_shift(const Tenseur &comp)
Corrects shift_auto in such a way that the total is equal to zero in the horizon,...
Definition bhole.C:403
Component of a tensorial field *** DEPRECATED : use class Scalar instead ***.
Definition cmp.h:446
void import(const Cmp &ci)
Assignment to another Cmp defined on a different mapping.
Definition cmp_import.C:73
void std_base_scal()
Sets the spectral bases of the Valeur va to the standard ones for a scalar.
Definition cmp.C:644
void annule(int l)
Sets the Cmp to zero in a given domain.
Definition cmp.C:348
void set_etat_zero()
Sets the logical state to ETATZERO (zero).
Definition cmp.C:289
void set_dzpuis(int)
Set a value to dzpuis.
Definition cmp.C:654
void import_symy(const Cmp &ci)
Assignment to another Cmp defined on a different mapping.
void set_val_inf(double val)
Sets the value of the Cmp to val at infinity.
Definition cmp_manip.C:126
void raccord(int n)
Performs the matching of the nucleus with respect to the first shell.
void import_asymy(const Cmp &ci)
Assignment to another Cmp defined on a different mapping.
Class for a star in a NS-BH binary system.
Definition et_bin_nsbh.h:79
const Tenseur & get_d_n_auto() const
Returns the gradient of {\tt n_auto} (Cartesian components with respect to {\tt ref_triad})
const Tenseur & get_confpsi_auto() const
Returns the part of the conformal factor $\Psi$ generated principaly by the star.
const Tenseur & get_n_auto() const
Returns the part of the lapse {\it N} generated principaly by the star.
const Tenseur & get_d_confpsi_auto() const
Returns the gradient of {\tt confpsi_auto} (Cartesian components with respect to {\tt ref_triad})
Affine radial mapping.
Definition map.h:2027
double integrale_surface(const Cmp &ci, double rayon) const
Performs the surface integration of ci on the sphere of radius rayon .
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
Coord y
y coordinate centered on the grid
Definition map.h:727
Coord r
r coordinate centered on the grid
Definition map.h:718
const Base_vect_spher & get_bvect_spher() const
Returns the orthonormal vectorial basis associated with the coordinates of the mapping.
Definition map.h:783
Coord x
x coordinate centered on the grid
Definition map.h:726
const Mg3d * get_mg() const
Gives the Mg3d on which the mapping is defined.
Definition map.h:765
int get_nzone() const
Returns the number of domains.
Definition grilles.h:448
Tensor handling *** DEPRECATED : use class Tensor instead ***.
Definition tenseur.h:301
const Base_vect * get_triad() const
Returns the vectorial basis (triad) on which the components are defined.
Definition tenseur.h:704
Cmp & set()
Read/write for a scalar (see also operator=(const Cmp&) ).
Definition tenseur.C:824
void sauve(FILE *) const
Save in a file.
Definition tenseur.C:1325
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 set_etat_zero()
Sets the logical state to ETATZERO (zero state).
Definition tenseur.C:645
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 set_triad(const Base_vect &new_triad)
Assigns a new vectorial basis (triad) of decomposition.
Definition tenseur.C:674
int get_etat() const
Returns the logical state.
Definition tenseur.h:707
Cmp pow(const Cmp &, int)
Power .
Definition cmp_math.C:348
int fread_be(int *aa, int size, int nb, FILE *fich)
Reads integer(s) from a binary file according to the big endian convention.
Definition fread_be.C:69
int fwrite_be(const int *aa, int size, int nb, FILE *fich)
Writes integer(s) into a binary file according to the big endian convention.
Definition fwrite_be.C:70
Lorene prototypes.
Definition app_hor.h:64