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