KADATH
eq_matching_order_array.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 "ope_eq.hpp"
22 #include "scalar.hpp"
23 #include "tensor_impl.hpp"
24 namespace Kadath {
25 Eq_matching_order_array::Eq_matching_order_array(const Domain* zedom, int dd, int bb, int other_dd, int other_bb, const Array<int>& ord,
26  Ope_eq* lhs, Ope_eq* rhs, int nused, Array<int>** pused) :
27  Equation(zedom, dd, 2, nused, pused), bound(bb), other_dom(other_dd), other_bound(other_bb), order(ord) {
28  parts[0] = lhs ;
29  parts[1] = rhs ;
30 }
31 
33 }
34 
35 void Eq_matching_order_array::export_val(int& conte, Term_eq** residus, Array<double>& sec, int& pos_res) const {
36 
37  assert (residus[conte]->get_type_data()==TERM_T) ;
38  assert (residus[conte+1]->get_type_data()==TERM_T) ;
39 
40 
41  int start = pos_res ;
42  dom->export_tau_boundary_array (residus[conte]->get_val_t(), ndom, bound, order, sec, pos_res, *n_cond, n_cmp_used, p_cmp_used) ;
43 
44  Array<double> auxi (pos_res - start) ;
45  auxi = 0. ;
46  int zero = 0 ;
48  (residus[conte+1]->get_val_t(), other_dom, other_bound, order, auxi, zero, *n_cond, n_cmp_used, p_cmp_used) ;
49 
50  for (int i =start ; i<pos_res ; i++)
51  sec.set(i) -= auxi(i-start) ;
52  conte +=2 ;
53 }
54 
55 void Eq_matching_order_array::export_der(int& conte, Term_eq** residus, Array<double>& sec, int& pos_res) const {
56 
57  assert (residus[conte]->get_type_data()==TERM_T) ;
58  assert (residus[conte+1]->get_type_data()==TERM_T) ;
59  int start = pos_res ;
60  dom->export_tau_boundary_array (residus[conte]->get_der_t(), ndom, bound, order, sec, pos_res, *n_cond, n_cmp_used, p_cmp_used) ;
61  Array<double> auxi (pos_res - start) ;
62  auxi = 0. ;
63  int zero = 0 ;
65  (residus[conte+1]->get_der_t(), other_dom, other_bound, order, auxi, zero, *n_cond, n_cmp_used, p_cmp_used) ;
66  for (int i =start ; i<pos_res ; i++)
67  sec.set(i) -= auxi(i-start) ;
68  conte +=2 ;
69 }
70 
73 }
74 
75 
77  if ((target==ndom) || (target==other_dom))
78  return true ;
79  else
80  return false ;
81 }
82 
83 }
reference set(const Index &pos)
Read/write of an element.
Definition: array.hpp:186
Abstract class that implements the fonctionnalities common to all the type of domains.
Definition: space.hpp:60
virtual Array< int > nbr_conditions_boundary_array(const Tensor &eq, int dom, int bound, const Array< int > &order, int n_cmp=-1, Array< int > **p_cmp=0x0) const
Computes number of discretized equations associated with a given tensorial equation on a boundary.
Definition: domain.cpp:1522
virtual void export_tau_boundary_array(const Tensor &eq, int dom, int bound, const Array< int > &order, Array< double > &res, int &pos_res, const Array< int > &ncond, int n_cmp=-1, Array< int > **p_cmp=0x0) const
Exports all the residual equations corresponding to one tensorial one on a given boundary It makes us...
Definition: domain.cpp:1557
Array< int > do_nbr_conditions(const Tensor &tt) const override
Computes the number of conditions associated with the equation.
void export_val(int &, Term_eq **, Array< double > &, int &) const override
Generates the discretized errors, from the various Term_eq computed by the equation.
const Array< int > & order
Orders of the equation wrt each variable.
void export_der(int &, Term_eq **, Array< double > &, int &) const override
Generates the discretized variations, from the various Term_eq computed by the equation.
int other_dom
Number of the Domain on the other side of the boundary.
int other_bound
Name of the boundary as seen from the other domain.
Eq_matching_order_array(const Domain *dom, int nd, int bb, int oz_nd, int oz_bb, const Array< int > &ord, Ope_eq *ope, Ope_eq *oz_ope, int n_cmp=-1, Array< int > **p_cmp=nullptr)
Constructor.
bool take_into_account(int) const override
Check whether the variation of the residual has to be taken into account when computing a given colum...
Class implementing an equation.
const Domain * dom
Pointer on the Domain where the equation is defined.
Array< int > ** p_cmp_used
Array of pointer on the indices of the used components.
int ndom
Number of the domain.
int n_cmp_used
Number of components used (by default the same thing as n_comp).
Array< int > * n_cond
Number of discretized equations, component by component.
MMPtr_array< Ope_eq > parts
Array of pointers on the various terms.
Abstract class that describes the various operators that can appear in the equations.
Definition: ope_eq.hpp:32
const Domain * get_domain(int i) const
returns a pointer on the domain.
Definition: space.hpp:1385
Tensor handling.
Definition: tensor.hpp:149
const Space & get_space() const
Returns the Space.
Definition: tensor.hpp:499
This class is intended to describe the manage objects appearing in the equations.
Definition: term_eq.hpp:62
Tensor const & get_val_t() const
Definition: term_eq.hpp:355