LORENE
ope_pois_tens_rr.C
1/*
2 * Methods of class Ope_pois_tens_rr
3 *
4 * (see file ope_elementary.h for documentation)
5 *
6 */
7
8/*
9 * Copyright (c) 2004 Jerome Novak
10 *
11 * This file is part of LORENE.
12 *
13 * LORENE is free software; you can redistribute it and/or modify
14 * it under the terms of the GNU General Public License version 2
15 * as published by the Free Software Foundation.
16 *
17 * LORENE is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
20 * GNU General Public License for more details.
21 *
22 * You should have received a copy of the GNU General Public License
23 * along with LORENE; if not, write to the Free Software
24 * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
25 *
26 */
27
28
29
30/*
31 * $Id: ope_pois_tens_rr.C,v 1.3 2016/12/05 16:18:12 j_novak Exp $
32 * $Log: ope_pois_tens_rr.C,v $
33 * Revision 1.3 2016/12/05 16:18:12 j_novak
34 * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
35 *
36 * Revision 1.2 2014/10/13 08:53:34 j_novak
37 * Lorene classes and functions now belong to the namespace Lorene.
38 *
39 * Revision 1.1 2004/12/23 16:30:15 j_novak
40 * New files and class for the solution of the rr component of the tensor Poisson
41 * equation.
42 *
43 *
44 *
45 * $Header: /cvsroot/Lorene/C++/Source/Ope_elementary/Ope_pois_vect_r/ope_pois_tens_rr.C,v 1.3 2016/12/05 16:18:12 j_novak Exp $
46 *
47 */
48
49#include"type_parite.h"
50#include "ope_elementary.h"
51
52namespace Lorene {
53Matrice ope_ptens_rr_mat(int , int , double , int, int ) ;
54Tbl sh_ptens_rr(int, int, double, int) ;
55Matrice cl_ptens_rr (const Matrice&, int, double, int, int) ;
56Matrice nondeg_ptens_rr (const Matrice&, int, double, int, int) ;
57Tbl val_solh (int, double, double, int) ;
58
59// Standard constructor :
60Ope_pois_tens_rr::Ope_pois_tens_rr (int nbr, int baser, double alf, double bet, int lq, int dz):
61 Ope_poisson(nbr, baser, alf, bet, lq, dz)
62{
63 assert (dzpuis == 4) ;
64 assert (l_quant > 1) ;
65}
66
67// Constructor by copy :
69{
70 assert (dzpuis == 4) ;
71 assert (l_quant > 1) ;
72}
73
74// Destructor :
76
77// True functions :
79 if (ope_mat != 0x0)
80 delete ope_mat ;
81
82 ope_mat = new Matrice
83 (ope_ptens_rr_mat(nr, l_quant, beta/alpha, dzpuis, base_r)) ;
84}
85
87 if (ope_mat == 0x0)
88 do_ope_mat() ;
89
90 if (ope_cl != 0x0)
91 delete ope_cl ;
92
93 ope_cl = new Matrice
94 (cl_ptens_rr(*ope_mat, l_quant, beta/alpha, dzpuis, base_r)) ;
95}
96
98 if (ope_cl == 0x0)
99 do_ope_cl() ;
100
101 if (non_dege != 0x0)
102 delete non_dege ;
103
104 non_dege = new Matrice
105 (nondeg_ptens_rr(*ope_cl, l_quant, beta/alpha, dzpuis, base_r)) ;
106}
107
109
110 assert (l_quant > 1) ;
111 int l1 = l_quant - 2 ;
112 int l2 = l_quant + 2 ;
113
114 Tbl valeurs1 (val_solh (l1, alpha, beta, base_r)) ;
115 Tbl valeurs2 (val_solh (l2, alpha, beta, base_r)) ;
116
117 assert (valeurs1.get_ndim() == valeurs2.get_ndim()) ;
118
119 if (valeurs1.get_ndim() == 2) {
120 // cas 2 sh
121 s_one_plus = valeurs1(0,0) ;
122 s_one_minus = valeurs1(0,1) ;
123 ds_one_plus = valeurs1(0,2) ;
124 ds_one_minus = valeurs1(0,3) ;
125
126 s_two_plus = valeurs2(1,0) ;
127 s_two_minus = valeurs2(1,1) ;
128 ds_two_plus = valeurs2(1,2) ;
129 ds_two_minus = valeurs2(1,3) ;
130 }
131 else {
132 // cas 1 sh :
133 Tbl& valeurs = (base_r == R_CHEBU) ? valeurs2 : valeurs1 ;
134 s_one_plus = valeurs(0) ;
135 s_one_minus = valeurs(1) ;
136 ds_one_plus = valeurs(2) ;
137 ds_one_minus = valeurs(3) ;
138 }
139
140 return sh_ptens_rr(nr, l_quant, beta/alpha, base_r) ;
141}
142
143}
Matrix handling.
Definition matrice.h:152
Matrice * ope_mat
Pointer on the matrix representation of the operator.
double ds_two_minus
Value of the derivative of the second homogeneous solution at the inner boundary.
double s_two_plus
Value of the second homogeneous solution at the outer boundary.
double s_one_minus
Value of the first homogeneous solution at the inner boundary.
double beta
Parameter of the associated mapping.
double ds_one_plus
Value of the derivative of the first homogeneous solution at the outer boundary.
double ds_one_minus
Value of the derivative of the first homogeneous solution at the inner boundary.
double alpha
Parameter of the associated mapping.
double s_two_minus
Value of the second homogeneous solution at the inner boundary.
int base_r
Radial basis of decomposition.
double s_one_plus
Value of the first homogeneous solution at the outer boundary.
Matrice * ope_cl
Pointer on the banded-matrix of the operator.
double ds_two_plus
Value of the derivative of the second homogeneous solution at the outer boundary.
Matrice * non_dege
Pointer on the non-degenerated matrix of the operator.
int nr
Number of radial points.
virtual void do_ope_cl() const
Computes the banded-matrix of the operator.
virtual Tbl get_solh() const
Computes the homogeneous solutions(s).
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.
virtual ~Ope_pois_tens_rr()
Destructor.
Ope_pois_tens_rr(int nbr, int baser, double alf, double bet, int lq, int dz)
Standard constructor.
int dzpuis
the associated dzpuis, if in the compactified domain.
int l_quant
quantum number
Ope_poisson(int nbr, int baser, double alf, double bet, int lq, int dz)
Standard constructor.
Definition ope_poisson.C:49
Basic array class.
Definition tbl.h:161
int get_ndim() const
Gives the number of dimensions (ie dim.ndim).
Definition tbl.h:400
#define R_CHEBU
base de Chebychev ordinaire (fin), dev. en 1/r
Lorene prototypes.
Definition app_hor.h:67