00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047
00048
00049
00050
00051 char sol_elliptic_sin_zec_C[] = "$Header: /cvsroot/Lorene/C++/Source/Non_class_members/PDE/sol_elliptic_sin_zec.C,v 1.6 2007/05/08 07:08:30 p_grandclement Exp $" ;
00052
00053
00054 #include <stdlib.h>
00055 #include <math.h>
00056 #include <gsl/gsl_sf_bessel.h>
00057
00058
00059 #include "tbl.h"
00060 #include "mtbl_cf.h"
00061 #include "map.h"
00062 #include "param_elliptic.h"
00063
00064
00065
00066
00067
00068
00069 Mtbl_cf elliptic_solver_sin_zec (const Param_elliptic& ope_var,
00070 const Mtbl_cf& source, double* amplis, double* phases) {
00071
00072
00073 int nz = source.get_mg()->get_nzone() ;
00074 assert (nz>1) ;
00075 assert (source.get_mg()->get_type_r(0) == RARE) ;
00076 assert (source.get_mg()->get_type_r(nz-1) == UNSURR) ;
00077 for (int l=1 ; l<nz-1 ; l++)
00078 assert(source.get_mg()->get_type_r(l) == FIN) ;
00079
00080
00081 int nr, nt, np ;
00082
00083
00084 Tbl *so ;
00085 Tbl *sol_hom ;
00086 Tbl *sol_part ;
00087
00088
00089
00090 Mtbl_cf solution_part(source.get_mg(), source.base) ;
00091 Mtbl_cf solution_hom_un(source.get_mg(), source.base) ;
00092 Mtbl_cf solution_hom_deux(source.get_mg(), source.base) ;
00093 Mtbl_cf resultat(source.get_mg(), source.base) ;
00094
00095 solution_part.annule_hard() ;
00096 solution_hom_un.annule_hard() ;
00097 solution_hom_deux.annule_hard() ;
00098 resultat.annule_hard() ;
00099
00100
00101 int conte = 0 ;
00102 for (int zone=0 ; zone<nz-1 ; zone++) {
00103 nr = source.get_mg()->get_nr(zone) ;
00104 nt = source.get_mg()->get_nt(zone) ;
00105 np = source.get_mg()->get_np(zone) ;
00106
00107 for (int k=0 ; k<np+1 ; k++)
00108 for (int j=0 ; j<nt ; j++) {
00109 if (ope_var.operateurs[conte] != 0x0) {
00110
00111 sol_hom = new Tbl(ope_var.operateurs[conte]->get_solh()) ;
00112
00113
00114 so = new Tbl(nr) ;
00115 so->set_etat_qcq() ;
00116 for (int i=0 ; i<nr ; i++)
00117 so->set(i) = source(zone, k, j, i) ;
00118
00119 sol_part = new Tbl(ope_var.operateurs[conte]->get_solp(*so)) ;
00120
00121
00122 for (int i=0 ; i<nr ; i++) {
00123 solution_part.set(zone, k, j, i) = (*sol_part)(i) ;
00124 if (sol_hom->get_ndim()==1)
00125 solution_hom_un.set(zone, k, j, i) = (*sol_hom)(i) ;
00126 else
00127 {
00128 solution_hom_un.set(zone, k, j, i) = (*sol_hom)(0,i) ;
00129 solution_hom_deux.set(zone, k, j, i) = (*sol_hom)(1,i) ;
00130 }
00131 }
00132
00133 delete so ;
00134 delete sol_hom ;
00135 delete sol_part ;
00136
00137 }
00138 conte ++ ;
00139 }
00140 }
00141
00142
00143
00144
00145
00146
00147
00148
00149 int start = 0 ;
00150 for (int k=0 ; k< source.get_mg()->get_np(0)+1; k++)
00151 for (int j=0 ; j<source.get_mg()->get_nt(0) ; j++) {
00152 if (ope_var.operateurs[start] != 0x0) {
00153
00154 int taille = 2*nz - 2 ;
00155 Matrice systeme (taille, taille) ;
00156 systeme.set_etat_qcq() ;
00157 for (int i=0 ; i<taille ; i++)
00158 for (int j2=0 ; j2<taille ; j2++)
00159 systeme.set(i,j2) = 0 ;
00160 Tbl sec_membre (taille) ;
00161 sec_membre.set_etat_qcq() ;
00162 for (int i=0 ; i<taille ; i++)
00163 sec_membre.set(i) = 0 ;
00164
00165
00166
00167 conte = start ;
00168
00169 systeme.set(0,0) = ope_var.G_plus(0) *
00170 ope_var.operateurs[conte]->val_sh_one_plus() ;
00171 systeme.set(1,0) =
00172 ope_var.dG_plus(0) * ope_var.operateurs[conte]->val_sh_one_plus() +
00173 ope_var.G_plus(0) * ope_var.operateurs[conte]->der_sh_one_plus() ;
00174
00175 sec_membre.set(0) -= ope_var.F_plus(0,k,j) +
00176 ope_var.G_plus(0) * ope_var.operateurs[conte]->val_sp_plus() ;
00177 sec_membre.set(1) -= ope_var.dF_plus(0,k,j) +
00178 ope_var.dG_plus(0) * ope_var.operateurs[conte]->val_sp_plus() +
00179 ope_var.G_plus(0) * ope_var.operateurs[conte]->der_sp_plus() ;
00180
00181
00182
00183
00184
00185 for (int l=1 ; l<nz-1 ; l++) {
00186
00187
00188 int np_prec = source.get_mg()->get_np(l-1) ;
00189 int nt_prec = source.get_mg()->get_nt(l-1) ;
00190 conte += (np_prec+1)*nt_prec ;
00191
00192 systeme.set(2*l-2, 2*l-1) = -ope_var.G_minus(l) *
00193 ope_var.operateurs[conte]->val_sh_one_minus() ;
00194 systeme.set(2*l-2, 2*l) = - ope_var.G_minus(l) *
00195 ope_var.operateurs[conte]->val_sh_two_minus() ;
00196 systeme.set(2*l-1, 2*l-1) =
00197 -ope_var.dG_minus(l)*ope_var.operateurs[conte]->val_sh_one_minus()-
00198 ope_var.G_minus(l)*ope_var.operateurs[conte]->der_sh_one_minus() ;
00199 systeme.set(2*l-1, 2*l) =
00200 -ope_var.dG_minus(l)*ope_var.operateurs[conte]->val_sh_two_minus()-
00201 ope_var.G_minus(l)*ope_var.operateurs[conte]->der_sh_two_minus() ;
00202
00203 sec_membre.set(2*l-2) += ope_var.F_minus(l,k,j) +
00204 ope_var.G_minus(l) * ope_var.operateurs[conte]->val_sp_minus() ;
00205 sec_membre.set(2*l-1) += ope_var.dF_minus(l,k,j) +
00206 ope_var.dG_minus(l) * ope_var.operateurs[conte]->val_sp_minus() +
00207 ope_var.G_minus(l) * ope_var.operateurs[conte]->der_sp_minus() ;
00208
00209
00210 systeme.set(2*l, 2*l-1) = ope_var.G_plus(l) *
00211 ope_var.operateurs[conte]->val_sh_one_plus() ;
00212 systeme.set(2*l, 2*l) = ope_var.G_plus(l) *
00213 ope_var.operateurs[conte]->val_sh_two_plus() ;
00214
00215 systeme.set(2*l+1, 2*l-1) =
00216 ope_var.dG_plus(l)*ope_var.operateurs[conte]->val_sh_one_plus()+
00217 ope_var.G_plus(l)*ope_var.operateurs[conte]->der_sh_one_plus() ;
00218 systeme.set(2*l+1, 2*l) =
00219 ope_var.dG_plus(l)*ope_var.operateurs[conte]->val_sh_two_plus()+
00220 ope_var.G_plus(l)*ope_var.operateurs[conte]->der_sh_two_plus() ;
00221
00222 sec_membre.set(2*l) -= ope_var.F_plus(l,k,j) +
00223 ope_var.G_plus(l) * ope_var.operateurs[conte]->val_sp_plus();
00224
00225 sec_membre.set(2*l+1) -= ope_var.dF_plus(l,k,j) +
00226 ope_var.dG_plus(l) * ope_var.operateurs[conte]->val_sp_plus() +
00227 ope_var.G_plus(l) * ope_var.operateurs[conte]->der_sp_plus() ;
00228
00229 }
00230
00231
00232 double val_sh = cos(phases[start])*ope_var.operateurs[conte]->val_sh_one_plus()
00233 + sin(phases[start])*ope_var.operateurs[conte]->val_sh_two_plus() ;
00234 double der_sh = cos(phases[start])*ope_var.operateurs[conte]->der_sh_one_plus()
00235 + sin(phases[start])*ope_var.operateurs[conte]->der_sh_two_plus() ;
00236
00237 systeme.set(taille-2, taille-1) = -ope_var.G_minus(nz-1) * val_sh ;
00238 systeme.set(taille-1, taille-1) = -ope_var.dG_minus(nz-1)*val_sh- ope_var.G_minus(nz-1)*der_sh ;
00239
00240 sec_membre.set(taille-2) += ope_var.F_minus(nz-1,k,j) ;
00241 sec_membre.set(taille-1) += ope_var.dF_minus(nz-1,k,j) ;
00242
00243
00244 if (taille > 2)
00245 systeme.set_band(2,2) ;
00246 else
00247 systeme.set_band(1,1) ;
00248
00249 systeme.set_lu() ;
00250 Tbl facteur (systeme.inverse(sec_membre)) ;
00251
00252 amplis[start] = facteur(taille-1) ;
00253
00254
00255
00256 nr = source.get_mg()->get_nr(0) ;
00257 for (int i=0 ; i<nr ; i++)
00258 resultat.set(0,k,j,i) = solution_part(0,k,j,i)
00259 +facteur(0)*solution_hom_un(0,k,j,i) ;
00260
00261
00262 for (int l=1 ; l<nz-1 ; l++) {
00263 nr = source.get_mg()->get_nr(l) ;
00264 for (int i=0 ; i<nr ; i++)
00265 resultat.set(l,k,j,i) = solution_part(l,k,j,i) +
00266 facteur(2*l-1)*solution_hom_un(l,k,j,i) +
00267 facteur(2*l)*solution_hom_deux(l,k,j,i) ;
00268 }
00269 }
00270 start ++ ;
00271 }
00272 return resultat;
00273 }
00274
00275