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