101#include "et_rot_mag.h"
102#include "utilitaires.h"
104#include "proto_f77.h"
105#include "graphique.h"
112 using namespace Unites ;
121 Param& par_poisson_At,
122 Param& par_poisson_Avect){
123 double relax_mag = 1.0 ;
125 int Z =
mp.get_mg()->get_nzone();
128 bool adapt(adapt_flag) ;
132 int nt =
mp.get_mg()->get_nt(
nzet-1) ;
133 for (
int l=0; l<Z; l++) assert(
mp.get_mg()->get_nt(l) == nt) ;
141 assert (mpr != 0x0) ;
142 for (
int j=0; j<nt; j++)
159 nphi.gradient_spher())());
161 nphi.gradient_spher())()) ;
177 Cmp npgrada(2*nphisr*(
A_phi.dsdr()+ATANT )) ;
179 BLAH -= grad3 + npgrada ;
192 for (
int j=0; j<nt; j++)
193 for (
int l=0; l<
nzet; l++)
194 for (
int i=0; i<
mp.get_mg()->get_nr(l); i++)
195 j_t.set(l,0,j,i) = ( (*
mp.r.c)(l,0,j,i) > Rsurf(j) ?
196 0. : tmp(l,0,j,i) ) ;
199 j_t.std_base_scal() ;
205 if (maxA_phi.
set(0) == 0) {
207 cout <<
"Initializing j_phi" << endl;
211 for (
int i=0; i<nd - 4;i++){
212 aini = aini + fabs(an_j(i));
214 aini = aini + fabs (an_j(3)) + fabs (an_j(5));
217 for (
int i=0; i<nd;i++){
218 bini = bini + fabs(bn_j(i));
224 j_phi.std_base_scal() ;
226 j_phi.std_base_scal() ;
231 j_phi.std_base_scal() ;
234 j_phi.std_base_scal() ;
239 j_phi.std_base_scal() ;
245 j_phi.std_base_scal() ;
249 cout <<
"ERROR" << endl;
250 cout <<
"initial_j = " << initial_j << endl;
259 double maxA_phi_surf = 0.;
260 for (
int j = 0; j < nt ; j++ ) {
262 double maxA_phi_surf_tmp = fabs(
A_phi(
nzet,0,j,0));
263 maxA_phi_surf=
max(maxA_phi_surf, maxA_phi_surf_tmp);
267 Cmp A_phi_scaled2 = (
A_phi - maxA_phi_surf)/ maxA_phi.
set(0);
270 for (
int j=0; j<nt; j++) {
271 for (
int l=0; l<Z; l++) {
272 for (
int i=0; i<
mp.get_mg()->get_nr(l); i++) {
273 if (A_phi_scaled2(l,0,j,i) < 0.) {
274 A_phi_scaled2.
set(l,0,j,i) = 0 ;
282 Cmp j_phi_ff = g_j (A_phi_scaled2, bn_j)/maxA_phi.
set(0);
290 + (
ener() +
press())*f_j(A_phi_scaled, an_j)/maxA_phi.
set(0)/g_si
294 j_phi.std_base_scal() ;
298 B_phi = N_j (A_phi_scaled2, bn_j);
299 B_phi.std_base_scal() ;
314 Tenseur source_tAphi(
mp, 1, CON,
mp.get_bvect_spher()) ;
323 source_tAphi.
set(0)=0 ;
324 source_tAphi.
set(1)=0 ;
330 Tenseur WORK_VECT(
mp, 1, CON,
mp.get_bvect_cart()) ;
332 for (
int i=0; i<3; i++) {
333 WORK_VECT.
set(i) = 0 ;
337 WORK_SCAL.
set() = 0 ;
339 double lambda_mag = 0. ;
342 if (source_tAphi.
get_etat() != ETATZERO) {
344 for (
int i=0; i<3; i++) {
345 if(source_tAphi(i).dz_nonzero()) {
346 assert( source_tAphi(i).get_dzpuis() == 4 ) ;
349 (source_tAphi.
set(i)).set_dzpuis(4) ;
354 source_tAphi.
poisson_vect(lambda_mag, par_poisson_Avect, AVECT, WORK_VECT,
357 Cmp A_phi_n(AVECT(2));
366 source_A_1t.
poisson(par_poisson_At, A_1t) ;
368 int L =
mp.get_mg()->get_nt(0);
386 for (
int p=0; p<
mp.get_mg()->get_np(0); p++) {
389 for(
int k=0;k<L;k++){
390 for(
int l=0;l<2*L;l++){
392 if(l==0) leg.
set(k,l)=1. ;
393 if(l==1) leg.
set(k,l)=
cos((*theta)(
l_surf()(p,k),p,k,0)) ;
394 if(l>=2) leg.
set(k,l) = double(2*l-1)/double(l)
396 * leg(k,l-1)-double(l-1)/double(l)*leg(k,l-2) ;
400 for(
int k=0;k<L;k++){
407 for(
int l=0;l<L;l++) MAT.
set(l,k) = leg(k,2*l)/
pow(Rsurf(k),2*l+1);
412 int* IPIV=
new int[L] ;
420 F77_dgesv(&L, &un, MAT.
t, &L, IPIV, VEC.
t, &L, &INFO) ;
424 for(
int k=0;k<L;k++) {VEC2.
set(k)=1. ; }
426 F77_dgesv(&L, &un, MAT_SAVE.
t, &L, IPIV, VEC2.
t, &L, &INFO) ;
430 for(
int nz=0;nz < Z; nz++){
431 for(
int i=0;i<
mp.get_mg()->get_nr(nz);i++){
432 for(
int k=0;k<L;k++){
433 psi.
set(nz,p,k,i) = 0. ;
434 psi2.
set(nz,p,k,i) = 0. ;
435 for(
int l=0;l<L;l++){
436 psi.
set(nz,p,k,i) += VEC(l)*leg(k,2*l) /
437 pow((*
mp.r.c)(nz,p,k,i),2*l+1);
438 psi2.
set(nz,p,k,i) += VEC2(l)*leg(k,2*l)/
439 pow((*
mp.r.c)(nz, p, k,i),2*l+1);
455 Cmp A_t_ext(A_1t + psi) ;
462 for (
int j=0; j<nt; j++)
463 for (
int l=0; l<Z; l++)
464 for (
int i=0; i<
mp.get_mg()->get_nr(l); i++)
465 A_0t.
set(l,0,j,i) = ( (*
mp.r.c)(l,0,j,i) > Rsurf(j) ?
466 A_1t(l,0,j,i) + psi(l,0,j,i) : tmp(l,0,j,i) ) ;
472 double Q_0 = -4*M_PI*(*asymp[1])(Z-1,0,0,0) ;
480 double Q_2 = -4*M_PI*(*asymp[1])(Z-1,0,0,0) ;
488 double C = (
Q-Q_0)/Q_2 ;
498 Cmp A_t_ext(A_0t + C*psi2) ;
505 for (
int j=0; j<nt; j++)
506 for (
int l=0; l<Z; l++)
507 for (
int i=0; i<
mp.get_mg()->get_nr(l); i++)
508 A_t_n.
set(l,0,j,i) = ( (*
mp.r.c)(l,0,j,i) > Rsurf(j) ?
509 A_0t(l,0,j,i) + C*psi2(l,0,j,i) :
510 A_0t(l,0,j,i) + C ) ;
520 A_t = relax_mag*A_t_n + (1.-relax_mag)*
A_t ;
521 A_phi = relax_mag*A_phi_n + (1. - relax_mag)*
A_phi ;
525 cout <<
"not implemented" << endl;
Component of a tensorial field *** DEPRECATED : use class Scalar instead ***.
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...
void div_r()
Division by r everywhere.
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.
void std_base_scal()
Sets the spectral bases of the Valeur va to the standard ones for a scalar.
void annule(int l)
Sets the Cmp to zero in a given domain.
const Cmp & laplacien(int zec_mult_r=4) const
Returns the Laplacian of *this.
int get_dzpuis() const
Returns dzpuis.
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.
Valeur ** asymptot(int n, const int flag=0) const
Asymptotic expansion at r = infinity.
void div_rsint()
Division by .
Cmp poisson() const
Solves the scalar Poisson equation with *this as a source.
Cmp j_phi
-component of the current 4-vector
Cmp A_phi
-component of the electromagnetic potential 1-form divided by .
Cmp A_t
t-component of the elecctromagnetic potential 1-form, divided by .
bool is_conduct() const
Tells if the star is made of conducting or isolating material.
Cmp j_t
t-component of the current 4-vector
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.
Cmp B_phi
-component of the magnetic field
double omega
Rotation angular velocity ([f_unit] ).
Tenseur & logn
Metric potential = logn_auto.
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 .
Tenseur bbb
Metric factor B.
Tenseur b_car
Square of the metric factor B.
Tenseur tnphi
Component of the shift vector.
int nzet
Number of domains of *mp occupied by the star.
Tenseur nnn
Total lapse function.
Map & mp
Mapping associated with the star.
Tenseur ener
Total energy density in the fluid frame.
Tenseur press
Fluid pressure.
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 .
Base class for pure radial mappings.
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...
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
double & set(int i)
Read/write of a particular element (index i) (1D case).
double * t
The array of double.
int get_dim(int i) const
Gives the i-th dimension (ie dim.dim[i]).
Tensor handling *** DEPRECATED : use class Tensor instead ***.
Cmp & set()
Read/write for a scalar (see also operator=(const Cmp&) ).
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
void set_std_base()
Set the standard spectal basis of decomposition for each component.
const Tenseur & gradient_spher() const
Returns the gradient of *this (Spherical coordinates) (scalar field only).
void change_triad(const Base_vect &new_triad)
Sets a new vectorial basis (triad) of decomposition and modifies the components accordingly.
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.
Values and coefficients of a (real-value) function.
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...
const Valeur & ssint() const
Returns of *this.
Tbl max(const Cmp &)
Maximum values of a Cmp in each domain.
Cmp pow(const Cmp &, int)
Power .
Cmp cos(const Cmp &)
Cosine.
Cmp abs(const Cmp &)
Absolute value.
Cmp log(const Cmp &)
Neperian logarithm.
Tenseur flat_scalar_prod_desal(const Tenseur &t1, const Tenseur &t2)
Same as flat_scalar_prod but with desaliasing.
virtual void adapt(const Cmp &ent, const Param &par, int nbr=0)=0
Adaptation of the mapping to a given scalar field.