LORENE
ope_helmholtz_plus.C
1/*
2 * Copyright (c) 2003 Philippe Grandclement
3 *
4 * This file is part of LORENE.
5 *
6 * LORENE is free software; you can redistribute it and/or modify
7 * it under the terms of the GNU General Public License version 2
8 * as published by the Free Software Foundation.
9 *
10 * LORENE is distributed in the hope that it will be useful,
11 * but WITHOUT ANY WARRANTY; without even the implied warranty of
12 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13 * GNU General Public License for more details.
14 *
15 * You should have received a copy of the GNU General Public License
16 * along with LORENE; if not, write to the Free Software
17 * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
18 *
19 */
20
21
22
23/*
24 * $Id: ope_helmholtz_plus.C,v 1.8 2016/12/05 16:18:11 j_novak Exp $
25 * $Log: ope_helmholtz_plus.C,v $
26 * Revision 1.8 2016/12/05 16:18:11 j_novak
27 * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
28 *
29 * Revision 1.7 2014/10/13 08:53:33 j_novak
30 * Lorene classes and functions now belong to the namespace Lorene.
31 *
32 * Revision 1.6 2014/10/06 15:13:15 j_novak
33 * Modified #include directives to use c++ syntax.
34 *
35 * Revision 1.5 2007/05/06 10:48:13 p_grandclement
36 * Modification of a few operators for the vorton project
37 *
38 * Revision 1.4 2004/01/15 09:15:38 p_grandclement
39 * Modification and addition of the Helmholtz operators
40 *
41 * Revision 1.3 2003/12/11 16:12:10 e_gourgoulhon
42 * Changed sqrt(2) to sqrt(double(2)).
43 *
44 * Revision 1.2 2003/12/11 15:57:27 p_grandclement
45 * include stdlib.h encore ...
46 *
47 * Revision 1.1 2003/12/11 14:48:50 p_grandclement
48 * Addition of ALL (and that is a lot !) the files needed for the general elliptic solver ... UNDER DEVELOPEMENT...
49 *
50 *
51 * $Header: /cvsroot/Lorene/C++/Source/Ope_elementary/ope_helmholtz_plus.C,v 1.8 2016/12/05 16:18:11 j_novak Exp $
52 *
53 */
54#include <cmath>
55#include <cstdlib>
56
57#include "proto.h"
58#include "ope_elementary.h"
59
60// Standard constructor :
61namespace Lorene {
62Ope_helmholtz_plus::Ope_helmholtz_plus (int nbr, int base, int lquant, double alf,
63 double bet, double mas):
64 Ope_elementary(nbr, base, alf, bet), lq(lquant), masse (mas) {
65}
66
67// Constructor by copy :
71
72// Destructor :
74
75// True functions :
77 if (ope_mat != 0x0)
78 delete ope_mat ;
79
80 ope_mat = new Matrice
81 (helmholtz_plus_mat(nr, lq, alpha, beta, masse, base_r)) ;
82}
83
85 if (ope_mat == 0x0)
86 do_ope_mat() ;
87
88 if (ope_cl != 0x0)
89 delete ope_cl ;
90
91 ope_cl = new Matrice
92 (cl_helmholtz_plus(*ope_mat, base_r)) ;
93}
94
96 if (ope_cl == 0x0)
97 do_ope_cl() ;
98
99 if (non_dege != 0x0)
100 delete non_dege ;
101
102 non_dege = new Matrice
103 (prepa_helmholtz_plus_nondege(*ope_cl, base_r)) ;
104}
105
107
108 if (non_dege == 0x0)
109 do_non_dege() ;
110
111 Tbl res(solp_helmholtz_plus (*ope_mat, *non_dege, so, alpha, beta, base_r));
112
113 Tbl valeurs (val_solp (res, alpha, base_r)) ;
114 sp_plus = valeurs(0) ;
115 sp_minus = valeurs(1) ;
116 dsp_plus = valeurs(2) ;
117 dsp_minus = valeurs(3) ;
118
119 return res ;
120}
121
123
124 Tbl res (solh_helmholtz_plus (nr, lq, alpha, beta, masse, base_r)) ;
125
126 // Un peu tricky...
127 if (res.get_ndim() == 1) {
128 Tbl val_lim (val_solp (res, alpha, base_r)) ;
129
130 s_one_plus = val_lim(0) ;
131 s_one_minus = val_lim(1) ;
132 ds_one_plus = val_lim(2) ;
133 ds_one_minus = val_lim(3) ;
134
135 }
136 else {
137 Tbl auxi (nr) ;
138 auxi.set_etat_qcq() ;
139 for (int i=0 ; i<nr ; i++)
140 auxi.set(i) = res(0,i) ;
141
142 Tbl val_one (val_solp (auxi, alpha, base_r)) ;
143
144 s_one_plus = val_one(0) ;
145 s_one_minus = val_one(1) ;
146 ds_one_plus = val_one(2) ;
147 ds_one_minus = val_one(3) ;
148
149 for (int i=0 ; i<nr ; i++)
150 auxi.set(i) = res(1,i) ;
151
152 Tbl val_two (val_solp (auxi, alpha, base_r)) ;
153
154 s_two_plus = val_two(0) ;
155 s_two_minus = val_two(1) ;
156 ds_two_plus = val_two(2) ;
157 ds_two_minus = val_two(3) ;
158 }
159 return res ;
160}
161
162
164
165 cout << "inc_l_quant not implemented for Helmholtz operator." << endl ;
166 abort() ;
167}
168}
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 dsp_plus
Value of the derivative of the particular solution at the outer boundary.
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.
Ope_elementary(int nbr, int baser, double alf, double eta)
Standard constructor, protected because the class is an abstract one.
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 sp_minus
Value of the particular solution at the inner boundary.
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.
double sp_plus
Value of the particular solution at the outer boundary.
int nr
Number of radial points.
double dsp_minus
Value of the derivative of the particular solution at the inner boundary.
int lq
The quantum number l.
virtual void do_non_dege() const
Computes the non-degenerated matrix of the operator.
virtual Tbl get_solh() const
Computes the homogeneous solutions(s).
virtual ~Ope_helmholtz_plus()
Destructor.
virtual void do_ope_cl() const
Computes the banded-matrix of the operator.
Ope_helmholtz_plus(int nbr, int baser, int lq, double alf, double bet, double mas)
Standard constructor.
double masse
The mass parameter m .
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 inc_l_quant()
Increases the quatum number l by one unit (CURRENTLY NOT IMPLEMENTED).
Basic array class.
Definition tbl.h:161
int get_ndim() const
Gives the number of dimensions (ie dim.ndim).
Definition tbl.h:400
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition tbl.C:364
double & set(int i)
Read/write of a particular element (index i) (1D case).
Definition tbl.h:281
Lorene prototypes.
Definition app_hor.h:67