File: map_af_poisson.C

package info (click to toggle)
lorene 0.0.0~cvs20161116%2Bdfsg-2
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid, trixie
  • size: 26,472 kB
  • sloc: cpp: 212,946; fortran: 21,645; makefile: 1,750; sh: 4
file content (219 lines) | stat: -rw-r--r-- 6,349 bytes parent folder | download | duplicates (2)
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
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
/*
 * Method of the class Map_af for the resolution of the scalar Poisson
 *  equation
 */

/*
 *   Copyright (c) 1999-2001 Eric Gourgoulhon
 *   Copyright (c) 1999-2001 Philippe Grandclement
 *
 *   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 as published by
 *   the Free Software Foundation; either version 2 of the License, or
 *   (at your option) any later version.
 *
 *   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_C[] = "$Header: /cvsroot/Lorene/C++/Source/Map/map_af_poisson.C,v 1.6 2014/10/13 08:53:02 j_novak Exp $" ;

/*
 * $Id: map_af_poisson.C,v 1.6 2014/10/13 08:53:02 j_novak Exp $
 * $Log: map_af_poisson.C,v $
 * Revision 1.6  2014/10/13 08:53:02  j_novak
 * Lorene classes and functions now belong to the namespace Lorene.
 *
 * Revision 1.5  2005/08/25 12:14:09  p_grandclement
 * Addition of a new method to solve the scalar Poisson equation, based on a multi-domain Tau-method
 *
 * Revision 1.4  2004/05/06 15:25:39  e_gourgoulhon
 * The case dzpuis=5 with null value in CED is well treated now.
 *
 * Revision 1.3  2004/02/20 10:55:23  j_novak
 * The versions dzpuis 5 -> 3 has been improved and polished. Should be
 * operational now...
 *
 * Revision 1.2  2004/02/06 10:53:52  j_novak
 * New dzpuis = 5 -> dzpuis = 3 case (not ready yet).
 *
 * Revision 1.1.1.1  2001/11/20 15:19:27  e_gourgoulhon
 * LORENE
 *
 * Revision 1.9  2000/05/22  13:46:48  phil
 * ajout du cas dzpuis = 3
 *
 * Revision 1.8  2000/02/09  14:44:24  eric
 * Traitement de dzpuis ameliore.
 *
 * Revision 1.7  1999/12/22  16:37:10  eric
 * Ajout de pot.set_dzpuis(0) a la fin.
 *
 * Revision 1.6  1999/12/22  15:11:03  eric
 * Remplacement du test source.get_mp() == this  par
 *  source.get_mp()->get_mg() == mg
 * (idem pour pot),
 * afin de permettre l'appel par Map_et::poisson.
 *
 * Revision 1.5  1999/12/21  13:02:37  eric
 * Changement de prototype de la routine poisson : la solution est
 *  desormais passee en argument (et non plus en valeur de retour)
 *  pour s'adapter au prototype general de la fonction virtuelle
 *   Map::poisson.
 *
 * Revision 1.4  1999/12/21  10:06:29  eric
 * Ajout de l'argument (muet) Param&.
 *
 * Revision 1.3  1999/12/07  16:48:50  phil
 * On fait ylm_i avant de quitter
 *
 * Revision 1.2  1999/12/02  16:12:22  eric
 * *** empty log message ***
 *
 * Revision 1.1  1999/12/02  14:30:07  eric
 * Initial revision
 *
 *
 * $Header: /cvsroot/Lorene/C++/Source/Map/map_af_poisson.C,v 1.6 2014/10/13 08:53:02 j_novak Exp $
 *
 */

// Header Lorene:
#include "map.h"
#include "cmp.h"

namespace Lorene {
Mtbl_cf sol_poisson(const Map_af&, const Mtbl_cf&, int, bool match = true) ;
Mtbl_cf sol_poisson_tau(const Map_af&, const Mtbl_cf&, int) ;
//*****************************************************************************

void Map_af::poisson(const Cmp& source, Param& , Cmp& pot) const {
    
    assert(source.get_etat() != ETATNONDEF) ; 
    assert(source.get_mp()->get_mg() == mg) ; 
    assert(pot.get_mp()->get_mg() == mg) ; 

    assert( source.check_dzpuis(2) || source.check_dzpuis(4) 
	    || source.check_dzpuis(3) || source.check_dzpuis(5) ) ; 
    
    bool match = true ;

    int dzpuis ; 
    
    if ( (source.dz_nonzero()) || (source.get_dzpuis() > 3)) { //##awkward??
	dzpuis = source.get_dzpuis() ; 
    }
    else{
	dzpuis = 4 ; 
    }

    match = !(dzpuis == 5) ;

    // Spherical harmonic expansion of the source
    // ------------------------------------------
    
    const Valeur& sourva = source.va ; 

    if (sourva.get_etat() == ETATZERO) {
	pot.set_etat_zero() ;
	return ;  
    }

    // Spectral coefficients of the source
    assert(sourva.get_etat() == ETATQCQ) ; 
    
    Valeur rho(sourva.get_mg()) ; 
    sourva.coef() ; 
    rho = *(sourva.c_cf) ;	// copy of the coefficients of the source
    
    rho.ylm() ;			// spherical harmonic transforms 
        
    // Call to the Mtbl_cf version
    // ---------------------------
    Mtbl_cf resu = sol_poisson(*this, *(rho.c_cf), dzpuis, match) ;
    
    // Final result returned as a Cmp
    // ------------------------------
    
    pot.set_etat_zero() ;  // to call Cmp::del_t().

    pot.set_etat_qcq() ; 
    
    pot.va = resu ;
    (pot.va).ylm_i() ; // On repasse en base standard.	    

    (dzpuis == 5) ? pot.set_dzpuis(3) : pot.set_dzpuis(0) ; 
    
}


		//----------------------
		// Tau version method
		//---------------------


void Map_af::poisson_tau(const Cmp& source, Param& , Cmp& pot) const {
    
    assert(source.get_etat() != ETATNONDEF) ; 
    assert(source.get_mp()->get_mg() == mg) ; 
    assert(pot.get_mp()->get_mg() == mg) ; 

    assert( source.check_dzpuis(2) || source.check_dzpuis(4) 
	    || source.check_dzpuis(3)) ; 

    int dzpuis ; 
    
    if ( (source.dz_nonzero()) || (source.get_dzpuis() > 3)) { //##awkward??
	dzpuis = source.get_dzpuis() ; 
    }
    else{
	dzpuis = 4 ; 
    }

    // Spherical harmonic expansion of the source
    // ------------------------------------------
    
    const Valeur& sourva = source.va ; 

    if (sourva.get_etat() == ETATZERO) {
	pot.set_etat_zero() ;
	return ;  
    }

    // Spectral coefficients of the source
    assert(sourva.get_etat() == ETATQCQ) ; 
    
    Valeur rho(sourva.get_mg()) ; 
    sourva.coef() ; 
    rho = *(sourva.c_cf) ;	// copy of the coefficients of the source
    
    rho.ylm() ;			// spherical harmonic transforms 
        
    // Call to the Mtbl_cf version
    // ---------------------------
    
    Mtbl_cf resu = sol_poisson_tau(*this, *(rho.c_cf), dzpuis) ;
    
    // Final result returned as a Cmp
    // ------------------------------
    
    pot.set_etat_zero() ;  // to call Cmp::del_t().

    pot.set_etat_qcq() ; 
    
    pot.va = resu ;
    (pot.va).ylm_i() ; // On repasse en base standard.	    
    pot.set_dzpuis(0) ;
}

}