21char sol_elliptic_no_zec_C[] =
"$Header: /cvsroot/Lorene/C++/Source/Non_class_members/PDE/sol_elliptic_no_zec.C,v 1.7 2014/10/13 08:53:30 j_novak Exp $" ;
65#include "param_elliptic.h"
73Mtbl_cf elliptic_solver_no_zec (
const Param_elliptic& ope_var,
const Mtbl_cf& source,
double valeur) {
77 assert (source.get_mg()->get_type_r(0) == RARE) ;
78 assert (source.get_mg()->get_type_r(nz-1) == UNSURR) ;
79 for (
int l=1 ; l<nz-1 ; l++)
80 assert(source.get_mg()->get_type_r(l) == FIN) ;
92 Mtbl_cf solution_part(source.get_mg(), source.base) ;
93 Mtbl_cf solution_hom_un(source.get_mg(), source.base) ;
94 Mtbl_cf solution_hom_deux(source.get_mg(), source.base) ;
95 Mtbl_cf resultat(source.get_mg(), source.base) ;
97 solution_part.annule_hard() ;
98 solution_hom_un.annule_hard() ;
99 solution_hom_deux.annule_hard() ;
100 resultat.annule_hard() ;
104 for (
int zone=0 ; zone<nz-1 ; zone++) {
105 nr = source.get_mg()->get_nr(zone) ;
106 nt = source.get_mg()->get_nt(zone) ;
107 np = source.get_mg()->get_np(zone) ;
109 for (
int k=0 ; k<np+1 ; k++)
110 for (
int j=0 ; j<nt ; j++) {
111 if (ope_var.operateurs[conte] != 0x0) {
113 sol_hom =
new Tbl(ope_var.operateurs[conte]->get_solh()) ;
118 for (
int i=0 ; i<nr ; i++)
119 so->set(i) = source(zone, k, j, i) ;
121 sol_part =
new Tbl(ope_var.operateurs[conte]->get_solp(*so)) ;
124 for (
int i=0 ; i<nr ; i++) {
125 solution_part.set(zone, k, j, i) = (*sol_part)(i) ;
126 if (sol_hom->get_ndim()==1)
127 solution_hom_un.set(zone, k, j, i) = (*sol_hom)(i) ;
130 solution_hom_un.set(zone, k, j, i) = (*sol_hom)(0,i) ;
131 solution_hom_deux.set(zone, k, j, i) = (*sol_hom)(1,i) ;
153 for (
int k=0 ; k<1 ; k++)
154 for (
int j=0 ; j<1 ; j++) {
155 if (ope_var.operateurs[start] != 0x0) {
157 int taille = 2*nz - 3 ;
158 Matrice systeme (taille, taille) ;
159 systeme.set_etat_qcq() ;
160 for (
int i=0 ; i<taille ; i++)
161 for (
int j2=0 ; j2<taille ; j2++)
162 systeme.set(i,j2) = 0 ;
163 Tbl sec_membre (taille) ;
164 sec_membre.set_etat_qcq() ;
165 for (
int i=0 ; i<taille ; i++)
166 sec_membre.set(i) = 0 ;
173 systeme.set(0,0) = ope_var.G_plus(0) *
174 ope_var.operateurs[conte]->val_sh_one_plus() ;
176 ope_var.dG_plus(0) * ope_var.operateurs[conte]->val_sh_one_plus() +
177 ope_var.G_plus(0) * ope_var.operateurs[conte]->der_sh_one_plus() ;
179 sec_membre.set(0) -= ope_var.F_plus(0,k,j) +
180 ope_var.G_plus(0) * ope_var.operateurs[conte]->val_sp_plus() ;
181 sec_membre.set(1) -= ope_var.dF_plus(0,k,j) +
182 ope_var.dG_plus(0) * ope_var.operateurs[conte]->val_sp_plus() +
183 ope_var.G_plus(0) * ope_var.operateurs[conte]->der_sp_plus() ;
190 for (
int l=1 ; l<nz-1 ; l++) {
193 int np_prec = source.get_mg()->get_np(l-1) ;
194 int nt_prec = source.get_mg()->get_nt(l-1) ;
195 conte += (np_prec+1)*nt_prec ;
197 systeme.set(2*l-2, 2*l-1) = -ope_var.G_minus(l) *
198 ope_var.operateurs[conte]->val_sh_one_minus() ;
199 systeme.set(2*l-2, 2*l) = - ope_var.G_minus(l) *
200 ope_var.operateurs[conte]->val_sh_two_minus() ;
201 systeme.set(2*l-1, 2*l-1) =
202 -ope_var.dG_minus(l)*ope_var.operateurs[conte]->val_sh_one_minus()-
203 ope_var.G_minus(l)*ope_var.operateurs[conte]->der_sh_one_minus() ;
204 systeme.set(2*l-1, 2*l) =
205 -ope_var.dG_minus(l)*ope_var.operateurs[conte]->val_sh_two_minus()-
206 ope_var.G_minus(l)*ope_var.operateurs[conte]->der_sh_two_minus() ;
208 sec_membre.set(2*l-2) += ope_var.F_minus(l,k,j) +
209 ope_var.G_minus(l) * ope_var.operateurs[conte]->val_sp_minus() ;
210 sec_membre.set(2*l-1) += ope_var.dF_minus(l,k,j) +
211 ope_var.dG_minus(l) * ope_var.operateurs[conte]->val_sp_minus() +
212 ope_var.G_minus(l) * ope_var.operateurs[conte]->der_sp_minus() ;
215 systeme.set(2*l, 2*l-1) = ope_var.G_plus(l) *
216 ope_var.operateurs[conte]->val_sh_one_plus() ;
217 systeme.set(2*l, 2*l) = ope_var.G_plus(l) *
218 ope_var.operateurs[conte]->val_sh_two_plus() ;
220 systeme.set(2*l+1, 2*l-1) =
221 ope_var.dG_plus(l)*ope_var.operateurs[conte]->val_sh_one_plus()+
222 ope_var.G_plus(l)*ope_var.operateurs[conte]->der_sh_one_plus() ;
223 systeme.set(2*l+1, 2*l) =
224 ope_var.dG_plus(l)*ope_var.operateurs[conte]->val_sh_two_plus()+
225 ope_var.G_plus(l)*ope_var.operateurs[conte]->der_sh_two_plus() ;
229 sec_membre.set(2*l) -= ope_var.F_plus(l,k,j) +
230 ope_var.G_plus(l) * ope_var.operateurs[conte]->val_sp_plus();
232 sec_membre.set(2*l+1) -= ope_var.dF_plus(l,k,j) +
233 ope_var.dG_plus(l) * ope_var.operateurs[conte]->val_sp_plus() +
234 ope_var.G_plus(l) * ope_var.operateurs[conte]->der_sp_plus() ;
237 sec_membre.set(2*l) -= -valeur + ope_var.F_plus(l,k,j) +
238 ope_var.G_plus(l) * ope_var.operateurs[conte]->val_sp_plus();
243 systeme.set_band(2,2) ;
245 systeme.set_band(1,1) ;
248 Tbl facteur (systeme.inverse(sec_membre)) ;
252 nr = source.get_mg()->get_nr(0) ;
253 for (
int i=0 ; i<nr ; i++)
254 resultat.set(0,k,j,i) = solution_part(0,k,j,i)
255 +facteur(0)*solution_hom_un(0,k,j,i) ;
258 for (
int l=1 ; l<nz-1 ; l++) {
259 nr = source.get_mg()->get_nr(l) ;
260 for (
int i=0 ; i<nr ; i++)
261 resultat.set(l,k,j,i) = solution_part(l,k,j,i) +
262 facteur(2*l-1)*solution_hom_un(l,k,j,i) +
263 facteur(2*l)*solution_hom_deux(l,k,j,i) ;
int get_nzone() const
Returns the number of domains.
const Mg3d * get_mg() const
Returns the Mg3d on which the Mtbl_cf is defined.