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