129#include "param_elliptic.h"
136 bool nullite = true ;
137 for (
int i=0; i<3; i++) {
138 assert(
cmp[i]->check_dzpuis(4)) ;
139 if (
cmp[i]->get_etat() != ETATZERO) nullite = false ;
141 assert ((method>=0) && (method<7)) ;
153 if (fabs(lambda+1) < 1.e-8)
169 if (fabs(lambda+1) < 1.e-8)
172 divf = (
potential(met_f) / (lambda + 1)) ;
175 int nz =
mp->get_mg()->get_nzone() ;
180 Scalar div_lnot0 = divf ;
182 Scalar source_r(*
mp) ;
184 if (div_lnot0.
get_etat() != ETATZERO) {
187 for (
int lz=0; lz<nz; lz++) {
188 int np =
mp->get_mg()->get_np(lz) ;
189 int nt =
mp->get_mg()->get_nt(lz) ;
190 int nr =
mp->get_mg()->get_nr(lz) ;
192 for (
int k=0; k<np+1; k++)
193 for (
int j=0; j<nt; j++) {
194 int l_q, m_q, base_r ;
195 if (nullite_plm(j, nt, k, np, va_div.
base) == 1) {
196 donne_lm(nz, lz, j, k, va_div.
base, m_q, l_q, base_r) ;
198 for (
int i=0; i<nr; i++)
199 va_div.
c_cf->
set(lz, k, j, i) = 0. ;
214 source_r += *(
cmp[0]) - lambda*divf.
dsdr() ;
216 if (source_r.
get_etat() != ETATZERO) {
222 for (
int lz=0; lz<nz; lz++)
231 Scalar source_eta = divf - f_r.
dsdr() ;
233 source_eta -= 2*f_r ;
242 Tenseur source_p(*
mp, 1, CON,
mp->get_bvect_spher() ) ;
244 for (
int i=0; i<3; i++) {
248 Tenseur vect_auxi (*
mp, 1, CON,
mp->get_bvect_cart()) ;
256 for (
int i=1; i<=3; i++)
257 resu.
set(i) = resu_p(i-1) ;
264 Tenseur source_p(*
mp, 1, CON,
mp->get_bvect_spher() ) ;
266 for (
int i=0; i<3; i++) {
276 for (
int i=1; i<=3; i++)
277 resu.
set(i) = resu_p(i-1) ;
284 if (fabs(lambda+1) < 1.e-8)
287 divf = (
potential(met_f) / (lambda + 1)) ;
290 int nz =
mp->get_mg()->get_nzone() ;
295 Scalar div_tmp = divf ;
297 Scalar source_r(*
mp) ;
299 if (div_tmp.
get_etat() != ETATZERO) {
301 for (
int lz=0; lz<nz; lz++) {
302 int np =
mp->get_mg()->get_np(lz) ;
303 int nt =
mp->get_mg()->get_nt(lz) ;
304 int nr =
mp->get_mg()->get_nr(lz) ;
306 for (
int k=0; k<np+1; k++)
307 for (
int j=0; j<nt; j++) {
308 int l_q, m_q, base_r ;
309 if (nullite_plm(j, nt, k, np, va_div.
base) == 1) {
310 donne_lm(nz, lz, j, k, va_div.
base, m_q, l_q, base_r) ;
312 for (
int i=0; i<nr; i++)
313 va_div.
c_cf->
set(lz, k, j, i) = 0. ;
328 source_r += *(
cmp[0]) - lambda*divf.
dsdr() ;
330 if (source_r.
get_etat() != ETATZERO) {
336 for (
int lz=0; lz<nz; lz++)
347 Scalar source_eta = *
cmp[1] ;
349 source_eta = source_eta.
dsdt() ;
353 Scalar dfrsr = f_r.
dsdr() ;
369 for (
int lz=0; lz<nz; lz++) {
370 bool ced = (
mp->get_mg()->get_type_r(lz) == UNSURR) ;
371 int np =
mp->get_mg()->get_np(lz) ;
372 int nt =
mp->get_mg()->get_nt(lz) ;
373 int nr =
mp->get_mg()->get_nr(lz) ;
375 for (
int k=0; k<np+1; k++)
376 for (
int j=0; j<nt; j++) {
377 int l_q, m_q, base_r ;
378 if (nullite_plm(j, nt, k, np, va_eta.
base) == 1) {
379 donne_lm(nz, lz, j, k, va_eta.
base, m_q, l_q, base_r) ;
381 for (
int i=0; i<nr; i++) {
384 -= (lambda + 2. / double(ced ? -l_q : (l_q+1) ))
385 * va_div.
c_cf->operator()(lz, k, j, i) ;
386 if (va_fsr.
c_cf->operator()(lz).
get_etat() != ETATZERO) {
388 += 2. / double(ced ? -l_q : (l_q+1) )
389 * va_dfr.
c_cf->operator()(lz, k, j, i) ;
391 -= 2.*( ced ? double(l_q+2)/double(l_q)
392 : double(l_q-1)/double(l_q+1) )
393 * va_fsr.
c_cf->
set(lz, k, j, i) ;
399 if (va_eta.
c != 0x0) {
409 for (
int lz=0; lz<nz; lz++)
420 if (fabs(lambda+1) < 1.e-8)
423 divf = (
potential(met_f) / (lambda + 1)) ;
426 int nz =
mp->get_mg()->get_nzone() ;
431 Scalar div_lnot0 = divf ;
433 Scalar source_r(*
mp) ;
435 if (div_lnot0.
get_etat() != ETATZERO) {
438 for (
int lz=0; lz<nz; lz++) {
439 int np =
mp->get_mg()->get_np(lz) ;
440 int nt =
mp->get_mg()->get_nt(lz) ;
441 int nr =
mp->get_mg()->get_nr(lz) ;
443 for (
int k=0; k<np+1; k++)
444 for (
int j=0; j<nt; j++) {
445 int l_q, m_q, base_r ;
446 if (nullite_plm(j, nt, k, np, va_div.
base) == 1) {
447 donne_lm(nz, lz, j, k, va_div.
base, m_q, l_q, base_r) ;
449 for (
int i=0; i<nr; i++)
450 va_div.
c_cf->
set(lz, k, j, i) = 0. ;
465 source_r += *(
cmp[0]) - lambda*divf.
dsdr() ;
467 if (source_r.
get_etat() != ETATZERO) {
474 for (
int lz=0; lz<nz; lz++)
482 Scalar source_eta = - *(
cmp[0]) ;
484 source_eta -= (lambda+2.)*divf ;
487 Scalar tmp = 2*f_r + f_r.
lapang() ;
490 tmp = (1.+lambda)*divf ;
493 source_eta = source_eta.
primr() ;
502 poisson_block(lambda, resu) ;
507 cout <<
"Vector::poisson : unexpected type of method !" << endl
508 <<
" method = " << method << endl ;
526 assert ((tspher != 0x0) || (tcart != 0x0)) ;
530 assert (tcart == 0x0) ;
531 met_f = &(
mp->flat_met_spher()) ;
535 assert (tspher == 0x0) ;
536 met_f = &(
mp->flat_met_cart()) ;
539 return (
poisson(lambda, *met_f, method) );
549 for (
int i=0; i<3; i++)
550 assert(
cmp[i]->check_dzpuis(4)) ;
552 assert ((method==0) || (method==2)) ;
561 int nitermax = par.
get_int(0) ;
577 if (fabs(lambda+1) < 1.e-8)
580 Scalar tmp =
potential(met_local) / (lambda + 1) ;
589 par_free.
add_int(nitermax, 0) ;
604 Tenseur source_p(*
mp, 1, CON,
mp->get_bvect_spher() ) ;
606 for (
int i=0; i<3; i++) {
611 Tenseur vect_auxi (*
mp, 1, CON,
mp->get_bvect_cart()) ;
613 for (
int i=0; i<3 ;i++){
614 vect_auxi.
set(i) = 0. ;
621 Tenseur resu_p(*
mp, 1, CON,
mp->get_bvect_cart() ) ;
623 source_p.
poisson_vect(lambda, par, resu_p, vect_auxi, scal_auxi) ;
626 for (
int i=1; i<=3; i++)
627 resu.
set(i) = resu_p(i-1) ;
633 cout <<
"Vector::poisson : unexpected type of method !" << endl
634 <<
" method = " << method << endl ;
Cartesian vectorial bases (triads).
Spherical orthonormal vectorial bases (triads).
Component of a tensorial field *** DEPRECATED : use class Scalar instead ***.
void annule_hard()
Sets the Cmp to zero in a hard way.
Flat metric for tensor calculation.
int get_etat() const
Returns the logical state.
Tbl & set(int l)
Read/write of the Tbl containing the coefficients in a given domain.
This class contains the parameters needed to call the general elliptic solver.
void set_poisson_vect_eta(int zone)
Sets the operator to be a regular elliptic operator to solve for the component of the vector Poisson...
void set_poisson_vect_r(int zone, bool only_l_zero=false)
Sets the operator to in all domains, for ; and to in all domains otherwise.
void add_double(const double &x, int position=0)
Adds the the address of a new double to the list.
Cmp & get_cmp_mod(int position=0) const
Returns the reference of a modifiable Cmp stored in the list.
const int & get_int(int position=0) const
Returns the reference of a int stored in the list.
const double & get_double(int position=0) const
Returns the reference of a double stored in the list.
void add_cmp_mod(Cmp &ti, int position=0)
Adds the address of a new modifiable Cmp to the list.
void add_int_mod(int &n, int position=0)
Adds the address of a new modifiable int to the list.
int & get_int_mod(int position=0) const
Returns the reference of a modifiable int stored in the list.
void add_int(const int &n, int position=0)
Adds the address of a new int to the list.
const Scalar & lapang() const
Returns the angular Laplacian of *this , where .
void div_r_dzpuis(int ced_mult_r)
Division by r everywhere but with the output flag dzpuis set to ced_mult_r .
const Scalar & dsdt() const
Returns of *this .
void div_sint()
Division by .
Scalar poisson() const
Solves the scalar Poisson equation with *this as a source.
virtual void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Scalar poisson_angu(double lambda=0) const
Solves the (generalized) angular Poisson equation with *this as source.
virtual void inc_dzpuis(int inc=1)
Increases by inc units the value of dzpuis and changes accordingly the values of the Scalar in the co...
virtual void set_etat_zero()
Sets the logical state to ETATZERO (zero).
const Scalar & dsdr() const
Returns of *this .
Valeur & set_spectral_va()
Returns va (read/write version).
void annule_hard()
Sets the Scalar to zero in a hard way.
void mult_sint()
Multiplication by .
int get_etat() const
Returns the logical state ETATNONDEF (undefined), ETATZERO (null) or ETATQCQ (ordinary).
void set_dzpuis(int)
Modifies the dzpuis flag.
Scalar sol_elliptic(Param_elliptic ¶ms) const
Resolution of a general elliptic equation, putting zero at infinity.
void mult_r_dzpuis(int ced_mult_r)
Multiplication by r everywhere but with the output flag dzpuis set to ced_mult_r .
const Vector & derive_con(const Metric &gam) const
Returns the "contravariant" derivative of *this with respect to some metric , by raising the index of...
virtual void dec_dzpuis(int dec=1)
Decreases by dec units the value of dzpuis and changes accordingly the values of the Scalar in the co...
Scalar primr(bool null_infty=true) const
Computes the radial primitive which vanishes for .
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 poisson_vect_oohara(double lambda, Param &par, Tenseur &shift, Tenseur &scal) const
Solves the vectorial Poisson equation .
void set_std_base()
Set the standard spectal basis of decomposition for each component.
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 : .
Values and coefficients of a (real-value) function.
void ylm()
Computes the coefficients of *this.
Mtbl * c
Values of the function at the points of the multi-grid.
Mtbl_cf * c_cf
Coefficients of the spectral expansion of the function.
void coef() const
Computes the coeffcients of *this.
void ylm_i()
Inverse of ylm().
Base_val base
Bases on which the spectral expansion is performed.
Vector_divfree poisson() const
Computes the solution of a vectorial Poisson equation with *this as a source:
Vector poisson(double lambda, int method=6) const
Solves the vector Poisson equation with *this as a source.
virtual const Scalar & mu() const
Gives the field such that the angular components of the vector are written:
Vector(const Map &map, int tipe, const Base_vect &triad_i)
Standard constructor.
const Vector_divfree & div_free(const Metric &) const
Returns the div-free vector in the Helmholtz decomposition.
const Scalar & potential(const Metric &) const
Returns the potential in the Helmholtz decomposition.
Scalar & set(int)
Read/write access to a component.
void set_vr_eta_mu(const Scalar &vr_i, const Scalar &eta_i, const Scalar &mu_i)
Defines the components through potentials and (see members p_eta and p_mu ), as well as the compon...
virtual void dec_dzpuis(int dec=1)
Decreases by dec units the value of dzpuis and changes accordingly the values in the compactified ext...
const Map *const mp
Mapping on which the numerical values at the grid points are defined.
Scalar ** cmp
Array of size n_comp of pointers onto the components.
virtual void set_etat_zero()
Sets the logical state of all components to ETATZERO (zero state).
const Base_vect * triad
Vectorial basis (triad) with respect to which the tensor components are defined.
virtual void poisson_angu(const Scalar &source, Param &par, Scalar &uu, double lambda=0) const =0
Computes the solution of the generalized angular Poisson equation.