binary_helical.C

00001 /*
00002  * Methods of Bin_star::helical
00003  *
00004  * (see file star.h for documentation)
00005  *
00006  */
00007 
00008 /*
00009  *   Copyright (c) 2006 Francois Limousin
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 char binary_helical_C[] = "$Header: /cvsroot/Lorene/C++/Source/Binary/binary_helical.C,v 1.4 2008/08/19 06:41:59 j_novak Exp $" ;
00030 
00031 /*
00032  * $Id: binary_helical.C,v 1.4 2008/08/19 06:41:59 j_novak Exp $
00033  * $Log: binary_helical.C,v $
00034  * Revision 1.4  2008/08/19 06:41:59  j_novak
00035  * Minor modifications to avoid warnings with gcc 4.3. Most of them concern
00036  * cast-type operations, and constant strings that must be defined as const char*
00037  *
00038  * Revision 1.3  2006/08/01 14:26:50  f_limousin
00039  * Small changes
00040  *
00041  * Revision 1.2  2006/06/05 17:05:57  f_limousin
00042  * *** empty log message ***
00043  *
00044  * Revision 1.1  2006/04/11 14:25:15  f_limousin
00045  * New version of the code : improvement of the computation of some
00046  * critical sources, estimation of the dirac gauge, helical symmetry...
00047  *
00048 
00049  *
00050  * $Header: /cvsroot/Lorene/C++/Source/Binary/binary_helical.C,v 1.4 2008/08/19 06:41:59 j_novak Exp $ *
00051  */
00052 
00053 
00054 // C headers
00055 #include <math.h>
00056 
00057 // Lorene headers
00058 #include "cmp.h"
00059 #include "tenseur.h"
00060 #include "metrique.h"
00061 #include "binary.h"
00062 #include "param.h"
00063 #include "graphique.h"
00064 #include "utilitaires.h"
00065 #include "tensor.h"
00066 #include "nbr_spx.h"
00067 #include "unites.h"
00068 
00069 void Binary::helical(){
00070 
00071     // Fundamental constants and units
00072     // -------------------------------
00073     using namespace Unites ;
00074 
00075  
00076     Sym_tensor lie_aij_1 (star1.mp, CON, star1.mp.get_bvect_cart()) ;
00077     Sym_tensor lie_aij_2 (star2.mp, CON, star2.mp.get_bvect_cart()) ;
00078 
00079     Scalar lie_K_1 (star1.mp) ;
00080     Scalar lie_K_2 (star2.mp) ;
00081 
00082     for (int ll=1; ll<=2; ll++) {
00083     
00084     Star_bin star_i (*et[ll-1]) ;
00085     
00086     Map& mp = star_i.mp ;
00087     const Mg3d* mg = mp.get_mg() ; 
00088     int nz = mg->get_nzone() ;      // total number of domains
00089     
00090     Metric_flat flat = star_i.flat ;
00091     Metric gtilde = star_i.gtilde ;
00092     Scalar nn = star_i.nn ;
00093     Scalar psi4 = star_i.psi4 ;
00094 
00095     // -------------------------------
00096     // AUXILIARY QUANTITIES
00097     // -------------------------------
00098 
00099     // Derivatives of N and logN
00100     //--------------------------
00101 
00102     const Vector dcov_logn_auto = star_i.logn_auto.derive_cov(flat) ;
00103     
00104     Tensor dcovdcov_logn_auto = (star_i.logn_auto.derive_cov(flat))
00105                                       .derive_cov(flat) ;
00106     dcovdcov_logn_auto.inc_dzpuis() ;
00107 
00108     // Derivatives of lnq, phi and Q
00109     //-------------------------------
00110 
00111     const Scalar phi (0.5 * (star_i.lnq - star_i.logn)) ;
00112     const Scalar phi_auto (0.5 * (star_i.lnq_auto - star_i.logn_auto)) ;
00113 
00114     const Vector dcov_phi_auto = phi_auto.derive_cov(flat) ;
00115     
00116     const Vector dcov_lnq = 2*star_i.dcov_phi + star_i.dcov_logn ;
00117     const Vector dcon_lnq = 2*star_i.dcon_phi + star_i.dcon_logn ;
00118     const Vector dcov_lnq_auto = star_i.lnq_auto.derive_cov(flat) ;
00119         Tensor dcovdcov_lnq_auto = dcov_lnq_auto.derive_cov(flat) ;
00120     dcovdcov_lnq_auto.inc_dzpuis() ;
00121 
00122     Scalar qq = exp(star_i.lnq) ;
00123     qq.std_spectral_base() ;
00124     const Vector& dcov_qq = qq.derive_cov(flat) ;
00125 
00126     Tensor dcovdcov_beta_auto = star_i.beta_auto.derive_cov(flat)
00127       .derive_cov(flat) ;
00128     dcovdcov_beta_auto.inc_dzpuis() ;
00129 
00130 
00131     // Derivatives of hij, gtilde... 
00132     //------------------------------
00133 
00134     Scalar psi2 (pow(star_i.psi4, 0.5)) ;
00135     psi2.std_spectral_base() ;
00136 
00137     const Tensor& dcov_hij = star_i.hij.derive_cov(flat) ;
00138     const Tensor& dcon_hij = star_i.hij.derive_con(flat) ;
00139     const Tensor& dcov_hij_auto = star_i.hij_auto.derive_cov(flat) ;
00140 
00141     const Sym_tensor& gtilde_cov = star_i.gtilde.cov() ;
00142     const Sym_tensor& gtilde_con = star_i.gtilde.con() ;
00143     const Tensor& dcov_gtilde = gtilde_cov.derive_cov(flat) ;
00144 
00145     // H^i and its derivatives ( = O in Dirac gauge)
00146     // ---------------------------------------------
00147 
00148     double lambda_dirac = 0. ;
00149 
00150     const Vector hdirac = lambda_dirac * star_i.hij.divergence(flat) ;
00151     const Vector hdirac_auto = lambda_dirac * 
00152         star_i.hij_auto.divergence(flat) ;
00153 
00154     Tensor dcov_hdirac = lambda_dirac * hdirac.derive_cov(flat) ;
00155     dcov_hdirac.inc_dzpuis() ;
00156     Tensor dcov_hdirac_auto = lambda_dirac * hdirac_auto.derive_cov(flat) ;
00157     dcov_hdirac_auto.inc_dzpuis() ;
00158     Tensor dcon_hdirac_auto = lambda_dirac * hdirac_auto.derive_con(flat) ;
00159     dcon_hdirac_auto.inc_dzpuis() ;
00160 
00161 
00162     // Function exp(-(r-r_0)^2/sigma^2)
00163     // --------------------------------
00164     
00165     double r0 = mp.val_r(nz-2, 1, 0, 0) ;
00166     double sigma = 1.*r0 ;
00167     double om = omega ;
00168     
00169     Scalar rr (mp) ;
00170     rr = mp.r ;
00171     
00172     Scalar ff (mp) ;
00173     ff = exp( -(rr - r0)*(rr - r0)/sigma/sigma ) ;
00174     for (int ii=0; ii<nz-1; ii++)
00175         ff.set_domain(ii) = 1. ;
00176     ff.set_outer_boundary(nz-1, 0) ;
00177     ff.std_spectral_base() ;
00178     
00179 
00180     // omdsdp
00181     // ---------------------------------
00182 
00183     Vector omdsdp (mp, CON, mp.get_bvect_cart()) ;
00184     Scalar yya (mp) ;
00185     yya = mp.ya ;
00186     Scalar xxa (mp) ;
00187     xxa = mp.xa ;
00188     Scalar zza (mp) ;
00189     zza = mp.za ;
00190 
00191     if (fabs(mp.get_rot_phi()) < 1e-10){ 
00192         omdsdp.set(1) = - om * yya * ff ;
00193         omdsdp.set(2) = om * xxa * ff ;
00194         omdsdp.set(3).annule_hard() ;
00195     }
00196     else{
00197         omdsdp.set(1) = om * yya * ff ;
00198         omdsdp.set(2) = - om * xxa * ff ;
00199         omdsdp.set(3).annule_hard() ;
00200     }
00201  
00202     omdsdp.set(1).set_outer_boundary(nz-1, 0) ;
00203     omdsdp.set(2).set_outer_boundary(nz-1, 0) ;     
00204     omdsdp.std_spectral_base() ;
00205 
00206 
00207     // Computation of helical A^{ij}
00208     // ------------------------------
00209 
00210     const Tensor& dbeta = star_i.beta_auto.derive_con(gtilde) ;
00211     Scalar div_beta = star_i.beta_auto.divergence(gtilde) ;
00212     
00213     Sym_tensor tkij_a (star_i.tkij_auto) ;
00214     for (int i=1; i<=3; i++) 
00215         for (int j=1; j<=i; j++) {
00216         tkij_a.set(i, j) = dbeta(i, j) + dbeta(j, i) - 
00217             double(2) /double(3) * div_beta * (gtilde.con())(i,j) ; 
00218         }
00219     
00220     tkij_a = tkij_a - star_i.hij_auto.derive_lie(omdsdp) ;  
00221     tkij_a = 0.5 * tkij_a / nn ;   
00222 
00223     Sym_tensor tkij_auto_cov = tkij_a.up_down(gtilde) ;
00224     
00225     Scalar a2_auto = contract(tkij_auto_cov, 0, 1, tkij_a, 0, 1,true) ; 
00226 
00227     // COMP 
00228     const Tensor& dbeta_comp = star_i.beta_comp.derive_con(gtilde) ;
00229         Scalar divbeta_comp = star_i.beta_comp.divergence(gtilde) ;
00230     
00231     Sym_tensor tkij_c (star_i.tkij_comp) ;
00232     for (int i=1; i<=3; i++) 
00233         for (int j=i; j<=3; j++) {      
00234         tkij_c.set(i, j) = dbeta_comp(i, j) + dbeta_comp(j, i) - 
00235             double(2) /double(3) * divbeta_comp * (gtilde.con())(i,j) ; 
00236         }
00237     
00238     tkij_c = tkij_c - star_i.hij_comp.derive_lie(omdsdp) ;  
00239     tkij_c = 0.5 * tkij_c / nn ;
00240     
00241     Scalar a2_comp = contract(tkij_auto_cov, 0, 1, tkij_c, 0, 1,true) ; 
00242 
00243 
00244 //  tkij_a = star_i.tkij_auto ;
00245 //  tkij_c = star_i.tkij_comp ;
00246 
00247 
00248     // Sources for N
00249     // ---------------
00250 
00251     Scalar source1(mp) ;
00252     Scalar source2(mp) ;
00253     Scalar source3(mp) ;
00254     Scalar source4(mp) ;
00255     Scalar source_tot(mp) ;
00256 
00257     source1 = qpig * star_i.psi4 % (star_i.ener_euler + star_i.s_euler) ; 
00258 
00259     source2 = star_i.psi4 % (a2_auto + a2_comp) ;
00260 
00261     source3 = - contract(dcov_logn_auto, 0, star_i.dcon_logn, 0, true) 
00262         - 2. * contract(contract(gtilde_con, 0, star_i.dcov_phi, 0), 
00263                 0, dcov_logn_auto, 0, true) ;
00264             
00265     source4 = - contract(star_i.hij, 0, 1, dcovdcov_logn_auto + 
00266                  dcov_logn_auto*star_i.dcov_logn, 0, 1) ;
00267 
00268     source_tot = source1 + source2 + source3 + source4 ;
00269     
00270     if (ll == 1) {
00271     lie_K_1 = nn / star_i.psi4 * (source_tot - star_i.logn_auto
00272                        .laplacian()) ;
00273     lie_K_1.dec_dzpuis(4) ;
00274     }
00275     if (ll == 2) {
00276     lie_K_2 = nn / star_i.psi4 * (source_tot - star_i.logn_auto
00277                        .laplacian()) ;
00278     lie_K_2.dec_dzpuis(4) ;
00279     }
00280 
00281 
00282     // Sources for hij
00283     // --------------
00284 
00285     Scalar source_tot_hij(mp) ;
00286     Tensor source_Sij(mp, 2, CON, mp.get_bvect_cart()) ;
00287     Tensor source_Rij(mp, 2, CON, mp.get_bvect_cart()) ;
00288     Tensor tens_temp(mp, 2, CON, mp.get_bvect_cart()) ;
00289     
00290     Tensor source_1 (mp, 2, CON, mp.get_bvect_cart()) ;
00291     Tensor source_2 (mp, 2, CON, mp.get_bvect_cart()) ;
00292     Tensor source_3a (mp, 2, CON, mp.get_bvect_cart()) ;
00293     Tensor source_3b (mp, 2, CON, mp.get_bvect_cart()) ;
00294     Tensor source_4 (mp, 2, CON, mp.get_bvect_cart()) ;
00295     Tensor source_5 (mp, 2, CON, mp.get_bvect_cart()) ;
00296     Tensor source_6 (mp, 2, CON, mp.get_bvect_cart()) ;
00297     
00298 
00299     source_1 = contract(dcon_hij, 1, dcov_lnq_auto, 0) ;
00300     
00301     source_2 = - contract(dcon_hij, 2, dcov_lnq_auto, 0) 
00302         - 2./3. * contract(hdirac, 0, dcov_lnq_auto, 0) * flat.con() ;
00303         
00304     // Lie derivative of A^{ij}
00305     // --------------------------
00306 
00307     Scalar decouple_logn = (star_i.logn_auto - 1.e-8)/
00308         (star_i.logn - 2.e-8) ;
00309 
00310     // Construction of Omega d/dphi
00311     // ----------------------------
00312         
00313     // Construction of D_k \Phi^i
00314     Itbl type (2) ;
00315     type.set(0) = CON ;
00316     type.set(1) = COV ;
00317 
00318     Tensor dcov_omdsdphi (mp, 2, type, mp.get_bvect_cart()) ;
00319     dcov_omdsdphi.set(1,1) = 0. ;
00320     dcov_omdsdphi.set(2,1) = om * ff ;
00321     dcov_omdsdphi.set(3,1) = 0. ;
00322     dcov_omdsdphi.set(2,2) = 0. ;
00323     dcov_omdsdphi.set(3,2) = 0. ;
00324     dcov_omdsdphi.set(3,3) = 0. ;
00325     dcov_omdsdphi.set(1,2) = -om *ff ;
00326     dcov_omdsdphi.set(1,3) = 0. ;
00327     dcov_omdsdphi.set(2,3) = 0. ;
00328     dcov_omdsdphi.std_spectral_base() ;
00329 
00330     source_3a = contract(tkij_a, 0, dcov_omdsdphi, 1) ;
00331     source_3a.inc_dzpuis() ;
00332 
00333     // Source 3b
00334     // ------------
00335 
00336     source_3b = - contract(omdsdp, 0, tkij_a.derive_cov(flat), 2) ; 
00337 
00338  
00339     // Source 4
00340     // ---------
00341 
00342     source_4 = - tkij_a.derive_lie(star_i.beta) ;
00343     source_4.inc_dzpuis() ;
00344     source_4 += - 2./3. * star_i.beta.divergence(flat) * tkij_a ;
00345 
00346     source_5 = dcon_hdirac_auto ;
00347         
00348     source_6 = - 2./3. * hdirac_auto.divergence(flat) * flat.con() ;
00349     source_6.inc_dzpuis() ;
00350 
00351     // Source terms for Sij
00352     //---------------------
00353         
00354     source_Sij = 8. * nn / psi4 * phi_auto.derive_con(gtilde) * 
00355         contract(gtilde_con, 0, star_i.dcov_phi, 0) ;
00356 
00357     source_Sij += 4. / psi4 * phi_auto.derive_con(gtilde) * 
00358         nn * contract(gtilde_con, 0, star_i.dcov_logn, 0) +
00359         4. / psi4 * nn * contract(gtilde_con, 0, star_i.dcov_logn, 0) *
00360         phi_auto.derive_con(gtilde) ;
00361 
00362     source_Sij += - nn / (3.*psi4) * gtilde_con * 
00363         ( 0.25 * contract(gtilde_con, 0, 1, contract(dcov_hij_auto, 0, 1,
00364                            dcov_gtilde, 0, 1), 0, 1)
00365           - 0.5 * contract(gtilde_con, 0, 1, contract(dcov_hij_auto, 0, 1,
00366                            dcov_gtilde, 0, 2), 0, 1)) ;
00367 
00368     source_Sij += - 8.*nn / (3.*psi4) * gtilde_con * 
00369      contract(dcov_phi_auto, 0, contract(gtilde_con, 0, star_i.dcov_phi, 0), 0) ;
00370 
00371     tens_temp = nn / (3.*psi4) * hdirac.divergence(flat)*star_i.hij_auto ;
00372     tens_temp.inc_dzpuis() ;
00373 
00374     source_Sij += tens_temp ;
00375        
00376     source_Sij += - 8./(3.*psi4) * contract(dcov_phi_auto, 0,
00377        nn*contract(gtilde_con, 0, star_i.dcov_logn, 0), 0) * gtilde_con ;
00378 
00379     source_Sij += 2.*nn* contract(gtilde_cov, 0, 1, tkij_a *
00380                 (tkij_a+tkij_c), 1, 3) ;
00381 
00382     source_Sij += - 2. * qpig * nn * ( psi4 * star_i.stress_euler 
00383               - 0.33333333333333333 * star_i.s_euler * gtilde_con ) ; 
00384         
00385     source_Sij += - 1./(psi4*psi2) * contract(gtilde_con, 1, 
00386             contract(gtilde_con, 1, qq*dcovdcov_lnq_auto + 
00387                  qq*dcov_lnq_auto*dcov_lnq, 0), 1) ;
00388 
00389     source_Sij += - 0.5/(psi4*psi2) * contract(contract(star_i.hij, 1,
00390                     dcov_hij_auto, 2), 1, dcov_qq, 0) -
00391         0.5/(psi4*psi2) * contract(contract(dcov_hij_auto, 2,
00392                     star_i.hij, 1), 1, dcov_qq, 0) ;
00393                     
00394     source_Sij += 0.5/(psi4*psi2) * contract(contract(star_i.hij, 0,
00395                      dcov_hij_auto, 2), 0, dcov_qq, 0) ;
00396 
00397     source_Sij += 1./(3.*psi4*psi2)*contract(gtilde_con, 0, 1,
00398        qq*dcovdcov_lnq_auto + qq*dcov_lnq_auto*dcov_lnq, 0, 1)
00399         *gtilde_con ;
00400 
00401     source_Sij += 1./(3.*psi4*psi2) * contract(hdirac, 0, 
00402                         dcov_qq, 0)*star_i.hij_auto ;
00403 
00404     // Source terms for Rij
00405     //---------------------
00406 
00407     source_Rij = contract(star_i.hij, 0, 1, dcov_hij_auto.derive_cov(flat),
00408                   2, 3) ;
00409     source_Rij.inc_dzpuis() ;
00410 
00411 
00412     source_Rij += - contract(star_i.hij_auto, 1, dcov_hdirac, 1) -
00413         contract(dcov_hdirac, 1, star_i.hij_auto, 1) ;
00414         
00415     source_Rij += contract(hdirac, 0, dcov_hij_auto, 2) ;
00416 
00417     source_Rij += - contract(contract(dcov_hij_auto, 1, dcov_hij, 2),
00418                  1, 3) ;
00419 
00420     source_Rij += - contract(gtilde_cov, 0, 1, contract(contract(
00421         gtilde_con, 0, dcov_hij_auto, 2), 0, dcov_hij, 2), 1, 3) ;
00422 
00423     source_Rij += contract(contract(contract(contract(gtilde_cov, 0, 
00424           dcov_hij_auto, 1), 2, gtilde_con, 1), 0, dcov_hij, 1), 0, 3) +
00425         contract(contract(contract(contract(gtilde_cov, 0, 
00426         dcov_hij_auto, 1), 0, dcov_hij, 1), 0, 3), 0, gtilde_con, 1) ;
00427 
00428     source_Rij += 0.5 * contract(gtilde_con*gtilde_con, 1, 3, 
00429          contract(dcov_hij_auto, 0, 1, dcov_gtilde, 0, 1), 0, 1) ;
00430 
00431     source_Rij = source_Rij * 0.5 ;
00432 
00433     for(int i=1; i<=3; i++) 
00434         for(int j=1; j<=i; j++) {
00435 
00436         source_tot_hij = source_1(i,j) + source_1(j,i) 
00437             + source_2(i,j) + 2.*psi4/nn * (
00438             source_4(i,j) - source_Sij(i,j)) 
00439             - 2.* source_Rij(i,j) +
00440             source_5(i,j) + source_5(j,i) + source_6(i,j) ;
00441         source_tot_hij.dec_dzpuis() ;
00442             
00443         source3 = 2.*psi4/nn * (source_3a(i,j) + source_3a(j,i)  
00444                  + source_3b(i,j)) ; 
00445 
00446         source_tot_hij = source_tot_hij + source3 ;
00447 
00448 
00449         if (ll == 1){
00450             if(i==1 && j==1) {
00451             Scalar lapl (star_i.hij_auto(1,1).laplacian()) ;
00452             lapl.dec_dzpuis() ;
00453             lie_aij_1.set(1,1) = nn / (2.*psi4) * 
00454                 (lapl-source_tot_hij) ;
00455             }
00456             if(i==2 && j==1) {
00457             Scalar lapl (star_i.hij_auto(2,1).laplacian()) ;
00458             lapl.dec_dzpuis() ;
00459             lie_aij_1.set(2,1) = nn / (2.*psi4) * 
00460                 (lapl-source_tot_hij) ;
00461             }
00462             if(i==3 && j==1) {
00463             Scalar lapl (star_i.hij_auto(3,1).laplacian()) ;
00464             lapl.dec_dzpuis() ;
00465             lie_aij_1.set(3,1) = nn / (2.*psi4) * 
00466                 (lapl-source_tot_hij) ;
00467             }
00468             if(i==2 && j==2) {
00469             Scalar lapl (star_i.hij_auto(2,2).laplacian()) ;
00470             lapl.dec_dzpuis() ;
00471             lie_aij_1.set(2,2) = nn / (2.*psi4) * 
00472                 (lapl-source_tot_hij) ;
00473             }
00474             if(i==3 && j==2) {
00475             Scalar lapl (star_i.hij_auto(3,2).laplacian()) ;
00476             lapl.dec_dzpuis() ;
00477             lie_aij_1.set(3,2) = nn / (2.*psi4) * 
00478                 (lapl-source_tot_hij) ;
00479             }
00480             if(i==3 && j==3) {
00481             Scalar lapl (star_i.hij_auto(3,3).laplacian()) ;
00482             lapl.dec_dzpuis() ;
00483             lie_aij_1.set(3,3) = nn / (2.*psi4) * 
00484                 (lapl-source_tot_hij) ;
00485             }
00486         }
00487 
00488         if (ll == 2){
00489             if(i==1 && j==1) {
00490             Scalar lapl (star_i.hij_auto(1,1).laplacian()) ;
00491             lapl.dec_dzpuis() ;
00492             lie_aij_2.set(1,1) = nn / (2.*psi4) * 
00493                 (lapl-source_tot_hij) ;
00494             }
00495             if(i==2 && j==1) {
00496             Scalar lapl (star_i.hij_auto(2,1).laplacian()) ;
00497             lapl.dec_dzpuis() ;
00498             lie_aij_2.set(2,1) = nn / (2.*psi4) * 
00499                 (lapl-source_tot_hij) ;
00500             }
00501             if(i==3 && j==1) {
00502             Scalar lapl (star_i.hij_auto(3,1).laplacian()) ;
00503             lapl.dec_dzpuis() ;
00504             lie_aij_2.set(3,1) = nn / (2.*psi4) * 
00505                 (lapl-source_tot_hij) ;
00506             }
00507             if(i==2 && j==2) {
00508             Scalar lapl (star_i.hij_auto(2,2).laplacian()) ;
00509             lapl.dec_dzpuis() ;
00510             lie_aij_2.set(2,2) = nn / (2.*psi4) * 
00511                 (lapl-source_tot_hij) ;
00512             }
00513             if(i==3 && j==2) {
00514             Scalar lapl (star_i.hij_auto(3,2).laplacian()) ;
00515             lapl.dec_dzpuis() ;
00516             lie_aij_2.set(3,2) = nn / (2.*psi4) * 
00517                 (lapl-source_tot_hij) ;
00518             }
00519             if(i==3 && j==3) {
00520             Scalar lapl (star_i.hij_auto(3,3).laplacian()) ;
00521             lapl.dec_dzpuis() ;
00522             lie_aij_2.set(3,3) = nn / (2.*psi4) * 
00523                 (lapl-source_tot_hij) ;
00524             }
00525         }
00526         }
00527     }
00528 
00529     lie_aij_1.dec_dzpuis(3) ;
00530     lie_aij_2.dec_dzpuis(3) ;
00531         
00532     int nz = star1.mp.get_mg()->get_nzone() ;
00533   
00534     // Construction of an auxiliar grid and mapping (Last domain is at lambda)
00535     double* bornes = new double [6] ;
00536     bornes[nz] = __infinity ;
00537     bornes[4] = M_PI / omega ;
00538     bornes[3] = M_PI / omega * 0.5 ;
00539     bornes[2] = M_PI / omega  * 0.2 ;
00540     bornes[1] = M_PI / omega  * 0.1 ;
00541     bornes[0] = 0 ;
00542     
00543     Map_af mapping (*(star1.mp.get_mg()), bornes) ;
00544     
00545     delete [] bornes ; 
00546  
00547     
00548     Sym_tensor lie_aij2_1 (star1.mp, CON, star1.mp.get_bvect_cart()) ;
00549     Sym_tensor lie_aij_tot_1 (star1.mp, CON, star1.mp.get_bvect_cart()) ;
00550     Sym_tensor lie_aij_tot (mapping, CON, mapping.get_bvect_cart()) ;
00551 
00552     Scalar lie_K2_1 (star1.mp) ;
00553     lie_K2_1.set_etat_qcq() ;
00554     Scalar lie_K_tot_1 (star1.mp) ;
00555 
00556 
00557     // Importation on the mapping 1
00558     // -------------------------------
00559 
00560     lie_K2_1.import(lie_K_2) ;
00561     lie_K_tot_1 = lie_K_1 + lie_K2_1 ;
00562     lie_K_tot_1.inc_dzpuis(2) ;
00563 
00564     lie_aij_2.change_triad(star1.mp.get_bvect_cart()) ;    
00565     for(int i=1; i<=3; i++) 
00566     for(int j=1; j<=i; j++) {
00567         lie_aij2_1.set(i,j).import(lie_aij_2(i,j)) ;
00568         lie_aij2_1.set(i,j).set_spectral_va().set_base(lie_aij_2(i,j).
00569                         get_spectral_va().get_base()) ;
00570     }
00571 
00572     lie_aij_tot_1 = lie_aij_1 + lie_aij2_1 ;
00573     lie_aij_tot_1.inc_dzpuis(2) ;
00574 
00575     
00576     Sym_tensor lie_kij_tot (star1.mp, CON, star1.mp.get_bvect_cart()) ;
00577     lie_kij_tot = lie_aij_tot_1/star1.psi4 + 1./3.*star1.gamma.con()*
00578     lie_K_tot_1 ;
00579 
00580 
00581     cout << " IN THE CENTER OF STAR 1 " << endl 
00582      << " ----------------------- " << endl ;
00583     /*
00584     cout << " components xx, xy, yy, xz, yz, zz" << endl ;
00585     for(int i=1; i<=3; i++) 
00586     for(int j=1; j<=i; j++) {
00587         Scalar resu(lie_kij_tot(i,j)*lie_kij_tot(i,j)) ;
00588         cout << "i = " << i << ", j = " << j << endl ; 
00589         cout << "norme de la diff " << endl 
00590          << norme(resu)/(nr*nt*np) << endl ;
00591 
00592         // Computation of the integral
00593         // -----------------------------
00594 
00595         Tbl integral (nz) ;
00596         integral.annule_hard()  ;
00597         Tbl integ (resu.integrale_domains()) ;
00598         for (int mm=0; mm<nz; mm++) 
00599         for (int pp=0; pp<=mm; pp++) 
00600             integral.set(mm) += integ(pp) ; 
00601         cout << sqrt(integral) / sqrt(omega) << endl ;  // To get 
00602         // dimensionless quantity
00603     }
00604     */
00605 
00606     cout << "L2 norm of L_k K^{ab} " << endl ;
00607     Scalar determinant (pow(star1.get_gamma().determinant(), 0.5)) ;
00608     determinant.std_spectral_base() ;
00609     Scalar resu(2.*contract(lie_kij_tot, 0, 1, 
00610              lie_kij_tot.up_down(star1.gamma), 0, 1)
00611         *determinant) ;
00612     Tbl integral (nz) ;
00613     integral.annule_hard() ;
00614     Tbl integ (resu.integrale_domains()) ;
00615     for (int mm=0; mm<nz; mm++) 
00616     for (int pp=0; pp<=mm; pp++) 
00617         integral.set(mm) += integ(pp) ; 
00618     cout << sqrt(integral) * sqrt(mass_adm()*ggrav) << endl ;  
00619     cout << sqrt(integral) / sqrt(omega) << endl ;  // To get 
00620     // dimensionless quantity
00621 
00622     
00623     cout << "omega = " << omega << endl ;
00624     cout << "mass_adm = " << mass_adm() << endl ;
00625     
00626 
00627     lie_kij_tot.dec_dzpuis(2) ;
00628     
00629     cout << "Position du centre de l'etoile x/lambda = "
00630      << -star1.get_mp().get_ori_x() * omega / M_PI << " in Lorene units"
00631      << endl ;
00632    
00633 
00634 
00635     // Importation on the mapping defined in the center of mass
00636     // -------------------------------------------------------------
00637 /*
00638     lie_aij_tot_1.change_triad(mapping.get_bvect_cart()) ;    
00639     for(int i=1; i<=3; i++) 
00640     for(int j=1; j<=i; j++) {
00641         lie_aij_tot.set(i,j).import(lie_aij_tot_1(i,j)) ;
00642         lie_aij_tot.set(i,j).set_spectral_va().set_base(lie_aij_tot_1(i,j).
00643                         get_spectral_va().get_base()) ;
00644     }
00645 
00646     lie_aij_tot.inc_dzpuis(2) ;
00647 
00648     cout << " IN THE CENTER OF MASS : " << endl 
00649      << " ----------------------- " << endl ;
00650 
00651     cout << " components xx, xy, yy, xz, yz, zz" << endl ;
00652     for(int i=1; i<=3; i++) 
00653     for(int j=1; j<=i; j++) {
00654         Scalar resu(lie_aij_tot(i,j)*lie_aij_tot(i,j)) ;
00655         cout << "i = " << i << ", j = " << j << endl ; 
00656         cout << "norme de la diff " << endl 
00657          << norme(resu)/(nr*nt*np) << endl ;
00658 
00659         Tbl integral (nz) ;
00660         integral.annule_hard()  ;
00661         Tbl integ (resu.integrale_domains()) ;
00662         for (int mm=0; mm<nz; mm++) 
00663         for (int pp=0; pp<=mm; pp++) 
00664             integral.set(mm) += integ(pp) ; 
00665         cout << sqrt(integral) / sqrt(omega) << endl ;  // To get 
00666         // dimensionless quantity
00667     }
00668 
00669     cout << "L2 norm of L_k K^{ab} " << endl ;
00670     resu = contract(lie_aij_tot, 0, 1, 
00671             lie_aij_tot.up_down(star1.gtilde), 0, 1) ;
00672     integral.annule_hard()  ;
00673     integ = resu.integrale_domains() ;
00674     for (int mm=0; mm<nz; mm++) 
00675     for (int pp=0; pp<=mm; pp++) 
00676         integral.set(mm) += integ(pp) ; 
00677     cout << sqrt(integral) / sqrt(omega) << endl ;  // To get 
00678     // dimensionless quantity
00679     */
00680 
00681     cout << "Omega M = " << omega * mass_adm()*ggrav << endl ;
00682 
00683 }

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