68 #include "param_elliptic.h" 76 Mtbl_cf elliptic_solver_no_zec (
const Param_elliptic& ope_var,
const Mtbl_cf& source,
double valeur) {
80 assert (source.get_mg()->get_type_r(0) == RARE) ;
81 assert (source.get_mg()->get_type_r(nz-1) == UNSURR) ;
82 for (
int l=1 ; l<nz-1 ; l++)
83 assert(source.get_mg()->get_type_r(l) == FIN) ;
95 Mtbl_cf solution_part(source.get_mg(), source.base) ;
96 Mtbl_cf solution_hom_un(source.get_mg(), source.base) ;
97 Mtbl_cf solution_hom_deux(source.get_mg(), source.base) ;
98 Mtbl_cf resultat(source.get_mg(), source.base) ;
100 solution_part.annule_hard() ;
101 solution_hom_un.annule_hard() ;
102 solution_hom_deux.annule_hard() ;
103 resultat.annule_hard() ;
107 for (
int zone=0 ; zone<nz-1 ; zone++) {
108 nr = source.get_mg()->get_nr(zone) ;
109 nt = source.get_mg()->get_nt(zone) ;
110 np = source.get_mg()->get_np(zone) ;
112 for (
int k=0 ; k<np+1 ; k++)
113 for (
int j=0 ; j<nt ; j++) {
114 if (ope_var.operateurs[conte] != 0x0) {
116 sol_hom =
new Tbl(ope_var.operateurs[conte]->get_solh()) ;
121 for (
int i=0 ; i<nr ; i++)
122 so->set(i) = source(zone, k, j, i) ;
124 sol_part =
new Tbl(ope_var.operateurs[conte]->get_solp(*so)) ;
127 for (
int i=0 ; i<nr ; i++) {
128 solution_part.set(zone, k, j, i) = (*sol_part)(i) ;
129 if (sol_hom->get_ndim()==1)
130 solution_hom_un.set(zone, k, j, i) = (*sol_hom)(i) ;
133 solution_hom_un.set(zone, k, j, i) = (*sol_hom)(0,i) ;
134 solution_hom_deux.set(zone, k, j, i) = (*sol_hom)(1,i) ;
156 for (
int k=0 ; k<1 ; k++)
157 for (
int j=0 ; j<1 ; j++) {
158 if (ope_var.operateurs[start] != 0x0) {
160 int taille = 2*nz - 3 ;
161 Matrice systeme (taille, taille) ;
162 systeme.set_etat_qcq() ;
163 for (
int i=0 ; i<taille ; i++)
164 for (
int j2=0 ; j2<taille ; j2++)
165 systeme.set(i,j2) = 0 ;
166 Tbl sec_membre (taille) ;
167 sec_membre.set_etat_qcq() ;
168 for (
int i=0 ; i<taille ; i++)
169 sec_membre.set(i) = 0 ;
176 systeme.set(0,0) = ope_var.G_plus(0) *
177 ope_var.operateurs[conte]->val_sh_one_plus() ;
179 ope_var.dG_plus(0) * ope_var.operateurs[conte]->val_sh_one_plus() +
180 ope_var.G_plus(0) * ope_var.operateurs[conte]->der_sh_one_plus() ;
182 sec_membre.set(0) -= ope_var.F_plus(0,k,j) +
183 ope_var.G_plus(0) * ope_var.operateurs[conte]->val_sp_plus() ;
184 sec_membre.set(1) -= ope_var.dF_plus(0,k,j) +
185 ope_var.dG_plus(0) * ope_var.operateurs[conte]->val_sp_plus() +
186 ope_var.G_plus(0) * ope_var.operateurs[conte]->der_sp_plus() ;
193 for (
int l=1 ; l<nz-1 ; l++) {
196 int np_prec = source.get_mg()->get_np(l-1) ;
197 int nt_prec = source.get_mg()->get_nt(l-1) ;
198 conte += (np_prec+1)*nt_prec ;
200 systeme.set(2*l-2, 2*l-1) = -ope_var.G_minus(l) *
201 ope_var.operateurs[conte]->val_sh_one_minus() ;
202 systeme.set(2*l-2, 2*l) = - ope_var.G_minus(l) *
203 ope_var.operateurs[conte]->val_sh_two_minus() ;
204 systeme.set(2*l-1, 2*l-1) =
205 -ope_var.dG_minus(l)*ope_var.operateurs[conte]->val_sh_one_minus()-
206 ope_var.G_minus(l)*ope_var.operateurs[conte]->der_sh_one_minus() ;
207 systeme.set(2*l-1, 2*l) =
208 -ope_var.dG_minus(l)*ope_var.operateurs[conte]->val_sh_two_minus()-
209 ope_var.G_minus(l)*ope_var.operateurs[conte]->der_sh_two_minus() ;
211 sec_membre.set(2*l-2) += ope_var.F_minus(l,k,j) +
212 ope_var.G_minus(l) * ope_var.operateurs[conte]->val_sp_minus() ;
213 sec_membre.set(2*l-1) += ope_var.dF_minus(l,k,j) +
214 ope_var.dG_minus(l) * ope_var.operateurs[conte]->val_sp_minus() +
215 ope_var.G_minus(l) * ope_var.operateurs[conte]->der_sp_minus() ;
218 systeme.set(2*l, 2*l-1) = ope_var.G_plus(l) *
219 ope_var.operateurs[conte]->val_sh_one_plus() ;
220 systeme.set(2*l, 2*l) = ope_var.G_plus(l) *
221 ope_var.operateurs[conte]->val_sh_two_plus() ;
223 systeme.set(2*l+1, 2*l-1) =
224 ope_var.dG_plus(l)*ope_var.operateurs[conte]->val_sh_one_plus()+
225 ope_var.G_plus(l)*ope_var.operateurs[conte]->der_sh_one_plus() ;
226 systeme.set(2*l+1, 2*l) =
227 ope_var.dG_plus(l)*ope_var.operateurs[conte]->val_sh_two_plus()+
228 ope_var.G_plus(l)*ope_var.operateurs[conte]->der_sh_two_plus() ;
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();
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() ;
240 sec_membre.set(2*l) -= -valeur + ope_var.F_plus(l,k,j) +
241 ope_var.G_plus(l) * ope_var.operateurs[conte]->val_sp_plus();
246 systeme.set_band(2,2) ;
248 systeme.set_band(1,1) ;
251 Tbl facteur (systeme.inverse(sec_membre)) ;
255 nr = source.get_mg()->get_nr(0) ;
256 for (
int i=0 ; i<nr ; i++)
257 resultat.set(0,k,j,i) = solution_part(0,k,j,i)
258 +facteur(0)*solution_hom_un(0,k,j,i) ;
261 for (
int l=1 ; l<nz-1 ; l++) {
262 nr = source.get_mg()->get_nr(l) ;
263 for (
int i=0 ; i<nr ; i++)
264 resultat.set(l,k,j,i) = solution_part(l,k,j,i) +
265 facteur(2*l-1)*solution_hom_un(l,k,j,i) +
266 facteur(2*l)*solution_hom_deux(l,k,j,i) ;
const Mg3d * get_mg() const
Returns the Mg3d on which the Mtbl_cf is defined.
int get_nzone() const
Returns the number of domains.