1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101
|
/*
* Resolution of the angular Poisson equation.
*
* (see file map.h for documentation)
*
*/
/*
* Copyright (c) 2003-2005 Eric Gourgoulhon & Jerome Novak
*
* This file is part of LORENE.
*
* LORENE is free software; you can redistribute it and/or modify
* it under the terms of the GNU General Public License version 2
* as published by the Free Software Foundation.
*
* LORENE is distributed in the hope that it will be useful,
* but WITHOUT ANY WARRANTY; without even the implied warranty of
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
* GNU General Public License for more details.
*
* You should have received a copy of the GNU General Public License
* along with LORENE; if not, write to the Free Software
* Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
*
*/
char map_af_poisson_angu_C[] = "$Header: /cvsroot/Lorene/C++/Source/Map/map_af_poisson_angu.C,v 1.4 2014/10/13 08:53:03 j_novak Exp $" ;
/*
* $Id: map_af_poisson_angu.C,v 1.4 2014/10/13 08:53:03 j_novak Exp $
* $Log: map_af_poisson_angu.C,v $
* Revision 1.4 2014/10/13 08:53:03 j_novak
* Lorene classes and functions now belong to the namespace Lorene.
*
* Revision 1.3 2005/04/04 21:31:31 e_gourgoulhon
* Added argument lambda to method poisson_angu
* to deal with the generalized angular Poisson equation:
* Lap_ang u + lambda u = source.
*
* Revision 1.2 2003/10/16 08:49:23 j_novak
* Added a flag to decide wether the output is in the Ylm or in the standard base.
*
* Revision 1.1 2003/10/15 21:11:26 e_gourgoulhon
* Added method poisson_angu.
*
*
*
* $Header: /cvsroot/Lorene/C++/Source/Map/map_af_poisson_angu.C,v 1.4 2014/10/13 08:53:03 j_novak Exp $
*
*/
// Lorene headers
#include "tensor.h"
#include "param.h"
namespace Lorene {
void Map_af::poisson_angu(const Scalar& source, Param& par, Scalar& uu,
double lambda) const {
assert(source.get_etat() != ETATNONDEF) ;
assert(&(source.get_mp()) == this ) ;
assert(&(uu.get_mp()) == this ) ;
// Spherical harmonic expansion of the source
// ------------------------------------------
const Valeur& sourva = source.get_spectral_va() ;
if (sourva.get_etat() == ETATZERO) {
uu.set_etat_zero() ;
return ;
}
// Spectral coefficients of the source
assert(sourva.get_etat() == ETATQCQ) ;
sourva.coef() ;
Valeur resu(mg) ;
resu = *(sourva.c_cf) ; // copy of the coefficients of the source
resu.ylm() ; // spherical harmonic transform
// Call to the Mtbl_cf version
// ---------------------------
(resu.c_cf)->poisson_angu(lambda) ;
if (par.get_n_int() == 0) resu.ylm_i() ; // Back to standard bases
//in the case of no flag present
// in the Param
// Final result returned as a Scalar
// ---------------------------------
uu = resu ;
uu.set_dzpuis( source.get_dzpuis() ) ; // dzpuis unchanged
}
}
|