WPILibC++ 2023.4.3
Eigen::CompleteOrthogonalDecomposition< _MatrixType > Class Template Reference

Complete orthogonal decomposition (COD) of a matrix. More...

#include </home/runner/work/allwpilib/allwpilib/wpimath/src/main/native/thirdparty/eigen/include/Eigen/src/QR/CompleteOrthogonalDecomposition.h>

Inheritance diagram for Eigen::CompleteOrthogonalDecomposition< _MatrixType >:
Eigen::SolverBase< CompleteOrthogonalDecomposition< _MatrixType > > Eigen::EigenBase< CompleteOrthogonalDecomposition< _MatrixType > >

Public Types

enum  { MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime , MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime }
 
typedef _MatrixType MatrixType
 
typedef SolverBase< CompleteOrthogonalDecompositionBase
 
typedef internal::plain_diag_type< MatrixType >::type HCoeffsType
 
typedef PermutationMatrix< ColsAtCompileTime, MaxColsAtCompileTimePermutationType
 
typedef internal::plain_row_type< MatrixType, Index >::type IntRowVectorType
 
typedef internal::plain_row_type< MatrixType >::type RowVectorType
 
typedef internal::plain_row_type< MatrixType, RealScalar >::type RealRowVectorType
 
typedef HouseholderSequence< MatrixType, typename internal::remove_all< typename HCoeffsType::ConjugateReturnType >::typeHouseholderSequenceType
 
typedef MatrixType::PlainObject PlainObject
 
- Public Types inherited from Eigen::SolverBase< CompleteOrthogonalDecomposition< _MatrixType > >
enum  
 
typedef EigenBase< CompleteOrthogonalDecomposition< _MatrixType > > Base
 
typedef internal::traits< CompleteOrthogonalDecomposition< _MatrixType > >::Scalar Scalar
 
typedef Scalar CoeffReturnType
 
typedef internal::add_const< Transpose< constDerived > >::type ConstTransposeReturnType
 
typedef internal::conditional< NumTraits< Scalar >::IsComplex, CwiseUnaryOp< internal::scalar_conjugate_op< Scalar >, ConstTransposeReturnType >, ConstTransposeReturnType >::type AdjointReturnType
 
- Public Types inherited from Eigen::EigenBase< CompleteOrthogonalDecomposition< _MatrixType > >
typedef Eigen::Index Index
 The interface type of indices. More...
 
typedef internal::traits< CompleteOrthogonalDecomposition< _MatrixType > >::StorageKind StorageKind
 

Public Member Functions

 CompleteOrthogonalDecomposition ()
 Default Constructor. More...
 
 CompleteOrthogonalDecomposition (Index rows, Index cols)
 Default Constructor with memory preallocation. More...
 
template<typename InputType >
 CompleteOrthogonalDecomposition (const EigenBase< InputType > &matrix)
 Constructs a complete orthogonal decomposition from a given matrix. More...
 
template<typename InputType >
 CompleteOrthogonalDecomposition (EigenBase< InputType > &matrix)
 Constructs a complete orthogonal decomposition from a given matrix. More...
 
HouseholderSequenceType householderQ (void) const
 
HouseholderSequenceType matrixQ (void) const
 
MatrixType matrixZ () const
 
const MatrixTypematrixQTZ () const
 
const MatrixTypematrixT () const
 
template<typename InputType >
CompleteOrthogonalDecompositioncompute (const EigenBase< InputType > &matrix)
 
const PermutationTypecolsPermutation () const
 
MatrixType::RealScalar absDeterminant () const
 
MatrixType::RealScalar logAbsDeterminant () const
 
Index rank () const
 
Index dimensionOfKernel () const
 
bool isInjective () const
 
bool isSurjective () const
 
bool isInvertible () const
 
const Inverse< CompleteOrthogonalDecompositionpseudoInverse () const
 
Index rows () const
 
Index cols () const
 
const HCoeffsTypehCoeffs () const
 
const HCoeffsTypezCoeffs () const
 
CompleteOrthogonalDecompositionsetThreshold (const RealScalar &threshold)
 Allows to prescribe a threshold to be used by certain methods, such as rank(), who need to determine when pivots are to be considered nonzero. More...
 
CompleteOrthogonalDecompositionsetThreshold (Default_t)
 Allows to come back to the default behavior, letting Eigen use its default formula for determining the threshold. More...
 
RealScalar threshold () const
 Returns the threshold that will be used by certain methods such as rank(). More...
 
Index nonzeroPivots () const
 
RealScalar maxPivot () const
 
ComputationInfo info () const
 Reports whether the complete orthogonal decomposition was successful. More...
 
template<typename RhsType , typename DstType >
void _solve_impl (const RhsType &rhs, DstType &dst) const
 
template<bool Conjugate, typename RhsType , typename DstType >
void _solve_impl_transposed (const RhsType &rhs, DstType &dst) const
 
- Public Member Functions inherited from Eigen::SolverBase< CompleteOrthogonalDecomposition< _MatrixType > >
 SolverBase ()
 Default constructor. More...
 
 ~SolverBase ()
 
const Solve< CompleteOrthogonalDecomposition< _MatrixType >, Rhs > solve (const MatrixBase< Rhs > &b) const
 
ConstTransposeReturnType transpose () const
 
AdjointReturnType adjoint () const
 
EIGEN_DEVICE_FUNC CompleteOrthogonalDecomposition< _MatrixType > & derived ()
 
EIGEN_DEVICE_FUNC const CompleteOrthogonalDecomposition< _MatrixType > & derived () const
 
- Public Member Functions inherited from Eigen::EigenBase< CompleteOrthogonalDecomposition< _MatrixType > >
EIGEN_DEVICE_FUNC CompleteOrthogonalDecomposition< _MatrixType > & derived ()
 
EIGEN_DEVICE_FUNC const CompleteOrthogonalDecomposition< _MatrixType > & derived () const
 
EIGEN_DEVICE_FUNC CompleteOrthogonalDecomposition< _MatrixType > & const_cast_derived () const
 
EIGEN_DEVICE_FUNC const CompleteOrthogonalDecomposition< _MatrixType > & const_derived () const
 
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR Index rows () const EIGEN_NOEXCEPT
 
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR Index cols () const EIGEN_NOEXCEPT
 
EIGEN_DEVICE_FUNC EIGEN_CONSTEXPR Index size () const EIGEN_NOEXCEPT
 
EIGEN_DEVICE_FUNC void evalTo (Dest &dst) const
 
EIGEN_DEVICE_FUNC void addTo (Dest &dst) const
 
EIGEN_DEVICE_FUNC void subTo (Dest &dst) const
 
EIGEN_DEVICE_FUNC void applyThisOnTheRight (Dest &dst) const
 
EIGEN_DEVICE_FUNC void applyThisOnTheLeft (Dest &dst) const
 

Protected Member Functions

template<bool Transpose_, typename Rhs >
void _check_solve_assertion (const Rhs &b) const
 
void computeInPlace ()
 Performs the complete orthogonal decomposition of the given matrix matrix. More...
 
template<bool Conjugate, typename Rhs >
void applyZOnTheLeftInPlace (Rhs &rhs) const
 Overwrites rhs with \( \mathbf{Z} * \mathbf{rhs} \) or \( \mathbf{\overline Z} * \mathbf{rhs} \) if Conjugate is set to true. More...
 
template<typename Rhs >
void applyZAdjointOnTheLeftInPlace (Rhs &rhs) const
 Overwrites rhs with \( \mathbf{Z}^* * \mathbf{rhs} \). More...
 
- Protected Member Functions inherited from Eigen::SolverBase< CompleteOrthogonalDecomposition< _MatrixType > >
void _check_solve_assertion (const Rhs &b) const
 

Static Protected Member Functions

static void check_template_parameters ()
 

Protected Attributes

ColPivHouseholderQR< MatrixTypem_cpqr
 
HCoeffsType m_zCoeffs
 
RowVectorType m_temp
 

Friends

template<typename Derived >
struct internal::solve_assertion
 

Detailed Description

template<typename _MatrixType>
class Eigen::CompleteOrthogonalDecomposition< _MatrixType >

Complete orthogonal decomposition (COD) of a matrix.

Parameters
MatrixTypethe type of the matrix of which we are computing the COD.

This class performs a rank-revealing complete orthogonal decomposition of a matrix A into matrices P, Q, T, and Z such that

\[ \mathbf{A} \, \mathbf{P} = \mathbf{Q} \, \begin{bmatrix} \mathbf{T} & \mathbf{0} \\ \mathbf{0} & \mathbf{0} \end{bmatrix} \, \mathbf{Z} \]

by using Householder transformations. Here, P is a permutation matrix, Q and Z are unitary matrices and T an upper triangular matrix of size rank-by-rank. A may be rank deficient.

This class supports the inplace decomposition mechanism.

See also
MatrixBase::completeOrthogonalDecomposition()

Member Typedef Documentation

◆ Base

template<typename _MatrixType >
typedef SolverBase<CompleteOrthogonalDecomposition> Eigen::CompleteOrthogonalDecomposition< _MatrixType >::Base

◆ HCoeffsType

template<typename _MatrixType >
typedef internal::plain_diag_type<MatrixType>::type Eigen::CompleteOrthogonalDecomposition< _MatrixType >::HCoeffsType

◆ HouseholderSequenceType

template<typename _MatrixType >
typedef HouseholderSequence< MatrixType, typename internal::remove_all< typename HCoeffsType::ConjugateReturnType>::type> Eigen::CompleteOrthogonalDecomposition< _MatrixType >::HouseholderSequenceType

◆ IntRowVectorType

template<typename _MatrixType >
typedef internal::plain_row_type<MatrixType,Index>::type Eigen::CompleteOrthogonalDecomposition< _MatrixType >::IntRowVectorType

◆ MatrixType

template<typename _MatrixType >
typedef _MatrixType Eigen::CompleteOrthogonalDecomposition< _MatrixType >::MatrixType

◆ PermutationType

template<typename _MatrixType >
typedef PermutationMatrix<ColsAtCompileTime, MaxColsAtCompileTime> Eigen::CompleteOrthogonalDecomposition< _MatrixType >::PermutationType

◆ PlainObject

template<typename _MatrixType >
typedef MatrixType::PlainObject Eigen::CompleteOrthogonalDecomposition< _MatrixType >::PlainObject

◆ RealRowVectorType

template<typename _MatrixType >
typedef internal::plain_row_type<MatrixType,RealScalar>::type Eigen::CompleteOrthogonalDecomposition< _MatrixType >::RealRowVectorType

◆ RowVectorType

template<typename _MatrixType >
typedef internal::plain_row_type<MatrixType>::type Eigen::CompleteOrthogonalDecomposition< _MatrixType >::RowVectorType

Member Enumeration Documentation

◆ anonymous enum

template<typename _MatrixType >
anonymous enum
Enumerator
MaxRowsAtCompileTime 
MaxColsAtCompileTime 

Constructor & Destructor Documentation

◆ CompleteOrthogonalDecomposition() [1/4]

template<typename _MatrixType >
Eigen::CompleteOrthogonalDecomposition< _MatrixType >::CompleteOrthogonalDecomposition ( )
inline

Default Constructor.

The default constructor is useful in cases in which the user intends to perform decompositions via CompleteOrthogonalDecomposition::compute(const* MatrixType&).

◆ CompleteOrthogonalDecomposition() [2/4]

template<typename _MatrixType >
Eigen::CompleteOrthogonalDecomposition< _MatrixType >::CompleteOrthogonalDecomposition ( Index  rows,
Index  cols 
)
inline

Default Constructor with memory preallocation.

Like the default constructor but with preallocation of the internal data according to the specified problem size.

See also
CompleteOrthogonalDecomposition()

◆ CompleteOrthogonalDecomposition() [3/4]

template<typename _MatrixType >
template<typename InputType >
Eigen::CompleteOrthogonalDecomposition< _MatrixType >::CompleteOrthogonalDecomposition ( const EigenBase< InputType > &  matrix)
inlineexplicit

Constructs a complete orthogonal decomposition from a given matrix.

This constructor computes the complete orthogonal decomposition of the matrix matrix by calling the method compute(). The default threshold for rank determination will be used. It is a short cut for:

matrix.cols());
cod.setThreshold(Default);
cod.compute(matrix);
Complete orthogonal decomposition (COD) of a matrix.
Definition: CompleteOrthogonalDecomposition.h:52
@ Default
Definition: Constants.h:362
See also
compute()

◆ CompleteOrthogonalDecomposition() [4/4]

template<typename _MatrixType >
template<typename InputType >
Eigen::CompleteOrthogonalDecomposition< _MatrixType >::CompleteOrthogonalDecomposition ( EigenBase< InputType > &  matrix)
inlineexplicit

Constructs a complete orthogonal decomposition from a given matrix.

This overloaded constructor is provided for inplace decomposition when MatrixType is a Eigen::Ref.

See also
CompleteOrthogonalDecomposition(const EigenBase&)

Member Function Documentation

◆ _check_solve_assertion()

template<typename _MatrixType >
template<bool Transpose_, typename Rhs >
void Eigen::CompleteOrthogonalDecomposition< _MatrixType >::_check_solve_assertion ( const Rhs &  b) const
inlineprotected

◆ _solve_impl()

template<typename _MatrixType >
template<typename RhsType , typename DstType >
void Eigen::CompleteOrthogonalDecomposition< _MatrixType >::_solve_impl ( const RhsType &  rhs,
DstType &  dst 
) const

◆ _solve_impl_transposed()

template<typename _MatrixType >
template<bool Conjugate, typename RhsType , typename DstType >
void Eigen::CompleteOrthogonalDecomposition< _MatrixType >::_solve_impl_transposed ( const RhsType &  rhs,
DstType &  dst 
) const

◆ absDeterminant()

template<typename MatrixType >
MatrixType::RealScalar Eigen::CompleteOrthogonalDecomposition< MatrixType >::absDeterminant
Returns
the absolute value of the determinant of the matrix of which *this is the complete orthogonal decomposition. It has only linear complexity (that is, O(n) where n is the dimension of the square matrix) as the complete orthogonal decomposition has already been computed.
Note
This is only for square matrices.
Warning
a determinant can be very big or small, so for matrices of large enough dimension, there is a risk of overflow/underflow. One way to work around that is to use logAbsDeterminant() instead.
See also
logAbsDeterminant(), MatrixBase::determinant()

◆ applyZAdjointOnTheLeftInPlace()

template<typename MatrixType >
template<typename Rhs >
void Eigen::CompleteOrthogonalDecomposition< MatrixType >::applyZAdjointOnTheLeftInPlace ( Rhs &  rhs) const
protected

Overwrites rhs with \( \mathbf{Z}^* * \mathbf{rhs} \).

◆ applyZOnTheLeftInPlace()

template<typename MatrixType >
template<bool Conjugate, typename Rhs >
void Eigen::CompleteOrthogonalDecomposition< MatrixType >::applyZOnTheLeftInPlace ( Rhs &  rhs) const
protected

Overwrites rhs with \( \mathbf{Z} * \mathbf{rhs} \) or \( \mathbf{\overline Z} * \mathbf{rhs} \) if Conjugate is set to true.

◆ check_template_parameters()

template<typename _MatrixType >
static void Eigen::CompleteOrthogonalDecomposition< _MatrixType >::check_template_parameters ( )
inlinestaticprotected

◆ cols()

template<typename _MatrixType >
Index Eigen::CompleteOrthogonalDecomposition< _MatrixType >::cols ( void  ) const
inline

◆ colsPermutation()

template<typename _MatrixType >
const PermutationType & Eigen::CompleteOrthogonalDecomposition< _MatrixType >::colsPermutation ( ) const
inline
Returns
a const reference to the column permutation matrix

◆ compute()

template<typename _MatrixType >
template<typename InputType >
CompleteOrthogonalDecomposition & Eigen::CompleteOrthogonalDecomposition< _MatrixType >::compute ( const EigenBase< InputType > &  matrix)
inline

◆ computeInPlace()

template<typename MatrixType >
void Eigen::CompleteOrthogonalDecomposition< MatrixType >::computeInPlace
protected

Performs the complete orthogonal decomposition of the given matrix matrix.

The result of the factorization is stored into *this, and a reference to *this is returned.

See also
class CompleteOrthogonalDecomposition, CompleteOrthogonalDecomposition(const MatrixType&)

◆ dimensionOfKernel()

template<typename _MatrixType >
Index Eigen::CompleteOrthogonalDecomposition< _MatrixType >::dimensionOfKernel ( ) const
inline
Returns
the dimension of the kernel of the matrix of which *this is the complete orthogonal decomposition.
Note
This method has to determine which pivots should be considered nonzero. For that, it uses the threshold value that you can control by calling setThreshold(const RealScalar&).

◆ hCoeffs()

template<typename _MatrixType >
const HCoeffsType & Eigen::CompleteOrthogonalDecomposition< _MatrixType >::hCoeffs ( ) const
inline
Returns
a const reference to the vector of Householder coefficients used to represent the factor Q.

For advanced uses only.

◆ householderQ()

Returns
the matrix Q as a sequence of householder transformations

◆ info()

template<typename _MatrixType >
ComputationInfo Eigen::CompleteOrthogonalDecomposition< _MatrixType >::info ( ) const
inline

Reports whether the complete orthogonal decomposition was successful.

Note
This function always returns Success. It is provided for compatibility with other factorization routines.
Returns
Success

◆ isInjective()

template<typename _MatrixType >
bool Eigen::CompleteOrthogonalDecomposition< _MatrixType >::isInjective ( ) const
inline
Returns
true if the matrix of which *this is the decomposition represents an injective linear map, i.e. has trivial kernel; false otherwise.
Note
This method has to determine which pivots should be considered nonzero. For that, it uses the threshold value that you can control by calling setThreshold(const RealScalar&).

◆ isInvertible()

template<typename _MatrixType >
bool Eigen::CompleteOrthogonalDecomposition< _MatrixType >::isInvertible ( ) const
inline
Returns
true if the matrix of which *this is the complete orthogonal decomposition is invertible.
Note
This method has to determine which pivots should be considered nonzero. For that, it uses the threshold value that you can control by calling setThreshold(const RealScalar&).

◆ isSurjective()

template<typename _MatrixType >
bool Eigen::CompleteOrthogonalDecomposition< _MatrixType >::isSurjective ( ) const
inline
Returns
true if the matrix of which *this is the decomposition represents a surjective linear map; false otherwise.
Note
This method has to determine which pivots should be considered nonzero. For that, it uses the threshold value that you can control by calling setThreshold(const RealScalar&).

◆ logAbsDeterminant()

template<typename MatrixType >
MatrixType::RealScalar Eigen::CompleteOrthogonalDecomposition< MatrixType >::logAbsDeterminant
Returns
the natural log of the absolute value of the determinant of the matrix of which *this is the complete orthogonal decomposition. It has only linear complexity (that is, O(n) where n is the dimension of the square matrix) as the complete orthogonal decomposition has already been computed.
Note
This is only for square matrices.
This method is useful to work around the risk of overflow/underflow that's inherent to determinant computation.
See also
absDeterminant(), MatrixBase::determinant()

◆ matrixQ()

template<typename _MatrixType >
HouseholderSequenceType Eigen::CompleteOrthogonalDecomposition< _MatrixType >::matrixQ ( void  ) const
inline

◆ matrixQTZ()

template<typename _MatrixType >
const MatrixType & Eigen::CompleteOrthogonalDecomposition< _MatrixType >::matrixQTZ ( ) const
inline
Returns
a reference to the matrix where the complete orthogonal decomposition is stored

◆ matrixT()

template<typename _MatrixType >
const MatrixType & Eigen::CompleteOrthogonalDecomposition< _MatrixType >::matrixT ( ) const
inline
Returns
a reference to the matrix where the complete orthogonal decomposition is stored.
Warning
The strict lower part and
cols() - rank()
Index cols() const
Definition: CompleteOrthogonalDecomposition.h:285
Index rank() const
Definition: CompleteOrthogonalDecomposition.h:235
right columns of this matrix contains internal values. Only the upper triangular part should be referenced. To get it, use
matrixT().template triangularView<Upper>()
const MatrixType & matrixT() const
Definition: CompleteOrthogonalDecomposition.h:183
For rank-deficient matrices, use
matrixR().topLeftCorner(rank(), rank()).template triangularView<Upper>()

◆ matrixZ()

template<typename _MatrixType >
MatrixType Eigen::CompleteOrthogonalDecomposition< _MatrixType >::matrixZ ( ) const
inline
Returns
the matrix Z.

◆ maxPivot()

template<typename _MatrixType >
RealScalar Eigen::CompleteOrthogonalDecomposition< _MatrixType >::maxPivot ( ) const
inline
Returns
the absolute value of the biggest pivot, i.e. the biggest diagonal coefficient of R.

◆ nonzeroPivots()

template<typename _MatrixType >
Index Eigen::CompleteOrthogonalDecomposition< _MatrixType >::nonzeroPivots ( ) const
inline
Returns
the number of nonzero pivots in the complete orthogonal decomposition. Here nonzero is meant in the exact sense, not in a fuzzy sense. So that notion isn't really intrinsically interesting, but it is still useful when implementing algorithms.
See also
rank()

◆ pseudoInverse()

template<typename _MatrixType >
const Inverse< CompleteOrthogonalDecomposition > Eigen::CompleteOrthogonalDecomposition< _MatrixType >::pseudoInverse ( ) const
inline
Returns
the pseudo-inverse of the matrix of which *this is the complete orthogonal decomposition.
Warning
: Do not compute this->pseudoInverse()*rhs to solve a linear systems. It is more efficient and numerically stable to call this->solve(rhs).

◆ rank()

template<typename _MatrixType >
Index Eigen::CompleteOrthogonalDecomposition< _MatrixType >::rank ( ) const
inline
Returns
the rank of the matrix of which *this is the complete orthogonal decomposition.
Note
This method has to determine which pivots should be considered nonzero. For that, it uses the threshold value that you can control by calling setThreshold(const RealScalar&).

◆ rows()

template<typename _MatrixType >
Index Eigen::CompleteOrthogonalDecomposition< _MatrixType >::rows ( void  ) const
inline

◆ setThreshold() [1/2]

template<typename _MatrixType >
CompleteOrthogonalDecomposition & Eigen::CompleteOrthogonalDecomposition< _MatrixType >::setThreshold ( const RealScalar &  threshold)
inline

Allows to prescribe a threshold to be used by certain methods, such as rank(), who need to determine when pivots are to be considered nonzero.

Most be called before calling compute().

When it needs to get the threshold value, Eigen calls threshold(). By default, this uses a formula to automatically determine a reasonable threshold. Once you have called the present method setThreshold(const RealScalar&), your value is used instead.

Parameters
thresholdThe new value to use as the threshold.

A pivot will be considered nonzero if its absolute value is strictly greater than \( \vert pivot \vert \leqslant threshold \times \vert maxpivot \vert \) where maxpivot is the biggest pivot.

If you want to come back to the default behavior, call setThreshold(Default_t)

◆ setThreshold() [2/2]

template<typename _MatrixType >
CompleteOrthogonalDecomposition & Eigen::CompleteOrthogonalDecomposition< _MatrixType >::setThreshold ( Default_t  )
inline

Allows to come back to the default behavior, letting Eigen use its default formula for determining the threshold.

You should pass the special object Eigen::Default as parameter here.

qr.setThreshold(Eigen::Default);

See the documentation of setThreshold(const RealScalar&).

◆ threshold()

template<typename _MatrixType >
RealScalar Eigen::CompleteOrthogonalDecomposition< _MatrixType >::threshold ( ) const
inline

Returns the threshold that will be used by certain methods such as rank().

See the documentation of setThreshold(const RealScalar&).

◆ zCoeffs()

template<typename _MatrixType >
const HCoeffsType & Eigen::CompleteOrthogonalDecomposition< _MatrixType >::zCoeffs ( ) const
inline
Returns
a const reference to the vector of Householder coefficients used to represent the factor Z.

For advanced uses only.

Friends And Related Function Documentation

◆ internal::solve_assertion

template<typename _MatrixType >
template<typename Derived >
friend struct internal::solve_assertion
friend

Member Data Documentation

◆ m_cpqr

template<typename _MatrixType >
ColPivHouseholderQR<MatrixType> Eigen::CompleteOrthogonalDecomposition< _MatrixType >::m_cpqr
protected

◆ m_temp

template<typename _MatrixType >
RowVectorType Eigen::CompleteOrthogonalDecomposition< _MatrixType >::m_temp
protected

◆ m_zCoeffs

template<typename _MatrixType >
HCoeffsType Eigen::CompleteOrthogonalDecomposition< _MatrixType >::m_zCoeffs
protected

The documentation for this class was generated from the following files: