File: param_elliptic_2d.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 (141 lines) | stat: -rw-r--r-- 4,154 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
/*
 *   Copyright (c) 2004 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 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 param_elliptic_2d_C[] = "$Header: /cvsroot/Lorene/C++/Source/Param_elliptic/param_elliptic_2d.C,v 1.3 2014/10/13 08:53:37 j_novak Exp $" ;

/*
 * $Id: param_elliptic_2d.C,v 1.3 2014/10/13 08:53:37 j_novak Exp $
 * $Log: param_elliptic_2d.C,v $
 * Revision 1.3  2014/10/13 08:53:37  j_novak
 * Lorene classes and functions now belong to the namespace Lorene.
 *
 * Revision 1.2  2014/10/06 15:13:15  j_novak
 * Modified #include directives to use c++ syntax.
 *
 * Revision 1.1  2004/08/24 09:14:49  p_grandclement
 * Addition of some new operators, like Poisson in 2d... It now requieres the
 * GSL library to work.
 *
 * Also, the way a variable change is stored by a Param_elliptic is changed and
 * no longer uses Change_var but rather 2 Scalars. The codes using that feature
 * will requiere some modification. (It should concern only the ones about monopoles)
 *
 * 
 * $Header: /cvsroot/Lorene/C++/Source/Param_elliptic/param_elliptic_2d.C,v 1.3 2014/10/13 08:53:37 j_novak Exp $
 *
 */

#include "headcpp.h"

#include <cmath>
#include <cstdlib>

#include "base_val.h" 
#include "map.h"
#include "ope_elementary.h"
#include "param_elliptic.h"
#include "change_var.h"
#include "scalar.h"


namespace Lorene {
void Param_elliptic::set_poisson_2d(const Scalar& source, bool indic) {

  int dzpuis = source.get_dzpuis() ;

  if (type_map != MAP_AFF) {
    cout << "set_poisson_2d only defined for an affine mapping..." << endl ;
    abort() ;
  }
  else {

    int nz = get_mp().get_mg()->get_nzone() ;
    
    int nr ;
    double alpha, beta ;
    int m_quant, l_quant, base_r_1d ;

    int conte = 0 ;
    for (int l=0 ; l<nz ; l++) {
      
      nr = get_mp().get_mg()->get_nr(l) ;
      alpha = get_alpha (l) ;
      beta = get_beta (l) ;

      for (int k=0 ; k<get_mp().get_mg()->get_np(l)+1 ; k++)
	for (int j=0 ; j<get_mp().get_mg()->get_nt(l) ; j++) {
	  if (operateurs[conte] != 0x0)	    
	    delete operateurs[conte] ;
	    source.get_spectral_va().base.give_quant_numbers(l, k, j, m_quant, l_quant, base_r_1d) ;
	    if (k!=1)
	      if ((indic) || ((!indic) && (l_quant !=0)))
		operateurs[conte] = new Ope_poisson_2d (nr, base_r_1d, alpha, beta, l_quant, dzpuis) ;
	      else
		operateurs[conte] = 0x0 ;
	    else
	      operateurs[conte] = 0x0 ;
	    conte ++ ;
	  }
    }
  }
}

void Param_elliptic::set_helmholtz_minus_2d(int zone, double masse, const Scalar& source) {

  int dzpuis = source.get_dzpuis() ;
  assert (masse > 0) ;

  if (type_map != MAP_AFF) {
    cout << "set_helmholtz_minus_2d only defined for an affine mapping..." << endl ;
    abort() ;
  }
  else {

    int nz = get_mp().get_mg()->get_nzone() ;
    if (zone == nz-1)  
      source.check_dzpuis(2) ; 
    int nr ;
    double alpha, beta ;
    int m_quant, l_quant, base_r_1d ;

    int conte = 0 ;
    for (int l=0 ; l<nz ; l++) {
      
      nr = get_mp().get_mg()->get_nr(l) ;
      alpha = get_alpha (l) ;
      beta = get_beta (l) ;

      for (int k=0 ; k<get_mp().get_mg()->get_np(l)+1 ; k++)
	for (int j=0 ; j<get_mp().get_mg()->get_nt(l) ; j++) {
	  if (l==zone) {
	    if (operateurs[conte] != 0x0) {	    
	      delete operateurs[conte] ;
	      source.get_spectral_va().base.give_quant_numbers 
		(l, k, j, m_quant, l_quant, base_r_1d) ;
	      operateurs[conte] = new Ope_helmholtz_minus_2d (nr, base_r_1d, alpha, beta, l_quant, masse, dzpuis) ;
	    }
	  } 
	  conte ++ ;
	}
    }
  }
}

}