00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028 char vector_poisson_C[] = "$Header: /cvsroot/Lorene/C++/Source/Tensor/vector_poisson.C,v 1.22 2005/02/14 13:01:50 j_novak Exp $" ;
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047
00048
00049
00050
00051
00052
00053
00054
00055
00056
00057
00058
00059
00060
00061
00062
00063
00064
00065
00066
00067
00068
00069
00070
00071
00072
00073
00074
00075
00076
00077
00078
00079
00080
00081
00082
00083
00084
00085
00086
00087
00088
00089
00090
00091
00092
00093
00094
00095
00096
00097
00098
00099
00100
00101
00102
00103
00104
00105
00106
00107
00108
00109
00110
00111
00112
00113 #include <stdlib.h>
00114 #include <math.h>
00115
00116
00117 #include "metric.h"
00118 #include "tenseur.h"
00119 #include "param.h"
00120 #include "param_elliptic.h"
00121 #include "proto.h"
00122
00123 Vector Vector::poisson(double lambda, const Metric_flat& met_f, int method)
00124 const {
00125
00126 bool nullite = true ;
00127 for (int i=0; i<3; i++) {
00128 assert(cmp[i]->check_dzpuis(4)) ;
00129 if (cmp[i]->get_etat() != ETATZERO) nullite = false ;
00130 }
00131 assert ((method>=0) && (method<7)) ;
00132
00133 Vector resu(*mp, CON, triad) ;
00134 if (nullite)
00135 resu.set_etat_zero() ;
00136 else {
00137
00138 switch (method) {
00139
00140 case 0 : {
00141
00142 Scalar poten(*mp) ;
00143 if (fabs(lambda+1) < 1.e-8)
00144 poten.set_etat_zero() ;
00145 else {
00146 poten = (potential(met_f) / (lambda + 1)).poisson() ;
00147 }
00148
00149 Vector grad = poten.derive_con(met_f) ;
00150 grad.dec_dzpuis(2) ;
00151
00152 return ( div_free(met_f).poisson() + grad) ;
00153 break ;
00154 }
00155
00156 case 1 : {
00157
00158 Scalar divf(*mp) ;
00159 if (fabs(lambda+1) < 1.e-8)
00160 divf.set_etat_zero() ;
00161 else {
00162 divf = (potential(met_f) / (lambda + 1)) ;
00163 }
00164
00165 int nz = mp->get_mg()->get_nzone() ;
00166
00167
00168
00169
00170 Scalar div_lnot0 = divf ;
00171 div_lnot0.div_r_dzpuis(4) ;
00172 Scalar source_r(*mp) ;
00173 Valeur& va_div = div_lnot0.set_spectral_va() ;
00174 if (div_lnot0.get_etat() != ETATZERO) {
00175 va_div.coef() ;
00176 va_div.ylm() ;
00177 for (int lz=0; lz<nz; lz++) {
00178 int np = mp->get_mg()->get_np(lz) ;
00179 int nt = mp->get_mg()->get_nt(lz) ;
00180 int nr = mp->get_mg()->get_nr(lz) ;
00181 if (va_div.c_cf->operator()(lz).get_etat() != ETATZERO)
00182 for (int k=0; k<np+1; k++)
00183 for (int j=0; j<nt; j++) {
00184 int l_q, m_q, base_r ;
00185 if (nullite_plm(j, nt, k, np, va_div.base) == 1) {
00186 donne_lm(nz, lz, j, k, va_div.base, m_q, l_q, base_r) ;
00187 if (l_q == 0)
00188 for (int i=0; i<nr; i++)
00189 va_div.c_cf->set(lz, k, j, i) = 0. ;
00190 }
00191 }
00192 }
00193 source_r.set_etat_qcq() ;
00194 source_r.set_spectral_va() = 2*(*va_div.c_cf) ;
00195 source_r.set_spectral_va().ylm_i() ;
00196 source_r.set_dzpuis(4) ;
00197 }
00198 else
00199 source_r.set_etat_zero() ;
00200
00201
00202
00203
00204 source_r += *(cmp[0]) - lambda*divf.dsdr() ;
00205 Scalar f_r(*mp) ;
00206 if (source_r.get_etat() != ETATZERO) {
00207
00208
00209
00210
00211 Param_elliptic param_fr(source_r) ;
00212 for (int lz=0; lz<nz; lz++)
00213 param_fr.set_poisson_vect_r(lz) ;
00214
00215 f_r = source_r.sol_elliptic(param_fr) ;
00216 }
00217 else
00218 f_r.set_etat_zero() ;
00219
00220 divf.dec_dzpuis(1) ;
00221 Scalar source_eta = divf - f_r.dsdr() ;
00222 source_eta.mult_r_dzpuis(0) ;
00223 source_eta -= 2*f_r ;
00224 resu.set_vr_eta_mu(f_r, source_eta.poisson_angu(), mu().poisson());
00225
00226 break ;
00227
00228 }
00229
00230 case 2 : {
00231
00232 Tenseur source_p(*mp, 1, CON, mp->get_bvect_spher() ) ;
00233 source_p.set_etat_qcq() ;
00234 for (int i=0; i<3; i++) {
00235 source_p.set(i) = Cmp(*cmp[i]) ;
00236 }
00237 source_p.change_triad(mp->get_bvect_cart()) ;
00238 Tenseur vect_auxi (*mp, 1, CON, mp->get_bvect_cart()) ;
00239 vect_auxi.set_etat_qcq() ;
00240 Tenseur scal_auxi (*mp) ;
00241 scal_auxi.set_etat_qcq() ;
00242
00243 Tenseur resu_p(source_p.poisson_vect(lambda, vect_auxi, scal_auxi)) ;
00244 resu_p.change_triad(mp->get_bvect_spher() ) ;
00245
00246 for (int i=1; i<=3; i++)
00247 resu.set(i) = resu_p(i-1) ;
00248
00249 break ;
00250 }
00251
00252 case 3 : {
00253
00254 Tenseur source_p(*mp, 1, CON, mp->get_bvect_spher() ) ;
00255 source_p.set_etat_qcq() ;
00256 for (int i=0; i<3; i++) {
00257 source_p.set(i) = Cmp(*cmp[i]) ;
00258 }
00259 source_p.change_triad(mp->get_bvect_cart()) ;
00260 Tenseur scal_auxi (*mp) ;
00261 scal_auxi.set_etat_qcq() ;
00262
00263 Tenseur resu_p(source_p.poisson_vect_oohara(lambda, scal_auxi)) ;
00264 resu_p.change_triad(mp->get_bvect_spher() ) ;
00265
00266 for (int i=1; i<=3; i++)
00267 resu.set(i) = resu_p(i-1) ;
00268
00269 break ;
00270 }
00271
00272 case 4 : {
00273 Scalar divf(*mp) ;
00274 if (fabs(lambda+1) < 1.e-8)
00275 divf.set_etat_zero() ;
00276 else {
00277 divf = (potential(met_f) / (lambda + 1)) ;
00278 }
00279
00280 int nz = mp->get_mg()->get_nzone() ;
00281
00282
00283
00284
00285 Scalar div_tmp = divf ;
00286 div_tmp.div_r_dzpuis(4) ;
00287 Scalar source_r(*mp) ;
00288 Valeur& va_div = div_tmp.set_spectral_va() ;
00289 if (div_tmp.get_etat() != ETATZERO) {
00290 va_div.ylm() ;
00291 for (int lz=0; lz<nz; lz++) {
00292 int np = mp->get_mg()->get_np(lz) ;
00293 int nt = mp->get_mg()->get_nt(lz) ;
00294 int nr = mp->get_mg()->get_nr(lz) ;
00295 if (va_div.c_cf->operator()(lz).get_etat() != ETATZERO)
00296 for (int k=0; k<np+1; k++)
00297 for (int j=0; j<nt; j++) {
00298 int l_q, m_q, base_r ;
00299 if (nullite_plm(j, nt, k, np, va_div.base) == 1) {
00300 donne_lm(nz, lz, j, k, va_div.base, m_q, l_q, base_r) ;
00301 if (l_q == 0)
00302 for (int i=0; i<nr; i++)
00303 va_div.c_cf->set(lz, k, j, i) = 0. ;
00304 }
00305 }
00306 }
00307 source_r.set_etat_qcq() ;
00308 source_r.set_spectral_va() = 2*(*va_div.c_cf) ;
00309 source_r.set_spectral_va().ylm_i() ;
00310 source_r.set_dzpuis(4) ;
00311 }
00312 else
00313 source_r.set_etat_zero() ;
00314
00315
00316
00317
00318 source_r += *(cmp[0]) - lambda*divf.dsdr() ;
00319 Scalar f_r(*mp) ;
00320 if (source_r.get_etat() != ETATZERO) {
00321
00322
00323
00324
00325 Param_elliptic param_fr(source_r) ;
00326 for (int lz=0; lz<nz; lz++)
00327 param_fr.set_poisson_vect_r(lz) ;
00328
00329 f_r = source_r.sol_elliptic(param_fr) ;
00330 }
00331 else
00332 f_r.set_etat_zero() ;
00333
00334
00335
00336
00337 Scalar source_eta = *cmp[1] ;
00338 source_eta.mult_sint() ;
00339 source_eta = source_eta.dsdt() ;
00340 source_eta.div_sint() ;
00341 source_eta = (source_eta + cmp[2]->stdsdp()).poisson_angu() ;
00342
00343 Scalar dfrsr = f_r.dsdr() ;
00344 dfrsr.div_r_dzpuis(4) ;
00345 Scalar frsr2 = f_r ;
00346 frsr2.div_r_dzpuis(2) ;
00347 frsr2.div_r_dzpuis(4) ;
00348
00349 Valeur& va_dfr = dfrsr.set_spectral_va() ;
00350 Valeur& va_fsr = frsr2.set_spectral_va() ;
00351 va_dfr.ylm() ;
00352 va_fsr.ylm() ;
00353
00354
00355
00356 Valeur& va_eta = source_eta.set_spectral_va() ;
00357 if (source_eta.get_etat() == ETATZERO) source_eta.annule_hard() ;
00358 va_eta.ylm() ;
00359 for (int lz=0; lz<nz; lz++) {
00360 bool ced = (mp->get_mg()->get_type_r(lz) == UNSURR) ;
00361 int np = mp->get_mg()->get_np(lz) ;
00362 int nt = mp->get_mg()->get_nt(lz) ;
00363 int nr = mp->get_mg()->get_nr(lz) ;
00364 if (va_eta.c_cf->operator()(lz).get_etat() != ETATZERO)
00365 for (int k=0; k<np+1; k++)
00366 for (int j=0; j<nt; j++) {
00367 int l_q, m_q, base_r ;
00368 if (nullite_plm(j, nt, k, np, va_eta.base) == 1) {
00369 donne_lm(nz, lz, j, k, va_eta.base, m_q, l_q, base_r) ;
00370 if (l_q > 0)
00371 for (int i=0; i<nr; i++) {
00372 if (va_div.c_cf->operator()(lz).get_etat() != ETATZERO)
00373 va_eta.c_cf->set(lz, k, j, i)
00374 -= (lambda + 2. / double(ced ? -l_q : (l_q+1) ))
00375 * va_div.c_cf->operator()(lz, k, j, i) ;
00376 if (va_fsr.c_cf->operator()(lz).get_etat() != ETATZERO) {
00377 va_eta.c_cf->set(lz, k, j, i)
00378 += 2. / double(ced ? -l_q : (l_q+1) )
00379 * va_dfr.c_cf->operator()(lz, k, j, i) ;
00380 va_eta.c_cf->set(lz, k, j, i)
00381 -= 2.*( ced ? double(l_q+2)/double(l_q)
00382 : double(l_q-1)/double(l_q+1) )
00383 * va_fsr.c_cf->set(lz, k, j, i) ;
00384 }
00385 }
00386 }
00387 }
00388 }
00389 if (va_eta.c != 0x0) {
00390 delete va_eta.c;
00391 va_eta.c = 0x0 ;
00392 }
00393 va_eta.ylm_i() ;
00394
00395
00396
00397
00398 Param_elliptic param_eta(source_eta) ;
00399 for (int lz=0; lz<nz; lz++)
00400 param_eta.set_poisson_vect_eta(lz) ;
00401
00402 resu.set_vr_eta_mu(f_r, source_eta.sol_elliptic(param_eta), mu().poisson()) ;
00403 break ;
00404
00405 }
00406
00407 case 5 : {
00408
00409 Scalar divf(*mp) ;
00410 if (fabs(lambda+1) < 1.e-8)
00411 divf.set_etat_zero() ;
00412 else {
00413 divf = (potential(met_f) / (lambda + 1)) ;
00414 }
00415
00416 int nz = mp->get_mg()->get_nzone() ;
00417
00418
00419
00420
00421 Scalar div_lnot0 = divf ;
00422 div_lnot0.div_r_dzpuis(4) ;
00423 Scalar source_r(*mp) ;
00424 Valeur& va_div = div_lnot0.set_spectral_va() ;
00425 if (div_lnot0.get_etat() != ETATZERO) {
00426 va_div.coef() ;
00427 va_div.ylm() ;
00428 for (int lz=0; lz<nz; lz++) {
00429 int np = mp->get_mg()->get_np(lz) ;
00430 int nt = mp->get_mg()->get_nt(lz) ;
00431 int nr = mp->get_mg()->get_nr(lz) ;
00432 if (va_div.c_cf->operator()(lz).get_etat() != ETATZERO)
00433 for (int k=0; k<np+1; k++)
00434 for (int j=0; j<nt; j++) {
00435 int l_q, m_q, base_r ;
00436 if (nullite_plm(j, nt, k, np, va_div.base) == 1) {
00437 donne_lm(nz, lz, j, k, va_div.base, m_q, l_q, base_r) ;
00438 if (l_q == 0)
00439 for (int i=0; i<nr; i++)
00440 va_div.c_cf->set(lz, k, j, i) = 0. ;
00441 }
00442 }
00443 }
00444 source_r.set_etat_qcq() ;
00445 source_r.set_spectral_va() = 2*(*va_div.c_cf) ;
00446 source_r.set_spectral_va().ylm_i() ;
00447 source_r.set_dzpuis(4) ;
00448 }
00449 else
00450 source_r.set_etat_zero() ;
00451
00452
00453
00454
00455 source_r += *(cmp[0]) - lambda*divf.dsdr() ;
00456 Scalar f_r(*mp) ;
00457 if (source_r.get_etat() != ETATZERO) {
00458
00459
00460
00461
00462
00463 Param_elliptic param_fr(source_r) ;
00464 for (int lz=0; lz<nz; lz++)
00465 param_fr.set_poisson_vect_r(lz) ;
00466
00467 f_r = source_r.sol_elliptic(param_fr) ;
00468 }
00469 else
00470 f_r.set_etat_zero() ;
00471
00472 Scalar source_eta = - *(cmp[0]) ;
00473 source_eta.mult_r_dzpuis(3) ;
00474 source_eta -= (lambda+2.)*divf ;
00475 source_eta.dec_dzpuis() ;
00476 f_r.set_spectral_va().ylm() ;
00477 Scalar tmp = 2*f_r + f_r.lapang() ;
00478 tmp.div_r_dzpuis(2) ;
00479 source_eta += tmp ;
00480 tmp = (1.+lambda)*divf ;
00481 tmp.mult_r_dzpuis(0) ;
00482 tmp += f_r ;
00483 source_eta = source_eta.primr() ;
00484 f_r.set_spectral_va().ylm_i() ;
00485 resu.set_vr_eta_mu(f_r, (tmp+source_eta).poisson_angu(), mu().poisson()) ;
00486 break ;
00487
00488 }
00489
00490 case 6 : {
00491
00492 poisson_block(lambda, resu) ;
00493 break ;
00494
00495 }
00496 default : {
00497 cout << "Vector::poisson : unexpected type of method !" << endl
00498 << " method = " << method << endl ;
00499 abort() ;
00500 break ;
00501 }
00502
00503 }
00504
00505 }
00506
00507 return resu ;
00508
00509 }
00510
00511 Vector Vector::poisson(double lambda, int method) const {
00512
00513 const Base_vect_spher* tspher = dynamic_cast<const Base_vect_spher*>(triad) ;
00514 const Base_vect_cart* tcart = dynamic_cast<const Base_vect_cart*>(triad) ;
00515
00516 assert ((tspher != 0x0) || (tcart != 0x0)) ;
00517 const Metric_flat* met_f = 0x0 ;
00518
00519 if (tspher != 0x0) {
00520 assert (tcart == 0x0) ;
00521 met_f = &(mp->flat_met_spher()) ;
00522 }
00523
00524 if (tcart != 0x0) {
00525 assert (tspher == 0x0) ;
00526 met_f = &(mp->flat_met_cart()) ;
00527 }
00528
00529 return ( poisson(lambda, *met_f, method) );
00530
00531 }
00532
00533
00534
00535
00536 Vector Vector::poisson(const double lambda, Param& par, int method) const {
00537
00538
00539 for (int i=0; i<3; i++)
00540 assert(cmp[i]->check_dzpuis(4)) ;
00541
00542 assert ((method==0) || (method==2)) ;
00543
00544 Vector resu(*mp, CON, triad) ;
00545
00546 switch (method) {
00547
00548 case 0 : {
00549
00550 Metric_flat met_local(*mp, *triad) ;
00551 int nitermax = par.get_int(0) ;
00552 int& niter = par.get_int_mod(0) ;
00553 double relax = par.get_double(0) ;
00554 double precis = par.get_double(1) ;
00555 Cmp& ss_phi = par.get_cmp_mod(0) ;
00556 Cmp& ss_khi = par.get_cmp_mod(1) ;
00557 Cmp& ss_mu = par.get_cmp_mod(2) ;
00558
00559 Param par_phi ;
00560 par_phi.add_int(nitermax, 0) ;
00561 par_phi.add_int_mod(niter, 0) ;
00562 par_phi.add_double(relax, 0) ;
00563 par_phi.add_double(precis, 1) ;
00564 par_phi.add_cmp_mod(ss_phi, 0) ;
00565
00566 Scalar poten(*mp) ;
00567 if (fabs(lambda+1) < 1.e-8)
00568 poten.set_etat_zero() ;
00569 else {
00570 Scalar tmp = potential(met_local) / (lambda + 1) ;
00571 tmp.inc_dzpuis(2) ;
00572 tmp.poisson(par_phi, poten) ;
00573 }
00574
00575 Vector grad = poten.derive_con(met_local) ;
00576 grad.dec_dzpuis(2) ;
00577
00578 Param par_free ;
00579 par_free.add_int(nitermax, 0) ;
00580 par_free.add_int_mod(niter, 0) ;
00581 par_free.add_double(relax, 0) ;
00582 par_free.add_double(precis, 1) ;
00583 par_free.add_cmp_mod(ss_khi, 0) ;
00584 par_free.add_cmp_mod(ss_mu, 1) ;
00585
00586 return div_free(met_local).poisson(par_free) + grad ;
00587 break ;
00588 }
00589
00590
00591
00592 case 2 : {
00593
00594 Tenseur source_p(*mp, 1, CON, mp->get_bvect_spher() ) ;
00595 source_p.set_etat_qcq() ;
00596 for (int i=0; i<3; i++) {
00597 source_p.set(i) = Cmp(*cmp[i]) ;
00598 }
00599 source_p.change_triad(mp->get_bvect_cart()) ;
00600
00601 Tenseur vect_auxi (*mp, 1, CON, mp->get_bvect_cart()) ;
00602 vect_auxi.set_etat_qcq() ;
00603 for (int i=0; i<3 ;i++){
00604 vect_auxi.set(i) = 0. ;
00605 }
00606 Tenseur scal_auxi (*mp) ;
00607 scal_auxi.set_etat_qcq() ;
00608 scal_auxi.set().annule_hard() ;
00609 scal_auxi.set_std_base() ;
00610
00611 Tenseur resu_p(*mp, 1, CON, mp->get_bvect_cart() ) ;
00612 resu_p.set_etat_qcq() ;
00613 source_p.poisson_vect(lambda, par, resu_p, vect_auxi, scal_auxi) ;
00614 resu_p.change_triad(mp->get_bvect_spher() ) ;
00615
00616 for (int i=1; i<=3; i++)
00617 resu.set(i) = resu_p(i-1) ;
00618
00619 break ;
00620 }
00621
00622 default : {
00623 cout << "Vector::poisson : unexpected type of method !" << endl
00624 << " method = " << method << endl ;
00625
00626
00627 abort() ;
00628 break ;
00629 }
00630
00631 }
00632 return resu ;
00633
00634 }
00635
00636
00637
00638