tslice_check_einstein.C

00001 /*
00002  *  Methods of class Time_slice to check Einstein equation solutions.
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 tslice_check_einstein_C[] = "$Header: /cvsroot/Lorene/C++/Source/Time_slice/tslice_check_einstein.C,v 1.8 2010/10/20 07:58:09 j_novak Exp $" ;
00029 
00030 /*
00031  * $Id: tslice_check_einstein.C,v 1.8 2010/10/20 07:58:09 j_novak Exp $
00032  * $Log: tslice_check_einstein.C,v $
00033  * Revision 1.8  2010/10/20 07:58:09  j_novak
00034  * Better implementation of the explicit time-integration. Not fully-tested yet.
00035  *
00036  * Revision 1.7  2007/11/06 12:10:56  j_novak
00037  * Corrected some mistakes.
00038  *
00039  * Revision 1.6  2007/11/06 11:53:34  j_novak
00040  * Use of contravariant version of the 3+1 equations.
00041  *
00042  * Revision 1.5  2007/06/28 14:40:36  j_novak
00043  * Dynamical check: the fields in the last domain are set to zero to avoid dzpuis problems
00044  *
00045  * Revision 1.4  2004/05/12 15:24:20  e_gourgoulhon
00046  * Reorganized the #include 's, taking into account that
00047  * time_slice.h contains now an #include "metric.h".
00048  *
00049  * Revision 1.3  2004/04/07 07:58:21  e_gourgoulhon
00050  * Constructor as Minkowski slice: added call to std_spectral_base()
00051  * after setting the lapse to 1.
00052  *
00053  * Revision 1.2  2004/04/05 12:38:45  j_novak
00054  * Minor modifs to prevent some warnings.
00055  *
00056  * Revision 1.1  2004/04/05 11:54:20  j_novak
00057  * First operational (but not tested!) version of checks of Eintein equations.
00058  *
00059  *
00060  * $Header: /cvsroot/Lorene/C++/Source/Time_slice/tslice_check_einstein.C,v 1.8 2010/10/20 07:58:09 j_novak Exp $
00061  *
00062  */
00063 
00064 // C headers
00065 #include <stdlib.h>
00066 #include <assert.h>
00067 
00068 // Lorene headers
00069 #include "time_slice.h"
00070 #include "unites.h"
00071 
00072 Tbl Time_slice::check_hamiltonian_constraint(const Scalar* energy_density,
00073                          ostream& ost, bool verb) const {
00074   using namespace Unites ;
00075 
00076   bool vacuum = ( energy_density == 0x0 ) ;
00077   
00078   Scalar field = trk() * trk() - contract( k_uu(), 0, 1, k_dd(), 0, 1 ) ;
00079   
00080   field.dec_dzpuis() ;  // dzpuis: 4 -> 3 
00081     
00082   field += gam().ricci_scal() ;
00083 
00084   const Scalar* matter ;
00085   if (vacuum) 
00086     matter = new Scalar (0*field) ;
00087   else
00088     matter = energy_density ;
00089 
00090   if (verb)  ost << endl ;
00091   const char* comment = 0x0 ;
00092   if (verb) comment = "Check of the Hamiltonian constraint" ;
00093   Tbl resu = maxabs(field - (4*qpig) * (*matter), comment, ost,verb ) ;
00094   if (verb) ost << endl ;
00095 
00096   if (vacuum) delete matter ;
00097 
00098   return resu ;
00099 
00100 }
00101     
00102 Tbl Time_slice::check_momentum_constraint(const Vector* momentum_density, 
00103                       ostream& ost, bool verb) const {
00104   using namespace Unites ;
00105   
00106   bool vacuum = ( momentum_density == 0x0 ) ;
00107   
00108   Vector field = k_uu().divergence(gam()) - trk().derive_con(gam()) ;
00109 
00110 
00111   const Vector* matter ;
00112   if (vacuum) 
00113     matter = new Vector (0*field) ;
00114   else {
00115     assert (momentum_density->get_index_type(0) == CON) ;
00116     matter = momentum_density ;
00117   }
00118   
00119   if (verb)  ost << endl ;
00120   const char* comment = 0x0 ;
00121   if (verb) comment = "Check of the momentum constraint" ;
00122   Tbl resu = maxabs(field - (2*qpig) * (*matter), comment, ost, verb ) ;
00123 
00124   if (verb)  ost << endl ;
00125 
00126   if (vacuum) delete matter ;
00127 
00128   return resu ;
00129 
00130 }
00131 
00132 Tbl Time_slice::check_dynamical_equations(const Sym_tensor* strain_tensor,
00133                       const Scalar* energy_density,
00134                       ostream& ost, bool verb) const {
00135 
00136   using namespace Unites ;
00137 
00138   bool vacuum = ( ( strain_tensor == 0x0 ) && ( energy_density == 0x0 ) ) ;
00139 
00140   Sym_tensor dyn_lhs = (k_dd_evol.time_derive(jtime, scheme_order)).up_down(gam()) ;
00141   int nz = dyn_lhs.get_mp().get_mg()->get_nzone() ;
00142   dyn_lhs.annule_domain(nz-1) ;
00143   dyn_lhs = dyn_lhs - k_dd().derive_lie(beta()).up_down(gam())  ;
00144   dyn_lhs.annule_domain(nz-1) ;
00145   
00146   const Sym_tensor* matter ;
00147   if (vacuum) 
00148     matter = new Sym_tensor(0 * dyn_lhs) ;
00149   else {
00150     const Scalar* ener = 0x0 ;
00151     const Sym_tensor* sij = 0x0 ;
00152     bool new_e = false ;
00153     bool new_s = false ;
00154     if (energy_density == 0x0) {
00155       ener = new Scalar ( 0 * nn() ) ;
00156       new_e = true ;
00157     }
00158     else {
00159       ener = energy_density ;
00160       if (strain_tensor == 0x0) {
00161     sij = new Sym_tensor(0 * dyn_lhs) ;
00162     new_s = true ;
00163       }
00164       else
00165     sij = strain_tensor ;
00166     }
00167     matter = new Sym_tensor( (sij->trace(gam()) - *ener)*gam().con() 
00168                   - 2*(*sij) ) ;
00169     if (new_e) delete ener ;
00170     if (new_s) delete sij ;
00171   }
00172 
00173   Sym_tensor dyn_rhs = nn()*( (gam().ricci().up_down(gam()) + trk()*k_uu() 
00174                              + qpig * (*matter))) ;
00175   dyn_rhs.annule_domain(nz-1) ;
00176   dyn_rhs = dyn_rhs - 2*nn()*contract(k_uu(), 1, k_dd().up(0, gam()), 1)  ;
00177   dyn_rhs.annule_domain(nz-1) ;
00178   dyn_rhs = dyn_rhs - nn().derive_con(gam()).derive_con(gam()) ;
00179   dyn_rhs.annule_domain(nz-1) ;
00180     
00181   if (verb) ost << endl ;
00182   const char* comment = 0x0 ;
00183   if (verb) comment = "Check of the dynamical equations" ;
00184   Tbl resu_tmp = maxabs(dyn_lhs - dyn_rhs, comment, ost, verb) ;
00185 
00186   if (verb) ost << endl ;
00187 
00188   if (vacuum) delete matter ;
00189 
00190   Tbl resu(4, 4, resu_tmp.get_dim(0)) ;
00191   resu.annule_hard() ; //## To initialize resu(0,0,...) which 
00192                        // does not correspond to a valid component
00193 
00194   for (int i=1; i<=3; i++) {
00195     for (int j=1; j<=i; j++) {
00196       Itbl idx(2) ;
00197       idx.set(0) = i ;
00198       idx.set(1) = j ;
00199       int pos = dyn_lhs.position(idx) ;
00200       assert (dyn_rhs.position(idx) == pos) ;
00201 
00202       for (int lz=0; lz<resu_tmp.get_dim(0); lz++)
00203     resu.set(i,j,lz) = resu_tmp(pos, lz) ;
00204     }
00205   }
00206 
00207   return resu ;
00208 
00209 }

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