00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028 char app_hor_finder_C[] = "$Header: /cvsroot/Lorene/C++/Source/App_hor/app_hor_finder.C,v 1.9 2012/01/02 13:52:57 j_novak Exp $" ;
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047
00048
00049
00050
00051
00052
00053
00054
00055
00056
00057
00058
00059
00060
00061
00062
00063
00064
00065
00066
00067
00068
00069
00070
00071
00072
00073
00074
00075
00076 #include <math.h>
00077 #include <assert.h>
00078
00079
00080 #include "app_hor.h"
00081 #include "graphique.h"
00082
00083
00084 bool ah_finder(const Metric& gamma, const Sym_tensor& k_dd, Valeur& h, Scalar& ex_fcn,
00085 double a_axis, double b_axis, double c_axis, bool verbose, bool print,
00086 double precis, double precis_exp, int step_max, int step_relax,
00087 double relax)
00088 {
00089
00090 bool ah_flag = false ;
00091
00092
00093 const Map& map = gamma.get_mp() ;
00094 const Mg3d* mg = map.get_mg() ;
00095 const Mg3d* g_angu = mg->get_angu() ;
00096 int nz = mg->get_nzone() ;
00097
00098 const Base_vect_spher& bspher = map.get_bvect_spher() ;
00099
00100 const Coord& rr = map.r ;
00101 const Coord& theta = map.tet ;
00102 const Coord& phi = map.phi ;
00103 const Coord& costh = map.cost ;
00104 const Coord& cosph = map.cosp ;
00105 const Coord& sinth = map.sint ;
00106 const Coord& sinph = map.sinp ;
00107
00108 double r_min = min(+rr)(0) ;
00109 double r_max = max(+rr)(nz-1) ;
00110
00111
00112
00113
00114 double aa = a_axis ;
00115 double bb = b_axis ;
00116 double cc = c_axis ;
00117
00118 Scalar ct(map) ;
00119 ct = costh ;
00120 ct.std_spectral_base() ;
00121
00122 Scalar st(map) ;
00123 st = sinth ;
00124 st.std_spectral_base() ;
00125
00126 Scalar cp(map) ;
00127 cp = cosph ;
00128 cp.std_spectral_base() ;
00129
00130 Scalar sp(map) ;
00131 sp = sinph ;
00132 sp.std_spectral_base() ;
00133
00134 Scalar h_tmp(map) ;
00135
00136 h_tmp = st*st*( cp*cp/aa/aa + sp*sp/bb/bb ) + ct*ct/cc/cc ;
00137 h_tmp = 1./sqrt( h_tmp ) ;
00138 h_tmp.std_spectral_base() ;
00139
00140 Valeur h_guess(g_angu) ;
00141 h_guess.annule_hard() ;
00142
00143 for (int l=0; l<nz; l++) {
00144
00145 int jmax = mg->get_nt(l) ;
00146 int kmax = mg->get_np(l) ;
00147
00148 for (int k=0; k<kmax; k++) {
00149 for (int j=0; j<jmax; j++) {
00150 h_guess.set(l,k,j,0) = h_tmp.val_grid_point(l,k,j,0) ;
00151 }
00152 }
00153 }
00154
00155 h_guess.std_base_scal() ;
00156
00157 h = h_guess ;
00158 h.std_base_scal() ;
00159
00160
00161
00162
00163 const Metric_flat& fmets = map.flat_met_spher() ;
00164
00165
00166
00167 double one_third = double(1) / double(3) ;
00168
00169 Scalar psi4 = gamma.determinant() / fmets.determinant() ;
00170 psi4 = pow( psi4, one_third ) ;
00171 psi4.std_spectral_base() ;
00172
00173
00174 Scalar ex_fcn_old(map) ;
00175 ex_fcn_old.set_etat_zero() ;
00176
00177
00178 Valeur ex_AH(g_angu) ;
00179 ex_AH.annule_hard() ;
00180
00181
00182 Vector s_unit(map, CON, bspher) ;
00183
00184 double relax_prev = double(1) - relax ;
00185 double diff_exfcn = 1. ;
00186 Tbl diff_h(nz) ;
00187 diff_h = 1. ;
00188
00189 bool no_AH_in_grid = false ;
00190
00191
00192
00193
00194
00195 for (int step=0 ;
00196 (max(diff_h) > precis) && (step < step_max) && (!no_AH_in_grid);
00197 step++) {
00198
00199
00200
00201
00202
00203
00204
00205 Scalar F(map) ;
00206 F.allocate_all() ;
00207
00208 for (int l=0; l<nz; l++) {
00209
00210 int imax = mg->get_nr(l) ;
00211 int jmax = mg->get_nt(l) ;
00212 int kmax = mg->get_np(l) ;
00213
00214 for (int k=0; k<kmax; k++) {
00215 for (int j=0; j<jmax; j++) {
00216 for (int i=0; i<imax; i++) {
00217
00218
00219 F.set_grid_point(l,k,j,i) = (+rr)(l,k,j,i) - h(l,k,j,0) ;
00220
00221 }
00222 }
00223 }
00224 }
00225
00226 F.std_spectral_base() ;
00227
00228
00229 Scalar dF_norm(map) ;
00230 dF_norm = contract( contract(gamma.con(), 0, F.derive_cov(gamma), 0),
00231 0, F.derive_cov(gamma), 0) ;
00232 dF_norm = sqrt( dF_norm ) ;
00233 dF_norm.std_spectral_base() ;
00234
00235 s_unit = F.derive_con(gamma) / dF_norm ;
00236
00237
00238 ex_fcn = s_unit.divergence(gamma) - k_dd.trace(gamma) +
00239 contract( s_unit, 0, contract(s_unit, 0, k_dd, 1), 0) ;
00240
00241
00242
00243
00244 Sym_tensor sou_1(map, CON, bspher) ;
00245 sou_1 = gamma.con() - fmets.con()/psi4 - s_unit*s_unit ;
00246
00247 Scalar source_tmp(map) ;
00248 source_tmp = contract( sou_1, 0, 1, F.derive_cov(fmets).derive_cov(fmets), 0, 1 ) ;
00249 source_tmp = source_tmp / dF_norm ;
00250
00251 Sym_tensor d_gam(map, COV, bspher) ;
00252 d_gam = contract( gamma.connect().get_delta(), 0, F.derive_cov(fmets), 0) ;
00253
00254 source_tmp -= contract( gamma.con() - s_unit*s_unit, 0, 1,
00255 d_gam, 0, 1) / dF_norm ;
00256
00257 source_tmp = psi4*dF_norm*( source_tmp - k_dd.trace(gamma) +
00258 contract(s_unit, 0, contract(s_unit, 0, k_dd, 1), 0) ) ;
00259
00260 source_tmp.std_spectral_base() ;
00261
00262
00263 Valeur sou_angu(g_angu) ;
00264
00265 sou_angu.annule_hard() ;
00266
00267 double h_min = min(h)(0) ;
00268 double h_max = max(h)(0) ;
00269 if ( (r_min < h_min) && (h_max < r_max) ) {
00270
00271 for (int l=0; l<nz; l++) {
00272
00273 int jmax = mg->get_nt(l) ;
00274 int kmax = mg->get_np(l) ;
00275 for (int k=0; k<kmax; k++) {
00276 for (int j=0; j<jmax; j++) {
00277 sou_angu.set(l,k,j,0) = source_tmp.val_point(h(l,k,j,0)
00278 ,(+theta)(l,k,j,0) ,(+phi)(l,k,j,0)) ;
00279 }
00280 }
00281 }
00282 sou_angu = h*h*sou_angu ;
00283 }
00284 else {
00285 no_AH_in_grid = true ;
00286 break ;
00287 }
00288 sou_angu.std_base_scal() ;
00289
00290
00291
00292
00293
00294
00295
00296
00297
00298 sou_angu.ylm() ;
00299
00300 Valeur h_new = sou_angu ;
00301
00302 h_new.c_cf->poisson_angu(-2.) ;
00303
00304 h_new.ylm_i() ;
00305
00306 if (h_new.c != 0x0)
00307 delete h_new.c ;
00308 h_new.c = 0x0 ;
00309 h_new.coef_i() ;
00310
00311
00312 diff_h = max(abs(h - h_new)) ;
00313
00314
00315
00316 if (step >= step_relax) {
00317 h_new = relax * h_new + relax_prev * h ;
00318 }
00319
00320
00321 h = h_new ;
00322
00323
00324 if (print)
00325 {
00326
00327 cout << "-------------------------------------" << endl ;
00328 cout << "App_hor iteration step: " << step << endl ;
00329 cout << " " << endl ;
00330
00331 cout << "Difference in h : " << diff_h << endl ;
00332
00333
00334 Tbl diff_exfcn_tbl = diffrel( ex_fcn, ex_fcn_old ) ;
00335 diff_exfcn = diff_exfcn_tbl(0) ;
00336 for (int l=1; l<nz; l++) {
00337 diff_exfcn += diff_exfcn_tbl(l) ;
00338 }
00339 diff_exfcn /= nz ;
00340 cout << "diff_exfcn : " << diff_exfcn << endl ;
00341
00342 ex_fcn_old = ex_fcn ;
00343
00344
00345 }
00346
00347 if ( (step == step_max-1) && (max(diff_h) > precis) ) {
00348
00349
00350
00351
00352 for (int l=0; l<nz; l++) {
00353
00354 int jmax = mg->get_nt(l) ;
00355 int kmax = mg->get_np(l) ;
00356
00357 for (int k=0; k<kmax; k++) {
00358 for (int j=0; j<jmax; j++) {
00359
00360 ex_AH.set(l,k,j,0) = ex_fcn.val_point(h(l,k,j,0),(+theta)(l,k,j,0)
00361 ,(+phi)(l,k,j,0)) ;
00362 }
00363 }
00364 }
00365
00366 if (verbose) {
00367 cout << " " << endl ;
00368 cout << "###############################################" << endl ;
00369 cout << "AH finder: maximum number of iteration reached!" << endl ;
00370 cout << " No convergence in the 2-surface h! " << endl ;
00371 cout << " max( difference in h ) > prescribed tolerance " << endl ;
00372 cout << " " << endl ;
00373 cout << " prescribed tolerance = " << precis << endl ;
00374 cout << " max( difference in h ) = " << max(diff_h) << endl ;
00375 cout << " max( expansion function on h ) = " << max(abs(ex_AH(0))) << endl ;
00376 cout << "###############################################" << endl ;
00377 cout << " " << endl ;
00378
00379 }
00380 }
00381
00382
00383 }
00384
00385
00386
00387
00388
00389
00390
00391 if (no_AH_in_grid) {
00392 if (print) {
00393 cout << " " << endl ;
00394 cout << "###############################################" << endl ;
00395 cout << " AH finder: no horizon found inside the grid!" << endl ;
00396 cout << " Grid: rmin= " << r_min << ", rmax= " << r_max << endl ;
00397 cout << "###############################################" << endl ;
00398 cout << " " << endl ;
00399 }
00400 }
00401 else {
00402 for (int l=0; l<nz; l++) {
00403
00404 int jmax = mg->get_nt(l) ;
00405 int kmax = mg->get_np(l) ;
00406
00407 for (int k=0; k<kmax; k++) {
00408 for (int j=0; j<jmax; j++) {
00409
00410 ex_AH.set(l,k,j,0) = ex_fcn.val_point(h(l,k,j,0),(+theta)(l,k,j,0)
00411 ,(+phi)(l,k,j,0)) ;
00412 }
00413 }
00414 }
00415
00416
00417
00418 if ( (max(diff_h) < precis) && (max(abs(ex_AH(0))) < precis_exp) ) {
00419
00420 ah_flag = true ;
00421
00422 if (verbose) {
00423 cout << " " << endl ;
00424 cout << "################################################" << endl ;
00425 cout << " AH finder: Apparent horizon found!!! " << endl ;
00426 cout << " Max error of the expansion function on h: " << endl ;
00427 cout << " max( expansion function on AH ) = " << max(abs(ex_AH(0))) << endl ;
00428 cout << "################################################" << endl ;
00429 cout << " " << endl ;
00430 }
00431
00432 }
00433
00434 if ( (max(diff_h) < precis) && (max(abs(ex_AH(0))) > precis_exp) ) {
00435
00436 if (print) {
00437 cout << " " << endl ;
00438 cout << "#############################################" << endl ;
00439 cout << " AH finder: convergence in the 2 surface h. " << endl ;
00440 cout << " But max error of the expansion function evaulated on h > precis_exp" << endl ;
00441 cout << " max( expansion function on AH ) = " << max(abs(ex_AH(0))) << endl ;
00442 cout << " Probably not an apparent horizon! " << endl ;
00443 cout << "#############################################" << endl ;
00444 cout << " " << endl ;
00445 }
00446
00447 }
00448 }
00449 return ah_flag ;
00450
00451 }
00452