bhole_equations_bin.C

00001 /*
00002  *   Copyright (c) 2001 Philippe Grandclement
00003  *
00004  *   This file is part of LORENE.
00005  *
00006  *   LORENE is free software; you can redistribute it and/or modify
00007  *   it under the terms of the GNU General Public License as published by
00008  *   the Free Software Foundation; either version 2 of the License, or
00009  *   (at your option) any later version.
00010  *
00011  *   LORENE is distributed in the hope that it will be useful,
00012  *   but WITHOUT ANY WARRANTY; without even the implied warranty of
00013  *   MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00014  *   GNU General Public License for more details.
00015  *
00016  *   You should have received a copy of the GNU General Public License
00017  *   along with LORENE; if not, write to the Free Software
00018  *   Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00019  *
00020  */
00021 
00022 
00023 char bhole_equations_bin_C[] = "$Header: /cvsroot/Lorene/C++/Source/Bhole_binaire/bhole_equations_bin.C,v 1.3 2006/04/27 09:12:32 p_grandclement Exp $" ;
00024 
00025 /*
00026  * $Id: bhole_equations_bin.C,v 1.3 2006/04/27 09:12:32 p_grandclement Exp $
00027  * $Log: bhole_equations_bin.C,v $
00028  * Revision 1.3  2006/04/27 09:12:32  p_grandclement
00029  * First try at irrotational black holes
00030  *
00031  * Revision 1.2  2002/10/16 14:36:33  j_novak
00032  * Reorganization of #include instructions of standard C++, in order to
00033  * use experimental version 3 of gcc.
00034  *
00035  * Revision 1.1.1.1  2001/11/20 15:19:28  e_gourgoulhon
00036  * LORENE
00037  *
00038  * Revision 2.6  2001/05/07  09:11:56  phil
00039  * *** empty log message ***
00040  *
00041  * Revision 2.5  2001/04/30  09:30:50  phil
00042  * filtre pour poisson vectoriel
00043  *
00044  * Revision 2.4  2001/04/26  12:04:34  phil
00045  * *** empty log message ***
00046  *
00047  * Revision 2.3  2001/04/25  15:08:48  phil
00048  * vire fait_tkij
00049  *
00050  * Revision 2.2  2001/04/02  12:16:11  phil
00051  * *** empty log message ***
00052  *
00053  * Revision 2.1  2001/03/22  10:41:36  phil
00054  * modification prototypage
00055  *
00056  * Revision 2.0  2001/03/01  08:18:04  phil
00057  * *** empty log message ***
00058  *
00059  *
00060  * $Header: /cvsroot/Lorene/C++/Source/Bhole_binaire/bhole_equations_bin.C,v 1.3 2006/04/27 09:12:32 p_grandclement Exp $
00061  *
00062  */
00063 
00064 //standard
00065 #include <stdlib.h>
00066 #include <math.h>
00067 
00068 // Lorene
00069 #include "nbr_spx.h"
00070 #include "tenseur.h"
00071 #include "bhole.h"
00072 #include "proto.h"
00073 #include "utilitaires.h"
00074 #include "graphique.h"
00075 
00076 // Resolution pour le lapse
00077 void Bhole_binaire::solve_lapse (double precision, double relax) {
00078     
00079     assert ((relax >0) && (relax<=1)) ;
00080     
00081     cout << "-----------------------------------------------" << endl ;
00082     cout << "Resolution LAPSE" << endl ;
00083     
00084     Tenseur lapse_un_old (hole1.n_auto) ;
00085     Tenseur lapse_deux_old (hole2.n_auto) ;
00086     
00087     Tenseur auxi_un (flat_scalar_prod(hole1.tkij_tot, hole1.tkij_auto)) ;
00088     Cmp kk_un (auxi_un(0, 0)+auxi_un(1, 1)+auxi_un(2, 2)) ;
00089     
00090     Tenseur auxi_deux (flat_scalar_prod(hole2.tkij_tot, hole2.tkij_auto)) ;
00091     Cmp kk_deux (auxi_deux(0, 0)+auxi_deux(1, 1)+auxi_deux(2, 2)) ;
00092     
00093     // Les sources
00094     
00095     Cmp source_un 
00096 (-2*flat_scalar_prod(hole1.grad_n_tot, hole1.psi_auto.gradient())()/hole1.psi_tot()
00097     +hole1.n_tot()*pow(hole1.psi_tot(), 4.)*kk_un) ;
00098     source_un.std_base_scal() ;
00099 
00100     Cmp source_deux   
00101 (-2*flat_scalar_prod(hole2.grad_n_tot, hole2.psi_auto.gradient())()/hole2.psi_tot()
00102     +hole2.n_tot()*pow(hole2.psi_tot(), 4.)*kk_deux) ;
00103     source_deux.std_base_scal() ;
00104        
00105     //On resout
00106     hole1.n_auto.set() = hole1.n_auto() - 1./2. ;
00107     hole2.n_auto.set() = hole2.n_auto() - 1./2. ; 
00108     
00109     dirichlet_binaire (source_un, source_deux, -1., -1., hole1.n_auto.set(),
00110                             hole2.n_auto.set(), 0, precision) ;
00111     
00112     hole1.n_auto.set() = hole1.n_auto() + 1./2. ;
00113     hole2.n_auto.set() = hole2.n_auto() + 1./2. ;
00114     
00115     hole1.n_auto.set().raccord(1) ;
00116     hole2.n_auto.set().raccord(1) ;
00117     
00118     // La relaxation :
00119     hole1.n_auto.set() = relax*hole1.n_auto() + (1-relax)*lapse_un_old() ;
00120     hole2.n_auto.set() = relax*hole2.n_auto() + (1-relax)*lapse_deux_old() ;
00121  
00122     hole1.fait_n_comp (hole2) ;
00123     hole2.fait_n_comp (hole1) ;
00124 }
00125 
00126 //Resolution sur Psi
00127 void Bhole_binaire::solve_psi (double precision, double relax) {
00128     
00129     assert ((relax>0) && (relax<=1)) ;
00130     
00131     cout << "-----------------------------------------------" << endl ;
00132     cout << "Resolution PSI" << endl ;
00133     
00134     Tenseur psi_un_old (hole1.psi_auto) ;
00135     Tenseur psi_deux_old (hole2.psi_auto) ;
00136     
00137     Tenseur auxi_un (flat_scalar_prod(hole1.tkij_tot, hole1.tkij_auto)) ;
00138     Cmp kk_un (auxi_un(0, 0)+auxi_un(1, 1)+auxi_un(2, 2)) ;
00139     
00140     Tenseur auxi_deux (flat_scalar_prod(hole2.tkij_tot, hole2.tkij_auto)) ;
00141     Cmp kk_deux (auxi_deux(0, 0)+auxi_deux(1, 1)+auxi_deux(2, 2)) ;
00142     
00143     // Les sources
00144     Cmp source_un (-pow(hole1.psi_tot(), 5.)*kk_un/8.) ;
00145     source_un.std_base_scal() ;
00146     
00147     Cmp source_deux (-pow(hole2.psi_tot(), 5.)*kk_deux/8.)  ;
00148     source_deux.std_base_scal() ;
00149        
00150     // Les valeurs limites :
00151     int np_un = hole1.mp.get_mg()->get_np(1) ;
00152     int nt_un = hole1.mp.get_mg()->get_nt(1) ;
00153     Valeur lim_un (hole1.mp.get_mg()->get_angu()) ;
00154     lim_un = 1 ;
00155     for (int k=0 ; k<np_un ; k++)
00156     for (int j=0 ; j<nt_un ; j++)
00157         lim_un.set(0, k, j, 0) = -0.5/hole1.rayon*hole1.psi_tot()(1, k, j, 0) ;
00158     lim_un.std_base_scal() ;
00159     
00160     int np_deux = hole2.mp.get_mg()->get_np(1) ;
00161     int nt_deux = hole2.mp.get_mg()->get_nt(1) ;
00162     Valeur lim_deux (hole2.mp.get_mg()->get_angu()) ;
00163     lim_deux = 1 ;
00164     for (int k=0 ; k<np_deux ; k++)
00165     for (int j=0 ; j<nt_deux ; j++)
00166         lim_deux.set(0, k, j, 0) = -0.5/hole2.rayon*hole2.psi_tot()(1, k, j, 0) ;
00167     lim_deux.std_base_scal() ;
00168  
00169     //On resout
00170     neumann_binaire (source_un, source_deux, lim_un, lim_deux, 
00171     hole1.psi_auto.set(), hole2.psi_auto.set(), 0, precision) ;
00172     
00173     hole1.psi_auto.set() = hole1.psi_auto() + 1./2. ;
00174     hole2.psi_auto.set() = hole2.psi_auto() + 1./2. ;
00175     
00176     hole1.psi_auto.set().raccord(1) ;
00177     hole2.psi_auto.set().raccord(1) ;
00178     
00179     // La relaxation :
00180     hole1.psi_auto.set() = relax*hole1.psi_auto() + (1-relax)*psi_un_old() ;
00181     hole2.psi_auto.set() = relax*hole2.psi_auto() + (1-relax)*psi_deux_old() ;
00182     
00183     hole1.fait_psi_comp (hole2) ;
00184     hole2.fait_psi_comp (hole1) ;
00185 }
00186 
00187 
00188 // Resolution sur shift a omega bloque.
00189 void Bhole_binaire::solve_shift (double precision, double relax) {
00190     
00191     cout << "------------------------------------------------" << endl ;
00192     cout << "Resolution shift : Omega = " << omega << endl ;
00193     
00194     // On determine les sources 
00195     Tenseur source_un (
00196     -6*flat_scalar_prod (hole1.taij_tot, hole1.psi_auto.gradient())/hole1.psi_tot
00197      +2*flat_scalar_prod (hole1.tkij_tot, hole1.n_auto.gradient())) ; 
00198     if (source_un.get_etat() == ETATZERO) {
00199     source_un.set_etat_qcq() ;
00200     for (int i=0 ; i<3 ; i++)
00201         source_un.set(i).set_etat_zero() ;
00202     }
00203     source_un.set_std_base() ;
00204     
00205     Tenseur source_deux (
00206     -6*flat_scalar_prod (hole2.taij_tot, hole2.psi_auto.gradient())/hole2.psi_tot
00207      +2*flat_scalar_prod (hole2.tkij_tot, hole2.n_auto.gradient())) ;
00208       if (source_deux.get_etat() == ETATZERO) {
00209     source_deux.set_etat_qcq() ;
00210     for (int i=0 ; i<3 ; i++)
00211         source_deux.set(i).set_etat_zero() ;
00212     }
00213     source_deux.set_std_base() ;
00214     
00215     // On filtre les hautes frequences.
00216     for (int i=0 ; i<3 ; i++) {
00217     if (source_un(i).get_etat() != ETATZERO)
00218         source_un.set(i).filtre(4) ;
00219     if (source_deux(i).get_etat() != ETATZERO)
00220         source_deux.set(i).filtre(4) ;
00221     }
00222     
00223     // Les alignemenents pour le signe des CL.
00224     double orientation_un = hole1.mp.get_rot_phi() ;
00225     assert ((orientation_un==0) || (orientation_un == M_PI)) ;
00226     
00227     double orientation_deux = hole2.mp.get_rot_phi() ;
00228     assert ((orientation_deux==0) || (orientation_deux == M_PI)) ;
00229     
00230     int aligne_un = (orientation_un == 0) ? 1 : -1 ;
00231     int aligne_deux = (orientation_deux == 0) ? 1 : -1 ;
00232     
00233     // On regarde si toutes les composantes sont nulles :
00234     int ind1 = 0 ;
00235     int ind2 = 0 ;
00236     for (int i=0 ; i<3 ; i++) {
00237     if (source_un(i).get_etat() == ETATQCQ)
00238         ind1 = 1 ;
00239     if (source_deux(i).get_etat() == ETATQCQ)
00240         ind2 = 1 ;
00241     }
00242     
00243     if (ind1==0)
00244     source_un.set_etat_zero() ;
00245     if (ind2==0)
00246     source_deux.set_etat_zero() ;
00247     
00248     // On determine les Cl en fonction de omega :
00249     int np_un = hole1.mp.get_mg()->get_np (1) ;
00250     int nt_un = hole1.mp.get_mg()->get_nt (1) ;
00251     
00252     Mtbl xa_mtbl_un (source_un.get_mp()->get_mg()) ;
00253     xa_mtbl_un.set_etat_qcq() ;
00254     Mtbl ya_mtbl_un (source_un.get_mp()->get_mg()) ;
00255     ya_mtbl_un.set_etat_qcq() ;
00256     
00257     xa_mtbl_un = source_un.get_mp()->xa ;
00258     ya_mtbl_un = source_un.get_mp()->ya ;
00259       
00260     Mtbl x_mtbl_un (source_un.get_mp()->get_mg()) ;
00261     x_mtbl_un.set_etat_qcq() ;
00262     Mtbl y_mtbl_un (source_un.get_mp()->get_mg()) ;
00263     y_mtbl_un.set_etat_qcq() ;
00264     
00265     x_mtbl_un = source_un.get_mp()->x ;
00266     y_mtbl_un = source_un.get_mp()->y ;
00267     
00268     // Les bases
00269     Base_val** bases_un = hole1.mp.get_mg()->std_base_vect_cart() ;
00270     Base_val** bases_deux = hole2.mp.get_mg()->std_base_vect_cart() ;
00271     
00272     Valeur lim_x_un (*hole1.mp.get_mg()->get_angu()) ;
00273     lim_x_un = 1 ; // Juste pour affecter dans espace des configs ;
00274     lim_x_un.set_etat_c_qcq() ;
00275     for (int k=0 ; k<np_un ; k++)
00276     for (int j=0 ; j<nt_un ; j++)
00277         lim_x_un.set(0, k, j, 0) = aligne_un*omega*ya_mtbl_un(0, 0, 0, 0) + aligne_un*hole1.omega_local*y_mtbl_un(1,k,j,0) ;
00278     lim_x_un.base = *bases_un[0] ;
00279     
00280     Valeur lim_y_un (*hole1.mp.get_mg()->get_angu()) ;
00281     lim_y_un = 1 ; // Juste pour affecter dans espace des configs ;
00282     lim_y_un.set_etat_c_qcq() ;
00283     for (int k=0 ; k<np_un ; k++)
00284     for (int j=0 ; j<nt_un ; j++)
00285         lim_y_un.set(0, k, j, 0) = -aligne_un*omega*xa_mtbl_un(0, 0, 0, 0) - aligne_un*hole1.omega_local*x_mtbl_un(1,k,j,0) ;;
00286     lim_y_un.base = *bases_un[1] ;
00287     
00288     Valeur lim_z_un (*hole1.mp.get_mg()->get_angu()) ;
00289     lim_z_un = 1 ;
00290      for (int k=0 ; k<np_un ; k++)
00291     for (int j=0 ; j<nt_un ; j++)
00292         lim_z_un.set(0, k, j, 0) = 0 ;
00293     lim_z_un.base = *bases_un[2] ;
00294     
00295     // On determine les Cl en fonction de omega :
00296     int np_deux = hole2.mp.get_mg()->get_np (1) ;
00297     int nt_deux = hole2.mp.get_mg()->get_nt (1) ;
00298     
00299     Mtbl xa_mtbl_deux (source_deux.get_mp()->get_mg()) ;
00300     xa_mtbl_deux.set_etat_qcq() ;
00301     Mtbl ya_mtbl_deux (source_deux.get_mp()->get_mg()) ;
00302     ya_mtbl_deux.set_etat_qcq() ;
00303     
00304     xa_mtbl_deux = source_deux.get_mp()->xa ;
00305     ya_mtbl_deux = source_deux.get_mp()->ya ;
00306 
00307     Mtbl x_mtbl_deux (source_deux.get_mp()->get_mg()) ;
00308     x_mtbl_deux.set_etat_qcq() ;
00309     Mtbl y_mtbl_deux (source_deux.get_mp()->get_mg()) ;
00310     y_mtbl_deux.set_etat_qcq() ;
00311     
00312     x_mtbl_deux = source_deux.get_mp()->x ;
00313     y_mtbl_deux = source_deux.get_mp()->y ;
00314     
00315     Valeur lim_x_deux (*hole2.mp.get_mg()->get_angu()) ;
00316     lim_x_deux = 1 ; // Juste pour affecter dans espace des configs ;
00317     lim_x_deux.set_etat_c_qcq() ;
00318     for (int k=0 ; k<np_deux ; k++)
00319     for (int j=0 ; j<nt_deux ; j++)
00320         lim_x_deux.set(0, k, j, 0) = aligne_deux*omega*ya_mtbl_deux(0, 0, 0, 0) + aligne_deux*hole2.omega_local*y_mtbl_deux(1,k,j,0) ;
00321     lim_x_deux.base = *bases_deux[0] ;
00322     
00323     Valeur lim_y_deux (*hole2.mp.get_mg()->get_angu()) ;
00324     lim_y_deux = 1 ; // Juste pour affecter dans espace des configs ;
00325     lim_y_deux.set_etat_c_qcq() ;
00326     for (int k=0 ; k<np_deux ; k++)
00327     for (int j=0 ; j<nt_deux ; j++)
00328        lim_y_deux.set(0, k, j, 0) = -aligne_deux*omega*xa_mtbl_deux(0, 0, 0, 0) - aligne_deux*hole2.omega_local*x_mtbl_deux(1,k,j,0) ;
00329     lim_y_deux.base = *bases_deux[1] ;
00330     
00331     Valeur lim_z_deux (*hole2.mp.get_mg()->get_angu()) ;
00332     lim_z_deux = 1 ;
00333     for (int k=0 ; k<np_deux ; k++)
00334     for (int j=0 ; j<nt_deux ; j++)
00335         lim_z_deux.set(0, k, j, 0) = 0 ;
00336     lim_z_deux.base = *bases_deux[2] ;
00337     
00338     for (int i=0 ; i<3 ; i++) {
00339     delete bases_un[i] ;
00340     delete bases_deux[i] ;
00341     }
00342     delete [] bases_un ;
00343     delete [] bases_deux ;
00344     
00345     // On resout le truc :
00346     Tenseur shift_un_old (hole1.shift_auto) ;
00347     Tenseur shift_deux_old (hole2.shift_auto) ;
00348     
00349     poisson_vect_binaire (1./3., source_un, source_deux, 
00350     lim_x_un, lim_y_un, lim_z_un, 
00351     lim_x_deux, lim_y_deux, lim_z_deux, 
00352     hole1.shift_auto, hole2.shift_auto, 0, precision) ;
00353     
00354     for (int i=0 ; i<3 ; i++) {
00355     hole1.shift_auto.set(i).raccord(1) ;
00356     hole2.shift_auto.set(i).raccord(1) ;
00357     }
00358     
00359     // On regularise les shift.
00360     hole1.shift_auto = relax*hole1.shift_auto + 
00361     (1-relax)*shift_un_old ;
00362     hole2.shift_auto = relax*hole2.shift_auto +
00363     (1-relax)*shift_deux_old ;
00364     
00365     double diff_un = regle (hole1.shift_auto, hole2.shift_auto, omega, hole1.omega_local) ;
00366     double diff_deux = regle (hole2.shift_auto, hole1.shift_auto, omega, hole2.omega_local) ;
00367     hole1.regul = diff_un ;
00368     hole2.regul = diff_deux ;
00369     
00370     cout << "Difference relatives : " << diff_un << " " << diff_deux << endl ;
00371 }

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