KADATH
eq_matching_import.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_import::Eq_matching_import(const Domain* zedom, int dd, int bb, Ope_eq* so, const Array<int>& ozers, int nused, Array<int>** pused) : Equation(zedom, dd, 1, nused, pused), bound(bb),
26  other_doms(ozers.get_size(1)), other_bounds(ozers.get_size(1)) {
27  parts[0] = so ;
28  for (int i=0 ; i<other_doms.get_size(0) ; i++) {
29  other_doms.set(i) = ozers(0, i) ;
30  other_bounds.set(i) = ozers(1, i) ;
31  }
32 }
33 
35 }
36 
37 void Eq_matching_import::export_val(int& conte, Term_eq** residus, Array<double>& sec, int& pos_res) const {
38 
39  assert (residus[conte]->get_type_data()==TERM_T) ;
40  dom->export_tau_boundary (residus[conte]->get_val_t(), ndom, bound, sec, pos_res, *n_cond, n_cmp_used, p_cmp_used) ;
41  conte ++ ;
42 }
43 
44 void Eq_matching_import::export_der(int& conte, Term_eq** residus, Array<double>& sec, int& pos_res) const {
45  assert (residus[conte]->get_type_data()==TERM_T) ;
46  dom->export_tau_boundary (residus[conte]->get_der_t(), ndom, bound, sec, pos_res, *n_cond, n_cmp_used, p_cmp_used) ;
47  conte ++ ;
48 }
49 
52 }
53 
54 bool Eq_matching_import::take_into_account (int target) const {
55 
56  bool res = (target==ndom) ? true : false ;
57  for (int i=0 ; i<other_doms.get_size(0) ; i++)
58  if (target==other_doms(i))
59  res = true ;
60  return res ;
61 }
62 }
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
Abstract class that implements the fonctionnalities common to all the type of domains.
Definition: space.hpp:60
virtual void export_tau_boundary(const Tensor &eq, int dom, int bound, 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 a tensorial one on a given boundary It makes use ...
Definition: domain.cpp:1545
virtual Array< int > nbr_conditions_boundary(const Tensor &eq, int dom, int bound, 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:1516
void export_der(int &, Term_eq **, Array< double > &, int &) const override
Generates the discretized variations, from the various Term_eq computed by the equation.
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...
Array< int > other_doms
Array containing the number of the domains being on the other side of the surface.
Array< int > other_bounds
Names of the boundary, as seen in the other domains.
Eq_matching_import(const Domain *dom, int nd, int bb, Ope_eq *eq, const Array< int > &ozers, int n_cmp=-1, Array< int > **p_cmp=nullptr)
Constructor.
void export_val(int &, Term_eq **, Array< double > &, int &) const override
Generates the discretized errors, from the various Term_eq computed by the equation.
int bound
Name of the boundary in the domain of the equation.
Array< int > do_nbr_conditions(const Tensor &tt) const override
Computes the number of conditions associated with the equation.
virtual ~Eq_matching_import()
Destructor.
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
Tensor handling.
Definition: tensor.hpp:149
This class is intended to describe the manage objects appearing in the equations.
Definition: term_eq.hpp:62