time_slice_conf.C

00001 /*
00002  *  Methods of class Time_slice_conf
00003  *
00004  *    (see file time_slice.h for documentation).
00005  *
00006  */
00007 
00008 /*
00009  *   Copyright (c) 2004 Eric Gourgoulhon, Jose Luis Jaramillo & Jerome Novak
00010  *
00011  *   This file is part of LORENE.
00012  *
00013  *   LORENE is free software; you can redistribute it and/or modify
00014  *   it under the terms of the GNU General Public License version 2
00015  *   as published by the Free Software Foundation.
00016  *
00017  *   LORENE is distributed in the hope that it will be useful,
00018  *   but WITHOUT ANY WARRANTY; without even the implied warranty of
00019  *   MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00020  *   GNU General Public License for more details.
00021  *
00022  *   You should have received a copy of the GNU General Public License
00023  *   along with LORENE; if not, write to the Free Software
00024  *   Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00025  *
00026  */
00027 
00028 char time_slice_conf_C[] = "$Header: /cvsroot/Lorene/C++/Source/Time_slice/time_slice_conf.C,v 1.16 2008/12/04 18:22:49 j_novak Exp $" ;
00029 
00030 /*
00031  * $Id: time_slice_conf.C,v 1.16 2008/12/04 18:22:49 j_novak Exp $
00032  * $Log: time_slice_conf.C,v $
00033  * Revision 1.16  2008/12/04 18:22:49  j_novak
00034  * Enhancement of the dzpuis treatment + various bug fixes.
00035  *
00036  * Revision 1.15  2008/12/02 15:02:22  j_novak
00037  * Implementation of the new constrained formalism, following Cordero et al. 2009
00038  * paper. The evolution eqs. are solved as a first-order system. Not tested yet!
00039  *
00040  * Revision 1.14  2004/06/24 20:36:54  e_gourgoulhon
00041  * Added method check_psi_dot.
00042  *
00043  * Revision 1.13  2004/05/31 09:08:18  e_gourgoulhon
00044  * Method sauve and constructor from binary file are now operational.
00045  *
00046  * Revision 1.12  2004/05/27 15:25:04  e_gourgoulhon
00047  * Added constructors from binary file, as well as corresponding
00048  * functions sauve and save.
00049  *
00050  * Revision 1.11  2004/05/12 15:24:20  e_gourgoulhon
00051  * Reorganized the #include 's, taking into account that
00052  * time_slice.h contains now an #include "metric.h".
00053  *
00054  * Revision 1.10  2004/05/10 09:10:05  e_gourgoulhon
00055  * Added "adm_mass_evol.downdate(jtime)" in methods set_*.
00056  *
00057  * Revision 1.9  2004/05/05 14:31:14  e_gourgoulhon
00058  * Method aa(): added *as a comment*  annulation of hh_point in the compactified
00059  * domain.
00060  *
00061  * Revision 1.8  2004/05/03 14:47:11  e_gourgoulhon
00062  * Corrected method aa().
00063  *
00064  * Revision 1.7  2004/04/08 16:43:26  e_gourgoulhon
00065  * Added methods set_*
00066  * Added test of determinant one in constructor and set_hh.
00067  *
00068  * Revision 1.6  2004/04/05 21:25:02  e_gourgoulhon
00069  * -- Added constructor as standard time slice of Minkowski spacetime.
00070  * -- Added some calls to Scalar::std_spectral_base() after
00071  *    non-arithmetical operations.
00072  *
00073  * Revision 1.5  2004/04/05 12:38:45  j_novak
00074  * Minor modifs to prevent some warnings.
00075  *
00076  * Revision 1.4  2004/04/01 16:09:02  j_novak
00077  * Trace of K_ij is now member of Time_slice (it was member of Time_slice_conf).
00078  * Added new methods for checking 3+1 Einstein equations (preliminary).
00079  *
00080  * Revision 1.3  2004/03/29 12:00:41  e_gourgoulhon
00081  * Many modifs.
00082  *
00083  * Revision 1.2  2004/03/28 21:32:23  e_gourgoulhon
00084  * Corrected error in method trk().
00085  *
00086  * Revision 1.1  2004/03/28 21:30:13  e_gourgoulhon
00087  * First version.
00088  *
00089  *
00090  * $Header: /cvsroot/Lorene/C++/Source/Time_slice/time_slice_conf.C,v 1.16 2008/12/04 18:22:49 j_novak Exp $
00091  *
00092  */
00093 
00094 // C headers
00095 #include <assert.h>
00096 
00097 // Lorene headers
00098 #include "time_slice.h"
00099 #include "utilitaires.h"
00100 
00101 
00102 
00103                 //--------------//
00104                 // Constructors //
00105                 //--------------//
00106 
00107 
00108 // Constructor from conformal decomposition
00109 // ----------------------------------------
00110 
00111 Time_slice_conf::Time_slice_conf(const Scalar& lapse_in, const Vector& shift_in,
00112             const Metric_flat& ff_in, const Scalar& psi_in, 
00113             const Sym_tensor& hh_in, const Sym_tensor& hata_in, 
00114             const Scalar& trk_in, int depth_in) 
00115                     : Time_slice(depth_in),
00116                       ff(ff_in),
00117                       psi_evol(psi_in, depth_in), 
00118                       npsi_evol(depth_in),
00119                       hh_evol(hh_in, depth_in), 
00120                       hata_evol(hata_in, depth_in),
00121               A_hata_evol(depth_in), B_hata_evol(depth_in){
00122 
00123     assert(hh_in.get_index_type(0) == CON) ; 
00124     assert(hh_in.get_index_type(1) == CON) ; 
00125     assert(hata_in.get_index_type(0) == CON) ; 
00126     assert(hata_in.get_index_type(1) == CON) ; 
00127 
00128     double time_init = the_time[jtime] ; 
00129 
00130     // Check whether det tgam^{ij} = det f^{ij} :
00131     // ----------------------------------------
00132     Sym_tensor tgam_in = ff_in.con() + hh_in ; 
00133     
00134     Scalar det_in = tgam_in(1, 1)*tgam_in(2, 2)*tgam_in(3, 3) 
00135         + tgam_in(1, 2)*tgam_in(2, 3)*tgam_in(3, 1)
00136         + tgam_in(1, 3)*tgam_in(2, 1)*tgam_in(3, 2) 
00137         - tgam_in(3, 1)*tgam_in(2, 2)*tgam_in(1, 3)
00138         - tgam_in(3, 2)*tgam_in(2, 3)*tgam_in(1, 1) 
00139         - tgam_in(3, 3)*tgam_in(2, 1)*tgam_in(1, 2) ;
00140     
00141     double diffdet = max(maxabs(det_in - 1. / ff.determinant(), 
00142         "Deviation of det tgam^{ij} from 1/f")) ;
00143     if ( diffdet > 1.e-13 ) {
00144         cerr << 
00145         "Time_slice_conf::Time_slice_conf : the input h^{ij} does not"
00146         << " ensure \n" << "  det tgam_{ij} = f  ! \n" 
00147         << "  error = " << diffdet << endl ; 
00148         abort() ; 
00149     }
00150           
00151     n_evol.update(lapse_in, jtime, time_init) ;     
00152     beta_evol.update(shift_in, jtime, time_init) ; 
00153     trk_evol.update(trk_in, jtime, time_init) ;
00154     A_hata() ;
00155     B_hata() ;
00156     
00157     set_der_0x0() ;  
00158     
00159 }
00160                  
00161 
00162 // Constructor from physical metric
00163 // --------------------------------                 
00164 
00165 Time_slice_conf::Time_slice_conf(const Scalar& lapse_in, const Vector& shift_in,
00166                const Sym_tensor& gamma_in, const Sym_tensor& kk_in,
00167                const Metric_flat& ff_in, int depth_in) 
00168                     : Time_slice(lapse_in, shift_in, gamma_in, kk_in, depth_in),
00169                       ff(ff_in),
00170                       psi_evol(depth_in), 
00171                       npsi_evol(depth_in),
00172                       hh_evol(depth_in), 
00173                       hata_evol(depth_in),
00174               A_hata_evol(depth_in), B_hata_evol(depth_in) {
00175                     
00176     set_der_0x0() ; // put here in order not to erase p_psi4
00177 
00178     double time_init = the_time[jtime] ; 
00179     
00180     assert( p_gamma != 0x0 ) ; 
00181     p_psi4 = new Scalar( pow( p_gamma->determinant() / ff.determinant(), 
00182                          0.3333333333333333) ) ;
00183     p_psi4->std_spectral_base() ;
00184 
00185     Scalar tmp = pow(*p_psi4, 0.25) ;
00186     tmp.std_spectral_base() ;
00187     psi_evol.update(tmp , jtime, time_init ) ; 
00188     
00189     hh_evol.update( (*p_psi4) * p_gamma->con() - ff.con(), 
00190                     jtime, time_init ) ; 
00191     
00192     hata_evol.update( tmp*tmp*(*p_psi4)*(*p_psi4) *( Time_slice::k_uu() 
00193                     - 0.3333333333333333 * trk_evol[jtime] * p_gamma->con() ), 
00194                     jtime, time_init ) ;                   
00195     A_hata() ;
00196     B_hata() ;
00197 }
00198 
00199 // Constructor as standard time slice of flat spacetime (Minkowski)
00200 // ----------------------------------------------------------------
00201 
00202 Time_slice_conf::Time_slice_conf(const Map& mp, const Base_vect& triad, 
00203                                  const Metric_flat& ff_in, int depth_in) 
00204     : Time_slice(mp, triad, depth_in),
00205       ff(ff_in), 
00206       psi_evol(depth_in), 
00207       npsi_evol(depth_in),
00208       hh_evol(depth_in), 
00209       hata_evol(depth_in),
00210       A_hata_evol(depth_in), B_hata_evol(depth_in) {
00211     
00212     double time_init = the_time[jtime] ; 
00213     
00214     // Psi identically one:
00215     Scalar tmp(mp) ; 
00216     tmp.set_etat_one() ; 
00217     tmp.std_spectral_base() ;
00218     psi_evol.update(tmp, jtime, time_init) ; 
00219     
00220     // N Psi identically one:
00221     npsi_evol.update(tmp, jtime, time_init) ; 
00222     
00223     // h^{ij} identically zero:
00224     Sym_tensor stmp(mp, CON, triad) ; 
00225     stmp.set_etat_zero() ; 
00226     hh_evol.update(stmp, jtime, time_init) ; 
00227     
00228     // \hat{A}^{ij} identically zero:
00229     hata_evol.update(stmp, jtime, time_init) ; 
00230 
00231     tmp.set_etat_zero() ;
00232     A_hata_evol.update(tmp, jtime, time_init) ;
00233     B_hata_evol.update(tmp, jtime, time_init) ;
00234 
00235     set_der_0x0() ; 
00236 
00237 }
00238 
00239 
00240 // Constructor from binary file             
00241 // ----------------------------
00242 
00243 Time_slice_conf::Time_slice_conf(const Map& mp, const Base_vect& triad, 
00244                         const Metric_flat& ff_in, FILE* fich, 
00245                         bool partial_read, int depth_in) 
00246         : Time_slice(mp, triad, fich, true, depth_in),
00247           ff(ff_in), 
00248           psi_evol(depth_in), 
00249           npsi_evol(depth_in),
00250           hh_evol(depth_in), 
00251           hata_evol(depth_in),
00252       A_hata_evol(depth_in), B_hata_evol(depth_in) {
00253 
00254     // Put here, not to destroy p_vec_X
00255     set_der_0x0() ; 
00256 
00257     // Reading of various fields
00258     // -------------------------
00259     
00260     int jmin = jtime - depth + 1 ; 
00261     int indicator ; 
00262 
00263     // Psi
00264     for (int j=jmin; j<=jtime; j++) {
00265         fread_be(&indicator, sizeof(int), 1, fich) ;
00266         if (indicator == 1) {
00267             Scalar psi_file(mp, *(mp.get_mg()), fich) ; 
00268             psi_evol.update(psi_file, j, the_time[j]) ; 
00269         }
00270     } 
00271     // hat{A}^{ij}
00272     for (int j=jmin; j<=jtime; j++) {
00273     fread_be(&indicator, sizeof(int), 1, fich) ;    
00274     if (indicator == 1) {
00275         Sym_tensor hat_A_file(mp, triad, fich) ; 
00276         hata_evol.update( hat_A_file, j, the_time[j] ) ; 
00277         }
00278     }
00279 
00280     //A and B...
00281     for (int j=jmin; j<=jtime; j++) {
00282         fread_be(&indicator, sizeof(int), 1, fich) ;
00283         if (indicator == 1) {
00284             Scalar A_file(mp, *(mp.get_mg()), fich) ; 
00285             A_hata_evol.update(A_file, j, the_time[j]) ; 
00286         }
00287     } 
00288     for (int j=jmin; j<=jtime; j++) {
00289         fread_be(&indicator, sizeof(int), 1, fich) ;
00290         if (indicator == 1) {
00291             Scalar B_file(mp, *(mp.get_mg()), fich) ; 
00292             B_hata_evol.update(B_file, j, the_time[j]) ; 
00293         }
00294     } 
00295 
00296     // Case of a full reading
00297     // -----------------------
00298     if (!partial_read) {
00299         cout << 
00300         "Time_slice constructor from file: the case of full reading\n"
00301         << " is not ready yet !" << endl ; 
00302         abort() ; 
00303     }
00304 
00305 } 
00306 
00307 
00308 
00309 
00310 // Copy constructor
00311 // ----------------
00312 
00313 Time_slice_conf::Time_slice_conf(const Time_slice_conf& tin) 
00314                     : Time_slice(tin), 
00315                       ff(tin.ff),
00316                       psi_evol(tin.psi_evol), 
00317                       npsi_evol(tin.npsi_evol),
00318                       hh_evol(tin.hh_evol), 
00319                       hata_evol(tin.hata_evol),
00320               A_hata_evol(tin.A_hata_evol),
00321               B_hata_evol(tin.B_hata_evol){
00322 
00323     set_der_0x0() ; 
00324                        
00325 }
00326                 //--------------//
00327                 //  Destructor  //
00328                 //--------------//
00329 
00330 Time_slice_conf::~Time_slice_conf(){
00331 
00332     Time_slice_conf::del_deriv() ; 
00333 
00334 }
00335 
00336                     //---------------------//
00337                     //  Memory management  //
00338                     //---------------------//
00339 
00340 void Time_slice_conf::del_deriv() const {
00341 
00342     if (p_tgamma != 0x0) delete p_tgamma ; 
00343     if (p_psi4 != 0x0) delete p_psi4 ; 
00344     if (p_ln_psi != 0x0) delete p_ln_psi ; 
00345     if (p_hdirac != 0x0) delete p_hdirac ; 
00346     if (p_vec_X != 0x0) delete p_vec_X ;
00347     
00348     set_der_0x0() ;
00349 
00350     Time_slice::del_deriv() ; 
00351 }
00352 
00353 
00354 void Time_slice_conf::set_der_0x0() const {
00355 
00356     p_tgamma = 0x0 ; 
00357     p_psi4 = 0x0 ; 
00358     p_ln_psi = 0x0 ; 
00359     p_hdirac = 0x0 ;
00360     p_vec_X = 0x0 ;
00361     
00362 }
00363 
00364 
00365                     //-----------------------//
00366                     // Mutators / assignment //
00367                     //-----------------------//
00368 
00369 void Time_slice_conf::operator=(const Time_slice_conf& tin) {
00370 
00371     Time_slice::operator=(tin) ; 
00372 
00373     psi_evol = tin.psi_evol ; 
00374     npsi_evol = tin.npsi_evol ; 
00375     hh_evol = tin.hh_evol ; 
00376     hata_evol = tin.hata_evol ; 
00377     A_hata_evol = tin.A_hata_evol ;
00378     B_hata_evol = tin.B_hata_evol ;
00379        
00380     del_deriv() ; 
00381     
00382 }
00383 
00384 void Time_slice_conf::operator=(const Time_slice& tin) {
00385 
00386     Time_slice::operator=(tin) ; 
00387 
00388     cerr << 
00389     "Time_slice_conf::operator=(const Time_slice& ) : not implemented yet !"
00390         << endl ;
00391     abort() ;       
00392     del_deriv() ; 
00393     
00394 }
00395 
00396 
00397 void Time_slice_conf::set_psi_del_npsi(const Scalar& psi_in) {
00398 
00399     psi_evol.update(psi_in, jtime, the_time[jtime]) ; 
00400 
00401     // Reset of quantities depending on Psi:
00402     if (p_psi4 != 0x0) {
00403         delete p_psi4 ; 
00404         p_psi4 = 0x0 ; 
00405     }
00406     if (p_ln_psi != 0x0) {
00407         delete p_ln_psi ; 
00408         p_ln_psi = 0x0 ; 
00409     }
00410     if (p_gamma != 0x0) {
00411         delete p_gamma ; 
00412         p_gamma = 0x0 ; 
00413     }
00414     npsi_evol.downdate(jtime) ; 
00415     gam_dd_evol.downdate(jtime) ; 
00416     gam_uu_evol.downdate(jtime) ; 
00417     adm_mass_evol.downdate(jtime) ;  
00418      
00419 }
00420 
00421 void Time_slice_conf::set_psi_del_n(const Scalar& psi_in) {
00422 
00423     psi_evol.update(psi_in, jtime, the_time[jtime]) ; 
00424 
00425     // Reset of quantities depending on Psi:
00426     if (p_psi4 != 0x0) {
00427         delete p_psi4 ; 
00428         p_psi4 = 0x0 ; 
00429     }
00430     if (p_ln_psi != 0x0) {
00431         delete p_ln_psi ; 
00432         p_ln_psi = 0x0 ; 
00433     }
00434     if (p_gamma != 0x0) {
00435         delete p_gamma ; 
00436         p_gamma = 0x0 ; 
00437     }
00438     n_evol.downdate(jtime) ; 
00439     gam_dd_evol.downdate(jtime) ; 
00440     gam_uu_evol.downdate(jtime) ; 
00441     adm_mass_evol.downdate(jtime) ;  
00442      
00443 }
00444 
00445 
00446 void Time_slice_conf::set_npsi_del_psi(const Scalar& npsi_in) {
00447 
00448     npsi_evol.update(npsi_in, jtime, the_time[jtime]) ; 
00449 
00450     // Reset of quantities depending on N Psi:
00451     psi_evol.downdate(jtime) ; 
00452     if (p_psi4 != 0x0) {
00453         delete p_psi4 ; 
00454         p_psi4 = 0x0 ; 
00455     }
00456     if (p_ln_psi != 0x0) {
00457         delete p_ln_psi ; 
00458         p_ln_psi = 0x0 ; 
00459     }
00460     if (p_gamma != 0x0) {
00461         delete p_gamma ; 
00462         p_gamma = 0x0 ; 
00463     }
00464     gam_dd_evol.downdate(jtime) ; 
00465     gam_uu_evol.downdate(jtime) ; 
00466     adm_mass_evol.downdate(jtime) ;  
00467      
00468 }
00469 
00470 
00471 void Time_slice_conf::set_npsi_del_n(const Scalar& npsi_in) {
00472 
00473     npsi_evol.update(npsi_in, jtime, the_time[jtime]) ; 
00474 
00475     // Reset of quantities depending on N Psi:
00476     n_evol.downdate(jtime) ; 
00477     adm_mass_evol.downdate(jtime) ;  
00478      
00479 }
00480 
00481 
00482 void Time_slice_conf::set_hh(const Sym_tensor& hh_in) {
00483 
00484     // Check whether det tgam^{ij} = det f^{ij} :
00485     // ----------------------------------------
00486     Sym_tensor tgam_in = ff.con() + hh_in ; 
00487     
00488     Scalar det_in = tgam_in(1, 1)*tgam_in(2, 2)*tgam_in(3, 3) 
00489         + tgam_in(1, 2)*tgam_in(2, 3)*tgam_in(3, 1)
00490         + tgam_in(1, 3)*tgam_in(2, 1)*tgam_in(3, 2) 
00491         - tgam_in(3, 1)*tgam_in(2, 2)*tgam_in(1, 3)
00492         - tgam_in(3, 2)*tgam_in(2, 3)*tgam_in(1, 1) 
00493         - tgam_in(3, 3)*tgam_in(2, 1)*tgam_in(1, 2) ;
00494     
00495     double diffdet = max(maxabs(det_in - 1. / ff.determinant(), 
00496         "Deviation of det tgam^{ij} from 1/f")) ;
00497     if ( diffdet > 1.e-13 ) {
00498         cerr << 
00499         "Time_slice_conf::set_hh : the input h^{ij} does not"
00500         << " ensure \n" << "  det tgam_{ij} = f  ! \n" 
00501         << "  error = " << diffdet << endl ; 
00502         abort() ; 
00503     }
00504           
00505     hh_evol.update(hh_in, jtime, the_time[jtime]) ; 
00506     
00507     // Reset of quantities depending on h^{ij}:
00508     if (p_tgamma != 0x0) {
00509         delete p_tgamma ;
00510         p_tgamma = 0x0 ; 
00511     } 
00512     if (p_hdirac != 0x0) {
00513         delete p_hdirac ; 
00514         p_hdirac = 0x0 ; 
00515     }
00516     if (p_gamma != 0x0) {
00517         delete p_gamma ; 
00518         p_gamma = 0x0 ;
00519     }
00520     gam_dd_evol.downdate(jtime) ; 
00521     gam_uu_evol.downdate(jtime) ; 
00522     adm_mass_evol.downdate(jtime) ;  
00523      
00524 }
00525 
00526 
00527 void Time_slice_conf::set_hata(const Sym_tensor& hata_in) {
00528 
00529     hata_evol.update(hata_in, jtime, the_time[jtime]) ; 
00530 
00531     // Reset of quantities depending on A^{ij}:
00532     A_hata_evol.downdate(jtime) ;
00533     B_hata_evol.downdate(jtime) ;
00534     k_dd_evol.downdate(jtime) ; 
00535     k_uu_evol.downdate(jtime) ; 
00536 
00537 }
00538 
00539 void Time_slice_conf::set_hata_TT(const Sym_tensor_tt& hata_tt) {
00540     
00541     Scalar tmp = hata_tt.compute_A(true) ;
00542     if (tmp.get_dzpuis() == 3)
00543     tmp.dec_dzpuis() ; // dzpuis 3->2
00544 
00545     A_hata_evol.update( tmp, jtime, the_time[jtime] ) ;
00546 
00547     tmp = hata_tt.compute_tilde_B_tt(true) ;
00548     if (tmp.get_dzpuis() == 3)
00549     tmp.dec_dzpuis() ; // dzpuis 3->2
00550 
00551     B_hata_evol.update( tmp, jtime, the_time[jtime] ) ;
00552 
00553     hata_evol.downdate(jtime) ;
00554     // Reset of quantities depending on A^{ij}:
00555     k_dd_evol.downdate(jtime) ; 
00556     k_uu_evol.downdate(jtime) ; 
00557 }
00558 
00559 void Time_slice_conf::set_hata_from_XAB(Param* par_bc, Param* par_mat) {
00560 
00561     assert (p_vec_X != 0x0) ;
00562     assert (A_hata_evol.is_known(jtime)) ;
00563     assert (B_hata_evol.is_known(jtime)) ;
00564 
00565     const Map& mp = p_vec_X->get_mp() ;
00566 
00567     Sym_tensor_tt hata_tt(mp, mp.get_bvect_spher(), ff) ;
00568     hata_tt.set_A_tildeB(A_hata_evol[jtime], B_hata_evol[jtime], par_bc, par_mat) ;
00569     hata_tt.inc_dzpuis(2) ;
00570 
00571     hata_evol.update( hata_tt + p_vec_X->ope_killing_conf(ff), jtime, the_time[jtime]) ;
00572 
00573     // Reset of quantities depending on A^{ij}:
00574     k_dd_evol.downdate(jtime) ; 
00575     k_uu_evol.downdate(jtime) ; 
00576     
00577 }
00578 
00579 
00580                 //-----------------------------------------------//
00581                 //  Update of fields from base class Time_slice  //
00582                 //-----------------------------------------------//
00583 
00584 const Scalar& Time_slice_conf::nn() const {
00585 
00586     if (!( n_evol.is_known(jtime) ) ) {
00587 
00588         assert( psi_evol.is_known(jtime) ) ; 
00589         assert( npsi_evol.is_known(jtime) ) ; 
00590         
00591         n_evol.update( npsi_evol[jtime] / psi_evol[jtime] , 
00592                         jtime, the_time[jtime] ) ; 
00593     }
00594 
00595     return n_evol[jtime] ;
00596 
00597 } 
00598 
00599 
00600 
00601 const Sym_tensor& Time_slice_conf::gam_dd() const {
00602 
00603     if (!( gam_dd_evol.is_known(jtime)) ) {
00604         gam_dd_evol.update( psi4() * tgam().cov(), jtime, the_time[jtime] ) ; 
00605     }
00606 
00607     return gam_dd_evol[jtime] ;
00608 
00609 }
00610 
00611 
00612 const Sym_tensor& Time_slice_conf::gam_uu() const {
00613 
00614     if (!( gam_uu_evol.is_known(jtime)) ) {
00615         gam_uu_evol.update( tgam().con() / psi4() , jtime, the_time[jtime] ) ; 
00616     }
00617 
00618     return gam_uu_evol[jtime] ;
00619 
00620 }
00621 
00622 
00623 const Sym_tensor& Time_slice_conf::k_dd() const {
00624 
00625     if ( ! (k_dd_evol.is_known(jtime)) ) {
00626        
00627         k_dd_evol.update( k_uu().up_down(gam()), jtime, the_time[jtime] ) ; 
00628         
00629     }
00630 
00631     return k_dd_evol[jtime] ;
00632 
00633 }
00634 
00635 
00636 const Sym_tensor& Time_slice_conf::k_uu() const {
00637 
00638     if ( ! (k_uu_evol.is_known(jtime)) ) {
00639        
00640         k_uu_evol.update( hata()/(psi4()*psi4()*psi()*psi()) 
00641               + 0.3333333333333333* trk()* gam().con(),
00642               jtime, the_time[jtime] ) ; 
00643     }
00644 
00645     return k_uu_evol[jtime] ;
00646 
00647 }
00648 
00649 
00650 
00651 
00652 
00653                 //-----------------------------------//
00654                 //  Update of fields from this class //
00655                 //-----------------------------------//
00656 
00657 const Scalar& Time_slice_conf::A_hata() const {
00658 
00659     if ( !( A_hata_evol.is_known(jtime) ) ) {
00660 
00661     assert( hata_evol.is_known(jtime) ) ;
00662     Scalar tmp = hata_evol[jtime].compute_A(true) ;
00663     if (tmp.get_dzpuis() == 3)
00664         tmp.dec_dzpuis() ; // dzpuis 3->2
00665 
00666     A_hata_evol.update( tmp, jtime, the_time[jtime] ) ;
00667     }
00668     return A_hata_evol[jtime] ;
00669 }
00670 
00671 const Scalar& Time_slice_conf::B_hata() const {
00672 
00673     if ( !( B_hata_evol.is_known(jtime) ) ) {
00674 
00675     assert( hata_evol.is_known(jtime) ) ;
00676     Scalar tmp = hata_evol[jtime].compute_tilde_B_tt(true) ;
00677     if (tmp.get_dzpuis() == 3)
00678         tmp.dec_dzpuis() ; // dzpuis 3->2
00679 
00680     B_hata_evol.update( tmp, jtime, the_time[jtime] ) ;
00681     }
00682     return A_hata_evol[jtime] ;
00683 }
00684 
00685 
00686 const Scalar& Time_slice_conf::psi() const {
00687 
00688     if (!( psi_evol.is_known(jtime) ) ) {
00689 
00690         assert( n_evol.is_known(jtime) ) ; 
00691         assert( npsi_evol.is_known(jtime) ) ; 
00692         
00693         psi_evol.update( npsi_evol[jtime] / n_evol[jtime] , jtime, the_time[jtime] ) ; 
00694     }
00695 
00696     return psi_evol[jtime] ;
00697 
00698 } 
00699 
00700 const Scalar& Time_slice_conf::psi4() const {
00701 
00702     if (p_psi4 == 0x0)  {
00703 
00704         p_psi4 = new Scalar( pow( psi(), 4.) ) ; 
00705         p_psi4->std_spectral_base() ;
00706     }
00707 
00708     return *p_psi4 ;
00709 
00710 } 
00711 
00712 const Scalar& Time_slice_conf::ln_psi() const {
00713 
00714     if (p_ln_psi == 0x0)  {
00715 
00716         p_ln_psi = new Scalar( log( psi() ) ) ; 
00717         p_ln_psi->std_spectral_base() ;
00718     }
00719 
00720     return *p_ln_psi ;
00721 
00722 } 
00723 
00724 
00725 const Scalar& Time_slice_conf::npsi() const {
00726 
00727     if (!( npsi_evol.is_known(jtime) ) ) {
00728         
00729         assert( n_evol.is_known(jtime) ) ; 
00730         assert( psi_evol.is_known(jtime) ) ; 
00731 
00732         npsi_evol.update( psi_evol[jtime] * n_evol[jtime], jtime, the_time[jtime] ) ; 
00733     }
00734 
00735     return npsi_evol[jtime] ;
00736 
00737 }
00738 
00739 
00740 const Metric& Time_slice_conf::tgam() const {
00741 
00742     if (p_tgamma == 0x0) {
00743         p_tgamma = new Metric( ff.con() + hh() ) ; 
00744     }
00745     
00746     return *p_tgamma ; 
00747 
00748 }
00749 
00750 
00751 const Sym_tensor& Time_slice_conf::hh(Param*, Param*) const {
00752 
00753     assert( hh_evol.is_known(jtime) ) ; 
00754     return hh_evol[jtime] ; 
00755 
00756 }
00757 
00758 Sym_tensor Time_slice_conf::aa() const {
00759 
00760     return hata()/(psi4()*psi()*psi()) ; 
00761 
00762 }
00763 
00764 
00765 const Sym_tensor& Time_slice_conf::hata() const {
00766 
00767     if( !(hata_evol.is_known(jtime)) ) {
00768 
00769         assert( hh_evol.is_known(jtime) ) ; 
00770 
00771         Sym_tensor hh_point = hh_evol.time_derive(jtime, scheme_order) ; 
00772         hh_point.inc_dzpuis(2) ; // dzpuis : 0 -> 2
00773     
00774         Sym_tensor resu = hh_point - hh().derive_lie(beta()) 
00775             - 0.6666666666666666 * beta().divergence(ff) * hh()
00776             + beta().ope_killing_conf(ff) ; 
00777 
00778         resu = psi4()*psi()*psi() * resu / (2*nn()) ;
00779         
00780         hata_evol.update(resu, jtime, the_time[jtime]) ;
00781         
00782     }
00783      
00784     return hata_evol[jtime] ; 
00785 
00786 }
00787 
00788 
00789 const Scalar& Time_slice_conf::trk() const {
00790 
00791     if( !(trk_evol.is_known(jtime)) ) {
00792 
00793         psi() ; 
00794         Scalar resu = -6*psi_evol.time_derive(jtime, scheme_order) / psi() ;
00795     resu.inc_dzpuis(2) ;
00796     resu += beta().divergence(ff) 
00797         + 6.*contract( beta(), 0, ln_psi().derive_cov(ff), 0 ) ;
00798         resu = resu / nn() ; 
00799         
00800         trk_evol.update(resu, jtime, the_time[jtime]) ;
00801     } 
00802     
00803     return trk_evol[jtime] ; 
00804 
00805 }
00806 
00807 
00808 const Vector& Time_slice_conf::hdirac() const {
00809 
00810     if (p_hdirac == 0x0) {
00811         p_hdirac = new Vector( hh().divergence(ff) ) ; 
00812     }
00813     
00814     return *p_hdirac ; 
00815 
00816 }
00817 
00818 const Vector& Time_slice_conf::vec_X(int methode_poisson) const {
00819 
00820     if (p_vec_X == 0x0) {
00821     assert ( hata_evol.is_known(jtime) ) ;
00822     Vector source = hata_evol[jtime].divergence(ff) ;
00823     source.inc_dzpuis() ; // dzpuis 3-> 4
00824         p_vec_X = new Vector( source.poisson( 1./3., methode_poisson ) ) ;
00825     }    
00826     return *p_vec_X ;     
00827 }
00828 
00829 void Time_slice_conf::compute_X_from_momentum_constraint
00830 (const Vector& hat_S, const Sym_tensor_tt& hata_tt, int iter_max, double precis,
00831  double relax, int meth_poisson) {
00832 
00833     // Some checks
00834     assert( hat_S.get_index_type(0) == CON ) ;
00835 #ifndef NDEBUG
00836     for (int i=1; i<=3; i++)
00837     for (int j=i; j<=3; j++)
00838         assert( hata_tt(i,j).get_dzpuis() == 2 ) ;
00839 #endif
00840     assert( hh_evol.is_known(jtime) ) ;
00841 
00842     // Initializations
00843     //----------------
00844     const Tensor_sym& delta = tgam().connect().get_delta() ;
00845     if (p_vec_X != 0x0) {
00846     delete p_vec_X ;
00847     p_vec_X = 0x0 ;
00848     }
00849 
00850     // Constant part of the source
00851     Vector source = hat_S - contract( delta, 1, 2, hata_tt, 0, 1 ) 
00852     - 2./3.*psi4()*psi()*psi()*contract( tgam().con(), 0, trk().derive_cov(ff), 0 );
00853 
00854     p_vec_X = new Vector( source.poisson( 1./3., meth_poisson ) ) ;
00855 
00856     // Iteration on the vector X
00857     //--------------------------
00858     for (int it=0; it<iter_max; it++) {
00859 
00860     Vector fuente = source 
00861         - contract( delta, 1, 2, p_vec_X->ope_killing_conf(ff), 0, 1 ) ;
00862 
00863     Vector X_new = fuente.poisson( 1./3., meth_poisson) ;
00864 
00865     // Control of the convergence
00866     double diff = 0. ;
00867     for (int i=1; i<=3; i++)
00868         diff += max( max( abs( X_new(i) - (*p_vec_X)(i) ) ) ) ;
00869 
00870     // Relaxation
00871     (*p_vec_X) = relax*X_new + ( 1. - relax )*(*p_vec_X) ;
00872 
00873     // If converged, gets out of the loop
00874     if (diff < precis) break ;
00875     }
00876 
00877     // Update of \hat{A}^{ij} and reset of related quantities
00878     hata_evol.update( hata_tt + p_vec_X->ope_killing_conf(ff), jtime, the_time[jtime] ) ;
00879 
00880     k_dd_evol.downdate(jtime) ; 
00881     k_uu_evol.downdate(jtime) ;     
00882 }
00883 
00884 void Time_slice_conf::set_AB_hata(const Scalar& A_in, const Scalar& B_in) {
00885 
00886     A_hata_evol.update(A_in, jtime, the_time[jtime]) ;
00887     B_hata_evol.update(B_in, jtime, the_time[jtime]) ;
00888 
00889     hata_evol.downdate(jtime) ;
00890     // Reset of quantities depending on A^{ij}:
00891     k_dd_evol.downdate(jtime) ; 
00892     k_uu_evol.downdate(jtime) ;     
00893 
00894 }
00895 
00896 
00897 void Time_slice_conf::check_psi_dot(Tbl& tlnpsi_dot, Tbl& tdiff, Tbl& tdiff_rel) const {
00898 
00899     // Computation of d/dt ln Psi :
00900      
00901     Scalar lnpsi_dot(psi().get_mp()) ; 
00902     if ( psi_evol.is_known(jtime-1) ) {       
00903         lnpsi_dot = psi_evol.time_derive(jtime, scheme_order) / psi() ; 
00904     }
00905     else {
00906         lnpsi_dot = npsi_evol.time_derive(jtime, scheme_order) / npsi()
00907                            - n_evol.time_derive(jtime, scheme_order) / nn()  ; 
00908     }
00909     
00910     tlnpsi_dot = max(abs(lnpsi_dot)) ; 
00911         
00912     // Error on the d/dt ln Psi relation : 
00913     
00914     Scalar diff = contract(beta(),0, ln_psi().derive_cov(ff),0) 
00915         + 0.1666666666666666 * ( beta().divergence(ff) - nn() * trk() ) ;
00916 
00917     diff.dec_dzpuis(2) ; 
00918     
00919     Tbl tref = max(abs(diff)) + tlnpsi_dot ; 
00920     
00921     diff -= lnpsi_dot ; 
00922     
00923     tdiff = max(abs(diff)) ;
00924             
00925     tdiff_rel = tdiff / tref  ; 
00926     
00927 }
00928 
00929 
00930                 //------------------//
00931                 //      output      //
00932                 //------------------//
00933 
00934 ostream& Time_slice_conf::operator>>(ostream& flux) const {
00935 
00936     Time_slice::operator>>(flux) ; 
00937 
00938     flux << "Triad on which the components of the flat metric are defined:\n" 
00939         << *(ff.get_triad()) << '\n' ;  
00940 
00941     if (psi_evol.is_known(jtime)) {
00942         maxabs( psi_evol[jtime], "Psi", flux) ;
00943     }
00944     if (npsi_evol.is_known(jtime)) {
00945         maxabs( npsi_evol[jtime], "N Psi", flux) ;
00946     }
00947     if (hh_evol.is_known(jtime)) {
00948         maxabs( hh_evol[jtime], "h^{ij}", flux) ;
00949     }
00950     if (hata_evol.is_known(jtime)) {
00951         maxabs( hata_evol[jtime], "hat{A}^{ij}", flux) ;
00952     }
00953 
00954     if (p_tgamma != 0x0) flux << 
00955         "Conformal metric tilde gamma is up to date" << endl ; 
00956     if (p_psi4 != 0x0) maxabs( *p_psi4, "Psi^4", flux) ; 
00957     if (p_ln_psi != 0x0) maxabs( *p_ln_psi, "ln(Psi)", flux) ; 
00958     if (p_hdirac != 0x0) maxabs( *p_hdirac, "H^i", flux) ; 
00959     if (p_vec_X != 0x0) maxabs( *p_vec_X, "X^i", flux) ; 
00960     
00961     return flux ; 
00962 
00963 }
00964 
00965 
00966 
00967 void Time_slice_conf::sauve(FILE* fich, bool partial_save) const {
00968     
00969     // Writing of quantities common to all derived classes of Time_slice
00970     // -----------------------------------------------------------------
00971     
00972     Time_slice::sauve(fich, true) ; 
00973     
00974     // Writing of quantities common to all derived classes of Time_slice_conf
00975     // ----------------------------------------------------------------------
00976     
00977     int jmin = jtime - depth + 1 ; 
00978 
00979     // Psi
00980     psi() ;     // forces the update at the current time step
00981     for (int j=jmin; j<=jtime; j++) {
00982         int indicator = (psi_evol.is_known(j)) ? 1 : 0 ; 
00983         fwrite_be(&indicator, sizeof(int), 1, fich) ;
00984         if (indicator == 1) psi_evol[j].sauve(fich) ;
00985     }
00986 
00987     // hat{A}
00988     hata() ;
00989     for (int j=jmin; j<=jtime; j++) {
00990     int indicator = ( hata_evol.is_known(j) ? 1 : 0 ) ;
00991     fwrite_be(&indicator, sizeof(int), 1, fich) ;
00992     if (indicator == 1) hata_evol[j].sauve(fich) ;
00993     }
00994 
00995     //A and B...
00996     A_hata() ;
00997     for (int j=jmin; j<=jtime; j++) {
00998         int indicator = (A_hata_evol.is_known(j)) ? 1 : 0 ; 
00999         fwrite_be(&indicator, sizeof(int), 1, fich) ;
01000         if (indicator == 1) A_hata_evol[j].sauve(fich) ;
01001     }
01002 
01003     B_hata() ;
01004     for (int j=jmin; j<=jtime; j++) {
01005         int indicator = (B_hata_evol.is_known(j)) ? 1 : 0 ; 
01006         fwrite_be(&indicator, sizeof(int), 1, fich) ;
01007         if (indicator == 1) B_hata_evol[j].sauve(fich) ;
01008     }   
01009 
01010     // Case of a complete save
01011     // -----------------------
01012     if (!partial_save) {
01013         cout << "Time_slice_conf::sauve: the full writing is not ready yet !" 
01014              << endl ; 
01015         abort() ; 
01016     }
01017     
01018 }

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