LORENE
sol_elliptic.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: sol_elliptic.C,v 1.11 2016/12/05 16:18:10 j_novak Exp $
25 * $Log: sol_elliptic.C,v $
26 * Revision 1.11 2016/12/05 16:18:10 j_novak
27 * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
28 *
29 * Revision 1.10 2014/10/13 08:53:30 j_novak
30 * Lorene classes and functions now belong to the namespace Lorene.
31 *
32 * Revision 1.9 2014/10/06 15:16:10 j_novak
33 * Modified #include directives to use c++ syntax.
34 *
35 * Revision 1.8 2008/07/10 11:20:33 p_grandclement
36 * mistake fixed in solh_helmholtz_minus
37 *
38 * Revision 1.7 2008/07/09 06:51:58 p_grandclement
39 * some corrections to helmholtz minus in the nucleus
40 *
41 * Revision 1.6 2005/01/25 15:44:20 j_novak
42 * Fixed a problem in the definition of nr at the end.
43 *
44 * Revision 1.5 2004/11/25 08:14:56 j_novak
45 * Modifs. comments
46 *
47 * Revision 1.4 2004/08/24 09:14:44 p_grandclement
48 * Addition of some new operators, like Poisson in 2d... It now requieres the
49 * GSL library to work.
50 *
51 * Also, the way a variable change is stored by a Param_elliptic is changed and
52 * no longer uses Change_var but rather 2 Scalars. The codes using that feature
53 * will requiere some modification. (It should concern only the ones about monopoles)
54 *
55 * Revision 1.3 2004/05/26 20:32:44 p_grandclement
56 * utilisation of val_r_jk
57 *
58 * Revision 1.2 2003/12/19 16:21:49 j_novak
59 * Shadow hunt
60 *
61 * Revision 1.1 2003/12/11 14:48:49 p_grandclement
62 * Addition of ALL (and that is a lot !) the files needed for the general elliptic solver ... UNDER DEVELOPEMENT...
63 *
64 *
65 * $Header: /cvsroot/Lorene/C++/Source/Non_class_members/PDE/sol_elliptic.C,v 1.11 2016/12/05 16:18:10 j_novak Exp $
66 *
67 */
68
69// Header C :
70#include <cstdlib>
71#include <cmath>
72
73// Headers Lorene :
74#include "param_elliptic.h"
75#include "tbl.h"
76#include "mtbl_cf.h"
77#include "map.h"
78
79
80 //----------------------------------------------
81 // Version Mtbl_cf
82 //----------------------------------------------
83
84
85
86namespace Lorene {
87Mtbl_cf elliptic_solver (const Param_elliptic& ope_var, const Mtbl_cf& source) {
88 // Verifications d'usage sur les zones
89 int nz = source.get_mg()->get_nzone() ;
90 assert (nz>1) ;
91 assert (source.get_mg()->get_type_r(0) == RARE) ;
92 assert (source.get_mg()->get_type_r(nz-1) == UNSURR) ;
93 for (int l=1 ; l<nz-1 ; l++)
94 assert(source.get_mg()->get_type_r(l) == FIN) ;
95
96 // donnees sur la zone
97 int nr, nt, np ;
98
99 //Rangement des valeurs intermediaires
100 Tbl *so ;
101 Tbl *sol_hom ;
102 Tbl *sol_part ;
103
104
105 // Rangement des solutions, avant raccordement
106 Mtbl_cf solution_part(source.get_mg(), source.base) ;
107 Mtbl_cf solution_hom_un(source.get_mg(), source.base) ;
108 Mtbl_cf solution_hom_deux(source.get_mg(), source.base) ;
109 Mtbl_cf resultat(source.get_mg(), source.base) ;
110
111 solution_part.annule_hard() ;
112 solution_hom_un.annule_hard() ;
113 solution_hom_deux.annule_hard() ;
114 resultat.annule_hard() ;
115
116 // Computation of the SP and SH's in every domain ...
117 int conte = 0 ;
118 for (int zone=0 ; zone<nz ; zone++) {
119
120 nr = source.get_mg()->get_nr(zone) ;
121 nt = source.get_mg()->get_nt(zone) ;
122 np = source.get_mg()->get_np(zone) ;
123
124 for (int k=0 ; k<np+1 ; k++)
125 for (int j=0 ; j<nt ; j++) {
126
127 if (ope_var.operateurs[conte] != 0x0) {
128
129 // Calcul de la SH
130 sol_hom = new Tbl(ope_var.operateurs[conte]->get_solh()) ;
131
132 //Calcul de la SP
133 so = new Tbl(nr) ;
134 so->set_etat_qcq() ;
135 for (int i=0 ; i<nr ; i++)
136 so->set(i) = source(zone, k, j, i) ;
137
138 sol_part = new Tbl(ope_var.operateurs[conte]->get_solp(*so)) ;
139
140 // Rangement dans les tableaux globaux ;
141 for (int i=0 ; i<nr ; i++) {
142 solution_part.set(zone, k, j, i) = (*sol_part)(i) ;
143 if (sol_hom->get_ndim()==1)
144 solution_hom_un.set(zone, k, j, i) = (*sol_hom)(i) ;
145 else
146 {
147 solution_hom_un.set(zone, k, j, i) = (*sol_hom)(0,i) ;
148 solution_hom_deux.set(zone, k, j, i) = (*sol_hom)(1,i) ;
149 }
150 }
151
152 delete so ;
153 delete sol_hom ;
154 delete sol_part ;
155
156 }
157 conte ++ ;
158 }
159 }
160
161
162 //-------------------------------------------------
163 // ON EST PARTI POUR LE RACCORD (Be carefull ....)
164 //-------------------------------------------------
165
166 // C'est pas simple toute cette sombre affaire...
167 // Que le cas meme nombre de points dans chaque domaines...
168
169 int start = 0 ;
170 for (int k=0 ; k<source.get_mg()->get_np(0)+1 ; k++)
171 for (int j=0 ; j<source.get_mg()->get_nt(0) ; j++) {
172 if (ope_var.operateurs[start] != 0x0) {
173
174 int taille = 2*nz - 2 ;
175 Matrice systeme (taille, taille) ;
176 systeme.set_etat_qcq() ;
177 for (int i=0 ; i<taille ; i++)
178 for (int j2=0 ; j2<taille ; j2++)
179 systeme.set(i,j2) = 0 ;
180 Tbl sec_membre (taille) ;
181 sec_membre.set_etat_qcq() ;
182 for (int i=0 ; i<taille ; i++)
183 sec_membre.set(i) = 0 ;
184
185 //---------
186 // Noyau :
187 //---------
188 conte = start ;
189
190 systeme.set(0,0) = ope_var.G_plus(0) *
191 ope_var.operateurs[conte]->val_sh_one_plus() ;
192 systeme.set(1,0) =
193 ope_var.dG_plus(0) * ope_var.operateurs[conte]->val_sh_one_plus() +
194 ope_var.G_plus(0) * ope_var.operateurs[conte]->der_sh_one_plus() ;
195
196 sec_membre.set(0) -= ope_var.F_plus(0,k,j) +
197 ope_var.G_plus(0) * ope_var.operateurs[conte]->val_sp_plus() ;
198 sec_membre.set(1) -= ope_var.dF_plus(0,k,j) +
199 ope_var.dG_plus(0) * ope_var.operateurs[conte]->val_sp_plus() +
200 ope_var.G_plus(0) * ope_var.operateurs[conte]->der_sp_plus() ;
201
202 //----------
203 // SHELLS :
204 //----------
205
206 for (int l=1 ; l<nz-1 ; l++) {
207
208 // On se met au bon endroit :
209 int np_prec = source.get_mg()->get_np(l-1) ;
210 int nt_prec = source.get_mg()->get_nt(l-1) ;
211 conte += (np_prec+1)*nt_prec ;
212
213 systeme.set(2*l-2, 2*l-1) = -ope_var.G_minus(l) *
214 ope_var.operateurs[conte]->val_sh_one_minus() ;
215 systeme.set(2*l-2, 2*l) = - ope_var.G_minus(l) *
216 ope_var.operateurs[conte]->val_sh_two_minus() ;
217 systeme.set(2*l-1, 2*l-1) =
218 -ope_var.dG_minus(l)*ope_var.operateurs[conte]->val_sh_one_minus()-
219 ope_var.G_minus(l)*ope_var.operateurs[conte]->der_sh_one_minus() ;
220 systeme.set(2*l-1, 2*l) =
221 -ope_var.dG_minus(l)*ope_var.operateurs[conte]->val_sh_two_minus()-
222 ope_var.G_minus(l)*ope_var.operateurs[conte]->der_sh_two_minus() ;
223
224 sec_membre.set(2*l-2) += ope_var.F_minus(l,k,j) +
225 ope_var.G_minus(l) * ope_var.operateurs[conte]->val_sp_minus() ;
226 sec_membre.set(2*l-1) += ope_var.dF_minus(l,k,j) +
227 ope_var.dG_minus(l) * ope_var.operateurs[conte]->val_sp_minus() +
228 ope_var.G_minus(l) * ope_var.operateurs[conte]->der_sp_minus() ;
229
230 // Valeurs en +1 :
231 systeme.set(2*l, 2*l-1) = ope_var.G_plus(l) *
232 ope_var.operateurs[conte]->val_sh_one_plus() ;
233 systeme.set(2*l, 2*l) = ope_var.G_plus(l) *
234 ope_var.operateurs[conte]->val_sh_two_plus() ;
235 systeme.set(2*l+1, 2*l-1) =
236 ope_var.dG_plus(l)*ope_var.operateurs[conte]->val_sh_one_plus()+
237 ope_var.G_plus(l)*ope_var.operateurs[conte]->der_sh_one_plus() ;
238 systeme.set(2*l+1, 2*l) =
239 ope_var.dG_plus(l)*ope_var.operateurs[conte]->val_sh_two_plus()+
240 ope_var.G_plus(l)*ope_var.operateurs[conte]->der_sh_two_plus() ;
241
242 sec_membre.set(2*l) -= ope_var.F_plus(l,k,j) +
243 ope_var.G_plus(l) * ope_var.operateurs[conte]->val_sp_plus();
244 sec_membre.set(2*l+1) -= ope_var.dF_plus(l,k,j) +
245 ope_var.dG_plus(l) * ope_var.operateurs[conte]->val_sp_plus() +
246 ope_var.G_plus(l) * ope_var.operateurs[conte]->der_sp_plus() ;
247 }
248
249 //-------
250 // ZEC :
251 //-------
252 int np_prec = source.get_mg()->get_np(nz-2) ;
253 int nt_prec = source.get_mg()->get_nt(nz-2) ;
254 conte += (np_prec+1)*nt_prec ;
255
256 systeme.set(taille-2, taille-1) = -ope_var.G_minus(nz-1) *
257 ope_var.operateurs[conte]->val_sh_one_minus() ;
258 systeme.set(taille-1, taille-1) =
259 -ope_var.dG_minus(nz-1)*ope_var.operateurs[conte]->val_sh_one_minus()-
260 ope_var.G_minus(nz-1)*ope_var.operateurs[conte]->der_sh_one_minus() ;
261
262 sec_membre.set(taille-2) += ope_var.F_minus(nz-1,k,j) +
263 ope_var.G_minus(nz-1)*ope_var.operateurs[conte]->val_sp_minus() ;
264 sec_membre.set(taille-1) += ope_var.dF_minus(nz-1,k,j) +
265 ope_var.dG_minus(nz-1) * ope_var.operateurs[conte]->val_sp_minus() +
266 ope_var.G_minus(nz-1) * ope_var.operateurs[conte]->der_sp_minus() ;
267
268 // On resout le systeme ...
269 if (taille > 2)
270 systeme.set_band(2,2) ;
271 else
272 systeme.set_band(1,1) ;
273
274 systeme.set_lu() ;
275 Tbl facteur (systeme.inverse(sec_membre)) ;
276
277 // On range tout ca :
278 // Noyau
279 nr = source.get_mg()->get_nr(0) ;
280 for (int i=0 ; i<nr ; i++)
281 resultat.set(0,k,j,i) = solution_part(0,k,j,i) +facteur(0)*solution_hom_un(0,k,j,i) ;
282
283 // Shells
284 for (int l=1 ; l<nz-1 ; l++) {
285 nr = source.get_mg()->get_nr(l) ;
286 for (int i=0 ; i<nr ; i++)
287 resultat.set(l,k,j,i) = solution_part(l,k,j,i) +
288 facteur(2*l-1)*solution_hom_un(l,k,j,i) +
289 facteur(2*l)*solution_hom_deux(l,k,j,i) ;
290 }
291
292 // Zec
293 nr = source.get_mg()->get_nr(nz-1) ;
294 for (int i=0 ; i<nr ; i++)
295 resultat.set(nz-1,k,j,i) = solution_part(nz-1,k,j,i) +
296 facteur(taille-1)*solution_hom_un(nz-1,k,j,i) ;
297 }
298 start ++ ;
299 }
300
301 return resultat;
302}
303
304}
Matrix handling.
Definition matrice.h:152
int get_nzone() const
Returns the number of domains.
Definition grilles.h:465
Coefficients storage for the multi-domain spectral method.
Definition mtbl_cf.h:196
const Mg3d * get_mg() const
Returns the Mg3d on which the Mtbl_cf is defined.
Definition mtbl_cf.h:463
This class contains the parameters needed to call the general elliptic solver.
Basic array class.
Definition tbl.h:161
Lorene prototypes.
Definition app_hor.h:67