polynomials_on_simplices.piecewise_polynomial.piecewise_polynomial module

Functionality for discontinuous Galerkin finite elements (piecewise polynomials) on a simplicial domain (triangulation) \(\mathcal{T}\).

The space of piecewise polynomials of degree r on \(\mathcal{T}, D\mathcal{P}_r (\mathcal{T}) = D\mathcal{P}_r (\mathcal{T}, \mathbb{R}^n) \subset L^2(\mathcal{T}, \mathbb{R}^n)\) is defined as

\[D\mathcal{P}_r (\mathcal{T}) = \{ v \in L^2(\mathcal{T}, \mathbb{R}^n) \big| v|_T \in \mathcal{P}_r(T, \mathbb{R}^n) \, \forall T \in \mathcal{T} \}.\]

Correspondingly the space of piecewise polynomials of degree r on \(\mathcal{T}\) which are zero on specified set of simplices or subsimplices of \(\mathcal{T}, D\mathcal{P}_{r, 0} (\mathcal{T}\) is defined as

\[D\mathcal{P}_{r, 0} (\mathcal{T}) = \{ v \in L^2_0(\mathcal{T}, \mathbb{R}^n) \big| v|_T \in \mathcal{P}_r(T, \mathbb{R}^n) \, \forall T \in \mathcal{T} \}.\]
class PiecewisePolynomialBase(coeff, triangles, vertices, r, tau=None, boundary_simplices=None, keep_boundary_dofs_last=False, support=None, bsp_tree=None)[source]

Bases: abc.ABC

Abstract base class for a piecewise polynomial function of degree r on a triangulation \(\mathcal{T}\), i.e. an element of \(D\mathcal{P}_r (\mathcal{T})\) or \(D\mathcal{P}_{r, 0} (\mathcal{T})\).

The domain dimension m = \(\dim \mathcal{T}\) and the target dimension n of the piecewise polynomial is given by the domain_dimension() and target_dimension() functions respectively.

The degree r of the piecewise polynomial is given by the degree() method.

Let \(N = \dim D\mathcal{P}_r (\mathcal{T}, \mathbb{R}^n)\). This value is given by the dimension() method.

Restriction to a simplex

On each simplex \(T \in \mathcal{T}\) the piecewise polynomial is given by a polynomial of degree r. This polynomial can be acquired using the restrict_to_simplex() method.

Basis

This class assumes that a basis \(\{ \bar{\varphi}_{\nu} \}_{\substack{\nu \in \mathbb{N}_0^m \\ | \nu | \leq r}}\) for \(\mathcal{P}_r (\Delta_c^m)\) has been chosen. From this bases \(\{ \varphi_{\nu, j} \}_{\substack{\nu \in \mathbb{N}_0^m \\ | \nu | \leq r}}\) for \(\mathcal{P}_r (T_j, \mathbb{R}^n)\) are constructed by

\[\varphi_{\nu, j}(x) = (\bar{\varphi}_{\nu} \circ \Phi_{T_j}^{-1})(x),\]

where \(\Phi_{T_j}\) is the unique affine map which maps the unit simplex \(\Delta_c^m\) to the simplex \(T_j\) (the i:th vertex of the unit simplex is mapped to the i:th vertex of \(T_j\)).

This class also assumes that a local-to-global map \(\tau, \tau : \{ 0, 1, \ldots, | \mathcal{T} | - 1 \} \times \mathbb{N}_0^m \to \mathbb{N}_0 \cup \{ -1 \}\) is available which maps local (simplex) basis functions to the index of the corresponding global basis function (or to -1 if the local basis function doesn’t correspond to a global basis function (which is the case if the piecewise polynomial has a fixed (prescribed) value on a (sub)simplex)). See generate_local_to_global_map(). Then a basis \(\{ \phi_i \}_{i = 1}^N\) is given by

\[\phi_i(x) = \sum_{(j, \nu) \in \operatorname{preim}_{\tau}(i)} \chi_{T_j}(x) \cdot \varphi_{\nu, j}(x),\]

And an arbitrary piecewise polynomial \(p(x)\) is given by specifying coefficients \(a_i \in \mathbb{R}^n, i = 1, 2, \ldots, N\) in front of each basis function.

\[p(x) = \sum_{i = 1}^N a_i \phi_i(x).\]

Hence the value of the piecewise polynomial function on a simplex \(T_j\) is given by

\[\begin{split}p(x) = \sum_{\substack{\nu \in \mathbb{N}_0^m \\ | \nu | \leq r \\ \tau(j, \nu) \neq -1}} a_{\tau(j, \nu)} \cdot \varphi_{\nu, j}(x), x \in T_j.\end{split}\]

The basis chosen for \(\mathcal{P}_r (\Delta_c^m)\) is given by the basis() method.

Algebraic structure

This class also defines the basic algebraic structures of the space of piecewise polynomials.

Ring structure:

Addition: \(+ : D\mathcal{P}_r (\mathcal{T}) \times D\mathcal{P}_r (\mathcal{T}) \to D\mathcal{P}_r (\mathcal{T}), (p_1 + p_2)(x) = p_1(x) + p_2(x)\).

Multiplication: \(\cdot : D\mathcal{P} (\mathcal{T}) \times D\mathcal{P} (\mathcal{T}) \to D\mathcal{P} (\mathcal{T}), (p_1 \cdot p_2)(x) = p_1(x) \cdot p_2(x)\) (only for n = 1).

Vector space structure:

Scalar multiplication: \(\cdot : \mathbb{R} \times D\mathcal{P}_r (\mathcal{T}) \to D\mathcal{P}_r (\mathcal{T}), (s \cdot p)(x) = s \cdot p(x)\).

Parameters:
  • coeff (List[Union[Scalar, Vector]]) – Coefficients for the piecewise polynomial in the \(\{ \phi_i \}_{i = 1}^N\) basis derived from the basis for \(\mathcal{P}_r (\Delta_c^m)\) used (see basis()) and the local-to-global map \(\tau\).
  • triangles – Triangles (or in general simplices) in the mesh \(\mathcal{T}\) (num_triangles by m + 1 array of indices).
  • vertices – Vertices in the mesh \(\mathcal{T}\) (num_vertices by m array of scalars).
  • r (int) – Degree of each polynomial in the piecewise polynomial.
  • tau (Callable \(\tau(j, \nu)\)) – Local-to-global map for mapping local basis functions to the index of the corresponding global basis function. Will be generated if not supplied.
  • boundary_simplices (List[List[int]]) – List of simplices or subsimplices on which the piecewise polynomial function should vanish (for an element of \(D\mathcal{P}_{r, 0} (\mathcal{T})\)) or which should be treated separately (if keep_boundary_dofs_last is set to True). Each simplex or subsimplex is specified as a list of vertex indices of the vertices that form the simplex.
  • keep_boundary_dofs_last (bool) – Whether or not to collect all global basis functions associated with any boundary simplex last in the enumeration of all basis functions. Enumerating basis functions associated with boundary simplices last is useful for handling \(D\mathcal{P}_{r, 0} (\mathcal{T})\) as a subset of \(D\mathcal{P}_r (\mathcal{T})\) in a practical way.
  • support (Optional[Set[int]]) – Indices of the triangles in the triangulation where the piecewise polynomial is supported. Will be generated if not supplied.
  • bsp_tree – Optional implementation detail. A binary space partitioning tree built around the triangulation for quicker lookup of triangle a point lies in. Will be generated if not supplied.
basis()[source]

Get basis for the space \(\mathcal{P}_r (\Delta_c^m)\) used to express the piecewise polynomial.

Returns:Unique identifier for the basis used.
Return type:str
degree()[source]

Get degree of each polynomial in the piecewise polynomial.

Returns:Polynomial degree.
Return type:int
degree_elevate(s)[source]

Express the piecewise polynomial using a higher degree polynomial basis.

Let \(p(x)\) be this piecewise polynomial. Let \(\{ \bar{\varphi}_{\nu, r} \}_{\substack{\nu \in \mathbb{N}_0^m \\ |\nu| \leq r}}\) be the polynomial basis for \(\mathcal{P}_r(\Delta_c^m)\) used for this piecewise polynomial, and let \(\{ \bar{\varphi}_{\nu, s} \}_{\substack{\nu \in \mathbb{N}_0^m \\ |\nu| \leq s}}, s \geq r\) be the corresponding higher degree polynomial basis for \(\mathcal{P}_s(\Delta_c^m)\). Then this function returns a piecewise polynomial \(q\) using the basis \(\{ \bar{\varphi}_{\nu, s} \}\) such that \(p(x) = q(x) \, \forall x \in \mathbb{R}^m\).

Parameters:s (int) – New degree for the basic polynomial basis the piecewise polynomial should use.
Returns:Elevation of this piecewise polynomial to the higher degree basis.
dimension()[source]

Get dimension of the space of piecewise polynomial the function belongs to.

Returns:Dimension of the function space.
Return type:int
domain_dimension()[source]

Get dimension of the domain of the piecewise polynomial.

Returns:Dimension of the domain of the piecewise polynomial.
Return type:int
evaluate_on_simplex(j, x)[source]

Evaluate the piecewise polynomial at a point \(x \in \mathbb{R}^m\) in the j:th simplex of the triangulation the piecewise polynomial is defined on.

In other words the polynomial at the j:th simplex is evaluated at the point x.

Parameters:
  • j (int) – Index of simplex where the piecewise polynomial should be evaluated.
  • x (float or length m Numpy array) – Point where the piecewise polynomial should be evaluated.
Returns:

Value of the piecewise polynomial.

Return type:

float or length n Numpy array.

static generate_local_to_global_map(triangles, r, boundary_simplices=None, keep_boundary_dofs_last=False)[source]

Generate a local-to-global map \(\tau\) for the space of piecewise polynomial functions of degree r on a triangulation \(\mathcal{T}, D \mathcal{P}_r (\mathcal{T})\) or \(D \mathcal{P}_{r, 0} (\mathcal{T})\). See generate_local_to_global_map().

multiply_with_constant(c)[source]

Multiplication of this piecewise polynomial with a constant scalar or a vector (only if n = 1), self * c.

Parameters:c (Union[float, Numpy array]) – Scalar or vector we should multiply this polynomial with.
Returns:Product of this piecewise polynomial with the constant.
Return type:Instance of self.__class__
restrict_to_simplex(i)[source]

Restriction of the piecewise polynomial to a specified simplex \(T_i \in \mathcal{T}\).

Parameters:i (int) – Index of the simplex we want to restrict the piecewise polynomial to (in \(0, 1, \ldots, | \mathcal{T} | - 1\)).
Returns:Polynomial which agrees with the piecewise polynomial on the simplex \(T_i\).
support()[source]

Get the indices of the simplices in \(\mathcal{T}\) where this piecewise polynomial is non-zero.

Returns:Indices of simplices whose intersection with the support of this piecewise polynomial is non-empty.
Return type:Set[int]
target_dimension()[source]

Get dimension of the target of the piecewise polynomial.

Returns:Dimension of the target of the piecewise polynomial.
Return type:int
generate_inverse_local_to_global_map(tau, num_triangles, num_dofs, r, m)[source]

Generate the inverse of a local-to-global map (i.e. a global-to-local map).

Let \(\tau : \{ 0, 1, \ldots, | \mathcal{T} | - 1 \} \times \mathbb{N}_0^m \to \mathbb{N}_0 \cup \{ -1 \}\) be a local-to-global map for the space of piecewise polynomials of degree r on a triangulation (see generate_local_to_global_map()). Then this function generates the inverse of \(\tau\), i.e. the map \(\tau^{-1} : \mathbb{N}_0 \to \{ 0, 1, \ldots, | \mathcal{T} | - 1 \} \times \mathbb{N}_0^m\) (\(\tau\) is not invertible where it has the value -1) such that

\[\tau(\tau^{-1} (i)) = i, \, i = 0, 1, \ldots, N - 1\]

where \(N\) is the number of degrees of freedom for the piecewise polynomial.

Parameters:
  • tau (Callable \(\tau(j, \nu)\)) – Local-to-global map that we want to invert.
  • num_triangles (int) – Number of triangles in the triangulation \(\mathcal{T}\) on which the piecewise polynomials associated with the local-to-global map are defined.
  • num_dofs (int) – Number of degrees of freedom (dimension) for the space of piecewise polynomials.
  • r (int) – Polynomial degree for the piecewise polynomials.
  • m (int) – Dimension of the domain of the piecewise polynomials.
Returns:

Inverse local-to-global map.

Return type:

Callable \(\tau^{-1}(i)\).

generate_inverse_vector_valued_local_to_global_map(tau, num_triangles, num_dofs, r, m, n)[source]

Generate the inverse of a local-to-global map (i.e. a global-to-local map).

Let \(\tau : \{ 0, 1, \ldots, | \mathcal{T} | - 1 \} \times \mathbb{N}_0^m \times \{ 0, 1, \ldots, n - 1 \} \to \mathbb{N}_0 \cup \{ -1 \}\) be a local-to-global map for the space of vector valued piecewise polynomials of degree r on a triangulation (see generate_vector_valued_local_to_global_map()). Then this function generates the inverse of \(\tau\), i.e. the map \(\tau^{-1} : \mathbb{N}_0 \to \{ 0, 1, \ldots, | \mathcal{T} | - 1 \} \times \mathbb{N}_0^m \times \{ 0, 1, \ldots, n - 1 \}\) (\(\tau\) is not invertible where it has the value -1) such that

\[\tau(\tau^{-1} (i)) = i, \, i = 0, 1, \ldots, N - 1\]

where \(N\) is the number of degrees of freedom for the vector valued piecewise polynomial.

Parameters:
  • tau (Callable \(\tau(j, \nu, k)\)) – Local-to-global map that we want to invert.
  • num_triangles (int) – Number of triangles in the triangulation \(\mathcal{T}\) on which the piecewise polynomials associated with the local-to-global map are defined.
  • num_dofs (int) – Number of degrees of freedom (dimension) for the space of piecewise polynomials.
  • r (int) – Polynomial degree for the piecewise polynomials.
  • m (int) – Dimension of the domain of the piecewise polynomials.
  • n (int) – Dimension of the target of the piecewise polynomials.
Returns:

Inverse local-to-global map.

Return type:

Callable \(\tau^{-1}(i)\).

generate_local_to_global_map(triangles, r, boundary_simplices=None, keep_boundary_dofs_last=False)[source]

Generate a local-to-global map \(\tau\) for the space of piecewise polynomial functions of degree r on a triangulation \(\mathcal{T}, D \mathcal{P}_r (\mathcal{T})\) or \(D \mathcal{P}_{r, 0} (\mathcal{T})\).

We have \(\tau : \{ 0, 1, \ldots, | \mathcal{T} | - 1 \} \times \mathbb{N}_0^m \to \mathbb{N}_0 \cup \{ -1 \}\).

A local-to-global map maps a local polynomial basis function on a simplex in the triangulation to a global basis index for a basis for the space of piecewise polynomials. The value -1 is used to indicate that a local basis function has no corresponding global basis function (see later).

Let \(d = \dim \mathcal{P}_r (\Delta_c^m)\), and let the triangles in \(\mathcal{T}\) be enumerated so that we have \(\mathcal{T} = \{ T_0, T_1, \ldots, T_{| \mathcal{T} | - 1} \}\). Given a basis \(\{ \bar{\varphi}_{\nu} \}_{\substack{\nu \in \mathbb{N}_0^m \\ | \nu | \leq r}}\) for \(\mathcal{P}_r (\Delta_c^m)\), the functions \(\{ \bar{\varphi}_{\nu} \circ \Phi_{T_j}^{-1} \}_{ \substack{\nu \in \mathbb{N}_0^m \\ | \nu | \leq r}}\) is a basis for \(\mathcal{P}_r(T_j)\), where \(\Phi_{T_j}, j = 0, 1, \ldots, | \mathcal{T} | - 1\) is the unique affine map which maps the unit simplex \(\Delta_c^m\) to the simplex \(T_j\). From these local bases and the local-to-global map \(\tau\) a basis \(\{ \phi_i \}_{i = 1}^N\) for \(D \mathcal{P}_r (\mathcal{T})\) is constructed, with

\[\begin{split}\phi_{\tau(j, \nu)}(x) = \begin{cases} \bar{\varphi}_{\nu} \circ \Phi_{T_j}^{-1}(x), & x \in T_j \\ 0, & \text{else} \end{cases}.\end{split}\]

Optionally a set of boundary simplices can be prescribed where the piecewise polynomial function should vanish. This is achieved by associating an invalid global index (-1) for all local basis functions supported on any of the boundary simplices). Alternatively local basis functions supported on boundary simplices can be associated with global indices placed last in the enumeration of global basis functions (by setting the keep_boundary_dofs_last to true).

Parameters:
  • triangles – Triangles (or in general simplices) in the mesh (num_triangles by n + 1 list of indices, where n is the dimension of each simplex).
  • r (int) – Polynomial degree for the piecewise polynomial functions.
  • boundary_simplices (List[List[int]]) – List of simplices or subsimplices on which the piecewise polynomial functions should vanish (for an element of \(D\mathcal{P}_{r, 0} (\mathcal{T})\)) or which should be treated separately (if keep_boundary_dofs_last is set to True). Each simplex or subsimplex is specified as a list of vertex indices of the vertices that form the simplex.
  • keep_boundary_dofs_last (bool) – Whether or not to collect all global basis functions associated with any boundary simplex last in the enumeration of all basis functions. Enumerating basis functions associated with boundary simplices last is useful for handling \(D\mathcal{P}_{r, 0} (\mathcal{T})\) as a subset of \(D\mathcal{P}_r (\mathcal{T})\) in a practical way.
Returns:

Tuple containing the local to global map \(\tau\) and the number of global basis functions. If keep_boundary_dofs_last is true then the number of global basis functions not supported on the boundary is also returned.

Return type:

Tuple[Callable \(\tau(j, \nu)\), int, Optional[int]].

generate_local_to_global_map_as_dictionaries(triangles, r, boundary_simplices=None, keep_boundary_dofs_last=False)[source]

Generate a local-to-global map for the space of piecewise polynomial functions of degree r on a triangulation \(\mathcal{T}, D \mathcal{P}_r (\mathcal{T})\) or \(D \mathcal{P}_{r, 0} (\mathcal{T})\). For details see generate_local_to_global_map(). This function differs from that function in that a list of dictionaries (one for each triangle in the mesh) with multi-indices (second function argument) as keys and global DOF indices (function value) as values is returned instead of a callable, which might be preferable in some cases.

Parameters:
  • triangles – Triangles (or in general simplices) in the mesh (num_triangles by n + 1 list of indices, where n is the dimension of each simplex).
  • r (int) – Polynomial degree for the piecewise polynomial functions.
  • boundary_simplices (List[List[int]]) – List of simplices or subsimplices on which the piecewise polynomial functions should vanish (for an element of \(D\mathcal{P}_{r, 0} (\mathcal{T})\)) or which should be treated separately (if keep_boundary_dofs_last is set to True). Each simplex or subsimplex is specified as a list of vertex indices of the vertices that form the simplex.
  • keep_boundary_dofs_last (bool) – Whether or not to collect all global basis functions associated with any boundary simplex last in the enumeration of all basis functions. Enumerating basis functions associated with boundary simplices last is useful for handling \(D\mathcal{P}_{r, 0} (\mathcal{T})\) as a subset of \(D\mathcal{P}_r (\mathcal{T})\) in a practical way.
Returns:

Tuple containing the local to global map dictionaries, one for each triangle, and the number of global basis functions. If keep_boundary_dofs_last is true then the number of global basis functions not supported on the boundary is also returned.

Return type:

Tuple[List[Dict[Tuple[int], int]], int, Optional[int]]

generate_vector_valued_local_to_global_map(triangles, r, n, boundary_simplices=None, keep_boundary_dofs_last=False, ordering='interleaved')[source]

Generate a local-to-global map \(\tau\) for the space of vector valued piecewise polynomial functions of degree r on a triangulation \(\mathcal{T}\) with values in \(\mathbb{R}^n\), D mathcal{P}_r (mathcal{T}, mathbb{R}^n)` or \(D \mathcal{P}_{r, 0} (\mathcal{T}, \mathbb{R}^n)\).

We have \(\tau : \{ 0, 1, \ldots, | \mathcal{T} | - 1 \} \times \mathbb{N}_0^m \times \{ 0, 1, \ldots, n - 1 \} \to \mathbb{N}_0 \cup \{ -1 \}\).

A local-to-global map maps a local polynomial basis function on a simplex in the triangulation to a global basis index for a basis for the space of vector valued piecewise polynomials. The value -1 is used to indicate that a local basis function has no corresponding global basis function (see later).

Let \(d = \dim \mathcal{P}_r (\Delta_c^m, \mathbb{R}^n)\), and let the triangles in \(\mathcal{T}\) be enumerated so that we have \(\mathcal{T} = \{ T_0, T_1, \ldots, T_{| \mathcal{T} | - 1} \}\). Given a basis \(\{ \bar{\varphi}_{\nu, i} \}_{\substack{\nu \in \mathbb{N}_0^m \\ | \nu | \leq r \\ i \in \{ 0, 1, \ldots, n - 1 \}}}\) for \(\mathcal{P}_r (\Delta_c^m, \mathbb{R}^n)\), the functions \(\{ \bar{\varphi}_{\nu, i} \circ \Phi_{T_j}^{-1} \}_{\substack{\nu \in \mathbb{N}_0^m \\ | \nu | \leq r \\ i \in \{ 0, 1, \ldots, n - 1 \}}}\) is a basis for \(\mathcal{P}_r(T_j, \mathbb{R}^n)\), where \(\Phi_{T_j}, j = 0, 1, \ldots, | \mathcal{T} | - 1\) is the unique affine map which maps the unit simplex \(\Delta_c^m\) to the simplex \(T_j\). From these local bases and the local-to-global map \(\tau\) a basis \(\{ \phi_i \}_{i = 1}^N\) for \(D \mathcal{P}_r (\mathcal{T}, \mathbb{R}^n)\) is constructed, with

\[\begin{split}\phi_{\tau(j, \nu, k)}(x) = \begin{cases} \bar{\varphi}_{\nu, k} \circ \Phi_{T_j}^{-1}(x), & x \in T_j \\ 0, & \text{else} \end{cases}.\end{split}\]

Optionally a set of boundary simplices can be prescribed where the piecewise polynomial function should vanish. This is achieved by associating an invalid global index (-1) for all local basis functions supported on any of the boundary simplices). Alternatively local basis functions supported on boundary simplices can be associated with global indices placed last in the enumeration of global basis functions (by setting the keep_boundary_dofs_last to true).

Parameters:
  • triangles – Triangles (or in general simplices) in the mesh (num_triangles by n + 1 list of indices, where n is the dimension of each simplex).
  • r (int) – Polynomial degree for the piecewise polynomial functions.
  • n (int) – Dimension of the target.
  • boundary_simplices (List[List[int]]) – List of simplices or subsimplices on which the piecewise polynomial functions should vanish (for an element of \(D\mathcal{P}_{r, 0} (\mathcal{T}, \mathbb{R}^n)\)) or which should be treated separately (if keep_boundary_dofs_last is set to True). Each simplex or subsimplex is specified as a list of vertex indices of the vertices that form the simplex.
  • keep_boundary_dofs_last (bool) – Whether or not to collect all global basis functions associated with any boundary simplex last in the enumeration of all basis functions. Enumerating basis functions associated with boundary simplices last is useful for handling \(D\mathcal{P}_{r, 0} (\mathcal{T}, \mathbb{R}^n)\) as a subset of \(D\mathcal{P}_r (\mathcal{T}, \mathbb{R}^n)\) in a practical way.
  • ordering (str) – How the vector valued basis functions are ordered. Can be “sequential” or “interleaved”. For sequential, sorting is first done on the index of the component that is non-zero, and then the non-zero component is sorted in the same way as the scalar valued basis functions. For “interleaved” basis functions are first sorted on their non-zero component in the same way as scalar valued basis functions, and then they are sorted on the index of the component that is non-zero.
Returns:

Tuple containing the local to global map \(\tau\) and the number of global basis functions. If keep_boundary_dofs_last is true then the number of global basis functions not supported on the boundary is also returned.

Return type:

Tuple[Callable \(\tau(j, \nu, k)\), int, Optional[int]].

piecewise_polynomials_equal(p1, p2, rel_tol=1e-09, abs_tol=1e-07)[source]

Check if two piecewise polynomials p1 and p2 are approximately equal by comparing their values on each simplex in the domain of the first piecewise polynomial.

For equality on a simplex the polynomials_equal_on_simplex() function is used.

Parameters:
  • p1 (Implementation of the PiecewisePolynomialBase interface) – First piecewise polynomial.
  • p2 (Callable p2(x)) – Second piecewise polynomial.
  • rel_tol (float) – Tolerance for the relative error. See math.isclose for details.
  • abs_tol (float) – Tolerance for the absolute error. See math.isclose for details.
Returns:

Whether or not the two piecewise polynomials are approximately equal.

Return type:

bool