isol_hole_compute_metric.C

00001 /*
00002  * Method of class Isol_Hole to compute metric data associated to a quasistationary single black hole spacetime slice.
00003  *
00004  * (see file isol_hole.h for documentation).
00005  *
00006  */
00007 
00008 /*
00009  *   Copyright (c) 2009 Nicolas Vasset
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 as published by
00015  *   the Free Software Foundation; either version 2 of the License, or
00016  *   (at your option) any later version.
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 // Headers C
00030 #include "math.h"
00031 
00032 // Headers Lorene
00033 #include "tenseur.h"
00034 #include "star.h"
00035 #include "param.h"
00036 #include "graphique.h"
00037 #include "unites.h"     
00038 #include "spheroid.h"
00039 #include "excision_surf.h"
00040 #include "excision_hor.h"
00041 #include "param_elliptic.h"
00042 #include "vector.h"
00043 #include "scalar.h"
00044 #include "diff.h"
00045 #include "proto.h"
00046 #include "isol_hole.h"
00047 
00048 
00049 void Isol_hole::compute_stat_metric(double precis, double relax, int mer_max, int mer_max2, bool isvoid){
00050     
00051     // Fundamental constants and units
00052     // -------------------------------
00053     using namespace Unites ;
00054 
00055     // Noticing the user that the iteration has started
00056 
00057     cout << "================================================================" << endl;
00058     cout << "STARTING THE MAIN ITERATION FOR COMPUTING METRIC FIELDS" << endl;
00059     cout << "        iteration parameters are the following:        " << endl;
00060     cout << "        convergence precision required:" << precis << endl;
00061     cout << "        max number of global steps    :" << mer_max << endl;
00062     cout << "        relaxation parameter          :" << relax   << endl;
00063     cout << "================================================================" << endl;
00064 
00065 
00066   // Construction of a multi-grid (Mg3d) and an affine mapping from the class mapping
00067   // --------------------------------------------------------------------------------
00068   const Map_af* map = dynamic_cast<const Map_af*>(&mp) ;
00069   const Mg3d* mgrid = (*map).get_mg();
00070     
00071   // Construct angular grid for h(theta,phi) 
00072   const Mg3d* g_angu = (*mgrid).get_angu_1dom() ;
00073   
00074   const int nz = (*mgrid).get_nzone();  // Number of domains
00075   int nt = (*mgrid).get_nt(1);  // Number of collocation points in theta in each domain
00076   const int np = (*mgrid).get_np(1); 
00077   const Coord& rr = (*map).r;
00078    Scalar rrr (*map) ; 
00079   rrr = rr ; 
00080   rrr.std_spectral_base();  
00081   assert((rrr.val_grid_point(1,0,0,0) - 1.) <= 1.e-9); // For now the code handles only horizons at r=1, corresponding to the first shell inner boundary. This test assures this is the case with our mapping.
00082  
00083   // Angular mapping defined as well
00084   //--------------------------------
00085   double r_limits2[] = {rrr.val_grid_point(1,0,0,0), rrr.val_grid_point(2,0,0,0)} ; 
00086   const Map_af map_2(*g_angu, r_limits2); //2D mapping; check if this is useful.
00087   const Metric_flat& mets = (*map).flat_met_spher() ;
00088  
00089     //----------------
00090     // Initializations
00091     // ---------------
00092 
00093 
00094 
00095   Scalar logn (*map) ; logn = log(lapse) ; 
00096   logn.std_spectral_base();
00097 
00098 
00099   Scalar logpsi(*map) ; logpsi = log(conf_fact) ;
00100   logpsi.std_spectral_base();
00101   Scalar psi4 (*map) ; psi4 = conf_fact*conf_fact*conf_fact*conf_fact ;
00102   Scalar npsi (*map) ;  npsi =conf_fact*lapse ;
00103   
00104 
00105   Scalar conf_fact_new(*map) ; conf_fact_new.annule_hard(); conf_fact_new.std_spectral_base(); 
00106   Scalar npsi_new(*map); npsi_new.annule_hard(); npsi_new.std_spectral_base();
00107 
00108   Vector shift_new (*map, CON, (*map).get_bvect_spher()); 
00109   for(int i=1; i<=3; i++){
00110     shift_new.set(i)=0;
00111   }
00112   shift_new.std_spectral_base();
00113  
00114   // Non conformally flat variables
00115 
00116  
00117   Sym_tensor gamtuu = mets.con() + hij; 
00118   Metric gamt(gamtuu);
00119   Metric gam(gamt.cov()*psi4) ;
00120   Sym_tensor gamma = gam.cov();
00121   
00122 
00123   // Extrinsic curvature variables
00124 
00125   Sym_tensor aa(*map, CON, (*map).get_bvect_spher());
00126   for (int iii= 1; iii<=3; iii++){ 
00127     for(int j=1; j<=3; j++){
00128       aa.set(iii,j)= 0;
00129     }
00130   }
00131   aa.std_spectral_base(); 
00132   Scalar aa_quad_scal(*map) ; aa_quad_scal = 0. ;
00133 
00134   Sym_tensor aa_hat(*map, CON, (*map).get_bvect_spher());
00135   for (int iii= 1; iii<=3; iii++){ 
00136     for(int j=1; j<=3; j++){
00137       aa_hat.set(iii,j)= 0;
00138     }
00139   }
00140  
00141   Sym_tensor kuu = aa/psi4 ;
00142   Sym_tensor kuu2 = aa_hat/(psi4*psi4*sqrt(psi4));
00143   Sym_tensor kdd = contract (gamma, 0, contract(gamma, 1, kuu, 0),1);
00144  
00145 
00146   // (2,1)-rank delta tensor: difference between ricci rotation coefficients. 
00147 
00148    Tensor delta = -0.5*contract( hij, 1, gamt.cov().derive_cov(mets), 2);
00149    Scalar tmp(*map);
00150 
00151    for (int i=1; i<=3; i++) {
00152      for (int j=1; j<=3; j++) {
00153        for (int k=1; k<=3; k++) {
00154      tmp = 0.;
00155      tmp =  -0.5 *(gamt.cov().derive_con(mets))(i,j,k);
00156      for (int l=1; l<=3; l++) {
00157        tmp += -0.5*( gamt.cov()(i,l)*(hij.derive_cov(mets))(k,l,j) + gamt.cov()(l,j)*(hij.derive_cov(mets))(k,l,i));
00158      }
00159          delta.set(k,i,j) += tmp ; 
00160        }
00161      }
00162    }
00163    
00164 
00165    // Conformal Rstar scalar(eq 61, Bonazzola et al. 2003) 
00166    Scalar Rstar = 
00167     0.25 * contract(gamt.con(), 0, 1,
00168             contract(hij.derive_cov(mets), 0, 1, gamt.cov().derive_cov(mets), 0, 1), 0, 1 ) 
00169     - 0.5  * contract(gamt.con(), 0, 1,
00170               contract(hij.derive_cov(mets), 0, 1, gamt.cov().derive_cov(mets), 0, 2), 0, 1 ) ; 
00171   
00172    Scalar norm(*map);
00173    Scalar norm3(*map);
00174 
00175    
00176    if (isvoid == false){
00177    
00178      cout <<"FAIL: case of non-void spacetime not treated yet" << endl;
00179    }
00180  
00181    else
00182      {
00183 
00184   // Parameters for the iteration 
00185 
00186 
00187   double diff_ent = 1 ; // initialization of difference marker between two iterations; 
00188 
00189   int util = 0; // Tool used to stop tensorial iteration at any wished step "util"
00190   
00194   
00195   
00196   for(int mer=0 ;(diff_ent > precis) && (mer<mer_max) ; mer++) {    
00197 
00198    //Global relaxation coefficient
00199 
00200     // Scalar variables linked to the norm of normal vector to horizon.
00201     norm = sqrt(1. + hij(1,1)); norm.std_spectral_base();
00202     norm3 = sqrt(1. + hij(3,3)); norm3.std_spectral_base();
00203     
00204 
00205 
00209 
00210 
00211 
00213       // Solving for (Psi-1) //
00215 
00216   
00217       // Setting of the boundary
00218 
00219       double   diric= 0.; 
00220       double   neum = 1.; 
00221 
00222       Vector ssalt = rrr.derive_cov(gam);
00223       Vector ssaltcon = ssalt.up_down(gam);
00224       Scalar ssnormalt = sqrt(contract (ssalt,0, ssaltcon, 0));
00225       ssnormalt.std_spectral_base();
00226       
00227       ssalt.annule_domain(nz-1);
00228       ssalt.annule_domain(0);
00229       ssaltcon.annule_domain(nz-1);
00230       ssaltcon.annule_domain(0);
00231       
00232       ssalt = ssalt/ssnormalt;
00233       ssaltcon = ssaltcon/ssnormalt;
00234       Vector ssconalt = ssaltcon*conf_fact*conf_fact;    // \tilde{s} in the notations of Gourgoulhon and Jaramillo, 2006
00235       ssconalt.std_spectral_base();
00236       ssconalt.annule_domain(nz-1);
00237       Scalar bound3bis =   -((1./conf_fact)*contract((contract(kdd,1,ssconalt,0)),0, ssconalt,0));
00238   
00239       bound3bis.annule_domain(nz-1);
00240       bound3bis += -conf_fact*ssconalt.divergence(gamt);
00241       bound3bis.annule_domain(nz-1);     
00242       bound3bis = 0.25*bound3bis;
00243       bound3bis += -contract(conf_fact.derive_cov(gamt),0,ssconalt,0) + conf_fact.dsdr();
00244       bound3bis.annule_domain(nz-1);
00245       bound3bis.std_spectral_base();
00246       bound3bis.set_spectral_va().ylm();    
00247 
00248       Mtbl_cf *boundd3bis = bound3bis.set_spectral_va().c_cf;
00249             
00250 
00252 
00253 
00254 
00255       // Computing the source
00256    
00257       Scalar source_conf_fact(*map) ; source_conf_fact=3. ; // Pour le fun... 
00258       source_conf_fact.std_spectral_base();                            
00259 
00260       Scalar d2logpsi = contract(conf_fact.derive_cov(mets).derive_cov(mets), 0, 1, hij, 0,1);  
00261       d2logpsi.inc_dzpuis(1);
00262  
00263       source_conf_fact = -(0.125* aa_quad_scal )/(psi4*conf_fact*conf_fact*conf_fact) +  conf_fact* 0.125* Rstar - d2logpsi; 
00264       
00265       source_conf_fact.std_spectral_base(); 
00266 
00267       if (source_conf_fact.get_etat() == ETATZERO) {
00268     source_conf_fact.annule_hard() ;
00269     source_conf_fact.set_dzpuis(4) ;
00270     source_conf_fact.std_spectral_base() ;
00271       }
00272       source_conf_fact.set_spectral_va().ylm();
00273     
00274       
00275 
00277 
00278       // System inversion   
00279 
00280       Param_elliptic source11(source_conf_fact);
00281       conf_fact_new = source_conf_fact.sol_elliptic_boundary(source11, *boundd3bis, diric , neum) + 1 ; // Resolution has been done for quantity Q-1, because our solver gives a vanishing solution at infinity! 
00282       
00283 
00285  
00286     // tests for resolution
00287 
00288     Scalar baba2 = (conf_fact_new-1).laplacian();
00289 //     cout << "psi+1-résolution" << endl;
00290 //     maxabs (baba2 - source_conf_fact);
00291    
00292     Scalar psinewbis = conf_fact_new -1. ; psinewbis.annule_domain(nz -1);
00293     psinewbis.std_spectral_base();
00294     psinewbis = psinewbis.dsdr();
00295     Scalar psinewfin2 (map_2) ;
00296     psinewfin2.allocate_all(); 
00297     psinewfin2.set_etat_qcq();  
00298     psinewfin2.std_spectral_base();
00299     
00300     for (int k=0; k<np; k++)
00301       for (int j=0; j<nt; j++) {
00302     
00303     psinewfin2.set_grid_point(0, k, j, 0) = psinewbis.val_grid_point(1, k,j,0) - bound3bis.val_grid_point(1, k, j, 0);
00304     
00305       }
00306     
00307     //  maxabs (psinewfin2);
00308     
00309 
00311 
00312     // Update during the loop
00313 
00314     conf_fact = conf_fact_new* (1-relax) + conf_fact* relax ;
00315     psi4 = conf_fact*conf_fact*conf_fact*conf_fact;
00316     logpsi = log(conf_fact) ; 
00317     logpsi.std_spectral_base();    
00318 
00319 
00320 
00321 
00325 
00326 
00327 
00329     // Solving for (N*Psi -1)/ 
00331 
00332 
00333     // Setting of the boundary
00334     assert (NorKappa == false) ; 
00335     Scalar bound(*map);
00336     bound = (boundNoK)*conf_fact -1;
00337     bound.annule_domain(nz -1);
00338     bound.std_spectral_base();
00339     bound.set_spectral_va().ylm();
00340     Mtbl_cf *boundd = bound.get_spectral_va().c_cf;
00341 
00342      diric =1; 
00343      neum = 0 ; 
00344 
00345 
00346     
00348 
00349     // Computing the source ...      
00350          Scalar d2lognpsi = contract(npsi.derive_cov(mets).derive_cov(mets), 0, 1, hij, 0,1);
00351        d2lognpsi.inc_dzpuis(1); //  dzpuis correction.
00352   
00353        Scalar source_npsi = npsi*(aa_quad_scal*(7./8.)/(psi4*psi4) + Rstar/8.) - d2lognpsi; // ?
00354     source_npsi.std_spectral_base();
00355     if (source_npsi.get_etat() == ETATZERO) {
00356       source_npsi.annule_hard() ;
00357       source_npsi.set_dzpuis(4) ;
00358       source_npsi.std_spectral_base() ;
00359     }
00360 
00361 
00363 
00364     // Inversion of the operator
00365     Param_elliptic source1 (source_npsi); 
00366     npsi_new = source_npsi.sol_elliptic_boundary(source1, *boundd, diric, neum) ;
00367 
00368     npsi_new = npsi_new +1;
00369 
00370   
00372 
00373     // Resolution tests in npsi
00374     Scalar baba = npsi_new.laplacian();
00375 //       cout << "resolution_npsi" << endl;
00376 //      maxabs (baba - source_npsi);
00377 
00378 
00379     //  cout << "bound_npsi" << endl;
00380         Scalar npsibound2 (map_2) ;
00381     npsibound2.allocate_all(); 
00382     npsibound2.set_etat_qcq();  
00383     npsibound2.std_spectral_base();      
00384     for (int k=0; k<np; k++)
00385       for (int j=0; j<nt; j++) {
00386              
00387     npsibound2.set_grid_point(0, k, j, 0) = npsi_new.val_grid_point(1, k,j,0) - bound.val_grid_point(1, k, j, 0) -1.;
00388              
00389       }
00390 
00391     //   maxabs (npsibound2);
00392     
00393 
00395 
00396     // Update during the loop
00397 
00398       
00399       npsi = npsi_new*(1-relax) + npsi* relax; 
00400        lapse = npsi/conf_fact; 
00401        logn = log(lapse);
00402        logn.std_spectral_base(); 
00403 
00404  
00408 
00410       //Resolution in Beta //
00412 
00413 
00414       // Setting of the boundary  
00415            
00416       bound = (boundNoK)/(conf_fact*conf_fact) ;
00417       bound.annule_domain(nz -1);
00418  
00419 
00420       Scalar hor_rot(*map); hor_rot.annule_hard(); hor_rot = Omega; // Rotation parameter for spacetime
00421       hor_rot.std_spectral_base() ; hor_rot.mult_rsint();
00422       hor_rot.annule_domain(nz -1);
00423 
00424       Vector limit = shift_new;
00425       Vector ephi(*map, CON, (*map).get_bvect_spher());
00426       ephi.set(1).annule_hard();
00427       ephi.set(2).annule_hard();
00428       ephi.set(3) = 1;
00429       ephi.std_spectral_base();
00430       ephi.annule_domain(nz -1);
00431 
00432       limit = bound*ssconalt + hor_rot*ephi;
00433       limit.std_spectral_base(); // Boundary is fixed by value of 3 components of a vector (rather than value of potentials)   
00434 
00435       Scalar Vrb = limit(1);  Vrb.set_spectral_va().ylm();
00436       Scalar mmub = limit.mu(); mmub.set_spectral_va().ylm(); 
00437       Scalar etab = limit.eta(); etab.set_spectral_va().ylm();
00438 
00439 
00440 
00442  
00443      // Computing the source
00444        Vector deltaA =  - 2*lapse*contract(delta, 1,2, aa, 0,1);
00445        Vector hijddb =  - contract (shift.derive_cov(mets).derive_cov(mets), 1,2, hij, 0,1) ;
00446        Vector hijddivb =  - 0.3333333333333* contract (shift.divergence(mets).derive_cov(mets),0, hij,1);
00447        hijddb.inc_dzpuis(); // dzpuis fixing patch... 
00448        hijddivb.inc_dzpuis(); 
00449        
00450        Vector sourcevect2(*map,CON, (*map).get_bvect_spher()); 
00451       sourcevect2= (2.* contract(aa, 1, lapse.derive_cov(mets),0)-12*lapse*contract(aa, 1, logpsi.derive_cov(mets), 0))  + deltaA + hijddb + hijddivb ; 
00452        
00453       sourcevect2.set(1).set_dzpuis(4);
00454       sourcevect2.set(2).set_dzpuis(4);
00455       sourcevect2.set(3).set_dzpuis(4);
00456       sourcevect2.std_spectral_base(); 
00457       if(sourcevect2.eta().get_etat() == ETATZERO)
00458     { sourcevect2.set(2).annule_hard();}
00459  
00460       
00461       double lam = (1./3.);    
00462     
00463    
00464 
00466 
00467       // System inversion
00468       sourcevect2.poisson_boundary2(lam, shift_new, Vrb, etab, mmub, 1., 0., 1. ,0. ,1. ,0.) ;   
00469 
00470 
00471 
00473     
00474     // resolution tests
00475     Vector source2 = contract(shift_new.derive_con(mets).derive_cov(mets), 1,2) + lam* contract(shift_new.derive_cov(mets), 0,1).derive_con(mets);
00476     source2.inc_dzpuis(1);
00477     //  maxabs (source2 - sourcevect2);   
00478 
00479     Scalar mufin = shift_new.mu();
00480     mufin.set_spectral_va().coef();
00481          
00482     Scalar mufin2 (map_2) ;
00483     mufin2.allocate_all(); 
00484     mufin2.set_etat_qcq();  
00485     mufin2.std_spectral_base();
00486          
00487     for (int k=0; k<np; k++)
00488       for (int j=0; j<nt; j++) {
00489              
00490     mufin2.set_grid_point(0, k, j, 0) = mufin.val_grid_point(1, k,j,0) - mmub.val_grid_point(1, k, j, 0);
00491              
00492       }
00493 
00494     //  maxabs (mufin2);
00495 
00496 
00497     Scalar brfin = shift_new(1);
00498     brfin.set_spectral_va().coef();
00499          
00500     Scalar brfin2 (map_2) ;
00501     brfin2.allocate_all(); 
00502     brfin2.set_etat_qcq();  
00503     brfin2.std_spectral_base();
00504          
00505     for (int k=0; k<np; k++)
00506       for (int j=0; j<nt; j++) {
00507              
00508     brfin2.set_grid_point(0, k, j, 0) = brfin.val_grid_point(1, k,j,0) - Vrb.val_grid_point(1, k, j, 0);
00509              
00510       }
00511     //  maxabs (brfin2);
00512 
00513 
00514 
00516 
00517     // Update during the loop
00518     for (int ii=1; ii <=3; ii++){
00519      shift.set(ii) = shift_new(ii)*(1-relax) + shift(ii)* relax;
00520     }
00521        
00522 
00523     diff_ent = max(maxabs(npsi_new - npsi )); // Convergence parameter (discutable relevance...) 
00524 
00525 
00526 
00530 
00531 
00532   
00534     // Tensor hij resolution        ///
00536 
00537     if (isCF == false){
00538    
00539     if (diff_ent <=5.e-3) { // No resolution until we are close to the result.
00540    
00541    //WARNING; parameter maybe to be changed according to convergence parameters in Poisson-Hole/Kerr1.9/pplncp.C
00542       util = util+1; // Loop marker for NCF equation.
00543      
00544 
00548 
00549      // Local convergence can be asked for hij equation. Allows to satisfy integrability conditions.
00550   double precis2 =  1.e5*precis ; // Here we ask for a local convergence; this can be improved later. WARNING; parameter maybe to be changed according to convergence parameters in Poisson-Hole/Kerr1.9/pplncp.C
00551 
00552   double diff_ent2 = 1 ; // Local convergence marker
00553   double relax2 = 0.5; // Local relaxation parameter. If not 1, the determinant condition won't be satisfied on a particular iteration.
00554 
00555 
00556   Sym_tensor sourcehij = hij; // Random initialization...
00557 
00558   for(int mer2=0 ;(diff_ent2 > precis2) && (mer2<mer_max2) ; mer2++) {    
00559 
00561 
00562     // Calculation of the source
00563 
00564     // The double Lie derivative term is taken care of in the subroutine. 
00565 
00566        secmembre_kerr(sourcehij);
00567        
00568       
00570     //System inversion (note that no boundary condition is imposed)
00571 
00572 
00573      Sym_tensor hij_new = hij;
00574 
00575      hij_new = boundfree_tensBC (sourcehij, shift , conf_fact, lapse, hij, precis2);
00576 
00577      cout << "maximum of convergence marker for the subiteration" << endl;
00578 
00579      diff_ent2 = max(maxabs(hij - hij_new));
00580       hij = relax2*hij_new + (1 - relax2)*hij;
00581 
00582      cout << "mer2, diffent2" << endl;
00583  
00584       cout << mer2 << endl;
00585       cout << diff_ent2 << endl;
00586 
00587 
00588   }
00589 
00591    // Resolution tests
00592 
00593 
00594        Sym_tensor gammatilde = mets.con() + hij;   
00595        Metric gammatilde2(gammatilde); Scalar detgam = gammatilde2.determinant();
00596  //       cout << "determinant of result" << endl;
00597 //       maxabs (detgam-1.);
00598    
00599 //      cout << "comment l'equation en hij est elle vérifiée?" << endl;
00600         
00601       Sym_tensor test =contract (hij.derive_cov(mets).derive_con(mets), 2,3);
00602       test.annule(nz-1, nz-1);
00603       test = test - hij.derive_lie(shift).derive_lie(shift)/ ((lapse/(conf_fact*conf_fact))*(lapse/(conf_fact*conf_fact)));        
00604       test.annule(nz-1, nz-1);
00605       Sym_tensor youps = test - sourcehij/((lapse/(conf_fact*conf_fact))*(lapse/(conf_fact*conf_fact)));
00606 //       maxabs (youps); 
00607 //       maxabs((youps).trace(mets));
00608 //       cout << " AAABBB" << endl;
00609 //       maxabs((youps).compute_A());
00610 //       maxabs((youps).compute_tilde_B());
00611 
00612     
00613     }
00614     }
00615  
00616 
00617 
00618 
00622 
00623     
00624 
00626      // Global variable update after an entire loop ////
00628 
00629    
00630        gamtuu = mets.con() + hij;
00631        gamt = gamtuu; // Métrique.
00632        gam = gamt.cov()*psi4;
00633        gamma = gam.cov();
00634        
00635 
00636        for (int i=1; i<=3; i++) {
00637      for (int j=1; j<=3; j++) {
00638        
00639        tmp = 0;       
00640        tmp = ((shift.derive_con(mets))(i,j) + (shift.derive_con(mets))(j,i) - (2./3.)*shift.divergence(mets) * mets.con()(i,j))*(1./(2.*lapse));       
00641        aa.set(i,j) = tmp ; 
00642        
00643      }
00644        }
00645 
00646        aa = aa -  (hij.derive_lie(shift) + (2./3.)*shift.divergence(mets)*hij)*(1./(2.*lapse)); //Non conformally flat correction; we assume here dhij/dt = 0.
00647   
00648        aa_hat = aa*psi4*sqrt(psi4); // Rescaling of traceless exrinsic curvature.
00649        aa_hat.std_spectral_base();
00650 
00651 
00652        hatA = aa_hat; // Probably obsolete very soon... replace aa_hat by hatA.
00653 
00654        Sym_tensor aaud = aa.up_down(gamt);
00655        Sym_tensor aaud_hat = aa_hat.up_down(gamt);
00656        aa_quad_scal =  contract(contract (aa_hat, 0, aaud_hat, 0), 0,1);
00657 
00658        delta = -0.5*contract( hij, 1, gamt.cov().derive_cov(mets), 2);
00659  
00660 
00661    for (int i=1; i<=3; i++) {
00662      for (int j=1; j<=3; j++) {
00663        for (int k=1; k<=3; k++) {
00664      tmp = 0.;
00665      tmp =  -0.5 *(gamt.cov().derive_con(mets))(i,j,k);
00666      for (int l=1; l<=3; l++) {
00667        tmp += -0.5*( gamt.cov()(i,l)*(hij.derive_cov(mets))(k,l,j) + gamt.cov()(l,j)*(hij.derive_cov(mets))(k,l,i));
00668      }
00669          delta.set(k,i,j) += tmp ; 
00670        }
00671      }
00672    } 
00673    
00674   Rstar = 
00675     0.25 * contract(gamt.con(), 0, 1,
00676             contract(gamt.con().derive_cov(mets), 0, 1, gamt.cov().derive_cov(mets), 0, 1), 0, 1 ) 
00677     - 0.5  * contract(gamt.con(), 0, 1,
00678               contract(gamt.con().derive_cov(mets), 0, 1, gamt.cov().derive_cov(mets), 0, 2), 0, 1 ) ;  
00679 
00680     kuu = aa/(psi4); 
00681     kdd =  contract (gamma, 0, contract(gamma, 1, kuu, 0),1);  
00682  
00683 
00684 
00686 
00687   // Convergence markers
00688     cout << "diffent" << endl; 
00689     cout<< diff_ent << endl;   
00690     cout <<"mer" << mer << endl;   
00691    
00692 
00693 
00694 
00698      
00699   //------------------------------------------------
00700   //     Check of Einstein equations (3+1 form)
00701   //------------------------------------------------
00702 
00703 
00704   // Lapse 
00705   //------
00706   Scalar lapse2(*map) ;
00707   lapse2 = lapse ;
00708   lapse2.std_spectral_base();
00709 
00710   // 3-metric
00711   //---------
00712   //   const Metric gam(mets.cov()*psi4) ;
00713    
00714   Sym_tensor gamt2(*map, COV, (*map).get_bvect_spher()); 
00715   for (int i=1; i<=3; i++)
00716     for (int j=1; j<=3; j++) 
00717       { gamt2.set(i,j)=gam.cov()(i,j); 
00718       }
00719   //Shift
00720   //-----
00721   Vector beta = shift ;
00722   
00723   // Extrinsic curvature
00724   //--------------------
00725 
00726   Scalar TrK3(*map);
00727   Sym_tensor k_uu = aa/(psi4) ;
00728   k_uu.dec_dzpuis(k_uu(1,1).get_dzpuis()); 
00729   Sym_tensor k_dd = k_uu.up_down(gam); // Another way of computing the same thing, just to be sure... 
00730 
00731   TrK3 = k_uu.trace(gam);
00732   //  TrK3.spectral_display("TraceKvraie", 1.e-10);
00733  
00734   // Hamiltonian constraint
00735   //-----------------------
00736   Scalar ham_constr = gam.ricci_scal() ;
00737   ham_constr.dec_dzpuis(3) ;
00738   ham_constr +=  TrK3*TrK3 - contract(k_uu, 0, 1, k_dd, 0, 1) ;
00739   // maxabs(ham_constr, "Hamiltonian constraint: ") ;
00740 
00741   ham_constr.set_spectral_va().ylm();
00742   //  ham_constr.spectral_display("ham_constr", 1.e-9);
00743   
00744 
00745  
00746   // Momentum constraint
00747   //-------------------
00748   Vector mom_constr = k_uu.divergence(gam)  - TrK3.derive_con(gam) ;
00749   mom_constr.dec_dzpuis(2) ;
00750 //   maxabs(mom_constr, "Momentum constraint: ") ;
00751 //   mom_constr(1).spectral_display("mom1", 1.e-9) ;
00752 //  mom_constr(2).spectral_display("mom2", 1.e-9) ;
00753 //  mom_constr(3).spectral_display("mom3", 1.e-9) ;
00754   
00755  
00756   // Evolution equations
00757   //--------------------
00758   Sym_tensor evol_eq = lapse2*gam.ricci() 
00759     - lapse2.derive_cov(gam).derive_cov(gam);
00760   evol_eq.dec_dzpuis() ;
00761   evol_eq += k_dd.derive_lie(beta) ;
00762   evol_eq.dec_dzpuis(2) ;
00763   evol_eq += lapse2*(TrK3*k_dd - 2*contract(k_dd, 1, k_dd.up(0, gam), 0) ) ;
00764 //   maxabs(evol_eq, "Evolution equations: ") ;
00765   
00766 //   evol_eq.trace(gam).spectral_display("evoltrace", 1.e-10);
00767 //   maxabs (evol_eq.trace(gam));
00768  
00769  
00770   //  evol_eq.spectral_display("evol", 1.e-10);
00771   
00772 
00773   }
00774 
00775      
00776 
00777     cout << "================================================================" << endl;
00778     cout << "                THE ITERATION HAS NOW CONVERGED" << endl;
00779     cout << "================================================================" << endl;
00780      }
00781   return; 
00782 }
00783 
00784 
00785 
00786 
00787 
00788 
00789 
00790 
00791 
00792 
00793 

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