LORENE
valeur.C
1/*
2 * Methods of class Valeur
3 *
4 * (see file valeur.h for documentation)
5 *
6 */
7
8/*
9 * Copyright (c) 1999-2000 Jean-Alain Marck
10 * Copyright (c) 1999-2003 Eric Gourgoulhon
11 * Copyright (c) 1999-2001 Philippe Grandclement
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: valeur.C,v 1.9 2016/12/05 16:18:20 j_novak Exp $
36 * $Log: valeur.C,v $
37 * Revision 1.9 2016/12/05 16:18:20 j_novak
38 * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
39 *
40 * Revision 1.8 2014/10/13 08:53:49 j_novak
41 * Lorene classes and functions now belong to the namespace Lorene.
42 *
43 * Revision 1.7 2014/10/06 15:13:22 j_novak
44 * Modified #include directives to use c++ syntax.
45 *
46 * Revision 1.6 2005/10/25 08:56:40 p_grandclement
47 * addition of std_spectral_base in the case of odd functions near the origin
48 *
49 * Revision 1.5 2004/11/23 12:45:00 f_limousin
50 * Add function filtre_tp(int nn, int nz1, int nz2).
51 *
52 * Revision 1.4 2003/10/19 19:52:56 e_gourgoulhon
53 * Added new method display_coef.
54 *
55 * Revision 1.3 2002/10/16 14:37:15 j_novak
56 * Reorganization of #include instructions of standard C++, in order to
57 * use experimental version 3 of gcc.
58 *
59 * Revision 1.2 2001/12/04 21:27:54 e_gourgoulhon
60 *
61 * All writing/reading to a binary file are now performed according to
62 * the big endian convention, whatever the system is big endian or
63 * small endian, thanks to the functions fwrite_be and fread_be
64 *
65 * Revision 1.1.1.1 2001/11/20 15:19:27 e_gourgoulhon
66 * LORENE
67 *
68 * Revision 2.37 2000/09/11 13:53:07 eric
69 * Ajout des membres p_mult_cp et p_mult_sp.
70 *
71 * Revision 2.36 2000/09/08 10:07:13 eric
72 * Ajout des methodes set_base_r, etc...
73 *
74 * Revision 2.35 1999/12/29 13:11:23 eric
75 * Ajout de la fonction val_point_jk.
76 *
77 * Revision 2.34 1999/12/20 16:35:23 eric
78 * Ajout de la fonction set_base.
79 *
80 * Revision 2.33 1999/12/10 16:09:47 eric
81 * Annulation des membres derives dans la fonction annule(int,int).
82 *
83 * Revision 2.32 1999/12/07 14:53:00 eric
84 * Changement ordre des arguments (phi,theta,xi) --> (xi,theta,phi)
85 * dans la routine val_point.
86 *
87 * Revision 2.31 1999/12/06 16:47:48 eric
88 * Ajout de la fonction val_point.
89 *
90 * Revision 2.30 1999/11/30 12:41:53 eric
91 * Le membre base est desormais un objet de type Base_val et non plus
92 * un pointeur vers une Base_val.
93 *
94 * Revision 2.29 1999/11/23 16:19:33 eric
95 * Suppression du constructeur par defaut.
96 *
97 * Revision 2.28 1999/11/23 14:32:16 novak
98 * Ajout des membres mult_ct et mult_st
99 *
100 * Revision 2.27 1999/11/22 15:41:33 eric
101 * Ajout de la fonction annule(int l).
102 *
103 * Revision 2.26 1999/11/19 11:51:01 eric
104 * Ajout du membre p_stdsdp.
105 *
106 * Revision 2.25 1999/11/19 10:53:17 eric
107 * *** empty log message ***
108 *
109 * Revision 2.24 1999/11/19 10:34:52 eric
110 * Corrections de diverses erreurs dans l'affectation (operator=):
111 * 1/ l'auto-affectation est desormais interdite
112 * 2/ l'affectation a des elements derives est desormais possible
113 *
114 * Revision 2.23 1999/11/19 09:32:15 eric
115 * Ajout du membre p_lapang.
116 *
117 * Revision 2.22 1999/11/16 13:28:32 novak
118 * Ajout de la gestion des pointeurs sur mult_x et scost
119 *
120 * Revision 2.21 1999/10/29 15:14:59 eric
121 * *** empty log message ***
122 *
123 * Revision 2.20 1999/10/28 13:24:15 phil
124 * copie passe avec ETATNONDEF
125 *
126 * Revision 2.19 1999/10/27 09:54:34 eric
127 * *** empty log message ***
128 *
129 * Revision 2.18 1999/10/22 08:42:53 eric
130 * *** empty log message ***
131 *
132 * Revision 2.17 1999/10/22 08:32:55 eric
133 * Anglicisation de operator<<
134 *
135 * Revision 2.16 1999/10/21 14:33:05 eric
136 * *** empty log message ***
137 *
138 * Revision 2.15 1999/10/21 14:21:23 eric
139 * Constructeur par lecture de fichier.
140 * Fonction sauve().
141 *
142 * Revision 2.14 1999/10/20 15:32:20 eric
143 * Ajout de la routine Valeur::std_base_scal().
144 * Ajout de operator=(const Mtbl_cf&).
145 *
146 * Revision 2.13 1999/10/19 15:30:33 eric
147 * Ajout de la fonction affiche_seuil.
148 *
149 * Revision 2.12 1999/10/18 15:16:17 eric
150 * *** empty log message ***
151 *
152 * Revision 2.11 1999/10/18 15:09:01 eric
153 * La fonction membre annule() est rebaptisee annule_hard().
154 * Introduction de la fonction membre annule(int, int).
155 *
156 * Revision 2.10 1999/10/18 13:43:29 eric
157 * Suppression de sxdsdx() et p_sxdsdx (car non-implementes).
158 *
159 * Revision 2.9 1999/10/13 15:52:35 eric
160 * Depoussierage.
161 *
162 * Revision 2.8 1999/04/12 13:05:40 phil
163 * *** empty log message ***
164 *
165 * Revision 2.7 1999/04/12 12:47:56 phil
166 * Correction du constructeur par recopie.
167 *
168 * Revision 2.6 1999/04/12 12:13:26 phil
169 * Correction de l'operateur Valeur = Valeur
170 *
171 * Revision 2.5 1999/03/01 15:10:42 eric
172 * *** empty log message ***
173 *
174 * Revision 2.4 1999/02/26 11:43:07 hyc
175 * *** empty log message ***
176 *
177 * Revision 2.3 1999/02/24 15:26:31 hyc
178 * *** empty log message ***
179 *
180 * Revision 2.2 1999/02/23 11:46:40 hyc
181 * *** empty log message ***
182 *
183 * $Header: /cvsroot/Lorene/C++/Source/Valeur/valeur.C,v 1.9 2016/12/05 16:18:20 j_novak Exp $
184 *
185 */
186
187// headers C
188#include <cassert>
189#include <cstdlib>
190
191// headers Lorene
192#include "valeur.h"
193#include "utilitaires.h"
194#include "proto.h"
195
196
197
198 //---------------//
199 // Constructeurs //
200 //---------------//
201
202namespace Lorene {
203Valeur::Valeur(const Mg3d& mgi) : mg(&mgi), base(mgi.get_nzone()) {
204
205 // C'est nouveau
206 nouveau() ;
207}
208
209Valeur::Valeur(const Mg3d* mgi) : mg(mgi), base(mgi->get_nzone()) {
210
211 // C'est nouveau
212 nouveau() ;
213}
214
215// Copie
216Valeur::Valeur(const Valeur& vc) : base(vc.base) {
217
218 // Tout par default
219 mg = vc.get_mg() ;
220 nouveau() ;
221
222 // Que faire ?
223 switch(vc.get_etat()) {
224 case ETATZERO:
225 set_etat_zero() ;
226 break ;
227
228 case ETATNONDEF:
230 break ;
231
232 case ETATQCQ:
233 assert((vc.c != 0x0) || (vc.c_cf != 0x0) ) ;
234 del_deriv() ;
235
236 if (vc.c != 0x0) {
237 if (c == 0x0) {
238 c = new Mtbl( *(vc.c) ) ;
239 }
240 else{
241 *c = *(vc.c) ;
242 }
243 }
244
245 if (vc.c_cf != 0x0) {
246 if (c_cf == 0x0) {
247 c_cf = new Mtbl_cf( *(vc.c_cf) ) ;
248 }
249 else{
250 *c_cf = *(vc.c_cf) ;
251 }
252 }
253
254 etat = ETATQCQ ;
255 break ;
256
257 default:
258 cout << "Etat pas possible" << endl ;
259 abort() ;
260 break ;
261 }
262}
263
264// Constructeur a partir d'une grille et d'un fichier
265// --------------------------------------------------
266Valeur::Valeur(const Mg3d& g, FILE* fd) : mg(&g), base(fd) {
267
268 // La multi-grille
269 Mg3d* mg_tmp = new Mg3d(fd) ; // la multi-grille d'origine
270 if (*mg != *mg_tmp) {
271 cout <<
272 "Valeur::Valeur(const Mg3d& g, FILE* fd): grid not consistent !"
273 << endl ;
274 abort() ;
275 }
276 delete mg_tmp ;
277
278 fread_be(&etat, sizeof(int), 1, fd) ; // L'etat
279
280 if (etat == ETATQCQ) {
281 c = new Mtbl(g, fd) ; // Les valeurs ds l'espace des conf.
282
283 // Tous les autres pointeurs sont mis a zero :
284 c_cf = 0x0 ;
285 set_der_0x0() ;
286 }
287 else {
288 c = 0x0 ;
289 c_cf = 0x0 ;
290 set_der_0x0() ;
291 }
292
293}
294
295 //--------------//
296 // Destructeurs //
297 //--------------//
298
299// Destructeur
301 del_t() ;
302}
303
304 //-------------//
305 // Affectation //
306 //-------------//
307
308// From double
309// -----------
310void Valeur::operator=(double x) {
311
312 // Cas x = 0
313 if (x == 0) {
314 set_etat_zero() ;
315 return ;
316 }
317
318 // Cas general
319 // les dependances
321
322 // Les bases
323 base.set_base_nondef() ;
324
325 *c = x ;
326}
327
328// From Valeur
329// -----------
330void Valeur::operator=(const Valeur& v) {
331
332 // Protections:
333 // -----------
334 assert(&v != this) ; // pour eviter l'auto-affectation
335 assert(mg == v.mg) ; // meme grille
336
337 // Copie de v :
338 // ------------
339 etat = v.etat ; // l'etat
340
341 base = v.base ; // Les bases
342
343 delete c ;
344 c = 0x0 ; // eventuellement provisoire...
345
346 delete c_cf ;
347 c_cf = 0x0 ; // eventuellement provisoire...
348
349 // La valeur eventuelle
350 switch(v.etat) {
351 case ETATNONDEF: {
352 set_etat_nondef() ;
353 break ; // valeur par defaut
354 }
355
356 case ETATZERO: {
357 set_etat_zero() ;
358 break ;
359 }
360
361 case ETATQCQ: {
362 if (v.c != 0x0) {
363 c = new Mtbl( *(v.c) ) ;
364 }
365
366 if (v.c_cf != 0x0) {
367 c_cf = new Mtbl_cf( *(v.c_cf) ) ;
368 }
369
370 etat = ETATQCQ ;
371
372 // On detruit les dependances (seulement lorsque tout est fini !)
373 del_deriv() ;
374
375 break ;
376 }
377
378 default: {
379 cout << "Unkwon state in Valeur::operator=(const Valeur&) !"
380 << endl ;
381 abort() ;
382 break ;
383 }
384 }
385
386}
387
388// From Mtbl
389// ---------
390void Valeur::operator=(const Mtbl& mt) {
391
392 // Protections
393 // -----------
394 assert(mt.get_etat() != ETATNONDEF) ;
395 assert(mg == mt.get_mg()) ; // meme grille
396 assert(&mt != c) ; // pour eviter l'autoaffectation
397
398 // Les bases
399 base.set_base_nondef() ;
400
401 delete c_cf ;
402 c_cf = 0x0 ;
403
404 // Suivant le cas...
405 switch(mt.get_etat()) {
406 case ETATZERO: {
407 set_etat_zero() ;
408 break ;
409 }
410
411 case ETATQCQ: {
412 delete c ;
413 c = new Mtbl(mt) ;
414
415 del_deriv() ; // On detruit les dependances...
416
417 etat = ETATQCQ ;
418 break ;
419 }
420
421 default: {
422 cout << "Unkwon state in Valeur::operator=(const Mtbl&) !"
423 << endl ;
424 abort() ;
425 break ;
426 }
427 }
428
429}
430
431// From Mtbl_cf
432// ------------
433void Valeur::operator=(const Mtbl_cf& mt) {
434
435 // Protections
436 // -----------
437 assert(mt.get_etat() != ETATNONDEF) ;
438 assert(mg == mt.get_mg()) ; // meme grille
439 assert(&mt != c_cf) ; // pour eviter l'autoaffectation
440
441 // Les bases
442 base = mt.base ;
443
444 delete c ;
445 c = 0x0 ;
446
447 // Suivant le cas...
448 switch(mt.get_etat()) {
449 case ETATZERO: {
450 set_etat_zero() ;
451 break ;
452 }
453
454 case ETATQCQ: {
455 delete c_cf ;
456 c_cf = new Mtbl_cf(mt) ;
457
458 del_deriv() ; // On detruit les dependances...
459
460 etat = ETATQCQ ;
461 break ;
462 }
463
464 default: {
465 cout << "Unkwon state in Valeur::operator=(const Mtbl_cf&) !"
466 << endl ;
467 abort() ;
468 break ;
469 }
470 }
471
472}
473
474 //------------//
475 // Sauvegarde //
476 //------------//
477
478void Valeur::sauve(FILE* fd) const {
479
480 base.sauve(fd) ; // la base
481 mg->sauve(fd) ; // la multi-grille
482 fwrite_be(&etat, sizeof(int), 1, fd) ; // l'etat
483
484 if (etat == ETATQCQ) {
485 if (c == 0x0) { // La sauvegarde s'effectue dans l'espace
486 coef_i() ; // des configurations
487 }
488 c->sauve(fd) ;
489 }
490
491}
492
493 //------------//
494 // Impression //
495 //------------//
496
497// Operator <<
498// -----------
499ostream& operator<<(ostream& o, const Valeur & vi) {
500
501 switch(vi.etat) {
502 case ETATNONDEF: {
503 o << "*** Valeur in UNDEFINED STATE" ;
504 break ;
505 }
506
507 case ETATZERO: {
508 o << "*** Valeur IDENTICALLY ZERO" ;
509 break ;
510 }
511
512 case ETATQCQ: {
513 if (vi.c != 0x0) {
514 o << "*** Valeur (configuration space) :" << endl ;
515 o << *(vi.c) << endl ;
516 }
517 if (vi.c_cf != 0x0) {
518 o << "*** Valeur (coefficients) :" << endl ;
519 o << *(vi.c_cf) << endl ;
520 }
521 break ;
522 }
523
524 default: {
525 cout << "operator<<(ostream&, const Valeur&) : unknown state !"
526 << endl ;
527 abort() ;
528 break ;
529 }
530
531 }
532
533 // Termine
534 return o ;
535}
536
537// display_coef
538// ------------
539void Valeur::display_coef(double thres, int precis, ostream& ost) const {
540
541 if (etat == ETATNONDEF) {
542 ost << " state: UNDEFINED" << endl ;
543 return ;
544 }
545
546 if (etat == ETATZERO) {
547 ost << " state: ZERO" << endl ;
548 return ;
549 }
550
551 coef() ; // the coefficients are required
552
553 c_cf->display(thres, precis, ost) ;
554
555}
556
557
558// affiche_seuil
559//---------------
560
561void Valeur::affiche_seuil(ostream& ost, int type, int precis,
562 double seuil) const {
563 ost << "*** Valeur " << endl ;
564
565
566 // Cas particuliers
567 //-----------------
568
569 if (etat == ETATNONDEF) {
570 ost << " state: UNDEFINED" << endl ;
571 return ;
572 }
573
574 if (etat == ETATZERO) {
575 ost << " state: ZERO" << endl ;
576 return ;
577 }
578
579 switch(type) {
580 case 0 : {
581 coef() ; // Mise a jour eventuelle des coefficients
582 ost << " Coefficients of the Valeur : " << endl ;
583 c_cf->affiche_seuil(ost, precis, seuil) ;
584 break ;
585 }
586
587 case 1 : {
588 coef_i() ; // Mise a jour eventuelle dans l'espace des conf.
589 ost << " Values of the Valeur at the collocation points: " << endl ;
590 c->affiche_seuil(ost, precis, seuil) ;
591 break ;
592 }
593
594 case 2 : {
595 coef() ; // Mise a jour eventuelle des coefficients
596 coef_i() ; // Mise a jour eventuelle dans l'espace des conf.
597 ost << " Coefficients of the Valeur : " << endl ;
598 c_cf->affiche_seuil(ost, precis, seuil) ;
599 ost << " Values of the Valeur at the collocation points: " << endl ;
600 c->affiche_seuil(ost, precis, seuil) ;
601 break ;
602 }
603
604 default : {
605 cout << "Unknown type in Valeur::affiche_seuil !" << endl ;
606 abort() ;
607 break ;
608 }
609 }
610
611
612}
613
614
615
616 //-----------------------//
617 // Gestion de la memoire //
618 //-----------------------//
619
621 // Les donnees
622 c = 0x0 ;
623 c_cf = 0x0 ;
624 set_der_0x0() ;
625 // L'etat
626 etat = ETATNONDEF ;
627}
628
630
631 delete c ;
632 c = 0x0 ;
633
634 delete c_cf ;
635 c_cf = 0x0 ;
636
637 del_deriv() ;
638
639}
640
642 // Destruction
643 delete p_dsdx ;
644 delete p_d2sdx2 ;
645 delete p_sx ;
646 delete p_sx2 ;
647 delete p_mult_x ;
648
649 delete p_dsdt ;
650 delete p_d2sdt2 ;
651 delete p_ssint ;
652 delete p_scost ;
653 delete p_mult_ct ;
654 delete p_mult_st ;
655
656 delete p_dsdp ;
657 delete p_stdsdp ;
658 delete p_d2sdp2 ;
659 delete p_mult_cp ;
660 delete p_mult_sp ;
661
662 delete p_lapang ;
663
664 // Pointeurs a 0x0
665 set_der_0x0() ;
666}
667
669 p_dsdx = 0x0 ;
670 p_d2sdx2 = 0x0 ;
671 p_sx = 0x0 ;
672 p_sx2 = 0x0 ;
673 p_mult_x = 0x0 ;
674
675 p_dsdt = 0x0 ;
676 p_d2sdt2 = 0x0 ;
677 p_ssint = 0x0 ;
678 p_scost = 0x0 ;
679 p_mult_ct = 0x0 ;
680 p_mult_st = 0x0 ;
681
682 p_dsdp = 0x0 ;
683 p_stdsdp = 0x0 ;
684 p_d2sdp2 = 0x0 ;
685 p_mult_cp = 0x0 ;
686 p_mult_sp = 0x0 ;
687
688 p_lapang = 0x0 ;
689}
690
691// ETATZERO
693 if (etat == ETATZERO) return ;
694 del_t() ;
695 etat = ETATZERO ;
696}
697// ETATNONDEF
699 if (etat == ETATNONDEF) return ;
700 del_t() ;
701 etat = ETATNONDEF ;
702}
703// ETATQCQ physique
705 // Detruit les dependances
706 del_deriv() ;
707 delete c_cf ; c_cf = 0x0 ;
708
709 if (c == 0x0) {
710 c = new Mtbl(mg) ;
711 }
712 etat = ETATQCQ ;
713}
714// ETATQCQ coefficients
716 // Detruit les dependances
717 del_deriv() ;
718 delete c ; c = 0x0 ;
719
720 if (c_cf == 0x0) {
721 c_cf = new Mtbl_cf(mg, base) ;
722 }
723 etat = ETATQCQ ;
724}
725// ZERO hard
727 // Detruit les dependances
728 del_deriv() ;
729
730 if (c == 0x0) {
731 c = new Mtbl(mg) ;
732 }
733 if (c_cf == 0x0) {
734 c_cf = new Mtbl_cf(mg, base) ;
735 }
736
737 c->annule_hard() ;
738 c_cf->annule_hard() ;
739
740 etat = ETATQCQ ;
741}
742
743
744// Sets the Valeur to zero in a given domain
745// -----------------------------------------
746
747void Valeur::annule(int l) {
748
749 annule(l, l) ;
750}
751
752
753
754// Sets the Valeur to zero in several domains
755// ------------------------------------------
756
757void Valeur::annule(int l_min, int l_max) {
758
759 // Cas particulier: annulation globale :
760 if ( (l_min == 0) && (l_max == mg->get_nzone()-1) ) {
761 set_etat_zero() ;
762 return ;
763 }
764
765 assert( etat != ETATNONDEF ) ;
766
767 if ( etat == ETATZERO ) {
768 return ; // rien n'a faire si c'est deja zero
769 }
770 else {
771 assert( etat == ETATQCQ ) ; // sinon...
772
773 if (c != 0x0) {
774 c->annule(l_min, l_max) ; // Annule le Mtbl
775 }
776
777 if (c_cf != 0x0) {
778 c_cf->annule(l_min, l_max) ; // Annule le Mtbl_cf
779 }
780
781 // Annulation des membres derives
782
783 if (p_dsdx != 0x0) p_dsdx->annule(l_min, l_max) ;
784 if (p_d2sdx2 != 0x0) p_d2sdx2->annule(l_min, l_max) ;
785 if (p_sx != 0x0) p_sx->annule(l_min, l_max) ;
786 if (p_sx2 != 0x0) p_sx2->annule(l_min, l_max) ;
787 if (p_mult_x != 0x0) p_mult_x->annule(l_min, l_max) ;
788
789 if (p_dsdt != 0x0) p_dsdt->annule(l_min, l_max) ;
790 if (p_d2sdt2 != 0x0) p_d2sdt2->annule(l_min, l_max) ;
791 if (p_ssint != 0x0) p_ssint->annule(l_min, l_max) ;
792 if (p_scost != 0x0) p_scost->annule(l_min, l_max) ;
793 if (p_mult_ct != 0x0) p_mult_ct->annule(l_min, l_max) ;
794 if (p_mult_st != 0x0) p_mult_st->annule(l_min, l_max) ;
795
796 if (p_dsdp != 0x0) p_dsdp->annule(l_min, l_max) ;
797 if (p_stdsdp != 0x0) p_stdsdp->annule(l_min, l_max) ;
798 if (p_d2sdp2 != 0x0) p_d2sdp2->annule(l_min, l_max) ;
799 if (p_mult_cp != 0x0) p_mult_cp->annule(l_min, l_max) ;
800 if (p_mult_sp != 0x0) p_mult_sp->annule(l_min, l_max) ;
801
802 if (p_lapang != 0x0) p_lapang->annule(l_min, l_max) ;
803
804 }
805
806}
807
808
809 //--------------------------------------//
810 // Spectral bases manipulation //
811 //--------------------------------------//
812
813void Valeur::set_base(const Base_val& base_i) {
814
815 base = base_i ;
816
817 if (c_cf != 0x0) {
818 if ( !(c_cf->base == base_i) ) {
819 delete c_cf ;
820 c_cf = 0x0 ;
821 }
822 }
823
824}
825
826
828
829 set_base( mg->std_base_scal() ) ;
830
831}
832
834
835 set_base( mg->std_base_scal_odd() ) ;
836
837}
838
839void Valeur::set_base_r(int l, int base_r) {
840
841 base.set_base_r(l, base_r) ;
842
843 if (c_cf != 0x0) {
844 if ( !(c_cf->base == base) ) {
845 delete c_cf ;
846 c_cf = 0x0 ;
847 }
848 }
849
850}
851
852void Valeur::set_base_t(int base_t) {
853
854 base.set_base_t(base_t) ;
855
856 if (c_cf != 0x0) {
857 if ( !(c_cf->base == base) ) {
858 delete c_cf ;
859 c_cf = 0x0 ;
860 }
861 }
862
863}
864
865void Valeur::set_base_p(int base_p) {
866
867 base.set_base_p(base_p) ;
868
869 if (c_cf != 0x0) {
870 if ( !(c_cf->base == base) ) {
871 delete c_cf ;
872 c_cf = 0x0 ;
873 }
874 }
875
876}
877
878
879
880
881 //-----------------------------------------------//
882 // Value at an arbitrary point //
883 //-----------------------------------------------//
884
885double Valeur::val_point(int l, double x, double theta, double phi) const {
886
887 assert(etat != ETATNONDEF) ;
888
889 double resu ;
890
891 if (etat == ETATZERO) {
892 resu = 0 ;
893 }
894 else{
895 assert(etat == ETATQCQ) ;
896 coef() ; // The coefficients are required
897 resu = c_cf->val_point(l, x, theta, phi) ; // Call to the Mtbl_cf version
898 }
899
900 return resu ;
901}
902
903double Valeur::val_point_jk(int l, double x, int j, int k) const {
904
905 assert(etat != ETATNONDEF) ;
906
907 double resu ;
908
909 if (etat == ETATZERO) {
910 resu = 0 ;
911 }
912 else{
913 assert(etat == ETATQCQ) ;
914 coef() ; // The coefficients are required
915 resu = c_cf->val_point_jk(l, x, j, k) ; // Call to the Mtbl_cf version
916 }
917
918 return resu ;
919}
920
921
922 //-----------------------------------------------//
923 // Filtres //
924 //-----------------------------------------------//
925
926
927void Valeur::filtre_tp(int nn, int nz1, int nz2) {
928
929 int nz = mg->get_nzone() ;
930 int nr = mg->get_nr(0) ;
931 int nt = mg->get_nt(0) ;
932 int np = mg->get_np(0) ;
933
934 if (etat != ETATZERO) {
935 assert (etat == ETATQCQ) ;
936 ylm() ;
937 for (int lz=nz1; lz<=nz2; lz++) {
938 if (c_cf->operator()(lz).get_etat() != ETATZERO)
939 for (int k=0; k<np+1; k++)
940 for (int j=0; j<nt; j++) {
941 int l_q, m_q, base_r ;
942 if (nullite_plm(j, nt, k, np, base) == 1) {
943 donne_lm(nz, lz, j, k, base, m_q, l_q,base_r) ;
944 if (l_q > nn)
945 for (int i=0; i<nr; i++)
946 c_cf->set(lz, k, j, i) = 0. ;
947 }
948 }
949 }
950 if (c != 0x0) {
951 delete c ;
952 c = 0x0 ;
953 }
954 }
955
956}
957}
Bases of the spectral expansions.
Definition base_val.h:325
Multi-domain grid.
Definition grilles.h:279
Coefficients storage for the multi-domain spectral method.
Definition mtbl_cf.h:196
Base_val base
Bases of the spectral expansions.
Definition mtbl_cf.h:210
int get_etat() const
Returns the logical state.
Definition mtbl_cf.h:466
const Mg3d * get_mg() const
Returns the Mg3d on which the Mtbl_cf is defined.
Definition mtbl_cf.h:463
Multi-domain array.
Definition mtbl.h:118
const Mg3d * get_mg() const
Gives the Mg3d on which the Mtbl is defined.
Definition mtbl.h:274
int get_etat() const
Gives the logical state.
Definition mtbl.h:277
void filtre_tp(int nn, int nz1, int nz2)
Sets the n lasts coefficients in to 0 in the domain nz1 to nz2 when expressed in spherical harmonics...
Definition valeur.C:927
void del_t()
Logical destructor.
Definition valeur.C:629
void set_etat_cf_qcq()
Sets the logical state to ETATQCQ (ordinary state) for values in the configuration space (Mtbl_cf c_c...
Definition valeur.C:715
double val_point(int l, double x, double theta, double phi) const
Computes the value of the field represented by *this at an arbitrary point, by means of the spectral ...
Definition valeur.C:885
~Valeur()
Destructor.
Definition valeur.C:300
void set_etat_c_qcq()
Sets the logical state to ETATQCQ (ordinary state) for values in the configuration space (Mtbl c ).
Definition valeur.C:704
double val_point_jk(int l, double x, int j, int k) const
Computes the value of the field represented by *this at an arbitrary point in , but collocation point...
Definition valeur.C:903
Valeur * p_dsdp
Pointer on .
Definition valeur.h:333
friend ostream & operator<<(ostream &, const Valeur &)
Display.
Definition valeur.C:499
Valeur * p_mult_sp
Pointer on .
Definition valeur.h:337
Valeur * p_mult_x
Pointer on .
Definition valeur.h:324
Valeur * p_scost
Pointer on .
Definition valeur.h:329
int get_etat() const
Returns the logical state.
Definition valeur.h:760
Valeur * p_sx2
Pointer on .
Definition valeur.h:323
Valeur * p_d2sdx2
Pointer on .
Definition valeur.h:321
void set_etat_zero()
Sets the logical state to ETATZERO (zero).
Definition valeur.C:692
void set_base(const Base_val &)
Sets the bases for spectral expansions (member base ).
Definition valeur.C:813
Valeur * p_dsdx
Pointer on .
Definition valeur.h:320
Valeur(const Mg3d &mgrid)
Constructor.
Definition valeur.C:203
void nouveau()
Memory allocation.
Definition valeur.C:620
void ylm()
Computes the coefficients of *this.
Definition valeur_ylm.C:141
Valeur * p_d2sdt2
Pointer on .
Definition valeur.h:327
Mtbl * c
Values of the function at the points of the multi-grid.
Definition valeur.h:309
Valeur * p_mult_cp
Pointer on .
Definition valeur.h:336
void display_coef(double threshold=1.e-7, int precision=4, ostream &ostr=cout) const
Displays the spectral coefficients and the associated basis functions.
Definition valeur.C:539
void annule(int l)
Sets the Valeur to zero in a given domain.
Definition valeur.C:747
const Mg3d * mg
Multi-grid Mgd3 on which this is defined.
Definition valeur.h:302
void operator=(const Valeur &a)
Assignement to another Valeur.
Definition valeur.C:330
void std_base_scal_odd()
Sets the bases for spectral expansions (member base ) to the standard odd ones for a scalar.
Definition valeur.C:833
Valeur * p_sx
Pointer on .
Definition valeur.h:322
void coef_i() const
Computes the physical value of *this.
void affiche_seuil(ostream &ostr, int type=0, int precision=4, double threshold=1.e-7) const
Prints only the values greater than a given threshold.
Definition valeur.C:561
Mtbl_cf * c_cf
Coefficients of the spectral expansion of the function.
Definition valeur.h:312
void coef() const
Computes the coeffcients of *this.
Valeur * p_ssint
Pointer on .
Definition valeur.h:328
void set_der_0x0()
Sets the pointers for derivatives to 0x0.
Definition valeur.C:668
Valeur * p_stdsdp
Pointer on .
Definition valeur.h:334
const Mg3d * get_mg() const
Returns the Mg3d on which the this is defined.
Definition valeur.h:763
void set_base_r(int l, int base_r)
Sets the expansion basis for r ( ) functions in a given domain.
Definition valeur.C:839
Valeur * p_dsdt
Pointer on .
Definition valeur.h:326
Valeur * p_lapang
Pointer on the angular Laplacian.
Definition valeur.h:339
Base_val base
Bases on which the spectral expansion is performed.
Definition valeur.h:315
void std_base_scal()
Sets the bases for spectral expansions (member base ) to the standard ones for a scalar.
Definition valeur.C:827
void annule_hard()
Sets the Valeur to zero in a hard way.
Definition valeur.C:726
Valeur * p_mult_ct
Pointer on .
Definition valeur.h:330
int etat
Logical state (ETATNONDEF , ETATQCQ or ETATZERO ).
Definition valeur.h:305
void set_base_t(int base_t)
Sets the expansion basis for functions in all domains.
Definition valeur.C:852
void set_base_p(int base_p)
Sets the expansion basis for functions in all domains.
Definition valeur.C:865
void set_etat_nondef()
Sets the logical state to ETATNONDEF (undefined).
Definition valeur.C:698
Valeur * p_d2sdp2
Pointer on .
Definition valeur.h:335
Valeur * p_mult_st
Pointer on .
Definition valeur.h:331
void del_deriv()
Logical destructor of the derivatives.
Definition valeur.C:641
void sauve(FILE *) const
Save in a file.
Definition valeur.C:478
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
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
Lorene prototypes.
Definition app_hor.h:67
Coord phi
coordinate centered on the grid
Definition map.h:732
Coord x
x coordinate centered on the grid
Definition map.h:738