KADATH
domain_shell_outer_adapted_affecte_tau_one_coef.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_outer_adapted::affecte_tau_one_coef_val_domain (Val_domain& so, int mlim, int cc, int& conte) const {
30 
31  int kmin = 2*mlim+2 ;
32 
33  so.is_zero = false ;
34  so.allocate_coef() ;
35  *so.cf=0. ;
36  Index pos_cf(nbr_coefs) ;
37 
38  bool found = false ;
39 
40  // True values
41  // Loop on phi :
42  for (int k=0 ; k<nbr_coefs(2)-1 ; k++)
43  if (k!=1) {
44  pos_cf.set(2) = k ;
45  // Loop on theta
46  int baset = (*so.get_base().bases_1d[1]) (k) ;
47  for (int j=0 ; j<nbr_coefs(1) ; j++) {
48  pos_cf.set(1) = j ;
49  bool true_tet = true ;
50  switch (baset) {
51  case COS_EVEN:
52  if ((j==0) && (k>=kmin))
53  true_tet = false ;
54  break ;
55  case COS_ODD:
56  if ((j==nbr_coefs(1)-1) || ((j==0) && (k>=kmin)))
57  true_tet = false ;
58  break ;
59  case SIN_EVEN:
60  if (((j==1)&&(k>=kmin+2))||(j==0) || (j==nbr_coefs(1)-1))
61  true_tet = false ;
62  break ;
63  case SIN_ODD:
64  if (((j==0)&&(k>=kmin+2))||(j==nbr_coefs(1)-1))
65  true_tet = false ;
66  break ;
67  default:
68  cerr << "Unknow theta basis in Domain_shell_outer_adapted::affecte_one_coef_val_domain" << endl ;
69  abort() ;
70  }
71 
72  if (true_tet)
73  for (int i=0 ; i<nbr_coefs(0) ; i++) {
74  pos_cf.set(0) = i ;
75  if (conte==cc) {
76  so.cf->set(pos_cf) = 1;
77  found = true ;
78  // regularity ??
79  if ((baset==COS_EVEN) || (baset==COS_ODD))
80  if (k>=kmin) {
81  pos_cf.set(1) = 0 ;
82  so.cf->set(pos_cf) = -1 ;
83  }
84 
85  if (baset==SIN_EVEN)
86  if (k>=kmin+2) {
87  pos_cf.set(1) = 1 ;
88  so.cf->set(pos_cf) = -j ;
89  }
90  if (baset==SIN_ODD)
91  if (k>=kmin+2) {
92  pos_cf.set(1) = 0 ;
93  so.cf->set(pos_cf) = -(2*j+1) ;
94  }
95 
96  }
97  else {
98  so.cf->set(pos_cf) = 0. ;
99  }
100  conte ++ ;
101  }
102  }
103  }
104  // If not found put to zero :
105  if (!found)
106  so.set_zero() ;
107 }
108 
109 void Domain_shell_outer_adapted::affecte_tau_one_coef (Tensor& tt, int dom, int cc, int& pos_cf) const {
110 
111  // Check right domain
112  assert (tt.get_space().get_domain(dom)==this) ;
113 
114  int val = tt.get_valence() ;
115  switch (val) {
116  case 0 :
117  affecte_tau_one_coef_val_domain (tt.set().set_domain(dom), 0, cc, pos_cf) ;
118  break ;
119  case 1 : {
120  bool found = false ;
121  // Cartesian basis
122  if (tt.get_basis().get_basis(dom)==CARTESIAN_BASIS) {
123  affecte_tau_one_coef_val_domain (tt.set(1).set_domain(dom), 0, cc, pos_cf) ;
124  affecte_tau_one_coef_val_domain (tt.set(2).set_domain(dom), 0, cc, pos_cf) ;
125  affecte_tau_one_coef_val_domain (tt.set(3).set_domain(dom), 0, cc, pos_cf) ;
126  found = true ;
127  }
128  // Spherical coordinates
129  if (tt.get_basis().get_basis(dom)==SPHERICAL_BASIS) {
130  affecte_tau_one_coef_val_domain (tt.set(1).set_domain(dom), 0, cc, pos_cf) ;
131  affecte_tau_one_coef_val_domain (tt.set(2).set_domain(dom), 1, cc, pos_cf) ;
132  affecte_tau_one_coef_val_domain (tt.set(3).set_domain(dom), 1, cc, pos_cf) ;
133  found = true ;
134  }
135 #ifndef REMOVE_ALL_CHECKS
136  if (!found) {
137  cerr << "Unknown type of vector Domain_shell_outer_adapted::affecte_tau_one_coef" << endl ;
138  abort() ;
139  }
140 #endif
141  }
142  break ;
143  case 2 : {
144  bool found = false ;
145  // Cartesian basis and symetric
146  if ((tt.get_basis().get_basis(dom)==CARTESIAN_BASIS) && (tt.get_n_comp()==6)) {
147  affecte_tau_one_coef_val_domain (tt.set(1,1).set_domain(dom), 0, cc, pos_cf) ;
148  affecte_tau_one_coef_val_domain (tt.set(1,2).set_domain(dom), 0, cc, pos_cf) ;
149  affecte_tau_one_coef_val_domain (tt.set(1,3).set_domain(dom), 0, cc, pos_cf) ;
150  affecte_tau_one_coef_val_domain (tt.set(2,2).set_domain(dom), 0, cc, pos_cf) ;
151  affecte_tau_one_coef_val_domain (tt.set(2,3).set_domain(dom), 0, cc, pos_cf) ;
152  affecte_tau_one_coef_val_domain (tt.set(3,3).set_domain(dom), 0, cc, pos_cf) ;
153  found = true ;
154  }
155  // Cartesian basis and not symetric
156  if ((tt.get_basis().get_basis(dom)==CARTESIAN_BASIS) && (tt.get_n_comp()==9)) {
157  affecte_tau_one_coef_val_domain (tt.set(1,1).set_domain(dom), 0, cc, pos_cf) ;
158  affecte_tau_one_coef_val_domain (tt.set(1,2).set_domain(dom), 0, cc, pos_cf) ;
159  affecte_tau_one_coef_val_domain (tt.set(1,3).set_domain(dom), 0, cc, pos_cf) ;
160  affecte_tau_one_coef_val_domain (tt.set(2,1).set_domain(dom), 0, cc, pos_cf) ;
161  affecte_tau_one_coef_val_domain (tt.set(2,2).set_domain(dom), 0, cc, pos_cf) ;
162  affecte_tau_one_coef_val_domain (tt.set(2,3).set_domain(dom), 0, cc, pos_cf) ;
163  affecte_tau_one_coef_val_domain (tt.set(3,1).set_domain(dom), 0, cc, pos_cf) ;
164  affecte_tau_one_coef_val_domain (tt.set(3,2).set_domain(dom), 0, cc, pos_cf) ;
165  affecte_tau_one_coef_val_domain (tt.set(3,3).set_domain(dom), 0, cc, pos_cf) ;
166  found = true ;
167  }
168  // Spherical coordinates and symetric
169  if ((tt.get_basis().get_basis(dom)==SPHERICAL_BASIS) && (tt.get_n_comp()==6)) {
170  affecte_tau_one_coef_val_domain (tt.set(1,1).set_domain(dom), 0,cc ,pos_cf) ;
171  affecte_tau_one_coef_val_domain (tt.set(1,2).set_domain(dom), 1, cc, pos_cf) ;
172  affecte_tau_one_coef_val_domain (tt.set(1,3).set_domain(dom), 1, cc, pos_cf) ;
173  affecte_tau_one_coef_val_domain (tt.set(2,2).set_domain(dom), 2, cc, pos_cf) ;
174  affecte_tau_one_coef_val_domain (tt.set(2,3).set_domain(dom), 2, cc, pos_cf) ;
175  affecte_tau_one_coef_val_domain (tt.set(3,3).set_domain(dom), 2, cc, pos_cf) ;
176  found = true ;
177  }
178  // Spherical coordinates and not symetric
179  if ((tt.get_basis().get_basis(dom)==SPHERICAL_BASIS) && (tt.get_n_comp()==9)) {
180  affecte_tau_one_coef_val_domain (tt.set(1,1).set_domain(dom), 0,cc ,pos_cf) ;
181  affecte_tau_one_coef_val_domain (tt.set(1,2).set_domain(dom), 1, cc, pos_cf) ;
182  affecte_tau_one_coef_val_domain (tt.set(1,3).set_domain(dom), 1, cc, pos_cf) ;
183  affecte_tau_one_coef_val_domain (tt.set(2,1).set_domain(dom), 1, cc, pos_cf) ;
184  affecte_tau_one_coef_val_domain (tt.set(2,2).set_domain(dom), 2, cc, pos_cf) ;
185  affecte_tau_one_coef_val_domain (tt.set(2,3).set_domain(dom), 2, cc, pos_cf) ;
186  affecte_tau_one_coef_val_domain (tt.set(3,1).set_domain(dom), 1, cc, pos_cf) ;
187  affecte_tau_one_coef_val_domain (tt.set(3,2).set_domain(dom), 2, cc, pos_cf) ;
188  affecte_tau_one_coef_val_domain (tt.set(3,3).set_domain(dom), 2, cc, pos_cf) ;
189  found = true ;
190  }
191 #ifndef REMOVE_ALL_CHECKS
192  if (!found) {
193  cerr << "Unknown type of 2-tensor Domain_shell_outer_adapted::affecte_tau_one_coef" << endl ;
194  abort() ;
195  }
196 #endif
197  }
198  break ;
199  default :
200  cerr << "Valence " << val << " not implemented in Domain_shell_outer_adapted::affecte_tau" << endl ;
201  break ;
202  }
203 }
204 }
205 
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
void affecte_tau_one_coef_val_domain(Val_domain &so, int mlim, int cf, int &pos_cf) const
Affects some coefficients to a Val_domain.
virtual void affecte_tau_one_coef(Tensor &, int, int, int &) const
Sets at most one coefficient of a Tensor to 1.
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 set_zero()
Sets the Val_domain to zero (logical state to zero and arrays destroyed).
Definition: val_domain.cpp:223
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
bool is_zero
Indicator used for null fields (for speed issues).
Definition: val_domain.hpp:74
const Base_spectral & get_base() const
Returns the basis of decomposition.
Definition: val_domain.hpp:122