LORENE
app_hor_finder.C
1/*
2 * Function ah_finder
3 *
4 * (see file app_hor.h for documentation)
5 */
6
7/*
8 * Copyright (c) 2005 Lap-Ming Lin & Jerome Novak
9 *
10 *
11 * This file is part of LORENE.
12 *
13 * LORENE is free software; you can redistribute it and/or modify
14 * it under the terms of the GNU General Public License version 2
15 * as published by the Free Software Foundation.
16 *
17 * LORENE is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
20 * GNU General Public License for more details.
21 *
22 * You should have received a copy of the GNU General Public License
23 * along with LORENE; if not, write to the Free Software
24 * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
25 *
26 */
27
28char app_hor_finder_C[] = "$Header: /cvsroot/Lorene/C++/Source/App_hor/app_hor_finder.C,v 1.11 2014/10/13 08:52:38 j_novak Exp $" ;
29
30/*
31 * $Id: app_hor_finder.C,v 1.11 2014/10/13 08:52:38 j_novak Exp $
32 * $Log: app_hor_finder.C,v $
33 * Revision 1.11 2014/10/13 08:52:38 j_novak
34 * Lorene classes and functions now belong to the namespace Lorene.
35 *
36 * Revision 1.10 2014/10/06 15:12:56 j_novak
37 * Modified #include directives to use c++ syntax.
38 *
39 * Revision 1.9 2012/01/02 13:52:57 j_novak
40 * New parameter 'verbose' to get less output if needed.
41 *
42 * Revision 1.8 2008/01/08 13:56:54 j_novak
43 * Minor modif.
44 *
45 * Revision 1.7 2007/10/23 12:26:08 j_novak
46 * Added a test for the case where there is no AH, h(theta,phi) is then going out of the grid
47 *
48 * Revision 1.6 2005/12/09 09:35:06 lm_lin
49 *
50 * Add more information to screen output if no convergence.
51 *
52 * Revision 1.5 2005/12/07 14:16:36 lm_lin
53 *
54 * Add option to turn off screen output if no horizon is found
55 * (for performance reason in hydrodynamics simulation).
56 *
57 * Revision 1.4 2005/12/07 11:11:09 lm_lin
58 *
59 * Add option to turn off screen output during iterations.
60 *
61 * Revision 1.3 2005/11/17 15:53:28 lm_lin
62 *
63 * A tiny fix.
64 *
65 * Revision 1.2 2005/11/17 14:20:43 lm_lin
66 *
67 * Check the expansion function evaluated on the apparent horizon after the
68 * iteration of the 2-surface converges.
69 *
70 * Revision 1.1 2005/10/13 08:51:15 j_novak
71 * New stuff for apparent horizon finder. For the moment, there is only an
72 * external function. A class should come soon...
73 *
74 *
75 *
76 * $Header: /cvsroot/Lorene/C++/Source/App_hor/app_hor_finder.C,v 1.11 2014/10/13 08:52:38 j_novak Exp $
77 *
78 */
79
80
81// C headers
82#include <cmath>
83#include <cassert>
84
85// Lorene headers
86#include "app_hor.h"
87#include "graphique.h"
88
89
90namespace Lorene {
91bool ah_finder(const Metric& gamma, const Sym_tensor& k_dd, Valeur& h, Scalar& ex_fcn,
92 double a_axis, double b_axis, double c_axis, bool verbose, bool print,
93 double precis, double precis_exp, int step_max, int step_relax,
94 double relax)
95{
96
97 bool ah_flag = false ;
98
99 //Get the mapping, grid, base vector...etc
100 const Map& map = gamma.get_mp() ;
101 const Mg3d* mg = map.get_mg() ;
102 const Mg3d* g_angu = mg->get_angu() ;
103 int nz = mg->get_nzone() ;
104
105 const Base_vect_spher& bspher = map.get_bvect_spher() ;
106
107 const Coord& rr = map.r ;
108 const Coord& theta = map.tet ;
109 const Coord& phi = map.phi ;
110 const Coord& costh = map.cost ;
111 const Coord& cosph = map.cosp ;
112 const Coord& sinth = map.sint ;
113 const Coord& sinph = map.sinp ;
114
115 double r_min = min(+rr)(0) ;
116 double r_max = max(+rr)(nz-1) ;
117
118 // Set up a triaxial ellipsoidal surface as the initial guess for h
119 //------------------------------------------------------------------
120
121 double aa = a_axis ;
122 double bb = b_axis ;
123 double cc = c_axis ;
124
125 Scalar ct(map) ;
126 ct = costh ;
127 ct.std_spectral_base() ;
128
129 Scalar st(map) ;
130 st = sinth ;
131 st.std_spectral_base() ;
132
133 Scalar cp(map) ;
134 cp = cosph ;
135 cp.std_spectral_base() ;
136
137 Scalar sp(map) ;
138 sp = sinph ;
139 sp.std_spectral_base() ;
140
141 Scalar h_tmp(map) ;
142
143 h_tmp = st*st*( cp*cp/aa/aa + sp*sp/bb/bb ) + ct*ct/cc/cc ;
144 h_tmp = 1./sqrt( h_tmp ) ;
145 h_tmp.std_spectral_base() ;
146
147 Valeur h_guess(g_angu) ;
148 h_guess.annule_hard() ;
149
150 for (int l=0; l<nz; l++) {
151
152 int jmax = mg->get_nt(l) ;
153 int kmax = mg->get_np(l) ;
154
155 for (int k=0; k<kmax; k++) {
156 for (int j=0; j<jmax; j++) {
157 h_guess.set(l,k,j,0) = h_tmp.val_grid_point(l,k,j,0) ;
158 }
159 }
160 }
161
162 h_guess.std_base_scal() ;
163
164 h = h_guess ; // initialize h to h_guess
165 h.std_base_scal() ;
166
167 //End setting initial guess for h
168 //------------------------------------
169
170 const Metric_flat& fmets = map.flat_met_spher() ;
171
172
173 // Define the conformal factor
174 double one_third = double(1) / double(3) ;
175
176 Scalar psi4 = gamma.determinant() / fmets.determinant() ;
177 psi4 = pow( psi4, one_third ) ;
178 psi4.std_spectral_base() ;
179
180 // The expansion function at the n-1 iteration step
181 Scalar ex_fcn_old(map) ;
182 ex_fcn_old.set_etat_zero() ;
183
184 // The expansion function evaulated on the 2 surface h
185 Valeur ex_AH(g_angu) ;
186 ex_AH.annule_hard() ;
187
188 // Normal unit vector of the level set surface F = r - h(\theta,phi)
189 Vector s_unit(map, CON, bspher) ;
190
191 double relax_prev = double(1) - relax ;
192 double diff_exfcn = 1. ;
193 Tbl diff_h(nz) ;
194 diff_h = 1. ;
195
196 bool no_AH_in_grid = false ;
197
198 //--------------------------------------------------------
199 // Start of iteration
200 //--------------------------------------------------------
201
202 for (int step=0 ;
203 (max(diff_h) > precis) && (step < step_max) && (!no_AH_in_grid);
204 step++) {
205
206
207 // ***To be fixed: the function "set_grid_point" does not delete the derived
208 // quantities of F.
209 // Temporary fix: Define the level set F inside the iteration loop...
210 //----------------------------------------------------------------
211
212 Scalar F(map) ; // level set function: F = r - h(theta,phi)
213 F.allocate_all() ;
214
215 for (int l=0; l<nz; l++) {
216
217 int imax = mg->get_nr(l) ;
218 int jmax = mg->get_nt(l) ;
219 int kmax = mg->get_np(l) ;
220
221 for (int k=0; k<kmax; k++) {
222 for (int j=0; j<jmax; j++) {
223 for (int i=0; i<imax; i++) {
224
225 // (+rr) converts rr to Mtbl
226 F.set_grid_point(l,k,j,i) = (+rr)(l,k,j,i) - h(l,k,j,0) ;
227
228 }
229 }
230 }
231 }
232
234
235 // Construct the unit normal vector s^i of the surface F
236 Scalar dF_norm(map) ;
237 dF_norm = contract( contract(gamma.con(), 0, F.derive_cov(gamma), 0),
238 0, F.derive_cov(gamma), 0) ;
239 dF_norm = sqrt( dF_norm ) ;
240 dF_norm.std_spectral_base() ;
241
242 s_unit = F.derive_con(gamma) / dF_norm ;
243
244 // The expansion function
245 ex_fcn = s_unit.divergence(gamma) - k_dd.trace(gamma) +
246 contract( s_unit, 0, contract(s_unit, 0, k_dd, 1), 0) ;
247
248 // Construct the source term for the angular Laplace equation
249 //---------------------------------------------------------
250
251 Sym_tensor sou_1(map, CON, bspher) ;
252 sou_1 = gamma.con() - fmets.con()/psi4 - s_unit*s_unit ;
253
254 Scalar source_tmp(map) ;
255 source_tmp = contract( sou_1, 0, 1, F.derive_cov(fmets).derive_cov(fmets), 0, 1 ) ;
256 source_tmp = source_tmp / dF_norm ;
257
258 Sym_tensor d_gam(map, COV, bspher) ;
259 d_gam = contract( gamma.connect().get_delta(), 0, F.derive_cov(fmets), 0) ;
260
261 source_tmp -= contract( gamma.con() - s_unit*s_unit, 0, 1,
262 d_gam, 0, 1) / dF_norm ;
263
264 source_tmp = psi4*dF_norm*( source_tmp - k_dd.trace(gamma) +
265 contract(s_unit, 0, contract(s_unit, 0, k_dd, 1), 0) ) ;
266
267 source_tmp.std_spectral_base() ;
268
269
270 Valeur sou_angu(g_angu) ; // source defined on the angular grid
271 // S(theta, phi) = S(h(theta,phi),theta,phi)
272 sou_angu.annule_hard() ;
273
274 double h_min = min(h)(0) ;
275 double h_max = max(h)(0) ;
276 if ( (r_min < h_min) && (h_max < r_max) ) {
277
278 for (int l=0; l<nz; l++) {
279
280 int jmax = mg->get_nt(l) ;
281 int kmax = mg->get_np(l) ;
282 for (int k=0; k<kmax; k++) {
283 for (int j=0; j<jmax; j++) {
284 sou_angu.set(l,k,j,0) = source_tmp.val_point(h(l,k,j,0)
285 ,(+theta)(l,k,j,0) ,(+phi)(l,k,j,0)) ;
286 }
287 }
288 }
289 sou_angu = h*h*sou_angu ; // Final source term: psi4*dF_norm*h^2*(source_tmp)
290 }
291 else {
292 no_AH_in_grid = true ;
293 break ;
294 }
295 sou_angu.std_base_scal() ;
296
297
298 // Done with the source term
299 //-----------------------------------------
300
301
302 // Start solving the equation L^2h - 2h = source
303 //-----------------------------------------------
304
305 sou_angu.ylm() ;
306
307 Valeur h_new = sou_angu ;
308
309 h_new.c_cf->poisson_angu(-2.) ;
310
311 h_new.ylm_i() ;
312
313 if (h_new.c != 0x0)
314 delete h_new.c ;
315 h_new.c = 0x0 ;
316 h_new.coef_i() ;
317
318 // Convergence condition:
319 diff_h = max(abs(h - h_new)) ;
320
321
322 // Relaxations
323 if (step >= step_relax) {
324 h_new = relax * h_new + relax_prev * h ;
325 }
326
327 // Recycling for the next step
328 h = h_new ;
329
330
331 if (print)
332 {
333
334 cout << "-------------------------------------" << endl ;
335 cout << "App_hor iteration step: " << step << endl ;
336 cout << " " << endl ;
337
338 cout << "Difference in h : " << diff_h << endl ;
339
340 // Check: calculate the difference between ex_fcn and ex_fcn_old
341 Tbl diff_exfcn_tbl = diffrel( ex_fcn, ex_fcn_old ) ;
342 diff_exfcn = diff_exfcn_tbl(0) ;
343 for (int l=1; l<nz; l++) {
344 diff_exfcn += diff_exfcn_tbl(l) ;
345 }
346 diff_exfcn /= nz ;
347 cout << "diff_exfcn : " << diff_exfcn << endl ;
348
349 ex_fcn_old = ex_fcn ; // recycling
350 // End check
351
352 }
353
354 if ( (step == step_max-1) && (max(diff_h) > precis) ) {
355
356
357 //Check: Evaluate the expansion function on the 2-surface
358
359 for (int l=0; l<nz; l++) {
360
361 int jmax = mg->get_nt(l) ;
362 int kmax = mg->get_np(l) ;
363
364 for (int k=0; k<kmax; k++) {
365 for (int j=0; j<jmax; j++) {
366
367 ex_AH.set(l,k,j,0) = ex_fcn.val_point(h(l,k,j,0),(+theta)(l,k,j,0)
368 ,(+phi)(l,k,j,0)) ;
369 }
370 }
371 }
372
373 if (verbose) {
374 cout << " " << endl ;
375 cout << "###############################################" << endl ;
376 cout << "AH finder: maximum number of iteration reached!" << endl ;
377 cout << " No convergence in the 2-surface h! " << endl ;
378 cout << " max( difference in h ) > prescribed tolerance " << endl ;
379 cout << " " << endl ;
380 cout << " prescribed tolerance = " << precis << endl ;
381 cout << " max( difference in h ) = " << max(diff_h) << endl ;
382 cout << " max( expansion function on h ) = " << max(abs(ex_AH(0))) << endl ;
383 cout << "###############################################" << endl ;
384 cout << " " << endl ;
385
386 }
387 }
388
389
390 } // End of iteration
391
392 //Done with the AH finder
393
394
395
396 //Check: Evaluate the expansion function on the 2-surface
397
398 if (no_AH_in_grid) {
399 if (print) {
400 cout << " " << endl ;
401 cout << "###############################################" << endl ;
402 cout << " AH finder: no horizon found inside the grid!" << endl ;
403 cout << " Grid: rmin= " << r_min << ", rmax= " << r_max << endl ;
404 cout << "###############################################" << endl ;
405 cout << " " << endl ;
406 }
407 }
408 else {
409 for (int l=0; l<nz; l++) {
410
411 int jmax = mg->get_nt(l) ;
412 int kmax = mg->get_np(l) ;
413
414 for (int k=0; k<kmax; k++) {
415 for (int j=0; j<jmax; j++) {
416
417 ex_AH.set(l,k,j,0) = ex_fcn.val_point(h(l,k,j,0),(+theta)(l,k,j,0)
418 ,(+phi)(l,k,j,0)) ;
419 }
420 }
421 }
422
423
424
425 if ( (max(diff_h) < precis) && (max(abs(ex_AH(0))) < precis_exp) ) {
426
427 ah_flag = true ;
428
429 if (verbose) {
430 cout << " " << endl ;
431 cout << "################################################" << endl ;
432 cout << " AH finder: Apparent horizon found!!! " << endl ;
433 cout << " Max error of the expansion function on h: " << endl ;
434 cout << " max( expansion function on AH ) = " << max(abs(ex_AH(0))) << endl ;
435 cout << "################################################" << endl ;
436 cout << " " << endl ;
437 }
438
439 }
440
441 if ( (max(diff_h) < precis) && (max(abs(ex_AH(0))) > precis_exp) ) {
442
443 if (print) {
444 cout << " " << endl ;
445 cout << "#############################################" << endl ;
446 cout << " AH finder: convergence in the 2 surface h. " << endl ;
447 cout << " But max error of the expansion function evaulated on h > precis_exp" << endl ;
448 cout << " max( expansion function on AH ) = " << max(abs(ex_AH(0))) << endl ;
449 cout << " Probably not an apparent horizon! " << endl ;
450 cout << "#############################################" << endl ;
451 cout << " " << endl ;
452 }
453
454 }
455 }
456 return ah_flag ;
457
458} // End ah_finder
459
460}
Spherical orthonormal vectorial bases (triads).
Definition base_vect.h:308
const Tensor_sym & get_delta() const
Returns the tensor which defines the connection with respect to the flat one: is the difference bet...
Definition connection.h:271
Active physical coordinates and mapping derivatives.
Definition coord.h:90
Base class for coordinate mappings.
Definition map.h:670
Coord cosp
Definition map.h:724
Coord sint
Definition map.h:721
Coord r
r coordinate centered on the grid
Definition map.h:718
Coord tet
coordinate centered on the grid
Definition map.h:719
Coord sinp
Definition map.h:723
const Base_vect_spher & get_bvect_spher() const
Returns the orthonormal vectorial basis associated with the coordinates of the mapping.
Definition map.h:783
Coord phi
coordinate centered on the grid
Definition map.h:720
Coord cost
Definition map.h:722
const Mg3d * get_mg() const
Gives the Mg3d on which the mapping is defined.
Definition map.h:765
const Metric_flat & flat_met_spher() const
Returns the flat metric associated with the spherical coordinates and with components expressed in th...
Definition map.C:321
Flat metric for tensor calculation.
Definition metric.h:261
virtual const Sym_tensor & con() const
Read-only access to the contravariant representation.
virtual const Scalar & determinant() const
Returns the determinant.
Metric for tensor calculation.
Definition metric.h:90
virtual const Sym_tensor & con() const
Read-only access to the contravariant representation.
Definition metric.C:290
const Map & get_mp() const
Returns the mapping.
Definition metric.h:202
virtual const Connection & connect() const
Returns the connection.
Definition metric.C:301
virtual const Scalar & determinant() const
Returns the determinant.
Definition metric.C:392
Multi-domain grid.
Definition grilles.h:273
const Mg3d * get_angu() const
Returns the pointer on the associated angular grid.
Definition mg3d.C:473
int get_np(int l) const
Returns the number of points in the azimuthal direction ( ) in domain no. l.
Definition grilles.h:462
int get_nt(int l) const
Returns the number of points in the co-latitude direction ( ) in domain no. l.
Definition grilles.h:457
int get_nzone() const
Returns the number of domains.
Definition grilles.h:448
int get_nr(int l) const
Returns the number of points in the radial direction ( ) in domain no. l.
Definition grilles.h:452
void poisson_angu(double lambda=0)
Resolution of the generalized angular Poisson equation.
Definition mtbl_cf_pde.C:83
Tensor field of valence 0 (or component of a tensorial field).
Definition scalar.h:387
const Vector & derive_cov(const Metric &gam) const
Returns the gradient (1-form = covariant vector) of *this
virtual void std_spectral_base()
Sets the spectral bases of the Valeur va to the standard ones for a scalar field.
Definition scalar.C:784
double val_grid_point(int l, int k, int j, int i) const
Returns the value of the field at a specified grid point.
Definition scalar.h:637
virtual void set_etat_zero()
Sets the logical state to ETATZERO (zero).
Definition scalar.C:324
virtual void allocate_all()
Sets the logical state to ETATQCQ (ordinary state) and performs the memory allocation of all the elem...
Definition scalar.C:365
double & set_grid_point(int l, int k, int j, int i)
Setting the value of the field at a given grid point.
Definition scalar.h:684
double val_point(double r, double theta, double phi) const
Computes the value of the field at an arbitrary point , by means of the spectral expansion.
Definition scalar.C:890
const Vector & derive_con(const Metric &gam) const
Returns the "contravariant" derivative of *this with respect to some metric , by raising the index of...
Class intended to describe valence-2 symmetric tensors.
Definition sym_tensor.h:223
Basic array class.
Definition tbl.h:161
Values and coefficients of a (real-value) function.
Definition valeur.h:287
void ylm()
Computes the coefficients of *this.
Definition valeur_ylm.C:138
Mtbl * c
Values of the function at the points of the multi-grid
Definition valeur.h:299
Tbl & set(int l)
Read/write of the value in a given domain (configuration space).
Definition valeur.h:363
void coef_i() const
Computes the physical value of *this.
Mtbl_cf * c_cf
Coefficients of the spectral expansion of the function.
Definition valeur.h:302
void ylm_i()
Inverse of ylm()
void std_base_scal()
Sets the bases for spectral expansions (member base ) to the standard ones for a scalar.
Definition valeur.C:824
void annule_hard()
Sets the Valeur to zero in a hard way.
Definition valeur.C:723
Tensor field of valence 1.
Definition vector.h:188
const Scalar & divergence(const Metric &) const
The divergence of this with respect to a Metric .
Definition vector.C:381
Cmp sqrt(const Cmp &)
Square root.
Definition cmp_math.C:220
Tbl diffrel(const Cmp &a, const Cmp &b)
Relative difference between two Cmp (norme version).
Definition cmp_math.C:504
Tbl min(const Cmp &)
Minimum values of a Cmp in each domain.
Definition cmp_math.C:458
Tbl max(const Cmp &)
Maximum values of a Cmp in each domain.
Definition cmp_math.C:435
Cmp pow(const Cmp &, int)
Power .
Definition cmp_math.C:348
Cmp abs(const Cmp &)
Absolute value.
Definition cmp_math.C:410
bool ah_finder(const Metric &gamma, const Sym_tensor &k_dd_in, Valeur &h, Scalar &ex_fcn, double a_axis, double b_axis, double c_axis, bool verbose=true, bool print=false, double precis=1.e-8, double precis_exp=1.e-6, int it_max=200, int it_relax=200, double relax_fac=1.)
Class for apparent horizon (under heavy development)
const Tensor & derive_cov(const Metric &gam) const
Returns the covariant derivative of this with respect to some metric .
Definition tensor.C:1002
Tensor trace(int ind1, int ind2) const
Trace on two different type indices.
Tenseur contract(const Tenseur &, int id1, int id2)
Self contraction of two indices of a Tenseur .
Lorene prototypes.
Definition app_hor.h:64