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