LORENE
map_af_elliptic_2d.C
1/*
2 * Method of the class Map_af for the resolution of the scalar 2D-Poisson
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
28
29
30/*
31 * $Id: map_af_elliptic_2d.C,v 1.3 2016/12/05 16:17:56 j_novak Exp $
32 * $Log: map_af_elliptic_2d.C,v $
33 * Revision 1.3 2016/12/05 16:17:56 j_novak
34 * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
35 *
36 * Revision 1.2 2014/10/13 08:53:02 j_novak
37 * Lorene classes and functions now belong to the namespace Lorene.
38 *
39 * Revision 1.1 2004/08/24 09:14:42 p_grandclement
40 * Addition of some new operators, like Poisson in 2d... It now requieres the
41 * GSL library to work.
42 *
43 * Also, the way a variable change is stored by a Param_elliptic is changed and
44 * no longer uses Change_var but rather 2 Scalars. The codes using that feature
45 * will requiere some modification. (It should concern only the ones about monopoles)
46 *
47 *
48 * $Header: /cvsroot/Lorene/C++/Source/Map/map_af_elliptic_2d.C,v 1.3 2016/12/05 16:17:56 j_novak Exp $
49 *
50 */
51
52// Header Lorene:
53#include "map.h"
54#include "scalar.h"
55#include "param_elliptic.h"
56
57//*****************************************************************************
58
59namespace Lorene {
60
62 const Scalar& source, Scalar& pot) const {
63
64 assert(source.get_etat() != ETATNONDEF) ;
65 assert(source.get_mp().get_mg() == mg) ;
66 assert(pot.get_mp().get_mg() == mg) ;
67
68 assert( source.check_dzpuis(2) || source.check_dzpuis(4)
69 || source.check_dzpuis(3)) ;
70
71 // Spherical harmonic expansion of the source
72 // ------------------------------------------
73
74 const Valeur& sourva = source.get_spectral_va() ;
75
76 if (sourva.get_etat() == ETATZERO) {
77 pot.set_etat_zero() ;
78 return ;
79 }
80
81 // Spectral coefficients of the source
82 assert(sourva.get_etat() == ETATQCQ) ;
83
84 Valeur rho(sourva.get_mg()) ;
85 sourva.coef() ;
86 rho = *(sourva.c_cf) ; // copy of the coefficients of the source
87
88 // On calcule les coefs radiaux de F et G
89 ope_var.var_F.set_spectral_va().coef() ;
90 ope_var.var_G.set_spectral_va().coef() ;
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 ;
99 pot.set_dzpuis(0) ;
100}
101
102
103}
void sol_elliptic_2d(Param_elliptic &, const Scalar &, Scalar &) const
General elliptic solver in a 2D case.
Coefficients storage for the multi-domain spectral method.
Definition mtbl_cf.h:196
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:393
virtual void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition scalar.C:359
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
virtual void set_etat_zero()
Sets the logical state to ETATZERO (zero).
Definition scalar.C:330
Valeur & set_spectral_va()
Returns va (read/write version).
Definition scalar.h:610
const Valeur & get_spectral_va() const
Returns va (read only version).
Definition scalar.h:607
int get_etat() const
Returns the logical state ETATNONDEF (undefined), ETATZERO (null) or ETATQCQ (ordinary).
Definition scalar.h:560
void set_dzpuis(int)
Modifies the dzpuis flag.
Definition scalar.C:814
Values and coefficients of a (real-value) function.
Definition valeur.h:297
int get_etat() const
Returns the logical state.
Definition valeur.h:760
Mtbl_cf * c_cf
Coefficients of the spectral expansion of the function.
Definition valeur.h:312
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:763
const Map & get_mp() const
Returns the mapping.
Definition tensor.h:874
Lorene prototypes.
Definition app_hor.h:67