KADATH
domain_shell_inner_adapted_affecte_tau.cpp
1 /*
2  Copyright 2017 Philippe Grandclement
3 
4  This file is part of Kadath.
5 
6  Kadath is free software: you can redistribute it and/or modify
7  it under the terms of the GNU General Public License as published by
8  the Free Software Foundation, either version 3 of the License, or
9  (at your option) any later version.
10 
11  Kadath is distributed in the hope that it will be useful,
12  but WITHOUT ANY WARRANTY; without even the implied warranty of
13  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14  GNU General Public License for more details.
15 
16  You should have received a copy of the GNU General Public License
17  along with Kadath. If not, see <http://www.gnu.org/licenses/>.
18 */
19 
20 #include "headcpp.hpp"
21 #include "adapted.hpp"
22 #include "point.hpp"
23 #include "array_math.hpp"
24 #include "scalar.hpp"
25 #include "tensor_impl.hpp"
26 #include "tensor.hpp"
27 
28 namespace Kadath {
29 void Domain_shell_inner_adapted::affecte_tau_val_domain (Val_domain& so, int mlim, const Array<double>& values, int& conte) const {
30 
31  int kmin = 2*mlim+2 ;
32 
33  so.allocate_coef() ;
34  *so.cf = 0. ;
35  Index pos_cf (nbr_coefs) ;
36 
37  // True values
38  // Loop on phi :
39  for (int k=0 ; k<nbr_coefs(2)-1 ; k++)
40  if (k!=1) {
41  pos_cf.set(2) = k ;
42  // Loop on theta
43  int baset = (*so.get_base().bases_1d[1]) (k) ;
44  for (int j=0 ; j<nbr_coefs(1) ; j++) {
45  pos_cf.set(1) = j ;
46  bool true_tet = true ;
47  switch (baset) {
48  case COS_EVEN:
49  if ((j==0) && (k>=kmin))
50  true_tet = false ;
51  break ;
52  case COS_ODD:
53  if ((j==nbr_coefs(1)-1) || ((j==0) && (k>=kmin)))
54  true_tet = false ;
55  break ;
56  case SIN_EVEN:
57  if (((j==1)&&(k>=kmin+2))||(j==0) || (j==nbr_coefs(1)-1))
58  true_tet = false ;
59  break ;
60  case SIN_ODD:
61  if (((j==0)&&(k>kmin+2)) || (j==nbr_coefs(1)-1))
62  true_tet = false ;
63  break ;
64  default:
65  cerr << "Unknow theta basis in Domain_shell_inner_adapted::affecte_tau_val_domain" << endl ;
66  abort() ;
67  }
68 
69  if (true_tet)
70  for (int i=0 ; i<nbr_coefs(0) ; i++) {
71  pos_cf.set(0) = i ;
72  so.cf->set(pos_cf) += values(conte);
73  conte ++ ;
74  }
75  }
76  }
77 
78  // Appropriate regularisation
79  // Loop on phi :
80  for (int k=0 ; k<nbr_coefs(2)-1 ; k++) {
81  pos_cf.set(2) = k ;
82  int baset = (*so.get_base().bases_1d[1]) (k) ;
83  // Loop on r :
84  for (int i=0 ; i<nbr_coefs(0) ; i++) {
85  pos_cf.set(0) = i ;
86  double sum = 0. ;
87  switch (baset) {
88  case COS_EVEN:
89  if (k>=kmin) {
90  for (int j=1 ; j<nbr_coefs(1) ; j++) {
91  pos_cf.set(1) = j ;
92  sum += (*so.cf)(pos_cf) ;
93  }
94  pos_cf.set(1) = 0 ;
95  so.cf->set(pos_cf) = -sum ;
96  }
97  break ;
98  case COS_ODD:
99  if (k>=kmin) {
100  for (int j=1 ; j<nbr_coefs(1) ; j++) {
101  pos_cf.set(1) = j ;
102  sum += (*so.cf)(pos_cf) ;
103  }
104  pos_cf.set(1) = 0 ;
105  so.cf->set(pos_cf) = -sum ;
106  }
107  break ;
108  case SIN_EVEN:
109  if (k>=kmin+2) {
110  for (int j=2 ; j<nbr_coefs(1) ; j++) {
111  pos_cf.set(1) = j ;
112  sum += j*(*so.cf)(pos_cf) ;
113  }
114  pos_cf.set(1) = 1 ;
115  so.cf->set(pos_cf) = -sum ;
116  }
117  break ;
118  case SIN_ODD:
119  if (k>=kmin+2) {
120  for (int j=1 ; j<nbr_coefs(1) ; j++) {
121  pos_cf.set(1) = j ;
122  sum += (2*j+1)*(*so.cf)(pos_cf) ;
123  }
124  pos_cf.set(1) = 0 ;
125  so.cf->set(pos_cf) = -sum ;
126  }
127  break ;
128  default:
129  cerr << "Unknow theta basis in Domain_shell_inner_adapted::affecte_tau_val_domain" << endl ;
130  abort() ;
131  }
132  }
133  }
134 }
135 
136 void Domain_shell_inner_adapted::affecte_tau (Tensor& tt, int dom, const Array<double>& cf, int& pos_cf) const {
137 
138  // Check right domain
139  assert (tt.get_space().get_domain(dom)==this) ;
140 
141  int val = tt.get_valence() ;
142  switch (val) {
143  case 0 :
144  affecte_tau_val_domain (tt.set().set_domain(dom), 0, cf, pos_cf) ;
145  break ;
146  case 1 : {
147  bool found = false ;
148  // Cartesian basis
149  if (tt.get_basis().get_basis(dom)==CARTESIAN_BASIS) {
150  affecte_tau_val_domain (tt.set(1).set_domain(dom), 0, cf, pos_cf) ;
151  affecte_tau_val_domain (tt.set(2).set_domain(dom), 0, cf, pos_cf) ;
152  affecte_tau_val_domain (tt.set(3).set_domain(dom), 0, cf, pos_cf) ;
153  found = true ;
154  }
155  // Spherical coordinates
156  if (tt.get_basis().get_basis(dom)==SPHERICAL_BASIS) {
157  affecte_tau_val_domain (tt.set(1).set_domain(dom), 0, cf, pos_cf) ;
158  affecte_tau_val_domain (tt.set(2).set_domain(dom), 1, cf, pos_cf) ;
159  affecte_tau_val_domain (tt.set(3).set_domain(dom), 1, cf, pos_cf) ;
160  found = true ;
161  }
162  if (!found) {
163  cerr << "Unknown type of vector Domain_shell_inner_adapted::affecte_tau" << endl ;
164  abort() ;
165  }
166  }
167  break ;
168  case 2 : {
169  bool found = false ;
170  // Cartesian basis and symetric
171  if ((tt.get_basis().get_basis(dom)==CARTESIAN_BASIS) && (tt.get_n_comp()==6)) {
172  affecte_tau_val_domain (tt.set(1,1).set_domain(dom), 0, cf, pos_cf) ;
173  affecte_tau_val_domain (tt.set(1,2).set_domain(dom), 0, cf, pos_cf) ;
174  affecte_tau_val_domain (tt.set(1,3).set_domain(dom), 0, cf, pos_cf) ;
175  affecte_tau_val_domain (tt.set(2,2).set_domain(dom), 0, cf, pos_cf) ;
176  affecte_tau_val_domain (tt.set(2,3).set_domain(dom), 0, cf, pos_cf) ;
177  affecte_tau_val_domain (tt.set(3,3).set_domain(dom), 0, cf, pos_cf) ;
178  found = true ;
179  }
180  // Cartesian basis and not symetric
181  if ((tt.get_basis().get_basis(dom)==CARTESIAN_BASIS) && (tt.get_n_comp()==9)) {
182  affecte_tau_val_domain (tt.set(1,1).set_domain(dom), 0, cf, pos_cf) ;
183  affecte_tau_val_domain (tt.set(1,2).set_domain(dom), 0, cf, pos_cf) ;
184  affecte_tau_val_domain (tt.set(1,3).set_domain(dom), 0, cf, pos_cf) ;
185  affecte_tau_val_domain (tt.set(2,1).set_domain(dom), 0, cf, pos_cf) ;
186  affecte_tau_val_domain (tt.set(2,2).set_domain(dom), 0, cf, pos_cf) ;
187  affecte_tau_val_domain (tt.set(2,3).set_domain(dom), 0, cf, pos_cf) ;
188  affecte_tau_val_domain (tt.set(3,1).set_domain(dom), 0, cf, pos_cf) ;
189  affecte_tau_val_domain (tt.set(3,2).set_domain(dom), 0, cf, pos_cf) ;
190  affecte_tau_val_domain (tt.set(3,3).set_domain(dom), 0, cf, pos_cf) ;
191  found = true ;
192  }
193  // Spherical coordinates and symetric
194  if ((tt.get_basis().get_basis(dom)==SPHERICAL_BASIS) && (tt.get_n_comp()==6)) {
195  affecte_tau_val_domain (tt.set(1,1).set_domain(dom), 0,cf ,pos_cf) ;
196  affecte_tau_val_domain (tt.set(1,2).set_domain(dom), 1, cf, pos_cf) ;
197  affecte_tau_val_domain (tt.set(1,3).set_domain(dom), 1, cf, pos_cf) ;
198  affecte_tau_val_domain (tt.set(2,2).set_domain(dom), 2, cf, pos_cf) ;
199  affecte_tau_val_domain (tt.set(2,3).set_domain(dom), 2, cf, pos_cf) ;
200  affecte_tau_val_domain (tt.set(3,3).set_domain(dom), 2, cf, pos_cf) ;
201  found = true ;
202  }
203  // Spherical coordinates and not symetric
204  if ((tt.get_basis().get_basis(dom)==SPHERICAL_BASIS) && (tt.get_n_comp()==9)) {
205  affecte_tau_val_domain (tt.set(1,1).set_domain(dom), 0,cf ,pos_cf) ;
206  affecte_tau_val_domain (tt.set(1,2).set_domain(dom), 1, cf, pos_cf) ;
207  affecte_tau_val_domain (tt.set(1,3).set_domain(dom), 1, cf, pos_cf) ;
208  affecte_tau_val_domain (tt.set(2,1).set_domain(dom), 1, cf, pos_cf) ;
209  affecte_tau_val_domain (tt.set(2,2).set_domain(dom), 2, cf, pos_cf) ;
210  affecte_tau_val_domain (tt.set(2,3).set_domain(dom), 2, cf, pos_cf) ;
211  affecte_tau_val_domain (tt.set(3,1).set_domain(dom), 1, cf, pos_cf) ;
212  affecte_tau_val_domain (tt.set(3,2).set_domain(dom), 2, cf, pos_cf) ;
213  affecte_tau_val_domain (tt.set(3,3).set_domain(dom), 2, cf, pos_cf) ;
214  found = true ;
215  }
216  if (!found) {
217  cerr << "Unknown type of 2-tensor Domain_shell_inner_adapted::affecte_tau" << endl ;
218  abort() ;
219  }
220  }
221  break ;
222  default :
223  cerr << "Valence " << val << " not implemented in Domain_shell_inner_adapted::affecte_tau" << endl ;
224  break ;
225  }
226 }
227 }
228 
reference set(const Index &pos)
Read/write of an element.
Definition: array.hpp:186
Bases_container bases_1d
Arrays containing the various basis of decomposition.
int get_basis(int nd) const
Read only the basis in a given domain.
Definition: base_tensor.hpp:93
virtual void affecte_tau(Tensor &, int, const Array< double > &, int &) const
Affects some coefficients to a Tensor.
void affecte_tau_val_domain(Val_domain &so, int mlim, const Array< double > &cf, int &pos_cf) const
Affects some coefficients to a Val_domain.
Dim_array nbr_coefs
Number of coefficients.
Definition: space.hpp:66
Class that gives the position inside a multi-dimensional Array.
Definition: index.hpp:38
int & set(int i)
Read/write of the position in a given dimension.
Definition: index.hpp:72
Val_domain & set_domain(int)
Read/write of a particular Val_domain.
Definition: scalar.hpp:555
const Domain * get_domain(int i) const
returns a pointer on the domain.
Definition: space.hpp:1385
Tensor handling.
Definition: tensor.hpp:149
Scalar & set(const Array< int > &ind)
Returns the value of a component (read/write version).
Definition: tensor_impl.hpp:91
const Base_tensor & get_basis() const
Returns the vectorial basis (triad) on which the components are defined.
Definition: tensor.hpp:504
int get_n_comp() const
Returns the number of stored components.
Definition: tensor.hpp:514
int get_valence() const
Returns the valence.
Definition: tensor.hpp:509
const Space & get_space() const
Returns the Space.
Definition: tensor.hpp:499
Class for storing the basis of decompositions of a field and its values on both the configuration and...
Definition: val_domain.hpp:69
void allocate_coef()
Allocates the values in the coefficient space and destroys the values in the configuration space.
Definition: val_domain.cpp:216
Array< double > * cf
Pointer on the Array of the values in the coefficients space.
Definition: val_domain.hpp:77
const Base_spectral & get_base() const
Returns the basis of decomposition.
Definition: val_domain.hpp:122