binhor_equations.C

00001 /*
00002  *   Copyright- (c) 2004-2005 Francois Limousin
00003  *                           Jose-Luis Jaramillo
00004  *
00005  *   This file is part of LORENE.
00006  *
00007  *   LORENE is free software; you can redistribute it and/or modify
00008  *   it under the terms of the GNU General Public License as published by
00009  *   the Free Software Foundation; either version 2 of the License, or
00010  *   (at your option) any later version.
00011  *
00012  *   LORENE is distributed in the hope that it will be useful,
00013  *   but WITHOUT ANY WARRANTY; without even the implied warranty of
00014  *   MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00015  *   GNU General Public License for more details.
00016  *
00017  *   You should have received a copy of the GNU General Public License
00018  *   along with LORENE; if not, write to the Free Software
00019  *   Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00020  *
00021  */
00022 
00023 
00024 char binhor_equations_C[] = "$Header: /cvsroot/Lorene/C++/Source/Bin_hor/binhor_equations.C,v 1.19 2008/02/06 18:20:02 f_limousin Exp $" ;
00025 
00026 /*
00027  * $Id: binhor_equations.C,v 1.19 2008/02/06 18:20:02 f_limousin Exp $
00028  * $Log: binhor_equations.C,v $
00029  * Revision 1.19  2008/02/06 18:20:02  f_limousin
00030  * Fixed an error in the triad
00031  *
00032  * Revision 1.18  2007/08/22 16:12:33  f_limousin
00033  * Changed the name of the function computing \tilde{\gamma}_{ij}
00034  *
00035  * Revision 1.17  2007/04/13 15:28:55  f_limousin
00036  * Lots of improvements, generalisation to an arbitrary state of
00037  * rotation, implementation of the spatial metric given by Samaya.
00038  *
00039  * Revision 1.16  2006/08/01 14:37:19  f_limousin
00040  * New version
00041  *
00042  * Revision 1.15  2006/06/29 08:51:00  f_limousin
00043  * *** empty log message ***
00044  *
00045  * Revision 1.14  2006/05/24 16:56:37  f_limousin
00046  * Many small modifs.
00047  *
00048  * Revision 1.13  2005/11/15 14:04:00  f_limousin
00049  * Minor change to control the resolution of the equation for psi.
00050  *
00051  * Revision 1.12  2005/10/23 16:39:43  f_limousin
00052  * Simplification of the equation in the case of a conformally
00053  * flat metric and maximal slicing
00054  *
00055  * Revision 1.11  2005/09/13 18:33:15  f_limousin
00056  * New function vv_bound_cart_bin(double) for computing binaries with
00057  * berlin condition for the shift vector.
00058  * Suppress all the symy and asymy in the importations.
00059  *
00060  * Revision 1.10  2005/07/11 08:21:57  f_limousin
00061  * Implementation of a new boundary condition for the lapse in the binary
00062  * case : boundary_nn_Dir_lapl().
00063  *
00064  * Revision 1.9  2005/06/09 16:12:04  f_limousin
00065  * Implementation of amg_mom_adm().
00066  *
00067  * Revision 1.8  2005/04/29 14:02:44  f_limousin
00068  * Important changes : manage the dependances between quantities (for
00069  * instance psi and psi4). New function write_global(ost).
00070  *
00071  * Revision 1.7  2005/03/10 17:21:52  f_limousin
00072  * Add the Berlin boundary condition for the shift.
00073  * Some changes to avoid warnings.
00074  *
00075  * Revision 1.6  2005/03/03 10:29:02  f_limousin
00076  * Delete omega in the parameters of the function boundary_beta_z().
00077  *
00078  * Revision 1.5  2005/02/24 17:25:23  f_limousin
00079  * The boundary conditions for psi, N and beta are now parameters in
00080  * par_init.d and par_coal.d.
00081  *
00082  * Revision 1.4  2005/02/11 18:21:38  f_limousin
00083  * Dirichlet_binaire and neumann_binaire are taking Scalars as arguments
00084  * instead of Cmps.
00085  *
00086  * Revision 1.3  2005/02/07 10:46:28  f_limousin
00087  * Many changes !! The sources are written differently to minimize the
00088  * numerical error, the boundary conditions are changed...
00089  *
00090  * Revision 1.2  2004/12/31 15:41:26  f_limousin
00091  * Correction of an error
00092  *
00093  * Revision 1.1  2004/12/29 16:11:34  f_limousin
00094  * First version
00095  *
00096  *
00097  * $Header: /cvsroot/Lorene/C++/Source/Bin_hor/binhor_equations.C,v 1.19 2008/02/06 18:20:02 f_limousin Exp $
00098  *
00099  */
00100 
00101 //standard
00102 #include <stdlib.h>
00103 #include <math.h>
00104 
00105 // Lorene
00106 #include "nbr_spx.h"
00107 #include "tensor.h"
00108 #include "tenseur.h"
00109 #include "isol_hor.h"
00110 #include "proto.h"
00111 #include "utilitaires.h"
00112 //#include "graphique.h"
00113 
00114 // Resolution for the lapse
00115 // ------------------------
00116 
00117 void Bin_hor::solve_lapse (double precision, double relax, int bound_nn,
00118                double lim_nn) {
00119     
00120     assert ((relax >0) && (relax<=1)) ;
00121     
00122     cout << "-----------------------------------------------" << endl ;
00123     cout << "Resolution LAPSE" << endl ;
00124     
00125     Scalar lapse_un_old (hole1.n_auto) ;
00126     Scalar lapse_deux_old (hole2.n_auto) ;
00127 
00128     Sym_tensor taa_un = hole1.aa.up_down(hole1.tgam) ;         
00129     Scalar aa_quad_un = contract(taa_un, 0, 1, hole1.aa_auto, 0, 1) ; 
00130 
00131     Sym_tensor taa_deux = hole2.aa.up_down(hole2.tgam) ;         
00132     Scalar aa_quad_deux = contract(taa_deux, 0, 1, hole2.aa_auto, 0, 1) ; 
00133        
00134     Tensor hdirac1 (contract((hole1.hh).derive_cov(hole1.ff),0,2)) ;
00135     Tensor hdirac2 (contract((hole2.hh).derive_cov(hole2.ff),0,2)) ;
00136 
00137     // Source 1
00138     // --------
00139  
00140     Scalar source_un (hole1.mp) ;
00141     
00142     // Conformally flat
00143     /*
00144     source_un = hole1.get_psi4()*hole1.nn*aa_quad_un
00145       -2.*contract(hole1.dn, 0, hole1.psi_auto
00146             .derive_con(hole1.ff), 0)/hole1.psi ;
00147     */
00148     
00149     source_un = hole1.get_psi4()*( hole1.nn*( aa_quad_un + 0.3333333333333333 * 
00150                         hole1.trK*hole1.trK*hole1.decouple)
00151                    - hole1.trK_point*hole1.decouple )     
00152 
00153        -2.*contract(contract(hole1.hh, 0, hole1.n_auto
00154            .derive_cov(hole1.ff), 0), 0, hole1.dpsi, 0)/hole1.psi
00155 
00156     -2.*contract(hole1.dn, 0, hole1.psi_auto
00157             .derive_con(hole1.ff), 0)/hole1.psi ;
00158 
00159     - contract(hdirac1, 0, hole1.n_auto.derive_cov(hole1.ff), 0) ; 
00160 
00161 
00162     Scalar tmp_un (hole1.mp) ;
00163     
00164     tmp_un = hole1.get_psi4()* contract(hole1.beta_auto, 0, hole1.trK.
00165                     derive_cov(hole1.ff), 0) 
00166     - contract( hole1.hh, 0, 1, hole1.n_auto.derive_cov(hole1.ff)
00167             .derive_cov(hole1.ff), 0, 1 ) ;
00168         
00169     tmp_un.inc_dzpuis() ; // dzpuis: 3 -> 4
00170 
00171     source_un += tmp_un ;
00172 
00173 
00174     // Source 2
00175     // ---------
00176 
00177     Scalar source_deux (hole2.mp) ;
00178 
00179     // Conformally flat
00180     /*
00181     source_deux = hole2.get_psi4()*hole2.nn*aa_quad_deux
00182       -2.*contract(hole2.dn, 0, hole2.psi_auto
00183             .derive_con(hole2.ff), 0)/hole2.psi ;
00184     */
00185     
00186     source_deux = hole2.get_psi4()*( hole2.nn*( aa_quad_deux + 0.3333333333333333
00187                       * hole2.trK*hole2.trK*hole2.decouple)
00188                    - hole2.trK_point*hole2.decouple ) 
00189 
00190     -2.*contract(contract(hole2.hh, 0, hole2.n_auto
00191            .derive_cov(hole2.ff), 0), 0, hole2.dpsi, 0)/hole2.psi
00192 
00193     -2.*contract(hole2.dn, 0, hole2.psi_auto
00194              .derive_con(hole2.ff), 0)/hole2.psi ;
00195 
00196      - contract(hdirac2, 0, hole2.n_auto.derive_cov(hole2.ff), 0) ; 
00197 
00198     Scalar tmp_deux (hole2.mp) ;
00199     
00200     tmp_deux = hole2.get_psi4()* contract(hole2.beta_auto, 0, hole2.trK.
00201                     derive_cov(hole2.ff), 0) 
00202     - contract( hole2.hh, 0, 1, hole2.n_auto.derive_cov(hole2.ff)
00203             .derive_cov(hole2.ff), 0, 1 ) ;
00204         
00205     tmp_deux.inc_dzpuis() ; // dzpuis: 3 -> 4
00206        
00207     source_deux += tmp_deux ;
00208     
00209     cout << "source lapse" << endl << norme(source_un) << endl ;
00210 
00211     // Boundary conditions and resolution
00212     // -----------------------------------
00213 
00214     Valeur lim_un (hole1.mp.get_mg()-> get_angu()) ;
00215     Valeur lim_deux (hole1.mp.get_mg()-> get_angu()) ;
00216 
00217     Scalar n_un_temp (hole1.n_auto) ;
00218     Scalar n_deux_temp (hole2.n_auto) ;
00219 
00220     switch (bound_nn) {
00221 
00222     case 0 : {
00223 
00224       lim_un = hole1.boundary_nn_Dir(lim_nn) ;
00225       lim_deux = hole2.boundary_nn_Dir(lim_nn) ;
00226        
00227         n_un_temp = n_un_temp - 1./2. ;
00228         n_deux_temp = n_deux_temp - 1./2. ;
00229  
00230         dirichlet_binaire (source_un, source_deux, lim_un, lim_deux, 
00231                    n_un_temp, n_deux_temp, 0, precision) ;
00232         break ;
00233     }
00234         
00235     case 1 : {
00236 
00237         lim_un = hole1.boundary_nn_Neu(lim_nn) ;
00238         lim_deux = hole2.boundary_nn_Neu(lim_nn) ;
00239 
00240         neumann_binaire (source_un, source_deux, lim_un, lim_deux, 
00241                    n_un_temp, n_deux_temp, 0, precision) ;
00242         break ;
00243     }
00244         
00245     default : {
00246         cout << "Unexpected type of boundary conditions for the lapse!" 
00247          << endl 
00248          << "  bound_nn = " << bound_nn << endl ; 
00249         abort() ;
00250         break ; 
00251     }
00252         
00253     } // End of switch  
00254 
00255     n_un_temp = n_un_temp + 1./2. ;
00256     n_deux_temp = n_deux_temp + 1./2. ;
00257    
00258     n_un_temp.raccord(3) ;
00259     n_deux_temp.raccord(3) ;
00260     
00261     // Check: has the Poisson equation been correctly solved ?
00262     // -----------------------------------------------------
00263     
00264     int nz = hole1.mp.get_mg()->get_nzone() ;
00265     cout << "lapse auto" << endl << norme (n_un_temp) << endl ;    
00266     Tbl tdiff_nn = diffrel(n_un_temp.laplacian(), source_un) ;
00267     
00268     cout << 
00269     "Relative error in the resolution of the equation for the lapse : "
00270      << endl ; 
00271     for (int l=0; l<nz; l++) {
00272     cout << tdiff_nn(l) << "  " ; 
00273     }
00274     cout << endl ;
00275 
00276     // Relaxation :
00277     // -------------
00278 
00279     n_un_temp = relax*n_un_temp + (1-relax)*lapse_un_old ;
00280     n_deux_temp = relax*n_deux_temp + (1-relax)*lapse_deux_old ;
00281  
00282     hole1.n_auto = n_un_temp ;
00283     hole2.n_auto = n_deux_temp ;
00284 
00285     hole1.n_comp_import(hole2) ;
00286     hole2.n_comp_import(hole1) ;
00287 
00288 }
00289 
00290 
00291 // Resolution for Psi
00292 // -------------------
00293 
00294 void Bin_hor::solve_psi (double precision, double relax, int bound_psi) {
00295     
00296     assert ((relax>0) && (relax<=1)) ;
00297     
00298     cout << "-----------------------------------------------" << endl ;
00299     cout << "Resolution PSI" << endl ;
00300     
00301     Scalar psi_un_old (hole1.psi_auto) ;
00302     Scalar psi_deux_old (hole2.psi_auto) ;
00303     
00304     Sym_tensor taa_un = hole1.aa.up_down(hole1.tgam) ;         
00305     Scalar aa_quad_un = contract(taa_un, 0, 1, hole1.aa_auto, 0, 1) ; 
00306 
00307     Sym_tensor taa_deux = hole2.aa.up_down(hole2.tgam) ;         
00308     Scalar aa_quad_deux = contract(taa_deux, 0, 1, hole2.aa_auto, 0, 1) ; 
00309     
00310     Tensor hdirac1 (contract((hole1.hh).derive_cov(hole1.ff),0,2)) ;
00311     Tensor hdirac2 (contract((hole2.hh).derive_cov(hole2.ff),0,2)) ;
00312 
00313     // Source 1
00314     // ---------
00315     
00316     Scalar source_un (hole1.mp) ;
00317     /*
00318     // Conformally flat source
00319     source_un.annule_hard() ;
00320     source_un.set_dzpuis(4) ;
00321     source_un += - hole1.psi*hole1.get_psi4()* 0.125* aa_quad_un ;
00322     source_un.std_spectral_base() ;
00323     */
00324     
00325     Scalar tmp_un (hole1.mp) ;
00326 
00327     tmp_un = 0.125* hole1.psi_auto * (hole1.tgam).ricci_scal() 
00328       - contract(hole1.hh, 0, 1, hole1.psi_auto.derive_cov(hole1.ff)
00329          .derive_cov(hole1.ff), 0, 1 ) ;
00330     tmp_un.inc_dzpuis() ; // dzpuis : 3 -> 4
00331     
00332     tmp_un -= contract(hdirac1, 0, hole1.psi_auto
00333             .derive_cov(hole1.ff), 0) ;  
00334 
00335     source_un = tmp_un - hole1.psi*hole1.get_psi4()* ( 0.125* aa_quad_un 
00336            - 8.33333333333333e-2* hole1.trK*hole1.trK*hole1.decouple ) ;
00337     source_un.std_spectral_base() ;
00338     
00339 
00340 
00341     // Source 2
00342     // ---------
00343    
00344     Scalar source_deux (hole2.mp) ;
00345     /*
00346     // Conformally flat source
00347     source_deux.annule_hard() ;
00348     source_deux.set_dzpuis(4) ;
00349     source_deux += - hole2.psi*hole2.get_psi4()* 0.125* aa_quad_deux ;
00350     source_deux.std_spectral_base() ;
00351     */
00352     
00353 
00354     Scalar tmp_deux (hole2.mp) ;
00355 
00356     tmp_deux = 0.125* hole2.psi_auto * (hole2.tgam).ricci_scal() 
00357       - contract(hole2.hh, 0, 1, hole2.psi_auto.derive_cov(hole2.ff)
00358          .derive_cov(hole2.ff), 0, 1 ) ;
00359     tmp_deux.inc_dzpuis() ; // dzpuis : 3 -> 4
00360     
00361     tmp_deux -= contract(hdirac2, 0, hole2.psi_auto
00362             .derive_cov(hole2.ff), 0) ;  
00363 
00364     source_deux = tmp_deux - hole2.psi*hole2.get_psi4()* ( 0.125* aa_quad_deux 
00365            - 8.33333333333333e-2* hole2.trK*hole2.trK*hole2.decouple ) ;
00366     source_deux.std_spectral_base() ;
00367     
00368 
00369     cout << "source psi" << endl << norme(source_un) << endl ;
00370 
00371     // Boundary conditions and resolution :
00372     // ------------------------------------
00373 
00374     Valeur lim_un (hole1.mp.get_mg()-> get_angu()) ;
00375     Valeur lim_deux (hole1.mp.get_mg()-> get_angu()) ;
00376 
00377     Scalar psi_un_temp (hole1.psi_auto) ;
00378     Scalar psi_deux_temp (hole2.psi_auto) ;
00379 
00380     switch (bound_psi) {
00381 
00382     case 0 : {
00383 
00384         lim_un = hole1.boundary_psi_app_hor() ;
00385         lim_deux = hole2.boundary_psi_app_hor() ;
00386 
00387         neumann_binaire (source_un, source_deux, lim_un, lim_deux, 
00388                  psi_un_temp, psi_deux_temp, 0, precision) ;
00389         break ;
00390     }
00391 
00392     default : {
00393         cout << "Unexpected type of boundary conditions for psi!" 
00394          << endl 
00395          << "  bound_psi = " << bound_psi << endl ; 
00396         abort() ;
00397         break ; 
00398     }
00399         
00400     } // End of switch  
00401 
00402     psi_un_temp = psi_un_temp + 1./2. ;
00403     psi_deux_temp = psi_deux_temp + 1./2. ;
00404      
00405     psi_un_temp.raccord(3) ;
00406     psi_deux_temp.raccord(3) ;
00407     
00408     // Check: has the Poisson equation been correctly solved ?
00409     // -----------------------------------------------------
00410     
00411     int nz = hole1.mp.get_mg()->get_nzone() ;
00412     cout << "psi auto" << endl << norme (psi_un_temp) << endl ;    
00413     Tbl tdiff_psi = diffrel(psi_un_temp.laplacian(), source_un) ;
00414     
00415     cout << 
00416     "Relative error in the resolution of the equation for psi : "
00417      << endl ; 
00418     for (int l=0; l<nz; l++) {
00419     cout << tdiff_psi(l) << "  " ; 
00420     }
00421     cout << endl ;
00422 
00423     // Relaxation :
00424     // -------------
00425 
00426     psi_un_temp = relax*psi_un_temp + (1-relax)*psi_un_old ;
00427     psi_deux_temp = relax*psi_deux_temp + (1-relax)*psi_deux_old ;
00428     
00429     hole1.psi_auto = psi_un_temp ;
00430     hole2.psi_auto = psi_deux_temp ;
00431 
00432     hole1.psi_comp_import(hole2) ;
00433     hole2.psi_comp_import(hole1) ;
00434 
00435     hole1.set_der_0x0() ;
00436     hole2.set_der_0x0() ;
00437 
00438     //set_hh_Samaya() ;
00439 
00440 }
00441 
00442 
00443 // Resolution for shift with omega fixed.
00444 // --------------------------------------
00445 
00446 void Bin_hor::solve_shift (double precision, double relax, int bound_beta,
00447                double omega_eff) {
00448     
00449     cout << "------------------------------------------------" << endl ;
00450     cout << "Resolution shift : Omega = " << omega << endl ;
00451     
00452     Sym_tensor taa_un = hole1.aa.up_down(hole1.tgam) ;         
00453     Scalar aa_quad_un = contract(taa_un, 0, 1, hole1.aa_auto, 0, 1) ; 
00454 
00455     Sym_tensor taa_deux = hole2.aa.up_down(hole2.tgam) ;         
00456     Scalar aa_quad_deux = contract(taa_deux, 0, 1, hole2.aa_auto, 0, 1) ; 
00457        
00458     Tensor hdirac1 (contract((hole1.hh).derive_cov(hole1.ff),0,2)) ;
00459     Tensor hdirac2 (contract((hole2.hh).derive_cov(hole2.ff),0,2)) ;
00460 
00461     // Source 1
00462     // ---------
00463 
00464     Vector source_un (hole1.mp, CON, hole1.mp.get_bvect_spher()) ;
00465     /*
00466     // Conformally flat source
00467     source_un = 2.* contract(hole1.aa_auto, 1, hole1.dn, 0)
00468       - 12.*hole1.nn*contract(hole1.aa_auto, 1, hole1.dpsi, 0)
00469       /hole1.psi;
00470     */
00471 
00472     
00473     Vector tmp_vect_un (hole1.mp, CON, hole1.mp.get_bvect_spher()) ;
00474     
00475     source_un = 2.* contract(hole1.aa_auto, 1, hole1.dn, 0)
00476       - 12.*hole1.nn*contract(hole1.aa_auto, 1, hole1.dpsi, 0)
00477       /hole1.psi;
00478 
00479      
00480     tmp_vect_un = 2./3.* hole1.trK.derive_con(hole1.tgam)
00481     * hole1.decouple ;
00482     tmp_vect_un.inc_dzpuis() ;
00483    
00484     source_un += 2.* hole1.nn * ( tmp_vect_un
00485             - contract(hole1.tgam.connect().get_delta(), 1, 2, 
00486                      hole1.aa_auto, 0, 1) ) ;
00487 
00488     Vector vtmp_un = contract(hole1.hh, 0, 1, 
00489                            hole1.beta_auto.derive_cov(hole1.ff)
00490                .derive_cov(hole1.ff), 1, 2)
00491       + 1./3.*contract(hole1.hh, 1, hole1.beta_auto
00492                    .divergence(hole1.ff).derive_cov(hole1.ff), 0) 
00493       - hdirac1.derive_lie(hole1.beta_auto) 
00494       + hole1.gamt_point.divergence(hole1.ff)*hole1.decouple ;      
00495     vtmp_un.inc_dzpuis() ; // dzpuis: 3 -> 4
00496 
00497     source_un -= vtmp_un ; 
00498         
00499     source_un += 2./3.* hole1.beta_auto.divergence(hole1.ff) 
00500     * hdirac1 ;
00501 
00502     source_un.std_spectral_base() ;
00503     
00504 
00505     // Source 2
00506     // ---------
00507 
00508     Vector source_deux (hole2.mp, CON, hole2.mp.get_bvect_spher()) ;
00509     /* 
00510     // Conformally flat source
00511     source_deux = 2.* contract(hole2.aa_auto, 1, hole2.dn, 0)
00512       - 12.*hole2.nn*contract(hole2.aa_auto, 1, hole2.dpsi, 0)
00513       /hole2.psi;
00514     */
00515 
00516     
00517     Vector tmp_vect_deux (hole2.mp, CON, hole2.mp.get_bvect_spher()) ;
00518     
00519     source_deux = 2.* contract(hole2.aa_auto, 1, hole2.dn, 0)
00520       - 12.*hole2.nn*contract(hole2.aa_auto, 1, hole2.dpsi, 0)
00521       /hole2.psi;
00522      
00523     tmp_vect_deux = 2./3.* hole2.trK.derive_con(hole2.tgam)
00524     * hole2.decouple ;
00525     tmp_vect_deux.inc_dzpuis() ;
00526    
00527     source_deux += 2.* hole2.nn * ( tmp_vect_deux
00528                - contract(hole2.tgam.connect().get_delta(), 1, 2, 
00529                      hole2.aa*hole2.decouple, 0, 1) ) ;
00530 
00531     Vector vtmp_deux = contract(hole2.hh, 0, 1, 
00532                            hole2.beta_auto.derive_cov(hole2.ff)
00533                .derive_cov(hole2.ff), 1, 2)
00534       + 1./3.*contract(hole2.hh, 1, hole2.beta_auto
00535                    .divergence(hole2.ff).derive_cov(hole2.ff), 0) 
00536       - hdirac2.derive_lie(hole2.beta_auto) 
00537       + hole2.gamt_point.divergence(hole2.ff)*hole2.decouple ;      
00538     vtmp_deux.inc_dzpuis() ; // dzpuis: 3 -> 4
00539 
00540     source_deux -= vtmp_deux ; 
00541         
00542     source_deux += 2./3.* hole2.beta_auto.divergence(hole2.ff) 
00543     * hdirac2 ;
00544 
00545     source_deux.std_spectral_base() ;
00546     
00547 
00548     Vector source_1 (source_un) ;
00549     Vector source_2 (source_deux) ;
00550     source_1.change_triad(hole1.mp.get_bvect_cart()) ;
00551     source_2.change_triad(hole2.mp.get_bvect_cart()) ;
00552 
00553     cout << "source shift_x" << endl << norme(source_1(1)) << endl ;
00554     cout << "source shift_y" << endl << norme(source_1(2)) << endl ;
00555     cout << "source shift_z" << endl << norme(source_1(3)) << endl ;
00556     
00557     // Filter for high frequencies.
00558     for (int i=1 ; i<=3 ; i++) {
00559     source_un.set(i).filtre(4) ;
00560     source_deux.set(i).filtre(4) ;
00561     }
00562     
00563     // Boundary conditions 
00564     // --------------------
00565 
00566     Valeur lim_x_un (hole1.mp.get_mg()-> get_angu()) ;
00567     Valeur lim_y_un (hole1.mp.get_mg()-> get_angu()) ;
00568     Valeur lim_z_un (hole1.mp.get_mg()-> get_angu()) ;
00569 
00570     Valeur lim_x_deux (hole2.mp.get_mg()-> get_angu()) ;
00571     Valeur lim_y_deux (hole2.mp.get_mg()-> get_angu()) ;
00572     Valeur lim_z_deux (hole2.mp.get_mg()-> get_angu()) ;
00573 
00574     switch (bound_beta) {
00575 
00576     case 0 : {
00577         
00578         lim_x_un = hole1.boundary_beta_x(omega, omega_eff) ;
00579         lim_y_un = hole1.boundary_beta_y(omega, omega_eff) ;
00580         lim_z_un = hole1.boundary_beta_z() ;
00581         
00582         lim_x_deux = hole2.boundary_beta_x(omega, omega_eff) ;
00583         lim_y_deux = hole2.boundary_beta_y(omega, omega_eff) ;
00584         lim_z_deux = hole2.boundary_beta_z() ;
00585         break ;
00586     }
00587 
00588     default : {
00589         cout << "Unexpected type of boundary conditions for beta!" 
00590          << endl 
00591          << "  bound_beta = " << bound_beta << endl ; 
00592         abort() ;
00593         break ; 
00594     }
00595         
00596     } // End of switch  
00597 
00598 
00599     // We solve :
00600     // -----------
00601 
00602     Vector beta_un_old (hole1.beta_auto) ;
00603     Vector beta_deux_old (hole2.beta_auto) ;
00604     Vector beta1 (hole1.beta_auto) ;
00605     Vector beta2 (hole2.beta_auto) ;
00606     
00607     poisson_vect_binaire (1./3., source_un, source_deux, 
00608     lim_x_un, lim_y_un, lim_z_un, 
00609     lim_x_deux, lim_y_deux, lim_z_deux, 
00610     beta1, beta2, 0, precision) ;
00611  
00612  
00613     beta1.change_triad(hole1.mp.get_bvect_cart()) ;
00614     beta2.change_triad(hole2.mp.get_bvect_cart()) ;
00615 
00616     for (int i=1 ; i<=3 ; i++) {
00617     beta1.set(i).raccord(3) ;
00618     beta2.set(i).raccord(3) ;
00619     }
00620 
00621     cout << "shift_auto x" << endl << norme(beta1(1)) << endl ;
00622     cout << "shift_auto y" << endl << norme(beta1(2)) << endl ;
00623     cout << "shift_auto z" << endl << norme(beta1(3)) << endl ;
00624 
00625     beta1.change_triad(hole1.mp.get_bvect_spher()) ;
00626     beta2.change_triad(hole2.mp.get_bvect_spher()) ;
00627 
00628     // Check: has the Poisson equation been correctly solved ?
00629     // -----------------------------------------------------
00630     
00631     int nz = hole1.mp.get_mg()->get_nzone() ;
00632     Vector lap_beta = (beta1.derive_con(hole1.ff)).divergence(hole1.ff) 
00633     + 1./3.* beta1.divergence(hole1.ff).derive_con(hole1.ff) ;
00634     source_un.dec_dzpuis() ;
00635 
00636     Tbl tdiff_beta_r = diffrel(lap_beta(1), source_un(1)) ; 
00637     Tbl tdiff_beta_t = diffrel(lap_beta(2), source_un(2)) ; 
00638     Tbl tdiff_beta_p = diffrel(lap_beta(3), source_un(3)) ; 
00639     
00640     cout << 
00641     "Relative error in the resolution of the equation for beta : "
00642      << endl ; 
00643     cout << "r component : " ;
00644     for (int l=0; l<nz; l++) {
00645     cout << tdiff_beta_r(l) << "  " ; 
00646     }
00647     cout << endl ;
00648     cout << "t component : " ;
00649     for (int l=0; l<nz; l++) {
00650     cout << tdiff_beta_t(l) << "  " ; 
00651     }
00652     cout << endl ;
00653     cout << "p component : " ;
00654     for (int l=0; l<nz; l++) {
00655     cout << tdiff_beta_p(l) << "  " ; 
00656     }
00657     cout << endl ;
00658     
00659     
00660     // Relaxation
00661     // -----------
00662 
00663     Vector beta1_new (hole1.mp, CON, hole1.mp.get_bvect_spher()) ;
00664     Vector beta2_new (hole2.mp, CON, hole2.mp.get_bvect_spher()) ;
00665     
00666 
00667     // Construction of Omega d/dphi
00668     // ----------------------------
00669     
00670     Vector omdsdp1 (hole1.mp, CON, hole1.mp.get_bvect_cart()) ;
00671     Scalar yya1 (hole1.mp) ;
00672     yya1 = hole1.mp.ya ;
00673     Scalar xxa1 (hole1.mp) ;
00674     xxa1 = hole1.mp.xa ;
00675    
00676     if (fabs(hole1.mp.get_rot_phi()) < 1e-10){ 
00677       omdsdp1.set(1) = - omega * yya1 ;
00678       omdsdp1.set(2) = omega * xxa1 ;
00679       omdsdp1.set(3).annule_hard() ;
00680     }
00681     else{
00682       omdsdp1.set(1) = omega * yya1 ;
00683       omdsdp1.set(2) = - omega * xxa1 ;
00684       omdsdp1.set(3).annule_hard() ;
00685     }
00686     
00687     omdsdp1.set(1).set_spectral_va()
00688       .set_base(*(hole1.mp.get_mg()->std_base_vect_cart()[0])) ;
00689     omdsdp1.set(2).set_spectral_va()
00690       .set_base(*(hole1.mp.get_mg()->std_base_vect_cart()[1])) ;
00691     omdsdp1.set(3).set_spectral_va()
00692       .set_base(*(hole1.mp.get_mg()->std_base_vect_cart()[2])) ;
00693     
00694     omdsdp1.annule_domain(nz-1) ;
00695     omdsdp1.change_triad(hole1.mp.get_bvect_spher()) ;
00696 
00697 
00698     Vector omdsdp2 (hole2.mp, CON, hole2.mp.get_bvect_cart()) ;
00699     Scalar yya2 (hole2.mp) ;
00700     yya2 = hole2.mp.ya ;
00701     Scalar xxa2 (hole2.mp) ;
00702     xxa2 = hole2.mp.xa ;
00703     
00704     if (fabs(hole2.mp.get_rot_phi()) < 1e-10){ 
00705       omdsdp2.set(1) = - omega * yya2 ;
00706       omdsdp2.set(2) = omega * xxa2 ;
00707       omdsdp2.set(3).annule_hard() ;
00708     }
00709     else{
00710       omdsdp2.set(1) = omega * yya2 ;
00711       omdsdp2.set(2) = - omega * xxa2 ;
00712       omdsdp2.set(3).annule_hard() ;
00713     }
00714     
00715     omdsdp2.set(1).set_spectral_va()
00716       .set_base(*(hole2.mp.get_mg()->std_base_vect_cart()[0])) ;
00717     omdsdp2.set(2).set_spectral_va()
00718       .set_base(*(hole2.mp.get_mg()->std_base_vect_cart()[1])) ;
00719     omdsdp2.set(3).set_spectral_va()
00720       .set_base(*(hole2.mp.get_mg()->std_base_vect_cart()[2])) ;
00721     
00722     omdsdp2.annule_domain(nz-1) ;
00723     omdsdp2.change_triad(hole2.mp.get_bvect_spher()) ;
00724 
00725     // New shift
00726     beta1_new = relax*(beta1+hole1.decouple*omdsdp1) + (1-relax)*beta_un_old ;
00727     beta2_new = relax*(beta2+hole2.decouple*omdsdp2) + (1-relax)*beta_deux_old ;
00728 
00729     hole1.beta_auto = beta1_new ;
00730     hole2.beta_auto = beta2_new ;
00731 
00732     hole1.beta_comp_import(hole2) ;
00733     hole2.beta_comp_import(hole1) ;
00734 
00735     // Regularisation of the shifts if necessary
00736     // -----------------------------------------
00737 
00738     int nnt = hole1.mp.get_mg()->get_nt(1) ;
00739     int nnp = hole1.mp.get_mg()->get_np(1) ;
00740     
00741     int check ;
00742     check = 0 ;
00743     for (int k=0; k<nnp; k++)
00744     for (int j=0; j<nnt; j++){
00745         if (fabs((hole1.n_auto+hole1.n_comp).val_grid_point(1, k, j , 0)) < 1e-8){
00746         check = 1 ;
00747         break ;
00748         }
00749     }
00750 
00751     if (check == 1){
00752     double diff_un = hole1.regularisation (hole1.beta_auto, 
00753                      hole2.beta_auto, omega) ;
00754     double diff_deux = hole2.regularisation (hole2.beta_auto, 
00755                        hole1.beta_auto, omega) ;
00756     hole1.regul = diff_un ;
00757     hole2.regul = diff_deux ;
00758     }
00759     
00760     else {
00761     hole1.regul = 0. ;
00762     hole2.regul = 0. ;
00763     }
00764     
00765 }

Generated on Tue Feb 7 01:35:15 2012 for LORENE by  doxygen 1.4.6