LORENE
et_rot_mag_mag_plus.C
1/*
2 * Computes magnetic fields and derived quantities for rotating equilibrium
3 *
4 * (see file et_rot_mag.h for documentation)
5 *
6 */
7
8/*
9 * Copyright (c) 2002 Emmanuel Marcq
10 * Copyright (c) 2002 Jerome Novak
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
30char et_rot_mag_mag_plus_C[] = "$Header: /cvsroot/Lorene/C++/Source/Etoile/et_rot_mag_mag_plus.C,v 1.3 2014/10/13 08:52:58 j_novak Exp $" ;
31
32/*
33 * $Id: et_rot_mag_mag_plus.C,v 1.3 2014/10/13 08:52:58 j_novak Exp $
34 * $Log: et_rot_mag_mag_plus.C,v $
35 * Revision 1.3 2014/10/13 08:52:58 j_novak
36 * Lorene classes and functions now belong to the namespace Lorene.
37 *
38 * Revision 1.2 2014/10/06 15:13:09 j_novak
39 * Modified #include directives to use c++ syntax.
40 *
41 * Revision 1.1 2012/08/12 17:48:35 p_cerda
42 * Magnetstar: New classes for magnetstar. Allowing for non-equatorial symmetry in Etoile et al. Adding B_phi in Et_rot_mag.
43 *
44 * Revision 1.14 2005/06/03 15:31:56 j_novak
45 * Better computation when more than one point in phi.
46 *
47 * Revision 1.13 2003/10/03 15:58:47 j_novak
48 * Cleaning of some headers
49 *
50 * Revision 1.12 2002/09/09 13:00:39 e_gourgoulhon
51 * Modification of declaration of Fortran 77 prototypes for
52 * a better portability (in particular on IBM AIX systems):
53 * All Fortran subroutine names are now written F77_* and are
54 * defined in the new file C++/Include/proto_f77.h.
55 *
56 * Revision 1.11 2002/06/05 15:15:59 j_novak
57 * The case of non-adapted mapping is treated.
58 * parmag.d and parrot.d have been merged.
59 *
60 * Revision 1.10 2002/06/03 13:23:16 j_novak
61 * The case when the mapping is not adapted is now treated
62 *
63 * Revision 1.9 2002/06/03 13:00:45 e_marcq
64 *
65 * conduc parameter read in parmag.d
66 *
67 * Revision 1.7 2002/05/20 10:31:59 j_novak
68 * *** empty log message ***
69 *
70 * Revision 1.6 2002/05/17 15:08:01 e_marcq
71 *
72 * Rotation progressive plug-in, units corrected, Q and a_j new member data
73 *
74 * Revision 1.5 2002/05/16 10:02:09 j_novak
75 * Errors in stress energy tensor corrected
76 *
77 * Revision 1.4 2002/05/15 09:54:00 j_novak
78 * First operational version
79 *
80 * Revision 1.3 2002/05/14 13:38:36 e_marcq
81 *
82 *
83 * Unit update, new outputs
84 *
85 * Revision 1.1 2002/05/10 09:26:52 j_novak
86 * Added new class Et_rot_mag for magnetized rotating neutron stars (under development)
87 *
88 *
89 * $Header: /cvsroot/Lorene/C++/Source/Etoile/et_rot_mag_mag_plus.C,v 1.3 2014/10/13 08:52:58 j_novak Exp $
90 *
91 */
92
93// Headers C
94#include <cstdlib>
95#include <cmath>
96
97// Headers Lorene
98#include "et_rot_mag.h"
99#include "utilitaires.h"
100#include "param.h"
101#include "proto_f77.h"
102#include "graphique.h"
103#include "unites.h"
104
105namespace Lorene {
106
107// Algo du papier de 1995
108
109 using namespace Unites ;
110
111void Et_rot_mag::magnet_comput_plus(const int adapt_flag,
112 const int initial_j,
113 const Tbl an_j,
114 Cmp (*f_j)(const Cmp&, const Tbl),
115 const Tbl bn_j,
116 Cmp (*g_j)(const Cmp&, const Tbl),
117 Cmp (*N_j)(const Cmp&, const Tbl),
118 Param& par_poisson_At,
119 Param& par_poisson_Avect){
120 double relax_mag = 1.0 ;
121
122 int Z = mp.get_mg()->get_nzone();
123
124 if(is_conduct()) {
125 bool adapt(adapt_flag) ;
126 /****************************************************************
127 * Assertion that all zones have same number of points in theta
128 ****************************************************************/
129 int nt = mp.get_mg()->get_nt(nzet-1) ;
130 for (int l=0; l<Z; l++) assert(mp.get_mg()->get_nt(l) == nt) ;
131
132 Tbl Rsurf(nt) ;
133 Rsurf.set_etat_qcq() ;
134 mp.r.fait() ;
135 mp.tet.fait() ;
136 Mtbl* theta = mp.tet.c ;
137 const Map_radial* mpr = dynamic_cast<const Map_radial*>(&mp) ;
138 assert (mpr != 0x0) ;
139 for (int j=0; j<nt; j++)
140 Rsurf.set(j) = mpr->val_r_jk(l_surf()(0,j), xi_surf()(0,j), j, 0) ;
141
142
143 // Calcul de A_0t dans l'etoile (conducteur parfait)
144
145 Cmp A_0t(- omega * A_phi) ;
146 A_0t.annule(nzet,Z-1) ;
147
148 Tenseur ATTENS(A_t) ;
149 Tenseur APTENS(A_phi) ;
150 Tenseur BMN(-logn) ;
151 BMN = BMN + log(bbb) ;
152 BMN.set_std_base() ;
153
154
156 nphi.gradient_spher())());
158 nphi.gradient_spher())()) ;
160 BMN.gradient_spher())()
162 BMN.gradient_spher())()) ;
163
164 Cmp ATANT(A_phi.srdsdt()); // Constrction par copie pour mapping
165
166 ATANT.va = ATANT.va.mult_ct().ssint() ;
167
168 Cmp ttnphi(tnphi()) ;
169 ttnphi.mult_rsint() ;
170 Cmp BLAH(- b_car()/(nnn()*nnn())*ttnphi*grad1) ;
171 BLAH -= (1+b_car()/(nnn()*nnn())*tnphi()*tnphi())*grad2 ;
172 Cmp nphisr(nphi()) ;
173 nphisr.div_r() ;
174 Cmp npgrada(2*nphisr*(A_phi.dsdr()+ATANT )) ;
175 npgrada.inc2_dzpuis() ;
176 BLAH -= grad3 + npgrada ;
177 Cmp gtt(-nnn()*nnn()+b_car()*tnphi()*tnphi()) ;
178 Cmp gtphi( - b_car()*ttnphi) ;
179
180 // Calcul de j_t grace a Maxwell-Gauss
181 Cmp tmp(((BLAH - A_0t.laplacien())/a_car() - gtphi*j_phi)
182 / gtt);
183 tmp.annule(nzet, Z-1) ;
184 if (adapt) {
185 j_t = tmp ;
186 }
187 else {
188 j_t.allocate_all() ;
189 for (int j=0; j<nt; j++)
190 for (int l=0; l<nzet; l++)
191 for (int i=0; i<mp.get_mg()->get_nr(l); i++)
192 j_t.set(l,0,j,i) = ( (*mp.r.c)(l,0,j,i) > Rsurf(j) ?
193 0. : tmp(l,0,j,i) ) ;
194 j_t.annule(nzet,Z-1) ;
195 }
197
198 // Calcul du courant j_phi
199
200 Tbl maxA_phi = max(abs(A_phi));
201
202 if (maxA_phi.set(0) == 0) {
203
204 cout << "Initializing j_phi" << endl;
205
206 double aini = 0.;
207 int nd = an_j.get_dim(0);
208 for (int i=0; i<nd - 4;i++){
209 aini = aini + fabs(an_j(i));
210 }
211 aini = aini + fabs (an_j(3)) + fabs (an_j(5));
212 double bini = 0.;
213 nd = bn_j.get_dim(0);
214 for (int i=0; i<nd;i++){
215 bini = bini + fabs(bn_j(i));
216 }
217
218 switch (initial_j) {
219 case 0 :
220 j_phi = (ener() + press())*aini + bini ;
222 j_phi.annule(nzet,Z-1) ;
224 break;
225
226 case 1 :
227 j_phi = (ener() + press())*aini + bini;
229 j_phi.mult_rsint();
230 j_phi.annule(nzet,Z-1) ;
232 break;
233
234 case 2 :
235 j_phi = (ener() + press())*aini + bini ;
237 j_phi.mult_cost();
238 j_phi.mult_rsint();
239 j_phi.set_dzpuis(2);
240 j_phi.mult_r();
241 j_phi.annule(nzet,Z-1) ;
243 break;
244
245 default :
246 cout << "ERROR" << endl;
247 cout << "initial_j = " << initial_j << endl;
248 abort();
249 }
250
251
252 }else{
253
254
255
256 double maxA_phi_surf = 0.;
257 for (int j = 0; j < nt ; j++ ) {
258 // double maxA_phi_surf_tmp = fabs(A_phi(nzet-1,0,j,mp.get_mg()->get_nr(nzet-1)-1));
259 double maxA_phi_surf_tmp = fabs(A_phi(nzet,0,j,0));
260 maxA_phi_surf= max(maxA_phi_surf, maxA_phi_surf_tmp);
261 }
262
263 Cmp A_phi_scaled = A_phi / maxA_phi.set(0);
264 Cmp A_phi_scaled2 = (A_phi - maxA_phi_surf)/ maxA_phi.set(0);
265
266 A_phi_scaled2.allocate_all() ;
267 for (int j=0; j<nt; j++) {
268 for (int l=0; l<Z; l++) {
269 for (int i=0; i<mp.get_mg()->get_nr(l); i++) {
270 if (A_phi_scaled2(l,0,j,i) < 0.) {
271 A_phi_scaled2.set(l,0,j,i) = 0 ;
272 }
273 }
274 }
275 }
276 A_phi_scaled2.std_base_scal() ;
277
278
279 Cmp j_phi_ff = g_j (A_phi_scaled2, bn_j)/maxA_phi.set(0);
280 //Cmp j_phi_ff = A_phi_scaled ;
281 j_phi_ff.std_base_scal() ;
282 j_phi_ff.div_rsint();
283 j_phi_ff.div_rsint();
284 j_phi_ff = j_phi_ff / b_car() / nnn() / nnn();
285
286 j_phi = omega * j_t
287 + (ener() + press())*f_j(A_phi_scaled, an_j)/maxA_phi.set(0)/g_si
288 + j_phi_ff;
289
290 j_phi.annule(nzet,Z-1) ;
292
293
294
295 B_phi = N_j (A_phi_scaled2, bn_j);
297 B_phi = B_phi / nnn();
298
299
300
301 }
302 // des_coupe_y(A_phi, 0., nzet, "Magnetic field") ;
303 // des_coupe_y(j_phi, 0., nzet, "Current",0x0,1.2,true,30) ;
304
305 // Resolution de Maxwell Ampere (-> A_phi)
306 // Calcul des termes sources avec A-t du pas precedent.
307
309 BMN.gradient_spher())());
310
311 Tenseur source_tAphi(mp, 1, CON, mp.get_bvect_spher()) ;
312
313 source_tAphi.set_etat_qcq() ;
314 Cmp tjphi(j_phi) ;
315 tjphi.mult_rsint() ;
316 Cmp tgrad1(grad1) ;
317 tgrad1.mult_rsint() ;
318 Cmp d_grad4(grad4) ;
319 d_grad4.div_rsint() ;
320 source_tAphi.set(0)=0 ;
321 source_tAphi.set(1)=0 ;
322 source_tAphi.set(2)= -b_car()*a_car()*(tjphi-tnphi()*j_t)
323 + b_car()/(nnn()*nnn())*(tgrad1+tnphi()*grad2)+d_grad4 ;
324
325 source_tAphi.change_triad(mp.get_bvect_cart());
326
327 Tenseur WORK_VECT(mp, 1, CON, mp.get_bvect_cart()) ;
328 WORK_VECT.set_etat_qcq() ;
329 for (int i=0; i<3; i++) {
330 WORK_VECT.set(i) = 0 ;
331 }
332 Tenseur WORK_SCAL(mp) ;
333 WORK_SCAL.set_etat_qcq() ;
334 WORK_SCAL.set() = 0 ;
335
336 double lambda_mag = 0. ; // No 3D version !
337
338 Tenseur AVECT(source_tAphi) ;
339 if (source_tAphi.get_etat() != ETATZERO) {
340
341 for (int i=0; i<3; i++) {
342 if(source_tAphi(i).dz_nonzero()) {
343 assert( source_tAphi(i).get_dzpuis() == 4 ) ;
344 }
345 else{
346 (source_tAphi.set(i)).set_dzpuis(4) ;
347 }
348 }
349
350 }
351 source_tAphi.poisson_vect(lambda_mag, par_poisson_Avect, AVECT, WORK_VECT,
352 WORK_SCAL) ;
354 Cmp A_phi_n(AVECT(2));
355 A_phi_n.mult_rsint() ;
356
357 // Resolution de Maxwell-Ampere : A_1
358
359 Cmp source_A_1t(-a_car()*(j_t*gtt + j_phi*gtphi) + BLAH);
360
361 Cmp A_1t(mp);
362 A_1t = 0 ;
363 source_A_1t.poisson(par_poisson_At, A_1t) ;
364
365 int L = mp.get_mg()->get_nt(0);
366
367 Tbl MAT(L,L) ;
368 Tbl MAT_PHI(L,L);
369 Tbl VEC(L) ;
370
371 MAT.set_etat_qcq() ;
372 VEC.set_etat_qcq() ;
373 MAT_PHI.set_etat_qcq() ;
374
375 Tbl leg(L,2*L) ;
376 leg.set_etat_qcq() ;
377
378 Cmp psi(mp);
379 Cmp psi2(mp);
380 psi.allocate_all() ;
381 psi2.allocate_all() ;
382
383 for (int p=0; p<mp.get_mg()->get_np(0); p++) {
384 // leg[k,l] : legendre_l(cos(theta_k))
385 // Construction par recurrence de degre 2
386 for(int k=0;k<L;k++){
387 for(int l=0;l<2*L;l++){
388
389 if(l==0) leg.set(k,l)=1. ;
390 if(l==1) leg.set(k,l)=cos((*theta)(l_surf()(p,k),p,k,0)) ;
391 if(l>=2) leg.set(k,l) = double(2*l-1)/double(l)
392 * cos((*theta)(l_surf()(p,k),p,k,0))
393 * leg(k,l-1)-double(l-1)/double(l)*leg(k,l-2) ;
394 }
395 }
396
397 for(int k=0;k<L;k++){
398
399 // Valeurs a la surface trouvees via va.val_point_jk(l,xisurf,k,p)
400
401 VEC.set(k) = A_0t.va.val_point_jk(l_surf()(p,k), xi_surf()(p,k), k, p)
402 -A_1t.va.val_point_jk(l_surf()(p,k), xi_surf()(p,k), k, p);
403
404 for(int l=0;l<L;l++) MAT.set(l,k) = leg(k,2*l)/pow(Rsurf(k),2*l+1);
405
406 }
407 // appel fortran :
408
409 int* IPIV=new int[L] ;
410 int INFO ;
411
412 Tbl MAT_SAVE(MAT) ;
413 Tbl VEC2(L) ;
414 VEC2.set_etat_qcq() ;
415 int un = 1 ;
416
417 F77_dgesv(&L, &un, MAT.t, &L, IPIV, VEC.t, &L, &INFO) ;
418
419 // coeffs a_l dans VEC
420
421 for(int k=0;k<L;k++) {VEC2.set(k)=1. ; }
422
423 F77_dgesv(&L, &un, MAT_SAVE.t, &L, IPIV, VEC2.t, &L, &INFO) ;
424
425 delete [] IPIV ;
426
427 for(int nz=0;nz < Z; nz++){
428 for(int i=0;i< mp.get_mg()->get_nr(nz);i++){
429 for(int k=0;k<L;k++){
430 psi.set(nz,p,k,i) = 0. ;
431 psi2.set(nz,p,k,i) = 0. ;
432 for(int l=0;l<L;l++){
433 psi.set(nz,p,k,i) += VEC(l)*leg(k,2*l) /
434 pow((*mp.r.c)(nz,p,k,i),2*l+1);
435 psi2.set(nz,p,k,i) += VEC2(l)*leg(k,2*l)/
436 pow((*mp.r.c)(nz, p, k,i),2*l+1);
437 }
438 }
439 }
440 }
441 }
442 psi.std_base_scal() ;
443 psi2.std_base_scal() ;
444
445 assert(psi.get_dzpuis() == 0) ;
446 int dif = A_1t.get_dzpuis() ;
447 if (dif > 0) {
448 for (int d=0; d<dif; d++) A_1t.dec_dzpuis() ;
449 }
450
451 if (adapt) {
452 Cmp A_t_ext(A_1t + psi) ;
453 A_t_ext.annule(0,nzet-1) ;
454 A_0t += A_t_ext ;
455 }
456 else {
457 tmp = A_0t ;
458 A_0t.allocate_all() ;
459 for (int j=0; j<nt; j++)
460 for (int l=0; l<Z; l++)
461 for (int i=0; i<mp.get_mg()->get_nr(l); i++)
462 A_0t.set(l,0,j,i) = ( (*mp.r.c)(l,0,j,i) > Rsurf(j) ?
463 A_1t(l,0,j,i) + psi(l,0,j,i) : tmp(l,0,j,i) ) ;
464 }
465 A_0t.std_base_scal() ;
466
467 Valeur** asymp = A_0t.asymptot(1) ;
468
469 double Q_0 = -4*M_PI*(*asymp[1])(Z-1,0,0,0) ; // utilise A_0t plutot que E
470 delete asymp[0] ;
471 delete asymp[1] ;
472
473 delete [] asymp ;
474
475 asymp = psi2.asymptot(1) ;
476
477 double Q_2 = -4*M_PI*(*asymp[1])(Z-1,0,0,0) ; // A_2t = psi2 a l'infini
478 delete asymp[0] ;
479 delete asymp[1] ;
480
481 delete [] asymp ;
482
483 // solution definitive de A_t:
484
485 double C = (Q-Q_0)/Q_2 ;
486
487 assert(psi2.get_dzpuis() == 0) ;
488 dif = A_0t.get_dzpuis() ;
489 if (dif > 0) {
490 for (int d=0; d<dif; d++) A_0t.dec_dzpuis() ;
491 }
492 Cmp A_t_n(mp) ;
493 if (adapt) {
494 A_t_n = A_0t + C ;
495 Cmp A_t_ext(A_0t + C*psi2) ;
496 A_t_ext.annule(0,nzet-1) ;
497 A_t_n.annule(nzet,Z-1) ;
498 A_t_n += A_t_ext ;
499 }
500 else {
501 A_t_n.allocate_all() ;
502 for (int j=0; j<nt; j++)
503 for (int l=0; l<Z; l++)
504 for (int i=0; i<mp.get_mg()->get_nr(l); i++)
505 A_t_n.set(l,0,j,i) = ( (*mp.r.c)(l,0,j,i) > Rsurf(j) ?
506 A_0t(l,0,j,i) + C*psi2(l,0,j,i) :
507 A_0t(l,0,j,i) + C ) ;
508 }
509 A_t_n.std_base_scal() ;
510
511 asymp = A_t_n.asymptot(1) ;
512
513 delete asymp[0] ;
514 delete asymp[1] ;
515
516 delete [] asymp ;
517 A_t = relax_mag*A_t_n + (1.-relax_mag)*A_t ;
518 A_phi = relax_mag*A_phi_n + (1. - relax_mag)*A_phi ;
519
520
521 }else{
522 cout << "not implemented" << endl;
523 abort();
524
525 }
526
527
528}
529
530
531
532
533
534
535
536
537
538
539
540
541
542}
Component of a tensorial field *** DEPRECATED : use class Scalar instead ***.
Definition cmp.h:446
void mult_rsint()
Multiplication by .
void allocate_all()
Sets the logical state to ETATQCQ (ordinary state) and performs the memory allocation of all the elem...
Definition cmp.C:323
void div_r()
Division by r everywhere.
Definition cmp_r_manip.C:78
void dec_dzpuis()
Decreases by 1 the value of dzpuis and changes accordingly the values of the Cmp in the external comp...
Valeur va
The numerical value of the Cmp
Definition cmp.h:464
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
const Cmp & laplacien(int zec_mult_r=4) const
Returns the Laplacian of *this.
Definition cmp_deriv.C:242
int get_dzpuis() const
Returns dzpuis.
Definition cmp.h:903
void mult_r()
Multiplication by r everywhere.
Definition cmp_r_manip.C:91
void inc2_dzpuis()
Increases by 2 the value of dzpuis and changes accordingly the values of the Cmp in the external comp...
Tbl & set(int l)
Read/write of the value in a given domain.
Definition cmp.h:724
void set_dzpuis(int)
Set a value to dzpuis.
Definition cmp.C:654
const Cmp & srdsdt() const
Returns of *this .
Definition cmp_deriv.C:105
Valeur ** asymptot(int n, const int flag=0) const
Asymptotic expansion at r = infinity.
void mult_cost()
Multiplication by \f\cos\theta\f$.
Cmp poisson() const
Solves the scalar Poisson equation with *this as a source.
Definition cmp_pde.C:94
const Cmp & dsdr() const
Returns of *this .
Definition cmp_deriv.C:84
void fait() const
Computes, at each point of the grid, the value of the coordinate or mapping derivative represented by...
Definition coord.C:116
Mtbl * c
The coordinate values at each grid point.
Definition coord.h:97
Cmp j_phi
-component of the current 4-vector
Definition et_rot_mag.h:159
Cmp A_phi
-component of the electromagnetic potential 1-form divided by .
Definition et_rot_mag.h:155
Cmp A_t
t-component of the elecctromagnetic potential 1-form, divided by .
Definition et_rot_mag.h:150
bool is_conduct() const
Tells if the star is made of conducting or isolating material.
Definition et_rot_mag.h:241
Cmp j_t
t-component of the current 4-vector
Definition et_rot_mag.h:158
virtual void magnet_comput_plus(const int adapt_flag, const int initial_j, const Tbl an_j, Cmp(*f_j)(const Cmp &x, const Tbl), const Tbl bn_j, Cmp(*g_j)(const Cmp &x, const Tbl), Cmp(*N_j)(const Cmp &x, const Tbl), Param &par_poisson_At, Param &par_poisson_Avect)
Computes the electromagnetic quantities solving the Maxwell equations (6) and (7) of [Bocquet,...
double Q
In the case of a perfect conductor, the requated baryonic charge.
Definition et_rot_mag.h:179
Cmp B_phi
-component of the magnetic field
Definition et_rot_mag.h:157
double omega
Rotation angular velocity ([f_unit] )
Definition etoile.h:1501
Tenseur & logn
Metric potential = logn_auto.
Definition etoile.h:1521
virtual const Itbl & l_surf() const
Description of the stellar surface: returns a 2-D Itbl containing the values of the domain index l on...
Tenseur nphi
Metric coefficient .
Definition etoile.h:1510
Tenseur bbb
Metric factor B.
Definition etoile.h:1504
Tenseur b_car
Square of the metric factor B.
Definition etoile.h:1507
Tenseur tnphi
Component of the shift vector.
Definition etoile.h:1515
int nzet
Number of domains of *mp occupied by the star.
Definition etoile.h:432
Tenseur nnn
Total lapse function.
Definition etoile.h:509
Map & mp
Mapping associated with the star.
Definition etoile.h:429
Tenseur ener
Total energy density in the fluid frame.
Definition etoile.h:460
Tenseur press
Fluid pressure.
Definition etoile.h:461
const Tbl & xi_surf() const
Description of the stellar surface: returns a 2-D Tbl containing the values of the radial coordinat...
Tenseur a_car
Total conformal factor .
Definition etoile.h:515
Base class for pure radial mappings.
Definition map.h:1536
virtual double val_r_jk(int l, double xi, int j, int k) const =0
Returns the value of the radial coordinate r for a given and a given collocation point in in a give...
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 r
r coordinate centered on the grid
Definition map.h:718
Coord tet
coordinate centered on the grid
Definition map.h:719
const Base_vect_spher & get_bvect_spher() const
Returns the orthonormal vectorial basis associated with the coordinates of the mapping.
Definition map.h:783
const Mg3d * get_mg() const
Gives the Mg3d on which the mapping is defined.
Definition map.h:765
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
Multi-domain array.
Definition mtbl.h:118
Parameter storage.
Definition param.h:125
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
double * t
The array of double.
Definition tbl.h:173
int get_dim(int i) const
Gives the i-th dimension (ie dim.dim[i])
Definition tbl.h:403
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 set_std_base()
Set the standard spectal basis of decomposition for each component.
Definition tenseur.C:1170
const Tenseur & gradient_spher() const
Returns the gradient of *this (Spherical coordinates) (scalar field only).
Definition tenseur.C:1548
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(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:707
Values and coefficients of a (real-value) function.
Definition valeur.h:287
const Valeur & mult_ct() const
Returns applied to *this.
double val_point_jk(int l, double x, int j, int k) const
Computes the value of the field represented by *this at an arbitrary point in , but collocation point...
Definition valeur.C:900
const Valeur & ssint() const
Returns of *this.
Tbl max(const Cmp &)
Maximum values of a Cmp in each domain.
Definition cmp_math.C:435
Cmp pow(const Cmp &, int)
Power .
Definition cmp_math.C:348
Cmp cos(const Cmp &)
Cosine.
Definition cmp_math.C:94
Cmp abs(const Cmp &)
Absolute value.
Definition cmp_math.C:410
Cmp log(const Cmp &)
Neperian logarithm.
Definition cmp_math.C:296
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.