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