Discontinuous Galerkin Library
#include "dg/algorithm.h"
Loading...
Searching...
No Matches
extrapolation.h
Go to the documentation of this file.
1#pragma once
2
3// #include <random>
4#include "blas.h"
5#include "topology/operator.h"
6
7namespace dg
8{
9
34template<class ContainerType0, class ContainerType1, class ContainerType2>
35std::vector<double> least_squares( const std::vector<ContainerType0>& bs, const ContainerType1 & b, const ContainerType2& weights)
36{
37 // it would be interesting to see how this algorithm fares against
38 // Gram-Schmidt/QR-factorization
39 // This implementation should have as many scalar dots as Gram-Schmidt
40 // namely (size^2+size)/2
41 // Solve B^T B a = B^T b
42 unsigned size = bs.size();
43 // B^T B is the "Gram matrix"
44 dg::SquareMatrix<double> op( size, 0.); // B^T B
45 std::vector<double> rhs( size, 0.);
46 std::vector<double> a(size,0.);
47 for( unsigned i=0; i<size; i++)
48 {
49 for( unsigned j=i; j<size; j++)
50 op(i,j) = dg::blas2::dot( bs[i], weights, bs[j]);
51 for( unsigned j=0; j<i; j++)
52 op(i,j) = op(j,i);
53 rhs[i] = dg::blas2::dot( bs[i], weights, b);
54 }
55 // possibly replace with Cholesky factorization?
56 std::vector<unsigned> p;
58 dg::lu_solve<double>( op, p, rhs);
59 return rhs;
60}
61
66template<class ContainerType0, class ContainerType1>
67std::vector<double> least_squares( const std::vector<ContainerType0>& bs, const ContainerType1 & b)
68{
69 return least_squares( bs, b, 1.);
70}
91template<class ContainerType0, class ContainerType1>
93{
95 using container_type = ContainerType0;
96 LeastSquaresExtrapolation( ){ m_counter = 0; }
102 LeastSquaresExtrapolation( unsigned max, const ContainerType0& copyable0, const ContainerType1& copyable1) {
103 set_max(max, copyable0, copyable1);
104 }
106 void set_max( unsigned max, const ContainerType0& copyable0,
107 const ContainerType1& copyable1)
108 {
109 m_counter = 0;
110 m_x.assign( max, copyable0);
111 m_y.assign( max, copyable1);
112 m_max = max;
113 }
116 unsigned get_max( ) const{
117 return m_counter;
118 }
119
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); // B^T b
131 for( unsigned i=0; i<size; i++)
132 rhs[i] = dg::blas1::dot( m_x[i], x);
133 // a = op_inv * rhs
134 dg::blas1::scal( y, beta);
135 for( unsigned i=0; i<size; i++)
136 {
137 for( unsigned j=0; j<size; j++)
138 opIi[j] = m_op_inv(i,j);
139 a[i] = dg::blas1::dot( rhs, opIi) ;
140 dg::blas1::axpby( alpha*a[i], m_y[i], 1., y);
141 }
142 }
148 void extrapolate( const ContainerType0& x, ContainerType1& y) const{
149 extrapolate( 1., x, 0., y);
150 }
151
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;
165 dg::SquareMatrix<double> op( size, 0.), op_inv( size, 0.); // B^T B
166 //i = 0
167 op(0,0) = dg::blas1::dot( x_new, x_new);
168 for( unsigned j=1; j<size; j++)
169 op(0,j) = op( j, 0) = dg::blas1::dot( x_new, m_x[j-1]);
170 // recursively fill in previous values
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);
174 // test if new value is linearly independent or zero
175 // maybe one can get a better control (with a tolerance value) over
176 // this test
177 try{
178 op_inv = dg::create::inverse( op);
179 }
180 catch ( std::runtime_error & e){
181 return;
182 }
183 m_op_inv = op_inv, m_op = op;
184 if( m_counter < m_max)
185 m_counter++;
186 //push out last value (keep track of what is oldest value
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());
189 blas1::copy( x_new, m_x[0]);
190 blas1::copy( y_new, m_y[0]);
191 }
192
193 private:
194 unsigned m_max, m_counter;
195 std::vector<ContainerType0> m_x;
196 std::vector<ContainerType1> m_y;
197 dg::SquareMatrix<double> m_op, m_op_inv;
198};
199
224template<class ContainerType>
226{
228 using container_type = ContainerType;
231 Extrapolation( ){ m_counter = 0; }
239 Extrapolation( unsigned max, const ContainerType& copyable) {
240 set_max(max, copyable);
241 }
243 void set_max( unsigned max, const ContainerType& copyable)
244 {
245 m_counter = 0;
246 m_x.assign( max, copyable);
247 m_t.assign( max, 0);
248 m_max = max;
249 }
257 unsigned get_max( ) const{
258 return m_counter;
259 }
260
261
274 bool exists( value_type t)const{
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)
278 return true;
279 return false;
280 }
281
292 template<class ContainerType0>
293 void extrapolate( value_type t, ContainerType0& new_x) const{
294 switch(m_counter)
295 {
296 case(0): dg::blas1::copy( 0, new_x);
297 break;
298 case(1): dg::blas1::copy( m_x[0], new_x);
299 break;
300 case(3): {
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]);
306 break;
307 }
308 default: {
309 value_type f0 = (t-m_t[1])/(m_t[0]-m_t[1]);
310 value_type f1 = (t-m_t[0])/(m_t[1]-m_t[0]);
311 dg::blas1::axpby( f0, m_x[0], f1, m_x[1], new_x);
312 }
313 }
314 }
315
329 template<class ContainerType0>
330 void derive( value_type t, ContainerType0& dot_x) const{
331 switch(m_counter)
332 {
333 case(0): dg::blas1::copy( 0, dot_x);
334 break;
335 case(1): dg::blas1::copy( 0, dot_x);
336 break;
337 case(3): {
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]);
343 break;
344 }
345 default: {
346 value_type f0 = 1./(m_t[0]-m_t[1]);
347 value_type f1 = 1./(m_t[1]-m_t[0]);
348 dg::blas1::axpby( f0, m_x[0], f1, m_x[1], dot_x);
349 }
350 }
351 }
352
359 template<class ContainerType0>
360 void extrapolate( ContainerType0& new_x) const{
361 value_type t = m_t[0] +1.;
362 extrapolate( t, new_x);
363 }
364
365
384 template<class MatrixType0, class ContainerType0, class ContainerType1>
385 void extrapolate_least_squares( MatrixType0&& A, const ContainerType0& b,
386 ContainerType0& new_x, const ContainerType1& weights)
387 {
388 if( m_counter < m_max)
389 {
390 extrapolate( new_x);
391 return;
392 }
393 // allocate m_b if not yet done
394 if( m_b.empty())
395 m_b.assign( m_max, b);
396 std::vector<const ContainerType*> x_ptrs = dg::asPointers( m_x);
397 // An attempt at the algorithm in https://arxiv.org/abs/2309.02156
398 //if( m_b.empty())
399 // m_b.assign( mm, b);
400 // if( m_zx.empty())
401 // m_zx.assign( mm, b);
402 // // compress the subspace via random linear combinations
403 // if( mm < m_max)
404 // {
405 // std::random_device rd{};
406 // std::mt19937 gen{rd()};
407 // std::normal_distribution<double> dist{0.0, 1.0};
408 // std::vector<double> z(m_max);
409 // for( unsigned i=0; i<mm; i++)
410 // {
411 // for( unsigned k=0; k<m_max; k++)
412 // z[k] = dist(gen);
413 // dg::blas2::gemv( 1., dg::asDenseMatrix(x_ptrs), z, 0., m_zx[i]);
414 // }
415 // x_ptrs = dg::asPointers( m_zx);
416 // }
417 // First compute bs
418 for( unsigned u=0; u<m_max; u++)
419 dg::apply( A, *x_ptrs[u], m_b[u]);
420
421 try{
422 std::vector<double> a = least_squares( m_b, b, weights);
423 dg::blas2::gemv( 1., dg::asDenseMatrix(x_ptrs), a, 0., new_x);
424 }
425 catch( std::runtime_error& err)
426 {
427 return extrapolate( new_x);
428 }
429
430 }
431
432
439 template<class ContainerType0>
440 void derive( ContainerType0& dot_x) const{
441 derive( m_t[0], dot_x);
442 }
443
450 template<class ContainerType0>
451 void update( value_type t_new, const ContainerType0& new_entry){
452 if( m_max == 0) return;
453 //check if entry is already there to avoid division by zero errors
454 for( unsigned i=0; i<m_counter; i++)
455 if( fabs(t_new - m_t[i]) <1e-14)
456 {
457 blas1::copy( new_entry, m_x[i]);
458 return;
459 }
460 if( m_counter < m_max) //don't update counter if Time entry was rejected
461 m_counter++;
462 //push out last value (keep track of what is oldest value
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());
465 m_t[0] = t_new;
466 blas1::copy( new_entry, m_x[0]);
467 }
474 template<class ContainerType0>
475 void update( const ContainerType0& new_entry){
476 value_type t_new = m_t[0] + 1;
477 update( t_new, new_entry);
478 }
479
484 const ContainerType& head()const{
485 return m_x[0];
486 }
488 ContainerType& tail(){
489 return m_x[m_max-1];
490 }
492 const ContainerType& tail()const{
493 return m_x[m_max-1];
494 }
495
496 private:
497 unsigned m_max, m_counter;
498 std::vector<value_type> m_t;
499 std::vector<ContainerType> m_x;
500 //only allocated if least squares extrapolate is used
501 std::vector<ContainerType> m_b;
502};
503
504
505}//namespace dg
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
@ y
y direction
@ x
x direction
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
std::vector< double > least_squares(const std::vector< ContainerType0 > &bs, const ContainerType1 &b, const ContainerType2 &weights)
Compute for given , weights and right hand side .
Definition extrapolation.h:35
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...
Extrapolate a polynomial passing through up to three points.
Definition extrapolation.h:226
void extrapolate_least_squares(MatrixType0 &&A, const ContainerType0 &b, ContainerType0 &new_x, const ContainerType1 &weights)
EXPERIMENTAL Perform a least squares extrapolation.
Definition extrapolation.h:385
void derive(value_type t, ContainerType0 &dot_x) const
Evaluate first derivative of interpolating polynomial
Definition extrapolation.h:330
void extrapolate(value_type t, ContainerType0 &new_x) const
Extrapolate value to given time.
Definition extrapolation.h:293
get_value_type< ContainerType > value_type
Definition extrapolation.h:227
bool exists(value_type t) const
Check if time exists in current points.
Definition extrapolation.h:274
unsigned get_max() const
Current extrapolation count.
Definition extrapolation.h:257
Extrapolation(unsigned max, const ContainerType &copyable)
Set maximum extrapolation order and allocate memory.
Definition extrapolation.h:239
void extrapolate(ContainerType0 &new_x) const
Extrapolate value (equidistant version)
Definition extrapolation.h:360
ContainerType & tail()
DEPRECATED write access to tail value ( the one that will be deleted in the next update,...
Definition extrapolation.h:488
void derive(ContainerType0 &dot_x) const
Evaluate first derivative of interpolating polynomial (equidistant version)
Definition extrapolation.h:440
const ContainerType & tail() const
DEPRECATED read access to tail value ( the one that will be deleted in the next update,...
Definition extrapolation.h:492
const ContainerType & head() const
return the current head (the one most recently inserted)
Definition extrapolation.h:484
ContainerType container_type
Definition extrapolation.h:228
void update(const ContainerType0 &new_entry)
insert a new entry
Definition extrapolation.h:475
Extrapolation()
Leave values uninitialized.
Definition extrapolation.h:231
void set_max(unsigned max, const ContainerType &copyable)
Set maximum extrapolation order and allocate memory.
Definition extrapolation.h:243
void update(value_type t_new, const ContainerType0 &new_entry)
insert a new entry, deleting the oldest entry or update existing entry
Definition extrapolation.h:451
Evaluate a least squares fit
Definition extrapolation.h:93
void set_max(unsigned max, const ContainerType0 &copyable0, const ContainerType1 &copyable1)
Set maximum number of vectors and allocate memory.
Definition extrapolation.h:106
void extrapolate(double alpha, const ContainerType0 &x, double beta, ContainerType1 &y) const
extrapolate value at a new unkown value
Definition extrapolation.h:127
LeastSquaresExtrapolation(unsigned max, const ContainerType0 &copyable0, const ContainerType1 &copyable1)
Set maximum number of vectors and allocate memory.
Definition extrapolation.h:102
void update(const ContainerType0 &x_new, const ContainerType1 &y_new)
insert a new entry / train the machine learning algorithm
Definition extrapolation.h:162
unsigned get_max() const
Definition extrapolation.h:116
LeastSquaresExtrapolation()
Definition extrapolation.h:96
get_value_type< ContainerType0 > value_type
Definition extrapolation.h:94
ContainerType0 container_type
Definition extrapolation.h:95
void extrapolate(const ContainerType0 &x, ContainerType1 &y) const
extrapolate value at a new unkown value
Definition extrapolation.h:148
Definition subroutines.h:124
Definition subroutines.h:22