21char param_elliptic_C[] =
"$Header: /cvsroot/Lorene/C++/Source/Param_elliptic/param_elliptic.C,v 1.21 2014/10/13 08:53:37 j_novak Exp $" ;
110#include "ope_elementary.h"
111#include "param_elliptic.h"
128 cout <<
"Unknown mapping in Param_elliptic" << endl ;
134double Param_elliptic::get_alpha(
int l)
const {
144 cout <<
"Unknown mapping in Param_elliptic" << endl ;
150double Param_elliptic::get_beta(
int l)
const {
154 return mp_af->get_beta()[l] ;
157 return mp_log->get_beta(l) ;
160 cout <<
"Unknown mapping in Param_elliptic" << endl ;
166int Param_elliptic::get_type(
int l)
const {
173 return mp_log->get_type(l) ;
176 cout <<
"Unknown mapping in Param_elliptic" << endl ;
221 if ((map_affine == 0x0) && (map_log == 0x0)) {
222 cout <<
"Param_elliptic not yet defined on this type of mapping" << endl ;
227 if (map_affine != 0x0) {
232 if (map_log != 0x0) {
237 int nz =
get_mp().get_mg()->get_nzone() ;
239 for (
int l=0 ; l<nz ; l++)
241 get_mp().get_mg()->get_nt(l) ;
245 for (
int l=0 ; l<nbr ; l++)
249 int base_r, m_quant, l_quant ;
252 for (
int l=0 ; l<nz ; l++) {
254 nr =
get_mp().get_mg()->get_nr(l) ;
256 for (
int k=0 ; k<
get_mp().get_mg()->get_np(l)+1 ; k++)
257 for (
int j=0 ; j<
get_mp().get_mg()->get_nt(l) ; j++) {
260 (l, k, j, m_quant, l_quant, base_r) ;
266 get_beta(l), l_quant, dzpuis) ;
269 if (
mp_log->get_type(l) == AFFINE)
272 get_beta(l), l_quant, dzpuis) ;
279 cout <<
"Unknown mapping in Param_elliptic::Param_elliptic"
289 var_F.annule_hard() ;
292 var_G.set_etat_qcq() ;
294 var_G.std_spectral_base() ;
313 for (
int l=0 ; l<
get_mp().get_mg()->get_nzone() ; l++)
316 for (
int i=0 ; i<nbr ; i++)
328 for (
int l=0 ; l<
get_mp().get_mg()->get_nzone() ; l++) {
330 np =
get_mp().get_mg()->get_np(l) ;
331 nt =
get_mp().get_mg()->get_nt(l) ;
333 for (
int k=0 ; k<np+1 ; k++)
334 for (
int j=0 ; j<nt ; j++) {
348 int nz =
get_mp().get_mg()->get_nzone() ;
351 int m_quant, base_r_1d, l_quant ;
357 for (
int l=0 ; l<nz ; l++) {
359 nr =
get_mp().get_mg()->get_nr(l) ;
361 for (
int k=0 ; k<
get_mp().get_mg()->get_np(l)+1 ; k++)
362 for (
int j=0 ; j<
get_mp().get_mg()->get_nt(l) ; j++) {
367 (l, k, j, m_quant, l_quant, base_r_1d) ;
369 get_beta(l) , masse) ;
378 if (
mp_log->get_type(zone) != AFFINE) {
379 cout <<
"Operator not define with LOG mapping..." << endl ;
383 for (
int l=0 ; l<nz ; l++) {
385 nr =
get_mp().get_mg()->get_nr(l) ;
387 for (
int k=0 ; k<
get_mp().get_mg()->get_np(l)+1 ; k++)
388 for (
int j=0 ; j<
get_mp().get_mg()->get_nt(l) ; j++) {
393 (l, k, j, m_quant, l_quant, base_r_1d) ;
395 get_alpha(l), get_beta(l), masse) ;
405 cout <<
"Unkown mapping in set_helmhotz_minus" << endl ;
416 int nz =
get_mp().get_mg()->get_nzone() ;
419 int m_quant, base_r_1d, l_quant ;
425 for (
int l=0 ; l<nz ; l++) {
427 nr =
get_mp().get_mg()->get_nr(l) ;
429 for (
int k=0 ; k<
get_mp().get_mg()->get_np(l)+1 ; k++)
430 for (
int j=0 ; j<
get_mp().get_mg()->get_nt(l) ; j++) {
433 int old_base =
operateurs[conte]->get_base_r() ;
438 (l, k, j, m_quant, l_quant, base_r_1d) ;
440 get_beta(l) , masse) ;
450 if (
mp_log->get_type(zone) != AFFINE) {
451 cout <<
"Operator not define with LOG mapping..." << endl ;
455 for (
int l=0 ; l<nz ; l++) {
457 nr =
get_mp().get_mg()->get_nr(l) ;
459 for (
int k=0 ; k<
get_mp().get_mg()->get_np(l)+1 ; k++)
460 for (
int j=0 ; j<
get_mp().get_mg()->get_nt(l) ; j++) {
463 int old_base =
operateurs[conte]->get_base_r() ;
468 (l, k, j, m_quant, l_quant, base_r_1d) ;
470 get_alpha(l), get_beta(l), masse) ;
481 cout <<
"Unkown mapping in set_helmhotz_minus" << endl ;
490 cout <<
"set_poisson_vect_r only defined for an affine mapping..." << endl ;
495 int nz =
get_mp().get_mg()->get_nzone() ;
500 for (
int l=0 ; l<nz ; l++) {
502 nr =
get_mp().get_mg()->get_nr(l) ;
504 for (
int k=0 ; k<
get_mp().get_mg()->get_np(l)+1 ; k++)
505 for (
int j=0 ; j<
get_mp().get_mg()->get_nt(l) ; j++) {
506 if ((
operateurs[conte] != 0x0) && (l==zone)) {
507 int old_base =
operateurs[conte]->get_base_r() ;
510 assert (op_pois !=0x0) ;
516 if ((!only_l_zero)||(lq_old == 0)) {
518 get_beta(l), lq_old, dzp) ;
535 int nz =
get_mp().get_mg()->get_nzone() ;
538 for (
int l=0 ; l<nz ; l++) {
540 bool ced = (
get_mp().get_mg()->get_type_r(l) == UNSURR ) ;
541 for (
int k=0 ; k<
get_mp().get_mg()->get_np(l)+1 ; k++)
542 for (
int j=0 ; j<
get_mp().get_mg()->get_nt(l) ; j++) {
543 if ((
operateurs[conte] != 0x0) && (l==zone)) {
546 assert (op_pois !=0x0) ;
563 cout <<
"set_poisson_tens_rr only defined for an affine mapping..."
569 int nz =
get_mp().get_mg()->get_nzone() ;
574 for (
int l=0 ; l<nz ; l++) {
576 nr =
get_mp().get_mg()->get_nr(l) ;
578 for (
int k=0 ; k<
get_mp().get_mg()->get_np(l)+1 ; k++)
579 for (
int j=0 ; j<
get_mp().get_mg()->get_nt(l) ; j++) {
580 if ((
operateurs[conte] != 0x0) && (l==zone)) {
581 int old_base =
operateurs[conte]->get_base_r() ;
584 assert (op_pois !=0x0) ;
591 (nr, old_base,get_alpha(l), get_beta(l), lq_old, dzp) ;
605 int nz =
get_mp().get_mg()->get_nzone() ;
613 for (
int l=0 ; l<nz ; l++) {
615 nr =
get_mp().get_mg()->get_nr(l) ;
617 for (
int k=0 ; k<
get_mp().get_mg()->get_np(l)+1 ; k++)
618 for (
int j=0 ; j<
get_mp().get_mg()->get_nt(l) ; j++) {
619 if ((
operateurs[conte] != 0x0) && (l==zone)) {
620 int old_base =
operateurs[conte]->get_base_r() ;
626 get_beta(l), a, b, c) ;
635 if (
mp_log->get_type(zone) != AFFINE) {
636 cout <<
"Operator not define with LOG mapping..." << endl ;
640 for (
int l=0 ; l<nz ; l++) {
642 nr =
get_mp().get_mg()->get_nr(l) ;
644 for (
int k=0 ; k<
get_mp().get_mg()->get_np(l)+1 ; k++)
645 for (
int j=0 ; j<
get_mp().get_mg()->get_nt(l) ; j++) {
646 if ((
operateurs[conte] != 0x0) && (l==zone)) {
647 int old_base =
operateurs[conte]->get_base_r() ;
653 get_beta(l), a, b, c) ;
663 cout <<
"Unkown mapping in set_sec_order_r2" << endl ;
671 if ((
type_map == MAP_AFF) || (
mp_log->get_type(zone) == AFFINE)) {
672 cout <<
"set_sec_order only defined for a log mapping" << endl ;
677 int nz =
get_mp().get_mg()->get_nzone() ;
682 for (
int l=0 ; l<nz ; l++) {
684 nr =
get_mp().get_mg()->get_nr(l) ;
686 for (
int k=0 ; k<
get_mp().get_mg()->get_np(l)+1 ; k++)
687 for (
int j=0 ; j<
get_mp().get_mg()->get_nt(l) ; j++) {
688 if ((
operateurs[conte] != 0x0) && (l==zone)) {
690 int old_base =
operateurs[conte]->get_base_r() ;
696 get_beta(l), a, b, c) ;
710 int nz =
get_mp().get_mg()->get_nzone() ;
713 int m_quant, base_r_1d, l_quant ;
720 for (
int l=0 ; l<nz ; l++) {
721 int dz = (l==nz-1) ? dzpuis : 0 ;
722 nr =
get_mp().get_mg()->get_nr(l) ;
724 for (
int k=0 ; k<
get_mp().get_mg()->get_np(l)+1 ; k++)
725 for (
int j=0 ; j<
get_mp().get_mg()->get_nt(l) ; j++) {
730 (l, k, j, m_quant, l_quant, base_r_1d) ;
732 get_beta(l), l_quant, dz) ;
741 if (
mp_log->get_type(zone) != AFFINE) {
742 cout <<
"Operator not define with LOG mapping..." << endl ;
746 for (
int l=0 ; l<nz ; l++) {
747 int dz = (l==nz-1) ? dzpuis : 0 ;
748 nr =
get_mp().get_mg()->get_nr(l) ;
750 for (
int k=0 ; k<
get_mp().get_mg()->get_np(l)+1 ; k++)
751 for (
int j=0 ; j<
get_mp().get_mg()->get_nt(l) ; j++) {
756 (l, k, j, m_quant, l_quant, base_r_1d) ;
758 get_alpha(l), get_beta(l), l_quant, dz) ;
768 cout <<
"Unkown mapping in set_ope_vorton" << endl ;
776 assert (so.
get_etat() != ETATNONDEF) ;
785 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.