116#include "ope_elementary.h"
117#include "param_elliptic.h"
134 cout <<
"Unknown mapping in Param_elliptic" << endl ;
140double Param_elliptic::get_alpha(
int l)
const {
150 cout <<
"Unknown mapping in Param_elliptic" << endl ;
156double Param_elliptic::get_beta(
int l)
const {
160 return mp_af->get_beta()[l] ;
163 return mp_log->get_beta(l) ;
166 cout <<
"Unknown mapping in Param_elliptic" << endl ;
172int Param_elliptic::get_type(
int l)
const {
179 return mp_log->get_type(l) ;
182 cout <<
"Unknown mapping in Param_elliptic" << endl ;
227 if ((map_affine == 0x0) && (map_log == 0x0)) {
228 cout <<
"Param_elliptic not yet defined on this type of mapping" << endl ;
233 if (map_affine != 0x0) {
238 if (map_log != 0x0) {
243 int nz =
get_mp().get_mg()->get_nzone() ;
245 for (
int l=0 ; l<nz ; l++)
247 get_mp().get_mg()->get_nt(l) ;
251 for (
int l=0 ; l<nbr ; l++)
255 int base_r, m_quant, l_quant ;
258 for (
int l=0 ; l<nz ; l++) {
260 nr =
get_mp().get_mg()->get_nr(l) ;
262 for (
int k=0 ; k<
get_mp().get_mg()->get_np(l)+1 ; k++)
263 for (
int j=0 ; j<
get_mp().get_mg()->get_nt(l) ; j++) {
266 (l, k, j, m_quant, l_quant, base_r) ;
272 get_beta(l), l_quant, dzpuis) ;
275 if (
mp_log->get_type(l) == AFFINE)
278 get_beta(l), l_quant, dzpuis) ;
285 cout <<
"Unknown mapping in Param_elliptic::Param_elliptic"
295 var_F.annule_hard() ;
298 var_G.set_etat_qcq() ;
300 var_G.std_spectral_base() ;
319 for (
int l=0 ; l<
get_mp().get_mg()->get_nzone() ; l++)
322 for (
int i=0 ; i<nbr ; i++)
334 for (
int l=0 ; l<
get_mp().get_mg()->get_nzone() ; l++) {
336 np =
get_mp().get_mg()->get_np(l) ;
337 nt =
get_mp().get_mg()->get_nt(l) ;
339 for (
int k=0 ; k<np+1 ; k++)
340 for (
int j=0 ; j<nt ; j++) {
354 int nz =
get_mp().get_mg()->get_nzone() ;
357 int m_quant, base_r_1d, l_quant ;
363 for (
int l=0 ; l<nz ; l++) {
365 nr =
get_mp().get_mg()->get_nr(l) ;
367 for (
int k=0 ; k<
get_mp().get_mg()->get_np(l)+1 ; k++)
368 for (
int j=0 ; j<
get_mp().get_mg()->get_nt(l) ; j++) {
373 (l, k, j, m_quant, l_quant, base_r_1d) ;
375 get_beta(l) , masse) ;
384 if (
mp_log->get_type(zone) != AFFINE) {
385 cout <<
"Operator not define with LOG mapping..." << endl ;
389 for (
int l=0 ; l<nz ; l++) {
391 nr =
get_mp().get_mg()->get_nr(l) ;
393 for (
int k=0 ; k<
get_mp().get_mg()->get_np(l)+1 ; k++)
394 for (
int j=0 ; j<
get_mp().get_mg()->get_nt(l) ; j++) {
399 (l, k, j, m_quant, l_quant, base_r_1d) ;
401 get_alpha(l), get_beta(l), masse) ;
411 cout <<
"Unkown mapping in set_helmhotz_minus" << endl ;
422 int nz =
get_mp().get_mg()->get_nzone() ;
425 int m_quant, base_r_1d, l_quant ;
431 for (
int l=0 ; l<nz ; l++) {
433 nr =
get_mp().get_mg()->get_nr(l) ;
435 for (
int k=0 ; k<
get_mp().get_mg()->get_np(l)+1 ; k++)
436 for (
int j=0 ; j<
get_mp().get_mg()->get_nt(l) ; j++) {
439 int old_base =
operateurs[conte]->get_base_r() ;
444 (l, k, j, m_quant, l_quant, base_r_1d) ;
446 get_beta(l) , masse) ;
456 if (
mp_log->get_type(zone) != AFFINE) {
457 cout <<
"Operator not define with LOG mapping..." << endl ;
461 for (
int l=0 ; l<nz ; l++) {
463 nr =
get_mp().get_mg()->get_nr(l) ;
465 for (
int k=0 ; k<
get_mp().get_mg()->get_np(l)+1 ; k++)
466 for (
int j=0 ; j<
get_mp().get_mg()->get_nt(l) ; j++) {
469 int old_base =
operateurs[conte]->get_base_r() ;
474 (l, k, j, m_quant, l_quant, base_r_1d) ;
476 get_alpha(l), get_beta(l), masse) ;
487 cout <<
"Unkown mapping in set_helmhotz_minus" << endl ;
496 cout <<
"set_poisson_vect_r only defined for an affine mapping..." << endl ;
501 int nz =
get_mp().get_mg()->get_nzone() ;
506 for (
int l=0 ; l<nz ; l++) {
508 nr =
get_mp().get_mg()->get_nr(l) ;
510 for (
int k=0 ; k<
get_mp().get_mg()->get_np(l)+1 ; k++)
511 for (
int j=0 ; j<
get_mp().get_mg()->get_nt(l) ; j++) {
512 if ((
operateurs[conte] != 0x0) && (l==zone)) {
513 int old_base =
operateurs[conte]->get_base_r() ;
516 assert (op_pois !=0x0) ;
522 if ((!only_l_zero)||(lq_old == 0)) {
524 get_beta(l), lq_old, dzp) ;
541 int nz =
get_mp().get_mg()->get_nzone() ;
544 for (
int l=0 ; l<nz ; l++) {
546 bool ced = (
get_mp().get_mg()->get_type_r(l) == UNSURR ) ;
547 for (
int k=0 ; k<
get_mp().get_mg()->get_np(l)+1 ; k++)
548 for (
int j=0 ; j<
get_mp().get_mg()->get_nt(l) ; j++) {
549 if ((
operateurs[conte] != 0x0) && (l==zone)) {
552 assert (op_pois !=0x0) ;
569 cout <<
"set_poisson_tens_rr only defined for an affine mapping..."
575 int nz =
get_mp().get_mg()->get_nzone() ;
580 for (
int l=0 ; l<nz ; l++) {
582 nr =
get_mp().get_mg()->get_nr(l) ;
584 for (
int k=0 ; k<
get_mp().get_mg()->get_np(l)+1 ; k++)
585 for (
int j=0 ; j<
get_mp().get_mg()->get_nt(l) ; j++) {
586 if ((
operateurs[conte] != 0x0) && (l==zone)) {
587 int old_base =
operateurs[conte]->get_base_r() ;
590 assert (op_pois !=0x0) ;
597 (nr, old_base,get_alpha(l), get_beta(l), lq_old, dzp) ;
611 int nz =
get_mp().get_mg()->get_nzone() ;
619 for (
int l=0 ; l<nz ; l++) {
621 nr =
get_mp().get_mg()->get_nr(l) ;
623 for (
int k=0 ; k<
get_mp().get_mg()->get_np(l)+1 ; k++)
624 for (
int j=0 ; j<
get_mp().get_mg()->get_nt(l) ; j++) {
625 if ((
operateurs[conte] != 0x0) && (l==zone)) {
626 int old_base =
operateurs[conte]->get_base_r() ;
632 get_beta(l), a, b, c) ;
641 if (
mp_log->get_type(zone) != AFFINE) {
642 cout <<
"Operator not define with LOG mapping..." << endl ;
646 for (
int l=0 ; l<nz ; l++) {
648 nr =
get_mp().get_mg()->get_nr(l) ;
650 for (
int k=0 ; k<
get_mp().get_mg()->get_np(l)+1 ; k++)
651 for (
int j=0 ; j<
get_mp().get_mg()->get_nt(l) ; j++) {
652 if ((
operateurs[conte] != 0x0) && (l==zone)) {
653 int old_base =
operateurs[conte]->get_base_r() ;
659 get_beta(l), a, b, c) ;
669 cout <<
"Unkown mapping in set_sec_order_r2" << endl ;
677 if ((
type_map == MAP_AFF) || (
mp_log->get_type(zone) == AFFINE)) {
678 cout <<
"set_sec_order only defined for a log mapping" << endl ;
683 int nz =
get_mp().get_mg()->get_nzone() ;
688 for (
int l=0 ; l<nz ; l++) {
690 nr =
get_mp().get_mg()->get_nr(l) ;
692 for (
int k=0 ; k<
get_mp().get_mg()->get_np(l)+1 ; k++)
693 for (
int j=0 ; j<
get_mp().get_mg()->get_nt(l) ; j++) {
694 if ((
operateurs[conte] != 0x0) && (l==zone)) {
696 int old_base =
operateurs[conte]->get_base_r() ;
702 get_beta(l), a, b, c) ;
716 int nz =
get_mp().get_mg()->get_nzone() ;
719 int m_quant, base_r_1d, l_quant ;
726 for (
int l=0 ; l<nz ; l++) {
727 int dz = (l==nz-1) ? dzpuis : 0 ;
728 nr =
get_mp().get_mg()->get_nr(l) ;
730 for (
int k=0 ; k<
get_mp().get_mg()->get_np(l)+1 ; k++)
731 for (
int j=0 ; j<
get_mp().get_mg()->get_nt(l) ; j++) {
736 (l, k, j, m_quant, l_quant, base_r_1d) ;
738 get_beta(l), l_quant, dz) ;
747 if (
mp_log->get_type(zone) != AFFINE) {
748 cout <<
"Operator not define with LOG mapping..." << endl ;
752 for (
int l=0 ; l<nz ; l++) {
753 int dz = (l==nz-1) ? dzpuis : 0 ;
754 nr =
get_mp().get_mg()->get_nr(l) ;
756 for (
int k=0 ; k<
get_mp().get_mg()->get_np(l)+1 ; k++)
757 for (
int j=0 ; j<
get_mp().get_mg()->get_nt(l) ; j++) {
762 (l, k, j, m_quant, l_quant, base_r_1d) ;
764 get_alpha(l), get_beta(l), l_quant, dz) ;
774 cout <<
"Unkown mapping in set_ope_vorton" << endl ;
782 assert (so.
get_etat() != ETATNONDEF) ;
791 assert (so.
get_etat() != ETATNONDEF) ;
Bases of the spectral expansions.
void give_quant_numbers(int, int, int, int &, int &, int &) const
Computes the various quantum numbers and 1d radial base.
const double * get_alpha() const
Returns the pointer on the array alpha.
Logarithmic radial mapping.
double get_alpha(int l) const
Returns in the domain l.
Base class for pure radial mappings.
Basic class for elementary elliptic operators.
Class for the Helmholtz operator ( ).
Class for the Helmholtz operator (m > 0).
Class for the operator of the rr component of the divergence-free tensor Poisson equation.
Class for the operator of the r component of the vector Poisson equation.
Class for the operator of the Poisson equation (i.e.
int get_dzpuis()
Returns the associated dzpuis, if in the compactified domain.
int get_lquant()
Returns the quantum number l.
virtual void dec_l_quant()
Decreases the quatum number l by one unit.
virtual void inc_l_quant()
Increases the quatum number l by one unit.
Class for operator of the type .
Class for operator of the type .
Class for the operator appearing for the vortons.
Scalar var_G
Multiplicative variable change that must be sphericaly symetric !
Itbl done_F
Stores what has been computed for F.
Tbl val_dG_minus
Values of the derivative of G at the inner boundaries of the various domains.
void set_sec_order_r2(int zone, double a, double b, double c)
Set the operator to in one domain (only in the shells).
void inc_l_quant(int zone)
Increases the quantum number l in the domain zone .
void set_helmholtz_minus(int zone, double mas, Scalar &so)
Set the operator to in one domain (not in the nucleus).
Tbl val_dF_minus
Values of the derivative of F at the inner boundaries of the various domains.
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...
const Map_log * mp_log
The mapping if log type.
Param_elliptic(const Scalar &)
Standard constructor from a Scalar.
const Map_radial & get_mp() const
Returns the mapping.
Tbl val_F_minus
Values of F at the inner boundaries of the various domains.
Tbl val_F_plus
Values of F at the outer boundaries of the various domains.
~Param_elliptic()
Destructor.
void set_sec_order(int zone, double a, double b, double c)
Set the operator to in one domain (only in the shells).
Tbl val_G_minus
Values of G at the inner boundaries of the various domains.
Tbl val_dG_plus
Values of the derivative of G at the outer boundaries of the various domains.
Ope_elementary ** operateurs
Array on the elementary operators.
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.
const Map_af * mp_af
The mapping, if affine.
void set_poisson_tens_rr(int zone)
Sets the operator to in all domains, for only.
void set_variable_F(const Scalar &)
Changes the variable function F.
void set_helmholtz_plus(int zone, double mas, Scalar &so)
Set the operator to in one domain (only in the shells).
void set_variable_G(const Scalar &)
Changes the variable function G.
int type_map
Type of mapping either MAP_AFF or MAP_LOG.
Tbl val_dF_plus
Values of the derivative of F at the outer boundaries of the various domains.
Tbl val_G_plus
Values of G at the outer boundaries of the various domains.
void set_ope_vorton(int zone, Scalar &so)
Set the operator to in one domain (not implemented in the nucleus).
Itbl done_G
Stores what has been computed for G.
Scalar var_F
Additive variable change function.
Tensor field of valence 0 (or component of a tensorial field).
int get_dzpuis() const
Returns dzpuis.
bool dz_nonzero() const
Returns true if the last domain is compactified and *this is not zero in this domain.
Valeur & set_spectral_va()
Returns va (read/write version).
const Valeur & get_spectral_va() const
Returns va (read only version).
int get_etat() const
Returns the logical state ETATNONDEF (undefined), ETATZERO (null) or ETATQCQ (ordinary).
const Base_val & get_base() const
Return the bases for spectral expansions (member base ).
void ylm()
Computes the coefficients of *this.
void coef() const
Computes the coeffcients of *this.
Base_val base
Bases on which the spectral expansion is performed.
#define R_CHEBI
base de Cheb. impaire (rare) seulement
const Map & get_mp() const
Returns the mapping.
const Mg3d * get_mg() const
Gives the Mg3d on which the mapping is defined.