bound_hor.C

00001 /*
00002  *  Method of class Isol_hor to compute boundary conditions
00003  *
00004  *    (see file isol_hor.h for documentation).
00005  *
00006  */
00007 
00008 /*
00009  *   Copyright (c) 2004  Jose Luis Jaramillo
00010  *                       Francois Limousin
00011  *
00012  *   This file is part of LORENE.
00013  *
00014  *   LORENE is free software; you can redistribute it and/or modify
00015  *   it under the terms of the GNU General Public License version 2
00016  *   as published by the Free Software Foundation.
00017  *
00018  *   LORENE is distributed in the hope that it will be useful,
00019  *   but WITHOUT ANY WARRANTY; without even the implied warranty of
00020  *   MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00021  *   GNU General Public License for more details.
00022  *
00023  *   You should have received a copy of the GNU General Public License
00024  *   along with LORENE; if not, write to the Free Software
00025  *   Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00026  *
00027  */
00028 
00029 char bound_hor_C[] = "$Header: /cvsroot/Lorene/C++/Source/Isol_hor/bound_hor.C,v 1.34 2008/08/27 11:22:25 j_novak Exp $" ;
00030 
00031 /*
00032  * $Id: bound_hor.C,v 1.34 2008/08/27 11:22:25 j_novak Exp $
00033  * $Log: bound_hor.C,v $
00034  * Revision 1.34  2008/08/27 11:22:25  j_novak
00035  * Minor modifications
00036  *
00037  * Revision 1.33  2008/08/19 06:42:00  j_novak
00038  * Minor modifications to avoid warnings with gcc 4.3. Most of them concern
00039  * cast-type operations, and constant strings that must be defined as const char*
00040  *
00041  * Revision 1.32  2007/04/13 15:28:35  f_limousin
00042  * Lots of improvements, generalisation to an arbitrary state of
00043  * rotation, implementation of the spatial metric given by Samaya.
00044  *
00045  * Revision 1.31  2006/08/01 14:35:48  f_limousin
00046  * Many small modifs
00047  *
00048  * Revision 1.30  2006/02/22 17:02:04  f_limousin
00049  * Removal of warnings
00050  *
00051  * Revision 1.29  2006/02/22 16:29:55  jl_jaramillo
00052  * corrections on the relaxation and boundary conditions
00053  *
00054  * Revision 1.27  2005/10/24 16:44:40  jl_jaramillo
00055  * Cook boundary condition ans minot bound of kss
00056  *
00057  * Revision 1.26  2005/10/21 16:20:55  jl_jaramillo
00058  * Version for the paper JaramL05
00059  *
00060  * Revision 1.25  2005/09/13 18:33:17  f_limousin
00061  * New function vv_bound_cart_bin(double) for computing binaries with
00062  * berlin condition for the shift vector.
00063  * Suppress all the symy and asymy in the importations.
00064  *
00065  * Revision 1.24  2005/09/12 12:33:54  f_limousin
00066  * Compilation Warning - Change of convention for the angular velocity
00067  * Add Berlin boundary condition in the case of binary horizons.
00068  *
00069  * Revision 1.23  2005/07/08 13:15:23  f_limousin
00070  * Improvements of boundary_vv_cart(), boundary_nn_lapl().
00071  * Add a fonction to compute the departure of axisymmetry.
00072  *
00073  * Revision 1.22  2005/06/09 08:05:32  f_limousin
00074  * Implement a new function sol_elliptic_boundary() and
00075  * Vector::poisson_boundary(...) which solve the vectorial poisson
00076  * equation (method 6) with an inner boundary condition.
00077  *
00078  * Revision 1.21  2005/05/12 14:48:07  f_limousin
00079  * New boundary condition for the lapse : boundary_nn_lapl().
00080  *
00081  * Revision 1.20  2005/04/29 14:04:17  f_limousin
00082  * Implementation of boundary_vv_x (y,z) for binary black holes.
00083  *
00084  * Revision 1.19  2005/04/19 16:40:51  jl_jaramillo
00085  * change of sign of ang_vel in vv_bound_cart. Convention of Phys. Rep.
00086  *
00087  * Revision 1.18  2005/04/08 12:16:52  f_limousin
00088  * Function set_psi(). And dependance in phi.
00089  *
00090  * Revision 1.17  2005/04/02 15:49:21  f_limousin
00091  * New choice (Lichnerowicz) for aaquad. New member data nz.
00092  *
00093  * Revision 1.16  2005/03/22 13:25:36  f_limousin
00094  * Small changes. The angular velocity and A^{ij} are computed
00095  * with a differnet sign.
00096  *
00097  * Revision 1.15  2005/03/10 10:19:42  f_limousin
00098  * Add the regularisation of the shift in the case of a single black hole
00099  * and lapse zero on the horizon.
00100  *
00101  * Revision 1.14  2005/03/03 10:00:55  f_limousin
00102  * The funtions beta_boost_x() and beta_boost_z() have been added.
00103  *
00104  * Revision 1.13  2005/02/24 17:21:04  f_limousin
00105  * Suppression of the function beta_bound_cart() and direct computation
00106  * of boundary_beta_x, y and z.
00107  *
00108  * Revision 1.12  2004/12/31 15:34:37  f_limousin
00109  * Modifications to avoid warnings
00110  *
00111  * Revision 1.11  2004/12/22 18:15:16  f_limousin
00112  * Many different changes.
00113  *
00114  * Revision 1.10  2004/11/24 19:30:58  jl_jaramillo
00115  * Berlin boundary conditions  vv_bound_cart
00116  *
00117  * Revision 1.9  2004/11/18 09:49:44  jl_jaramillo
00118  * Some new conditions for the shift (Neumann + Dirichlet)
00119  *
00120  * Revision 1.8  2004/11/05 10:52:26  f_limousin
00121  * Replace double aa by double cc in argument of boundary_beta_x
00122  * boundary_beta_y and boundary_beta_z to avoid warnings.
00123  *
00124  * Revision 1.7  2004/10/29 15:42:14  jl_jaramillo
00125  * Static shift boundary conbdition
00126  *
00127  * Revision 1.6  2004/10/01 16:46:51  f_limousin
00128  * Added a pure Dirichlet boundary condition
00129  *
00130  * Revision 1.5  2004/09/28 16:06:41  f_limousin
00131  * Correction of an error when taking the bases of the boundary
00132  * condition for the shift.
00133  *
00134  * Revision 1.4  2004/09/17 13:36:23  f_limousin
00135  * Add some new boundary conditions
00136  *
00137  * Revision 1.2  2004/09/09 16:53:49  f_limousin
00138  * Add the two lines $Id: bound_hor.C,v 1.34 2008/08/27 11:22:25 j_novak Exp $Log: for CVS.
00139  *
00140  *
00141  * $Header: /cvsroot/Lorene/C++/Source/Isol_hor/bound_hor.C,v 1.34 2008/08/27 11:22:25 j_novak Exp $
00142  *
00143  */
00144 
00145 // C++ headers
00146 #include "headcpp.h"
00147 
00148 // C headers
00149 #include <stdlib.h>
00150 #include <assert.h>
00151 
00152 // Lorene headers
00153 #include "time_slice.h"
00154 #include "isol_hor.h"
00155 #include "metric.h"
00156 #include "evolution.h"
00157 #include "unites.h"
00158 #include "graphique.h"
00159 #include "utilitaires.h"
00160 #include "param.h"
00161 
00162 
00163 // Dirichlet boundary condition for Psi 
00164 //-------------------------------------
00165 // ONE HAS TO GUARANTEE THAT BETA IS NOT ZERO, BUT IT IS PROPORTIONAL TO THE RADIAL VECTOR
00166 
00167 const Valeur Isol_hor::boundary_psi_Dir_evol() const{
00168 
00169     Scalar tmp = - 6 * contract(beta(), 0, psi().derive_cov(ff), 0) ;
00170   tmp = tmp / (contract(beta().derive_cov(ff), 0, 1) - nn() * trk() ) - 1 ;
00171 
00172   // We have substracted 1, since we solve for zero condition at infinity 
00173   //and then we add 1 to the solution  
00174 
00175   Valeur psi_bound (mp.get_mg()->get_angu() )  ;
00176   
00177   int nnp = mp.get_mg()->get_np(1) ;
00178   int nnt = mp.get_mg()->get_nt(1) ;
00179             
00180   psi_bound = 1 ;
00181             
00182   for (int k=0 ; k<nnp ; k++)
00183     for (int j=0 ; j<nnt ; j++)
00184       psi_bound.set(0, k, j, 0) = tmp.val_grid_point(1, k, j, 0) ;
00185 
00186   psi_bound.std_base_scal() ;
00187 
00188   return psi_bound ;
00189 
00190 }
00191 
00192 
00193 // Neumann boundary condition for Psi 
00194 //-------------------------------------
00195 
00196 const Valeur Isol_hor::boundary_psi_Neu_evol()const {
00197 
00198   // Introduce 2-trace gamma tilde dot 
00199   Scalar tmp = - 1./ 6. * psi() * (beta().divergence(ff) - nn() * trk() ) 
00200     - beta()(2)* psi().derive_cov(ff)(2) - beta()(3)* psi().derive_cov(ff)(3) ;
00201 
00202   tmp = tmp / beta()(1) ;
00203 
00204   // in this case you don't have to substract any value
00205  
00206   Valeur psi_bound (mp.get_mg()->get_angu() )  ;
00207   
00208   int nnp = mp.get_mg()->get_np(1) ;
00209   int nnt = mp.get_mg()->get_nt(1) ;
00210             
00211   psi_bound = 1 ;
00212             
00213   for (int k=0 ; k<nnp ; k++)
00214     for (int j=0 ; j<nnt ; j++)
00215       psi_bound.set(0, k, j, 0) = tmp.val_grid_point(1, k, j, 0) ;
00216 
00217   psi_bound.std_base_scal() ;
00218   
00219   return psi_bound ;
00220 
00221 }
00222 
00223 
00224 const Valeur Isol_hor::boundary_psi_Dir_spat()const {
00225 
00226     Scalar tmp = psi() * psi() * psi() * trk() 
00227       - contract(k_dd(), 0, 1, tradial_vect_hor() * tradial_vect_hor(), 0, 1) 
00228       / psi()
00229       - 4.* contract(tradial_vect_hor(), 0, psi().derive_cov(ff), 0) ;
00230 
00231     // rho = 1 is the standart case.
00232     double rho = 1. ;
00233     tmp += (tradial_vect_hor().divergence(ff)) * (rho - 1.)*psi() ;
00234 
00235     tmp = tmp / (rho * tradial_vect_hor().divergence(ff)) - 1. ;
00236 
00237     tmp.std_spectral_base() ;
00238     tmp.inc_dzpuis(2) ;
00239 
00240   // We have substracted 1, since we solve for zero condition at infinity 
00241   //and then we add 1 to the solution  
00242  
00243   Valeur psi_bound (mp.get_mg()->get_angu() )  ;
00244   
00245   int nnp = mp.get_mg()->get_np(1) ;
00246   int nnt = mp.get_mg()->get_nt(1) ;
00247             
00248   psi_bound = 1 ;
00249             
00250   for (int k=0 ; k<nnp ; k++)
00251     for (int j=0 ; j<nnt ; j++)
00252       psi_bound.set(0, k, j, 0) = tmp.val_grid_point(1, k, j, 0)  ;
00253 
00254   psi_bound.std_base_scal() ;
00255   
00256   return psi_bound ;
00257 
00258 }
00259 
00260 const Valeur Isol_hor::boundary_psi_Neu_spat()const {
00261 
00262     Scalar tmp = psi() * psi() * psi() * trk() 
00263       - contract(k_dd(), 0, 1, tradial_vect_hor() * tradial_vect_hor(), 0, 1) 
00264       / psi()
00265       - psi() * tradial_vect_hor().divergence(ff) 
00266       - 4 * ( tradial_vect_hor()(2) * psi().derive_cov(ff)(2) 
00267           + tradial_vect_hor()(3) * psi().derive_cov(ff)(3) ) ;
00268 
00269   tmp = tmp / (4 * tradial_vect_hor()(1)) ;
00270 
00271   // in this case you don't have to substract any value
00272  
00273   Valeur psi_bound (mp.get_mg()->get_angu() )  ;
00274   
00275   int nnp = mp.get_mg()->get_np(1) ;
00276   int nnt = mp.get_mg()->get_nt(1) ;
00277             
00278   psi_bound = 1 ;
00279             
00280   for (int k=0 ; k<nnp ; k++)
00281     for (int j=0 ; j<nnt ; j++)
00282       psi_bound.set(0, k, j, 0) = tmp.val_grid_point(1, k, j, 0) ;
00283 
00284   psi_bound.std_base_scal() ;
00285   
00286   return psi_bound ;
00287 
00288 }
00289 
00290 const Valeur Isol_hor::boundary_psi_app_hor()const {
00291 
00292     int nnp = mp.get_mg()->get_np(1) ;
00293     int nnt = mp.get_mg()->get_nt(1) ;
00294 
00295     Valeur psi_bound (mp.get_mg()->get_angu()) ;
00296 
00297     psi_bound = 1 ; // Juste pour affecter dans espace des configs ;
00298 
00299 //    if (psi_comp_evol.is_known(jtime)) {
00300 //    for (int k=0 ; k<nnp ; k++)
00301 //  for (int j=0 ; j<nnt ; j++)
00302 //      psi_bound.set(0, k, j, 0) = - 0.5/radius*(psi_auto().val_grid_point(1, k, j, 0) + psi_comp().val_grid_point(1, k, j, 0)) ;
00303 //    }
00304 //    else {
00305     for (int k=0 ; k<nnp ; k++)
00306     for (int j=0 ; j<nnt ; j++)
00307         psi_bound.set(0, k, j, 0) = - 0.5/radius*psi().val_grid_point(1, k, j, 0) ;
00308 //    }
00309 
00310 
00311     psi_bound.std_base_scal() ;
00312 
00313     return psi_bound ;
00314 
00315 }
00316 
00317 const Valeur Isol_hor::boundary_psi_Dir()const {
00318 
00319   Scalar rho (mp) ;
00320   rho = 5. ;
00321   rho.std_spectral_base() ;
00322 
00323 
00324   Scalar tmp(mp) ;
00325   //  tmp = nn()+1. - 1 ;
00326 
00327   
00328   //Scalar b_tilde_tmp = b_tilde() ;
00329   //b_tilde_tmp.set_domain(0) = 1. ;
00330   //tmp = pow(nn()/b_tilde_tmp, 0.5) ;
00331   
00332 
00333   tmp = 1.5 ;
00334   tmp.std_spectral_base() ;
00335 
00336   //tmp = 1/ (2* nn()) ;
00337 
00338   tmp = (tmp + rho * psi())/(1 + rho) ; 
00339 
00340   tmp = tmp - 1. ;
00341 
00342 
00343 
00344   Valeur psi_bound (mp.get_mg()->get_angu()) ;
00345     
00346   psi_bound = 1 ;  // Juste pour affecter dans espace des configs ;
00347 
00348   int nnp = mp.get_mg()->get_np(1) ;
00349   int nnt = mp.get_mg()->get_nt(1) ;
00350   
00351   for (int k=0 ; k<nnp ; k++)
00352     for (int j=0 ; j<nnt ; j++)
00353       psi_bound.set(0, k, j, 0) = tmp.val_grid_point(1, k, j, 0) ;
00354   
00355   psi_bound.std_base_scal() ;
00356   
00357   return  psi_bound ;
00358 
00359 }
00360 
00361 // Dirichlet boundary condition on nn using the extrinsic curvature
00362 // (No time evolution taken into account! Make this)
00363 //--------------------------------------------------------------------------
00364 const Valeur Isol_hor::boundary_nn_Dir_kk()const {
00365 
00366   Scalar tmp(mp) ;
00367   double rho = 0. ;
00368 
00369   //  Scalar kk_rr = 0.8*psi().derive_cov(mp.flat_met_spher())(1) / psi() ;
00370   Scalar kk_rr = contract( gam().radial_vect() * gam().radial_vect(), 0, 1
00371                , k_dd(), 0, 1 ) ;
00372 
00373   Scalar k_kerr (mp) ;
00374   k_kerr = 0.1 ;//1.*kappa_hor() ;
00375   k_kerr.std_spectral_base() ;
00376   k_kerr.inc_dzpuis(2) ;
00377 
00378   Scalar temp (rho*nn()) ;
00379   temp.inc_dzpuis(2) ;
00380 
00381   tmp = contract(gam().radial_vect(), 0, nn().derive_cov(ff), 0) + temp
00382       - k_kerr ;
00383 
00384   tmp = tmp / (kk_rr + rho) - 1;
00385   
00386   Scalar diN (contract(gam().radial_vect(), 0, nn().derive_cov(ff), 0)) ;
00387   cout << "k_kerr = " << k_kerr.val_grid_point(1, 0, 0, 0) << endl ;
00388   cout << "D^i N  = " << diN.val_grid_point(1, 0, 0, 0) << endl ;
00389   cout << "kss = " << kk_rr.val_grid_point(1, 0, 0, 0) << endl ;
00390 
00391   // We have substracted 1, since we solve for zero condition at infinity 
00392   //and then we add 1 to the solution  
00393 
00394   int nnp = mp.get_mg()->get_np(1) ;
00395   int nnt = mp.get_mg()->get_nt(1) ;
00396 
00397   Valeur nn_bound (mp.get_mg()->get_angu()) ;
00398     
00399   nn_bound = 1 ;  // Juste pour affecter dans espace des configs ;
00400   
00401 
00402   for (int k=0 ; k<nnp ; k++)
00403     for (int j=0 ; j<nnt ; j++)
00404       nn_bound.set(0, k, j, 0) = tmp.val_grid_point(1, k, j, 0) ;
00405   
00406   nn_bound.std_base_scal() ;
00407   
00408   return  nn_bound ;
00409 
00410 }
00411 
00412 
00413 const Valeur Isol_hor::boundary_nn_Neu_kk(int step) const {
00414   
00415   const Vector& dnnn = nn().derive_cov(ff) ;
00416   double rho = 5. ;
00417   
00418   step = 100 ;  // Truc bidon pour eviter warning
00419 
00420   Scalar kk_rr = contract( gam().radial_vect() * gam().radial_vect(), 0, 1
00421                , k_dd(), 0, 1 ) ; 
00422  
00423   Scalar k_kerr (mp) ;
00424   k_kerr = kappa_hor() ; 
00425   //k_kerr = 0.06 ;
00426   k_kerr.std_spectral_base() ;
00427   k_kerr.inc_dzpuis(2) ;
00428   
00429   Scalar k_hor (mp) ;
00430   k_hor =  kappa_hor() ; 
00431   k_hor.std_spectral_base() ;
00432 
00433   //  Vector beta_tilde_d (beta().down(0, met_gamt)) ;
00434   //Scalar naass = 1./2. * contract( met_gamt.radial_vect() * met_gamt.radial_vect(), 0, 1, beta_tilde_d.ope_killing_conf(met_gamt) , 0, 1 ) ;
00435   
00436   Scalar naass = contract( met_gamt.radial_vect() * met_gamt.radial_vect(), 
00437                0, 1, aa_auto().up_down(met_gamt), 0, 1) ;
00438 
00439   Scalar traceK = 1./3. * nn() * trk() * 
00440                contract( met_gamt.radial_vect() * met_gamt.radial_vect(), 0, 1
00441                , met_gamt.cov() , 0, 1 ) ;
00442 
00443   Scalar sdN (contract(gam().radial_vect(), 0, nn().derive_cov(ff), 0)) ;
00444 
00445 
00446 
00447   Scalar tmp =  psi() * psi() * ( k_kerr +  naass + 1.* traceK) 
00448     -  met_gamt.radial_vect()(2) * dnnn(2) 
00449     - met_gamt.radial_vect()(3) * dnnn(3) 
00450     + ( rho - 1 ) * met_gamt.radial_vect()(1)  * dnnn(1)  ;
00451 
00452 
00453   tmp = tmp / (rho * met_gamt.radial_vect()(1))  ;
00454 
00455 
00456 
00457 
00458   Scalar rhs ( sdN - nn() * kk_rr) ;
00459   cout << "kappa_pres = " << k_kerr.val_grid_point(1, 0, 0, 0) << endl ;
00460   cout << "kappa_hor = " << k_hor.val_grid_point(1, 0, 0, 0) << endl ;
00461   cout << "sDN  = " << sdN.val_grid_point(1, 0, 0, 0) << endl ;
00462 
00463   // in this case you don't have to substract any value
00464  
00465   int nnp = mp.get_mg()->get_np(1) ;
00466   int nnt = mp.get_mg()->get_nt(1) ;
00467   
00468   Valeur nn_bound (mp.get_mg()->get_angu()) ;
00469   
00470   nn_bound = 1 ;  // Juste pour affecter dans espace des configs ; 
00471   
00472   for (int k=0 ; k<nnp ; k++)
00473     for (int j=0 ; j<nnt ; j++)
00474       nn_bound.set(0, k, j, 0) = tmp.val_grid_point(1, k, j, 0) ;
00475   
00476   nn_bound.std_base_scal() ;
00477   
00478   return  nn_bound ;
00479   
00480 }
00481 
00482 const Valeur Isol_hor::boundary_nn_Neu_Cook() const {
00483   
00484   const Vector& dnnn = nn().derive_cov(ff) ;
00485   double rho = 1. ;
00486 
00487 
00488   Scalar kk_rr = contract( gam().radial_vect() * gam().radial_vect(), 0, 1
00489                , k_dd(), 0, 1 ) ; 
00490 
00491   Sym_tensor qq_uu = gam_uu() - gam().radial_vect() * gam().radial_vect() ;
00492 
00493   // Free function L_theta = h
00494   //--------------------------
00495   Scalar L_theta (mp) ;
00496   L_theta = .0;
00497   L_theta.std_spectral_base() ;
00498   L_theta.inc_dzpuis(4) ;
00499 
00500   //Divergence of Omega
00501   //-------------------
00502   Vector hom1 = nn().derive_cov(met_gamt)/nn()  ;
00503   hom1 = contract(qq_uu.down(1, gam()), 0, hom1, 0) ;
00504 
00505   Vector hom2 = -contract( k_dd(), 1, gam().radial_vect() , 0) ;
00506   hom2 = contract(qq_uu.down(1, gam()), 0, hom2, 0) ;
00507 
00508   Vector hom = hom1 + hom2 ;
00509   
00510   Scalar div_omega = 1.*contract( qq_uu, 0, 1,  (1.*hom1 + 1.*hom2).derive_cov(gam()), 0, 1) ;
00511   div_omega.inc_dzpuis() ;
00512 
00513   //---------------------
00514 
00515   
00516   // Two-dimensional Ricci scalar
00517   //-----------------------------
00518 
00519   Scalar rr (mp) ;
00520   rr = mp.r ;
00521 
00522   Scalar Ricci_conf(mp) ;
00523   Ricci_conf = 2 / rr / rr ;
00524   Ricci_conf.std_spectral_base() ;
00525 
00526   Scalar Ricci(mp) ;
00527   Scalar log_psi (log(psi())) ;
00528   log_psi.std_spectral_base() ;
00529   Ricci = pow(psi(), -4) * (Ricci_conf - 4*log_psi.lapang())/rr /rr ;
00530   Ricci.std_spectral_base() ;
00531   Ricci.inc_dzpuis(4) ;
00532   //-------------------------------
00533 
00534  Scalar theta_k = -1/(2*nn()) * (gam().radial_vect().divergence(gam()) -
00535                   kk_rr + trk() ) ;
00536   
00537   int nnp = mp.get_mg()->get_np(1) ;
00538   int nnt = mp.get_mg()->get_nt(1) ;
00539   /*  
00540   for (int k=0 ; k<nnp ; k++)
00541     for (int j=0 ; j<nnt ; j++){
00542        cout << theta_k.val_grid_point(1, k, j, 0) << endl ;
00543        cout << nn().val_grid_point(1, k, j, 0) << endl ;
00544     }
00545   */
00546 
00547 
00548 
00549   //Source
00550   //------
00551   Scalar  source = div_omega + contract( qq_uu, 0, 1,  hom * hom , 0, 1) - 0.5 * Ricci - L_theta;
00552   source = source / theta_k ;
00553 
00554   Scalar tmp = ( source + nn() * kk_rr + rho * contract(gam().radial_vect(), 0,
00555                             nn().derive_cov(gam()), 0))/(1+rho) 
00556     - gam().radial_vect()(2) * dnnn(2) - gam().radial_vect()(3) * dnnn(3)  ;
00557   
00558   tmp = tmp / gam().radial_vect()(1) ;
00559 
00560   // in this case you don't have to substract any value
00561  
00562   Valeur nn_bound (mp.get_mg()->get_angu()) ;
00563   
00564   nn_bound = 1 ;  // Juste pour affecter dans espace des configs ; 
00565   
00566   
00567   for (int k=0 ; k<nnp ; k++)
00568     for (int j=0 ; j<nnt ; j++)
00569       nn_bound.set(0, k, j, 0) = tmp.val_grid_point(1, k, j, 0) ;
00570 
00571   nn_bound.std_base_scal() ;
00572   
00573   return  nn_bound ;
00574   
00575 }
00576 
00577 
00578 
00579 const Valeur Isol_hor::boundary_nn_Dir_eff(double cc)const {
00580 
00581   Scalar tmp(mp) ;
00582 
00583   tmp = - nn().derive_cov(ff)(1) ;
00584 
00585   // rho = 1 is the standart case
00586   double rho = 1. ;
00587   tmp.dec_dzpuis(2) ;
00588   tmp += cc * (rho -1)*nn() ;
00589   tmp = tmp / (rho*cc) ;
00590 
00591   tmp = tmp - 1. ;
00592   
00593   // We have substracted 1, since we solve for zero condition at infinity 
00594   //and then we add 1 to the solution  
00595 
00596   int nnp = mp.get_mg()->get_np(1) ;
00597   int nnt = mp.get_mg()->get_nt(1) ;
00598 
00599   Valeur nn_bound (mp.get_mg()->get_angu()) ;
00600     
00601   nn_bound = 1 ;   // Juste pour affecter dans espace des configs ;
00602 
00603   for (int k=0 ; k<nnp ; k++)
00604     for (int j=0 ; j<nnt ; j++)
00605       nn_bound.set(0, k, j, 0) = tmp.val_grid_point(1, k, j, 0) ;
00606   
00607   nn_bound.std_base_scal() ;
00608   
00609   return  nn_bound ;
00610 
00611 }
00612 
00613 
00614 
00615 const Valeur Isol_hor::boundary_nn_Neu_eff(double cc)const  {
00616   
00617   Scalar tmp = - cc * nn() ;
00618   //  Scalar tmp = - nn()/psi()*psi().dsdr() ;
00619 
00620   // in this case you don't have to substract any value
00621  
00622   int nnp = mp.get_mg()->get_np(1) ;
00623   int nnt = mp.get_mg()->get_nt(1) ;
00624 
00625   Valeur nn_bound (mp.get_mg()->get_angu()) ;
00626     
00627   nn_bound = 1 ;   // Juste pour affecter dans espace des configs ;
00628   
00629   for (int k=0 ; k<nnp ; k++)
00630     for (int j=0 ; j<nnt ; j++)
00631       nn_bound.set(0, k, j, 0) = tmp.val_grid_point(1, k, j, 0) ;
00632   
00633   nn_bound.std_base_scal() ;
00634   
00635   return  nn_bound ;
00636 
00637 }
00638 
00639 
00640 const Valeur Isol_hor::boundary_nn_Dir(double cc)const {
00641 
00642   Scalar rho(mp);
00643   rho = 0. ;  // 0 is the standard case
00644 
00645   Scalar tmp(mp) ;
00646   tmp = cc ; 
00647   
00648   /*
00649   if (b_tilde().val_grid_point(1, 0, 0, 0) < 0.08)
00650     tmp = 0.25   ;
00651   else {
00652     cout << "OK" << endl ;
00653     //   des_profile(nn(), 0, 20, M_PI/2, M_PI) ;
00654     rho = 5. ;
00655     tmp =   b_tilde()*psi()*psi() ;
00656   }
00657   */
00658   
00659   //tmp = 1./(2*psi()) ;
00660   //  tmp = - psi() * nn().dsdr() / (psi().dsdr())  ;
00661 
00662   // We  have substracted 1, since we solve for zero condition at infinity 
00663   //and then we add 1 to the solution  
00664 
00665   tmp = (tmp + rho * nn())/(1 + rho) ;
00666 
00667   tmp = tmp - 1 ;
00668 
00669   int nnp = mp.get_mg()->get_np(1) ;
00670   int nnt = mp.get_mg()->get_nt(1) ;
00671 
00672   Valeur nn_bound (mp.get_mg()->get_angu()) ;
00673     
00674   nn_bound = 1 ;  // Juste pour affecter dans espace des configs ;
00675   
00676   for (int k=0 ; k<nnp ; k++)
00677     for (int j=0 ; j<nnt ; j++)
00678       nn_bound.set(0, k, j, 0) = tmp.val_grid_point(1, k, j, 0) ;
00679   
00680   nn_bound.std_base_scal() ;
00681   
00682   return  nn_bound ;
00683 
00684 }
00685 
00686 const Valeur Isol_hor::boundary_nn_Dir_lapl(int mer) const {
00687 
00688   int nnp = mp.get_mg()->get_np(1) ;
00689   int nnt = mp.get_mg()->get_nt(1) ;
00690    
00691 
00692   //Preliminaries
00693   //-------------
00694   
00695   //Metrics on S^2
00696   //--------------
00697 
00698   Sym_tensor q_uu = gam_uu() - gam().radial_vect() * gam().radial_vect() ;
00699   Sym_tensor q_dd = q_uu.up_down(gam()) ;
00700 
00701   Scalar det_q = q_dd(2,2) * q_dd(3,3) - q_dd(2,3) * q_dd(3,2) ;
00702   Scalar square_q (pow(det_q, 0.5)) ;
00703   square_q.std_spectral_base() ;
00704 
00705   Sym_tensor qhat_uu = square_q * q_uu ;
00706 
00707   Sym_tensor ffr_uu = ff.con() - ff.radial_vect() * ff.radial_vect() ;
00708   Sym_tensor ffr_dd = ff.cov() - ff.radial_vect().down(0, ff) * ff.radial_vect().down(0,ff) ;
00709 
00710   Sym_tensor h_uu = qhat_uu - ffr_uu ;
00711 
00712   //------------------
00713 
00714 
00715   // Source of the "round" laplacian: 
00716   //---------------------------------
00717 
00718   //Source 1: "Divergence" term
00719   //---------------------------
00720   Vector kqs = contract(k_dd(), 1, gam().radial_vect(), 0 ) ;
00721   kqs = contract( q_uu.down(0, gam()), 1, kqs, 0) ;
00722   Scalar div_kqs = contract( q_uu, 0, 1, kqs.derive_cov(gam()), 0, 1) ;
00723 
00724   Scalar source1 = div_kqs ;
00725   source1 *= square_q ;
00726 
00727 
00728   // Source 2: correction term
00729   //--------------------------
00730     
00731   Vector corr = - contract(h_uu, 1, nn().derive_cov(ff), 0) / nn() ;
00732   Scalar source2 = contract( ffr_dd, 0, 1, corr.derive_con(ff), 0, 1 ) ;
00733 
00734 
00735   // Source 3: (physical) divergence of Omega
00736   //-----------------------------------------
00737 
00738   Scalar div_omega(mp) ;
00739   
00740   //Source from (L_l\theta_k=0)
00741   
00742   Scalar L_theta (mp) ;
00743   L_theta = 0. ;
00744   L_theta.std_spectral_base() ;
00745   
00746   Scalar kk_rr = contract( gam().radial_vect() * gam().radial_vect(), 0, 1
00747                , k_dd(), 0, 1 ) ; 
00748      
00749  //   Scalar kappa (mp) ;
00750  //   kappa = kappa_hor() ;
00751  //   kappa.std_spectral_base() ;
00752  //   kappa.set_dzpuis(2) ;
00753   
00754  
00755   Scalar kappa = contract(gam().radial_vect(), 0, nn().derive_cov(gam()), 0) - nn() * kk_rr ;
00756   Scalar theta_k = -1/(2*nn()) * (gam().radial_vect().divergence(gam()) -
00757                   kk_rr  + trk() ) ;
00758   
00759   Sym_tensor qqq = gam_uu() - gam().radial_vect() * gam().radial_vect() ;
00760   
00761   Vector hom = nn().derive_cov(met_gamt) - contract( k_dd(), 1, gam().radial_vect() , 0) ;
00762   hom = contract(qqq.down(1, gam()), 0, hom, 0) ;
00763   
00764   Scalar rr(mp) ;
00765   rr = mp.r ;
00766   
00767   Scalar Ricci_conf = 2 / rr /rr ;
00768   Ricci_conf.std_spectral_base() ;
00769   
00770   Scalar log_psi (log(psi())) ;
00771   log_psi.std_spectral_base() ;
00772   Scalar Ricci = pow(psi(), -4) * (Ricci_conf - 4*log_psi.lapang())/rr /rr ;
00773   Ricci.std_spectral_base() ;
00774   Ricci.inc_dzpuis(4) ;
00775   
00776   
00777   div_omega =  L_theta - contract(qqq, 0, 1, hom * hom, 0, 1) + 0.5 * Ricci
00778     + theta_k * kappa ;
00779   div_omega.dec_dzpuis() ;
00780     
00781   //End of Source from (L_l\theta_k=0)
00782   //----------------------------------
00783 
00784   div_omega = 0. ;
00785   //  div_omega = 0.01*log(1/(2*psi())) ;
00786   div_omega.std_spectral_base() ;
00787   div_omega.inc_dzpuis(3) ;
00788 
00789  
00790   //Construction of source3 (correction of div_omega by the square root of the determinant)
00791   Scalar source3 =  div_omega ;
00792   source3 *= square_q ;
00793 
00794 
00795   //"Correction" to the div_omega term (no Y_00 component must be present)
00796   
00797   double corr_const = mp.integrale_surface(source3, radius  + 1e-15)/(4. * M_PI) ;
00798   cout << "Constant part of div_omega = " << corr_const <<endl ;
00799 
00800   Scalar corr_div_omega (mp) ;
00801   corr_div_omega = corr_const ;
00802   corr_div_omega.std_spectral_base() ;
00803   corr_div_omega.set_dzpuis(3) ;
00804 
00805   source3 -=  corr_div_omega ;
00806 
00807  
00808 
00809   //Source
00810   //------
00811   
00812   Scalar source =  source1 + source2 + source3 ;
00813 
00814 
00815 
00816   //Resolution of round angular laplacian
00817   //-------------------------------------
00818   
00819   Scalar source_tmp = source ;
00820 
00821   Scalar logn (mp) ;
00822   logn = source.poisson_angu() ;
00823 
00824   double cc = 0.2 ; // Integration constant
00825 
00826   Scalar lapse (mp) ;
00827   lapse = (exp(logn)*cc) ;
00828   lapse.std_spectral_base() ;
00829   
00830 
00831   //Test of Divergence of Omega
00832   //---------------------------
00833   
00834   ofstream err ;
00835   err.open ("err_laplBC.d", ofstream::app ) ;
00836   
00837   hom = nn().derive_cov(gam())/nn() 
00838     - contract( k_dd(), 1, gam().radial_vect() , 0) ;
00839   hom = contract(q_uu.down(1, gam()), 0, hom, 0) ;
00840   
00841   Scalar div_omega_test = contract( q_uu, 0, 1,  hom.derive_cov(gam()), 0, 1) ;
00842   
00843   
00844   Scalar err_lapl = div_omega_test - div_omega ;
00845   
00846   double max_err = err_lapl.val_grid_point(1, 0, 0, 0) ;
00847   double min_err = err_lapl.val_grid_point(1, 0, 0, 0) ;
00848   
00849   for (int k=0 ; k<nnp ; k++)
00850     for (int j=0 ; j<nnt ; j++){
00851       if (err_lapl.val_grid_point(1, k, j, 0) > max_err)
00852     max_err = err_lapl.val_grid_point(1, k, j, 0) ;
00853       if (err_lapl.val_grid_point(1, k, j, 0) < min_err)
00854     min_err = err_lapl.val_grid_point(1, k, j, 0) ;
00855     }
00856   
00857    err << mer << " " << max_err << " " << min_err << endl ;
00858 
00859    err.close() ;
00860 
00861 
00862 
00863   //Construction of the Valeur
00864   //--------------------------
00865 
00866   lapse = lapse - 1. ;
00867 
00868   //  int nnp = mp.get_mg()->get_np(1) ;
00869   //  int nnt = mp.get_mg()->get_nt(1) ;
00870   
00871   Valeur nn_bound (mp.get_mg()->get_angu()) ;
00872   
00873   nn_bound = 1 ;  // Juste pour affecter dans espace des configs ;
00874   
00875   for (int k=0 ; k<nnp ; k++)
00876     for (int j=0 ; j<nnt ; j++)
00877       nn_bound.set(0, k, j, 0) = lapse.val_grid_point(1, k, j, 0) ;
00878   
00879   nn_bound.std_base_scal() ;
00880   
00881     return  nn_bound ;
00882 
00883 }
00884 
00885 
00886 // Component r of boundary value of beta (using expression in terms 
00887 // of radial vector)
00888 // --------------------------------------
00889 
00890 const Valeur Isol_hor:: boundary_beta_r()const {
00891 
00892   Scalar tmp (mp) ;
00893 
00894   tmp = nn() * radial_vect_hor()(1) ;
00895 
00896   Base_val base_tmp (radial_vect_hor()(1).get_spectral_va().base) ;
00897 
00898   int nnp = mp.get_mg()->get_np(1) ;
00899   int nnt = mp.get_mg()->get_nt(1) ;
00900 
00901   Valeur bnd_beta_r (mp.get_mg()->get_angu()) ;
00902     
00903   bnd_beta_r = 1. ;  // Juste pour affecter dans espace des configs ;
00904   
00905   for (int k=0 ; k<nnp ; k++)
00906     for (int j=0 ; j<nnt ; j++)
00907       bnd_beta_r.set(1, k, j, 0) = tmp.val_grid_point(1, k, j, 0) ;
00908 
00909   bnd_beta_r.set_base_r(0, base_tmp.get_base_r(0)) ;
00910   for (int l=0 ; l<(*mp.get_mg()).get_nzone()-1 ; l++)
00911       bnd_beta_r.set_base_r(l, base_tmp.get_base_r(l)) ;
00912   bnd_beta_r.set_base_r((*mp.get_mg()).get_nzone()-1, base_tmp.get_base_r((*mp.get_mg()).get_nzone()-1)) ;
00913   bnd_beta_r.set_base_t(tmp.get_spectral_va().get_base().get_base_t(1)) ;
00914   bnd_beta_r.set_base_p(tmp.get_spectral_va().get_base().get_base_p(1)) ;
00915 
00916 //  bnd_beta_r.set_base(*(mp.get_mg()->std_base_vect_spher()[0])) ;
00917 
00918   cout << "norme de lim_vr" << endl << norme(bnd_beta_r) << endl ;
00919   cout << "bases" << endl << bnd_beta_r.base << endl ;
00920   
00921   return  bnd_beta_r ;
00922 
00923 
00924 }
00925 
00926 
00927 // Component theta of boundary value of beta (using expression in terms 
00928 // of radial vector)
00929 // ------------------------------------------
00930 
00931 const Valeur Isol_hor::boundary_beta_theta()const {
00932   
00933   Scalar tmp(mp) ;  
00934   
00935   tmp = nn() * radial_vect_hor()(2) ;
00936   Base_val base_tmp (radial_vect_hor()(2).get_spectral_va().base) ;
00937 
00938 
00939   int nnp = mp.get_mg()->get_np(1) ;
00940   int nnt = mp.get_mg()->get_nt(1) ;
00941 
00942   Valeur bnd_beta_theta (mp.get_mg()->get_angu()) ;
00943     
00944   bnd_beta_theta = 1. ;   // Juste pour affecter dans espace des configs ;
00945   
00946   for (int k=0 ; k<nnp ; k++)
00947     for (int j=0 ; j<nnt ; j++)
00948       bnd_beta_theta.set(1, k, j, 0) = tmp.val_grid_point(1, k, j, 0) ;
00949   
00950 //  bnd_beta_theta.set_base(*(mp.get_mg()->std_base_vect_spher()[1])) ;
00951 
00952   bnd_beta_theta.set_base_r(0, base_tmp.get_base_r(0)) ;
00953   for (int l=0 ; l<(*mp.get_mg()).get_nzone()-1 ; l++)
00954       bnd_beta_theta.set_base_r(l, base_tmp.get_base_r(l)) ;
00955   bnd_beta_theta.set_base_r((*mp.get_mg()).get_nzone()-1, base_tmp.get_base_r((*mp.get_mg()).get_nzone()-1)) ;
00956   bnd_beta_theta.set_base_t(base_tmp.get_base_t(1)) ;
00957   bnd_beta_theta.set_base_p(base_tmp.get_base_p(1)) ;
00958 
00959   cout << "bases" << endl << bnd_beta_theta.base << endl ;
00960 
00961 
00962   return  bnd_beta_theta ;
00963 
00964 
00965 }
00966  
00967 // Component phi of boundary value of beta (using expression in terms 
00968 // of radial vector) 
00969 // -----------------------------------------------------------
00970 
00971 const Valeur Isol_hor::boundary_beta_phi(double om)const {
00972 
00973   Scalar tmp (mp) ;
00974 
00975   Scalar ang_vel(mp) ;
00976   ang_vel = om ;
00977   ang_vel.std_spectral_base() ;
00978   ang_vel.mult_rsint() ;
00979 
00980   tmp = nn() * radial_vect_hor()(3)  -  ang_vel ;
00981   Base_val base_tmp (ang_vel.get_spectral_va().base) ;
00982 
00983   int nnp = mp.get_mg()->get_np(1) ;
00984   int nnt = mp.get_mg()->get_nt(1) ;
00985 
00986   Valeur bnd_beta_phi (mp.get_mg()->get_angu()) ;
00987     
00988   bnd_beta_phi = 1. ; // Juste pour affecter dans espace des configs ;
00989   
00990   for (int k=0 ; k<nnp ; k++)
00991     for (int j=0 ; j<nnt ; j++)
00992       bnd_beta_phi.set(1, k, j, 0) = tmp.val_grid_point(1, k, j, 0) ;
00993     
00994   for (int l=0 ; l<(*mp.get_mg()).get_nzone()-1 ; l++)
00995 
00996 //  bnd_beta_phi.set_base(*(mp.get_mg()->std_base_vect_spher()[2])) ;
00997 
00998   bnd_beta_phi.set_base_r(0, base_tmp.get_base_r(0)) ;
00999   for (int l=0 ; l<(*mp.get_mg()).get_nzone()-1 ; l++)
01000       bnd_beta_phi.set_base_r(l, base_tmp.get_base_r(l)) ;
01001   bnd_beta_phi.set_base_r((*mp.get_mg()).get_nzone()-1, base_tmp.get_base_r((*mp.get_mg()).get_nzone()-1)) ;
01002   bnd_beta_phi.set_base_t(base_tmp.get_base_t(1)) ;
01003   bnd_beta_phi.set_base_p(base_tmp.get_base_p(1)) ;
01004 
01005 //  cout << "bound beta_phi" << endl << bnd_beta_phi << endl ;
01006 
01007   return  bnd_beta_phi ;
01008 
01009 }
01010 
01011 // Component x of boundary value of beta
01012 //--------------------------------------
01013 
01014 const Valeur Isol_hor:: boundary_beta_x(double om)const {
01015 
01016     // Les alignemenents pour le signe des CL.
01017     double orientation = mp.get_rot_phi() ;
01018     assert ((orientation == 0) || (orientation == M_PI)) ;
01019     int aligne = (orientation == 0) ? 1 : -1 ;
01020    
01021     int nnp = mp.get_mg()->get_np(1) ;
01022     int nnt = mp.get_mg()->get_nt(1) ;
01023 
01024     Vector tmp_vect = nn() * gam().radial_vect() ;
01025     tmp_vect.change_triad(mp.get_bvect_cart() ) ;
01026 
01027     //Isol_hor boundary conditions
01028   
01029     Valeur lim_x (mp.get_mg()->get_angu()) ;
01030     
01031     lim_x = 1 ;  // Juste pour affecter dans espace des configs ;
01032   
01033     Mtbl ya_mtbl (mp.get_mg()) ;
01034     ya_mtbl.set_etat_qcq() ;
01035     ya_mtbl = mp.ya ;
01036 
01037     Scalar cosp (mp) ;
01038     cosp = mp.cosp ;
01039     Scalar cost (mp) ;
01040     cost = mp.cost ;
01041     Scalar sinp (mp) ;
01042     sinp = mp.sinp ;
01043     Scalar sint (mp) ;
01044     sint = mp.sint ;
01045 
01046     Scalar dep_phi (mp) ;
01047     dep_phi = 0.0*sint*cosp ;
01048 
01049     for (int k=0 ; k<nnp ; k++)
01050     for (int j=0 ; j<nnt ; j++)
01051       lim_x.set(0, k, j, 0) = aligne * om * ya_mtbl(1, k, j, 0) * (1 + 
01052                   dep_phi.val_grid_point(1, k, j, 0))
01053         + tmp_vect(1).val_grid_point(1, k, j, 0) ;
01054     
01055   lim_x.set_base(*(mp.get_mg()->std_base_vect_cart()[0])) ;
01056 
01057   return  lim_x ;
01058 }
01059 
01060 
01061 // Component y of boundary value of beta 
01062 //--------------------------------------
01063 
01064 const Valeur Isol_hor:: boundary_beta_y(double om)const {
01065 
01066     // Les alignemenents pour le signe des CL.
01067     double orientation = mp.get_rot_phi() ;
01068     assert ((orientation == 0) || (orientation == M_PI)) ;
01069     int aligne = (orientation == 0) ? 1 : -1 ;
01070     
01071 
01072     int nnp = mp.get_mg()->get_np(1) ;
01073     int nnt = mp.get_mg()->get_nt(1) ;
01074 
01075     Vector tmp_vect = nn() * gam().radial_vect() ;
01076     tmp_vect.change_triad(mp.get_bvect_cart() ) ;
01077 
01078     //Isol_hor boundary conditions
01079   
01080     Valeur lim_y (mp.get_mg()->get_angu()) ;
01081     
01082     lim_y = 1 ;  // Juste pour affecter dans espace des configs ;
01083   
01084     Mtbl xa_mtbl (mp.get_mg()) ;
01085     xa_mtbl.set_etat_qcq() ;
01086     xa_mtbl = mp.xa ;
01087 
01088     Scalar cosp (mp) ;
01089     cosp = mp.cosp ;
01090     Scalar cost (mp) ;
01091     cost = mp.cost ;
01092     Scalar sinp (mp) ;
01093     sinp = mp.sinp ;
01094     Scalar sint (mp) ;
01095     sint = mp.sint ;
01096 
01097     Scalar dep_phi (mp) ;
01098     dep_phi = 0.0*sint*cosp ;
01099     
01100     for (int k=0 ; k<nnp ; k++)
01101     for (int j=0 ; j<nnt ; j++)
01102       lim_y.set(0, k, j, 0) = - aligne * om * xa_mtbl(1, k, j, 0) *(1 +
01103                     dep_phi.val_grid_point(1, k, j, 0))
01104         + tmp_vect(2).val_grid_point(1, k, j, 0) ;
01105     
01106   lim_y.set_base(*(mp.get_mg()->std_base_vect_cart()[1])) ;
01107 
01108   return  lim_y ;
01109 }
01110 
01111 // Component z of boundary value of beta 
01112 //--------------------------------------
01113 
01114 const Valeur Isol_hor:: boundary_beta_z()const {
01115     
01116     int nnp = mp.get_mg()->get_np(1) ;
01117     int nnt = mp.get_mg()->get_nt(1) ;
01118 
01119     Vector tmp_vect = nn() * gam().radial_vect() ;
01120     tmp_vect.change_triad(mp.get_bvect_cart() ) ;
01121 
01122     //Isol_hor boundary conditions
01123   
01124     Valeur lim_z (mp.get_mg()->get_angu()) ;
01125     
01126     lim_z = 1 ;  // Juste pour affecter dans espace des configs ;
01127   
01128     for (int k=0 ; k<nnp ; k++)
01129     for (int j=0 ; j<nnt ; j++)
01130         lim_z.set(0, k, j, 0) = tmp_vect(3).val_grid_point(1, k, j, 0) ;
01131     
01132   lim_z.set_base(*(mp.get_mg()->std_base_vect_cart()[2])) ;
01133 
01134   return  lim_z ;
01135 }
01136 
01137 const Valeur Isol_hor::beta_boost_x() const {
01138  
01139     Valeur lim_x (mp.get_mg()->get_angu()) ;
01140     lim_x = boost_x ;  // Juste pour affecter dans espace des configs ;
01141   
01142     lim_x.set_base(*(mp.get_mg()->std_base_vect_cart()[0])) ;
01143 
01144     return  lim_x ;
01145 }
01146    
01147 
01148 const Valeur Isol_hor::beta_boost_z() const {
01149  
01150     Valeur lim_z (mp.get_mg()->get_angu()) ;
01151     lim_z = boost_z ;  // Juste pour affecter dans espace des configs ;
01152       
01153     lim_z.set_base(*(mp.get_mg()->std_base_vect_cart()[2])) ;
01154 
01155   return  lim_z ;
01156 }
01157       
01158 
01159 // Neumann boundary condition for b_tilde 
01160 //---------------------------------------
01161 
01162 const Valeur Isol_hor::boundary_b_tilde_Neu()const {
01163   
01164   // Introduce 2-trace gamma tilde dot
01165 
01166   Vector s_tilde = met_gamt.radial_vect() ;
01167 
01168   Scalar hh_tilde = contract(s_tilde.derive_cov(met_gamt), 0, 1) ;
01169 
01170   //des_profile(hh_tilde, 1.00001, 10, M_PI/2., 0., "H_tilde") ;
01171 
01172   Scalar tmp (mp) ;
01173 
01174   tmp = + b_tilde() * hh_tilde - 2 * ( s_tilde(2) * b_tilde().derive_cov(ff)(2)
01175                      + s_tilde(3) * b_tilde().derive_cov(ff)(3) ) ;
01176   
01177   Scalar constant (mp) ;
01178   constant = 0. ;
01179   constant.std_spectral_base() ;
01180   constant.inc_dzpuis(2) ;
01181 
01182   tmp += constant ;
01183   tmp = tmp / (2 *  s_tilde(1) ) ;
01184 
01185   // in this case you don't have to substract any value
01186  
01187   Valeur b_tilde_bound (mp.get_mg()->get_angu() )  ;
01188   
01189   int nnp = mp.get_mg()->get_np(1) ;
01190   int nnt = mp.get_mg()->get_nt(1) ;
01191             
01192   b_tilde_bound = 1 ;
01193             
01194   for (int k=0 ; k<nnp ; k++)
01195     for (int j=0 ; j<nnt ; j++)
01196       b_tilde_bound.set(0, k, j, 0) = tmp.val_grid_point(1, k, j, 0) ;
01197 
01198   b_tilde_bound.std_base_scal() ;
01199   
01200   return b_tilde_bound ;
01201 
01202 }
01203 
01204 
01205 const Valeur Isol_hor::boundary_b_tilde_Dir()const {
01206 
01207   Vector s_tilde = met_gamt.radial_vect() ;
01208 
01209   Scalar hh_tilde = contract(s_tilde.derive_cov(met_gamt), 0, 1) ;
01210 
01211   Scalar tmp = 2 * contract (s_tilde, 0, b_tilde().derive_cov(ff) , 0) ; 
01212 
01213   Scalar constant (mp) ;
01214   constant = -1. ;
01215   constant.std_spectral_base() ;
01216   constant.inc_dzpuis(2) ;
01217 
01218   tmp -= constant ;
01219     
01220   tmp = tmp / hh_tilde   ;
01221 
01222   //  des_profile(tmp, 1.00001, 10, M_PI/2., 0., "tmp") ;
01223 
01224   // We have substracted 1, since we solve for zero condition at infinity 
01225   //and then we add 1 to the solution  
01226 
01227   Valeur b_tilde_bound (mp.get_mg()->get_angu() )  ;
01228   
01229   int nnp = mp.get_mg()->get_np(1) ;
01230   int nnt = mp.get_mg()->get_nt(1) ;
01231             
01232   b_tilde_bound = 1 ;
01233             
01234   for (int k=0 ; k<nnp ; k++)
01235     for (int j=0 ; j<nnt ; j++)
01236       b_tilde_bound.set(0, k, j, 0) = tmp.val_grid_point(1, k, j, 0) ;
01237 
01238   b_tilde_bound.std_base_scal() ;
01239 
01240   return b_tilde_bound ;
01241 
01242 }
01243 
01244 const Vector Isol_hor::vv_bound_cart(double om) const{
01245 
01246     // Preliminaries
01247     //--------------
01248 
01249     // HH_tilde
01250     Vector s_tilde =  met_gamt.radial_vect() ;
01251 
01252     Scalar hh_tilde = contract(s_tilde.derive_cov(met_gamt), 0, 1) ;
01253     hh_tilde.dec_dzpuis(2) ;
01254 
01255     // Tangential part of the shift
01256     Vector tmp_vect = b_tilde() * s_tilde ;
01257 
01258     Vector VV =  b_tilde() * s_tilde - beta() ;
01259 
01260     // "Acceleration" term V^a \tilde{D}_a ln M
01261     Scalar accel = 2 * contract( VV, 0, contract( s_tilde, 0, s_tilde.down(0,
01262                                        met_gamt).derive_cov(met_gamt), 1), 0) ;
01263 
01264 
01265     // Divergence of V^a
01266     Sym_tensor qq_spher = met_gamt.con() - s_tilde * s_tilde ;
01267     Scalar div_VV = contract( qq_spher.down(0, met_gamt), 0, 1, VV.derive_cov(met_gamt), 0, 1) ;
01268 
01269 
01270     Scalar cosp (mp) ;
01271     cosp = mp.cosp ;
01272     Scalar cost (mp) ;
01273     cost = mp.cost ;
01274     Scalar sinp (mp) ;
01275     sinp = mp.sinp ;
01276     Scalar sint (mp) ;
01277     sint = mp.sint ;
01278 
01279     Scalar dep_phi (mp) ;
01280     dep_phi = 0.0*sint*cosp ;
01281 
01282     // Les alignemenents pour le signe des CL.
01283     double orientation = mp.get_rot_phi() ;
01284     assert ((orientation == 0) || (orientation == M_PI)) ;
01285     int aligne = (orientation == 0) ? 1 : -1 ;
01286 
01287     Vector angular (mp, CON, mp.get_bvect_cart()) ;
01288     Scalar yya (mp) ;
01289     yya = mp.ya ;
01290     Scalar xxa (mp) ;
01291     xxa = mp.xa ;
01292 
01293     angular.set(1) = - aligne * om * yya * (1 + dep_phi) ;
01294     angular.set(2) = aligne * om * xxa * (1 + dep_phi) ;
01295     angular.set(3).annule_hard() ;
01296 
01297 
01298     angular.set(1).set_spectral_va()
01299     .set_base(*(mp.get_mg()->std_base_vect_cart()[0])) ;
01300     angular.set(2).set_spectral_va()
01301     .set_base(*(mp.get_mg()->std_base_vect_cart()[1])) ;
01302     angular.set(3).set_spectral_va()
01303     .set_base(*(mp.get_mg()->std_base_vect_cart()[2])) ;
01304 
01305     angular.change_triad(mp.get_bvect_spher()) ;
01306 
01307 /*
01308   Scalar ang_vel (mp) ;
01309   ang_vel = om * (1 + dep_phi) ;
01310   ang_vel.std_spectral_base() ;
01311   ang_vel.mult_rsint() ;
01312 */
01313 
01314     Scalar kss (mp) ;
01315     kss = - 0.2 ;
01316     //    kss.std_spectral_base() ;
01317     //    kss.inc_dzpuis(2) ;
01318   
01319     
01320     //kss from from L_lN=0 condition
01321     //------------------------------
01322     kss = - 0.15 / nn() ;
01323     kss.inc_dzpuis(2) ;
01324     kss += contract(gam().radial_vect(), 0, nn().derive_cov(ff), 0) / nn() ;
01325     
01326 
01327     cout << "kappa_hor = " << kappa_hor() << endl ;
01328 
01329     /*
01330     // Apparent horizon condition
01331     //---------------------------
01332     kss = trk() ;
01333     kss -= (4* contract(psi().derive_cov(met_gamt), 0, met_gamt.radial_vect(), 
01334             0)/psi() +
01335         contract(met_gamt.radial_vect().derive_cov(met_gamt), 0, 1)) /
01336       (psi() * psi()) ;
01337     */
01338   
01339 
01340     // beta^r component
01341     //-----------------
01342     double rho = 5. ; // rho>1 ; rho=1 "pure Dirichlet" version
01343     Scalar beta_r_corr = (rho - 1) * b_tilde() * hh_tilde;
01344     beta_r_corr.inc_dzpuis(2) ;
01345     Scalar nn_KK = nn() * trk() ;
01346     nn_KK.inc_dzpuis(2) ;
01347     
01348     beta_r_corr.set_dzpuis(2) ;
01349     nn_KK.set_dzpuis(2) ;
01350     accel.set_dzpuis(2) ;
01351     div_VV.set_dzpuis(2) ;
01352     
01353     Scalar b_tilde_new (mp) ;
01354     b_tilde_new = 2 * contract(s_tilde, 0, b_tilde().derive_cov(ff), 0)
01355     + beta_r_corr
01356     - 3 * nn() * kss + nn_KK + accel + div_VV  ;
01357     
01358     b_tilde_new = b_tilde_new / (hh_tilde * rho) ;
01359     
01360     b_tilde_new.dec_dzpuis(2) ;
01361     
01362     tmp_vect.set(1) = b_tilde_new * s_tilde(1) + angular(1) ;
01363     tmp_vect.set(2) = b_tilde_new * s_tilde(2) + angular(2) ;
01364     tmp_vect.set(3) = b_tilde_new * s_tilde(3) + angular(3) ;
01365     
01366     tmp_vect.change_triad(mp.get_bvect_cart() ) ;
01367     
01368     return tmp_vect ;    
01369 }
01370 
01371 
01372 const Vector Isol_hor::vv_bound_cart_bin(double, int ) const{
01373   /*
01374   // Les alignemenents pour le signe des CL.
01375   double orientation = mp.get_rot_phi() ;
01376   assert ((orientation == 0) || (orientation == M_PI)) ;
01377   int aligne = (orientation == 0) ? 1 : -1 ;
01378   
01379   Vector angular (mp, CON, mp.get_bvect_cart()) ;
01380   Scalar yya (mp) ;
01381   yya = mp.ya ;
01382   Scalar xxa (mp) ;
01383   xxa = mp.xa ;
01384   
01385   angular.set(1) = aligne * om * yya ;
01386   angular.set(2) = - aligne * om * xxa ;
01387   angular.set(3).annule_hard() ;
01388   
01389   angular.set(1).set_spectral_va()
01390     .set_base(*(mp.get_mg()->std_base_vect_cart()[0])) ;
01391   angular.set(2).set_spectral_va()
01392     .set_base(*(mp.get_mg()->std_base_vect_cart()[1])) ;
01393   angular.set(3).set_spectral_va()
01394     .set_base(*(mp.get_mg()->std_base_vect_cart()[2])) ;
01395   
01396   angular.change_triad(mp.get_bvect_spher()) ;
01397   
01398   // HH_tilde
01399   Vector s_tilde =  met_gamt.radial_vect() ;
01400   
01401   Scalar hh_tilde = contract(s_tilde.derive_cov(met_gamt), 0, 1) ;
01402   hh_tilde.dec_dzpuis(2) ;
01403   
01404   Scalar btilde = b_tilde() - contract(angular, 0, s_tilde.up_down(met_gamt), 0) ;
01405   
01406   // Tangential part of the shift
01407   Vector tmp_vect = btilde * s_tilde ;
01408   
01409   // Value of kss
01410   // --------------
01411   
01412   Scalar kss (mp) ;
01413     kss = - 0.26 ;
01414   kss.std_spectral_base() ;
01415   kss.inc_dzpuis(2) ;
01416   
01417   // Apparent horizon boundary condition
01418   // -----------------------------------
01419 
01420   
01421   // kss frome a fich
01422 
01423   // Construction of the binary
01424   // --------------------------
01425   
01426   char* name_fich = "/home/francois/resu/bin_hor/Test/b11_9x9x8/bin.dat" ;
01427   
01428   FILE* fich_s = fopen(name_fich, "r") ;
01429   Mg3d grid_s (fich_s) ;
01430   Map_af map_un_s (grid_s, fich_s) ;
01431   Map_af map_deux_s (grid_s, fich_s) ;
01432   Bin_hor bin (map_un_s, map_deux_s, fich_s) ;
01433   fclose(fich_s) ;
01434   
01435   // Inititialisation of fields :
01436   // ---------------------------- 
01437   
01438   bin.set(1).n_comp_import (bin(2)) ;
01439   bin.set(1).psi_comp_import (bin(2)) ;
01440   bin.set(2).n_comp_import (bin(1)) ;
01441   bin.set(2).psi_comp_import (bin(1)) ;
01442   bin.decouple() ;
01443   bin.extrinsic_curvature() ;
01444   
01445   kss = contract(bin(jj).get_k_dd(), 0, 1, bin(jj).get_gam().radial_vect()*
01446          bin(jj).get_gam().radial_vect(), 0, 1) ;
01447 
01448 
01449   cout << "kappa_hor = " << kappa_hor() << endl ;
01450   
01451   // beta^r component
01452   //-----------------
01453   double rho = 5. ; // rho=0 "pure Dirichlet" version
01454   Scalar beta_r_corr = rho * btilde * hh_tilde;
01455   beta_r_corr.inc_dzpuis(2) ;
01456   Scalar nn_KK = nn() * trk() ;
01457   nn_KK.inc_dzpuis(2) ;
01458   
01459   beta_r_corr.set_dzpuis(2) ;
01460   nn_KK.set_dzpuis(2) ;
01461   
01462   Scalar b_tilde_new (mp) ;
01463   b_tilde_new = 2 * contract(s_tilde, 0, btilde.derive_cov(ff), 0)
01464     + beta_r_corr - 3 * nn() * kss + nn_KK ;
01465   
01466   b_tilde_new = b_tilde_new / (hh_tilde * (rho+1)) ;
01467   
01468   b_tilde_new.dec_dzpuis(2) ;
01469   
01470   tmp_vect.set(1) = b_tilde_new * s_tilde(1) + angular(1) ;
01471   tmp_vect.set(2) = b_tilde_new * s_tilde(2) + angular(2) ;
01472   tmp_vect.set(3) = b_tilde_new * s_tilde(3) + angular(3) ; 
01473   
01474   tmp_vect.change_triad(mp.get_bvect_cart() ) ;
01475   
01476   return tmp_vect ;
01477   */
01478     Vector pipo(mp, CON, mp.get_bvect_cart()) ;
01479     return pipo ;
01480 }
01481 
01482 
01483 // Component x of boundary value of V^i 
01484 //-------------------------------------
01485 
01486 const Valeur Isol_hor:: boundary_vv_x(double om)const {
01487   
01488   int nnp = mp.get_mg()->get_np(1) ;
01489   int nnt = mp.get_mg()->get_nt(1) ;
01490 
01491   //Isol_hor boundary conditions
01492   
01493   Valeur lim_x (mp.get_mg()->get_angu()) ;
01494     
01495   lim_x = 1 ;  // Juste pour affecter dans espace des configs ;
01496  
01497   Scalar vv_x = vv_bound_cart(om)(1) ;
01498 
01499   for (int k=0 ; k<nnp ; k++)
01500     for (int j=0 ; j<nnt ; j++)
01501       lim_x.set(0, k, j, 0) = vv_x.val_grid_point(1, k, j, 0) ;
01502   
01503   lim_x.set_base(*(mp.get_mg()->std_base_vect_cart()[0])) ;
01504 
01505   return  lim_x ;
01506 
01507 
01508 }
01509 
01510 
01511 // Component y of boundary value of V^i
01512 //--------------------------------------
01513 
01514 const Valeur Isol_hor:: boundary_vv_y(double om)const {
01515  
01516   int nnp = mp.get_mg()->get_np(1) ;
01517   int nnt = mp.get_mg()->get_nt(1) ;
01518 
01519   // Isol_hor boundary conditions
01520   
01521   Valeur lim_y (mp.get_mg()->get_angu()) ;
01522     
01523   lim_y = 1 ;  // Juste pour affecter dans espace des configs ;
01524  
01525   Scalar vv_y = vv_bound_cart(om)(2) ;
01526 
01527   for (int k=0 ; k<nnp ; k++)
01528       for (int j=0 ; j<nnt ; j++)
01529       lim_y.set(0, k, j, 0) = vv_y.val_grid_point(1, k, j, 0) ;
01530 
01531   lim_y.set_base(*(mp.get_mg()->std_base_vect_cart()[1])) ;
01532 
01533   return  lim_y ;
01534 }
01535 
01536 
01537 // Component z of boundary value of V^i
01538 //-------------------------------------
01539 
01540 const Valeur Isol_hor:: boundary_vv_z(double om)const {
01541 
01542   int nnp = mp.get_mg()->get_np(1) ;
01543   int nnt = mp.get_mg()->get_nt(1) ;
01544 
01545   // Isol_hor boundary conditions
01546  
01547   Valeur lim_z (mp.get_mg()->get_angu()) ;
01548     
01549   lim_z = 1 ;   // Juste pour affecter dans espace des configs ;
01550   
01551   Scalar vv_z = vv_bound_cart(om)(3) ;
01552 
01553   for (int k=0 ; k<nnp ; k++)
01554     for (int j=0 ; j<nnt ; j++)
01555       lim_z.set(0, k, j, 0) = vv_z.val_grid_point(1, k, j, 0) ;
01556  
01557   lim_z.set_base(*(mp.get_mg()->std_base_vect_cart()[2])) ;
01558 
01559   return  lim_z ;
01560 }
01561 
01562 // Component x of boundary value of V^i
01563 //-------------------------------------
01564 
01565 const Valeur Isol_hor:: boundary_vv_x_bin(double om, int jj)const {
01566 
01567     int nnp = mp.get_mg()->get_np(1) ;
01568     int nnt = mp.get_mg()->get_nt(1) ;
01569 
01570     //Isol_hor boundary conditions
01571 
01572     Valeur lim_x (mp.get_mg()->get_angu()) ;
01573 
01574     lim_x = 1 ;  // Juste pour affecter dans espace des configs ;
01575 
01576     Scalar vv_x = vv_bound_cart_bin(om, jj)(1) ;
01577 
01578     for (int k=0 ; k<nnp ; k++)
01579     for (int j=0 ; j<nnt ; j++)
01580         lim_x.set(0, k, j, 0) = vv_x.val_grid_point(1, k, j, 0) ;
01581 
01582     lim_x.set_base(*(mp.get_mg()->std_base_vect_cart()[0])) ;
01583 
01584     return  lim_x ;
01585 
01586 
01587 }
01588 
01589 // Component y of boundary value of V^i
01590 //--------------------------------------
01591 
01592 const Valeur Isol_hor:: boundary_vv_y_bin(double om, int jj)const {
01593 
01594     int nnp = mp.get_mg()->get_np(1) ;
01595     int nnt = mp.get_mg()->get_nt(1) ;
01596 
01597     // Isol_hor boundary conditions
01598 
01599     Valeur lim_y (mp.get_mg()->get_angu()) ;
01600 
01601     lim_y = 1 ;  // Juste pour affecter dans espace des configs ;
01602 
01603     Scalar vv_y = vv_bound_cart_bin(om, jj)(2) ;
01604 
01605     for (int k=0 ; k<nnp ; k++)
01606     for (int j=0 ; j<nnt ; j++)
01607         lim_y.set(0, k, j, 0) = vv_y.val_grid_point(1, k, j, 0) ;
01608 
01609     lim_y.set_base(*(mp.get_mg()->std_base_vect_cart()[1])) ;
01610 
01611     return  lim_y ;
01612 }
01613 
01614 
01615 // Component z of boundary value of V^i
01616 //-------------------------------------
01617 
01618 const Valeur Isol_hor:: boundary_vv_z_bin(double om, int jj)const {
01619 
01620     int nnp = mp.get_mg()->get_np(1) ;
01621     int nnt = mp.get_mg()->get_nt(1) ;
01622 
01623     // Isol_hor boundary conditions
01624 
01625     Valeur lim_z (mp.get_mg()->get_angu()) ;
01626 
01627     lim_z = 1 ;   // Juste pour affecter dans espace des configs ;
01628 
01629     Scalar vv_z = vv_bound_cart_bin(om, jj)(3) ;
01630 
01631     for (int k=0 ; k<nnp ; k++)
01632     for (int j=0 ; j<nnt ; j++)
01633         lim_z.set(0, k, j, 0) = vv_z.val_grid_point(1, k, j, 0) ;
01634 
01635     lim_z.set_base(*(mp.get_mg()->std_base_vect_cart()[2])) ;
01636 
01637     return  lim_z ;
01638 }

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