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