File: tenseur_pde_regu.C

package info (click to toggle)
lorene 0.0.0~cvs20161116%2Bdfsg-1
  • links: PTS, VCS
  • area: main
  • in suites: bookworm, bullseye, buster, stretch
  • size: 26,444 kB
  • ctags: 13,953
  • sloc: cpp: 212,946; fortran: 21,645; makefile: 1,750; sh: 4
file content (192 lines) | stat: -rw-r--r-- 5,776 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
/*
 *  Method of class Tenseur to solve a vector Poisson equation
 *  by regularizing its source.
 *
 *  (see file tenseur.h for documentation).
 *
 */

/*
 *   Copyright (c) 2000-2001 Keisuke Taniguchi
 *   Copyright (c) 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 tenseur_pde_regu_C[] = "$Header: /cvsroot/Lorene/C++/Source/Tenseur/tenseur_pde_regu.C,v 1.4 2014/10/13 08:53:42 j_novak Exp $" ;

/*
 * $Id: tenseur_pde_regu.C,v 1.4 2014/10/13 08:53:42 j_novak Exp $
 * $Log: tenseur_pde_regu.C,v $
 * Revision 1.4  2014/10/13 08:53:42  j_novak
 * Lorene classes and functions now belong to the namespace Lorene.
 *
 * Revision 1.3  2003/10/03 15:58:51  j_novak
 * Cleaning of some headers
 *
 * Revision 1.2  2002/08/07 16:14:11  j_novak
 * class Tenseur can now also handle tensor densities, this should be transparent to older codes
 *
 * Revision 1.1.1.1  2001/11/20 15:19:30  e_gourgoulhon
 * LORENE
 *
 * Revision 2.1  2001/01/15  11:01:34  phil
 * vire version sans parametres
 *
 * Revision 2.0  2000/10/06  15:34:03  keisuke
 * Initial revision.
 *
 *
 * $Header: /cvsroot/Lorene/C++/Source/Tenseur/tenseur_pde_regu.C,v 1.4 2014/10/13 08:53:42 j_novak Exp $
 *
 */

// Header Lorene:
#include "param.h"
#include "tenseur.h"

		    //-----------------------------------//
		    //      Vectorial Poisson equation	 //
		    //-----------------------------------//

// Version avec parametres
// -----------------------
namespace Lorene {
void Tenseur::poisson_vect_regu(int k_div, int nzet, double unsgam1,
				double lambda, Param& para, Tenseur& shift,
				Tenseur& vecteur, Tenseur& scalaire) const {
    assert (lambda != -1) ;
    
    // Verifications d'usage ...
    assert (valence == 1) ;
    assert (shift.get_valence() == 1) ;
    assert (shift.get_type_indice(0) == type_indice(0)) ;
    assert (vecteur.get_valence() == 1) ;
    assert (vecteur.get_type_indice(0) == type_indice(0)) ;
    assert (scalaire.get_valence() == 0) ;
    assert (etat != ETATNONDEF) ;

    // Nothing to do if the source is zero
    if (etat == ETATZERO) {

	shift.set_etat_zero() ; 

	vecteur.set_etat_qcq() ;
	for (int i=0; i<3; i++) {
	    vecteur.set(i) = 0 ; 
	}

	scalaire.set_etat_qcq() ;
	scalaire.set() = 0 ;  

	return ; 
    }

    for (int i=0 ; i<3 ; i++)
	assert ((*this)(i).check_dzpuis(4)) ;

    Tenseur vecteur_regu(*mp, 1, CON, mp->get_bvect_cart(), metric, poids) ;
    Tenseur vecteur_div(*mp, 1, CON, mp->get_bvect_cart(), metric, poids) ;
    Tenseur dvect_div(*mp, 1, CON, mp->get_bvect_cart(), metric, poids) ;
    Tenseur souvect_regu(*mp, 1, CON, mp->get_bvect_cart(), metric, poids) ;
    Tenseur souvect_div(*mp, 1, CON, mp->get_bvect_cart(), metric, poids) ;

    vecteur_regu.set_etat_qcq() ;
    vecteur_div.set_etat_qcq() ;
    dvect_div.set_etat_qcq() ;
    souvect_regu.set_etat_qcq() ;
    souvect_div.set_etat_qcq() ;

    // On construit le tableau contenant le terme P_i ...

    // Apply only to x and y components because poisson_regular does not
    // work for z component due to the symmetry.
    for (int i=0 ; i<2 ; i++) {
	Param* par = mp->donne_para_poisson_vect(para, i) ; 

	(*this)(i).poisson_regular(k_div, nzet, unsgam1, *par,
				   vecteur.set(i),
				   vecteur_regu.set(i), vecteur_div.set(i),
				   dvect_div,
				   souvect_regu.set(i), souvect_div.set(i)) ;

	delete par ; 
    }

    Param* par = mp->donne_para_poisson_vect(para, 2) ;

    (*this)(2).poisson(*par, vecteur.set(2)) ;

    delete par ;

    vecteur.set_triad( *triad ) ; 
    
    // Equation de Poisson scalaire :
    Tenseur source_scal (-skxk(*this)) ;
      
    assert (source_scal().check_dzpuis(3)) ; 

    par = mp->donne_para_poisson_vect(para, 3) ; 

    Tenseur scalaire_regu(*mp, metric, poids) ;
    Tenseur scalaire_div(*mp, metric, poids) ;
    Tenseur dscal_div(*mp, 1, CON, mp->get_bvect_cart(), metric, poids) ;
    Cmp souscal_regu(mp) ;
    Cmp souscal_div(mp) ;

    scalaire_regu.set_etat_qcq() ;
    scalaire_div.set_etat_qcq() ;
    dscal_div.set_etat_qcq() ;

    souscal_regu.std_base_scal() ;
    souscal_div.std_base_scal() ;

    source_scal().poisson_regular(k_div, nzet, unsgam1, *par,
				  scalaire.set(),
				  scalaire_regu.set(), scalaire_div.set(),
				  dscal_div, souscal_regu, souscal_div) ;
    
    delete par ; 

    // On construit le tableau contenant le terme d xsi / d x_i ...
    Tenseur auxiliaire(scalaire) ;
    Tenseur dxsi (auxiliaire.gradient()) ;
    dxsi.dec2_dzpuis() ;
 
    // On construit le tableau contenant le terme x_k d P_k / d x_i
    Tenseur dp (skxk(vecteur.gradient())) ;
    dp.dec_dzpuis() ;
    
    // Il ne reste plus qu'a tout ranger dans P :
    // The final computation is done component by component because
    // d_khi and x_d_w are covariant comp. whereas w_shift is
    // contravariant

    shift.set_etat_qcq() ; 

    for (int i=0 ; i<3 ; i++)
	shift.set(i) = (lambda+2)/2/(lambda+1) * vecteur(i) 
			    - (lambda/2/(lambda+1)) * (dxsi(i) + dp(i)) ;   
			    
    shift.set_triad( *(vecteur.triad) ) ; 

}


}