LORENE
map_af_elliptic_pseudo_1d.C
1/*
2 * Method of the class Map_af for the resolution of the scalar pseudo_1d
3 * equation
4 */
5
6/*
7 * Copyright (c) 2004 Philippe Grandclement
8 *
9 * This file is part of LORENE.
10 *
11 * LORENE is free software; you can redistribute it and/or modify
12 * it under the terms of the GNU General Public License as published by
13 * the Free Software Foundation; either version 2 of the License, or
14 * (at your option) any later version.
15 *
16 * LORENE is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
19 * GNU General Public License for more details.
20 *
21 * You should have received a copy of the GNU General Public License
22 * along with LORENE; if not, write to the Free Software
23 * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
24 *
25 */
26
27
28char map_af_elliptic_pseudo_1d_C[] = "$Header: /cvsroot/Lorene/C++/Source/Map/map_af_elliptic_pseudo_1d.C,v 1.2 2014/10/13 08:53:02 j_novak Exp $" ;
29
30/*
31 * $Id: map_af_elliptic_pseudo_1d.C,v 1.2 2014/10/13 08:53:02 j_novak Exp $
32 * $Log: map_af_elliptic_pseudo_1d.C,v $
33 * Revision 1.2 2014/10/13 08:53:02 j_novak
34 * Lorene classes and functions now belong to the namespace Lorene.
35 *
36 * Revision 1.1 2004/08/24 09:14:42 p_grandclement
37 * Addition of some new operators, like Poisson in 2d... It now requieres the
38 * GSL library to work.
39 *
40 * Also, the way a variable change is stored by a Param_elliptic is changed and
41 * no longer uses Change_var but rather 2 Scalars. The codes using that feature
42 * will requiere some modification. (It should concern only the ones about monopoles)
43 *
44 *
45 * $Header: /cvsroot/Lorene/C++/Source/Map/map_af_elliptic_pseudo_1d.C,v 1.2 2014/10/13 08:53:02 j_novak Exp $
46 *
47 */
48
49// Header Lorene:
50#include "valeur.h"
51#include "map.h"
52#include "scalar.h"
53#include "param_elliptic.h"
54
55//*****************************************************************************
56
57namespace Lorene {
58
60 const Scalar& source, Scalar& pot) const {
61
62 assert(source.get_etat() != ETATNONDEF) ;
63 assert(source.get_mp().get_mg() == mg) ;
64 assert(pot.get_mp().get_mg() == mg) ;
65
66 assert( source.check_dzpuis(2)) ;
67
68 // Spherical harmonic expansion of the source
69 // ------------------------------------------
70
71 const Valeur& sourva = source.get_spectral_va() ;
72
73 if (sourva.get_etat() == ETATZERO) {
74 pot.set_etat_zero() ;
75 return ;
76 }
77
78 // Spectral coefficients of the source
79 assert(sourva.get_etat() == ETATQCQ) ;
80
81 Valeur rho(sourva.get_mg()) ;
82 sourva.coef() ;
83 rho = *(sourva.c_cf) ; // copy of the coefficients of the source
84 rho.val_propre_1d() ;
85
86 // On calcule les coefs radiaux de F et G
87 ope_var.var_F.set_spectral_va().coef() ;
89 ope_var.var_G.set_spectral_va().coef() ;
90
91
92 // Call to the Mtbl_cf version
93 // ---------------------------
94 Mtbl_cf resu = elliptic_solver (ope_var, *(rho.c_cf)) ;
95
96 pot.set_etat_zero() ;
97 pot.set_etat_qcq() ;
98 pot.set_spectral_va() = resu ;
100 pot.set_dzpuis(0) ;
101
102}
103
104
105}
void sol_elliptic_pseudo_1d(Param_elliptic &, const Scalar &, Scalar &) const
General elliptic solver in a pseudo 1d case.
const Mg3d * get_mg() const
Gives the Mg3d on which the mapping is defined.
Definition map.h:765
const Mg3d * mg
Pointer on the multi-grid Mgd3 on which this is defined
Definition map.h:676
Coefficients storage for the multi-domain spectral method.
Definition mtbl_cf.h:186
This class contains the parameters needed to call the general elliptic solver.
Scalar var_G
Multiplicative variable change that must be sphericaly symetric !
Scalar var_F
Additive variable change function.
Tensor field of valence 0 (or component of a tensorial field).
Definition scalar.h:387
virtual void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition scalar.C:353
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:873
virtual void set_etat_zero()
Sets the logical state to ETATZERO (zero).
Definition scalar.C:324
Valeur & set_spectral_va()
Returns va (read/write version)
Definition scalar.h:604
const Valeur & get_spectral_va() const
Returns va (read only version)
Definition scalar.h:601
int get_etat() const
Returns the logical state ETATNONDEF (undefined), ETATZERO (null) or ETATQCQ (ordinary).
Definition scalar.h:554
void set_dzpuis(int)
Modifies the dzpuis flag.
Definition scalar.C:808
Values and coefficients of a (real-value) function.
Definition valeur.h:287
void val_propre_1d_i()
Inverse transformation of val_propre_1d.
int get_etat() const
Returns the logical state.
Definition valeur.h:726
Mtbl_cf * c_cf
Coefficients of the spectral expansion of the function.
Definition valeur.h:302
void coef() const
Computes the coeffcients of *this.
const Mg3d * get_mg() const
Returns the Mg3d on which the this is defined.
Definition valeur.h:729
void val_propre_1d()
Set the basis to the eigenvalues of .
const Map & get_mp() const
Returns the mapping.
Definition tensor.h:861
Lorene prototypes.
Definition app_hor.h:64