• Home
  • Line#
  • Scopes#
  • Navigate#
  • Raw
  • Download
1 /*
2  [auto_generated]
3  boost/numeric/odeint/external/openmp/openmp_nested_algebra.hpp
4 
5  [begin_description]
6  Nested parallelized algebra for OpenMP.
7  [end_description]
8 
9  Copyright 2013 Karsten Ahnert
10  Copyright 2013 Mario Mulansky
11  Copyright 2013 Pascal Germroth
12 
13  Distributed under the Boost Software License, Version 1.0.
14  (See accompanying file LICENSE_1_0.txt or
15  copy at http://www.boost.org/LICENSE_1_0.txt)
16  */
17 
18 
19 #ifndef BOOST_NUMERIC_ODEINT_EXTERNAL_OPENMP_OPENMP_NESTED_ALGEBRA_HPP_INCLUDED
20 #define BOOST_NUMERIC_ODEINT_EXTERNAL_OPENMP_OPENMP_NESTED_ALGEBRA_HPP_INCLUDED
21 
22 #include <boost/assert.hpp>
23 #include <boost/range.hpp>
24 #include <boost/numeric/odeint/algebra/norm_result_type.hpp>
25 #include <boost/numeric/odeint/util/n_ary_helper.hpp>
26 
27 namespace boost {
28 namespace numeric {
29 namespace odeint {
30 
31 /** \brief OpenMP-parallelized algebra, wrapping another, non-parallelized algebra.
32  *
33  * NestedState must be a model of Random Access Range, where the elements are sub-states
34  * which will be processed in parallel.
35  */
36 template< class InnerAlgebra >
37 struct openmp_nested_algebra
38 {
39 
40 #if __cplusplus >= 201103L // C++11 supports _Pragma
41 
42 #define BOOST_ODEINT_GEN_LOCAL(z, n, unused) \
43     BOOST_ASSERT_MSG( len == boost::size(s ## n), "All nested state ranges must have the same size." ); \
44     typename boost::range_iterator<S ## n>::type beg ## n = boost::begin(s ## n);
45 #define BOOST_ODEINT_GEN_BODY(n) \
46     const size_t len = boost::size(s0); \
47     BOOST_PP_REPEAT(n, BOOST_ODEINT_GEN_LOCAL, ~) \
48     _Pragma("omp parallel for schedule(runtime)") \
49     for( size_t i = 0 ; i < len ; i++ ) \
50         InnerAlgebra::for_each##n( \
51             BOOST_PP_ENUM_BINARY_PARAMS(n, beg, [i] BOOST_PP_INTERCEPT) , \
52             op \
53         );
BOOST_ODEINT_GEN_FOR_EACHboost::numeric::odeint::openmp_nested_algebra54 BOOST_ODEINT_GEN_FOR_EACH(BOOST_ODEINT_GEN_BODY)
55 #undef BOOST_ODEINT_GEN_BODY
56 #undef BOOST_ODEINT_GEN_LOCAL
57 
58 #else
59 
60     template< class S0 , class Op > static void for_each1 ( S0 &s0 , Op op ) {
61         const size_t len = boost::size(s0);
62         typename boost::range_iterator<S0>::type beg0 = boost::begin(s0);
63         #pragma omp parallel for schedule(runtime)
64         for( size_t i = 0 ; i < len ; i++ ) InnerAlgebra::for_each1( beg0 [i] , op );
65     }
66     template< class S0 , class S1 , class Op > static void for_each2 ( S0 &s0 , S1 &s1 , Op op ) {
67         const size_t len = boost::size(s0);
68         typename boost::range_iterator<S0>::type beg0 = boost::begin(s0);
69         typename boost::range_iterator<S1>::type beg1 = boost::begin(s1);
70         #pragma omp parallel for schedule(runtime)
71         for( size_t i = 0 ; i < len ; i++ ) InnerAlgebra::for_each2( beg0 [i] , beg1 [i] , op );
72     }
73     template< class S0 , class S1 , class S2 , class Op > static void for_each3 ( S0 &s0 , S1 &s1 , S2 &s2 , Op op ) {
74         const size_t len = boost::size(s0);
75         typename boost::range_iterator<S0>::type beg0 = boost::begin(s0);
76         typename boost::range_iterator<S1>::type beg1 = boost::begin(s1);
77         typename boost::range_iterator<S2>::type beg2 = boost::begin(s2);
78         #pragma omp parallel for schedule(runtime)
79         for( size_t i = 0 ; i < len ; i++ ) InnerAlgebra::for_each3( beg0 [i] , beg1 [i] , beg2 [i] , op );
80     }
81     template< class S0 , class S1 , class S2 , class S3 , class Op > static void for_each4 ( S0 &s0 , S1 &s1 , S2 &s2 , S3 &s3 , Op op ) {
82         const size_t len = boost::size(s0);
83         typename boost::range_iterator<S0>::type beg0 = boost::begin(s0);
84         typename boost::range_iterator<S1>::type beg1 = boost::begin(s1);
85         typename boost::range_iterator<S2>::type beg2 = boost::begin(s2);
86         typename boost::range_iterator<S3>::type beg3 = boost::begin(s3);
87         #pragma omp parallel for schedule(runtime)
88         for( size_t i = 0 ; i < len ; i++ ) InnerAlgebra::for_each4( beg0 [i] , beg1 [i] , beg2 [i] , beg3 [i] , op );
89     }
90     template< class S0 , class S1 , class S2 , class S3 , class S4 , class Op > static void for_each5 ( S0 &s0 , S1 &s1 , S2 &s2 , S3 &s3 , S4 &s4 , Op op ) {
91         const size_t len = boost::size(s0);
92         typename boost::range_iterator<S0>::type beg0 = boost::begin(s0);
93         typename boost::range_iterator<S1>::type beg1 = boost::begin(s1);
94         typename boost::range_iterator<S2>::type beg2 = boost::begin(s2);
95         typename boost::range_iterator<S3>::type beg3 = boost::begin(s3);
96         typename boost::range_iterator<S4>::type beg4 = boost::begin(s4);
97         #pragma omp parallel for schedule(runtime)
98         for( size_t i = 0 ; i < len ; i++ ) InnerAlgebra::for_each5( beg0 [i] , beg1 [i] , beg2 [i] , beg3 [i] , beg4 [i] , op );
99     }
100     template< class S0 , class S1 , class S2 , class S3 , class S4 , class S5 , class Op > static void for_each6 ( S0 &s0 , S1 &s1 , S2 &s2 , S3 &s3 , S4 &s4 , S5 &s5 , Op op ) {
101         const size_t len = boost::size(s0);
102         typename boost::range_iterator<S0>::type beg0 = boost::begin(s0);
103         typename boost::range_iterator<S1>::type beg1 = boost::begin(s1);
104         typename boost::range_iterator<S2>::type beg2 = boost::begin(s2);
105         typename boost::range_iterator<S3>::type beg3 = boost::begin(s3);
106         typename boost::range_iterator<S4>::type beg4 = boost::begin(s4);
107         typename boost::range_iterator<S5>::type beg5 = boost::begin(s5);
108         #pragma omp parallel for schedule(runtime)
109         for( size_t i = 0 ; i < len ; i++ ) InnerAlgebra::for_each6( beg0 [i] , beg1 [i] , beg2 [i] , beg3 [i] , beg4 [i] , beg5 [i] , op );
110     }
111     template< class S0 , class S1 , class S2 , class S3 , class S4 , class S5 , class S6 , class Op > static void for_each7 ( S0 &s0 , S1 &s1 , S2 &s2 , S3 &s3 , S4 &s4 , S5 &s5 , S6 &s6 , Op op ) {
112         const size_t len = boost::size(s0);
113         typename boost::range_iterator<S0>::type beg0 = boost::begin(s0);
114         typename boost::range_iterator<S1>::type beg1 = boost::begin(s1);
115         typename boost::range_iterator<S2>::type beg2 = boost::begin(s2);
116         typename boost::range_iterator<S3>::type beg3 = boost::begin(s3);
117         typename boost::range_iterator<S4>::type beg4 = boost::begin(s4);
118         typename boost::range_iterator<S5>::type beg5 = boost::begin(s5);
119         typename boost::range_iterator<S6>::type beg6 = boost::begin(s6);
120         #pragma omp parallel for schedule(runtime)
121         for( size_t i = 0 ; i < len ; i++ ) InnerAlgebra::for_each7( beg0 [i] , beg1 [i] , beg2 [i] , beg3 [i] , beg4 [i] , beg5 [i] , beg6 [i] , op );
122     }
123     template< class S0 , class S1 , class S2 , class S3 , class S4 , class S5 , class S6 , class S7 , class Op > static void for_each8 ( S0 &s0 , S1 &s1 , S2 &s2 , S3 &s3 , S4 &s4 , S5 &s5 , S6 &s6 , S7 &s7 , Op op ) {
124         const size_t len = boost::size(s0);
125         typename boost::range_iterator<S0>::type beg0 = boost::begin(s0);
126         typename boost::range_iterator<S1>::type beg1 = boost::begin(s1);
127         typename boost::range_iterator<S2>::type beg2 = boost::begin(s2);
128         typename boost::range_iterator<S3>::type beg3 = boost::begin(s3);
129         typename boost::range_iterator<S4>::type beg4 = boost::begin(s4);
130         typename boost::range_iterator<S5>::type beg5 = boost::begin(s5);
131         typename boost::range_iterator<S6>::type beg6 = boost::begin(s6);
132         typename boost::range_iterator<S7>::type beg7 = boost::begin(s7);
133         #pragma omp parallel for schedule(runtime)
134         for( size_t i = 0 ; i < len ; i++ ) InnerAlgebra::for_each8( beg0 [i] , beg1 [i] , beg2 [i] , beg3 [i] , beg4 [i] , beg5 [i] , beg6 [i] , beg7 [i] , op );
135     }
136     template< class S0 , class S1 , class S2 , class S3 , class S4 , class S5 , class S6 , class S7 , class S8 , class Op > static void for_each9 ( S0 &s0 , S1 &s1 , S2 &s2 , S3 &s3 , S4 &s4 , S5 &s5 , S6 &s6 , S7 &s7 , S8 &s8 , Op op ) {
137         const size_t len = boost::size(s0);
138         typename boost::range_iterator<S0>::type beg0 = boost::begin(s0);
139         typename boost::range_iterator<S1>::type beg1 = boost::begin(s1);
140         typename boost::range_iterator<S2>::type beg2 = boost::begin(s2);
141         typename boost::range_iterator<S3>::type beg3 = boost::begin(s3);
142         typename boost::range_iterator<S4>::type beg4 = boost::begin(s4);
143         typename boost::range_iterator<S5>::type beg5 = boost::begin(s5);
144         typename boost::range_iterator<S6>::type beg6 = boost::begin(s6);
145         typename boost::range_iterator<S7>::type beg7 = boost::begin(s7);
146         typename boost::range_iterator<S8>::type beg8 = boost::begin(s8);
147         #pragma omp parallel for schedule(runtime)
148         for( size_t i = 0 ; i < len ; i++ ) InnerAlgebra::for_each9( beg0 [i] , beg1 [i] , beg2 [i] , beg3 [i] , beg4 [i] , beg5 [i] , beg6 [i] , beg7 [i] , beg8 [i] , op );
149     }
150     template< class S0 , class S1 , class S2 , class S3 , class S4 , class S5 , class S6 , class S7 , class S8 , class S9 , class Op > static void for_each10 ( S0 &s0 , S1 &s1 , S2 &s2 , S3 &s3 , S4 &s4 , S5 &s5 , S6 &s6 , S7 &s7 , S8 &s8 , S9 &s9 , Op op ) {
151         const size_t len = boost::size(s0);
152         typename boost::range_iterator<S0>::type beg0 = boost::begin(s0);
153         typename boost::range_iterator<S1>::type beg1 = boost::begin(s1);
154         typename boost::range_iterator<S2>::type beg2 = boost::begin(s2);
155         typename boost::range_iterator<S3>::type beg3 = boost::begin(s3);
156         typename boost::range_iterator<S4>::type beg4 = boost::begin(s4);
157         typename boost::range_iterator<S5>::type beg5 = boost::begin(s5);
158         typename boost::range_iterator<S6>::type beg6 = boost::begin(s6);
159         typename boost::range_iterator<S7>::type beg7 = boost::begin(s7);
160         typename boost::range_iterator<S8>::type beg8 = boost::begin(s8);
161         typename boost::range_iterator<S9>::type beg9 = boost::begin(s9);
162         #pragma omp parallel for schedule(runtime)
163         for( size_t i = 0 ; i < len ; i++ ) InnerAlgebra::for_each10( beg0 [i] , beg1 [i] , beg2 [i] , beg3 [i] , beg4 [i] , beg5 [i] , beg6 [i] , beg7 [i] , beg8 [i] , beg9 [i] , op );
164     }
165     template< class S0 , class S1 , class S2 , class S3 , class S4 , class S5 , class S6 , class S7 , class S8 , class S9 , class S10 , class Op > static void for_each11 ( S0 &s0 , S1 &s1 , S2 &s2 , S3 &s3 , S4 &s4 , S5 &s5 , S6 &s6 , S7 &s7 , S8 &s8 , S9 &s9 , S10 &s10 , Op op ) {
166         const size_t len = boost::size(s0);
167         typename boost::range_iterator<S0>::type beg0 = boost::begin(s0);
168         typename boost::range_iterator<S1>::type beg1 = boost::begin(s1);
169         typename boost::range_iterator<S2>::type beg2 = boost::begin(s2);
170         typename boost::range_iterator<S3>::type beg3 = boost::begin(s3);
171         typename boost::range_iterator<S4>::type beg4 = boost::begin(s4);
172         typename boost::range_iterator<S5>::type beg5 = boost::begin(s5);
173         typename boost::range_iterator<S6>::type beg6 = boost::begin(s6);
174         typename boost::range_iterator<S7>::type beg7 = boost::begin(s7);
175         typename boost::range_iterator<S8>::type beg8 = boost::begin(s8);
176         typename boost::range_iterator<S9>::type beg9 = boost::begin(s9);
177         typename boost::range_iterator<S10>::type beg10 = boost::begin(s10);
178         #pragma omp parallel for schedule(runtime)
179         for( size_t i = 0 ; i < len ; i++ ) InnerAlgebra::for_each11( beg0 [i] , beg1 [i] , beg2 [i] , beg3 [i] , beg4 [i] , beg5 [i] , beg6 [i] , beg7 [i] , beg8 [i] , beg9 [i] , beg10 [i] , op );
180     }
181     template< class S0 , class S1 , class S2 , class S3 , class S4 , class S5 , class S6 , class S7 , class S8 , class S9 , class S10 , class S11 , class Op > static void for_each12 ( S0 &s0 , S1 &s1 , S2 &s2 , S3 &s3 , S4 &s4 , S5 &s5 , S6 &s6 , S7 &s7 , S8 &s8 , S9 &s9 , S10 &s10 , S11 &s11 , Op op ) {
182         const size_t len = boost::size(s0);
183         typename boost::range_iterator<S0>::type beg0 = boost::begin(s0);
184         typename boost::range_iterator<S1>::type beg1 = boost::begin(s1);
185         typename boost::range_iterator<S2>::type beg2 = boost::begin(s2);
186         typename boost::range_iterator<S3>::type beg3 = boost::begin(s3);
187         typename boost::range_iterator<S4>::type beg4 = boost::begin(s4);
188         typename boost::range_iterator<S5>::type beg5 = boost::begin(s5);
189         typename boost::range_iterator<S6>::type beg6 = boost::begin(s6);
190         typename boost::range_iterator<S7>::type beg7 = boost::begin(s7);
191         typename boost::range_iterator<S8>::type beg8 = boost::begin(s8);
192         typename boost::range_iterator<S9>::type beg9 = boost::begin(s9);
193         typename boost::range_iterator<S10>::type beg10 = boost::begin(s10);
194         typename boost::range_iterator<S11>::type beg11 = boost::begin(s11);
195         #pragma omp parallel for schedule(runtime)
196         for( size_t i = 0 ; i < len ; i++ ) InnerAlgebra::for_each12( beg0 [i] , beg1 [i] , beg2 [i] , beg3 [i] , beg4 [i] , beg5 [i] , beg6 [i] , beg7 [i] , beg8 [i] , beg9 [i] , beg10 [i] , beg11 [i] , op );
197     }
198     template< class S0 , class S1 , class S2 , class S3 , class S4 , class S5 , class S6 , class S7 , class S8 , class S9 , class S10 , class S11 , class S12 , class Op > static void for_each13 ( S0 &s0 , S1 &s1 , S2 &s2 , S3 &s3 , S4 &s4 , S5 &s5 , S6 &s6 , S7 &s7 , S8 &s8 , S9 &s9 , S10 &s10 , S11 &s11 , S12 &s12 , Op op ) {
199         const size_t len = boost::size(s0);
200         typename boost::range_iterator<S0>::type beg0 = boost::begin(s0);
201         typename boost::range_iterator<S1>::type beg1 = boost::begin(s1);
202         typename boost::range_iterator<S2>::type beg2 = boost::begin(s2);
203         typename boost::range_iterator<S3>::type beg3 = boost::begin(s3);
204         typename boost::range_iterator<S4>::type beg4 = boost::begin(s4);
205         typename boost::range_iterator<S5>::type beg5 = boost::begin(s5);
206         typename boost::range_iterator<S6>::type beg6 = boost::begin(s6);
207         typename boost::range_iterator<S7>::type beg7 = boost::begin(s7);
208         typename boost::range_iterator<S8>::type beg8 = boost::begin(s8);
209         typename boost::range_iterator<S9>::type beg9 = boost::begin(s9);
210         typename boost::range_iterator<S10>::type beg10 = boost::begin(s10);
211         typename boost::range_iterator<S11>::type beg11 = boost::begin(s11);
212         typename boost::range_iterator<S12>::type beg12 = boost::begin(s12);
213         #pragma omp parallel for schedule(runtime)
214         for( size_t i = 0 ; i < len ; i++ ) InnerAlgebra::for_each13( beg0 [i] , beg1 [i] , beg2 [i] , beg3 [i] , beg4 [i] , beg5 [i] , beg6 [i] , beg7 [i] , beg8 [i] , beg9 [i] , beg10 [i] , beg11 [i] , beg12 [i] , op );
215     }
216     template< class S0 , class S1 , class S2 , class S3 , class S4 , class S5 , class S6 , class S7 , class S8 , class S9 , class S10 , class S11 , class S12 , class S13 , class Op > static void for_each14 ( S0 &s0 , S1 &s1 , S2 &s2 , S3 &s3 , S4 &s4 , S5 &s5 , S6 &s6 , S7 &s7 , S8 &s8 , S9 &s9 , S10 &s10 , S11 &s11 , S12 &s12 , S13 &s13 , Op op ) {
217         const size_t len = boost::size(s0);
218         typename boost::range_iterator<S0>::type beg0 = boost::begin(s0);
219         typename boost::range_iterator<S1>::type beg1 = boost::begin(s1);
220         typename boost::range_iterator<S2>::type beg2 = boost::begin(s2);
221         typename boost::range_iterator<S3>::type beg3 = boost::begin(s3);
222         typename boost::range_iterator<S4>::type beg4 = boost::begin(s4);
223         typename boost::range_iterator<S5>::type beg5 = boost::begin(s5);
224         typename boost::range_iterator<S6>::type beg6 = boost::begin(s6);
225         typename boost::range_iterator<S7>::type beg7 = boost::begin(s7);
226         typename boost::range_iterator<S8>::type beg8 = boost::begin(s8);
227         typename boost::range_iterator<S9>::type beg9 = boost::begin(s9);
228         typename boost::range_iterator<S10>::type beg10 = boost::begin(s10);
229         typename boost::range_iterator<S11>::type beg11 = boost::begin(s11);
230         typename boost::range_iterator<S12>::type beg12 = boost::begin(s12);
231         typename boost::range_iterator<S13>::type beg13 = boost::begin(s13);
232         #pragma omp parallel for schedule(runtime)
233         for( size_t i = 0 ; i < len ; i++ ) InnerAlgebra::for_each14( beg0 [i] , beg1 [i] , beg2 [i] , beg3 [i] , beg4 [i] , beg5 [i] , beg6 [i] , beg7 [i] , beg8 [i] , beg9 [i] , beg10 [i] , beg11 [i] , beg12 [i] , beg13 [i] , op );
234     }
235     template< class S0 , class S1 , class S2 , class S3 , class S4 , class S5 , class S6 , class S7 , class S8 , class S9 , class S10 , class S11 , class S12 , class S13 , class S14 , class Op > static void for_each15 ( S0 &s0 , S1 &s1 , S2 &s2 , S3 &s3 , S4 &s4 , S5 &s5 , S6 &s6 , S7 &s7 , S8 &s8 , S9 &s9 , S10 &s10 , S11 &s11 , S12 &s12 , S13 &s13 , S14 &s14 , Op op ) {
236         const size_t len = boost::size(s0);
237         typename boost::range_iterator<S0>::type beg0 = boost::begin(s0);
238         typename boost::range_iterator<S1>::type beg1 = boost::begin(s1);
239         typename boost::range_iterator<S2>::type beg2 = boost::begin(s2);
240         typename boost::range_iterator<S3>::type beg3 = boost::begin(s3);
241         typename boost::range_iterator<S4>::type beg4 = boost::begin(s4);
242         typename boost::range_iterator<S5>::type beg5 = boost::begin(s5);
243         typename boost::range_iterator<S6>::type beg6 = boost::begin(s6);
244         typename boost::range_iterator<S7>::type beg7 = boost::begin(s7);
245         typename boost::range_iterator<S8>::type beg8 = boost::begin(s8);
246         typename boost::range_iterator<S9>::type beg9 = boost::begin(s9);
247         typename boost::range_iterator<S10>::type beg10 = boost::begin(s10);
248         typename boost::range_iterator<S11>::type beg11 = boost::begin(s11);
249         typename boost::range_iterator<S12>::type beg12 = boost::begin(s12);
250         typename boost::range_iterator<S13>::type beg13 = boost::begin(s13);
251         typename boost::range_iterator<S14>::type beg14 = boost::begin(s14);
252         #pragma omp parallel for schedule(runtime)
253         for( size_t i = 0 ; i < len ; i++ ) InnerAlgebra::for_each15( beg0 [i] , beg1 [i] , beg2 [i] , beg3 [i] , beg4 [i] , beg5 [i] , beg6 [i] , beg7 [i] , beg8 [i] , beg9 [i] , beg10 [i] , beg11 [i] , beg12 [i] , beg13 [i] , beg14 [i] , op );
254     }
255 
256 #endif
257 
258 
259     template< class NestedState >
260     static typename norm_result_type< typename NestedState::value_type >::type norm_inf( const NestedState &s )
261     {
262         typedef typename boost::range_iterator<const NestedState>::type iterator;
263         typedef typename std::iterator_traits<iterator>::value_type value_type;
264         typedef typename norm_result_type<value_type>::type result_type;
265         result_type init = static_cast< result_type >( 0 );
266         const size_t len = boost::size(s);
267         iterator beg = boost::begin(s);
268 #       pragma omp parallel for reduction(max: init) schedule(dynamic)
269         for( size_t i = 0 ; i < len ; i++ )
270             init = (std::max)( init , InnerAlgebra::norm_inf( beg[i] ) );
271         return init;
272     }
273 
274 };
275 
276 
277 }
278 }
279 }
280 
281 #endif
282