1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211
|
/*
[auto_generated]
libs/numeric/odeint/test/multi_array.cpp
[begin_description]
tba.
[end_description]
Copyright 2009-2012 Karsten Ahnert
Copyright 2009-2012 Mario Mulansky
Distributed under the Boost Software License, Version 1.0.
(See accompanying file LICENSE_1_0.txt or
copy at http://www.boost.org/LICENSE_1_0.txt)
*/
#include <boost/config.hpp>
#ifdef BOOST_MSVC
#pragma warning(disable:4996)
#endif
#define BOOST_TEST_MODULE odeint_multi_array
#include <boost/test/unit_test.hpp>
#include <boost/numeric/odeint/algebra/multi_array_algebra.hpp>
#include <boost/numeric/odeint/util/multi_array_adaption.hpp>
#include <boost/numeric/odeint/stepper/runge_kutta4.hpp>
#include <boost/numeric/odeint/integrate/integrate_const.hpp>
#include <boost/numeric/odeint/stepper/runge_kutta_dopri5.hpp>
#include <boost/numeric/odeint/stepper/generation.hpp>
using namespace boost::unit_test;
using namespace boost::numeric::odeint;
typedef boost::multi_array< double , 1 > vector_type;
typedef boost::multi_array< double , 2 > matrix_type;
typedef boost::multi_array< double , 3 > tensor_type;
BOOST_AUTO_TEST_SUITE( multi_array_test )
BOOST_AUTO_TEST_CASE( test_resizeable )
{
BOOST_CHECK( is_resizeable< vector_type >::value );
BOOST_CHECK( is_resizeable< matrix_type >::value );
BOOST_CHECK( is_resizeable< tensor_type >::value );
}
BOOST_AUTO_TEST_CASE( test_same_size_vector1 )
{
vector_type v1( boost::extents[12] );
vector_type v2( boost::extents[12] );
BOOST_CHECK( same_size( v1 , v2 ) );
}
BOOST_AUTO_TEST_CASE( test_same_size_vector2 )
{
vector_type v1( boost::extents[12] );
vector_type v2( boost::extents[13] );
BOOST_CHECK( !same_size( v1 , v2 ) );
}
BOOST_AUTO_TEST_CASE( test_same_size_vector3 )
{
vector_type v1( boost::extents[12] );
vector_type v2( boost::extents[vector_type::extent_range(-1,11)] );
BOOST_CHECK( !same_size( v1 , v2 ) );
}
BOOST_AUTO_TEST_CASE( test_same_size_vector4 )
{
vector_type v1( boost::extents[12] );
vector_type v2( boost::extents[vector_type::extent_range(-1,10)] );
BOOST_CHECK( !same_size( v1 , v2 ) );
}
BOOST_AUTO_TEST_CASE( test_same_size_vector5 )
{
vector_type v1( boost::extents[vector_type::extent_range(-1,10)] );
vector_type v2( boost::extents[vector_type::extent_range(-1,10)] );
BOOST_CHECK( same_size( v1 , v2 ) );
}
BOOST_AUTO_TEST_CASE( test_same_size_matrix1 )
{
matrix_type m1( boost::extents[12][4] );
matrix_type m2( boost::extents[12][4] );
BOOST_CHECK( same_size( m1 , m2 ) );
}
BOOST_AUTO_TEST_CASE( test_same_size_matrix2 )
{
matrix_type m1( boost::extents[12][4] );
matrix_type m2( boost::extents[12][3] );
BOOST_CHECK( !same_size( m1 , m2 ) );
}
BOOST_AUTO_TEST_CASE( test_same_size_matrix3 )
{
matrix_type m1( boost::extents[12][matrix_type::extent_range(-1,2)] );
matrix_type m2( boost::extents[12][3] );
BOOST_CHECK( !same_size( m1 , m2 ) );
}
BOOST_AUTO_TEST_CASE( test_same_size_matrix4 )
{
matrix_type m1( boost::extents[12][matrix_type::extent_range(-1,1)] );
matrix_type m2( boost::extents[12][3] );
BOOST_CHECK( !same_size( m1 , m2 ) );
}
BOOST_AUTO_TEST_CASE( test_same_size_tensor1 )
{
tensor_type t1( boost::extents[ tensor_type::extent_range( -2 , 9 ) ]
[ tensor_type::extent_range( 5 , 11 ) ]
[ tensor_type::extent_range( -1 , 2 ) ] );
tensor_type t2( boost::extents[ tensor_type::extent_range( 2 , 13 ) ]
[ tensor_type::extent_range( -1 , 5 ) ]
[ tensor_type::extent_range( 1 , 4 ) ] );
BOOST_CHECK( !same_size( t1 , t2 ) );
}
BOOST_AUTO_TEST_CASE( test_same_size_tensor2 )
{
tensor_type t1( boost::extents[ tensor_type::extent_range( -2 , 9 ) ]
[ tensor_type::extent_range( 5 , 11 ) ]
[ tensor_type::extent_range( -1 , 2 ) ] );
tensor_type t2( boost::extents[ tensor_type::extent_range( -2 , 9 ) ]
[ tensor_type::extent_range( 5 , 11 ) ]
[ tensor_type::extent_range( -1 , 2 ) ] );
BOOST_CHECK( same_size( t1 , t2 ) );
}
// Tests for tensor_type
BOOST_AUTO_TEST_CASE( test_resize_vector1 )
{
vector_type v1( boost::extents[4] );
vector_type v2;
resize( v2 , v1 );
BOOST_CHECK( same_size( v1 , v2 ) );
}
BOOST_AUTO_TEST_CASE( test_resize_vector2 )
{
vector_type v1( boost::extents[ vector_type::extent_range( -2 , 9 ) ] );
vector_type v2;
resize( v2 , v1 );
BOOST_CHECK( same_size( v1 , v2 ) );
}
BOOST_AUTO_TEST_CASE( test_resize_vector3 )
{
vector_type v1( boost::extents[ vector_type::extent_range( -2 , 9 ) ] );
vector_type v2( boost::extents[ vector_type::extent_range( 2 , 9 ) ] );
BOOST_CHECK( !same_size( v1 , v2 ) );
resize( v2 , v1 );
BOOST_CHECK( same_size( v1 , v2 ) );
}
struct mult4
{
void operator()( double &a , double const &b ) const
{
a = b * 4.0;
}
};
BOOST_AUTO_TEST_CASE( test_for_each2_vector )
{
vector_type v1( boost::extents[ vector_type::extent_range( -2 , 9 ) ] );
vector_type v2( boost::extents[ vector_type::extent_range( 2 , 13 ) ] );
for( int i=-2 ; i<9 ; ++i )
v1[i] = i * 2;
multi_array_algebra::for_each2( v2 , v1 , mult4() );
for( int i=2 ; i<13 ; ++i )
BOOST_CHECK_EQUAL( v2[i] , v1[i-4]*4.0 );
}
struct vector_ode
{
const static size_t n = 128;
void operator()( const vector_type &x , vector_type &dxdt , double t ) const
{
dxdt[-1] = x[n] - 2.0 * x[-1] + x[0];
for( size_t i=0 ; i<n ; ++i )
dxdt[i] = x[i-1] - 2.0 * x[i] + x[i+1];
dxdt[n] = x[-1] - 2.0 * x[n] + x[n-1];
}
};
BOOST_AUTO_TEST_CASE( test_rk4_vector )
{
vector_type v1( boost::extents[ vector_type::extent_range( -1 , vector_ode::n+1 ) ] );
integrate_const( runge_kutta4< vector_type >() , vector_ode() , v1 , 0.0 , 1.0 , 0.01 );
}
BOOST_AUTO_TEST_CASE( test_dopri5_vector )
{
vector_type v1( boost::extents[ vector_type::extent_range( -1 , vector_ode::n+1 ) ] );
integrate_const( make_dense_output( 1.0e-6 , 1.0e-6 , runge_kutta_dopri5< vector_type >() ) , vector_ode() , v1 , 0.0 , 1.0 , 0.01 );
}
BOOST_AUTO_TEST_SUITE_END()
|