LORENE
param_elliptic_2d.C
1 /*
2  * Copyright (c) 2004 Philippe Grandclement
3  *
4  * This file is part of LORENE.
5  *
6  * LORENE is free software; you can redistribute it and/or modify
7  * it under the terms of the GNU General Public License version 2
8  * as published by the Free Software Foundation.
9  *
10  * LORENE is distributed in the hope that it will be useful,
11  * but WITHOUT ANY WARRANTY; without even the implied warranty of
12  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13  * GNU General Public License for more details.
14  *
15  * You should have received a copy of the GNU General Public License
16  * along with LORENE; if not, write to the Free Software
17  * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
18  *
19  */
20 
21 
22 
23 /*
24  * $Id: param_elliptic_2d.C,v 1.5 2018/11/16 14:34:36 j_novak Exp $
25  * $Log: param_elliptic_2d.C,v $
26  * Revision 1.5 2018/11/16 14:34:36 j_novak
27  * Changed minor points to avoid some compilation warnings.
28  *
29  * Revision 1.4 2016/12/05 16:18:14 j_novak
30  * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
31  *
32  * Revision 1.3 2014/10/13 08:53:37 j_novak
33  * Lorene classes and functions now belong to the namespace Lorene.
34  *
35  * Revision 1.2 2014/10/06 15:13:15 j_novak
36  * Modified #include directives to use c++ syntax.
37  *
38  * Revision 1.1 2004/08/24 09:14:49 p_grandclement
39  * Addition of some new operators, like Poisson in 2d... It now requieres the
40  * GSL library to work.
41  *
42  * Also, the way a variable change is stored by a Param_elliptic is changed and
43  * no longer uses Change_var but rather 2 Scalars. The codes using that feature
44  * will requiere some modification. (It should concern only the ones about monopoles)
45  *
46  *
47  * $Header: /cvsroot/Lorene/C++/Source/Param_elliptic/param_elliptic_2d.C,v 1.5 2018/11/16 14:34:36 j_novak Exp $
48  *
49  */
50 
51 #include "headcpp.h"
52 
53 #include <cmath>
54 #include <cstdlib>
55 
56 #include "base_val.h"
57 #include "map.h"
58 #include "ope_elementary.h"
59 #include "param_elliptic.h"
60 #include "change_var.h"
61 #include "scalar.h"
62 
63 
64 namespace Lorene {
65 void Param_elliptic::set_poisson_2d(const Scalar& source, bool indic) {
66 
67  int dzpuis = source.get_dzpuis() ;
68 
69  if (type_map != MAP_AFF) {
70  cout << "set_poisson_2d only defined for an affine mapping..." << endl ;
71  abort() ;
72  }
73  else {
74 
75  int nz = get_mp().get_mg()->get_nzone() ;
76 
77  int nr ;
78  double alpha, beta ;
79  int m_quant, l_quant, base_r_1d ;
80 
81  int conte = 0 ;
82  for (int l=0 ; l<nz ; l++) {
83 
84  nr = get_mp().get_mg()->get_nr(l) ;
85  alpha = get_alpha (l) ;
86  beta = get_beta (l) ;
87 
88  for (int k=0 ; k<get_mp().get_mg()->get_np(l)+1 ; k++)
89  for (int j=0 ; j<get_mp().get_mg()->get_nt(l) ; j++) {
90  if (operateurs[conte] != 0x0)
91  delete operateurs[conte] ;
92  source.get_spectral_va().base.give_quant_numbers(l, k, j, m_quant, l_quant, base_r_1d) ;
93  if (k!=1)
94  if ((indic) || ((!indic) && (l_quant !=0)))
95  operateurs[conte] = new Ope_poisson_2d (nr, base_r_1d, alpha, beta, l_quant, dzpuis) ;
96  else
97  operateurs[conte] = 0x0 ;
98  else
99  operateurs[conte] = 0x0 ;
100  conte ++ ;
101  }
102  }
103  }
104 }
105 
106 void Param_elliptic::set_helmholtz_minus_2d(int zone, double masse, const Scalar& source) {
107 
108  int dzpuis = source.get_dzpuis() ;
109  assert (masse > 0) ;
110 
111  if (type_map != MAP_AFF) {
112  cout << "set_helmholtz_minus_2d only defined for an affine mapping..." << endl ;
113  abort() ;
114  }
115  else {
116 
117  int nz = get_mp().get_mg()->get_nzone() ;
118  if (zone == nz-1)
119  source.check_dzpuis(2) ;
120  int nr ;
121  double alpha, beta ;
122  int m_quant, l_quant, base_r_1d ;
123 
124  int conte = 0 ;
125  for (int l=0 ; l<nz ; l++) {
126 
127  nr = get_mp().get_mg()->get_nr(l) ;
128  alpha = get_alpha (l) ;
129  beta = get_beta (l) ;
130 
131  for (int k=0 ; k<get_mp().get_mg()->get_np(l)+1 ; k++)
132  for (int j=0 ; j<get_mp().get_mg()->get_nt(l) ; j++) {
133  if (l==zone) {
134  if (operateurs[conte] != 0x0) {
135  delete operateurs[conte] ;
137  (l, k, j, m_quant, l_quant, base_r_1d) ;
138  operateurs[conte] = new Ope_helmholtz_minus_2d (nr, base_r_1d, alpha, beta, l_quant, masse, dzpuis) ;
139  }
140  }
141  conte ++ ;
142  }
143  }
144  }
145 }
146 
147 }
int get_np(int l) const
Returns the number of points in the azimuthal direction ( ) in domain no. l.
Definition: grilles.h:479
void give_quant_numbers(int, int, int, int &, int &, int &) const
Computes the various quantum numbers and 1d radial base.
Lorene prototypes.
Definition: app_hor.h:67
const Mg3d * get_mg() const
Gives the Mg3d on which the mapping is defined.
Definition: map.h:777
Tensor field of valence 0 (or component of a tensorial field).
Definition: scalar.h:393
Base_val base
Bases on which the spectral expansion is performed.
Definition: valeur.h:315
int get_dzpuis() const
Returns dzpuis.
Definition: scalar.h:563
int type_map
Type of mapping either MAP_AFF or MAP_LOG.
Class for the operator of the Poisson equation in 2D.
void set_helmholtz_minus_2d(int zone, double mas, const Scalar &)
Set the 2D Helmholtz operator (with minus sign)
const Map_radial & get_mp() const
Returns the mapping.
int get_nzone() const
Returns the number of domains.
Definition: grilles.h:465
Ope_elementary ** operateurs
Array on the elementary operators.
void set_poisson_2d(const Scalar &, bool indic=false)
Set everything to do a 2d-Poisson, with or without l=0 (not put by default...)
int get_nr(int l) const
Returns the number of points in the radial direction ( ) in domain no. l.
Definition: grilles.h:469
Class for the operator of the Helmholtz equation in 2D.
int get_nt(int l) const
Returns the number of points in the co-latitude direction ( ) in domain no. l.
Definition: grilles.h:474
bool check_dzpuis(int dzi) const
Returns false if the last domain is compactified and *this is not zero in this domain and dzpuis is n...
Definition: scalar.C:879
const Valeur & get_spectral_va() const
Returns va (read only version)
Definition: scalar.h:607