33 #ifndef BGEOT_POLY_H__
34 #define BGEOT_POLY_H__
65 std::vector<short_type> v;
71 typedef std::vector<short_type>::iterator iterator;
72 typedef std::vector<short_type>::const_iterator const_iterator;
73 typedef std::vector<short_type>::reverse_iterator reverse_iterator;
74 typedef std::vector<short_type>::const_reverse_iterator const_reverse_iterator;
78 iterator begin() { dirty();
return v.begin(); }
79 const_iterator begin()
const {
return v.begin(); }
80 iterator end() { dirty();
return v.end(); }
81 const_iterator end()
const {
return v.end(); }
83 reverse_iterator rbegin() { dirty();
return v.rbegin(); }
84 const_reverse_iterator rbegin()
const {
return v.rbegin(); }
85 reverse_iterator rend() { dirty();
return v.rend(); }
86 const_reverse_iterator rend()
const {
return v.rend(); }
88 size_type size()
const {
return v.size(); }
180 template<
typename T>
class polynomial :
public std::vector<T> {
187 typedef typename std::vector<T>::iterator iterator;
188 typedef typename std::vector<T>::const_iterator const_iterator;
189 typedef typename std::vector<T>::reverse_iterator reverse_iterator;
190 typedef typename std::vector<T>::const_reverse_iterator const_reverse_iterator;
233 {
polynomial res = *
this; res /= e;
return res; }
245 {
return(this->
real_degree()==0) && (this->size()==0 || (*this)[0]==T(0)); }
246 template <
typename ITER> T horner(power_index &mi,
short_type k,
251 template <
typename ITER> T
eval(
const ITER &it)
const;
255 { n = 0; d = 0; (*this)[0] = 0.0; }
264 : std::vector<T>(
alpha(nn,dd))
265 { n = nn; d = dd; std::fill(this->begin(), this->end(), T(0)); }
269 : std::vector<T>(
alpha(nn,dd)) {
271 std::fill(this->begin(), this->end(), T(0));
277 {
polynomial res = *
this; res *= Q;
return res; }
281 if (dim() != Q.
dim())
return false;
282 const_iterator it1 = this->begin(), ite1 = this->end();
283 const_iterator it2 = Q.begin(), ite2 = Q.end();
284 for ( ; it1 != ite1 && it2 != ite2; ++it1, ++it2)
285 if (*it1 != *it2)
return false;
286 for ( ; it1 != ite1; ++it1)
if (*it1 != T(0))
return false;
287 for ( ; it2 != ite2; ++it2)
if (*it2 != T(0))
return false;
293 {
polynomial res = *
this; res *= e;
return res; }
296 const_reverse_iterator it = this->rbegin(), ite = this->rend();
298 for ( ; it != ite; ++it, --l) {
if (*it != T(0))
break; }
305 this->resize(
alpha(n,dd));
306 if (dd > d) std::fill(this->begin() +
alpha(n,d), this->end(), T(0));
313 GMM_ASSERT2(n == power.size(),
"dimensions mismatch");
314 if (i >= this->size()) { change_degree(power.
degree()); }
315 ((*this)[i]) += coeff;
320 GMM_ASSERT2(Q.
dim() == dim(),
"dimensions mismatch");
323 iterator it = this->begin();
324 const_iterator itq = Q.begin(), ite = Q.end();
325 for ( ; itq != ite; ++itq, ++it) *it += *itq;
331 GMM_ASSERT2(Q.
dim() == dim() && dim() != 0,
"dimensions mismatch");
334 iterator it = this->begin();
335 const_iterator itq = Q.begin(), ite = Q.end();
336 for ( ; itq != ite; ++itq, ++it) *it -= *itq;
343 iterator itq = Q.begin(), ite = Q.end();
344 for ( ; itq != ite; ++itq) *itq = -(*itq);
350 GMM_ASSERT2(Q.
dim() == dim(),
"dimensions mismatch");
353 change_degree(0); (*this)[0] = T(0);
356 if (dim() > 0) miq[dim()-1] = Q.
degree();
357 const_reverse_iterator itq = Q.rbegin(), ite = Q.rend();
358 for ( ; itq != ite; ++itq, --miq) {
360 reverse_iterator ita = aux.rbegin(), itae = aux.rend();
361 std::fill(mia.begin(), mia.end(), 0);
362 if (dim() > 0) mia[dim()-1] = aux.
degree();
363 for ( ; ita != itae; ++ita, --mia)
365 power_index::iterator mita = mia.begin(), mitq = miq.begin();
366 power_index::iterator mit = mitot.begin(), mite = mia.end();
367 for ( ; mita != mite; ++mita, ++mitq, ++mit)
374 add_monomial((*itq) * (*ita), mitot);
386 change_degree(0); n =
short_type(n+Q.
dim()); (*this)[0] = T(0);
390 const_reverse_iterator itq = Q.rbegin(), ite = Q.rend();
391 for ( ; itq != ite; ++itq, --miq)
393 reverse_iterator ita = aux.rbegin(), itae = aux.rend();
394 std::fill(mia.begin(), mia.end(), 0);
396 for ( ; ita != itae; ++ita, --mia)
398 std::copy(mia.begin(), mia.end(), mitot.begin());
399 std::copy(miq.begin(), miq.end(), mitot.begin() + aux.
dim());
400 add_monomial((*itq) * (*ita), mitot);
408 iterator it = this->begin(), ite = this->end();
409 for ( ; it != ite; ++it) (*it) *= e;
415 iterator it = this->begin(), ite = this->end();
416 for ( ; it != ite; ++it) (*it) /= e;
422 GMM_ASSERT2(k < n,
"index out of range");
424 iterator it = this->begin(), ite = this->end();
426 for ( ; it != ite; ++it) {
427 if ((*it) != T(0) && mi[k] > 0)
428 { mi[k]--; (*this)[mi.
global_index()] = (*it) * T(mi[k] + 1); mi[k]++; }
435 template<
typename T>
template<
typename ITER>
441 T v = (*(it+k-1)), res = T(0);
452 template<
typename T>
template<
typename ITER>
455 unsigned deg = degree();
456 const_iterator P = this->begin();
457 if (deg == 0)
return P[0];
460 for (
size_type i=0; i < dim(); ++i) s += it[i]*P[i+1];
467 if (deg == 2)
return P[0] + x*(P[1] + x*(P[2]));
468 if (deg == 3)
return P[0] + x*(P[1] + x*(P[2] + x*(P[3])));
469 if (deg == 4)
return P[0] + x*(P[1] + x*(P[2] + x*(P[3] + x*(P[4]))));
470 if (deg == 5)
return P[0] + x*(P[1] + x*(P[2] + x*(P[3] + x*(P[4] + x*(P[5])))));
471 if (deg == 6)
return P[0] + x*(P[1] + x*(P[2] + x*(P[3] + x*(P[4] + x*(P[5] + x*(P[6]))))));
476 if (deg == 2)
return P[0] + x*(P[1] + x*(P[3])) + y*(P[2] + x*(P[4]) + y*(P[5]));
477 if (deg == 3)
return P[0] + x*(P[1] + x*(P[3] + x*(P[6]))) + y*(P[2] + x*(P[4] + x*(P[7])) + y*(P[5] + x*(P[8]) + y*(P[9])));
478 if (deg == 4)
return P[0] + x*(P[1] + x*(P[3] + x*(P[6] + x*(P[10])))) + y*(P[2] + x*(P[4] + x*(P[7] + x*(P[11]))) + y*(P[5] + x*(P[8] + x*(P[12])) + y*(P[9] + x*(P[13]) + y*(P[14]))));
479 if (deg == 5)
return P[0] + x*(P[1] + x*(P[3] + x*(P[6] + x*(P[10] + x*(P[15]))))) + y*(P[2] + x*(P[4] + x*(P[7] + x*(P[11] + x*(P[16])))) + y*(P[5] + x*(P[8] + x*(P[12] + x*(P[17]))) + y*(P[9] + x*(P[13] + x*(P[18])) + y*(P[14] + x*(P[19]) + y*(P[20])))));
480 if (deg == 6)
return P[0] + x*(P[1] + x*(P[3] + x*(P[6] + x*(P[10] + x*(P[15] + x*(P[21])))))) + y*(P[2] + x*(P[4] + x*(P[7] + x*(P[11] + x*(P[16] + x*(P[22]))))) + y*(P[5] + x*(P[8] + x*(P[12] + x*(P[17] + x*(P[23])))) + y*(P[9] + x*(P[13] + x*(P[18] + x*(P[24]))) + y*(P[14] + x*(P[19] + x*(P[25])) + y*(P[20] + x*(P[26]) + y*(P[27]))))));
486 if (deg == 2)
return P[0] + x*(P[1] + x*(P[4])) + y*(P[2] + x*(P[5]) + y*(P[7])) + z*(P[3] + x*(P[6]) + y*(P[8]) + z*(P[9]));
487 if (deg == 3)
return P[0] + x*(P[1] + x*(P[4] + x*(P[10]))) + y*(P[2] + x*(P[5] + x*(P[11])) + y*(P[7] + x*(P[13]) + y*(P[16]))) + z*(P[3] + x*(P[6] + x*(P[12])) + y*(P[8] + x*(P[14]) + y*(P[17])) + z*(P[9] + x*(P[15]) + y*(P[18]) + z*(P[19])));
488 if (deg == 4)
return P[0] + x*(P[1] + x*(P[4] + x*(P[10] + x*(P[20])))) + y*(P[2] + x*(P[5] + x*(P[11] + x*(P[21]))) + y*(P[7] + x*(P[13] + x*(P[23])) + y*(P[16] + x*(P[26]) + y*(P[30])))) + z*(P[3] + x*(P[6] + x*(P[12] + x*(P[22]))) + y*(P[8] + x*(P[14] + x*(P[24])) + y*(P[17] + x*(P[27]) + y*(P[31]))) + z*(P[9] + x*(P[15] + x*(P[25])) + y*(P[18] + x*(P[28]) + y*(P[32])) + z*(P[19] + x*(P[29]) + y*(P[33]) + z*(P[34]))));
489 if (deg == 5)
return P[0] + x*(P[1] + x*(P[4] + x*(P[10] + x*(P[20] + x*(P[35]))))) + y*(P[2] + x*(P[5] + x*(P[11] + x*(P[21] + x*(P[36])))) + y*(P[7] + x*(P[13] + x*(P[23] + x*(P[38]))) + y*(P[16] + x*(P[26] + x*(P[41])) + y*(P[30] + x*(P[45]) + y*(P[50]))))) + z*(P[3] + x*(P[6] + x*(P[12] + x*(P[22] + x*(P[37])))) + y*(P[8] + x*(P[14] + x*(P[24] + x*(P[39]))) + y*(P[17] + x*(P[27] + x*(P[42])) + y*(P[31] + x*(P[46]) + y*(P[51])))) + z*(P[9] + x*(P[15] + x*(P[25] + x*(P[40]))) + y*(P[18] + x*(P[28] + x*(P[43])) + y*(P[32] + x*(P[47]) + y*(P[52]))) + z*(P[19] + x*(P[29] + x*(P[44])) + y*(P[33] + x*(P[48]) + y*(P[53])) + z*(P[34] + x*(P[49]) + y*(P[54]) + z*(P[55])))));
490 if (deg == 6)
return P[0] + x*(P[1] + x*(P[4] + x*(P[10] + x*(P[20] + x*(P[35] + x*(P[56])))))) + y*(P[2] + x*(P[5] + x*(P[11] + x*(P[21] + x*(P[36] + x*(P[57]))))) + y*(P[7] + x*(P[13] + x*(P[23] + x*(P[38] + x*(P[59])))) + y*(P[16] + x*(P[26] + x*(P[41] + x*(P[62]))) + y*(P[30] + x*(P[45] + x*(P[66])) + y*(P[50] + x*(P[71]) + y*(P[77])))))) + z*(P[3] + x*(P[6] + x*(P[12] + x*(P[22] + x*(P[37] + x*(P[58]))))) + y*(P[8] + x*(P[14] + x*(P[24] + x*(P[39] + x*(P[60])))) + y*(P[17] + x*(P[27] + x*(P[42] + x*(P[63]))) + y*(P[31] + x*(P[46] + x*(P[67])) + y*(P[51] + x*(P[72]) + y*(P[78]))))) + z*(P[9] + x*(P[15] + x*(P[25] + x*(P[40] + x*(P[61])))) + y*(P[18] + x*(P[28] + x*(P[43] + x*(P[64]))) + y*(P[32] + x*(P[47] + x*(P[68])) + y*(P[52] + x*(P[73]) + y*(P[79])))) + z*(P[19] + x*(P[29] + x*(P[44] + x*(P[65]))) + y*(P[33] + x*(P[48] + x*(P[69])) + y*(P[53] + x*(P[74]) + y*(P[80]))) + z*(P[34] + x*(P[49] + x*(P[70])) + y*(P[54] + x*(P[75]) + y*(P[81])) + z*(P[55] + x*(P[76]) + y*(P[82]) + z*(P[83]))))));
534 return horner(mi, dim(), 0, it);
537 template<
typename ITER>
538 typename std::iterator_traits<ITER>::value_type
540 typename std::iterator_traits<ITER>::value_type res
541 =
typename std::iterator_traits<ITER>::value_type(1);
542 power_index::const_iterator mit = mi.begin(), mite = mi.end();
543 for ( ; mit != mite; ++mit, ++it)
554 typename polynomial<T>::const_iterator it = P.begin(), ite = P.end();
556 if (it != ite && *it != T(0))
557 { o << *it; first =
false; ++it; ++n; ++mi; }
558 for ( ; it != ite ; ++it, ++mi ) {
560 bool first_var =
true;
561 if (!first) {
if (*it < T(0)) o <<
" - ";
else o <<
" + "; }
562 else if (*it < T(0)) o <<
"-";
563 if (gmm::abs(*it)!=T(1)) { o << gmm::abs(*it); first_var =
false; }
566 if (!first_var) o <<
"*";
568 if (P.
dim() <= 7) o <<
"xyzwvut"[j];
570 if (mi[j] > 1) o <<
"^" << mi[j];
575 if (n == 0) o <<
"0";
590 GMM_ASSERT2(S.
dim() == 1 && subs_dim < P.
dim(),
591 "wrong arguments for polynomial substitution");
594 std::vector< polynomial<T> > Spow(1);
596 for (
size_type k=0; k < P.size(); ++k, ++pi) {
597 if (P[k] == T(0))
continue;
598 while (pi[subs_dim] >= Spow.size())
599 Spow.push_back(S*Spow.back());
610 template<
typename U,
typename T>
611 polynomial<T> operator *(T c,
const polynomial<T> &p)
612 { polynomial<T> q = p; q *= c;
return q; }
614 typedef polynomial<opt_long_scalar_type> base_poly;
618 inline base_poly null_poly(
short_type n) {
return base_poly(n, 0); }
620 { base_poly res=base_poly(n, 0); res.one();
return res; }
622 {
return base_poly(n, 1, k); }
635 template<
typename T>
class rational_fraction :
public std::vector<T> {
638 polynomial<T> numerator_, denominator_;
642 const polynomial<T> &numerator()
const {
return numerator_; }
643 const polynomial<T> &denominator()
const {
return denominator_; }
645 short_type dim()
const {
return numerator_.dim(); }
648 rational_fraction &operator +=(
const rational_fraction &Q) {
649 numerator_ = numerator_*Q.denominator() + Q.numerator()*denominator_;
650 denominator_ *= Q.denominator();
654 rational_fraction &operator -=(
const rational_fraction &Q) {
655 numerator_ = numerator_*Q.denominator() - Q.numerator()*denominator_;
656 denominator_ *= Q.denominator();
660 rational_fraction
operator +(
const rational_fraction &Q)
const
661 { rational_fraction R = *
this; R += Q;
return R; }
663 rational_fraction
operator -(
const rational_fraction &Q)
const
664 { rational_fraction R = *
this; R -= Q;
return R; }
666 rational_fraction
operator +(
const polynomial<T> &Q)
const
667 { rational_fraction R(numerator_+Q*denominator_, denominator_);
return R; }
669 rational_fraction
operator -(
const polynomial<T> &Q)
const
670 { rational_fraction R(numerator_-Q*denominator_, denominator_);
return R; }
672 { rational_fraction R(-numerator_, denominator_);
return R; }
674 rational_fraction &operator *=(
const rational_fraction &Q)
675 { numerator_*=Q.numerator(); denominator_*=Q.denominator();
return *
this; }
677 rational_fraction &operator /=(
const rational_fraction &Q)
678 { numerator_*=Q.denominator(); denominator_*=Q.numerator();
return *
this; }
680 rational_fraction operator *(
const rational_fraction &Q)
const
681 { rational_fraction R = *
this; R *= Q;
return R; }
683 rational_fraction operator /(
const rational_fraction &Q)
const
684 { rational_fraction R = *
this; R /= Q;
return R; }
686 rational_fraction &operator *=(
const T &e)
687 { numerator_ *= e;
return *
this; }
689 rational_fraction operator *(
const T &e)
const
690 { rational_fraction R = *
this; R *= e;
return R; }
692 rational_fraction &operator /=(
const T &e)
693 { denominator_ *= e;
return *
this; }
695 rational_fraction operator /(
const T &e)
const
696 { rational_fraction res = *
this; res /= e;
return res; }
699 {
return (numerator_==Q.numerator() && denominator_==Q.denominator()); }
701 bool operator !=(
const rational_fraction &Q)
const
705 polynomial<T> der_num = numerator_; der_num.derivative(k);
706 polynomial<T> der_den = denominator_; der_den.derivative(k);
707 if (der_den.is_zero()) {
708 if (der_num.is_zero()) this->
clear();
709 else numerator_ = der_num;
711 numerator_ = der_num * denominator_ - der_den * numerator_;
712 denominator_ = denominator_ * denominator_;
716 void one() { numerator_.one(); denominator_.one(); }
717 void clear() { numerator_.clear(); denominator_.one(); }
718 template <
typename ITER> T eval(
const ITER &it)
const {
719 typedef typename gmm::number_traits<T>::magnitude_type R;
720 T a = numerator_.eval(it), b = denominator_.eval(it);
722 std::vector<T> p(it, it+dim());
725 else gmm::scale(p, R(0.9999999));
726 a = numerator_.eval(p.begin());
727 b = denominator_.eval(p.begin());
729 if (a != T(0)) a /= b;
734 : numerator_(1,0), denominator_(1,0) {
clear(); }
737 : numerator_(dim_,0), denominator_(dim_,0) {
clear(); }
739 rational_fraction(
const polynomial<T> &numer)
740 : numerator_(numer), denominator_(numer.dim(),0) { denominator_.one(); }
742 rational_fraction(
const polynomial<T> &numer,
const polynomial<T> &denom)
743 : numerator_(numer), denominator_(denom)
744 { GMM_ASSERT1(numer.dim() == denom.dim(),
"Dimensions mismatch"); }
750 const rational_fraction<T> &Q) {
751 rational_fraction<T> R(P*Q.denominator()+Q.numerator(), Q.denominator());
757 const rational_fraction<T> &Q) {
758 rational_fraction<T> R(P*Q.denominator()-Q.numerator(), Q.denominator());
762 template<
typename T> std::ostream &
operator <<
763 (std::ostream &o,
const rational_fraction<T>& P) {
764 o <<
"[" << P.numerator() <<
"]/[" << P.denominator() <<
"]";
768 typedef rational_fraction<opt_long_scalar_type> base_rational_fraction;
defines and typedefs for namespace bgeot
This class deals with plain polynomials with several variables.
short_type real_degree() const
gives the degree of the polynomial, considering only non-zero coefficients
polynomial(short_type dim_, short_type degree_, short_type k)
Constructor for the polynomial 'x' (k=0), 'y' (k=1), 'z' (k=2) etc.
polynomial & operator+=(const polynomial &Q)
Add Q to P. P contains the result.
short_type dim() const
Gives the dimension (number of variables)
polynomial & operator-=(const polynomial &Q)
Subtract Q from P. P contains the result.
void add_monomial(const T &coeff, const power_index &power)
Add to the polynomial a monomial of coefficient a and correpsonding to the power index pi.
polynomial operator*(const polynomial &Q) const
Multiply P with Q.
polynomial(short_type dim_, short_type degree_)
Constructor.
T eval(const ITER &it) const
Evaluate the polynomial.
void change_degree(short_type dd)
Change the degree of the polynomial to d.
void derivative(short_type k)
Derivative of P with respect to the variable k. P contains the result.
polynomial operator+(const polynomial &Q) const
Add Q to P.
bool operator!=(const polynomial &Q) const
operator !=.
polynomial & operator*=(const polynomial &Q)
Multiply P with Q. P contains the result.
short_type degree() const
Gives the degree of the polynomial.
polynomial & operator/=(const T &e)
Divide P with the scalar a. P contains the result.
polynomial operator/(const T &e) const
Divide P with the scalar a.
void direct_product(const polynomial &Q)
Product of P and Q considering that variables of Q come after variables of P.
bool operator==(const polynomial &Q) const
operator ==.
Vector of integer (16 bits type) which represent the powers of a monomial.
power_index()
Constructor.
const power_index & operator++()
Gives the next power index.
short_type degree() const
Gives the degree.
size_type global_index() const
Gives the global number of the index (i.e.
const power_index & operator--()
Gives the next previous index.
Stores interdependent getfem objects.
number_traits< typename linalg_traits< V >::value_type >::magnitude_type vect_norm2(const V &v)
Euclidean norm of a vector.
void fill_random(L &l)
fill a vector or matrix with random value (uniform [-1,1]).
void clear(L &l)
clear (fill with zeros) a vector or matrix.
gmm::uint16_type short_type
used as the common short type integer in the library
base_poly read_base_poly(short_type n, std::istream &f)
read a base_poly on the stream ist.
rational_fraction< T > operator-(const polynomial< T > &P, const rational_fraction< T > &Q)
Subtract Q from P.
polynomial< T > poly_substitute_var(const polynomial< T > &P, const polynomial< T > &S, size_type subs_dim)
polynomial variable substitution
rational_fraction< T > operator+(const polynomial< T > &P, const rational_fraction< T > &Q)
Add Q to P.
std::ostream & operator<<(std::ostream &o, const convex_structure &cv)
Print the details of the convex structure cvs to the output stream o.
bool operator==(const pconvex_structure &p1, const pconvex_structure &p2)
Stored objects must be compared by keys, because there is a possibility that they are duplicated in s...
size_t size_type
used as the common size type in the library
size_type alpha(short_type n, short_type d)
Return the value of which is the number of monomials of a polynomial of variables and degree .