LORENE
tensor.C
1 /*
2  * Methods of class Tensor
3  *
4  * (see file tensor.h for documentation)
5  *
6  */
7 
8 /*
9  * Copyright (c) 2003-2004 Eric Gourgoulhon & Jerome Novak
10  *
11  * Copyright (c) 1999-2001 Philippe Grandclement (for preceding class Tenseur)
12  *
13  * This file is part of LORENE.
14  *
15  * LORENE is free software; you can redistribute it and/or modify
16  * it under the terms of the GNU General Public License as published by
17  * the Free Software Foundation; either version 2 of the License, or
18  * (at your option) any later version.
19  *
20  * LORENE is distributed in the hope that it will be useful,
21  * but WITHOUT ANY WARRANTY; without even the implied warranty of
22  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
23  * GNU General Public License for more details.
24  *
25  * You should have received a copy of the GNU General Public License
26  * along with LORENE; if not, write to the Free Software
27  * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
28  *
29  */
30 
31 
32 
33 
34 /*
35  * $Id: tensor.C,v 1.48 2025/03/31 08:29:50 j_novak Exp $
36  * $Log: tensor.C,v $
37  * Revision 1.48 2025/03/31 08:29:50 j_novak
38  * Methods to copyt Tensors and Scalars to smaller grids...
39  *
40  * Revision 1.47 2025/03/28 09:12:50 j_novak
41  * News methods to copy Tensors and Scalars to larger grids. To be used for de-aliasing.
42  *
43  * Revision 1.46 2023/08/31 08:27:26 g_servignat
44  * Added the possibility to filter in the r direction within the ylm filter. An order filtering of 0 means no filtering (for all 3 directions).
45  *
46  * Revision 1.45 2023/08/28 09:53:33 g_servignat
47  * Added ylm filter for Tensor and Scalar in theta + phi directions
48  *
49  * Revision 1.44 2016/12/05 16:18:17 j_novak
50  * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
51  *
52  * Revision 1.43 2014/10/13 08:53:44 j_novak
53  * Lorene classes and functions now belong to the namespace Lorene.
54  *
55  * Revision 1.42 2014/10/06 15:13:20 j_novak
56  * Modified #include directives to use c++ syntax.
57  *
58  * Revision 1.41 2013/06/05 15:43:49 j_novak
59  * Suppression of leg_spectral_base()
60  *
61  * Revision 1.40 2013/01/11 15:44:54 j_novak
62  * Addition of Legendre bases (part 2).
63  *
64  * Revision 1.39 2007/12/21 16:07:08 j_novak
65  * Methods to filter Tensor, Vector and Sym_tensor objects.
66  *
67  * Revision 1.38 2005/10/25 08:56:38 p_grandclement
68  * addition of std_spectral_base in the case of odd functions near the origin
69  *
70  * Revision 1.37 2005/05/18 11:45:44 j_novak
71  * Added del_deriv() calls at the end of inc/dec_dzpuis.
72  *
73  * Revision 1.36 2004/07/08 12:21:53 j_novak
74  * Replaced tensor::annule_extern_c2 with tensor::annule_extern_cn for a
75  * more general transition.
76  *
77  * Revision 1.35 2004/06/17 06:55:07 e_gourgoulhon
78  * Added method annule_extern_c2.
79  *
80  * Revision 1.34 2004/03/04 09:47:51 e_gourgoulhon
81  * Method annule_domain(int, int): added call to virtual function
82  * del_deriv() at the end !
83  *
84  * Revision 1.33 2004/02/27 21:14:27 e_gourgoulhon
85  * Modif of method derive_con to create proper type of the result.
86  *
87  * Revision 1.32 2004/02/19 22:10:14 e_gourgoulhon
88  * Added argument "comment" in method spectral_display.
89  *
90  * Revision 1.31 2004/02/05 15:03:47 e_gourgoulhon
91  * Corrected bug in method derive_con().
92  *
93  * Revision 1.30 2004/01/27 13:05:11 j_novak
94  * Removed the method Tensor::mult_r_ced()
95  *
96  * Revision 1.29 2004/01/19 16:32:13 e_gourgoulhon
97  * Added operator()(int, int, int, int) and set(int, int, int, int)
98  * for direct access to components of valence 4 tensors.
99  *
100  * Revision 1.28 2004/01/04 20:55:23 e_gourgoulhon
101  * Method spectral_display(): added printing of type of class (through typeid).
102  *
103  * Revision 1.27 2003/12/30 23:09:47 e_gourgoulhon
104  * Change in methods derive_cov() and divergence() to take into account
105  * the change of name: Metric::get_connect() --> Metric::connect().
106  *
107  * Revision 1.26 2003/12/27 15:01:50 e_gourgoulhon
108  * Method derive_con(): the position of the derivation index has
109  * been changed from the first one to the last one.
110  *
111  * Revision 1.25 2003/11/03 22:34:41 e_gourgoulhon
112  * Method dec_dzpuis: changed the name of argument dec --> decrem
113  * (in order not to shadow some globally defined dec).
114  *
115  * Revision 1.24 2003/10/29 11:02:54 e_gourgoulhon
116  * Functions dec_dzpuis and inc_dzpuis have now an integer argument to
117  * specify by which amount dzpuis is to be increased.
118  * Accordingly methods dec2_dzpuis and inc2_dzpuis have been suppressed
119  *
120  * Revision 1.23 2003/10/27 10:49:48 e_gourgoulhon
121  * Slightly modified operator<<.
122  *
123  * Revision 1.22 2003/10/19 19:57:00 e_gourgoulhon
124  * -- Added new method spectral_display
125  * -- slightly rearranged the operator<<
126  *
127  * Revision 1.21 2003/10/16 15:27:46 e_gourgoulhon
128  * Name of method annule(int ) changed to annule_domain(int ).
129  *
130  * Revision 1.20 2003/10/16 14:21:36 j_novak
131  * The calculation of the divergence of a Tensor is now possible.
132  *
133  * Revision 1.19 2003/10/13 13:52:40 j_novak
134  * Better managment of derived quantities.
135  *
136  * Revision 1.18 2003/10/11 16:47:10 e_gourgoulhon
137  * Suppressed the call to Ibtl::set_etat_qcq() after the construction
138  * of the Itbl's, thanks to the new property of the Itbl class.
139  *
140  * Revision 1.17 2003/10/11 14:48:40 e_gourgoulhon
141  * Line 344: suppressed assert(resu == -1)
142  * and added a break command to quit the loop.
143  *
144  * Revision 1.16 2003/10/08 14:24:09 j_novak
145  * replaced mult_r_zec with mult_r_ced
146  *
147  * Revision 1.15 2003/10/07 15:25:38 e_gourgoulhon
148  * Added call to del_derive_met in del_deriv().
149  *
150  * Revision 1.14 2003/10/07 09:10:00 j_novak
151  * Use of ::contract instead of up()
152  *
153  * Revision 1.13 2003/10/06 20:51:43 e_gourgoulhon
154  * In methods set: changed name "indices" to "idx" to avoid shadowing
155  * of class member.
156  *
157  * Revision 1.12 2003/10/06 16:17:31 j_novak
158  * Calculation of contravariant derivative and Ricci scalar.
159  *
160  * Revision 1.11 2003/10/06 13:58:48 j_novak
161  * The memory management has been improved.
162  * Implementation of the covariant derivative with respect to the exact Tensor
163  * type.
164  *
165  * Revision 1.10 2003/10/05 21:11:22 e_gourgoulhon
166  * - Added method std_spectral_base().
167  * - Removed method change_triad() from this file.
168  *
169  * Revision 1.9 2003/10/03 15:09:38 j_novak
170  * Improved display
171  *
172  * Revision 1.8 2003/10/03 11:21:48 j_novak
173  * More methods for the class Metric
174  *
175  * Revision 1.7 2003/10/01 11:56:31 e_gourgoulhon
176  * Corrected error: '=' replaced by '==' in two assert tests.
177  *
178  * Revision 1.6 2003/09/30 08:38:23 j_novak
179  * added a header typeinfo
180  *
181  * Revision 1.5 2003/09/29 12:52:57 j_novak
182  * Methods for changing the triad are implemented.
183  *
184  * Revision 1.4 2003/09/26 14:33:53 j_novak
185  * Arithmetic functions for the class Tensor
186  *
187  * Revision 1.3 2003/09/24 15:10:54 j_novak
188  * Suppression of the etat flag in class Tensor (still present in Scalar)
189  *
190  * Revision 1.2 2003/09/23 08:52:17 e_gourgoulhon
191  * new version
192  *
193  * Revision 1.1 2003/09/22 12:52:51 e_gourgoulhon
194  * First version: not ready yet !
195  *
196  *
197  * $Header: /cvsroot/Lorene/C++/Source/Tensor/tensor.C,v 1.48 2025/03/31 08:29:50 j_novak Exp $
198  *
199  */
200 
201 // Headers Lorene
202 #include "metric.h"
203 #include "utilitaires.h"
204 
205  //--------------//
206  // Constructors //
207  //--------------//
208 
209 // Standard constructor
210 // --------------------
211 namespace Lorene {
212 Tensor::Tensor(const Map& map, int val, const Itbl& tipe,
213  const Base_vect& triad_i)
214  : mp(&map), valence(val), triad(&triad_i), type_indice(tipe),
215  n_comp(int(pow(3., val))){
216 
217  // Des verifs :
218  assert (valence >= 0) ;
219  assert (tipe.get_ndim() == 1) ;
220  assert (valence == tipe.get_dim(0)) ;
221  for (int i=0 ; i<valence ; i++)
222  assert ((tipe(i) == COV) || (tipe(i) == CON)) ;
223 
224  cmp = new Scalar*[n_comp] ;
225 
226  for (int i=0 ; i<n_comp ; i++)
227  cmp[i] = new Scalar(map) ;
228 
229  set_der_0x0() ;
230 
231 }
232 
233 // Standard constructor with the triad passed as a pointer
234 // -------------------------------------------------------
235 Tensor::Tensor(const Map& map, int val, const Itbl& tipe,
236  const Base_vect* triad_i)
237  : mp(&map), valence(val), triad(triad_i), type_indice(tipe),
238  n_comp(int(pow(3., val))){
239 
240  // Des verifs :
241  assert (valence >= 0) ;
242  assert (tipe.get_ndim() == 1) ;
243  assert (valence == tipe.get_dim(0)) ;
244  for (int i=0 ; i<valence ; i++)
245  assert ((tipe(i) == COV) || (tipe(i) == CON)) ;
246 
247  cmp = new Scalar*[n_comp] ;
248 
249  for (int i=0 ; i<n_comp ; i++)
250  cmp[i] = new Scalar(map) ;
251 
252  set_der_0x0() ;
253 
254 }
255 
256 
257 
258 
259 // Standard constructor when all the indices are of the same type
260 // --------------------------------------------------------------
261 Tensor::Tensor(const Map& map, int val, int tipe, const Base_vect& triad_i)
262  : mp(&map), valence(val), triad(&triad_i), type_indice(val),
263  n_comp(int(pow(3., val))){
264 
265  // Des verifs :
266  assert (valence >= 0) ;
267  assert ((tipe == COV) || (tipe == CON)) ;
268 
269  type_indice = tipe ;
270 
271  cmp = new Scalar*[n_comp] ;
272 
273  for (int i=0 ; i<n_comp ; i++)
274  cmp[i] = new Scalar(map) ;
275 
276  set_der_0x0() ;
277 }
278 
279 // Copy constructor
280 // ----------------
281 Tensor::Tensor (const Tensor& source) :
282  mp(source.mp), valence(source.valence), triad(source.triad),
283  type_indice(source.type_indice) {
284 
285  n_comp = int(pow(3., valence)) ;
286 
287  cmp = new Scalar*[n_comp] ;
288  for (int i=0 ; i<n_comp ; i++) {
289 
290  // the following writing takes into account the case where
291  // source belongs to a derived class of Tensor with a different
292  // storage of components :
293 
294  int place_source = source.position(indices(i)) ;
295  cmp[i] = new Scalar(*source.cmp[place_source]) ;
296  }
297 
298  set_der_0x0() ;
299 
300 }
301 
302 
303 // Constructor from a file
304 // -----------------------
305 Tensor::Tensor(const Map& mapping, const Base_vect& triad_i, FILE* fd)
306  : mp(&mapping), triad(&triad_i), type_indice(fd){
307 
308  fread_be(&valence, sizeof(int), 1, fd) ;
309 
310  if (valence != 0) {
311  Base_vect* triad_fich = Base_vect::bvect_from_file(fd) ;
312  assert( *triad_fich == *triad) ;
313  delete triad_fich ;
314  }
315  else{
316  triad = 0x0 ;
317  }
318 
319  fread_be(&n_comp, sizeof(int), 1, fd) ;
320 
321  cmp = new Scalar*[n_comp] ;
322  for (int i=0 ; i<n_comp ; i++)
323  cmp[i] = new Scalar(*mp, *(mp->get_mg()), fd) ;
324 
325  set_der_0x0() ;
326 }
327 
328 
329 // Constructor for a scalar field: to be used by the derived
330 // class {\tt Scalar}
331 //-----------------------------------------------------------
332 Tensor::Tensor(const Map& map) : mp(&map), valence(0), triad(0x0),
333  type_indice(0), n_comp(1) {
334 
335  cmp = new Scalar*[n_comp] ;
336  cmp[0] = 0x0 ;
337 
338  set_der_0x0() ;
339 }
340 
341 
342 // Constructor used by the derived classes
343 // ---------------------------------------
344 Tensor::Tensor (const Map& map, int val, const Itbl& tipe, int compo,
345  const Base_vect& triad_i) :
346  mp(&map), valence(val), triad(&triad_i), type_indice(tipe), n_comp(compo)
347 {
348 
349  // Des verifs :
350  assert (valence >= 0) ;
351  assert (tipe.get_ndim() == 1) ;
352  assert (n_comp > 0) ;
353  assert (valence == tipe.get_dim(0)) ;
354  for (int i=0 ; i<valence ; i++)
355  assert ((tipe(i) == COV) || (tipe(i) == CON)) ;
356 
357  cmp = new Scalar*[n_comp] ;
358 
359  for (int i=0 ; i<n_comp ; i++)
360  cmp[i] = new Scalar(map) ;
361 
362  set_der_0x0() ;
363 
364 }
365 
366 // Constructor used by the derived classes when all the indices are of
367 // the same type.
368 // -------------------------------------------------------------------
369 Tensor::Tensor (const Map& map, int val, int tipe, int compo,
370  const Base_vect& triad_i) :
371  mp(&map), valence(val), triad(&triad_i), type_indice(val), n_comp(compo)
372 {
373 
374  // Des verifs :
375  assert (valence >= 0) ;
376  assert (n_comp >= 0) ;
377  assert ((tipe == COV) || (tipe == CON)) ;
378 
379  type_indice = tipe ;
380 
381  cmp = new Scalar*[n_comp] ;
382 
383  for (int i=0 ; i<n_comp ; i++)
384  cmp[i] = new Scalar(map) ;
385 
386  set_der_0x0() ;
387 
388 }
389 
390  //--------------//
391  // Destructor //
392  //--------------//
393 
394 
396 
398 
399  for (int i=0 ; i<n_comp ; i++) {
400  if (cmp[i] != 0x0)
401  delete cmp[i] ;
402  }
403  delete [] cmp ;
404 }
405 
406 
407 
408 void Tensor::del_deriv() const {
409 
410  for (int i=0; i<N_MET_MAX; i++)
411  del_derive_met(i) ;
412 
413  set_der_0x0() ;
414 
415 }
416 
417 void Tensor::set_der_0x0() const {
418 
419  for (int i=0; i<N_MET_MAX; i++)
420  set_der_met_0x0(i) ;
421 
422 }
423 
424 void Tensor::del_derive_met(int j) const {
425 
426  assert( (j>=0) && (j<N_MET_MAX) ) ;
427 
428  if (met_depend[j] != 0x0) {
429  for (int i=0 ; i<N_TENSOR_DEPEND ; i++)
430  if (met_depend[j]->tensor_depend[i] == this)
431  met_depend[j]->tensor_depend[i] = 0x0 ;
432  if (p_derive_cov[j] != 0x0)
433  delete p_derive_cov[j] ;
434  if (p_derive_con[j] != 0x0)
435  delete p_derive_con[j] ;
436  if (p_divergence[j] != 0x0)
437  delete p_divergence[j] ;
438 
439  set_der_met_0x0(j) ;
440  }
441 }
442 
443 void Tensor::set_der_met_0x0(int i) const {
444 
445  assert( (i>=0) && (i<N_MET_MAX) ) ;
446  met_depend[i] = 0x0 ;
447  p_derive_cov[i] = 0x0 ;
448  p_derive_con[i] = 0x0 ;
449  p_divergence[i] = 0x0 ;
450 
451 }
452 
453 int Tensor::get_place_met(const Metric& metre) const {
454  int resu = -1 ;
455  for (int i=0; i<N_MET_MAX; i++)
456  if (met_depend[i] == &metre) {
457  resu = i ;
458  break ;
459  }
460  return resu ;
461 }
462 
463 void Tensor::set_dependance (const Metric& met) const {
464 
465  int nmet = 0 ;
466  bool deja = false ;
467  for (int i=0; i<N_MET_MAX; i++) {
468  if (met_depend[i] == &met) deja = true ;
469  if ((!deja) && (met_depend[i] != 0x0)) nmet++ ;
470  }
471  if (nmet == N_MET_MAX) {
472  cout << "Too many metrics in Tensor::set_dependances" << endl ;
473  abort() ;
474  }
475  if (!deja) {
476  int conte = 0 ;
477  while ((conte < N_TENSOR_DEPEND) && (met.tensor_depend[conte] != 0x0))
478  conte ++ ;
479 
480  if (conte == N_TENSOR_DEPEND) {
481  cout << "Too many dependancies in Tensor::set_dependances " << endl ;
482  abort() ;
483  }
484  else {
485  met.tensor_depend[conte] = this ;
486  met_depend[nmet] = &met ;
487  }
488  }
489 }
490 
492 
493  del_deriv() ;
494  for (int i=0 ; i<n_comp ; i++) {
495  cmp[i]->set_etat_qcq() ;
496  }
497 }
498 
500 
501  del_deriv() ;
502  for (int i=0 ; i<n_comp ; i++) {
503  cmp[i]->set_etat_nondef() ;
504  }
505 }
506 
508 
509  del_deriv() ;
510  for (int i=0 ; i<n_comp ; i++) {
511  cmp[i]->set_etat_zero() ;
512  }
513 }
514 
515 
516 // Allocates everything
517 // --------------------
519 
520  del_deriv() ;
521  for (int i=0 ; i<n_comp ; i++) {
522  cmp[i]->allocate_all() ;
523  }
524 
525 }
526 
527 
528 
529 void Tensor::set_triad(const Base_vect& bi) {
530 
531  triad = &bi ;
532 
533 }
534 
535 int Tensor::position (const Itbl& idx) const {
536 
537  assert (idx.get_ndim() == 1) ;
538  assert (idx.get_dim(0) == valence) ;
539 
540  for (int i=0 ; i<valence ; i++)
541  assert ((idx(i)>=1) && (idx(i)<=3)) ;
542  int res = 0 ;
543  for (int i=0 ; i<valence ; i++)
544  res = 3*res+(idx(i)-1) ;
545 
546  return res;
547 }
548 
549 Itbl Tensor::indices (int place) const {
550 
551  assert ((place >= 0) && (place < n_comp)) ;
552 
553  Itbl res(valence) ;
554 
555  for (int i=valence-1 ; i>=0 ; i--) {
556  res.set(i) = div(place, 3).rem ;
557  place = int((place-res(i))/3) ;
558  res.set(i)++ ;
559  }
560  return res ;
561 }
562 
563 void Tensor::operator=(const Tensor& t) {
564 
565  assert (valence == t.valence) ;
566 
567  triad = t.triad ;
568 
569  for (int i=0 ; i<valence ; i++)
570  assert(t.type_indice(i) == type_indice(i)) ;
571 
572  for (int i=0 ; i<n_comp ; i++) {
573  int place_t = t.position(indices(i)) ;
574  *cmp[i] = *t.cmp[place_t] ;
575  }
576 
577  del_deriv() ;
578 
579 }
580 
581 void Tensor::operator+=(const Tensor& t) {
582 
583  assert (valence == t.valence) ;
584  assert (triad == t.triad) ;
585  for (int i=0 ; i<valence ; i++)
586  assert(t.type_indice(i) == type_indice(i)) ;
587 
588  for (int i=0 ; i<n_comp ; i++) {
589  int place_t = t.position(indices(i)) ;
590  *cmp[i] += *t.cmp[place_t] ;
591  }
592 
593  del_deriv() ;
594 
595 }
596 
597 void Tensor::operator-=(const Tensor& t) {
598 
599  assert (valence == t.valence) ;
600  assert (triad == t.triad) ;
601  for (int i=0 ; i<valence ; i++)
602  assert(t.type_indice(i) == type_indice(i)) ;
603 
604  for (int i=0 ; i<n_comp ; i++) {
605  int place_t = t.position(indices(i)) ;
606  *cmp[i] -= *t.cmp[place_t] ;
607  }
608 
609  del_deriv() ;
610 
611 }
612 
613 
614 
615 // Affectation d'un tenseur d'ordre 2 :
616 Scalar& Tensor::set(int ind1, int ind2) {
617 
618  assert (valence == 2) ;
619 
620  Itbl ind (valence) ;
621  ind.set(0) = ind1 ;
622  ind.set(1) = ind2 ;
623 
624  int place = position(ind) ;
625 
626  del_deriv() ;
627  return *cmp[place] ;
628 }
629 
630 // Affectation d'un tenseur d'ordre 3 :
631 Scalar& Tensor::set(int ind1, int ind2, int ind3) {
632 
633  assert (valence == 3) ;
634 
635  Itbl idx(valence) ;
636  idx.set(0) = ind1 ;
637  idx.set(1) = ind2 ;
638  idx.set(2) = ind3 ;
639  int place = position(idx) ;
640  del_deriv() ;
641 
642  return *cmp[place] ;
643 }
644 
645 
646 // Affectation d'un tenseur d'ordre 4 :
647 Scalar& Tensor::set(int ind1, int ind2, int ind3, int ind4) {
648 
649  assert (valence == 4) ;
650 
651  Itbl idx(valence) ;
652  idx.set(0) = ind1 ;
653  idx.set(1) = ind2 ;
654  idx.set(2) = ind3 ;
655  idx.set(3) = ind4 ;
656  int place = position(idx) ;
657  del_deriv() ;
658 
659  return *cmp[place] ;
660 }
661 
662 
663 // Affectation cas general
664 Scalar& Tensor::set(const Itbl& idx) {
665 
666  assert (idx.get_ndim() == 1) ;
667  assert (idx.get_dim(0) == valence) ;
668 
669  int place = position(idx) ;
670 
671  del_deriv() ;
672  return *cmp[place] ;
673 }
674 
675 // Annulation dans des domaines
677 
678  annule(l, l) ;
679 }
680 
681 void Tensor::annule(int l_min, int l_max) {
682 
683  // Cas particulier: annulation globale :
684  if ( (l_min == 0) && (l_max == mp->get_mg()->get_nzone()-1) ) {
685  set_etat_zero() ;
686  return ;
687  }
688 
689  // Annulation des composantes:
690  for (int i=0 ; i<n_comp ; i++) {
691  cmp[i]->annule(l_min, l_max) ;
692  }
693 
694  // The derived members are no longer up to date:
695  del_deriv() ;
696 
697 }
698 
699 
700 void Tensor::annule_extern_cn(int lrac, int deg) {
701 
702  // Not applicable in the nucleus nor the CED:
703  assert( mp->get_mg()->get_type_r(lrac) == FIN ) ;
704 
705  int nz = mp->get_mg()->get_nzone() ;
706 #ifndef NDEBUG
707  if ((2*deg+1) >= mp->get_mg()->get_nr(lrac))
708  cout << "Tensor::annule_extern_cn : \n"
709  << "WARNING!! \n"
710  << "The number of coefficients in r is too low \n"
711  << "to do a clean matching..." << endl ;
712 #endif
713  // Boundary of domain lrac
714  double r_min = mp->val_r(lrac, -1., 0., 0.) ;
715  double r_max = mp->val_r(lrac, 1., 0., 0.) ;
716 
717  //Definition of binomial coefficients array
718  Itbl binom(deg+1, deg+1) ;
719  binom.annule_hard() ;
720  binom.set(0,0) = 1 ;
721  for (int n=1; n<=deg; n++) {
722  binom.set(n,0) = 1 ;
723  for (int k=1; k<=n; k++)
724  binom.set(n,k) = binom(n-1, k) + binom(n-1, k-1) ;
725  }
726 
727  // Coefficient of the second polynomial factor
728  Tbl coef(deg+1) ;
729  coef.set_etat_qcq() ;
730  coef.set(deg) = 1 ;
731  int sg = -1 ;
732  for (int i=deg-1; i>=0; i--) {
733 
734  coef.set(i) = double(r_max*(i+1)*coef(i+1)
735  + sg*binom(deg,i)*(2*deg+1)*pow(r_min,deg-i))
736  / double(deg+i+1) ;
737  sg *= -1 ;
738  }
739 
740  // Normalization to have 1 at r_min
741  double aa = coef(deg) ;
742  for (int i = deg-1; i>=0; i--)
743  aa = r_min*aa + coef(i) ;
744  aa *= pow(r_min - r_max, deg+1) ;
745  aa = 1/aa ;
746 
747  Mtbl mr = mp->r ;
748  Tbl rr = mr(lrac) ;
749 
750  Tbl poly(rr) ;
751  poly = coef(deg) ;
752  for (int i=deg-1; i>=0; i--)
753  poly = rr*poly + coef(i) ;
754  poly *= aa*pow((rr-r_max), deg+1) ;
755 
756  Scalar rac(*mp) ;
757  rac.allocate_all() ;
758  for (int l=0; l<lrac; l++) rac.set_domain(l) = 1 ;
759  rac.set_domain(lrac) = poly ;
760  rac.annule(lrac+1,nz-1) ;
761  rac.std_spectral_base() ;
762 
763  for (int ic=0; ic<n_comp; ic++) *(cmp[ic]) *= rac ;
764 
765  del_deriv() ;
766 }
767 
768 
769 
770 const Scalar& Tensor::operator()(int indice1, int indice2) const {
771 
772  assert(valence == 2) ;
773 
774  Itbl idx(2) ;
775  idx.set(0) = indice1 ;
776  idx.set(1) = indice2 ;
777  return *cmp[position(idx)] ;
778 
779 }
780 
781 const Scalar& Tensor::operator()(int indice1, int indice2, int indice3) const {
782 
783  assert(valence == 3) ;
784 
785  Itbl idx(3) ;
786  idx.set(0) = indice1 ;
787  idx.set(1) = indice2 ;
788  idx.set(2) = indice3 ;
789  return *cmp[position(idx)] ;
790 }
791 
792 
793 const Scalar& Tensor::operator()(int indice1, int indice2, int indice3,
794  int indice4) const {
795 
796  assert(valence == 4) ;
797 
798  Itbl idx(4) ;
799  idx.set(0) = indice1 ;
800  idx.set(1) = indice2 ;
801  idx.set(2) = indice3 ;
802  idx.set(3) = indice4 ;
803  return *cmp[position(idx)] ;
804 }
805 
806 
807 
808 const Scalar& Tensor::operator()(const Itbl& ind) const {
809 
810  assert (ind.get_ndim() == 1) ;
811  assert (ind.get_dim(0) == valence) ;
812  return *cmp[position(ind)] ;
813 
814 }
815 
816 
817 // Gestion de la CED :
818 void Tensor::dec_dzpuis(int decrem) {
819 
820  for (int i=0 ; i<n_comp ; i++)
821  cmp[i]->dec_dzpuis(decrem) ;
822 
823  del_deriv() ;
824 }
825 
826 void Tensor::inc_dzpuis(int inc) {
827 
828  for (int i=0 ; i<n_comp ; i++)
829  cmp[i]->inc_dzpuis(inc) ;
830 
831  del_deriv() ;
832 }
833 
834 
835 // Le cout :
836 ostream& operator<<(ostream& flux, const Tensor &source ) {
837 
838  flux << '\n' ;
839  flux << "Lorene class : " << typeid(source).name()
840  << " Valence : " << source.valence << '\n' ;
841 
842  if (source.get_triad() != 0x0) {
843  flux << "Vectorial basis (triad) on which the components are defined :"
844  << '\n' ;
845  flux << *(source.get_triad()) << '\n' ;
846  }
847 
848  if (source.valence != 0) {
849  flux << "Type of the indices : " ;
850  for (int i=0 ; i<source.valence ; i++) {
851  flux << "index " << i << " : " ;
852  if (source.type_indice(i) == CON)
853  flux << " contravariant." << '\n' ;
854  else
855  flux << " covariant." << '\n' ;
856  if ( i < source.valence-1 ) flux << " " ;
857  }
858  flux << '\n' ;
859  }
860 
861  for (int i=0 ; i<source.n_comp ; i++) {
862 
863  if (source.valence == 0) {
864  flux <<
865  "===================== Scalar field ========================= \n" ;
866  }
867  else {
868  flux << "================ Component " ;
869  Itbl num_indices = source.indices(i) ;
870  for (int j=0 ; j<source.valence ; j++) {
871  flux << " " << num_indices(j) ;
872  }
873  flux << " ================ \n" ;
874  }
875  flux << '\n' ;
876 
877  flux << *source.cmp[i] << '\n' ;
878  }
879 
880  return flux ;
881 }
882 
883 
884 void Tensor::spectral_display(const char* comment,
885  double thres, int precis, ostream& ost) const {
886 
887  if (comment != 0x0) {
888  ost << comment << " : " << endl ;
889  }
890 
891  ost << "Lorene class : " << typeid(*this).name()
892  << " Valence : " << valence << '\n' ;
893 
894  for (int ic=0; ic<n_comp; ic++) {
895 
896  if (valence == 0) {
897  ost <<
898  "===================== Scalar field ========================= \n" ;
899  }
900  else {
901  ost << "================ Component " ;
902  Itbl num_indices = indices(ic) ;
903  for (int j=0 ; j<valence ; j++) {
904  ost << " " << num_indices(j) ;
905  }
906  ost << " ================ \n" ;
907  }
908  ost << '\n' ;
909 
910  cmp[ic]->spectral_display(0x0, thres, precis, ost) ;
911  ost << '\n' ;
912  }
913 }
914 
915 
916 void Tensor::sauve(FILE* fd) const {
917 
918  type_indice.sauve(fd) ; // type des composantes
919  fwrite_be(&valence, sizeof(int), 1, fd) ; // la valence
920 
921  if (valence != 0) {
922  triad->sauve(fd) ; // Vectorial basis
923  }
924 
925  fwrite_be(&n_comp, sizeof(int), 1, fd) ; // nbre composantes
926  for (int i=0 ; i<n_comp ; i++)
927  cmp[i]->sauve(fd) ;
928 
929 }
930 
931 
932 
933 
934 
935 // Sets the standard spectal bases of decomposition for each component
937 
938  switch (valence) {
939 
940  case 0 : {
941  cmp[0]->std_spectral_base() ;
942  break ;
943  }
944 
945  case 1 : {
946  cout <<
947  "Tensor::std_spectral_base: should not be called on a Tensor"
948  << " of valence 1 but on a Vector !" << endl ;
949  abort() ;
950  break ;
951  }
952 
953  case 2 : {
954 
955  Base_val** bases = 0x0 ;
956  if( triad->identify() == (mp->get_bvect_cart()).identify() ) {
957  bases = mp->get_mg()->std_base_vect_cart() ;
958  }
959  else {
960  assert( triad->identify() == (mp->get_bvect_spher()).identify()) ;
961  bases = mp->get_mg()->std_base_vect_spher() ;
962  }
963 
964  Itbl ind(2) ;
965  for (int i=0 ; i<n_comp ; i++) {
966  ind = indices(i) ;
967  cmp[i]->set_spectral_base( (*bases[ind(0)-1]) *
968  (*bases[ind(1)-1]) ) ;
969  }
970 
971  for (int i=0 ; i<3 ; i++) {
972  delete bases[i] ;
973  }
974  delete [] bases ;
975  break ;
976 
977  }
978 
979 
980  default : {
981 
982  cout << "Tensor::std_spectral_base: the case valence = " << valence
983  << " is not treated yet !" << endl ;
984  abort() ;
985  break ;
986  }
987  }
988 }
989 
990 // Sets the standard spectal bases of decomposition for each component (odd in the nucleus)
991 
993 
994  switch (valence) {
995 
996  case 0 : {
997  cmp[0]->std_spectral_base_odd() ;
998  break ;
999  }
1000 
1001  default : {
1002 
1003  cout << "Tensor::std_spectral_base_odd: the case valence = " << valence
1004  << " is not treated yet !" << endl ;
1005  abort() ;
1006  break ;
1007  }
1008  }
1009 }
1010 
1011 
1012 const Tensor& Tensor::derive_cov(const Metric& metre) const {
1013 
1014  set_dependance(metre) ;
1015  int j = get_place_met(metre) ;
1016  assert ((j>=0) && (j<N_MET_MAX)) ;
1017  if (p_derive_cov[j] == 0x0) {
1018  p_derive_cov[j] = metre.connect().p_derive_cov(*this) ;
1019  }
1020  return *p_derive_cov[j] ;
1021 }
1022 
1023 
1024 const Tensor& Tensor::derive_con(const Metric& metre) const {
1025 
1026  set_dependance(metre) ;
1027  int j = get_place_met(metre) ;
1028  assert ((j>=0) && (j<N_MET_MAX)) ;
1029  if (p_derive_con[j] == 0x0) {
1030 
1031  if (valence == 0) {
1032  p_derive_con[j] =
1033  new Vector( contract(derive_cov(metre), 0, metre.con(), 0) ) ;
1034  }
1035  else {
1036  const Tensor_sym* tsym = dynamic_cast<const Tensor_sym*>(this) ;
1037 
1038  if (tsym != 0x0) { // symmetric case, preserved by derive_con
1039  const Tensor& dercov = derive_cov(metre) ;
1040  Itbl type_ind = dercov.get_index_type() ;
1041  type_ind.set(valence) = CON ;
1042  p_derive_con[j] = new Tensor_sym(*mp, valence+1, type_ind, *triad,
1043  tsym->sym_index1(), tsym->sym_index2()) ;
1044 
1045  *(p_derive_con[j]) = contract(dercov, valence, metre.con(), 0) ;
1046  // valence is the number of the last index of derive_cov
1047  // (the "derivation" index)
1048  }
1049  else { // general case, no symmetry
1050 
1051  p_derive_con[j] = new Tensor( contract(derive_cov(metre), valence,
1052  metre.con(), 0) ) ;
1053  // valence is the number of the last index of derive_cov
1054  // (the "derivation" index)
1055  }
1056 
1057  }
1058 
1059  }
1060 
1061  return *p_derive_con[j] ;
1062 
1063 }
1064 
1065 const Tensor& Tensor::divergence(const Metric& metre) const {
1066 
1067  set_dependance(metre) ;
1068  int j = get_place_met(metre) ;
1069  assert ((j>=0) && (j<N_MET_MAX)) ;
1070  if (p_divergence[j] == 0x0) {
1071  p_divergence[j] = metre.connect().p_divergence(*this) ;
1072  }
1073  return *p_divergence[j] ;
1074 }
1075 
1076 void Tensor::exponential_filter_r(int lzmin, int lzmax, int p,
1077  double alpha) {
1078  if( triad->identify() == (mp->get_bvect_cart()).identify() )
1079  for (int i=0; i<n_comp; i++)
1080  cmp[i]->exponential_filter_r(lzmin, lzmax, p, alpha) ;
1081  else {
1082  cout << "Tensor::exponential_filter_r : " << endl ;
1083  cout << "Only Cartesian triad is implemented!" << endl ;
1084  cout << "Exiting..." << endl ;
1085  abort() ;
1086  }
1087 }
1088 
1089 void Tensor::exponential_filter_ylm(int lzmin, int lzmax, int p,
1090  double alpha) {
1091  if( triad->identify() == (mp->get_bvect_cart()).identify() )
1092  for (int i=0; i<n_comp; i++)
1093  cmp[i]->exponential_filter_ylm(lzmin, lzmax, p, alpha) ;
1094  else {
1095  cout << "Tensor::exponential_filter_ylm : " << endl ;
1096  cout << "Only Cartesian triad is implemented!" << endl ;
1097  cout << "Exiting..." << endl ;
1098  abort() ;
1099  }
1100 }
1101 
1102 void Tensor::exponential_filter_ylm_phi(int lzmin, int lzmax, int p_r, int p_tet, int p_phi,
1103  double alpha) {
1104  if( triad->identify() == (mp->get_bvect_cart()).identify() )
1105  for (int i=0; i<n_comp; i++)
1106  cmp[i]->exponential_filter_ylm_phi(lzmin, lzmax, p_r, p_tet, p_phi, alpha) ;
1107  else {
1108  cout << "Tensor::exponential_filter_ylm : " << endl ;
1109  cout << "Only Cartesian triad is implemented!" << endl ;
1110  cout << "Exiting..." << endl ;
1111  abort() ;
1112  }
1113 }
1114 
1115 
1116 void Tensor::copy_coefs_from_smaller_grid(const Tensor& tens_small_grid) {
1117 
1118  for (int i=0; i<n_comp; i++) {
1119  cmp[i]->copy_coefs_from_small_grid(*(tens_small_grid.cmp[i])) ;
1120  }
1121 }
1122 
1123 void Tensor::copy_coefs_from_larger_grid(const Tensor& tens_large_grid) {
1124 
1125  for (int i=0; i<n_comp; i++) {
1126  cmp[i]->copy_coefs_from_large_grid(*(tens_large_grid.cmp[i])) ;
1127  }
1128 }
1129 
1130 
1131 
1132 
1133 } // namespace Lorene
void annule_domain(int l)
Sets the Tensor to zero in a given domain.
Definition: tensor.C:676
Base_val ** std_base_vect_cart() const
Returns the standard spectral bases for the Cartesian components of a vector.
Metric for tensor calculation.
Definition: metric.h:90
virtual void set_etat_qcq()
Sets the logical state of all components to ETATQCQ (ordinary state).
Definition: tensor.C:491
virtual const Sym_tensor & con() const
Read-only access to the contravariant representation.
Definition: metric.C:293
virtual void operator=(const Tensor &)
Assignment to a Tensor.
Definition: tensor.C:563
int & set(int i)
Read/write of a particular element (index i ) (1D case)
Definition: itbl.h:247
void copy_coefs_from_small_grid(const Scalar &scal_small_grid)
Copies the content of the argument Scalar to this defined on a larger grid, with similar mappings...
Definition: scalar_manip.C:505
int n_comp
Number of stored components, depending on the symmetry.
Definition: tensor.h:324
virtual void set_etat_zero()
Sets the logical state to ETATZERO (zero).
Definition: scalar.C:330
virtual void annule(int l_min, int l_max)
Sets the Scalar to zero in several domains.
Definition: scalar.C:397
int sym_index1() const
Number of the first symmetric index (0<= id_sym1 < valence )
Definition: tensor.h:1190
Multi-domain array.
Definition: mtbl.h:118
const Base_vect_spher & get_bvect_spher() const
Returns the orthonormal vectorial basis associated with the coordinates of the mapping.
Definition: map.h:810
int get_ndim() const
Gives the number of dimensions (ie dim.ndim )
Definition: itbl.h:323
virtual void del_derive_met(int) const
Logical destructor of the derivatives depending on the i-th element of met_depend ...
Definition: tensor.C:424
Lorene prototypes.
Definition: app_hor.h:67
void set_dependance(const Metric &) const
To be used to describe the fact that the derivatives members have been calculated with met ...
Definition: tensor.C:463
const Mg3d * get_mg() const
Gives the Mg3d on which the mapping is defined.
Definition: map.h:792
void annule_hard()
Sets the Itbl to zero in a hard way.
Definition: itbl.C:275
double & set(int i)
Read/write of a particular element (index i) (1D case)
Definition: tbl.h:301
Tensor field of valence 0 (or component of a tensorial field).
Definition: scalar.h:399
virtual void set_etat_nondef()
Sets the logical state of all components to ETATNONDEF (undefined state).
Definition: tensor.C:499
virtual void sauve(FILE *) const
Save in a file.
Definition: base_vect.C:153
Base class for coordinate mappings.
Definition: map.h:697
void copy_coefs_from_smaller_grid(const Tensor &tens_small_grid)
Copies the content of the argument Tensor to this defined on a larger grid, with similar mappings...
Definition: tensor.C:1116
virtual Tensor * p_derive_cov(const Tensor &tens) const
Computes the covariant derivative of a tensor (with respect to the current connection).
Definition: connection.C:310
int sym_index2() const
Number of the second symmetric index (id_sym1 < id_sym2 < valence )
Definition: tensor.h:1195
Basic integer array class.
Definition: itbl.h:122
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
const Base_vect * triad
Vectorial basis (triad) with respect to which the tensor components are defined.
Definition: tensor.h:315
virtual void sauve(FILE *) const
Save in a binary file.
Definition: tensor.C:916
virtual void std_spectral_base_odd()
Sets the standard odd spectal bases of decomposition for each component.
Definition: tensor.C:992
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
virtual void set_etat_nondef()
Sets the logical state to ETATNONDEF (undefined).
Definition: scalar.C:350
virtual void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition: scalar.C:359
int get_place_met(const Metric &) const
Returns the position of the pointer on metre in the array met_depend .
Definition: tensor.C:453
Tensor field of valence 1.
Definition: vector.h:188
Vectorial bases (triads) with respect to which the tensorial components are defined.
Definition: base_vect.h:105
Tbl & set_domain(int l)
Read/write of the value in a given domain.
Definition: scalar.h:627
virtual void exponential_filter_ylm_phi(int lzmin, int lzmax, int p_r, int p_tet, int p_phi, double alpha=-16.)
Applies exponential filters to all components (see Scalar::exponential_filter_ylm_phi )...
Definition: tensor.C:1102
void copy_coefs_from_larger_grid(const Tensor &tens_large_grid)
Copies the content of the argument Tensor to this defined on a smaller grid, with similar mappings...
Definition: tensor.C:1123
static Base_vect * bvect_from_file(FILE *)
Construction of a vectorial basis from a file (see sauve(FILE* ) ).
Tensor * p_derive_cov[N_MET_MAX]
Array of pointers on the covariant derivatives of this with respect to various metrics.
Definition: tensor.h:347
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 set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition: tbl.C:364
virtual void exponential_filter_r(int lzmin, int lzmax, int p, double alpha=-16.)
Applies exponential filters to all components (see Scalar::exponential_filter_r ).
Definition: tensor.C:1076
virtual double val_r(int l, double xi, double theta, double pphi) const =0
Returns the value of the radial coordinate r for a given in a given domain.
void operator-=(const Tensor &)
-= Tensor
Definition: tensor.C:597
Itbl type_indice
1D array of integers (class Itbl ) of size valence containing the type of each index: COV for a cova...
Definition: tensor.h:322
const Base_vect * get_triad() const
Returns the vectorial basis (triad) on which the components are defined.
Definition: tensor.h:907
virtual int position(const Itbl &ind) const
Returns the position in the array cmp of a component given by its indices.
Definition: tensor.C:535
int get_index_type(int i) const
Gives the type (covariant or contravariant) of the index number i .
Definition: tensor.h:927
virtual void allocate_all()
Performs the memory allocation of all the elements, down to the double arrays of the Tbl s...
Definition: tensor.C:518
virtual const Connection & connect() const
Returns the connection.
Definition: metric.C:304
Scalar ** cmp
Array of size n_comp of pointers onto the components.
Definition: tensor.h:327
virtual void spectral_display(const char *comment=0x0, double threshold=1.e-7, int precision=4, ostream &ostr=cout) const
Displays the spectral coefficients and the associated basis functions.
Definition: scalar.C:747
void sauve(FILE *) const
Save in a file.
Definition: itbl.C:229
Base_val ** std_base_vect_spher() const
Returns the standard spectral bases for the spherical components of a vector.
int get_nzone() const
Returns the number of domains.
Definition: grilles.h:467
virtual int identify() const =0
Returns a number to identify the sub-classe of Base_vect the object belongs to.
virtual void std_spectral_base_odd()
Sets the spectral bases of the Valeur va to the standard odd ones for a scalar field.
Definition: scalar.C:797
int fwrite_be(const int *aa, int size, int nb, FILE *fich)
Writes integer(s) into a binary file according to the big endian convention.
Definition: fwrite_be.C:73
Tensor * p_derive_con[N_MET_MAX]
Array of pointers on the contravariant derivatives of this with respect to various metrics...
Definition: tensor.h:355
Cmp pow(const Cmp &, int)
Power .
Definition: cmp_math.C:351
Tensor handling.
Definition: tensor.h:300
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:1012
const Tensor & divergence(const Metric &gam) const
Computes the divergence of this with respect to some metric .
Definition: tensor.C:1065
virtual void del_deriv() const
Deletes the derived quantities.
Definition: tensor.C:408
virtual void inc_dzpuis(int inc=1)
Increases by inc units the value of dzpuis and changes accordingly the values in the compactified ext...
Definition: tensor.C:826
virtual void spectral_display(const char *comment=0x0, double threshold=1.e-7, int precision=4, ostream &ostr=cout) const
Displays the spectral coefficients and the associated basis functions of each component.
Definition: tensor.C:884
void set_spectral_base(const Base_val &)
Sets the spectral bases of the Valeur va
Definition: scalar.C:803
int fread_be(int *aa, int size, int nb, FILE *fich)
Reads integer(s) from a binary file according to the big endian convention.
Definition: fread_be.C:72
virtual ~Tensor()
Destructor.
Definition: tensor.C:395
void copy_coefs_from_large_grid(const Scalar &scal_large_grid)
Copies the content of the argument Scalar to this defined on a smaller grid, with similar mappings...
Definition: scalar_manip.C:590
int valence
Valence of the tensor (0 = scalar, 1 = vector, etc...)
Definition: tensor.h:310
int get_nr(int l) const
Returns the number of points in the radial direction ( ) in domain no. l.
Definition: grilles.h:471
Tensor * p_divergence[N_MET_MAX]
Array of pointers on the divergence of this with respect to various metrics.
Definition: tensor.h:363
void operator+=(const Tensor &)
+= Tensor
Definition: tensor.C:581
const Tensor & derive_con(const Metric &gam) const
Returns the "contravariant" derivative of this with respect to some metric , by raising the last inde...
Definition: tensor.C:1024
Bases of the spectral expansions.
Definition: base_val.h:325
virtual void exponential_filter_ylm(int lzmin, int lzmax, int p, double alpha=-16.)
Applies exponential filters to all components (see Scalar::exponential_filter_ylm )...
Definition: tensor.C:1089
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:818
void set_der_met_0x0(int) const
Sets all the i-th components of met_depend , p_derive_cov , etc...
Definition: tensor.C:443
const Metric * met_depend[N_MET_MAX]
Array on the Metric &#39;s which were used to compute derived quantities, like p_derive_cov ...
Definition: tensor.h:339
Scalar & set(const Itbl &ind)
Returns the value of a component (read/write version).
Definition: tensor.C:664
void annule_extern_cn(int l_0, int deg)
Performs a smooth (C^n) transition in a given domain to zero.
Definition: tensor.C:700
const Tensor * tensor_depend[N_TENSOR_DEPEND]
Pointer on the dependancies, that means the array contains pointers on all the Tensor whom derivative...
Definition: metric.h:139
int get_dim(int i) const
Gives the i th dimension (ie {tt dim.dim[i] )
Definition: itbl.h:326
const Scalar & operator()(const Itbl &ind) const
Returns the value of a component (read-only version).
Definition: tensor.C:808
void set_triad(const Base_vect &new_triad)
Assigns a new vectorial basis (triad) of decomposition.
Definition: tensor.C:529
virtual void set_etat_zero()
Sets the logical state of all components to ETATZERO (zero state).
Definition: tensor.C:507
Basic array class.
Definition: tbl.h:164
Symmetric tensors (with respect to two of their arguments).
Definition: tensor.h:1078
void set_der_0x0() const
Sets the pointers on derived quantities to 0x0.
Definition: tensor.C:417
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
virtual void annule(int l_min, int l_max)
Sets the Tensor to zero in several domains.
Definition: tensor.C:681
virtual void std_spectral_base()
Sets the standard spectal bases of decomposition for each component.
Definition: tensor.C:936
virtual Itbl indices(int pos) const
Returns the indices of a component given by its position in the array cmp .
Definition: tensor.C:549
Coord r
r coordinate centered on the grid
Definition: map.h:745
Tensor(const Map &map, int val, const Itbl &tipe, const Base_vect &triad_i)
Standard constructor.
Definition: tensor.C:212
virtual Tensor * p_divergence(const Tensor &tens) const
Computes the divergence of a tensor (with respect to the current connection).
Definition: connection.C:466