34template<
class ContainerType0,
class ContainerType1,
class ContainerType2>
35std::vector<double>
least_squares(
const std::vector<ContainerType0>& bs,
const ContainerType1 & b,
const ContainerType2& weights)
42 unsigned size = bs.size();
45 std::vector<double> rhs( size, 0.);
46 std::vector<double> a(size,0.);
47 for(
unsigned i=0; i<size; i++)
49 for(
unsigned j=i; j<size; j++)
51 for(
unsigned j=0; j<i; j++)
56 std::vector<unsigned> p;
66template<
class ContainerType0,
class ContainerType1>
67std::vector<double>
least_squares(
const std::vector<ContainerType0>& bs,
const ContainerType1 & b)
91template<
class ContainerType0,
class ContainerType1>
103 set_max(max, copyable0, copyable1);
106 void set_max(
unsigned max,
const ContainerType0& copyable0,
107 const ContainerType1& copyable1)
110 m_x.assign( max, copyable0);
111 m_y.assign( max, copyable1);
127 void extrapolate(
double alpha,
const ContainerType0& x,
double beta,
128 ContainerType1&
y)
const{
129 unsigned size = m_counter;
130 thrust::host_vector<double> rhs( size, 0.), a(rhs), opIi(rhs);
131 for(
unsigned i=0; i<size; i++)
135 for(
unsigned i=0; i<size; i++)
137 for(
unsigned j=0; j<size; j++)
138 opIi[j] = m_op_inv(i,j);
162 void update(
const ContainerType0& x_new,
const ContainerType1& y_new){
163 if( m_max == 0)
return;
164 unsigned size = m_counter < m_max ? m_counter + 1 : m_max;
168 for(
unsigned j=1; j<size; j++)
171 for(
unsigned i=1; i<size; i++)
172 for(
unsigned j=1; j<size; j++)
173 op(i,j) = m_op(i-1,j-1);
180 catch ( std::runtime_error & e){
183 m_op_inv = op_inv, m_op = op;
184 if( m_counter < m_max)
187 std::rotate( m_x.rbegin(), m_x.rbegin()+1, m_x.rend());
188 std::rotate( m_y.rbegin(), m_y.rbegin()+1, m_y.rend());
194 unsigned m_max, m_counter;
195 std::vector<ContainerType0> m_x;
196 std::vector<ContainerType1> m_y;
224template<
class ContainerType>
243 void set_max(
unsigned max,
const ContainerType& copyable)
246 m_x.assign( max, copyable);
275 if( m_max == 0)
return false;
276 for(
unsigned i=0; i<m_counter; i++)
277 if( fabs(t - m_t[i]) <1e-14)
292 template<
class ContainerType0>
301 value_type f0 = (t-m_t[1])*(t-m_t[2])/(m_t[0]-m_t[1])/(m_t[0]-m_t[2]);
302 value_type f1 =-(t-m_t[0])*(t-m_t[2])/(m_t[0]-m_t[1])/(m_t[1]-m_t[2]);
303 value_type f2 = (t-m_t[0])*(t-m_t[1])/(m_t[2]-m_t[0])/(m_t[2]-m_t[1]);
305 f0, m_x[0], f1, m_x[1], f2, m_x[2]);
329 template<
class ContainerType0>
338 value_type f0 =-(-2.*t+m_t[1]+m_t[2])/(m_t[0]-m_t[1])/(m_t[0]-m_t[2]);
339 value_type f1 = (-2.*t+m_t[0]+m_t[2])/(m_t[0]-m_t[1])/(m_t[1]-m_t[2]);
340 value_type f2 =-(-2.*t+m_t[0]+m_t[1])/(m_t[2]-m_t[0])/(m_t[2]-m_t[1]);
342 f0, m_x[0], f1, m_x[1], f2, m_x[2]);
359 template<
class ContainerType0>
384 template<
class MatrixType0,
class ContainerType0,
class ContainerType1>
386 ContainerType0& new_x,
const ContainerType1& weights)
388 if( m_counter < m_max)
395 m_b.assign( m_max, b);
418 for(
unsigned u=0; u<m_max; u++)
425 catch( std::runtime_error& err)
439 template<
class ContainerType0>
440 void derive( ContainerType0& dot_x)
const{
450 template<
class ContainerType0>
452 if( m_max == 0)
return;
454 for(
unsigned i=0; i<m_counter; i++)
455 if( fabs(t_new - m_t[i]) <1e-14)
460 if( m_counter < m_max)
463 std::rotate( m_x.rbegin(), m_x.rbegin()+1, m_x.rend());
464 std::rotate( m_t.rbegin(), m_t.rbegin()+1, m_t.rend());
474 template<
class ContainerType0>
475 void update(
const ContainerType0& new_entry){
477 update( t_new, new_entry);
484 const ContainerType&
head()
const{
492 const ContainerType&
tail()
const{
497 unsigned m_max, m_counter;
498 std::vector<value_type> m_t;
499 std::vector<ContainerType> m_x;
501 std::vector<ContainerType> m_b;
A square nxn matrix.
Definition operator.h:31
void copy(const ContainerTypeIn &source, ContainerTypeOut &target)
Definition blas1.h:243
void axpby(value_type alpha, const ContainerType1 &x, value_type1 beta, ContainerType &y)
Definition blas1.h:306
auto dot(const ContainerType1 &x, const ContainerType2 &y)
Binary reproducible Euclidean dot product between two vectors
Definition blas1.h:153
void evaluate(ContainerType &y, BinarySubroutine f, Functor g, const ContainerType0 &x0, const ContainerTypes &...xs)
Definition blas1.h:612
void scal(ContainerType &x, value_type alpha)
Definition blas1.h:263
void gemv(get_value_type< ContainerType1 > alpha, MatrixType &&M, const ContainerType1 &x, get_value_type< ContainerType1 > beta, ContainerType2 &y)
Alias for blas2::symv ;.
Definition blas2.h:339
auto dot(const ContainerType1 &x, const MatrixType &m, const ContainerType2 &y)
; Binary reproducible general dot product
Definition blas2.h:94
void apply(get_value_type< ContainerType1 > alpha, MatrixType &&M, const ContainerType1 &x, get_value_type< ContainerType1 > beta, ContainerType2 &y)
Alias for dg::blas2::symv ;.
Definition blas2.h:490
std::vector< const ContainerType * > asPointers(const std::vector< ContainerType > &in)
Convert a vector of vectors to a vector of pointers.
Definition densematrix.h:110
auto asDenseMatrix(const std::vector< const ContainerType * > &in)
Lightweight DenseMatrix for dg::blas2::gemv.
Definition densematrix.h:75
typename TensorTraits< std::decay_t< Vector > >::value_type get_value_type
Definition tensor_traits.h:45
dg::SquareMatrix< T > inverse(const dg::SquareMatrix< T > &in)
Invert a square matrix.
Definition operator.h:514
void lu_solve(const dg::SquareMatrix< T > &lu, const std::vector< unsigned > &p, std::vector< T > &b)
Solve the linear system with the LU decomposition.
Definition operator.h:682
T lu_pivot(dg::SquareMatrix< T > &m, std::vector< unsigned > &p)
LU Decomposition with partial pivoting.
Definition operator.h:441
This is the namespace for all functions and classes defined and used by the discontinuous Galerkin li...
Definition subroutines.h:124
Definition subroutines.h:22