LORENE
vector_poisson.C
1 /*
2  * Methods for solving vector Poisson equation.
3  *
4  * (see file vector.h for documentation).
5  *
6  */
7 
8 /*
9  * Copyright (c) 2003 Eric Gourgoulhon & Jerome Novak
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: vector_poisson.C,v 1.26 2026/03/04 15:47:04 j_novak Exp $
32  * $Log: vector_poisson.C,v $
33  * Revision 1.26 2026/03/04 15:47:04 j_novak
34  * Added 'verbose' flag (true by default) to limit the output of poisson solvers.
35  *
36  * Revision 1.25 2016/12/05 16:18:18 j_novak
37  * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
38  *
39  * Revision 1.24 2014/10/13 08:53:45 j_novak
40  * Lorene classes and functions now belong to the namespace Lorene.
41  *
42  * Revision 1.23 2014/10/06 15:13:21 j_novak
43  * Modified #include directives to use c++ syntax.
44  *
45  * Revision 1.22 2005/02/14 13:01:50 j_novak
46  * p_eta and p_mu are members of the class Vector. Most of associated functions
47  * have been moved from the class Vector_divfree to the class Vector.
48  *
49  * Revision 1.21 2005/01/10 15:36:09 j_novak
50  * In method 5: added transformation back from the Ylm base.
51  *
52  * Revision 1.20 2004/12/23 10:23:06 j_novak
53  * Improved method 5 in the case when some components are zero.
54  * Changed Vector_divfree::poisson to deduce eta from the equation.
55  *
56  * Revision 1.19 2004/08/24 09:14:50 p_grandclement
57  * Addition of some new operators, like Poisson in 2d... It now requieres the
58  * GSL library to work.
59  *
60  * Also, the way a variable change is stored by a Param_elliptic is changed and
61  * no longer uses Change_var but rather 2 Scalars. The codes using that feature
62  * will requiere some modification. (It should concern only the ones about monopoles)
63  *
64  * Revision 1.18 2004/07/27 09:40:05 j_novak
65  * Yet another method for solving vector Poisson eq. with spherical coordinates.
66  *
67  * Revision 1.17 2004/06/14 15:15:58 j_novak
68  * New method (No.4) to solve the vector Poisson eq. The output is continuous
69  * for all (spherical) components.
70  *
71  * Revision 1.16 2004/05/26 07:30:59 j_novak
72  * Version 1.15 was not good.
73  *
74  * Revision 1.15 2004/05/25 15:15:47 f_limousin
75  * Function Vector::poisson with parameters now returns a Vector (the
76  * result) instead of void.
77  *
78  * Revision 1.12 2004/03/26 17:05:24 j_novak
79  * Added new method n.3 using Tenseur::poisson_vect_oohara
80  *
81  * Revision 1.11 2004/03/11 08:48:45 f_limousin
82  * Implement method Vector::poisson with parameters, only with method
83  * 2 yet.
84  *
85  * Revision 1.10 2004/03/10 16:38:38 e_gourgoulhon
86  * Modified the prototype of poisson with param. to let it
87  * agree with declaration in vector.h.
88  *
89  * Revision 1.9 2004/03/03 09:07:03 j_novak
90  * In Vector::poisson(double, int), the flat metric is taken from the mapping.
91  *
92  * Revision 1.8 2004/02/24 17:00:25 j_novak
93  * Added a forgotten term.
94  *
95  * Revision 1.7 2004/02/24 09:46:20 j_novak
96  * Correction to cope with SGI compiler's warnings.
97  *
98  * Revision 1.6 2004/02/22 15:47:46 j_novak
99  * Added 2 more methods to solve the vector poisson equation. Method 1 is not
100  * tested yet.
101  *
102  * Revision 1.5 2004/02/20 10:53:41 j_novak
103  * Minor modifs.
104  *
105  * Revision 1.4 2004/02/16 17:40:14 j_novak
106  * Added a version of poisson with the flat metric as argument (avoids
107  * unnecessary calculations by decompose_div)
108  *
109  * Revision 1.3 2003/10/29 11:04:34 e_gourgoulhon
110  * dec2_dpzuis() replaced by dec_dzpuis(2).
111  * inc2_dpzuis() replaced by inc_dzpuis(2).
112  *
113  * Revision 1.2 2003/10/22 13:08:06 j_novak
114  * Better handling of dzpuis flags
115  *
116  * Revision 1.1 2003/10/20 15:15:42 j_novak
117  * New method Vector::poisson().
118  *
119  *
120  * $Headers: $
121  *
122  */
123 
124 // Lorene headers
125 #include "metric.h"
126 #include "tenseur.h"
127 #include "param.h"
128 #include "param_elliptic.h"
129 #include "proto.h"
130 
131 namespace Lorene {
132 Vector Vector::poisson(double lambda, const Metric_flat& met_f, int method)
133  const {
134 
135  bool nullite = true ;
136  for (int i=0; i<3; i++) {
137  assert(cmp[i]->check_dzpuis(4)) ;
138  if (cmp[i]->get_etat() != ETATZERO) nullite = false ;
139  }
140  assert ((method>=0) && (method<7)) ;
141 
142  Vector resu(*mp, CON, triad) ;
143  if (nullite)
144  resu.set_etat_zero() ;
145  else {
146 
147  switch (method) {
148 
149  case 0 : {
150 
151  Scalar poten(*mp) ;
152  if (fabs(lambda+1) < 1.e-8)
153  poten.set_etat_zero() ;
154  else {
155  poten = (potential(met_f) / (lambda + 1)).poisson() ;
156  }
157 
158  Vector grad = poten.derive_con(met_f) ;
159  grad.dec_dzpuis(2) ;
160 
161  return ( div_free(met_f).poisson() + grad) ;
162  break ;
163  }
164 
165  case 1 : {
166 
167  Scalar divf(*mp) ;
168  if (fabs(lambda+1) < 1.e-8)
169  divf.set_etat_zero() ;
170  else {
171  divf = (potential(met_f) / (lambda + 1)) ;
172  }
173 
174  int nz = mp->get_mg()->get_nzone() ;
175 
176  //-----------------------------------
177  // Removal of the l=0 part of div(F)
178  //-----------------------------------
179  Scalar div_lnot0 = divf ;
180  div_lnot0.div_r_dzpuis(4) ;
181  Scalar source_r(*mp) ;
182  Valeur& va_div = div_lnot0.set_spectral_va() ;
183  if (div_lnot0.get_etat() != ETATZERO) {
184  va_div.coef() ;
185  va_div.ylm() ;
186  for (int lz=0; lz<nz; lz++) {
187  int np = mp->get_mg()->get_np(lz) ;
188  int nt = mp->get_mg()->get_nt(lz) ;
189  int nr = mp->get_mg()->get_nr(lz) ;
190  if (va_div.c_cf->operator()(lz).get_etat() != ETATZERO)
191  for (int k=0; k<np+1; k++)
192  for (int j=0; j<nt; j++) {
193  int l_q, m_q, base_r ;
194  if (nullite_plm(j, nt, k, np, va_div.base) == 1) {
195  donne_lm(nz, lz, j, k, va_div.base, m_q, l_q, base_r) ;
196  if (l_q == 0)
197  for (int i=0; i<nr; i++)
198  va_div.c_cf->set(lz, k, j, i) = 0. ;
199  }
200  }
201  }
202  source_r.set_etat_qcq() ;
203  source_r.set_spectral_va() = 2*(*va_div.c_cf) ; //2*div(F)
204  source_r.set_spectral_va().ylm_i() ;
205  source_r.set_dzpuis(4) ;
206  }
207  else
208  source_r.set_etat_zero() ;
209 
210  //------------------------
211  // Other source terms ....
212  //------------------------
213  source_r += *(cmp[0]) - lambda*divf.dsdr() ;
214  Scalar f_r(*mp) ;
215  if (source_r.get_etat() != ETATZERO) {
216 
217  //------------------------
218  // The elliptic operator
219  //------------------------
220  Param_elliptic param_fr(source_r) ;
221  for (int lz=0; lz<nz; lz++)
222  param_fr.set_poisson_vect_r(lz) ;
223 
224  f_r = source_r.sol_elliptic(param_fr) ;
225  }
226  else
227  f_r.set_etat_zero() ;
228 
229  divf.dec_dzpuis(1) ;
230  Scalar source_eta = divf - f_r.dsdr() ;
231  source_eta.mult_r_dzpuis(0) ;
232  source_eta -= 2*f_r ;
233  resu.set_vr_eta_mu(f_r, source_eta.poisson_angu(), mu().poisson());
234 
235  break ;
236 
237  }
238 
239  case 2 : {
240 
241  Tenseur source_p(*mp, 1, CON, mp->get_bvect_spher() ) ;
242  source_p.set_etat_qcq() ;
243  for (int i=0; i<3; i++) {
244  source_p.set(i) = Cmp(*cmp[i]) ;
245  }
246  source_p.change_triad(mp->get_bvect_cart()) ;
247  Tenseur vect_auxi (*mp, 1, CON, mp->get_bvect_cart()) ;
248  vect_auxi.set_etat_qcq() ;
249  Tenseur scal_auxi (*mp) ;
250  scal_auxi.set_etat_qcq() ;
251 
252  Tenseur resu_p(source_p.poisson_vect(lambda, vect_auxi, scal_auxi)) ;
253  resu_p.change_triad(mp->get_bvect_spher() ) ;
254 
255  for (int i=1; i<=3; i++)
256  resu.set(i) = resu_p(i-1) ;
257 
258  break ;
259  }
260 
261  case 3 : {
262 
263  Tenseur source_p(*mp, 1, CON, mp->get_bvect_spher() ) ;
264  source_p.set_etat_qcq() ;
265  for (int i=0; i<3; i++) {
266  source_p.set(i) = Cmp(*cmp[i]) ;
267  }
268  source_p.change_triad(mp->get_bvect_cart()) ;
269  Tenseur scal_auxi (*mp) ;
270  scal_auxi.set_etat_qcq() ;
271 
272  Tenseur resu_p(source_p.poisson_vect_oohara(lambda, scal_auxi)) ;
273  resu_p.change_triad(mp->get_bvect_spher() ) ;
274 
275  for (int i=1; i<=3; i++)
276  resu.set(i) = resu_p(i-1) ;
277 
278  break ;
279  }
280 
281  case 4 : {
282  Scalar divf(*mp) ;
283  if (fabs(lambda+1) < 1.e-8)
284  divf.set_etat_zero() ;
285  else {
286  divf = (potential(met_f) / (lambda + 1)) ;
287  }
288 
289  int nz = mp->get_mg()->get_nzone() ;
290 
291  //-----------------------------------
292  // Removal of the l=0 part of div(F)
293  //-----------------------------------
294  Scalar div_tmp = divf ;
295  div_tmp.div_r_dzpuis(4) ;
296  Scalar source_r(*mp) ;
297  Valeur& va_div = div_tmp.set_spectral_va() ;
298  if (div_tmp.get_etat() != ETATZERO) {
299  va_div.ylm() ;
300  for (int lz=0; lz<nz; lz++) {
301  int np = mp->get_mg()->get_np(lz) ;
302  int nt = mp->get_mg()->get_nt(lz) ;
303  int nr = mp->get_mg()->get_nr(lz) ;
304  if (va_div.c_cf->operator()(lz).get_etat() != ETATZERO)
305  for (int k=0; k<np+1; k++)
306  for (int j=0; j<nt; j++) {
307  int l_q, m_q, base_r ;
308  if (nullite_plm(j, nt, k, np, va_div.base) == 1) {
309  donne_lm(nz, lz, j, k, va_div.base, m_q, l_q, base_r) ;
310  if (l_q == 0)
311  for (int i=0; i<nr; i++)
312  va_div.c_cf->set(lz, k, j, i) = 0. ;
313  }
314  }
315  }
316  source_r.set_etat_qcq() ;
317  source_r.set_spectral_va() = 2*(*va_div.c_cf) ; //2*div(F)
318  source_r.set_spectral_va().ylm_i() ;
319  source_r.set_dzpuis(4) ;
320  }
321  else
322  source_r.set_etat_zero() ;
323 
324  //------------------------
325  // Other source terms ....
326  //------------------------
327  source_r += *(cmp[0]) - lambda*divf.dsdr() ;
328  Scalar f_r(*mp) ;
329  if (source_r.get_etat() != ETATZERO) {
330 
331  //------------------------
332  // The elliptic operator
333  //------------------------
334  Param_elliptic param_fr(source_r) ;
335  for (int lz=0; lz<nz; lz++)
336  param_fr.set_poisson_vect_r(lz) ;
337 
338  f_r = source_r.sol_elliptic(param_fr) ;
339  }
340  else
341  f_r.set_etat_zero() ;
342 
343  //--------------------------
344  // Equation for eta
345  //--------------------------
346  Scalar source_eta = *cmp[1] ;
347  source_eta.mult_sint() ;
348  source_eta = source_eta.dsdt() ;
349  source_eta.div_sint() ;
350  source_eta = (source_eta + cmp[2]->stdsdp()).poisson_angu() ;
351 
352  Scalar dfrsr = f_r.dsdr() ;
353  dfrsr.div_r_dzpuis(4) ;
354  Scalar frsr2 = f_r ;
355  frsr2.div_r_dzpuis(2) ;
356  frsr2.div_r_dzpuis(4) ;
357 
358  Valeur& va_dfr = dfrsr.set_spectral_va() ;
359  Valeur& va_fsr = frsr2.set_spectral_va() ;
360  va_dfr.ylm() ;
361  va_fsr.ylm() ;
362 
363  //Since the operator depends on the domain, the source
364  //must be transformed accordingly.
365  Valeur& va_eta = source_eta.set_spectral_va() ;
366  if (source_eta.get_etat() == ETATZERO) source_eta.annule_hard() ;
367  va_eta.ylm() ;
368  for (int lz=0; lz<nz; lz++) {
369  bool ced = (mp->get_mg()->get_type_r(lz) == UNSURR) ;
370  int np = mp->get_mg()->get_np(lz) ;
371  int nt = mp->get_mg()->get_nt(lz) ;
372  int nr = mp->get_mg()->get_nr(lz) ;
373  if (va_eta.c_cf->operator()(lz).get_etat() != ETATZERO)
374  for (int k=0; k<np+1; k++)
375  for (int j=0; j<nt; j++) {
376  int l_q, m_q, base_r ;
377  if (nullite_plm(j, nt, k, np, va_eta.base) == 1) {
378  donne_lm(nz, lz, j, k, va_eta.base, m_q, l_q, base_r) ;
379  if (l_q > 0)
380  for (int i=0; i<nr; i++) {
381  if (va_div.c_cf->operator()(lz).get_etat() != ETATZERO)
382  va_eta.c_cf->set(lz, k, j, i)
383  -= (lambda + 2. / double(ced ? -l_q : (l_q+1) ))
384  * va_div.c_cf->operator()(lz, k, j, i) ;
385  if (va_fsr.c_cf->operator()(lz).get_etat() != ETATZERO) {
386  va_eta.c_cf->set(lz, k, j, i)
387  += 2. / double(ced ? -l_q : (l_q+1) )
388  * va_dfr.c_cf->operator()(lz, k, j, i) ;
389  va_eta.c_cf->set(lz, k, j, i)
390  -= 2.*( ced ? double(l_q+2)/double(l_q)
391  : double(l_q-1)/double(l_q+1) )
392  * va_fsr.c_cf->set(lz, k, j, i) ;
393  }
394  } //Loop on r
395  } //nullite_plm
396  } //Loop on theta
397  } // Loop on domains
398  if (va_eta.c != 0x0) {
399  delete va_eta.c;
400  va_eta.c = 0x0 ;
401  }
402  va_eta.ylm_i() ;
403 
404  //------------------------
405  // The elliptic operator
406  //------------------------
407  Param_elliptic param_eta(source_eta) ;
408  for (int lz=0; lz<nz; lz++)
409  param_eta.set_poisson_vect_eta(lz) ;
410 
411  resu.set_vr_eta_mu(f_r, source_eta.sol_elliptic(param_eta), mu().poisson()) ;
412  break ;
413 
414  }
415 
416  case 5 : {
417 
418  Scalar divf(*mp) ;
419  if (fabs(lambda+1) < 1.e-8)
420  divf.set_etat_zero() ;
421  else {
422  divf = (potential(met_f) / (lambda + 1)) ;
423  }
424 
425  int nz = mp->get_mg()->get_nzone() ;
426 
427  //-----------------------------------
428  // Removal of the l=0 part of div(F)
429  //-----------------------------------
430  Scalar div_lnot0 = divf ;
431  div_lnot0.div_r_dzpuis(4) ;
432  Scalar source_r(*mp) ;
433  Valeur& va_div = div_lnot0.set_spectral_va() ;
434  if (div_lnot0.get_etat() != ETATZERO) {
435  va_div.coef() ;
436  va_div.ylm() ;
437  for (int lz=0; lz<nz; lz++) {
438  int np = mp->get_mg()->get_np(lz) ;
439  int nt = mp->get_mg()->get_nt(lz) ;
440  int nr = mp->get_mg()->get_nr(lz) ;
441  if (va_div.c_cf->operator()(lz).get_etat() != ETATZERO)
442  for (int k=0; k<np+1; k++)
443  for (int j=0; j<nt; j++) {
444  int l_q, m_q, base_r ;
445  if (nullite_plm(j, nt, k, np, va_div.base) == 1) {
446  donne_lm(nz, lz, j, k, va_div.base, m_q, l_q, base_r) ;
447  if (l_q == 0)
448  for (int i=0; i<nr; i++)
449  va_div.c_cf->set(lz, k, j, i) = 0. ;
450  }
451  }
452  }
453  source_r.set_etat_qcq() ;
454  source_r.set_spectral_va() = 2*(*va_div.c_cf) ; //2*div(F)
455  source_r.set_spectral_va().ylm_i() ;
456  source_r.set_dzpuis(4) ;
457  }
458  else
459  source_r.set_etat_zero() ;
460 
461  //------------------------
462  // Other source terms ....
463  //------------------------
464  source_r += *(cmp[0]) - lambda*divf.dsdr() ;
465  Scalar f_r(*mp) ;
466  if (source_r.get_etat() != ETATZERO) {
467 
468  //------------------------
469  // The elliptic operator
470  //------------------------
471 
472  Param_elliptic param_fr(source_r) ;
473  for (int lz=0; lz<nz; lz++)
474  param_fr.set_poisson_vect_r(lz) ;
475 
476  f_r = source_r.sol_elliptic(param_fr) ;
477  }
478  else
479  f_r.set_etat_zero() ;
480 
481  Scalar source_eta = - *(cmp[0]) ;
482  source_eta.mult_r_dzpuis(3) ;
483  source_eta -= (lambda+2.)*divf ;
484  source_eta.dec_dzpuis() ;
485  f_r.set_spectral_va().ylm() ;
486  Scalar tmp = 2*f_r + f_r.lapang() ;
487  tmp.div_r_dzpuis(2) ;
488  source_eta += tmp ;
489  tmp = (1.+lambda)*divf ;
490  tmp.mult_r_dzpuis(0) ;
491  tmp += f_r ;
492  source_eta = source_eta.primr() ;
493  f_r.set_spectral_va().ylm_i() ;
494  resu.set_vr_eta_mu(f_r, (tmp+source_eta).poisson_angu(), mu().poisson()) ;
495  break ;
496 
497  }
498 
499  case 6 : {
500 
501  poisson_block(lambda, resu) ;
502  break ;
503 
504  }
505  default : {
506  cerr << "Vector::poisson : unexpected type of method !" << endl
507  << " method = " << method << endl ;
508  abort() ;
509  break ;
510  }
511 
512  } // End of switch
513 
514  } // End of non-null case
515 
516  return resu ;
517 
518 }
519 
520 Vector Vector::poisson(double lambda, int method) const {
521 
522  const Base_vect_spher* tspher = dynamic_cast<const Base_vect_spher*>(triad) ;
523  const Base_vect_cart* tcart = dynamic_cast<const Base_vect_cart*>(triad) ;
524 
525  assert ((tspher != 0x0) || (tcart != 0x0)) ;
526  const Metric_flat* met_f = 0x0 ;
527 
528  if (tspher != 0x0) {
529  assert (tcart == 0x0) ;
530  met_f = &(mp->flat_met_spher()) ;
531  }
532 
533  if (tcart != 0x0) {
534  assert (tspher == 0x0) ;
535  met_f = &(mp->flat_met_cart()) ;
536  }
537 
538  return ( poisson(lambda, *met_f, method) );
539 
540 }
541 
542 // Version with parameters
543 // -----------------------
544 
545 Vector Vector::poisson(const double lambda, Param& par, int method) const {
546 
547 
548  for (int i=0; i<3; i++)
549  assert(cmp[i]->check_dzpuis(4)) ;
550 
551  assert ((method==0) || (method==2)) ;
552 
553  Vector resu(*mp, CON, triad) ;
554 
555  switch (method) {
556 
557  case 0 : {
558 
559  Metric_flat met_local(*mp, *triad) ;
560  int nitermax = par.get_int(0) ;
561  int& niter = par.get_int_mod(0) ;
562  double relax = par.get_double(0) ;
563  double precis = par.get_double(1) ;
564  Cmp& ss_phi = par.get_cmp_mod(0) ;
565  Cmp& ss_khi = par.get_cmp_mod(1) ;
566  Cmp& ss_mu = par.get_cmp_mod(2) ;
567 
568  Param par_phi ;
569  par_phi.add_int(nitermax, 0) ;
570  par_phi.add_int_mod(niter, 0) ;
571  par_phi.add_double(relax, 0) ;
572  par_phi.add_double(precis, 1) ;
573  par_phi.add_cmp_mod(ss_phi, 0) ;
574 
575  Scalar poten(*mp) ;
576  if (fabs(lambda+1) < 1.e-8)
577  poten.set_etat_zero() ;
578  else {
579  Scalar tmp = potential(met_local) / (lambda + 1) ;
580  tmp.inc_dzpuis(2) ;
581  tmp.poisson(par_phi, poten) ;
582  }
583 
584  Vector grad = poten.derive_con(met_local) ;
585  grad.dec_dzpuis(2) ;
586 
587  Param par_free ;
588  par_free.add_int(nitermax, 0) ;
589  par_free.add_int_mod(niter, 0) ;
590  par_free.add_double(relax, 0) ;
591  par_free.add_double(precis, 1) ;
592  par_free.add_cmp_mod(ss_khi, 0) ;
593  par_free.add_cmp_mod(ss_mu, 1) ;
594 
595  return div_free(met_local).poisson(par_free) + grad ;
596  break ;
597  }
598 
599 
600 
601  case 2 : {
602 
603  Tenseur source_p(*mp, 1, CON, mp->get_bvect_spher() ) ;
604  source_p.set_etat_qcq() ;
605  for (int i=0; i<3; i++) {
606  source_p.set(i) = Cmp(*cmp[i]) ;
607  }
608  source_p.change_triad(mp->get_bvect_cart()) ;
609 
610  Tenseur vect_auxi (*mp, 1, CON, mp->get_bvect_cart()) ;
611  vect_auxi.set_etat_qcq() ;
612  for (int i=0; i<3 ;i++){
613  vect_auxi.set(i) = 0. ;
614  }
615  Tenseur scal_auxi (*mp) ;
616  scal_auxi.set_etat_qcq() ;
617  scal_auxi.set().annule_hard() ;
618  scal_auxi.set_std_base() ;
619 
620  Tenseur resu_p(*mp, 1, CON, mp->get_bvect_cart() ) ;
621  resu_p.set_etat_qcq() ;
622  source_p.poisson_vect(lambda, par, resu_p, vect_auxi, scal_auxi) ;
623  resu_p.change_triad(mp->get_bvect_spher() ) ;
624 
625  for (int i=1; i<=3; i++)
626  resu.set(i) = resu_p(i-1) ;
627 
628  break ;
629  }
630 
631  default : {
632  cerr << "Vector::poisson : unexpected type of method !" << endl
633  << " method = " << method << endl ;
634 
635 
636  abort() ;
637  break ;
638  }
639 
640  }// End of switch
641  return resu ;
642 
643 }
644 
645 
646 
647 
648 }
void set_poisson_vect_eta(int zone)
Sets the operator to be a regular elliptic operator to solve for the component of the vector Poisson...
Mtbl_cf * c_cf
Coefficients of the spectral expansion of the function.
Definition: valeur.h:312
Component of a tensorial field *** DEPRECATED : use class Scalar instead ***.
Definition: cmp.h:449
void add_int(const int &n, int position=0)
Adds the address of a new int to the list.
Definition: param.C:249
void ylm_i()
Inverse of ylm()
Definition: valeur_ylm_i.C:134
const Scalar & lapang() const
Returns the angular Laplacian of *this , where .
Definition: scalar_deriv.C:460
int get_np(int l) const
Returns the number of points in the azimuthal direction ( ) in domain no. l.
Definition: grilles.h:481
void coef() const
Computes the coeffcients of *this.
Definition: valeur_coef.C:151
virtual void set_etat_zero()
Sets the logical state to ETATZERO (zero).
Definition: scalar.C:330
const Vector_divfree & div_free(const Metric &) const
Returns the div-free vector in the Helmholtz decomposition.
Definition: vector.C:510
void set_std_base()
Set the standard spectal basis of decomposition for each component.
Definition: tenseur.C:1186
const Base_vect_spher & get_bvect_spher() const
Returns the orthonormal vectorial basis associated with the coordinates of the mapping.
Definition: map.h:809
Lorene prototypes.
Definition: app_hor.h:67
Tbl & set(int l)
Read/write of the Tbl containing the coefficients in a given domain.
Definition: mtbl_cf.h:304
void ylm()
Computes the coefficients of *this.
Definition: valeur_ylm.C:141
const Scalar & dsdt() const
Returns of *this .
Definition: scalar_deriv.C:207
const Mg3d * get_mg() const
Gives the Mg3d on which the mapping is defined.
Definition: map.h:791
Flat metric for tensor calculation.
Definition: metric.h:261
Tensor field of valence 0 (or component of a tensorial field).
Definition: scalar.h:402
Scalar poisson() const
Solves the scalar Poisson equation with *this as a source.
Definition: scalar_pde.C:142
Scalar primr(bool null_infty=true) const
Computes the radial primitive which vanishes for .
Definition: scalar_integ.C:106
void mult_sint()
Multiplication by .
void annule_hard()
Sets the Cmp to zero in a hard way.
Definition: cmp.C:341
const Base_vect * triad
Vectorial basis (triad) with respect to which the tensor components are defined.
Definition: tensor.h:315
Values and coefficients of a (real-value) function.
Definition: valeur.h:297
int get_etat() const
Returns the logical state ETATNONDEF (undefined), ETATZERO (null) or ETATQCQ (ordinary).
Definition: scalar.h:569
virtual void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition: scalar.C:359
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:401
Tensor field of valence 1.
Definition: vector.h:188
Scalar poisson_angu(double lambda=0) const
Solves the (generalized) angular Poisson equation with *this as source.
Definition: scalar_pde.C:192
const Metric_flat & flat_met_cart() const
Returns the flat metric associated with the Cartesian coordinates and with components expressed in th...
Definition: map.C:337
void set_dzpuis(int)
Modifies the dzpuis flag.
Definition: scalar.C:814
virtual void dec_dzpuis(int dec=1)
Decreases by dec units the value of dzpuis and changes accordingly the values in the compactified ext...
Definition: tensor.C:818
void annule_hard()
Sets the Scalar to zero in a hard way.
Definition: scalar.C:386
Cmp & get_cmp_mod(int position=0) const
Returns the reference of a modifiable Cmp stored in the list.
Definition: param.C:1052
Cmp & set()
Read/write for a scalar (see also operator=(const Cmp&) ).
Definition: tenseur.C:840
const int & get_int(int position=0) const
Returns the reference of a int stored in the list.
Definition: param.C:295
void change_triad(const Base_vect &new_triad)
Sets a new vectorial basis (triad) of decomposition and modifies the components accordingly.
Definition: tenseur.C:684
virtual void inc_dzpuis(int inc=1)
Increases by inc units the value of dzpuis and changes accordingly the values of the Scalar in the co...
Vector poisson(double lambda, int method=6) const
Solves the vector Poisson equation with *this as a source.
Base_val base
Bases on which the spectral expansion is performed.
Definition: valeur.h:315
void set_vr_eta_mu(const Scalar &vr_i, const Scalar &eta_i, const Scalar &mu_i)
Defines the components through potentials and (see members p_eta and p_mu ), as well as the compon...
Definition: vector_etamu.C:198
Scalar ** cmp
Array of size n_comp of pointers onto the components.
Definition: tensor.h:327
Scalar sol_elliptic(Param_elliptic &params) const
Resolution of a general elliptic equation, putting zero at infinity.
Definition: scalar_pde.C:226
Mtbl * c
Values of the function at the points of the multi-grid.
Definition: valeur.h:309
Parameter storage.
Definition: param.h:125
int get_nzone() const
Returns the number of domains.
Definition: grilles.h:467
int & get_int_mod(int position=0) const
Returns the reference of a modifiable int stored in the list.
Definition: param.C:433
const Scalar & stdsdp() const
Returns of *this .
Definition: scalar_deriv.C:237
Vector_divfree poisson() const
Computes the solution of a vectorial Poisson equation with *this as a source: .
This class contains the parameters needed to call the general elliptic solver.
Cartesian vectorial bases (triads).
Definition: base_vect.h:201
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:471
const Base_vect_cart & get_bvect_cart() const
Returns the Cartesian basis associated with the coordinates (x,y,z) of the mapping, i.e.
Definition: map.h:817
const Scalar & potential(const Metric &) const
Returns the potential in the Helmholtz decomposition.
Definition: vector.C:498
void add_double(const double &x, int position=0)
Adds the the address of a new double to the list.
Definition: param.C:318
const Scalar & dsdr() const
Returns of *this .
Definition: scalar_deriv.C:112
virtual const Scalar & mu() const
Gives the field such that the angular components of the vector are written: .
Definition: vector_etamu.C:107
void div_sint()
Division by .
virtual void set_etat_zero()
Sets the logical state of all components to ETATZERO (zero state).
Definition: tensor.C:507
const double & get_double(int position=0) const
Returns the reference of a double stored in the list.
Definition: param.C:364
int get_nt(int l) const
Returns the number of points in the co-latitude direction ( ) in domain no. l.
Definition: grilles.h:476
Valeur & set_spectral_va()
Returns va (read/write version)
Definition: scalar.h:619
Scalar & set(int)
Read/write access to a component.
Definition: vector.C:302
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition: tenseur.C:652
void add_cmp_mod(Cmp &ti, int position=0)
Adds the address of a new modifiable Cmp to the list.
Definition: param.C:1007
int get_type_r(int l) const
Returns the type of sampling in the radial direction in domain no.
Definition: grilles.h:493
const Map *const mp
Mapping on which the numerical values at the grid points are defined.
Definition: tensor.h:307
void div_r_dzpuis(int ced_mult_r)
Division by r everywhere but with the output flag dzpuis set to ced_mult_r .
virtual void dec_dzpuis(int dec=1)
Decreases by dec units the value of dzpuis and changes accordingly the values of the Scalar in the co...
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:327
void mult_r_dzpuis(int ced_mult_r)
Multiplication by r everywhere but with the output flag dzpuis set to ced_mult_r ...
Tensor handling *** DEPRECATED : use class Tensor instead ***.
Definition: tenseur.h:307
void add_int_mod(int &n, int position=0)
Adds the address of a new modifiable int to the list.
Definition: param.C:388
void set_poisson_vect_r(int zone, bool only_l_zero=false)
Sets the operator to in all domains, for ; and to in all domains otherwise.