LORENE
ope_elementary.h
1 /*
2  * Definition of Lorene classes Ope_elementary
3  *
4  */
5 
6 /*
7  * Copyright (c) 2003 Philippe Grandclement
8  *
9  * This file is part of LORENE.
10  *
11  * LORENE is free software; you can redistribute it and/or modify
12  * it under the terms of the GNU General Public License version 2
13  * as published by the Free Software Foundation.
14  *
15  * LORENE is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
18  * GNU General Public License for more details.
19  *
20  * You should have received a copy of the GNU General Public License
21  * along with LORENE; if not, write to the Free Software
22  * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
23  *
24  */
25 
26 #ifndef __OPE_ELEMENTARY_H_
27 #define __OPE_ELEMENTARY_H_
28 
29 /*
30  * $Id: ope_elementary.h,v 1.12 2014/10/13 08:52:36 j_novak Exp $
31  * $Log: ope_elementary.h,v $
32  * Revision 1.12 2014/10/13 08:52:36 j_novak
33  * Lorene classes and functions now belong to the namespace Lorene.
34  *
35  * Revision 1.11 2007/05/06 10:48:08 p_grandclement
36  * Modification of a few operators for the vorton project
37  *
38  * Revision 1.10 2007/04/24 09:04:11 p_grandclement
39  * Addition of an operator for the vortons
40  *
41  * Revision 1.9 2004/12/23 16:30:14 j_novak
42  * New files and class for the solution of the rr component of the tensor Poisson
43  * equation.
44  *
45  * Revision 1.8 2004/08/24 09:14:40 p_grandclement
46  * Addition of some new operators, like Poisson in 2d... It now requieres the
47  * GSL library to work.
48  *
49  * Also, the way a variable change is stored by a Param_elliptic is changed and
50  * no longer uses Change_var but rather 2 Scalars. The codes using that feature
51  * will requiere some modification. (It should concern only the ones about monopoles)
52  *
53  * Revision 1.7 2004/06/22 08:49:57 p_grandclement
54  * Addition of everything needed for using the logarithmic mapping
55  *
56  * Revision 1.6 2004/06/14 15:07:10 j_novak
57  * New methods for the construction of the elliptic operator appearing in
58  * the vector Poisson equation (acting on eta).
59  *
60  * Revision 1.5 2004/05/10 15:28:21 j_novak
61  * First version of functions for the solution of the r-component of the
62  * vector Poisson equation.
63  *
64  * Revision 1.4 2004/03/23 14:54:45 j_novak
65  * More documentation
66  *
67  * Revision 1.3 2004/03/09 09:15:11 p_grandclement
68  * Correction pour le raccord...
69  *
70  * Revision 1.2 2004/03/05 09:18:48 p_grandclement
71  * Addition of operator sec_order_r2
72  *
73  * Revision 1.1 2003/12/11 14:57:00 p_grandclement
74  * I had forgotten the .h (sorry folks...)
75  *
76  * Revision 1.2 2001/12/11 06:44:41 e_gourgoulhon
77  *
78  * $Header: /cvsroot/Lorene/C++/Include/ope_elementary.h,v 1.12 2014/10/13 08:52:36 j_novak Exp $
79  *
80  */
81 
82 #include "matrice.h"
83 
84 namespace Lorene {
103 
104  protected:
105 
106  int nr ;
107  int base_r ;
108  double alpha ;
109  double beta ;
110 
114  mutable Matrice* ope_mat ;
118  mutable Matrice* ope_cl ;
122  mutable Matrice* non_dege ;
123 
127  mutable double s_one_plus ;
131  mutable double s_one_minus ;
136  mutable double ds_one_plus ;
141  mutable double ds_one_minus ;
142 
146  mutable double s_two_plus ;
150  mutable double s_two_minus ;
155  mutable double ds_two_plus ;
160  mutable double ds_two_minus ;
161 
165  mutable double sp_minus ;
169  mutable double sp_plus ;
173  mutable double dsp_minus ;
177  mutable double dsp_plus ;
178 
179  // Constructors destructor
180  protected:
189  explicit Ope_elementary (int nbr , int baser , double alf, double eta) ;
190  Ope_elementary (const Ope_elementary&) ;
191 
192  public:
193  virtual ~Ope_elementary() ;
194 
195  public:
200  double val_sh_one_minus() const {return s_one_minus ;} ;
205  double val_sh_one_plus() const {return s_one_plus ;} ;
211  double der_sh_one_minus() const {return ds_one_minus ;} ;
217  double der_sh_one_plus() const {return ds_one_plus ;} ;
218 
223  double val_sh_two_minus() const {return s_two_minus ;} ;
228  double val_sh_two_plus() const {return s_two_plus ;} ;
234  double der_sh_two_minus() const {return ds_two_minus ;} ;
240  double der_sh_two_plus() const {return ds_two_plus ;} ;
241 
245  double val_sp_minus() const {return sp_minus ;} ;
249  double val_sp_plus() const {return sp_plus ;} ;
254  double der_sp_minus() const {return dsp_minus ;} ;
259  double der_sp_plus() const {return dsp_plus ;} ;
260 
262  double get_alpha() const {return alpha ;} ;
263 
265  double get_beta() const {return beta ;} ;
266 
268  int get_base_r() const {return base_r ;} ;
269 
272  if (ope_mat ==0x0)
273  do_ope_mat() ;
274  return *ope_mat ;
275 }
276 
279  if (ope_cl ==0x0)
280  do_ope_cl() ;
281  return *ope_cl ;
282 }
283 
286  if (non_dege ==0x0)
287  do_non_dege() ;
288  return *non_dege ;
289 }
290  private:
294  virtual void do_ope_mat() const = 0 ;
298  virtual void do_ope_cl() const = 0 ;
302  virtual void do_non_dege() const = 0 ;
303 
304  public:
308  virtual Tbl get_solp(const Tbl& so) const = 0 ;
312  virtual Tbl get_solh() const = 0 ;
316  virtual void inc_l_quant() = 0 ;
317 } ;
318 
325 class Ope_poisson : public Ope_elementary {
326 
327  protected:
328  int l_quant ;
329  int dzpuis ;
330 
331  public:
342  Ope_poisson (int nbr, int baser, double alf, double bet, int lq, int dz) ;
343  Ope_poisson (const Ope_poisson&) ;
344  virtual ~Ope_poisson() ;
345 
347  int get_dzpuis() {return dzpuis ;} ;
348 
350  int get_lquant() {return l_quant;} ;
351 
352  private:
356  virtual void do_ope_mat() const ;
360  virtual void do_ope_cl() const ;
364  virtual void do_non_dege() const ;
365 
366  public:
370  virtual Tbl get_solp(const Tbl& so) const ;
374  virtual Tbl get_solh() const ;
378  virtual void inc_l_quant() ;
382  virtual void dec_l_quant() ;
383 } ;
384 
391 
392  protected:
393  int lq ;
394  double masse ;
395 
396  public:
407  Ope_helmholtz_minus (int nbr, int baser, int lq, double alf, double bet,
408  double mas) ;
410  virtual ~Ope_helmholtz_minus() ;
411 
412  private:
416  virtual void do_ope_mat() const ;
420  virtual void do_ope_cl() const ;
424  virtual void do_non_dege() const ;
425 
426  public:
430  virtual Tbl get_solp(const Tbl& so) const ;
434  virtual Tbl get_solh() const ;
438  virtual void inc_l_quant() ;
439 } ;
440 
447 
448  protected:
449  int lq ;
450  double masse ;
451 
452  public:
463  Ope_helmholtz_plus (int nbr, int baser, int lq, double alf, double bet,
464  double mas) ;
466  virtual ~Ope_helmholtz_plus() ;
467 
468  private:
472  virtual void do_ope_mat() const ;
476  virtual void do_ope_cl() const ;
480  virtual void do_non_dege() const ;
481 
482  public:
486  virtual Tbl get_solp(const Tbl& so) const ;
490  virtual Tbl get_solh() const ;
494  virtual void inc_l_quant() ;
495 } ;
496 
504 
505  protected:
506 
507  double a_param ;
508  double b_param ;
509  double c_param ;
510 
511  public:
524  Ope_sec_order_r2 (int nbr, int baser, double alf, double bet,
525  double a, double b, double c) ;
526 
528  virtual ~Ope_sec_order_r2() ;
529 
530  private:
534  virtual void do_ope_mat() const ;
538  virtual void do_ope_cl() const ;
542  virtual void do_non_dege() const ;
543 
544  public:
548  virtual Tbl get_solp(const Tbl& so) const ;
552  virtual Tbl get_solh() const ;
556  virtual void inc_l_quant() ;
557 } ;
558 
566 
567  protected:
568 
569  double a_param ;
570  double b_param ;
571  double c_param ;
572 
573  public:
586  Ope_sec_order (int nbr, int baser, double alf, double bet,
587  double a, double b, double c) ;
588 
589  Ope_sec_order (const Ope_sec_order&) ;
590  virtual ~Ope_sec_order() ;
591 
592  private:
596  virtual void do_ope_mat() const ;
600  virtual void do_ope_cl() const ;
604  virtual void do_non_dege() const ;
605 
606  public:
607 
611  virtual Tbl get_solp(const Tbl& so) const ;
615  virtual Tbl get_solh() const ;
619  virtual void inc_l_quant() ;
620 } ;
621 
633 class Ope_pois_vect_r : public Ope_poisson {
634 
635  public:
646  Ope_pois_vect_r (int nbr, int baser, double alf, double bet, int lq, int dz) ;
647  Ope_pois_vect_r (const Ope_pois_vect_r&) ;
648  virtual ~Ope_pois_vect_r() ;
649 
650  private:
654  virtual void do_ope_mat() const ;
658  virtual void do_ope_cl() const ;
662  virtual void do_non_dege() const ;
663 
664  public:
668  virtual Tbl get_solh() const ;
669 } ;
670 
680 
681  public:
692  Ope_pois_tens_rr (int nbr, int baser, double alf, double bet, int lq, int dz) ;
694  virtual ~Ope_pois_tens_rr() ;
695 
696  private:
700  virtual void do_ope_mat() const ;
704  virtual void do_ope_cl() const ;
708  virtual void do_non_dege() const ;
709 
710  public:
714  virtual Tbl get_solh() const ;
715 } ;
716 
724 
725  protected:
726  int l_quant ;
727  int dzpuis ;
728 
729  public:
740  Ope_poisson_2d (int nbr, int baser, double alf, double bet, int lq, int dz) ;
741  Ope_poisson_2d (const Ope_poisson_2d&) ;
742  virtual ~Ope_poisson_2d() ;
743 
745  int get_dzpuis() {return dzpuis ;} ;
746 
748  int get_lquant() {return l_quant;} ;
749 
750  private:
754  virtual void do_ope_mat() const ;
758  virtual void do_ope_cl() const ;
762  virtual void do_non_dege() const ;
763 
764  public:
768  virtual Tbl get_solp(const Tbl& so) const ;
772  virtual Tbl get_solh() const ;
776  virtual void inc_l_quant() ;
780  virtual void dec_l_quant() ;
781 
782 } ;
783 
790 
791  protected:
792  int l_quant ;
793  double masse ;
794  int dzpuis ;
795 
796  public:
808  Ope_helmholtz_minus_2d (int nbr, int baser, double alf, double bet, int lq, double masse, int dz) ;
810  virtual ~Ope_helmholtz_minus_2d() ;
811 
813  int get_dzpuis() {return dzpuis ;} ;
814 
816  int get_lquant() {return l_quant;} ;
817 
819  double get_masse() {return masse;} ;
820 
821  private:
825  virtual void do_ope_mat() const ;
829  virtual void do_ope_cl() const ;
833  virtual void do_non_dege() const ;
834 
835  public:
839  virtual Tbl get_solp(const Tbl& so) const ;
843  virtual Tbl get_solh() const ;
847  virtual void inc_l_quant() ;
851  virtual void dec_l_quant() ;
852 } ;
853 
860 
861  protected:
862  int l_quant ;
863 
864  public:
874  Ope_poisson_pseudo_1d (int nbr, int baser, double alf, double bet, int lq) ;
876  virtual ~Ope_poisson_pseudo_1d() ;
877 
879  int get_lquant() {return l_quant;} ;
880 
881  private:
885  virtual void do_ope_mat() const ;
889  virtual void do_ope_cl() const ;
893  virtual void do_non_dege() const ;
894 
895  public:
899  virtual Tbl get_solp(const Tbl& so) const ;
903  virtual Tbl get_solh() const ;
907  virtual void inc_l_quant() ;
911  virtual void dec_l_quant() ;
912 
913 } ;
914 
915 
923 
924  protected:
925  int l_quant ;
926  double masse ;
927  int dzpuis ;
928 
929  public:
941  Ope_helmholtz_minus_pseudo_1d (int nbr, int baser, double alf, double bet,
942  int lq, double masse, int dz) ;
944  virtual ~Ope_helmholtz_minus_pseudo_1d() ;
945 
947  int get_dzpuis() {return dzpuis ;} ;
948 
950  int get_lquant() {return l_quant;} ;
951 
953  double get_masse() {return masse;} ;
954 
955  private:
959  virtual void do_ope_mat() const ;
963  virtual void do_ope_cl() const ;
967  virtual void do_non_dege() const ;
968 
969  public:
973  virtual Tbl get_solp(const Tbl& so) const ;
977  virtual Tbl get_solh() const ;
981  virtual void inc_l_quant() ;
985  virtual void dec_l_quant() ;
986 } ;
987 
994 class Ope_vorton : public Ope_elementary {
995 
996  protected:
997  int l_quant ;
998  int dzpuis ;
999 
1000  public:
1011  Ope_vorton (int nbr, int baser, double alf, double bet, int lq, int dz) ;
1012  Ope_vorton (const Ope_vorton&) ;
1013  virtual ~Ope_vorton() ;
1014 
1016  int get_dzpuis() {return dzpuis ;} ;
1017 
1019  int get_lquant() {return l_quant;} ;
1020 
1021  private:
1025  virtual void do_ope_mat() const ;
1029  virtual void do_ope_cl() const ;
1033  virtual void do_non_dege() const ;
1034 
1035  public:
1039  virtual Tbl get_solp(const Tbl& so) const ;
1043  virtual Tbl get_solh() const ;
1047  virtual void inc_l_quant() ;
1051  virtual void dec_l_quant() ;
1052 } ;
1053 
1054 }
1055 #endif
1056 
double get_beta() const
Returns beta}.
virtual Tbl get_solp(const Tbl &so) const
Computes the particular solution, given the source so .
virtual void do_ope_cl() const
Computes the banded-matrix of the operator.
double alpha
Parameter of the associated mapping.
virtual void do_non_dege() const
Computes the non-degenerated matrix of the operator.
int dzpuis
the associated dzpuis, if in the compactified domain.
virtual ~Ope_helmholtz_minus()
Destructor.
virtual Tbl get_solh() const
Computes the homogeneous solutions(s).
virtual void dec_l_quant()
Decreases the quatum number l by one unit.
Definition: ope_poisson.C:158
Ope_poisson_pseudo_1d(int nbr, int baser, double alf, double bet, int lq)
Standard constructor.
Matrice get_ope_cl()
Returns the banded matrix representation.
virtual void inc_l_quant()
Increases the quatum number l by one unit.
virtual Tbl get_solp(const Tbl &so) const
Computes the particular solution, given the source so .
double val_sh_two_minus() const
Returns the value of the second homogeneous solution at the inner boundary.
virtual void do_ope_cl() const
Computes the banded-matrix of the operator.
virtual void do_ope_cl() const
Computes the banded-matrix of the operator.
Definition: ope_poisson.C:75
double val_sh_one_minus() const
Returns the value of the first homogeneous solution at the inner boundary.
double get_masse()
Returns the mass term.
virtual void do_ope_cl() const =0
Computes the banded-matrix of the operator.
Ope_helmholtz_minus_pseudo_1d(int nbr, int baser, double alf, double bet, int lq, double masse, int dz)
Standard constructor.
Matrice * ope_cl
Pointer on the banded-matrix of the operator.
int get_lquant()
Returns the quantum number l.
int get_dzpuis()
Returns the associated dzpuis, if in the compactified domain.
virtual void do_ope_cl() const
Computes the banded-matrix of the operator.
virtual ~Ope_pois_tens_rr()
Destructor.
double s_one_minus
Value of the first homogeneous solution at the inner boundary.
virtual void do_non_dege() const
Computes the non-degenerated matrix of the operator.
virtual void dec_l_quant()
Decreases the quatum number l by one unit.
double ds_two_minus
Value of the derivative of the second homogeneous solution at the inner boundary. ...
int l_quant
quantum number
virtual ~Ope_poisson_pseudo_1d()
Destructor.
virtual void do_ope_cl() const
Computes the banded-matrix of the operator.
double beta
Parameter of the associated mapping.
int dzpuis
the associated dzpuis, if in the compactified domain.
virtual void do_non_dege() const
Computes the non-degenerated matrix of the operator.
virtual void do_ope_cl() const
Computes the banded-matrix of the operator.
Lorene prototypes.
Definition: app_hor.h:67
virtual ~Ope_sec_order_r2()
Destructor.
Matrice * ope_mat
Pointer on the matrix representation of the operator.
double dsp_minus
Value of the derivative of the particular solution at the inner boundary.
virtual void do_non_dege() const
Computes the non-degenerated matrix of the operator.
virtual ~Ope_sec_order()
Destructor.
Definition: ope_sec_order.C:48
virtual Tbl get_solh() const
Computes the homogeneous solutions(s).
int get_lquant()
Returns the quantum number l.
int dzpuis
the associated dzpuis, if in the compactified domain.
int get_dzpuis()
Returns the associated dzpuis, if in the compactified domain.
virtual void do_non_dege() const
Computes the non-degenerated matrix of the operator.
double ds_two_plus
Value of the derivative of the second homogeneous solution at the outer boundary. ...
Class for the operator of the modified Helmholtz equation in pseudo-1d.
int get_lquant()
Returns the quantum number l.
virtual Tbl get_solh() const
Computes the homogeneous solutions(s).
virtual Tbl get_solp(const Tbl &so) const
Computes the particular solution, given the source so .
Definition: ope_poisson.C:97
Ope_vorton(int nbr, int baser, double alf, double bet, int lq, int dz)
Standard constructor.
Definition: ope_vorton.C:35
virtual Tbl get_solp(const Tbl &so) const
Computes the particular solution, given the source so .
Ope_helmholtz_plus(int nbr, int baser, int lq, double alf, double bet, double mas)
Standard constructor.
Matrice get_non_dege()
Returns the non degenerate matrix representation.
virtual Tbl get_solh() const
Computes the homogeneous solutions(s).
virtual Tbl get_solp(const Tbl &so) const
Computes the particular solution, given the source so .
virtual void inc_l_quant()
Increases the quatum number l by one unit (CURRENTLY NOT IMPLEMENTED)
virtual void do_ope_mat() const
Computes the matrix of the operator.
virtual void dec_l_quant()
Decreases the quatum number l by one unit.
virtual Tbl get_solp(const Tbl &so) const
Computes the particular solution, given the source so .
virtual void inc_l_quant()
Increases the quatum number l by one unit (CURRENTLY NOT IMPLEMENTED)
int l_quant
quantum number
virtual ~Ope_helmholtz_minus_2d()
Destructor.
virtual void do_ope_cl() const
Computes the banded-matrix of the operator.
double a_param
The parameter a .
int dzpuis
the associated dzpuis, if in the compactified domain.
virtual void inc_l_quant()
Increases the quatum number l by one unit (CURRENTLY NOT IMPLEMENTED)
Class for operator of the type .
double sp_minus
Value of the particular solution at the inner boundary.
virtual Tbl get_solh() const
Computes the homogeneous solutions(s).
virtual ~Ope_poisson()
Destructor.
Definition: ope_poisson.C:64
int get_dzpuis()
Returns the associated dzpuis, if in the compactified domain.
virtual void inc_l_quant()
Increases the quatum number l by one unit.
Definition: ope_poisson.C:139
Matrice get_ope_mat()
Returns the matrix representation.
virtual void do_ope_mat() const
Computes the matrix of the operator.
double get_alpha() const
Returns alpha .
virtual void do_ope_mat() const
Computes the matrix of the operator.
Definition: ope_poisson.C:67
int l_quant
quantum number
int base_r
Radial basis of decomposition.
double sp_plus
Value of the particular solution at the outer boundary.
virtual void do_non_dege() const
Computes the non-degenerated matrix of the operator.
double der_sp_plus() const
Returns the value of the derivative particular solution at the outer boundary.
Class for the operator appearing for the vortons.
virtual void do_non_dege() const
Computes the non-degenerated matrix of the operator.
double masse
The mass term.
virtual void do_non_dege() const
Computes the non-degenerated matrix of the operator.
virtual ~Ope_poisson_2d()
Destructor.
Class for the Helmholtz operator (m > 0).
int lq
The quantum number l.
double b_param
The parameter .
Ope_pois_vect_r(int nbr, int baser, double alf, double bet, int lq, int dz)
Standard constructor.
virtual Tbl get_solp(const Tbl &so) const
Computes the particular solution, given the source so .
virtual void inc_l_quant()
Increases the quatum number l by one unit.
Ope_helmholtz_minus_2d(int nbr, int baser, double alf, double bet, int lq, double masse, int dz)
Standard constructor.
Ope_elementary(int nbr, int baser, double alf, double eta)
Standard constructor, protected because the class is an abstract one.
double ds_one_plus
Value of the derivative of the first homogeneous solution at the outer boundary.
Class for operator of the type .
virtual void do_ope_mat() const
Computes the matrix of the operator.
double dsp_plus
Value of the derivative of the particular solution at the outer boundary.
int get_dzpuis()
Returns the associated dzpuis, if in the compactified domain.
double masse
The mass parameter m .
virtual void do_ope_mat() const
Computes the matrix of the operator.
virtual Tbl get_solh() const
Computes the homogeneous solutions(s).
double c_param
The parameter .
double der_sh_two_plus() const
Returns the value of the derivative of the second homogeneous solution at the outer boundary...
virtual void inc_l_quant()
Increases the quatum number l by one unit.
Ope_helmholtz_minus(int nbr, int baser, int lq, double alf, double bet, double mas)
Standard constructor.
virtual void dec_l_quant()
Decreases the quatum number l by one unit.
Definition: ope_vorton.C:54
double der_sh_two_minus() const
Returns the value of the derivative of the second homogeneous solution at the inner boundary...
double get_masse()
Returns the mass term.
virtual Tbl get_solp(const Tbl &so) const
Computes the particular solution, given the source so .
Matrix handling.
Definition: matrice.h:152
virtual Tbl get_solh() const
Computes the homogeneous solutions(s).
double a_param
The parameter .
int get_lquant()
Returns the quantum number l.
Class for the operator of the Poisson equation in 2D.
virtual void do_non_dege() const
Computes the non-degenerated matrix of the operator.
virtual void do_non_dege() const =0
Computes the non-degenerated matrix of the operator.
virtual void inc_l_quant()=0
Increases the quatum number l by one unit.
double val_sh_one_plus() const
Returns the value of the first homogeneous solution at the outer boundary.
Ope_pois_tens_rr(int nbr, int baser, double alf, double bet, int lq, int dz)
Standard constructor.
double b_param
The parameter b .
double s_two_minus
Value of the second homogeneous solution at the inner boundary.
Ope_sec_order(int nbr, int baser, double alf, double bet, double a, double b, double c)
Standard constructor.
Definition: ope_sec_order.C:35
Class for the operator of the rr component of the divergence-free tensor Poisson equation...
int get_base_r() const
Returns base_r}.
double der_sp_minus() const
Returns the value of the derivative particular solution at the inner boundary.
virtual void do_ope_cl() const
Computes the banded-matrix of the operator.
double val_sp_plus() const
Returns the value of the particular solution at the outer boundary.
Ope_sec_order_r2(int nbr, int baser, double alf, double bet, double a, double b, double c)
Standard constructor.
Class for the operator of the Poisson equation in pseudo 1d.
virtual Tbl get_solp(const Tbl &so) const
Computes the particular solution, given the source so .
Ope_poisson_2d(int nbr, int baser, double alf, double bet, int lq, int dz)
Standard constructor.
virtual void do_ope_mat() const
Computes the matrix of the operator.
Class for the operator of the Poisson equation (i.e.
double val_sp_minus() const
Returns the value of the particular solution at the inner boundary.
virtual void do_ope_cl() const
Computes the banded-matrix of the operator.
virtual Tbl get_solh() const
Computes the homogeneous solutions(s).
Definition: ope_poisson.C:113
int dzpuis
the associated dzpuis, if in the compactified domain.
virtual Tbl get_solh() const
Computes the homogeneous solutions(s).
Class for the operator of the r component of the vector Poisson equation.
virtual void do_ope_mat() const
Computes the matrix of the operator.
virtual Tbl get_solp(const Tbl &so) const
Computes the particular solution, given the source so .
virtual void dec_l_quant()
Decreases the quatum number l by one unit.
double s_one_plus
Value of the first homogeneous solution at the outer boundary.
Basic class for elementary elliptic operators.
virtual ~Ope_helmholtz_plus()
Destructor.
int lq
The quantum number l.
virtual void inc_l_quant()
Increases the quatum number by one unit (CURRENTLY NOT IMPLEMENTED)
Definition: ope_sec_order.C:50
virtual Tbl get_solp(const Tbl &so) const =0
Computes the particular solution, given the source so .
virtual void do_ope_mat() const
Computes the matrix of the operator.
virtual void do_non_dege() const
Computes the non-degenerated matrix of the operator.
double der_sh_one_plus() const
Returns the value of the derivative of the first homogeneous solution at the outer boundary...
virtual void do_ope_mat() const
Computes the matrix of the operator.
virtual void do_ope_mat() const =0
Computes the matrix of the operator.
virtual void do_ope_cl() const
Computes the banded-matrix of the operator.
Class for the Helmholtz operator ( ).
virtual Tbl get_solh() const
Computes the homogeneous solutions(s).
int nr
Number of radial points.
double ds_one_minus
Value of the derivative of the first homogeneous solution at the inner boundary.
virtual Tbl get_solh() const =0
Computes the homogeneous solutions(s).
Class for the operator of the Helmholtz equation in 2D.
virtual void do_non_dege() const
Computes the non-degenerated matrix of the operator.
virtual void do_ope_mat() const
Computes the matrix of the operator.
double masse
The mass parameter m .
double der_sh_one_minus() const
Returns the value of the derivative of the first homogeneous solution at the inner boundary...
Basic array class.
Definition: tbl.h:164
virtual ~Ope_pois_vect_r()
Destructor.
virtual void inc_l_quant()
Increases the quatum number l by one unit.
virtual void do_ope_cl() const
Computes the banded-matrix of the operator.
virtual void do_ope_mat() const
Computes the matrix of the operator.
virtual void do_ope_mat() const
Computes the matrix of the operator.
double c_param
The parameter c .
double s_two_plus
Value of the second homogeneous solution at the outer boundary.
virtual void dec_l_quant()
Decreases the quatum number l by one unit.
double val_sh_two_plus() const
Returns the value of the second homogeneous solution at the outer boundary.
int get_lquant()
Returns the quantum number l.
int get_dzpuis()
Returns the associated dzpuis, if in the compactified domain.
virtual void do_ope_cl() const
Computes the banded-matrix of the operator.
virtual void do_non_dege() const
Computes the non-degenerated matrix of the operator.
Definition: ope_poisson.C:86
virtual Tbl get_solh() const
Computes the homogeneous solutions(s).
virtual ~Ope_elementary()
Destructor.
virtual void inc_l_quant()
Increases the quatum number l by one unit.
Definition: ope_vorton.C:48
virtual ~Ope_vorton()
Destructor.
Definition: ope_vorton.C:46
int get_lquant()
Returns the quantum number l.
Ope_poisson(int nbr, int baser, double alf, double bet, int lq, int dz)
Standard constructor.
Definition: ope_poisson.C:49
virtual Tbl get_solh() const
Computes the homogeneous solutions(s).
Matrice * non_dege
Pointer on the non-degenerated matrix of the operator.