LORENE
sol_elliptic_boundary.C
1/*
2 * Copyright (c) 2003 Philippe Grandclement
3 * Jose Luis Jaramillo
4 * Francois Limousin
5 *
6 * This file is part of LORENE.
7 *
8 * LORENE is free software; you can redistribute it and/or modify
9 * it under the terms of the GNU General Public License version 2
10 * as published by the Free Software Foundation.
11 *
12 * LORENE is distributed in the hope that it will be useful,
13 * but WITHOUT ANY WARRANTY; without even the implied warranty of
14 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15 * GNU General Public License for more details.
16 *
17 * You should have received a copy of the GNU General Public License
18 * along with LORENE; if not, write to the Free Software
19 * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
20 *
21 */
22
23
24
25/*
26 * $Id: sol_elliptic_boundary.C,v 1.5 2016/12/05 16:18:10 j_novak Exp $
27 * $Log: sol_elliptic_boundary.C,v $
28 * Revision 1.5 2016/12/05 16:18:10 j_novak
29 * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
30 *
31 * Revision 1.4 2014/10/13 08:53:30 j_novak
32 * Lorene classes and functions now belong to the namespace Lorene.
33 *
34 * Revision 1.3 2014/10/06 15:16:10 j_novak
35 * Modified #include directives to use c++ syntax.
36 *
37 * Revision 1.2 2008/08/20 15:03:55 n_vasset
38 * Correction on how the boundary condition is imposed
39 *
40 * Revision 1.1 2005/06/09 08:01:05 f_limousin
41 * Implement a new function sol_elliptic_boundary() and
42 * Vector::poisson_boundary(...) which solve the vectorial poisson
43 * equation (method 6) with an inner boundary condition.
44 *
45 *
46 * $Header: /cvsroot/Lorene/C++/Source/Non_class_members/PDE/sol_elliptic_boundary.C,v 1.5 2016/12/05 16:18:10 j_novak Exp $
47 *
48 */
49
50// Header C :
51#include <cstdlib>
52#include <cmath>
53
54// Headers Lorene :
55#include "param_elliptic.h"
56#include "tbl.h"
57#include "mtbl_cf.h"
58#include "map.h"
59
60
61 //----------------------------------------------
62 // Version Mtbl_cf
63 //----------------------------------------------
64
65
66
67namespace Lorene {
68Mtbl_cf elliptic_solver_boundary (const Param_elliptic& ope_var, const Mtbl_cf& source,
69 const Mtbl_cf& bound, double fact_dir, double fact_neu ) {
70 // Verifications d'usage sur les zones
71 int nz = source.get_mg()->get_nzone() ;
72 assert (nz>1) ;
73 assert (source.get_mg()->get_type_r(0) == RARE) ;
74 assert (source.get_mg()->get_type_r(nz-1) == UNSURR) ;
75 for (int l=1 ; l<nz-1 ; l++)
76 assert(source.get_mg()->get_type_r(l) == FIN) ;
77
78 // donnees sur la zone
79 int nr, nt, np ;
80
81 //Rangement des valeurs intermediaires
82 Tbl *so ;
83 Tbl *sol_hom ;
84 Tbl *sol_part ;
85
86
87 // Rangement des solutions, avant raccordement
88 Mtbl_cf solution_part(source.get_mg(), source.base) ;
89 Mtbl_cf solution_hom_un(source.get_mg(), source.base) ;
90 Mtbl_cf solution_hom_deux(source.get_mg(), source.base) ;
91 Mtbl_cf resultat(source.get_mg(), source.base) ;
92
93 solution_part.annule_hard() ;
94 solution_hom_un.annule_hard() ;
95 solution_hom_deux.annule_hard() ;
96 resultat.annule_hard() ;
97
98 // Computation of the SP and SH's in every domain ...
99 int conte = 0 ;
100 for (int zone=0 ; zone<nz ; zone++) {
101
102 nr = source.get_mg()->get_nr(zone) ;
103 nt = source.get_mg()->get_nt(zone) ;
104 np = source.get_mg()->get_np(zone) ;
105
106 for (int k=0 ; k<np+1 ; k++)
107 for (int j=0 ; j<nt ; j++) {
108
109 if (ope_var.operateurs[conte] != 0x0) {
110
111 // Calcul de la SH
112 sol_hom = new Tbl(ope_var.operateurs[conte]->get_solh()) ;
113
114 //Calcul de la SP
115 so = new Tbl(nr) ;
116 so->set_etat_qcq() ;
117 for (int i=0 ; i<nr ; i++)
118 so->set(i) = source(zone, k, j, i) ;
119
120 sol_part = new Tbl(ope_var.operateurs[conte]->get_solp(*so)) ;
121
122 // Rangement dans les tableaux globaux ;
123 for (int i=0 ; i<nr ; i++) {
124 solution_part.set(zone, k, j, i) = (*sol_part)(i) ;
125 if (sol_hom->get_ndim()==1)
126 solution_hom_un.set(zone, k, j, i) = (*sol_hom)(i) ;
127 else
128 {
129 solution_hom_un.set(zone, k, j, i) = (*sol_hom)(0,i) ;
130 solution_hom_deux.set(zone, k, j, i) = (*sol_hom)(1,i) ;
131 }
132 }
133
134 delete so ;
135 delete sol_hom ;
136 delete sol_part ;
137
138 }
139 conte ++ ;
140 }
141 }
142
143
144 //-------------------------------------------------
145 // ON EST PARTI POUR LE RACCORD (Be careful ....)
146 //-------------------------------------------------
147
148 // C'est pas simple toute cette sombre affaire...
149 // Que de cas meme nombre de points dans chaque domaines...
150
151 int start = 0 ;
152 for (int k=0 ; k<source.get_mg()->get_np(0)+1 ; k++)
153 for (int j=0 ; j<source.get_mg()->get_nt(0) ; j++) {
154 if (ope_var.operateurs[start] != 0x0) {
155
156 int taille = 2*nz - 3 ;
157 Matrice systeme (taille, taille) ;
158 systeme.set_etat_qcq() ;
159 for (int i=0 ; i<taille ; i++)
160 for (int j2=0 ; j2<taille ; j2++)
161 systeme.set(i,j2) = 0 ;
162 Tbl sec_membre (taille) ;
163 sec_membre.set_etat_qcq() ;
164 for (int i=0 ; i<taille ; i++)
165 sec_membre.set(i) = 0 ;
166
167 //----------
168 // BOUNDARY :
169 //----------
170 conte = start ;
171
172 // Boundary value at an angular point
173
174 // Setting the right hand side term
175
176 sec_membre.set(0) -= bound.val_in_bound_jk(1, j, k)/sqrt(2.) ; //Relevant value is in the first shell and only in the first coefficient
177
178
179 //--------------
180 // FIRST SHELL :
181 //--------------
182
183 int l_1=1 ;
184
185 // On se met au bon endroit :
186 int np_prec_1 = source.get_mg()->get_np(l_1-1) ;
187 int nt_prec_1 = source.get_mg()->get_nt(l_1-1) ;
188 conte += (np_prec_1+1)*nt_prec_1 ;
189
190 systeme.set(0, 0) = fact_dir * (-ope_var.G_minus(l_1) *
191 ope_var.operateurs[conte]->val_sh_one_minus() )
192 + fact_neu *
193 ( -ope_var.dG_minus(l_1)*ope_var.operateurs[conte]->val_sh_one_minus()-
194 ope_var.G_minus(l_1)*ope_var.operateurs[conte]->der_sh_one_minus() );
195
196 systeme.set(0, 1) = fact_dir * (- ope_var.G_minus(l_1) *
197 ope_var.operateurs[conte]->val_sh_two_minus() ) + fact_neu *
198 (-ope_var.dG_minus(l_1)*ope_var.operateurs[conte]->val_sh_two_minus()-
199 ope_var.G_minus(l_1)*ope_var.operateurs[conte]->der_sh_two_minus() ) ;
200
201 //Completing the right hand side term
202 sec_membre.set(0) += fact_dir * (ope_var.F_minus(l_1,k,j) +
203 ope_var.G_minus(l_1) * ope_var.operateurs[conte]->val_sp_minus() ) +
204 fact_neu * ( ope_var.dF_minus(l_1,k,j) +
205 ope_var.dG_minus(l_1) * ope_var.operateurs[conte]->val_sp_minus() +
206 ope_var.G_minus(l_1) * ope_var.operateurs[conte]->der_sp_minus() ) ;
207
208
209 // Valeurs en +1 :
210 systeme.set(2*l_1-1, 2*l_1-2) = ope_var.G_plus(l_1) *
211 ope_var.operateurs[conte]->val_sh_one_plus() ;
212 systeme.set(2*l_1-1, 2*l_1-1) = ope_var.G_plus(l_1) *
213 ope_var.operateurs[conte]->val_sh_two_plus() ;
214 systeme.set(2*l_1, 2*l_1-2) =
215 ope_var.dG_plus(l_1)*ope_var.operateurs[conte]->val_sh_one_plus()+
216 ope_var.G_plus(l_1)*ope_var.operateurs[conte]->der_sh_one_plus() ;
217 systeme.set(2*l_1, 2*l_1-1) =
218 ope_var.dG_plus(l_1)*ope_var.operateurs[conte]->val_sh_two_plus()+
219 ope_var.G_plus(l_1)*ope_var.operateurs[conte]->der_sh_two_plus() ;
220
221 sec_membre.set(2*l_1-1) -= ope_var.F_plus(l_1,k,j) +
222 ope_var.G_plus(l_1) * ope_var.operateurs[conte]->val_sp_plus();
223 sec_membre.set(2*l_1) -= ope_var.dF_plus(l_1,k,j) +
224 ope_var.dG_plus(l_1) * ope_var.operateurs[conte]->val_sp_plus() +
225 ope_var.G_plus(l_1) * ope_var.operateurs[conte]->der_sp_plus() ;
226
227
228
229 //----------
230 // SHELLS :
231 //----------
232// assert (nz-1 > 2) ; // At least two shells
233 for (int l=2 ; l<nz-1 ; l++) {
234
235 // On se met au bon endroit :
236 int np_prec = source.get_mg()->get_np(l-1) ;
237 int nt_prec = source.get_mg()->get_nt(l-1) ;
238 conte += (np_prec+1)*nt_prec ;
239
240 systeme.set(2*l-3, 2*l-2) = -ope_var.G_minus(l) *
241 ope_var.operateurs[conte]->val_sh_one_minus() ;
242 systeme.set(2*l-3, 2*l-1) = - ope_var.G_minus(l) *
243 ope_var.operateurs[conte]->val_sh_two_minus() ;
244 systeme.set(2*l-2, 2*l-2) =
245 -ope_var.dG_minus(l)*ope_var.operateurs[conte]->val_sh_one_minus()-
246 ope_var.G_minus(l)*ope_var.operateurs[conte]->der_sh_one_minus() ;
247 systeme.set(2*l-2, 2*l-1) =
248 -ope_var.dG_minus(l)*ope_var.operateurs[conte]->val_sh_two_minus()-
249 ope_var.G_minus(l)*ope_var.operateurs[conte]->der_sh_two_minus() ;
250
251 sec_membre.set(2*l-3) += ope_var.F_minus(l,k,j) +
252 ope_var.G_minus(l) * ope_var.operateurs[conte]->val_sp_minus() ;
253 sec_membre.set(2*l-2) += ope_var.dF_minus(l,k,j) +
254 ope_var.dG_minus(l) * ope_var.operateurs[conte]->val_sp_minus() +
255 ope_var.G_minus(l) * ope_var.operateurs[conte]->der_sp_minus() ;
256
257 // Valeurs en +1 :
258 systeme.set(2*l-1, 2*l-2) = ope_var.G_plus(l) *
259 ope_var.operateurs[conte]->val_sh_one_plus() ;
260 systeme.set(2*l-1, 2*l-1) = ope_var.G_plus(l) *
261 ope_var.operateurs[conte]->val_sh_two_plus() ;
262 systeme.set(2*l, 2*l-2) =
263 ope_var.dG_plus(l)*ope_var.operateurs[conte]->val_sh_one_plus()+
264 ope_var.G_plus(l)*ope_var.operateurs[conte]->der_sh_one_plus() ;
265 systeme.set(2*l, 2*l-1) =
266 ope_var.dG_plus(l)*ope_var.operateurs[conte]->val_sh_two_plus()+
267 ope_var.G_plus(l)*ope_var.operateurs[conte]->der_sh_two_plus() ;
268
269 sec_membre.set(2*l-1) -= ope_var.F_plus(l,k,j) +
270 ope_var.G_plus(l) * ope_var.operateurs[conte]->val_sp_plus();
271 sec_membre.set(2*l) -= ope_var.dF_plus(l,k,j) +
272 ope_var.dG_plus(l) * ope_var.operateurs[conte]->val_sp_plus() +
273 ope_var.G_plus(l) * ope_var.operateurs[conte]->der_sp_plus() ;
274 }
275
276 //-------
277 // ZEC :
278 //-------
279 int np_prec = source.get_mg()->get_np(nz-2) ;
280 int nt_prec = source.get_mg()->get_nt(nz-2) ;
281 conte += (np_prec+1)*nt_prec ;
282
283 systeme.set(taille-2, taille-1) = -ope_var.G_minus(nz-1) *
284 ope_var.operateurs[conte]->val_sh_one_minus() ;
285 systeme.set(taille-1, taille-1) =
286 -ope_var.dG_minus(nz-1)*ope_var.operateurs[conte]->val_sh_one_minus()-
287 ope_var.G_minus(nz-1)*ope_var.operateurs[conte]->der_sh_one_minus() ;
288
289 sec_membre.set(taille-2) += ope_var.F_minus(nz-1,k,j) +
290 ope_var.G_minus(nz-1)*ope_var.operateurs[conte]->val_sp_minus() ;
291 sec_membre.set(taille-1) += ope_var.dF_minus(nz-1,k,j) +
292 ope_var.dG_minus(nz-1) * ope_var.operateurs[conte]->val_sp_minus() +
293 ope_var.G_minus(nz-1) * ope_var.operateurs[conte]->der_sp_minus() ;
294
295 // On resout le systeme ...
296 if (taille > 2)
297 systeme.set_band(2,2) ;
298 else
299 systeme.set_band(1,1) ;
300
301 systeme.set_lu() ;
302 Tbl facteur (systeme.inverse(sec_membre)) ;
303
304 // On range tout ca :
305
306 // Shells
307 for (int l=1 ; l<nz-1 ; l++) {
308 nr = source.get_mg()->get_nr(l) ;
309 for (int i=0 ; i<nr ; i++)
310 resultat.set(l,k,j,i) = solution_part(l,k,j,i) +
311 facteur(2*l-2)*solution_hom_un(l,k,j,i) +
312 facteur(2*l-1)*solution_hom_deux(l,k,j,i) ;
313 }
314
315 // Zec
316 nr = source.get_mg()->get_nr(nz-1) ;
317 for (int i=0 ; i<nr ; i++)
318 resultat.set(nz-1,k,j,i) = solution_part(nz-1,k,j,i) +
319 facteur(taille-1)*solution_hom_un(nz-1,k,j,i) ;
320 }
321 start ++ ;
322 }
323
324 return resultat;
325}
326
327}
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
Cmp sqrt(const Cmp &)
Square root.
Definition cmp_math.C:223
Lorene prototypes.
Definition app_hor.h:67