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
21char param_elliptic_2d_C[] = "$Header: /cvsroot/Lorene/C++/Source/Param_elliptic/param_elliptic_2d.C,v 1.3 2014/10/13 08:53:37 j_novak Exp $" ;
22
23/*
24 * $Id: param_elliptic_2d.C,v 1.3 2014/10/13 08:53:37 j_novak Exp $
25 * $Log: param_elliptic_2d.C,v $
26 * Revision 1.3 2014/10/13 08:53:37 j_novak
27 * Lorene classes and functions now belong to the namespace Lorene.
28 *
29 * Revision 1.2 2014/10/06 15:13:15 j_novak
30 * Modified #include directives to use c++ syntax.
31 *
32 * Revision 1.1 2004/08/24 09:14:49 p_grandclement
33 * Addition of some new operators, like Poisson in 2d... It now requieres the
34 * GSL library to work.
35 *
36 * Also, the way a variable change is stored by a Param_elliptic is changed and
37 * no longer uses Change_var but rather 2 Scalars. The codes using that feature
38 * will requiere some modification. (It should concern only the ones about monopoles)
39 *
40 *
41 * $Header: /cvsroot/Lorene/C++/Source/Param_elliptic/param_elliptic_2d.C,v 1.3 2014/10/13 08:53:37 j_novak Exp $
42 *
43 */
44
45#include "headcpp.h"
46
47#include <cmath>
48#include <cstdlib>
49
50#include "base_val.h"
51#include "map.h"
52#include "ope_elementary.h"
53#include "param_elliptic.h"
54#include "change_var.h"
55#include "scalar.h"
56
57
58namespace Lorene {
60
61 int dzpuis = source.get_dzpuis() ;
62
63 if (type_map != MAP_AFF) {
64 cout << "set_poisson_2d only defined for an affine mapping..." << endl ;
65 abort() ;
66 }
67 else {
68
69 int nz = get_mp().get_mg()->get_nzone() ;
70
71 int nr ;
72 double alpha, beta ;
73 int m_quant, l_quant, base_r_1d ;
74
75 int conte = 0 ;
76 for (int l=0 ; l<nz ; l++) {
77
78 nr = get_mp().get_mg()->get_nr(l) ;
79 alpha = get_alpha (l) ;
80 beta = get_beta (l) ;
81
82 for (int k=0 ; k<get_mp().get_mg()->get_np(l)+1 ; k++)
83 for (int j=0 ; j<get_mp().get_mg()->get_nt(l) ; j++) {
84 if (operateurs[conte] != 0x0)
85 delete operateurs[conte] ;
86 source.get_spectral_va().base.give_quant_numbers(l, k, j, m_quant, l_quant, base_r_1d) ;
87 if (k!=1)
88 if ((indic) || ((!indic) && (l_quant !=0)))
89 operateurs[conte] = new Ope_poisson_2d (nr, base_r_1d, alpha, beta, l_quant, dzpuis) ;
90 else
91 operateurs[conte] = 0x0 ;
92 else
93 operateurs[conte] = 0x0 ;
94 conte ++ ;
95 }
96 }
97 }
98}
99
101
102 int dzpuis = source.get_dzpuis() ;
103 assert (masse > 0) ;
104
105 if (type_map != MAP_AFF) {
106 cout << "set_helmholtz_minus_2d only defined for an affine mapping..." << endl ;
107 abort() ;
108 }
109 else {
110
111 int nz = get_mp().get_mg()->get_nzone() ;
112 if (zone == nz-1)
113 source.check_dzpuis(2) ;
114 int nr ;
115 double alpha, beta ;
116 int m_quant, l_quant, base_r_1d ;
117
118 int conte = 0 ;
119 for (int l=0 ; l<nz ; l++) {
120
121 nr = get_mp().get_mg()->get_nr(l) ;
122 alpha = get_alpha (l) ;
123 beta = get_beta (l) ;
124
125 for (int k=0 ; k<get_mp().get_mg()->get_np(l)+1 ; k++)
126 for (int j=0 ; j<get_mp().get_mg()->get_nt(l) ; j++) {
127 if (l==zone) {
128 if (operateurs[conte] != 0x0) {
129 delete operateurs[conte] ;
130 source.get_spectral_va().base.give_quant_numbers
131 (l, k, j, m_quant, l_quant, base_r_1d) ;
132 operateurs[conte] = new Ope_helmholtz_minus_2d (nr, base_r_1d, alpha, beta, l_quant, masse, dzpuis) ;
133 }
134 }
135 conte ++ ;
136 }
137 }
138 }
139}
140
141}
Time evolution with partial storage (*** under development ***).
Definition evolution.h:371
Class for the operator of the Helmholtz equation in 2D.
Class for the operator of the Poisson equation in 2D.
const Map_radial & get_mp() const
Returns the mapping.
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 type_map
Type of mapping either MAP_AFF or MAP_LOG.
void set_helmholtz_minus_2d(int zone, double mas, const Scalar &)
Set the 2D Helmholtz operator (with minus sign)
Tensor field of valence 0 (or component of a tensorial field).
Definition scalar.h:387
Lorene prototypes.
Definition app_hor.h:64