00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026 #ifndef __OPE_ELEMENTARY_H_
00027 #define __OPE_ELEMENTARY_H_
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047
00048
00049
00050
00051
00052
00053
00054
00055
00056
00057
00058
00059
00060
00061
00062
00063
00064
00065
00066
00067
00068
00069
00070
00071
00072
00073
00074
00075
00076
00077
00078
00079 #include "matrice.h"
00080
00098 class Ope_elementary {
00099
00100 protected:
00101
00102 int nr ;
00103 int base_r ;
00104 double alpha ;
00105 double beta ;
00106
00110 mutable Matrice* ope_mat ;
00114 mutable Matrice* ope_cl ;
00118 mutable Matrice* non_dege ;
00119
00123 mutable double s_one_plus ;
00127 mutable double s_one_minus ;
00132 mutable double ds_one_plus ;
00137 mutable double ds_one_minus ;
00138
00142 mutable double s_two_plus ;
00146 mutable double s_two_minus ;
00151 mutable double ds_two_plus ;
00156 mutable double ds_two_minus ;
00157
00161 mutable double sp_minus ;
00165 mutable double sp_plus ;
00169 mutable double dsp_minus ;
00173 mutable double dsp_plus ;
00174
00175
00176 protected:
00185 explicit Ope_elementary (int nbr , int baser , double alf, double eta) ;
00186 Ope_elementary (const Ope_elementary&) ;
00187
00188 public:
00189 virtual ~Ope_elementary() ;
00190
00191 public:
00196 double val_sh_one_minus() const {return s_one_minus ;} ;
00201 double val_sh_one_plus() const {return s_one_plus ;} ;
00207 double der_sh_one_minus() const {return ds_one_minus ;} ;
00213 double der_sh_one_plus() const {return ds_one_plus ;} ;
00214
00219 double val_sh_two_minus() const {return s_two_minus ;} ;
00224 double val_sh_two_plus() const {return s_two_plus ;} ;
00230 double der_sh_two_minus() const {return ds_two_minus ;} ;
00236 double der_sh_two_plus() const {return ds_two_plus ;} ;
00237
00241 double val_sp_minus() const {return sp_minus ;} ;
00245 double val_sp_plus() const {return sp_plus ;} ;
00250 double der_sp_minus() const {return dsp_minus ;} ;
00255 double der_sp_plus() const {return dsp_plus ;} ;
00256
00258 double get_alpha() const {return alpha ;} ;
00259
00261 double get_beta() const {return beta ;} ;
00262
00264 int get_base_r() const {return base_r ;} ;
00265
00267 Matrice get_ope_mat() {
00268 if (ope_mat ==0x0)
00269 do_ope_mat() ;
00270 return *ope_mat ;
00271 }
00272
00274 Matrice get_ope_cl() {
00275 if (ope_cl ==0x0)
00276 do_ope_cl() ;
00277 return *ope_cl ;
00278 }
00279
00281 Matrice get_non_dege() {
00282 if (non_dege ==0x0)
00283 do_non_dege() ;
00284 return *non_dege ;
00285 }
00286 private:
00290 virtual void do_ope_mat() const = 0 ;
00294 virtual void do_ope_cl() const = 0 ;
00298 virtual void do_non_dege() const = 0 ;
00299
00300 public:
00304 virtual Tbl get_solp(const Tbl& so) const = 0 ;
00308 virtual Tbl get_solh() const = 0 ;
00312 virtual void inc_l_quant() = 0 ;
00313 } ;
00314
00321 class Ope_poisson : public Ope_elementary {
00322
00323 protected:
00324 int l_quant ;
00325 int dzpuis ;
00326
00327 public:
00338 Ope_poisson (int nbr, int baser, double alf, double bet, int lq, int dz) ;
00339 Ope_poisson (const Ope_poisson&) ;
00340 virtual ~Ope_poisson() ;
00341
00343 int get_dzpuis() {return dzpuis ;} ;
00344
00346 int get_lquant() {return l_quant;} ;
00347
00348 private:
00352 virtual void do_ope_mat() const ;
00356 virtual void do_ope_cl() const ;
00360 virtual void do_non_dege() const ;
00361
00362 public:
00366 virtual Tbl get_solp(const Tbl& so) const ;
00370 virtual Tbl get_solh() const ;
00374 virtual void inc_l_quant() ;
00378 virtual void dec_l_quant() ;
00379 } ;
00380
00386 class Ope_helmholtz_minus : public Ope_elementary {
00387
00388 protected:
00389 int lq ;
00390 double masse ;
00391
00392 public:
00403 Ope_helmholtz_minus (int nbr, int baser, int lq, double alf, double bet,
00404 double mas) ;
00405 Ope_helmholtz_minus (const Ope_helmholtz_minus&) ;
00406 virtual ~Ope_helmholtz_minus() ;
00407
00408 private:
00412 virtual void do_ope_mat() const ;
00416 virtual void do_ope_cl() const ;
00420 virtual void do_non_dege() const ;
00421
00422 public:
00426 virtual Tbl get_solp(const Tbl& so) const ;
00430 virtual Tbl get_solh() const ;
00434 virtual void inc_l_quant() ;
00435 } ;
00436
00442 class Ope_helmholtz_plus : public Ope_elementary {
00443
00444 protected:
00445 int lq ;
00446 double masse ;
00447
00448 public:
00459 Ope_helmholtz_plus (int nbr, int baser, int lq, double alf, double bet,
00460 double mas) ;
00461 Ope_helmholtz_plus (const Ope_helmholtz_plus&) ;
00462 virtual ~Ope_helmholtz_plus() ;
00463
00464 private:
00468 virtual void do_ope_mat() const ;
00472 virtual void do_ope_cl() const ;
00476 virtual void do_non_dege() const ;
00477
00478 public:
00482 virtual Tbl get_solp(const Tbl& so) const ;
00486 virtual Tbl get_solh() const ;
00490 virtual void inc_l_quant() ;
00491 } ;
00492
00499 class Ope_sec_order_r2 : public Ope_elementary {
00500
00501 protected:
00502
00503 double a_param ;
00504 double b_param ;
00505 double c_param ;
00506
00507 public:
00520 Ope_sec_order_r2 (int nbr, int baser, double alf, double bet,
00521 double a, double b, double c) ;
00522
00523 Ope_sec_order_r2 (const Ope_sec_order_r2&) ;
00524 virtual ~Ope_sec_order_r2() ;
00525
00526 private:
00530 virtual void do_ope_mat() const ;
00534 virtual void do_ope_cl() const ;
00538 virtual void do_non_dege() const ;
00539
00540 public:
00544 virtual Tbl get_solp(const Tbl& so) const ;
00548 virtual Tbl get_solh() const ;
00552 virtual void inc_l_quant() ;
00553 } ;
00554
00561 class Ope_sec_order : public Ope_elementary {
00562
00563 protected:
00564
00565 double a_param ;
00566 double b_param ;
00567 double c_param ;
00568
00569 public:
00582 Ope_sec_order (int nbr, int baser, double alf, double bet,
00583 double a, double b, double c) ;
00584
00585 Ope_sec_order (const Ope_sec_order&) ;
00586 virtual ~Ope_sec_order() ;
00587
00588 private:
00592 virtual void do_ope_mat() const ;
00596 virtual void do_ope_cl() const ;
00600 virtual void do_non_dege() const ;
00601
00602 public:
00603
00607 virtual Tbl get_solp(const Tbl& so) const ;
00611 virtual Tbl get_solh() const ;
00615 virtual void inc_l_quant() ;
00616 } ;
00617
00629 class Ope_pois_vect_r : public Ope_poisson {
00630
00631 public:
00642 Ope_pois_vect_r (int nbr, int baser, double alf, double bet, int lq, int dz) ;
00643 Ope_pois_vect_r (const Ope_pois_vect_r&) ;
00644 virtual ~Ope_pois_vect_r() ;
00645
00646 private:
00650 virtual void do_ope_mat() const ;
00654 virtual void do_ope_cl() const ;
00658 virtual void do_non_dege() const ;
00659
00660 public:
00664 virtual Tbl get_solh() const ;
00665 } ;
00666
00675 class Ope_pois_tens_rr : public Ope_poisson {
00676
00677 public:
00688 Ope_pois_tens_rr (int nbr, int baser, double alf, double bet, int lq, int dz) ;
00689 Ope_pois_tens_rr (const Ope_pois_tens_rr&) ;
00690 virtual ~Ope_pois_tens_rr() ;
00691
00692 private:
00696 virtual void do_ope_mat() const ;
00700 virtual void do_ope_cl() const ;
00704 virtual void do_non_dege() const ;
00705
00706 public:
00710 virtual Tbl get_solh() const ;
00711 } ;
00712
00719 class Ope_poisson_2d : public Ope_elementary {
00720
00721 protected:
00722 int l_quant ;
00723 int dzpuis ;
00724
00725 public:
00736 Ope_poisson_2d (int nbr, int baser, double alf, double bet, int lq, int dz) ;
00737 Ope_poisson_2d (const Ope_poisson_2d&) ;
00738 virtual ~Ope_poisson_2d() ;
00739
00741 int get_dzpuis() {return dzpuis ;} ;
00742
00744 int get_lquant() {return l_quant;} ;
00745
00746 private:
00750 virtual void do_ope_mat() const ;
00754 virtual void do_ope_cl() const ;
00758 virtual void do_non_dege() const ;
00759
00760 public:
00764 virtual Tbl get_solp(const Tbl& so) const ;
00768 virtual Tbl get_solh() const ;
00772 virtual void inc_l_quant() ;
00776 virtual void dec_l_quant() ;
00777
00778 } ;
00779
00785 class Ope_helmholtz_minus_2d : public Ope_elementary {
00786
00787 protected:
00788 int l_quant ;
00789 double masse ;
00790 int dzpuis ;
00791
00792 public:
00804 Ope_helmholtz_minus_2d (int nbr, int baser, double alf, double bet, int lq, double masse, int dz) ;
00805 Ope_helmholtz_minus_2d (const Ope_helmholtz_minus_2d&) ;
00806 virtual ~Ope_helmholtz_minus_2d() ;
00807
00809 int get_dzpuis() {return dzpuis ;} ;
00810
00812 int get_lquant() {return l_quant;} ;
00813
00815 double get_masse() {return masse;} ;
00816
00817 private:
00821 virtual void do_ope_mat() const ;
00825 virtual void do_ope_cl() const ;
00829 virtual void do_non_dege() const ;
00830
00831 public:
00835 virtual Tbl get_solp(const Tbl& so) const ;
00839 virtual Tbl get_solh() const ;
00843 virtual void inc_l_quant() ;
00847 virtual void dec_l_quant() ;
00848 } ;
00849
00855 class Ope_poisson_pseudo_1d : public Ope_elementary {
00856
00857 protected:
00858 int l_quant ;
00859
00860 public:
00870 Ope_poisson_pseudo_1d (int nbr, int baser, double alf, double bet, int lq) ;
00871 Ope_poisson_pseudo_1d (const Ope_poisson_pseudo_1d&) ;
00872 virtual ~Ope_poisson_pseudo_1d() ;
00873
00875 int get_lquant() {return l_quant;} ;
00876
00877 private:
00881 virtual void do_ope_mat() const ;
00885 virtual void do_ope_cl() const ;
00889 virtual void do_non_dege() const ;
00890
00891 public:
00895 virtual Tbl get_solp(const Tbl& so) const ;
00899 virtual Tbl get_solh() const ;
00903 virtual void inc_l_quant() ;
00907 virtual void dec_l_quant() ;
00908
00909 } ;
00910
00911
00918 class Ope_helmholtz_minus_pseudo_1d : public Ope_elementary {
00919
00920 protected:
00921 int l_quant ;
00922 double masse ;
00923 int dzpuis ;
00924
00925 public:
00937 Ope_helmholtz_minus_pseudo_1d (int nbr, int baser, double alf, double bet,
00938 int lq, double masse, int dz) ;
00939 Ope_helmholtz_minus_pseudo_1d (const Ope_helmholtz_minus_pseudo_1d&) ;
00940 virtual ~Ope_helmholtz_minus_pseudo_1d() ;
00941
00943 int get_dzpuis() {return dzpuis ;} ;
00944
00946 int get_lquant() {return l_quant;} ;
00947
00949 double get_masse() {return masse;} ;
00950
00951 private:
00955 virtual void do_ope_mat() const ;
00959 virtual void do_ope_cl() const ;
00963 virtual void do_non_dege() const ;
00964
00965 public:
00969 virtual Tbl get_solp(const Tbl& so) const ;
00973 virtual Tbl get_solh() const ;
00977 virtual void inc_l_quant() ;
00981 virtual void dec_l_quant() ;
00982 } ;
00983
00990 class Ope_vorton : public Ope_elementary {
00991
00992 protected:
00993 int l_quant ;
00994 int dzpuis ;
00995
00996 public:
01007 Ope_vorton (int nbr, int baser, double alf, double bet, int lq, int dz) ;
01008 Ope_vorton (const Ope_vorton&) ;
01009 virtual ~Ope_vorton() ;
01010
01012 int get_dzpuis() {return dzpuis ;} ;
01013
01015 int get_lquant() {return l_quant;} ;
01016
01017 private:
01021 virtual void do_ope_mat() const ;
01025 virtual void do_ope_cl() const ;
01029 virtual void do_non_dege() const ;
01030
01031 public:
01035 virtual Tbl get_solp(const Tbl& so) const ;
01039 virtual Tbl get_solh() const ;
01043 virtual void inc_l_quant() ;
01047 virtual void dec_l_quant() ;
01048 } ;
01049
01050 # endif
01051