LORENE
ope_poisson.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_poisson.C,v 1.4 2016/12/05 16:18:11 j_novak Exp $
25 * $Log: ope_poisson.C,v $
26 * Revision 1.4 2016/12/05 16:18:11 j_novak
27 * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
28 *
29 * Revision 1.3 2014/10/13 08:53:33 j_novak
30 * Lorene classes and functions now belong to the namespace Lorene.
31 *
32 * Revision 1.2 2004/06/14 15:07:11 j_novak
33 * New methods for the construction of the elliptic operator appearing in
34 * the vector Poisson equation (acting on eta).
35 *
36 * Revision 1.1 2003/12/11 14:48:50 p_grandclement
37 * Addition of ALL (and that is a lot !) the files needed for the general elliptic solver ... UNDER DEVELOPEMENT...
38 *
39 *
40 * $Header: /cvsroot/Lorene/C++/Source/Ope_elementary/ope_poisson.C,v 1.4 2016/12/05 16:18:11 j_novak Exp $
41 *
42 */
43
44#include "proto.h"
45#include "ope_elementary.h"
46
47// Standard constructor :
48namespace Lorene {
49Ope_poisson::Ope_poisson (int nbr, int baser, double alf, double bet, int lq, int dz):
50 Ope_elementary(nbr, baser, alf, bet), l_quant (lq),
51 dzpuis (dz) {
52
53 assert ((dzpuis==2) || (dzpuis==3) || (dzpuis==4)) ;
54}
55
56// Constructor by copy :
58 l_quant (so.l_quant), dzpuis (so.dzpuis) {
59
60 assert ((dzpuis==2) || (dzpuis==3) || (dzpuis==4)) ;
61}
62
63// Destructor :
65
66// True functions :
68 if (ope_mat != 0x0)
69 delete ope_mat ;
70
71 ope_mat = new Matrice
72 (laplacien_mat(nr, l_quant, beta/alpha, dzpuis, base_r)) ;
73}
74
76 if (ope_mat == 0x0)
77 do_ope_mat() ;
78
79 if (ope_cl != 0x0)
80 delete ope_cl ;
81
82 ope_cl = new Matrice
83 (combinaison(*ope_mat, l_quant, beta/alpha, dzpuis, base_r)) ;
84}
85
87 if (ope_cl == 0x0)
88 do_ope_cl() ;
89
90 if (non_dege != 0x0)
91 delete non_dege ;
92
93 non_dege = new Matrice
94 (prepa_nondege(*ope_cl, l_quant, beta/alpha, dzpuis, base_r)) ;
95}
96
97Tbl Ope_poisson::get_solp (const Tbl& so) const {
98
99 if (non_dege == 0x0)
100 do_non_dege() ;
101
102 Tbl res(solp(*ope_mat, *non_dege, alpha, beta, so, dzpuis, base_r)) ;
103
104 Tbl valeurs (val_solp (res, alpha, base_r)) ;
105 sp_plus = valeurs(0) ;
106 sp_minus = valeurs(1) ;
107 dsp_plus = valeurs(2) ;
108 dsp_minus = valeurs(3) ;
109
110 return res ;
111}
112
114
115 Tbl valeurs (val_solh (l_quant, alpha, beta, base_r)) ;
116 if (valeurs.get_ndim() == 2) {
117 // cas 2 sh
118 s_one_plus = valeurs(0,0) ;
119 s_one_minus = valeurs(0,1) ;
120 ds_one_plus = valeurs(0,2) ;
121 ds_one_minus = valeurs(0,3) ;
122
123 s_two_plus = valeurs(1,0) ;
124 s_two_minus = valeurs(1,1) ;
125 ds_two_plus = valeurs(1,2) ;
126 ds_two_minus = valeurs(1,3) ;
127 }
128 else {
129 // cas 1 sh :
130 s_one_plus = valeurs(0) ;
131 s_one_minus = valeurs(1) ;
132 ds_one_plus = valeurs(2) ;
133 ds_one_minus = valeurs(3) ;
134 }
135
136 return solh(nr, l_quant, beta/alpha, base_r) ;
137}
138
140
141 if (ope_mat != 0x0) {
142 delete ope_mat ;
143 ope_mat = 0x0 ;
144 }
145
146 if (ope_cl != 0x0) {
147 delete ope_cl ;
148 ope_cl = 0x0 ;
149 }
150
151 if (non_dege != 0x0) {
152 delete non_dege ;
153 non_dege = 0x0 ;
154 }
155 l_quant ++ ;
156}
157
159
160 assert(l_quant > 0) ;
161
162 if (ope_mat != 0x0) {
163 delete ope_mat ;
164 ope_mat = 0x0 ;
165 }
166
167 if (ope_cl != 0x0) {
168 delete ope_cl ;
169 ope_cl = 0x0 ;
170 }
171
172 if (non_dege != 0x0) {
173 delete non_dege ;
174 non_dege = 0x0 ;
175 }
176 l_quant -- ;
177}
178}
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.
virtual Tbl get_solh() const
Computes the homogeneous solutions(s).
virtual ~Ope_poisson()
Destructor.
Definition ope_poisson.C:64
int dzpuis
the associated dzpuis, if in the compactified domain.
virtual void dec_l_quant()
Decreases the quatum number l by one unit.
virtual void do_non_dege() const
Computes the non-degenerated matrix of the operator.
Definition ope_poisson.C:86
virtual void do_ope_cl() const
Computes the banded-matrix of the operator.
Definition ope_poisson.C:75
virtual void do_ope_mat() const
Computes the matrix of the operator.
Definition ope_poisson.C:67
virtual void inc_l_quant()
Increases the quatum number l by one unit.
int l_quant
quantum number
virtual Tbl get_solp(const Tbl &so) const
Computes the particular solution, given the source so .
Definition ope_poisson.C:97
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
Lorene prototypes.
Definition app_hor.h:67