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
21char ope_helmholtz_plus_C[] = "$Header: /cvsroot/Lorene/C++/Source/Ope_elementary/ope_helmholtz_plus.C,v 1.7 2014/10/13 08:53:33 j_novak Exp $" ;
22
23/*
24 * $Id: ope_helmholtz_plus.C,v 1.7 2014/10/13 08:53:33 j_novak Exp $
25 * $Log: ope_helmholtz_plus.C,v $
26 * Revision 1.7 2014/10/13 08:53:33 j_novak
27 * Lorene classes and functions now belong to the namespace Lorene.
28 *
29 * Revision 1.6 2014/10/06 15:13:15 j_novak
30 * Modified #include directives to use c++ syntax.
31 *
32 * Revision 1.5 2007/05/06 10:48:13 p_grandclement
33 * Modification of a few operators for the vorton project
34 *
35 * Revision 1.4 2004/01/15 09:15:38 p_grandclement
36 * Modification and addition of the Helmholtz operators
37 *
38 * Revision 1.3 2003/12/11 16:12:10 e_gourgoulhon
39 * Changed sqrt(2) to sqrt(double(2)).
40 *
41 * Revision 1.2 2003/12/11 15:57:27 p_grandclement
42 * include stdlib.h encore ...
43 *
44 * Revision 1.1 2003/12/11 14:48:50 p_grandclement
45 * Addition of ALL (and that is a lot !) the files needed for the general elliptic solver ... UNDER DEVELOPEMENT...
46 *
47 *
48 * $Header: /cvsroot/Lorene/C++/Source/Ope_elementary/ope_helmholtz_plus.C,v 1.7 2014/10/13 08:53:33 j_novak Exp $
49 *
50 */
51#include <cmath>
52#include <cstdlib>
53
54#include "proto.h"
55#include "ope_elementary.h"
56
57// Standard constructor :
58namespace Lorene {
60 double bet, double mas):
61 Ope_elementary(nbr, base, alf, bet), lq(lquant), masse (mas) {
62}
63
64// Constructor by copy :
68
69// Destructor :
71
72// True functions :
74 if (ope_mat != 0x0)
75 delete ope_mat ;
76
77 ope_mat = new Matrice
78 (helmholtz_plus_mat(nr, lq, alpha, beta, masse, base_r)) ;
79}
80
82 if (ope_mat == 0x0)
83 do_ope_mat() ;
84
85 if (ope_cl != 0x0)
86 delete ope_cl ;
87
88 ope_cl = new Matrice
89 (cl_helmholtz_plus(*ope_mat, base_r)) ;
90}
91
93 if (ope_cl == 0x0)
94 do_ope_cl() ;
95
96 if (non_dege != 0x0)
97 delete non_dege ;
98
99 non_dege = new Matrice
100 (prepa_helmholtz_plus_nondege(*ope_cl, base_r)) ;
101}
102
104
105 if (non_dege == 0x0)
106 do_non_dege() ;
107
108 Tbl res(solp_helmholtz_plus (*ope_mat, *non_dege, so, alpha, beta, base_r));
109
110 Tbl valeurs (val_solp (res, alpha, base_r)) ;
111 sp_plus = valeurs(0) ;
112 sp_minus = valeurs(1) ;
113 dsp_plus = valeurs(2) ;
114 dsp_minus = valeurs(3) ;
115
116 return res ;
117}
118
120
121 Tbl res (solh_helmholtz_plus (nr, lq, alpha, beta, masse, base_r)) ;
122
123 // Un peu tricky...
124 if (res.get_ndim() == 1) {
125 Tbl val_lim (val_solp (res, alpha, base_r)) ;
126
127 s_one_plus = val_lim(0) ;
128 s_one_minus = val_lim(1) ;
129 ds_one_plus = val_lim(2) ;
130 ds_one_minus = val_lim(3) ;
131
132 }
133 else {
134 Tbl auxi (nr) ;
135 auxi.set_etat_qcq() ;
136 for (int i=0 ; i<nr ; i++)
137 auxi.set(i) = res(0,i) ;
138
139 Tbl val_one (val_solp (auxi, alpha, base_r)) ;
140
141 s_one_plus = val_one(0) ;
142 s_one_minus = val_one(1) ;
143 ds_one_plus = val_one(2) ;
144 ds_one_minus = val_one(3) ;
145
146 for (int i=0 ; i<nr ; i++)
147 auxi.set(i) = res(1,i) ;
148
149 Tbl val_two (val_solp (auxi, alpha, base_r)) ;
150
151 s_two_plus = val_two(0) ;
152 s_two_minus = val_two(1) ;
153 ds_two_plus = val_two(2) ;
154 ds_two_minus = val_two(3) ;
155 }
156 return res ;
157}
158
159
161
162 cout << "inc_l_quant not implemented for Helmholtz operator." << endl ;
163 abort() ;
164}
165}
Time evolution with partial storage (*** under development ***).
Definition evolution.h:371
Matrix handling.
Definition matrice.h:152
Basic class for elementary elliptic operators.
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.
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.
Class for the Helmholtz operator (m > 0).
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
Lorene prototypes.
Definition app_hor.h:64