excision_hor.C

00001 /*
00002  *  Definition of methods for the class Spheroid and its subclass App_hor
00003  *
00004  */
00005 
00006 /*
00007  *   Copyright (c) 2009  Jose-Luis Jaramillo & Jerome Novak & Nicolas Vasset
00008  *
00009  *   This file is part of LORENE.
00010  *
00011  *   LORENE is free software; you can redistribute it and/or modify
00012  *   it under the terms of the GNU General Public License version 
00013  *   as published by the Free Software Foundation.
00014  *
00015  *   LORENE is distributed in the hope that it will be useful,
00016  *   but WITHOUT ANY WARRANTY; without even the implied warranty of
00017  *   MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00018  *   GNU General Public License for more details.
00019  *
00020  *   You should have received a copy of the GNU General Public License
00021  *   along with LORENE; if not, write to the Free Software
00022  *   Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00023  *
00024  */
00025 
00026 char excision_hor_C[] = "$Header: /cvsroot/Lorene/C++/Source/App_hor/excision_hor.C,v 1.2 2009/10/23 13:19:15 j_novak Exp $" ;
00027 
00028 /*
00029  * $Header: /cvsroot/Lorene/C++/Source/App_hor/excision_hor.C,v 1.2 2009/10/23 13:19:15 j_novak Exp $
00030  *
00031  */
00032 
00033 // C headers
00034 #include <math.h>
00035 #include <assert.h>
00036 
00037 // Lorene headers
00038 #include "excision_hor.h"
00039 
00040 //---------------//
00041 //  Constructors //
00042 //--------------// 
00043  
00044 
00045 Excision_hor::Excision_hor(const Scalar& h_in, const Metric& gij, const Sym_tensor& Kij2, const Scalar& ppsi, const Scalar& nn, const Vector& beta, const Sym_tensor& Tij2, double timestep, int int_nos):
00046   sph(h_in, gij, Kij2),
00047   conf_fact(ppsi),
00048   lapse(nn),
00049   shift(beta),
00050   gamij (gij),
00051   Kij(Kij2),
00052   delta_t(timestep),
00053   no_of_steps(int_nos),
00054   Tij(Tij2)
00055 {
00056   
00057  set_der_0x0() ;
00058 
00059 }
00060 
00061 
00062 
00063 
00064 
00065 //Copy constructor//
00066 
00067 Excision_hor::Excision_hor(const Excision_hor &exc_in) :sph(exc_in.sph),
00068                                 conf_fact(exc_in.conf_fact),
00069                                 lapse(exc_in.lapse),                                        
00070                                 shift(exc_in.shift),
00071                                 gamij (exc_in.gamij),
00072                                 Kij (exc_in.Kij),
00073                                 delta_t(exc_in.delta_t),
00074                             no_of_steps(exc_in.no_of_steps),
00075                             Tij(exc_in.Tij)
00076                                 
00077 {
00078   set_der_0x0() ; 
00079   
00080 }
00081 //------------//
00082 //Destructor //
00083 //-----------//
00084 
00085 Excision_hor::~Excision_hor()
00086 {
00087   del_deriv() ;
00088 }
00089 
00090 // -----------------//
00091 // Memory management//
00092 //------------------//
00093 void Excision_hor::del_deriv() const {
00094   
00095  if (p_get_BC_conf_fact != 0x0) delete p_get_BC_conf_fact ;
00096  if (p_get_BC_bmN != 0x0) delete p_get_BC_bmN ;
00097  if (p_get_BC_bpN != 0x0) delete p_get_BC_bpN ;
00098  if (p_get_BC_shift != 0x0) delete p_get_BC_shift ;
00099   set_der_0x0() ;
00100 }
00101 
00102 void Excision_hor::set_der_0x0() const {
00103   p_get_BC_conf_fact = 0x0 ;
00104   p_get_BC_bmN = 0x0 ;
00105   p_get_BC_bpN = 0x0 ;
00106   p_get_BC_shift = 0x0 ;
00107 
00108 } 
00109 
00110 
00111  
00112 //---------//
00113 //Accessors//
00114 //---------//
00115 
00116 
00117 
00118 // Source for the Neumann BC on the conformal factor
00119 const Scalar& Excision_hor::get_BC_conf_fact() const{
00120   if (p_get_BC_conf_fact == 0x0){
00121     Sym_tensor gamconfcov = gamij.cov()/pow(conf_fact, 4); 
00122     gamconfcov.std_spectral_base();
00123     Metric gamconf(gamconfcov);
00124     Vector tilde_s = gamconf.radial_vect();
00125     Scalar bound_psi =   -((1./conf_fact)*contract((contract(Kij,1,tilde_s,0)),0, tilde_s,0));    
00126     bound_psi += -conf_fact*tilde_s.divergence(gamconf);
00127     bound_psi = 0.25*bound_psi;
00128     bound_psi += -contract(conf_fact.derive_cov(gamconf),0,tilde_s,0) + conf_fact.dsdr();
00129     bound_psi.std_spectral_base();
00130     bound_psi.set_spectral_va().ylm();    
00131     p_get_BC_conf_fact = new Scalar(bound_psi);
00132 
00133 }
00134   return *p_get_BC_conf_fact ;
00135 }
00136 
00137 
00138 
00139 // Case 0: Source of Dirichlet BC for (b-N), based on an entropy prescription.
00140 // WARNING: the argument value has to be carefully fixed w.r.t initial data for (attempted) continuity.
00141 
00142 //  Case 1: Source of Dirichlet BC for (b-N), from a component of projected Einstein Equations.
00143 //  Requires a 2d poisson solver for a non-flat metric.
00144 
00145 
00146 const Scalar& Excision_hor::get_BC_bmN(int choice_bmN, double value) const{
00147   if (p_get_BC_bmN == 0x0){
00148 
00149     switch(choice_bmN){
00150 
00151     case 0 : {
00152 
00153   Scalar thetaminus = sph.theta_minus();
00154       Scalar theta_minus3 (lapse.get_mp()); 
00155  
00156   theta_minus3.allocate_all();
00157   theta_minus3.std_spectral_base();
00158 
00159   int nz = (*lapse.get_mp().get_mg()).get_nzone();
00160   int nr = (*lapse.get_mp().get_mg()).get_nr(1);
00161   int nt = (*lapse.get_mp().get_mg()).get_nt(1);
00162   int np = (*lapse.get_mp().get_mg()).get_np(1);
00163 
00164 
00165   for (int f= 0; f<nz; f++)
00166     for (int k=0; k<np; k++)
00167       for (int j=0; j<nt; j++) {
00168     for (int l=0; l<nr; l++) {
00169         
00170       theta_minus3.set_grid_point(f,k,j,l) = thetaminus.val_grid_point(0,k,j,0);
00171     
00172     }
00173       }
00174   if (nz >2){
00175      theta_minus3.annule_domain(0);
00176     theta_minus3.annule_domain(nz - 1);
00177    }
00178 
00179 
00180     Scalar bound_bmN(lapse.get_mp()); 
00181     bound_bmN = - value*theta_minus3; bound_bmN.std_spectral_base();
00182     bound_bmN.set_spectral_va().ylm();
00183     p_get_BC_bmN = new Scalar(bound_bmN); 
00184     }
00185 
00186     case 1 : {
00187 
00188     Scalar bound_bmN(lapse.get_mp());
00189     bound_bmN.allocate_all();
00190     bound_bmN.std_spectral_base();
00191 
00192 // Radial vector for the full 3-metric.
00193     Vector sss = gamij.radial_vect();
00194     Vector sss_down = sss.up_down(gamij);
00195     Scalar bb = contract (shift,0, sss_down,0);
00196     Scalar bmN3 = bb - lapse; bmN3.set_spectral_va().ylm();
00197     Scalar bpN3 = bb + lapse; bpN3.set_spectral_va().ylm();  
00198     
00199     int nt = (*lapse.get_mp().get_mg()).get_nt(1);
00200     int np = (*lapse.get_mp().get_mg()).get_np(1);
00201     
00202     Scalar bmN(sph.get_hsurf().get_mp());
00203     bmN.allocate_all();
00204     bmN.std_spectral_base();
00205     bmN.set_spectral_va().ylm();
00206      Scalar bpN(sph.get_hsurf().get_mp());
00207     bpN.allocate_all();
00208     bpN.std_spectral_base();
00209     bpN.set_spectral_va().ylm();
00210 
00211     for (int k=0; k<np; k++)
00212       for (int j=0; j<nt; j++) {
00213 
00214           bmN.set_grid_point(0,k,j,0) = bmN3.val_grid_point(1,k,j,0);
00215       bpN.set_grid_point(0,k,j,0) = bpN3.val_grid_point(1,k,j,0);
00216     }
00217 
00218     Scalar bmN_new(bmN.get_mp());
00219     bmN_new.allocate_all();
00220     bmN_new.std_spectral_base();
00221      
00222     double diff_ent = 1.; 
00223     double precis = 1.e-9;
00224     int mer_max = 200;
00225     double relax = 0.;
00226     for(int mer=0 ;(diff_ent > precis) && (mer<mer_max) ; mer++) {    
00227   
00228       // Calculation of some source terms.
00229      
00230       Scalar hsurf = sph.get_hsurf();
00231       hsurf.set_spectral_va().ylm();
00232       const Metric_flat& fmets = hsurf.get_mp().flat_met_spher() ;
00233  
00234     Scalar shear_up = sph.shear(); shear_up.up_down(sph.get_qab());
00235  
00236     Scalar B_source = 0.5*contract(contract(sph.shear(),0, shear_up, 0),0,1) + 4.*M_PI*Tij.trace(sph.get_qab()); // Redo the matter terms.
00237     Scalar A_source = 0.5*sph.get_ricci() - contract(sph.derive_cov2d(sph.get_ll()), 0, 1) - contract(sph.get_ll(),0, sph.get_ll().up_down(sph.get_qab()),0) - 8.*M_PI*Tij.trace(sph.get_qab()); // Redo the matter terms.
00238    
00239     Scalar op_bmN_add = - 2.*contract(sph.derive_cov2d(bmN),0, sph.get_ll(),0) + A_source*bmN;
00240 
00241     Scalar source_bmN = B_source*bpN - op_bmN_add;  
00242     source_bmN.set_spectral_va().ylm();
00243 
00244     Scalar sqrtqh2 = sph.sqrt_q()*hsurf*hsurf;
00245     sqrtqh2.set_spectral_va().ylm();
00246    
00247     source_bmN = sqrtqh2*source_bmN;
00248      
00249      // Conformal decomposition of the 2-metric
00250         
00251     Sym_tensor qab_con = sph.get_qab().con();
00252     qab_con = qab_con/(hsurf*hsurf); // Renormalization due to the triad still not built-in spheroid class
00253     //This is provisory work.
00254 
00255     // h^ab as q^ab = (f^ab + h^ab) / sqrt_q
00256     Sym_tensor hab =(qab_con*sqrtqh2 - fmets.con()) / (hsurf*hsurf) ;
00257     // for the sake of clarity
00258     hab.set(1,1) = 1. ;
00259     hab.set(1,2) = 0. ;
00260     hab.set(1,3) = 0. ;
00261     hab.std_spectral_base() ;
00262     //end
00263     // Complete source for the angular laplacian.
00264     Scalar d_bmN = sph.derive_cov2dflat(bmN);
00265     d_bmN.set_spectral_va().ylm();
00266       Scalar d2_bmN = sph.derive_cov2dflat(d_bmN);
00267     d2_bmN.set_spectral_va().ylm();
00268  
00269     Scalar source_add = - hsurf*hsurf*contract(hab, 0,1, d2_bmN, 0,1) + sqrtqh2*contract(contract(qab_con,0,1,sph.delta(),1,2),0,d_bmN,0) ;   
00270     source_add.set_spectral_va().ylm();
00271     source_bmN = source_bmN + source_add;
00272     //
00273     
00274     // System inversion
00275     bmN_new = source_bmN.poisson_angu(0.);
00276     
00277     // Actualisation of the principal variable, convergence parameter.
00278     diff_ent = max(maxabs(bmN - bmN_new));
00279     
00280     bmN = relax*bmN + (1. - relax)*bmN_new;   
00281      
00282     }
00283       bound_bmN = bmN;
00284     bound_bmN.set_spectral_va().ylm();
00285           p_get_BC_bmN = new Scalar(bound_bmN); 
00286     }
00287 
00288 }
00289   }
00290   return *p_get_BC_bmN ;
00291 
00292 }
00293 
00294 
00295 // Case 0:  Arbitrary Dirichlet BC for (b+N), fixed by a parabolic driver towards a constant value.
00296 //  Case 1: Source of Dirichlet BC for (b+N), from a component of projected Einstein Equations.
00297 const Scalar& Excision_hor::get_BC_bpN(int choice_bpN, double c_bpn_lap, double c_bpn_fin, Scalar *bpN_fin) const{
00298   if (p_get_BC_bpN == 0x0){
00299 
00300     switch(choice_bpN) {
00301  
00302     case 0 : {
00303 
00304       Vector sss = gamij.radial_vect();
00305      Vector sss_down = sss.up_down(gamij);
00306       Scalar bb = contract (shift,0, sss_down,0);
00307   Scalar bpN = bb + lapse;
00308   Scalar ff = lapse*(c_bpn_lap*bpN.lapang() + c_bpn_fin*(bpN- *bpN_fin));
00309   ff.std_spectral_base();
00310  
00311 
00312   // Definition of k_1
00313   Scalar k_1 =delta_t*ff;
00314    
00315   // Intermediate value of b-N, for Runge-Kutta 2nd order scheme
00316   Scalar bpN_int = bpN + k_1; bpN_int.std_spectral_base();
00317 
00318   // Recalculation of ff with intermediate values. 
00319   Scalar ff_int =  lapse*(c_bpn_lap*bpN_int.lapang() + c_bpn_fin*bpN_int);
00320  
00321   // Definition of k_2
00322   Scalar k_2 = delta_t*ff_int; 
00323   k_2.std_spectral_base();
00324  
00325   // Result of RK2 evolution
00326   Scalar bound_bpN =  bpN + k_2;
00327   bound_bpN.std_spectral_base();
00328   bound_bpN.set_spectral_va().ylm();
00329 
00330         p_get_BC_bpN = new Scalar(bound_bpN); 
00331     
00332 
00333     }
00334 
00335     case 1 : {
00336 
00337     
00338   int nz = (*lapse.get_mp().get_mg()).get_nzone();
00339   int nr = (*lapse.get_mp().get_mg()).get_nr(1);
00340   int nt = (*lapse.get_mp().get_mg()).get_nt(1);
00341   int np = (*lapse.get_mp().get_mg()).get_np(1);
00342 
00343     
00344     
00345   Scalar bmN3 = get_BC_bmN(0, 1.); // change the argument.
00346     
00347     Scalar bmN(sph.get_hsurf().get_mp());
00348     bmN.allocate_all();
00349     bmN.std_spectral_base();
00350     bmN.set_spectral_va().ylm();
00351  
00352     
00353 
00354     for (int k=0; k<np; k++)
00355       for (int j=0; j<nt; j++) {
00356 
00357         
00358       bmN.set_grid_point(0,k,j,0) = bmN3.val_grid_point(1,k,j,0);
00359     
00360     }
00361 
00362 
00363     Scalar bound_bpN(lapse.get_mp());
00364     bound_bpN.allocate_all();
00365     bound_bpN.std_spectral_base();
00366  
00367     // Definition of source terms in relation (6) of Jaramillo et al. 2007.
00368     // All is done on the spheroid of radius r=1.
00369  
00370     Scalar shear_up = sph.shear(); shear_up.up_down(sph.get_qab());
00371  
00372     Scalar B_source = 0.5*contract(contract(sph.shear(),0, shear_up, 0),0,1) + 4.*M_PI*Tij.trace(sph.get_qab()); // Redo the matter terms.
00373     Scalar A_source = 0.5*sph.get_ricci() - contract(sph.derive_cov2d(sph.get_ll()), 0, 1) - contract(sph.get_ll(),0, sph.get_ll().up_down(sph.get_qab()),0) - 8.*M_PI*Tij.trace(sph.get_qab()); // Redo the matter terms.
00374   
00375     // Curved 2d Laplacian of (b -N).
00376    
00377     Sym_tensor interlap = sph.derive_cov2d(sph.derive_cov2d(bmN)); 
00378     interlap.up(0,sph.get_qab()); 
00379     Sym_tensor lap_bmN = contract(interlap,0,1);
00380  
00381     Scalar op_bmN = lap_bmN - 2.*contract(sph.derive_cov2d(bmN),0, sph.get_ll(),0) + A_source*bmN;
00382           
00383    Scalar  bound_bpN2 = op_bmN/B_source;
00384     bound_bpN2.std_spectral_base();
00385     bound_bpN2.set_spectral_va().ylm();
00386 
00387 
00388   for (int f= 0; f<nz; f++)
00389     for (int k=0; k<np; k++)
00390       for (int j=0; j<nt; j++) {
00391     for (int l=0; l<nr; l++) {
00392         
00393       bound_bpN.set_grid_point(f,k,j,l) = bound_bpN2.val_grid_point(0,k,j,0);
00394     
00395     }
00396       }
00397   if (nz >2){
00398      bound_bpN.annule_domain(0);
00399     bound_bpN.annule_domain(nz - 1);
00400    }
00401  
00402  
00403 
00404         p_get_BC_bpN = new Scalar(bound_bpN); 
00405     
00406  
00407 }
00408 
00409 }
00410 }
00411   return *p_get_BC_bpN ;
00412 }
00413 
00414 
00415 
00416 
00417 // Source for the Dirichlet BC on the shift
00418 // The tangential shift is fixed using a parabolic driver based on the conformal Killing equation in the dynamical case.
00419 
00420 const Vector& Excision_hor::get_BC_shift( double c_V_lap ) const{
00421   if (p_get_BC_shift == 0x0){
00422 
00423     // Radial vector for the full 3-metric.
00424      Vector sss = gamij.radial_vect();
00425      Vector sss_down = sss.up_down(gamij);
00426      
00427 //     // Boundary value for the radial part of the shift: parabolic driver for (b-N)
00428      //  Scalar bound = lapse ; 
00429      Scalar bb = 0.5*(*p_get_BC_bpN + *p_get_BC_bmN) ; // TO CHANGE: additional function?-> put choice-bb
00430   
00431 
00432   // Tangent part of the shift, with parabolic driver
00433   
00434   
00435   Vector V_par = shift - bb*sss;
00436   Sym_tensor q_upup = gamij.con() - sss*sss;
00437 
00438   
00439   // Calculation of the conformal 2d laplacian of V
00440   Tensor q_updown = q_upup.down(1, gamij); 
00441   Tensor dd_V = V_par.derive_con(gamij);
00442   dd_V = contract(q_updown, 1, contract(q_updown,1 ,dd_V, 0), 1);
00443   Vector lap_V = contract(q_updown, 1, contract(dd_V.derive_cov(gamij),1,2), 0);
00444   
00445   // 3d interpolation of the Ricci scalar on the surface.
00446   
00447   Scalar ricci2 = sph.get_ricci();
00448   
00449      // Start Mapping interpolation
00450  
00451       Scalar ricci3 (lapse.get_mp()); 
00452  
00453   ricci3.allocate_all();
00454   ricci3.std_spectral_base();
00455 
00456   int nz = (*lapse.get_mp().get_mg()).get_nzone();
00457   int nr = (*lapse.get_mp().get_mg()).get_nr(1);
00458   int nt = (*lapse.get_mp().get_mg()).get_nt(1);
00459   int np = (*lapse.get_mp().get_mg()).get_np(1);
00460 
00461 
00462   for (int f= 0; f<nz; f++)
00463     for (int k=0; k<np; k++)
00464       for (int j=0; j<nt; j++) {
00465     for (int l=0; l<nr; l++) {
00466         
00467       ricci3.set_grid_point(f,k,j,l) = ricci2.val_grid_point(0,k,j,0);
00468 
00469     
00470     }
00471       }
00472   if (nz >2){
00473      ricci3.annule_domain(0);
00474     ricci3.annule_domain(nz - 1);
00475 
00476   }
00477     // End Mapping interpolation
00478   
00479     // Construction of the Ricci COV tensor on the sphere
00480  
00481   Sym_tensor ricci_t = gamij.cov() - sss_down*sss_down;
00482   ricci_t = 0.5*ricci3*ricci_t;
00483   ricci_t.std_spectral_base();
00484  
00485   Tensor ricci_t_updown = contract(q_upup,0, ricci_t,0); 
00486   
00487   // Calculation of ff 
00488 
00489   Vector ffV = c_V_lap*lapse*(lap_V + contract(ricci_t_updown,1, V_par,0));
00490   ffV.std_spectral_base();
00491 
00492 
00493   // Definition of k_1
00494   Vector k_1V =delta_t*ffV;
00495    
00496   // Intermediate value of Npsi, for Runge-Kutta 2nd order scheme
00497   if (nz >2){
00498     k_1V.annule_domain(nz-1);
00499   }                             // Patch to avoid dzpuis problems if existent.
00500   Vector V_par_int = V_par + k_1V;// V_par_int.std_spectral_base();
00501 
00502   // Recalculation of ff with intermediate values.
00503 
00504   Sym_tensor dd_V_int = V_par_int.derive_con(gamij);
00505   dd_V_int = contract(q_updown, 1, contract(q_updown,1 ,dd_V_int, 0), 1);
00506   Vector lap_V_int = contract(q_updown, 1, contract(dd_V_int.derive_cov(gamij),1,2), 0);
00507  
00508   Vector ffV_int =  c_V_lap*lapse*(lap_V_int + contract(ricci_t_updown,1, V_par_int,0));
00509  
00510   // Definition of k_2
00511   Vector k_2V = delta_t*ffV_int; 
00512   //  k_2.std_spectral_base();
00513  
00514   // Result of RK2 evolution
00515   if (nz >2){
00516     k_2V.annule_domain(nz-1);
00517   }
00518   Vector bound_V = V_par + k_2V;
00519   //  bound_V.std_spectral_base();
00520 
00521        // Construction of the total shift boundary condition
00522        Vector bound_shift = bb*sss + bound_V;
00523        bound_shift.std_spectral_base();
00524        p_get_BC_shift = new Vector(bound_shift);
00525 }
00526   return *p_get_BC_shift ;
00527 }
00528  
00529 
00530 
00531  void Excision_hor::sauve(FILE* ) const {
00532 
00533    cout << "c'est pas fait!" << endl ;
00534    return ; 
00535 
00536  }

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