KADATH
solver.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 "system_of_eqs.hpp"
21 #include "scalar.hpp"
22 #include "tensor_impl.hpp"
23 #include "metric.hpp"
24 #include "exceptions.hpp"
25 namespace Kadath {
27 
28 
29  Array<double> sec (sec_member()) ;
30  Array<double> errors (neq_int + neq) ;
31  errors = 0 ;
32 
33  int pos = 0 ;
34  for (int i=0 ; i<neq_int ; i++) {
35  errors.set(i) = fabs(sec(pos)) ;
36  pos ++ ;
37  }
38 
39  for (int i=0 ; i<neq ; i++) {
40  double max = 0 ;
41  for (int j=0 ; j<eq[i]->get_n_cond_tot() ; j++) {
42  if (fabs(sec(pos)) > max)
43  max = fabs(sec(pos)) ;
44  pos ++ ;
45  }
46  errors.set(neq_int+i) = max ;
47  }
48  return errors ;
49 }
50 
52 {
53 
54  if (met!=0x0)
55  for (int d=dom_min ; d<=dom_max ; d++)
56  met->update(d) ;
57  for (int i=0 ; i<ndef ; i++)
58  def[i]->compute_res() ;
59 
60  int conte = 0 ;
61  for (int i=0 ; i<neq ; i++)
62  eq[i]->apply(conte, results.set_data()) ;
63 
64  // Need to assert the size :
65  if (nbr_conditions==-1) {
66  nbr_conditions = 0 ;
67  for (int i=0 ; i<neq_int ; i++)
68  nbr_conditions ++ ;
69  for (int i=0 ; i<neq ; i++)
70  nbr_conditions += eq[i]->get_n_cond_tot() ;
71  }
72 
73 }
74 
76 
77  this->vars_to_terms() ;
79  // Computation of the second member itself :
81  res = 0 ;
82  int conte = 0 ;
83  int pos_res = 0 ;
84  for (int i=0 ; i<neq_int ; i++) {
85  res.set(pos_res) = eq_int[i]->get_val() ;
86  pos_res ++ ;
87  }
88 
89  for (int i=0 ; i<neq ; i++)
90  eq[i]->export_val(conte, results.set_data(), res, pos_res) ;
91  return res ;
92 }
93 
95 
96  xx_to_ders(xx) ;
97  if (met!=0x0)
98  for (int d=dom_min ; d<=dom_max ; d++)
99  met->update(d) ;
100  // Delete the definitions
101  for (int i=0 ; i<ndef ; i++)
102  def[i]->compute_res() ;
103 
104  int conte = 0 ;
105  for (int i=0 ; i<neq ; i++)
106  eq[i]->apply(conte, results.set_data()) ;
107 
108  if (nbr_conditions==-1) {
109  cerr << "Number of conditions unknown ; call sec_member first" << endl ;
110  abort() ;
111  }
112 
114  conte = 0 ;
115  int pos_res = 0 ;
116  res = 0 ;
117  for (int i=0 ; i<neq_int ; i++) {
118  res.set(pos_res) = eq_int[i]->get_der() ;
119  pos_res ++ ;
120  }
121 
122  for (int i=0 ; i<neq ; i++)
123  eq[i]->export_der(conte, results.set_data(), res, pos_res) ;
124  return res ;
125 }
126 
128 
129  assert ((cc>=0) && (cc<nbr_unknowns)) ;
130 
131  // Affecte nterms derivatives :
132  int conte = 0 ;
133  int zedom = -1 ;
134  bool is_var_double = false ;
135  Array<int> zedoms (2) ;
136  zedoms = - 1 ;
137  // Variable Domains :
138  espace.affecte_coef_to_variable_domains(conte, cc, zedoms) ;
139  if (zedoms(0)!=-1)
141  else {
142  for (int i=0 ; i<nterm_cst ; i++)
143  cst[i]->set_der_zero() ;
144  for (int i=0 ; i<nterm ; i++)
145  term[i]->set_der_zero() ;
146  }
147 
148  // Double
149  for (int i=0 ; i<nvar_double ; i++) {
150  if (conte==cc) {
151  for (int dd=dom_min ; dd<=dom_max ; dd++)
152  term_double[i*ndom+(dd-dom_min)]->set_der_d(1.) ;
153  is_var_double = true ;
154  }
155  else
156  for (int dd=dom_min ; dd<=dom_max ; dd++)
157  term_double[i*ndom+(dd-dom_min)]->set_der_d(0.) ;
158  conte ++ ;
159  }
160 
161 
162  // Fields
163  for (int i=0 ; i<nterm ; i++) {
164  int dom = term[i]->get_dom() ;
165  Tensor auxi (term[i]->get_val_t(), false) ;
166  try {
167  espace.get_domain(dom)->affecte_tau_one_coef(auxi, dom, cc, conte);
168  }
169 
170  catch(Unknown_base_error & e) {
171  std::cerr << "Error in System_of_eqs[rank=" << mpi_proc_rank << "]::do_col_J(cc = " << cc
172  /*<< ")\n"*/;
173 // std::cerr << "affecte_tau_one_coef raised the following exception : " << e.what() << "\n";
174 // std::cerr << "while calling espace.get_domain(dom)->affecte_tau_one_coef(auxi, dom, cc, conte)\n";
175  std::cerr << " with auxi=term[i]->get_val_t(), where i="<< i ;
176  std::cerr << ", dom=" << dom << " and conte=" << conte << std::endl;
177  throw e;
178  }
179  for (int j=0 ; j<auxi.get_n_comp() ; j++) {
180  // Si la base n'est pas affectee on la met
181  if ((!auxi(auxi.indices(j))(dom).check_if_zero()) && (!auxi(auxi.indices(j))(dom).get_base().is_def()))
182  auxi.set(auxi.indices(j)).set_domain(dom).set_base() = term[i]->get_val_t()(auxi.indices(j))(dom).get_base() ;
183  if ((zedom==-1) && (!auxi(auxi.indices(j))(dom).check_if_zero()))
184  zedom = dom ;
185  }
186 
187  term[i]->set_der_t(auxi + term[i]->get_der_t()) ;
188  }
189 
190  // Delete the metric derivative terms :
191  if (met!=0x0)
192  for (int d=dom_min ; d<=dom_max ; d++)
193  met->update(d) ;
194 
195  // Delete the definitions
196  for (int i=0 ; i<ndef ; i++)
197  def[i]->compute_res() ;
198 
199  conte = 0 ;
200  for (int i=0 ; i<neq ; i++) {
201  if ((is_var_double) || (eq[i]->take_into_account(zedom)) || (eq[i]->take_into_account(zedoms(0))) || (eq[i]->take_into_account(zedoms(1)))) {
202  eq[i]->apply(conte, results.set_data()) ;//NOT THREAD SAFE (but both eq and results are own by System_of_eqs)
203  }
204  else
205  conte += eq[i]->n_ope ;
206  }
207 
208  if (nbr_conditions==-1) {
209  cerr << "Number of conditions unknown ; call sec_member first" << endl ;
210  abort() ;
211  }
212 
214  conte = 0 ;
215  int pos_res = 0 ;
216  res = 0 ;
217  for (int i=0 ; i<neq_int ; i++) {
218  res.set(pos_res) = eq_int[i]->get_der() ;
219  pos_res ++ ;
220  }
221 
222  for (int i=0 ; i<neq ; i++) {
223  if ((is_var_double) || (eq[i]->take_into_account(zedom))|| (eq[i]->take_into_account(zedoms(0))) || (eq[i]->take_into_account(zedoms(1)))) {
224  eq[i]->export_der(conte, results.set_data(), res, pos_res) ;
225  }
226  else
227  {
228  pos_res += eq[i]->n_cond_tot ;
229  conte += eq[i]->n_ope ;
230  }
231  }
232 
233  return res ;
234 }
235 
237 
238  for (int i=0 ; i<nterm ; i++) {
239  int dom = term[i]->get_dom() ;
240  bool todo = false ;
241  for (int d=0 ; d<zedoms.get_size(0) ; d++)
242  if (zedoms(d)==dom)
243  todo = true ;
244  if (todo)
245  espace.get_domain(dom)->update_term_eq (term[i]) ;
246  else
247  term[i]->set_der_zero() ;
248  }
249 
250  for (int i=0 ; i<nterm_cst ; i++)
251  if (cst[i]->get_type_data() == TERM_T) {
252  int dom = cst[i]->get_dom() ;
253  bool todo = false ;
254  for (int d=0 ; d<zedoms.get_size(0) ; d++)
255  if (zedoms(d)==dom)
256  todo = true ;
257  if (todo)
258  espace.get_domain(dom)->update_term_eq (cst[i]) ;
259  else
260  cst[i]->set_der_zero() ;
261  }
262 }}
reference set(const Index &pos)
Read/write of an element.
Definition: array.hpp:186
int get_size(int i) const
Returns the size of a given dimension.
Definition: array.hpp:331
virtual void update_term_eq(Term_eq *so) const
Update the value of a field, after the shape of the Domain has been changed by the system.
Definition: domain.cpp:1749
virtual void affecte_tau_one_coef(Tensor &so, int dom, int cc, int &pos_cf) const
Sets at most one coefficient of a Tensor to 1.
Definition: domain.cpp:1575
virtual void update()
Updates the derived quantities (Christoffels etc..) This is done only for the ones that are needed,...
Definition: metric.cpp:99
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
virtual void affecte_coef_to_variable_domains(int &conte, int cc, Array< int > &doms) const
The variation of the functions describing the shape of the Domain are affected from the unknowns of t...
Definition: space.hpp:1401
MMPtr_array< Ope_def > def
Pointers on the definition (i.e. on the Ope_def that is needed to compute the result).
Array< double > sec_member()
Computes the second member of the Newton-Raphson equations.
Definition: solver.cpp:75
int nvar_double
Number of unknowns that are numbers (i.e. not fields)
MMPtr_array< Term_eq > results
Pointers on the residual of the various equations.
int nterm
Number of Term_eq corresponding to the unknown fields.
Array< double > do_col_J(int i)
Computes one column of the Jacobian.
Definition: solver.cpp:127
int mpi_proc_rank
Sylvain's stuff.
MMPtr_array< Term_eq > term_double
Pointers on the Term_eq corresponding to the unknowns that are numbers.
int nbr_unknowns
Number of unknowns (basically the number of coefficients of all the unknown fields,...
int ndom
Number of domains used.
int neq_int
Number of integral equations (i.e. which are doubles)
void update_terms_from_variable_domains(const Array< int > &zedoms)
Updates the variations of the Term_eq that comes from the fact that some Domains are variable (i....
Definition: solver.cpp:236
MMPtr_array< Equation > eq
Pointers onto the field equations.
virtual void compute_nbr_of_conditions()
Sylvain's stuff.
Definition: solver.cpp:51
MMPtr_array< Term_eq > cst
Pointers on the Term_eq coming from the constants passed by the user.
Array< double > check_equations()
Computes the residual of all the equations.
Definition: solver.cpp:26
void vars_to_terms()
Copies the various unknowns (doubles and Tensors) into their Term_eq counterparts.
int ndef
Number of definitions.
int nterm_cst
Number of Term_eq coming from the constants passed by the user.
int nbr_conditions
Total number of conditions (the number of coefficients of all the equations, once regularities are ta...
MMPtr_array< Term_eq > term
Pointers on the Term_eq corresponding to the unknown fields.
int neq
Number of field equations.
int dom_max
Highest domain number.
const Space & espace
Associated Space.
Array< double > do_JX(const Array< double > &xx)
Computes the product where is the Jacobian of the system.
Definition: solver.cpp:94
MMPtr_array< Eq_int > eq_int
Pointers onto the integral equations.
void xx_to_ders(const Array< double > &vder)
Sets the values the variation of the fields.
int dom_min
Smallest domain number.
Metric * met
Pointer on the associated Metric, if defined.
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
int get_n_comp() const
Returns the number of stored components.
Definition: tensor.hpp:514
virtual Array< int > indices(int pos) const
Gives the values of the indices corresponding to a location in the array used for storage of the comp...
Definition: tensor.hpp:484
Base_spectral & set_base()
Sets the basis of decomposition.
Definition: val_domain.hpp:126