21 char sol_elliptic_fixe_der_zeroC[] =
"$Header: /cvsroot/Lorene/C++/Source/Non_class_members/PDE/sol_elliptic_fixe_der_zero.C,v 1.5 2014/10/13 08:53:30 j_novak Exp $" ;
36 #include "param_elliptic.h" 46 Mtbl_cf elliptic_solver_fixe_der_zero (
double valeur,
47 const Param_elliptic& ope_var,
48 const Mtbl_cf& source) {
52 assert (source.get_mg()->get_type_r(0) == RARE) ;
53 assert (source.get_mg()->get_type_r(nz-1) == UNSURR) ;
54 for (
int l=1 ; l<nz-1 ; l++)
55 assert(source.get_mg()->get_type_r(l) == FIN) ;
67 Mtbl_cf solution_part(source.get_mg(), source.base) ;
68 Mtbl_cf solution_hom_un(source.get_mg(), source.base) ;
69 Mtbl_cf solution_hom_deux(source.get_mg(), source.base) ;
70 Mtbl_cf resultat(source.get_mg(), source.base) ;
72 solution_part.annule_hard() ;
73 solution_hom_un.annule_hard() ;
74 solution_hom_deux.annule_hard() ;
75 resultat.annule_hard() ;
79 for (
int zone=0 ; zone<nz ; zone++) {
80 nr = source.get_mg()->get_nr(zone) ;
81 nt = source.get_mg()->get_nt(zone) ;
82 np = source.get_mg()->get_np(zone) ;
84 for (
int k=0 ; k<np+1 ; k++)
85 for (
int j=0 ; j<nt ; j++) {
86 if (ope_var.operateurs[conte] != 0x0) {
89 sol_hom =
new Tbl(ope_var.operateurs[conte]->get_solh()) ;
94 for (
int i=0 ; i<nr ; i++)
95 so->set(i) = source(zone, k, j, i) ;
97 sol_part =
new Tbl(ope_var.operateurs[conte]->get_solp(*so)) ;
100 for (
int i=0 ; i<nr ; i++) {
101 solution_part.set(zone, k, j, i) = (*sol_part)(i) ;
102 if (sol_hom->get_ndim()==1)
103 solution_hom_un.set(zone, k, j, i) = (*sol_hom)(i) ;
106 solution_hom_un.set(zone, k, j, i) = (*sol_hom)(0,i) ;
107 solution_hom_deux.set(zone, k, j, i) = (*sol_hom)(1,i) ;
129 for (
int k=0 ; k<source.get_mg()->get_np(0)+1 ; k++)
130 for (
int j=0 ; j<source.get_mg()->get_nt(0) ; j++) {
131 if (ope_var.operateurs[start] != 0x0) {
133 int taille = 2*nz - 2 ;
134 Matrice systeme (taille, taille) ;
135 systeme.set_etat_qcq() ;
136 for (
int i=0 ; i<taille ; i++)
137 for (
int j2=0 ; j2<taille ; j2++)
138 systeme.set(i,j2) = 0 ;
139 Tbl sec_membre (taille) ;
140 sec_membre.set_etat_qcq() ;
141 for (
int i=0 ; i<taille ; i++)
142 sec_membre.set(i) = 0 ;
149 systeme.set(0,0) = ope_var.G_plus(0) *
150 ope_var.operateurs[conte]->val_sh_one_plus() ;
154 ope_var.dG_minus(0) * ope_var.operateurs[conte]->val_sh_one_minus() +
155 ope_var.G_minus(0) * ope_var.operateurs[conte]->der_sh_one_minus() ;
157 sec_membre.set(0) -= ope_var.F_plus(0,k,j) +
158 ope_var.G_plus(0) * ope_var.operateurs[conte]->val_sp_plus() ;
160 if ((k==0) && (j==0))
161 sec_membre.set(1) -= -valeur +
162 ope_var.dF_minus(0,k,j) +
163 ope_var.dG_minus(0) * ope_var.operateurs[conte]->val_sp_minus() +
164 ope_var.G_minus(0) * ope_var.operateurs[conte]->der_sp_minus() ;
170 for (
int l=1 ; l<nz-1 ; l++) {
173 int np_prec = source.get_mg()->get_np(l-1) ;
174 int nt_prec = source.get_mg()->get_nt(l-1) ;
175 conte += (np_prec+1)*nt_prec ;
178 systeme.set(2*l-2, 2*l-1) = -ope_var.G_minus(l) *
179 ope_var.operateurs[conte]->val_sh_one_minus() ;
180 systeme.set(2*l-2, 2*l) = - ope_var.G_minus(l) *
181 ope_var.operateurs[conte]->val_sh_two_minus() ;
182 if ((l!=1) || (k!=0) || (j!=0)) {
183 systeme.set(2*l-1, 2*l-1) =
184 -ope_var.dG_minus(l)*ope_var.operateurs[conte]->val_sh_one_minus()-
185 ope_var.G_minus(l)*ope_var.operateurs[conte]->der_sh_one_minus() ;
186 systeme.set(2*l-1, 2*l) =
187 -ope_var.dG_minus(l)*ope_var.operateurs[conte]->val_sh_two_minus()-
188 ope_var.G_minus(l)*ope_var.operateurs[conte]->der_sh_two_minus() ;
190 sec_membre.set(2*l-2) += ope_var.F_minus(l,k,j) +
191 ope_var.G_minus(l) * ope_var.operateurs[conte]->val_sp_minus() ;
192 if ((l!=1) || (k!=0) || (j!=0)) {
193 sec_membre.set(2*l-1) += ope_var.dF_minus(l,k,j) +
194 ope_var.dG_minus(l) * ope_var.operateurs[conte]->val_sp_minus() +
195 ope_var.G_minus(l) * ope_var.operateurs[conte]->der_sp_minus() ;
200 systeme.set(2*l, 2*l-1) = ope_var.G_plus(l) *
201 ope_var.operateurs[conte]->val_sh_one_plus() ;
202 systeme.set(2*l, 2*l) = ope_var.G_plus(l) *
203 ope_var.operateurs[conte]->val_sh_two_plus() ;
204 systeme.set(2*l+1, 2*l-1) =
205 ope_var.dG_plus(l)*ope_var.operateurs[conte]->val_sh_one_plus()+
206 ope_var.G_plus(l)*ope_var.operateurs[conte]->der_sh_one_plus() ;
207 systeme.set(2*l+1, 2*l) =
208 ope_var.dG_plus(l)*ope_var.operateurs[conte]->val_sh_two_plus()+
209 ope_var.G_plus(l)*ope_var.operateurs[conte]->der_sh_two_plus() ;
211 sec_membre.set(2*l) -= ope_var.F_plus(l,k,j) +
212 ope_var.G_plus(l) * ope_var.operateurs[conte]->val_sp_plus();
213 sec_membre.set(2*l+1) -= ope_var.dF_plus(l,k,j) +
214 ope_var.dG_plus(l) * ope_var.operateurs[conte]->val_sp_plus() +
215 ope_var.G_plus(l) * ope_var.operateurs[conte]->der_sp_plus() ;
221 int np_prec = source.get_mg()->get_np(nz-2) ;
222 int nt_prec = source.get_mg()->get_nt(nz-2) ;
223 conte += (np_prec+1)*nt_prec ;
225 systeme.set(taille-2, taille-1) = -ope_var.G_minus(nz-1) *
226 ope_var.operateurs[conte]->val_sh_one_minus() ;
227 systeme.set(taille-1, taille-1) =
228 -ope_var.dG_minus(nz-1)*ope_var.operateurs[conte]->val_sh_one_minus()-
229 ope_var.G_minus(nz-1)*ope_var.operateurs[conte]->der_sh_one_minus() ;
231 sec_membre.set(taille-2) += ope_var.F_minus(nz-1,k,j) +
232 ope_var.G_minus(nz-1)*ope_var.operateurs[conte]->val_sp_minus() ;
233 sec_membre.set(taille-1) += ope_var.dF_minus(nz-1,k,j) +
234 ope_var.dG_minus(nz-1) * ope_var.operateurs[conte]->val_sp_minus() +
235 ope_var.G_minus(nz-1) * ope_var.operateurs[conte]->der_sp_minus() ;
239 systeme.set_band(2,2) ;
241 systeme.set_band(1,1) ;
244 Tbl facteur (systeme.inverse(sec_membre)) ;
248 nr = source.get_mg()->get_nr(0) ;
249 for (
int i=0 ; i<nr ; i++)
250 resultat.set(0,k,j,i) = solution_part(0,k,j,i)
251 +facteur(0)*solution_hom_un(0,k,j,i) ;
254 for (
int l=1 ; l<nz-1 ; l++) {
255 nr = source.get_mg()->get_nr(l) ;
256 for (
int i=0 ; i<nr ; i++)
257 resultat.set(l,k,j,i) = solution_part(l,k,j,i) +
258 facteur(2*l-1)*solution_hom_un(l,k,j,i) +
259 facteur(2*l)*solution_hom_deux(l,k,j,i) ;
263 nr = source.get_mg()->get_nr(nz-1) ;
264 for (
int i=0 ; i<nr ; i++)
265 resultat.set(nz-1,k,j,i) = solution_part(nz-1,k,j,i) +
266 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.