20 #ifndef __SYSTEM_OF_EQS_HPP_
21 #define __SYSTEM_OF_EQS_HPP_
23 #include "headcpp.hpp"
28 #include "profiled_object.hpp"
31 #include "term_eq.hpp"
34 #include "matrice.hpp"
35 #include "list_comp.hpp"
40 constexpr
int VARMAX{1000};
42 #define GMRES_SPARSE 0
43 #define BICGSTAB_SPARSE 1
72 enum :
bool { DO_NOT_TRANSPOSE =
false, TRANSPOSE =
true};
118 MMPtr_array<Term_eq>
cst ;
127 MMPtr_array<Ope_def>
def ;
154 MMPtr_array<Equation>
eq ;
284 virtual void add_var (
const char* name,
double&
var) ;
296 virtual void add_cst (
const char* name,
double cst) ;
308 virtual void add_def (
const char* name) ;
314 virtual void add_def (
int dd,
const char* name) ;
346 bool isvar_double (
const char* target,
int& which)
const ;
355 bool isvar (
const char* target,
int& which,
int& valence,
char*& name_ind,
Array<int>*& type_ind)
const ;
364 bool iscst (
const char* target,
int& which,
int& valence,
char*& name_ind,
Array<int>*& type_ind)
const ;
374 bool isdef (
int dd,
const char* target,
int& which,
int& valence,
char*& name_ind,
Array<int>*& type_ind)
const ;
381 bool isdef_glob (
int dd,
const char* target,
int& which)
const ;
389 bool is_ope_minus (
const char* input,
char* output)
const ;
395 bool isdouble (
const char* input,
double& output)
const ;
402 bool ismet (
const char* input,
char*& name_ind,
int& type_ind)
const ;
407 bool ismet (
const char* input)
const ;
447 bool is_ope_bin (
const char* input,
char* p1,
char* p2,
char symb)
const ;
454 bool is_ope_uni (
const char* input,
char* p1,
const char* nameope)
const ;
462 bool is_ope_uni (
const char* input,
char* p1,
char* p2,
const char* nameope)
const ;
470 bool is_ope_deriv (
const char* input,
char* p1,
int& typeder,
char& nameind)
const ;
478 bool is_ope_deriv_flat (
const char* input,
char* p1,
int& typeder,
char& nameind)
const ;
493 bool is_ope_pow (
const char* input,
char* p1,
int& expo)
const ;
500 bool is_ope_partial (
const char* input,
char* p1,
char& nameind)
const ;
508 bool is_ope_der_var (
int dd,
const char* input,
char* p1,
int& which)
const ;
564 virtual void add_eq_bc (
int dom,
int bb,
const char*
eq,
int n_cmp = -1,
Array<int>** p_cmp=
nullptr) ;
748 virtual void add_eq_vel_pot (
int dom,
int order,
const char*
eq,
const char* const_part) ;
825 virtual void add_eq_mode (
int dom,
int bb,
const char*
eq,
const Index& pos_cf,
double val) ;
905 bool transpose = DO_NOT_TRANSPOSE);
910 bool transpose = DO_NOT_TRANSPOSE, std::vector<std::vector<std::size_t>> *dm =
nullptr);
919 template<Computational_model computational_model = default_computational_model>
920 bool do_newton(
double prec,
double &error);
944 template<Computational_model computational_model = default_computational_model>
959 template<Computational_model computational_model = default_computational_model>
965 template<Computational_model computational_model = default_computational_model>
971 template<Computational_model computational_model = default_computational_model>
977 template<Computational_model computational_model = default_computational_model>
983 template<Computational_model computational_model = default_computational_model>
989 template<Computational_model computational_model = default_computational_model>
990 void update_fields(
double lambda, vector<double>
const& old_var_double, vector<Tensor>
const& old_var_fields, vector<double>
const& p_var_double, vector<Tensor>
const& p_var_fields);
995 template<Computational_model computational_model = default_computational_model>
996 void compute_old_and_var(
Array<double> const& xx, vector<double>& old_var_double, vector<Tensor>& old_var_fields, vector<double>& p_var_double, vector<Tensor>& p_var_fields);
1001 template<Computational_model computational_model = default_computational_model>
1007 template<Computational_model computational_model = default_computational_model>
1013 template<Computational_model computational_model = default_computational_model>
1019 template<Computational_model computational_model = default_computational_model>
1040 friend class Metric_relax ;
1613 Eq_matching_exception(
const Domain*
dom,
int nd,
int bb,
int oz_nd,
int oz_bb,
Ope_eq* ope,
Ope_eq* oz_ope,
const Param& par,
Ope_eq* ope_exc,
int n_cmp = -1,
Array<int>** p_cmp=
nullptr) ;
1653 void sl_init_ (
int*,
int*,
int*);
1654 void blacs_gridexit_ (
int*);
1655 void blacs_gridinfo_ (
int*,
int*,
int*,
int*,
int*);
1656 int numroc_ (
int*,
int*,
int*,
int*,
int*);
1657 void descinit_ (
int*,
int*,
int*,
int*,
int* ,
int*,
int*,
int*,
int*,
int*);
1658 void pdgesv_ (
int*,
int*,
double*,
int*,
int*,
int*,
int*,
double*,
int*,
int*,
int*,
int*);
1659 int blacs_pnum_ (
int*,
int*,
int*);
1660 void Cpdgemr2d (
int,
int,
double*,
int,
int,
int*,
double*,
int,
int,
int*,
int);
1664 inline void split_two_d (
int target,
int& low,
int& up) {
1665 int start = int(sqrt(target));
1669 bool endloop = (low*up==target) ?
true :
false;
1678 up = int(target/low);
1687 bool System_of_eqs::do_newton<Computational_model::sequential>(
double prec,
double &error);
1689 bool System_of_eqs::do_newton<Computational_model::mpi_parallel>(
double,
double &);
1691 bool System_of_eqs::do_newton<Computational_model::gpu_sequential>(
double,
double &);
1693 bool System_of_eqs::do_newton<Computational_model::gpu_mpi_parallel>(
double,
double &);
1694 template<>
bool System_of_eqs::do_newton_with_linesearch<Computational_model::mpi_parallel>(
double ,
double &,
int ,
double );
1698 template<Computational_model computational_model>
1701 cerr <<
"Error: System_of_eq::do_newton is not implemented for the "
1702 << computational_model_name(computational_model)
1703 <<
" computational model." << endl;
1707 template<Computational_model computational_model>
1710 cerr <<
"Error: System_of_eq::do_newton is not implemented for the "
1711 << computational_model_name(computational_model)
1712 <<
" computational model." << endl;
1716 template<Computational_model computational_model>
1717 void check_size_VS_unknowns(
int n)
1719 cerr <<
"Error: System_of_eq::check_size_VS_unknowns is not implemented for the "
1720 << computational_model_name(computational_model)
1721 <<
" computational model." << endl;
1724 template<Computational_model computational_model>
1725 void check_bsize(
int bsize)
1727 cerr <<
"Error: System_of_eq::check_bsize is not implemented for the "
1728 << computational_model_name(computational_model)
1729 <<
" computational model." << endl;
1732 template<Computational_model computational_model>
1733 void compute_matloc(Array<double>& matloc_in,
int nn,
int bsize)
1735 cerr <<
"Error: System_of_eq::compute_matloc is not implemented for the "
1736 << computational_model_name(computational_model)
1737 <<
" computational model." << endl;
1740 template<Computational_model computational_model>
1741 void translate_second_member(Array<double>& secloc, Array<double>
const& second,
int nn,
int bsize,
int nprow,
int myrow,
int mycol)
1743 cerr <<
"Error: System_of_eq::translate_second_member is not implemented for the "
1744 << computational_model_name(computational_model)
1745 <<
" computational model." << endl;
1748 template<Computational_model computational_model>
1749 void get_global_solution(Array<double>& auxi, Array<double>
const& secloc,
int nn,
int bsize,
int nprow,
int myrow,
int mycol)
1751 cerr <<
"Error: System_of_eq::get_global_solution is not implemented for the "
1752 << computational_model_name(computational_model)
1753 <<
" computational model." << endl;
1756 template<Computational_model computational_model>
1757 void update_fields(
double lambda, vector<double>
const& old_var_double, vector<Tensor>
const& old_var_fields, vector<double>
const& p_var_double, vector<Tensor>
const& p_var_fields)
1759 cerr <<
"Error: System_of_eq::update_fields is not implemented for the "
1760 << computational_model_name(computational_model)
1761 <<
" computational model." << endl;
1764 template<Computational_model computational_model>
1765 void compute_old_and_var(Array<double>
const& xx, vector<double>& old_var_double, vector<Tensor>& old_var_fields, vector<double>& p_var_double, vector<Tensor>& p_var_fields)
1767 cerr <<
"Error: System_of_eq::compute_old_and_var is not implemented for the "
1768 << computational_model_name(computational_model)
1769 <<
" computational model." << endl;
1772 template<Computational_model computational_model>
1773 double compute_f(Array<double>
const& second)
1775 cerr <<
"Error: System_of_eq::compute_f is not implemented for the "
1776 << computational_model_name(computational_model)
1777 <<
" computational model." << endl;
1781 template<Computational_model computational_model>
1782 void compute_p(Array<double>& xx, Array<double>
const& second,
int nn)
1784 cerr <<
"Error: System_of_eq::compute_p is not implemented for the "
1785 << computational_model_name(computational_model)
1786 <<
" computational model." << endl;
1789 template<Computational_model computational_model>
1790 void check_positive(
double delta)
1792 cerr <<
"Error: System_of_eq::check_positive is not implemented for the "
1793 << computational_model_name(computational_model)
1794 <<
" computational model." << endl;
1797 template<Computational_model computational_model>
1798 void check_negative(
double delta)
1800 cerr <<
"Error: System_of_eq::check_negative is not implemented for the "
1801 << computational_model_name(computational_model)
1802 <<
" computational model." << endl;
Class for storing the basis of decompositions of a field.
Abstract class that implements the fonctionnalities common to all the type of domains.
Class for enforcing boundary condition.
virtual ~Eq_bc_exception()
Destructor.
Eq_bc_exception(const Domain *dom, int nd, int bound, Ope_eq *ope, Ope_eq *ope_constant)
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.
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...
void export_der(int &, Term_eq **, Array< double > &, int &) const override
Generates the discretized variations, from the various Term_eq computed by the equation.
Array< int > do_nbr_conditions(const Tensor &tt) const override
Computes the number of conditions associated with the equation.
int bound
Order of the equation.
Class for a boundary condition.
Eq_bc_order_array(const Domain *dom, int nd, int bb, const Array< int > &ord, Ope_eq *ope, int n_cmp=-1, Array< int > **p_cmp=nullptr)
Constructor.
virtual ~Eq_bc_order_array()
Destructor.
const Array< int > & order
Orders of the equation wrt each variable.
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 > do_nbr_conditions(const Tensor &tt) const override
Computes the number of conditions associated with the equation.
void export_der(int &, Term_eq **, Array< double > &, int &) const override
Generates the discretized variations, from the various Term_eq computed by 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.
Class for an equation representing a boundary condition on some surface.
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 > 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.
Eq_bc(const Domain *dom, int nd, int bb, Ope_eq *ope, int n_cmp=-1, Array< int > **p_cmp=nullptr)
Constructor.
virtual ~Eq_bc()
Destructor.
void export_der(int &, Term_eq **, Array< double > &, int &) const override
Generates the discretized variations, from the various Term_eq computed by the equation.
Equation for describing a first integral equation (i.e.
Array< int > do_nbr_conditions(const Tensor &tt) const override
Computes the number of conditions associated with the equation.
Eq_first_integral(const System_of_eqs *syst, const Domain *dom, int dommin, int dommax, const char *integ_part, const char *const_part)
Constructor.
int dom_min
Index of the first Domain.
void export_val(int &, Term_eq **, Array< double > &, int &) const override
Generates the discretized errors, from the various Term_eq computed by the equation.
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 dom_max
Index of the last Domain.
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...
virtual ~Eq_first_integral()
Destructor.
Class for a zeroth order equation in a Domain.
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 > do_nbr_conditions(const Tensor &tt) const override
Computes the number of conditions associated with the equation.
virtual ~Eq_full()
Destructor.
Eq_full(const Domain *dom, int nd, Ope_eq *ope, 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.
void export_der(int &, Term_eq **, Array< double > &, int &) const override
Generates the discretized variations, from the various Term_eq computed by the equation.
Class for bulk equations that are solved stricly inside a given domain.
void export_der(int &, Term_eq **, Array< double > &, int &) const override
Generates the discretized variations, from the various Term_eq computed by the equation.
Array< int > do_nbr_conditions(const Tensor &tt) const override
Computes the number of conditions associated with 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...
virtual ~Eq_inside()
Destructor.
Eq_inside(const Domain *dom, int nd, Ope_eq *ope, 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.
Class implementing an integral equation.
Ope_eq ** parts
Array of pointers on the various terms.
double get_der() const
Return the variation of the equation.
int n_ope
Number of terms.
Eq_int(int nop)
Constructor just sets n_ope.
void set_part(int i, Ope_eq *ope)
Sets one of the Ope_eq needed by the equation.
double get_val() const
Return the value of the equation.
Equation for a matching condition, except for one coefficient where an alternative condition is enfor...
virtual ~Eq_matching_exception()
Destructor.
Eq_matching_exception(const Domain *dom, int nd, int bb, int oz_nd, int oz_bb, Ope_eq *ope, Ope_eq *oz_ope, const Param &par, Ope_eq *ope_exc, int n_cmp=-1, Array< int > **p_cmp=nullptr)
Constructor.
Array< int > do_nbr_conditions(const Tensor &tt) const override
Computes the number of conditions associated with the equation.
const Param & parameters
Parameters needed for describing the exception.
void export_val(int &, Term_eq **, Array< double > &, int &) const override
Generates the discretized errors, 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...
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.
Class for an equation representing the matching of quantities accross a boundary using the "import" r...
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 for an equation representing the matching of quantities accross a boundary.
void export_val(int &, Term_eq **, Array< double > &, int &) const override
Generates the discretized errors, from the various Term_eq computed by the equation.
Array< int > do_nbr_conditions(const Tensor &tt) const override
Computes the number of conditions associated with the equation.
virtual ~Eq_matching_non_std()
Destructor.
Array< int > other_bounds
Names of the boundary, as seen in the other domains.
void do_which_points(const Base_spectral &base, int start)
Computes the collocation points used.
Array< int > other_doms
Array containing the number of the domains being on the other side of the surface.
int bound
Name of the boundary in the domain of the equation.
Eq_matching_non_std(const Domain *dom, int nd, int bb, const Array< int > &ozers, int n_cmp=-1, Array< int > **p_cmp=nullptr)
Constructor.
void apply(int &, Term_eq **) override
Computes the terms involved in computing the residual of the equations.
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...
Index ** which_points
Lists the collocation points on the boundary (probably...)
Class for an equation representing the matching of quantities accross a boundary.
void export_der(int &, Term_eq **, Array< double > &, int &) const override
Generates the discretized variations, from the various Term_eq computed by the equation.
Eq_matching_one_side(const Domain *dom, int nd, int bb, int oz_nd, int oz_bb, Ope_eq *ope, Ope_eq *oz_ope, int n_cmp=-1, Array< int > **p_cmp=nullptr)
Constructor.
virtual ~Eq_matching_one_side()
Destructor.
void export_val(int &, Term_eq **, Array< double > &, int &) const override
Generates the discretized errors, 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...
int other_dom
Number of the other domain.
Array< int > do_nbr_conditions(const Tensor &tt) const override
Computes the number of conditions associated with the equation.
int other_bound
Name of the boundary in the other domain.
int bound
Name of the boundary in the domain of the equation.
Class for a matching condition.
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.
virtual ~Eq_matching_order_array()
Destructor.
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 for an equation representing the matching of quantities accross a boundary.
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...
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 other_dom
Number of the other domain.
virtual ~Eq_matching()
Destructor.
Array< int > do_nbr_conditions(const Tensor &tt) const override
Computes the number of conditions associated with the equation.
Eq_matching(const Domain *dom, int nd, int bb, int oz_nd, int oz_bb, Ope_eq *ope, Ope_eq *oz_ope, int n_cmp=-1, Array< int > **p_cmp=nullptr)
Constructor.
int bound
Name of the boundary in the domain of the equation.
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_bound
Name of the boundary in the other domain.
Class for a first order equation in a Domain.
virtual ~Eq_one_side()
Destructor.
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...
void export_val(int &, Term_eq **, Array< double > &, int &) const override
Generates the discretized errors, from the various Term_eq computed by the equation.
void export_der(int &, Term_eq **, Array< double > &, int &) const override
Generates the discretized variations, from the various Term_eq computed by the equation.
Array< int > do_nbr_conditions(const Tensor &tt) const override
Computes the number of conditions associated with the equation.
Eq_one_side(const Domain *dom, int nd, Ope_eq *ope, int n_cmp=-1, Array< int > **p_cmp=nullptr)
Constructor.
Class for an equation in a Domain which order is passed, for each variable.
Array< int > do_nbr_conditions(const Tensor &tt) const override
Computes the number of conditions associated with the equation.
Eq_order_array(const Domain *dom, int nd, const Array< int > &ord, Ope_eq *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...
virtual ~Eq_order_array()
Destructor.
void export_val(int &, Term_eq **, Array< double > &, int &) const override
Generates the discretized errors, from the various Term_eq computed by the equation.
void export_der(int &, Term_eq **, Array< double > &, int &) const override
Generates the discretized variations, from the various Term_eq computed by the equation.
const Array< int > & order
Orders of the equation wrt each variable.
Class for bulk equation which order is passed as a parameter.
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...
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 order
Order of the equation.
virtual ~Eq_order()
Destructor.
void export_val(int &, Term_eq **, Array< double > &, int &) const override
Generates the discretized errors, from the various Term_eq computed by the equation.
Array< int > do_nbr_conditions(const Tensor &tt) const override
Computes the number of conditions associated with the equation.
Eq_order(const Domain *dom, int nd, int ord, Ope_eq *ope, int n_cmp=-1, Array< int > **p_cmp=nullptr)
Constructor.
Class for the velocity potential in irrotational binray neutron stars.
int order
Order of 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...
virtual ~Eq_vel_pot()
Destructor.
Eq_vel_pot(const Domain *dom, int nd, int ord, Ope_eq *ope, Ope_eq *ope_constant)
Constructor.
void export_der(int &, Term_eq **, Array< double > &, int &) const override
Generates the discretized variations, from the various Term_eq computed by 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.
Array< int > do_nbr_conditions(const Tensor &tt) const override
Computes the number of conditions associated with the equation.
Class implementing an equation.
bool called
Indicator checking whther the result has been computed already once.
const Domain * dom
Pointer on the Domain where the equation is defined.
virtual Array< int > do_nbr_conditions(const Tensor &tt) const =0
Computes the number of conditions associated with the equation.
Array< int > ** p_cmp_used
Array of pointer on the indices of the used components.
virtual void apply(int &conte, Term_eq **res)
Computes the terms involved in computing the residual of the equations.
virtual ~Equation()
Destructor.
int ndom
Number of the domain.
int n_cmp_used
Number of components used (by default the same thing as n_comp).
int n_ope
Number of terms involved in the equation (one for bulk, two or more fot matching.....
virtual void export_val(int &conte, Term_eq **residuals, Array< double > &sec, int &pos_sec) const =0
Generates the discretized errors, from the various Term_eq computed by the equation.
int n_comp
Number of components of the residual (1 for a scalar, 6 for a symmetric rank-2 tensor etc).
Array< int > * n_cond
Number of discretized equations, component by component.
int get_n_cond_tot() const
virtual void export_der(int &conte, Term_eq **residuals, Array< double > &sec, int &pos_sec) const =0
Generates the discretized variations, from the various Term_eq computed by the equation.
virtual bool take_into_account(int target) const =0
Check whether the variation of the residual has to be taken into account when computing a given colum...
MMPtr_array< Ope_eq > parts
Array of pointers on the various terms.
Equation(const Domain *dom, int nd, int nope, int n_cmp=-1, Array< int > **p_cmp=nullptr)
Constructor.
int n_cond_tot
Total number of discretized equations (essentially the number of all coefficients of the residual).
Class that gives the position inside a multi-dimensional Array.
Class for storing a list of tensorial components.
Class to manage asymptotically anti de Sitter metrics.
Class to manage anti de Sitter metrics.
Class to deal with a metric with a conformaly flat metric.
Class to deal with a metric with a conformal decomposition (mainly used for AADS spacetimes) The true...
Class to deal with a metric which determinant is 1.
Class to deal with arbitrary type of metric but that is constant.
Class to deal with a conformal metric in the Dirac gauge.
Class to deal with a conformal metric in the Dirac gauge.
Class that deals with flat metric assuming a axisymmetric setting (nothing depends on ).
Class that deals with flat metric.
Class to deal with arbitrary type of metric.
Class to deal with a metric independant of with a conformal decomposition and constant.
Class to deal with a metric independant of with a conformal decomposition (mainly used for AADS spac...
Class to deal with a metric independant of and constant.
Class to deal with a metric independant of .
Purely abstract class for metric handling.
Operator for a global definition (i.e.
Abstract class that describes the various operators that can appear in the equations.
The class Point is used to store the coordinates of a point.
The Space_spheric_adapted class fills the space with one shell adapted on the inside,...
Spacetime intended for binary black hole configurations in full general relativity (see constructor f...
Spacetime intended for binary black hole configurations (see constructor for details about the domain...
Spacetime intended for fake binary neutron star configurations (see constructor for details about the...
Spacetime intended for binary neutron stars configurations (see constructor for details about the dom...
The Space_bispheric class fills the space with bispherical coordinates (see the constructor for more ...
The Space_critic ; set the space with two critic domains, separated at .
The Space_polar_adapted class fills the space with one polar nucleus, one polar shell adapted on the ...
The Space_polar_periodic class fills the space with one polar nucleus and several polar shells,...
The Space_polar class fills the space with one polar nucleus, several polar shells and a compactified...
The Space_spheric_adapted class fills the space with one nucleus, one shell adapted on the outside,...
The Space_spheric class fills the space with one nucleus, several shells and a compactified domain,...
The Space class is an ensemble of domains describing the whole space of the computation.
Class used to describe and solve a system of equations.
MMPtr_array< Ope_def > def
Pointers on the definition (i.e. on the Ope_def that is needed to compute the result).
bool isvar(const char *target, int &which, int &valence, char *&name_ind, Array< int > *&type_ind) const
Check if a string is an unknown field.
bool do_newton_with_linesearch(double precision, double &error, int ntrymax=10, double stepmax=1.0)
Does one step of the Newton-Raphson iteration with a linesearch algorithm.
virtual void vars_to_terms_impl()
Sylvain's stuff.
MMPtr_array< Param > paruser
Parameters used by the user defined operators (single argument).
virtual void compute_matrix_adjacent(Array< double > &matrix, int n, int first_col=0, int n_col=ALL_COLUMNS, int num_proc=1, bool transpose=DO_NOT_TRANSPOSE, std::vector< std::vector< std::size_t >> *dm=nullptr)
Sylvain's stuff.
Array< double > sec_member()
Computes the second member of the Newton-Raphson equations.
virtual void add_eq_inside(int dom, const char *eq, int n_cmp=-1, Array< int > **p_cmp=nullptr)
Addition of an equation to be solved inside a domain (assumed to be second order).
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 nopeuser_bin
Number of operators defined by the user (with two arguments).
bool isvar_double(const char *target, int &which) const
Check if a string is an unknown (number).
virtual void add_var(const char *name, double &var)
Addition of a variable (number case)
const Metric * get_met() const
Returns a pointer on the Metric.
int nterm
Number of Term_eq corresponding to the unknown fields.
Term_eq * give_term_double(int which, int dd) const
Returns a pointer on a Term_eq corresponding to an unknown number.
int get_mpi_world_size() const
Returns the total number of MPI processes.
virtual void add_eq_full(int dom, const char *eq, int n_cmp=-1, Array< int > **p_cmp=nullptr)
Addition of an equation to be solved inside a domain (assumed to be zeroth order i....
void init_proc_data()
Sylvain's stuff.
const Space & get_space() const
Returns the space.
Array< double > do_col_J(int i)
Computes one column of the Jacobian.
virtual void add_eq_order(int dom, int order, const char *eq, int n_cmp=-1, Array< int > **p_cmp=nullptr)
Addition of an equation to be solved inside a domain (of arbitrary order).
bool is_ope_der_var(int dd, const char *input, char *p1, int &which) const
Checks if a string represents the derivative wrt a numerical coordinate of a given Domain (like "a,...
~System_of_eqs() override
Destructor.
MMPtr_array< double > var_double
Pointer on the unknowns that are numbers (i.e. not fields)
bool isricci_tensor(const char *input, char *&name_ind, Array< int > *&type_ind) const
Checks if a string represents the Ricci tensor.
MMPtr_array< Ope_def_global > def_glob
Pointers on the global definitions.
Array< int > assoc_var
Array giving the correspondance with the var pointers.
virtual void add_def(const char *name)
Addition of a definition.
void check_positive(double delta)
Tests the positivity of   Used by do_newton_with_linesearch ; only implemented in parallel version.
static constexpr int ALL_COLUMNS
Dummy variable for the purpose of better readability.
std::ostream * output_stream
Default output stream for log messages.
virtual void add_cst(const char *name, double cst)
Addition of a constant (number case)
MMPtr_array< char > names_def_glob
Names of the global definitions.
bool is_ope_minus(const char *input, char *output) const
Checks if a string contains the operator minus.
void add_eq_bc_exception(int dom, int bound, const char *eq, const char *const_part)
Addition of an boundary equation with an exception for .
MMPtr_array< char > names_var
Names of the unknown fields.
bool is_ope_uni(const char *input, char *p1, const char *nameope) const
Checks if a string represents an operator of the type "ope(a)".
int mpi_proc_rank
Sylvain's stuff.
Term_eq * give_cst_hard(double xx, int dd) const
Returns a pointer on a Term_eq corresponding to a constant generated on the fly.
bool is_ope_deriv_background(const char *input, char *p1, int &typeder, char &nameind) const
Checks if a string represents the covariant derivative wrt a background metric.
virtual void add_eq_matching_import(int dom, int bb, const char *eq, int n_cmp=-1, Array< int > **p_cmp=nullptr)
Addition of an equation describing a matching condition between domains using the ("import" setting) ...
virtual void add_eq_mode(int dom, int bb, const char *eq, const Index &pos_cf, double val)
Addition of an equation prescribing the value of one coefficient of a scalar field,...
virtual void add_eq_matching(int dom, int bb, const char *eq, int n_cmp=-1, Array< int > **p_cmp=nullptr)
Addition of an equation describing a matching condition between two domains (standard setting)
void compute_matloc(Array< double > &matloc_in, int nn, int bsize)
Computes the local part of the Jacobian Used by do_newton_with_linesearch ; only implemented in paral...
MMPtr_array< Term_eq > cst_hard
Pointers on the Term_eq coming from the constants generated on the fly (when encoutering things like ...
Term_eq * give_term(int which, int dd) const
Returns a pointer on a Term_eq corresponding to an unknown field.
Ope_def * give_def(int i) const
Returns a pointer on a definition (better to use give_val_def if one wants to access the result of so...
bool isdouble(const char *input, double &output) const
Checks if a string is a double.
MMPtr_array< Term_eq > term_double
Pointers on the Term_eq corresponding to the unknowns that are numbers.
bool ischristo(const char *input, char *&name_ind, Array< int > *&type_ind) const
Checks if a string represents the Christoffel symbols.
int nbr_unknowns
Number of unknowns (basically the number of coefficients of all the unknown fields,...
MMPtr_array< char > names_opeuser_bin
Names of the user defined operators (with two arguments).
int ndom
Number of domains used.
int ndef_glob
Number of global definitions (the one that require the knowledge of the whole space to give the resul...
bool isriemann(const char *input, char *&name_ind, Array< int > *&type_ind) const
Checks if a string represents the Riemann tensor.
Tensor give_val_def(const char *name) const
Gives the result of a definition.
Ope_def_global * give_def_glob(int i) const
Returns a pointer on a global definition.
bool is_ope_pow(const char *input, char *p1, int &expo) const
Checks if a string represents the power of something (like "a^2").
int neq_int
Number of integral equations (i.e. which are doubles)
int get_mpi_proc_rank() const
Returns the MPI rank.
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....
bool iscst(const char *target, int &which, int &valence, char *&name_ind, Array< int > *&type_ind) const
Check if a string is a constant (can required indices manipulation and/or inner contraction).
MMPtr_array< char > names_opeuser
Names of the user defined operators (single argument).
virtual void add_eq_bc(int dom, int bb, const char *eq, int n_cmp=-1, Array< int > **p_cmp=nullptr)
Addition of an equation describing a boundary condition.
bool is_ope_partial(const char *input, char *p1, char &nameind) const
Checks if a string represents the partial derivative (like "partial_i a")
MMPtr_array< Equation > eq
Pointers onto the field equations.
bool is_ope_bin(const char *input, char *p1, char *p2, char symb) const
Checks if a string represents an operator of the type "a + b".
virtual void add_eq_point(int dom, const char *eq, const Point &MM)
Addition of an equation saying that the value of a field must be zero at one point (arbitrary).
bool ismet(const char *input, char *&name_ind, int &type_ind) const
Checks if a string is a metric.
MMPtr_array< Index > which_coef
Stores the "true" coefficients on some boundaries (probably deprecated).
int nopeuser
Number of operators defined by the user (single argument).
MMPtr_array< Param > paruser_bin
Parameters used by the user defined operators (with two arguments).
MMPtr_array< Tensor > var
Pointer on the unknown fields.
virtual void compute_nbr_of_conditions()
Sylvain's stuff.
bool do_newton(double prec, double &error)
Does one step of the Newton-Raphson iteration.
virtual void add_eq_matching_one_side(int dom, int bb, const char *eq, int n_cmp=-1, Array< int > **p_cmp=nullptr)
Addition of an equation describing a matching condition between two domains (specialized function for...
std::ostream & get_output_stream() const
Returns the default output stream reference.
Array< int > assoc_var_double
Array giving the correspondance with the var_double pointers.
bool isdef_glob(int dd, const char *target, int &which) const
Check if a string is a global definition.
unsigned get_niter() const
Returns the current iteration number.
MMPtr_array< Term_eq > cst
Pointers on the Term_eq coming from the constants passed by the user.
void compute_p(Array< double > &xx, Array< double > const &second, int nn)
Inner routine for the linesearch algorithm Used by do_newton_with_linesearch ; only implemented in pa...
void xx_to_vars(const Array< double > &val, int &conte)
Sets the values the of the fields.
Array< double > check_equations()
Computes the residual of all the equations.
void check_negative(double delta)
Tests the positivity of   Used by do_newton_with_linesearch ; only implemented in parallel version.
Ope_eq * give_ope(int dom, const char *name, int bb=0) const
Function that reads a string and returns a pointer on the generated Ope_eq.
void check_size_VS_unknowns(int n)
Tests the value of the number of unknowns.
void get_global_solution(Array< double > &auxi, Array< double > const &secloc, int nn, int bsize, int nprow, int myrow, int mycol)
Solves the linear problem in Newton-Raphson.
MMPtr_array< char > names_cst
Names of the constants passed by the user.
int get_nbr_unknowns() const
Returns the number of unknowns.
void compute_old_and_var(Array< double > const &xx, vector< double > &old_var_double, vector< Tensor > &old_var_fields, vector< double > &p_var_double, vector< Tensor > &p_var_fields)
Update the fields when some domains have been modified.
virtual void add_ope(const char *name, Term_eq(*pope)(const Term_eq &, Param *), Param *par)
Addition of a user defined operator (one argument version)
virtual void add_eq_val(int dom, const char *eq, const Index &pos)
Addition of an equation saying that the value of a field must be zero at one collocation point.
bool isricci_scalar(const char *input, char *&name_ind, Array< int > *&type_ind) const
Checks if a string represents the Ricci scalar.
static std::size_t default_block_size
Defines the sub-matrix size in the scalapack 2D cyclic block decomposition.
void vars_to_terms()
Copies the various unknowns (doubles and Tensors) into their Term_eq counterparts.
unsigned niter
Counter toward the number of times the do_newton method has been called.
int ndef
Number of definitions.
virtual void add_eq_vel_pot(int dom, int order, const char *eq, const char *const_part)
Addition of an equation for the velocity potential of irrotational binaries.
double compute_f(Array< double > const &second)
Inner routine for the linesearch algorithm Used by do_newton_with_linesearch ; only implemented in pa...
MMPtr_array< char > names_var_double
Names of the unknowns that are numbers (i.e. not fields)
int nterm_cst
Number of Term_eq coming from the constants passed by the user.
void do_arnoldi(int n, Array< double > &qi, Matrice &Hmat)
Performs the Arnoldi algorithm (under developpement)
virtual void add_def_global(const char *name)
Addition of a global definition.
virtual void add_eq_matching_exception(int dom, int bb, const char *eq, const Param &par, const char *eq_exception, int n_cmp=-1, Array< int > **p_cmp=nullptr)
Addition of a matching condition, except for one coefficient where an alternative condition is enforc...
char * name_met
Name by which the metric is recognized.
int get_dom_max() const
Returns the highest index of the domains.
bool is_ope_deriv(const char *input, char *p1, int &typeder, char &nameind) const
Checks if a string represents the covariant derivative.
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.
System_of_eqs(const Space &so, int i)
Constructor, nothing is done.
bool isdef(int dd, const char *target, int &which, int &valence, char *&name_ind, Array< int > *&type_ind) const
Check if a string is a definition (can required indices manipulation and/or inner contraction).
int neq
Number of field equations.
int dom_max
Highest domain number.
void newton_update_vars(Array< double > const &xx)
Update the values of var and var_double from the solution of the linear system of the last Newton ite...
const Space & espace
Associated Space.
void check_bsize(int bsize)
Tests the not too many processors are used.
MMPtr_array< char > names_def
Names of the definitions.
int ncst
Number of constants passed by the user.
int ncst_hard
Number of constants generated on the fly (when encoutering things like "2.2" etc.....
virtual void add_eq_one_side(int dom, const char *eq, int n_cmp=-1, Array< int > **p_cmp=nullptr)
Addition of an equation to be solved inside a domain (assumed to be first order).
Array< double > val_cst_hard
Values of the constants generated on the fly (when encoutering things like "2.2" etc....
int get_dom_min() const
Returns the smallest index of the domains.
int nterm_double
Number of Term_eq corresponding to the unknowns that are numbers.
Array< double > do_JX(const Array< double > &xx)
Computes the product where is the Jacobian of the system.
MMPtr_array< Eq_int > eq_int
Pointers onto the integral equations.
Term_eq * give_cst(int which, int dd) const
Returns a pointer on a Term_eq corresponding to a constant.
int mpi_world_size
Sylvain's stuff.
void xx_to_ders(const Array< double > &vder)
Sets the values the variation of the fields.
int dom_min
Smallest domain number.
bool is_ope_deriv_flat(const char *input, char *p1, int &typeder, char &nameind) const
Checks if a string represents the flat covariant derivative.
virtual void add_eq_val_mode(int dom, const char *eq, const Index &pos_cf, double val)
Addition of an equation prescribing the value of one coefficient of a scalar field.
Metric * met
Pointer on the associated Metric, if defined.
virtual void compute_matrix_cyclic(Array< double > &matrix, int n, int first_col=0, int n_col=ALL_COLUMNS, int num_proc=1, bool transpose=DO_NOT_TRANSPOSE)
Compute some columns of the jacobian of the system of equations.
virtual void add_eq_matching_non_std(int dom, int bb, const char *eq, int n_cmp=-1, Array< int > **p_cmp=nullptr)
Addition of an equation describing a matching condition between domains.
void update_gmres(const Array< double > &)
Perfoems one step of the GMRES method (under developpement)
void update_fields(double lambda, vector< double > const &old_var_double, vector< Tensor > const &old_var_fields, vector< double > const &p_var_double, vector< Tensor > const &p_var_fields)
Update the fields after a Newton-Raphson iteration.
System_of_eqs(const System_of_eqs &)=delete
Constructor by copy.
int nvar
Number of unknown fields.
Term_eq(* opeuser_bin[VARMAX])(const Term_eq &, const Term_eq &, Param *)
Pointers on the functions used by the user defined operators (with two arguments).
System_of_eqs(const Space &so)
Standard constructor nothing is done.
void set_output_stream(std::ostream &new_output_stream)
Sets a new output stream reference.
Term_eq(* opeuser[VARMAX])(const Term_eq &, Param *)
Pointers on the functions used by the user defined operators (single argument).
Output_data current_output_data
Data related to the last newton iterations.
static constexpr std::size_t nb_core_per_node
Sylvain's stuff.
virtual void add_eq_first_integral(int dom_min, int dom_max, const char *integ_part, const char *const_part)
Addition of an equation representing a first integral.
void translate_second_member(Array< double > &secloc, Array< double > const &second, int nn, int bsize, int nprow, int myrow, int mycol)
Distributes the second member of Newton-Raphson accross the various processors.
int get_nbr_conditions() const
Returns the number of conditions.
This class is intended to describe the manage objects appearing in the equations.
Duration t_newton_update
Sylvain's stuff.
unsigned n_iter
Sylvain's stuff.
int problem_size
Sylvain's stuff.
Duration t_load_matrix
Sylvain's stuff.
Duration t_inv_matrix
Sylvain's stuff.
double current_error
Sylvain's stuff.
Duration t_trans_matrix
Sylvain's stuff.