74 #include "param_elliptic.h" 87 Mtbl_cf elliptic_solver (
const Param_elliptic& ope_var,
const Mtbl_cf& source) {
91 assert (source.get_mg()->get_type_r(0) == RARE) ;
92 assert (source.get_mg()->get_type_r(nz-1) == UNSURR) ;
93 for (
int l=1 ; l<nz-1 ; l++)
94 assert(source.get_mg()->get_type_r(l) == FIN) ;
106 Mtbl_cf solution_part(source.get_mg(), source.base) ;
107 Mtbl_cf solution_hom_un(source.get_mg(), source.base) ;
108 Mtbl_cf solution_hom_deux(source.get_mg(), source.base) ;
109 Mtbl_cf resultat(source.get_mg(), source.base) ;
111 solution_part.annule_hard() ;
112 solution_hom_un.annule_hard() ;
113 solution_hom_deux.annule_hard() ;
114 resultat.annule_hard() ;
118 for (
int zone=0 ; zone<nz ; zone++) {
120 nr = source.get_mg()->get_nr(zone) ;
121 nt = source.get_mg()->get_nt(zone) ;
122 np = source.get_mg()->get_np(zone) ;
124 for (
int k=0 ; k<np+1 ; k++)
125 for (
int j=0 ; j<nt ; j++) {
127 if (ope_var.operateurs[conte] != 0x0) {
130 sol_hom =
new Tbl(ope_var.operateurs[conte]->get_solh()) ;
135 for (
int i=0 ; i<nr ; i++)
136 so->set(i) = source(zone, k, j, i) ;
138 sol_part =
new Tbl(ope_var.operateurs[conte]->get_solp(*so)) ;
141 for (
int i=0 ; i<nr ; i++) {
142 solution_part.set(zone, k, j, i) = (*sol_part)(i) ;
143 if (sol_hom->get_ndim()==1)
144 solution_hom_un.set(zone, k, j, i) = (*sol_hom)(i) ;
147 solution_hom_un.set(zone, k, j, i) = (*sol_hom)(0,i) ;
148 solution_hom_deux.set(zone, k, j, i) = (*sol_hom)(1,i) ;
170 for (
int k=0 ; k<source.get_mg()->get_np(0)+1 ; k++)
171 for (
int j=0 ; j<source.get_mg()->get_nt(0) ; j++) {
172 if (ope_var.operateurs[start] != 0x0) {
174 int taille = 2*nz - 2 ;
175 Matrice systeme (taille, taille) ;
176 systeme.set_etat_qcq() ;
177 for (
int i=0 ; i<taille ; i++)
178 for (
int j2=0 ; j2<taille ; j2++)
179 systeme.set(i,j2) = 0 ;
180 Tbl sec_membre (taille) ;
181 sec_membre.set_etat_qcq() ;
182 for (
int i=0 ; i<taille ; i++)
183 sec_membre.set(i) = 0 ;
190 systeme.set(0,0) = ope_var.G_plus(0) *
191 ope_var.operateurs[conte]->val_sh_one_plus() ;
193 ope_var.dG_plus(0) * ope_var.operateurs[conte]->val_sh_one_plus() +
194 ope_var.G_plus(0) * ope_var.operateurs[conte]->der_sh_one_plus() ;
196 sec_membre.set(0) -= ope_var.F_plus(0,k,j) +
197 ope_var.G_plus(0) * ope_var.operateurs[conte]->val_sp_plus() ;
198 sec_membre.set(1) -= ope_var.dF_plus(0,k,j) +
199 ope_var.dG_plus(0) * ope_var.operateurs[conte]->val_sp_plus() +
200 ope_var.G_plus(0) * ope_var.operateurs[conte]->der_sp_plus() ;
206 for (
int l=1 ; l<nz-1 ; l++) {
209 int np_prec = source.get_mg()->get_np(l-1) ;
210 int nt_prec = source.get_mg()->get_nt(l-1) ;
211 conte += (np_prec+1)*nt_prec ;
213 systeme.set(2*l-2, 2*l-1) = -ope_var.G_minus(l) *
214 ope_var.operateurs[conte]->val_sh_one_minus() ;
215 systeme.set(2*l-2, 2*l) = - ope_var.G_minus(l) *
216 ope_var.operateurs[conte]->val_sh_two_minus() ;
217 systeme.set(2*l-1, 2*l-1) =
218 -ope_var.dG_minus(l)*ope_var.operateurs[conte]->val_sh_one_minus()-
219 ope_var.G_minus(l)*ope_var.operateurs[conte]->der_sh_one_minus() ;
220 systeme.set(2*l-1, 2*l) =
221 -ope_var.dG_minus(l)*ope_var.operateurs[conte]->val_sh_two_minus()-
222 ope_var.G_minus(l)*ope_var.operateurs[conte]->der_sh_two_minus() ;
224 sec_membre.set(2*l-2) += ope_var.F_minus(l,k,j) +
225 ope_var.G_minus(l) * ope_var.operateurs[conte]->val_sp_minus() ;
226 sec_membre.set(2*l-1) += ope_var.dF_minus(l,k,j) +
227 ope_var.dG_minus(l) * ope_var.operateurs[conte]->val_sp_minus() +
228 ope_var.G_minus(l) * ope_var.operateurs[conte]->der_sp_minus() ;
231 systeme.set(2*l, 2*l-1) = ope_var.G_plus(l) *
232 ope_var.operateurs[conte]->val_sh_one_plus() ;
233 systeme.set(2*l, 2*l) = ope_var.G_plus(l) *
234 ope_var.operateurs[conte]->val_sh_two_plus() ;
235 systeme.set(2*l+1, 2*l-1) =
236 ope_var.dG_plus(l)*ope_var.operateurs[conte]->val_sh_one_plus()+
237 ope_var.G_plus(l)*ope_var.operateurs[conte]->der_sh_one_plus() ;
238 systeme.set(2*l+1, 2*l) =
239 ope_var.dG_plus(l)*ope_var.operateurs[conte]->val_sh_two_plus()+
240 ope_var.G_plus(l)*ope_var.operateurs[conte]->der_sh_two_plus() ;
242 sec_membre.set(2*l) -= ope_var.F_plus(l,k,j) +
243 ope_var.G_plus(l) * ope_var.operateurs[conte]->val_sp_plus();
244 sec_membre.set(2*l+1) -= ope_var.dF_plus(l,k,j) +
245 ope_var.dG_plus(l) * ope_var.operateurs[conte]->val_sp_plus() +
246 ope_var.G_plus(l) * ope_var.operateurs[conte]->der_sp_plus() ;
252 int np_prec = source.get_mg()->get_np(nz-2) ;
253 int nt_prec = source.get_mg()->get_nt(nz-2) ;
254 conte += (np_prec+1)*nt_prec ;
256 systeme.set(taille-2, taille-1) = -ope_var.G_minus(nz-1) *
257 ope_var.operateurs[conte]->val_sh_one_minus() ;
258 systeme.set(taille-1, taille-1) =
259 -ope_var.dG_minus(nz-1)*ope_var.operateurs[conte]->val_sh_one_minus()-
260 ope_var.G_minus(nz-1)*ope_var.operateurs[conte]->der_sh_one_minus() ;
262 sec_membre.set(taille-2) += ope_var.F_minus(nz-1,k,j) +
263 ope_var.G_minus(nz-1)*ope_var.operateurs[conte]->val_sp_minus() ;
264 sec_membre.set(taille-1) += ope_var.dF_minus(nz-1,k,j) +
265 ope_var.dG_minus(nz-1) * ope_var.operateurs[conte]->val_sp_minus() +
266 ope_var.G_minus(nz-1) * ope_var.operateurs[conte]->der_sp_minus() ;
270 systeme.set_band(2,2) ;
272 systeme.set_band(1,1) ;
275 Tbl facteur (systeme.inverse(sec_membre)) ;
279 nr = source.get_mg()->get_nr(0) ;
280 for (
int i=0 ; i<nr ; i++)
281 resultat.set(0,k,j,i) = solution_part(0,k,j,i) +facteur(0)*solution_hom_un(0,k,j,i) ;
284 for (
int l=1 ; l<nz-1 ; l++) {
285 nr = source.get_mg()->get_nr(l) ;
286 for (
int i=0 ; i<nr ; i++)
287 resultat.set(l,k,j,i) = solution_part(l,k,j,i) +
288 facteur(2*l-1)*solution_hom_un(l,k,j,i) +
289 facteur(2*l)*solution_hom_deux(l,k,j,i) ;
293 nr = source.get_mg()->get_nr(nz-1) ;
294 for (
int i=0 ; i<nr ; i++)
295 resultat.set(nz-1,k,j,i) = solution_part(nz-1,k,j,i) +
296 facteur(taille-1)*solution_hom_un(nz-1,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.