ALP User Documentation 0.7.0
Algebraic Programming User Documentation
sparse_nn_single_inference.hpp
Go to the documentation of this file.
1
2/*
3 * Copyright 2021 Huawei Technologies Co., Ltd.
4 *
5 * Licensed under the Apache License, Version 2.0 (the "License");
6 * you may not use this file except in compliance with the License.
7 * You may obtain a copy of the License at
8 *
9 * http://www.apache.org/licenses/LICENSE-2.0
10 *
11 * Unless required by applicable law or agreed to in writing, software
12 * distributed under the License is distributed on an "AS IS" BASIS,
13 * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
14 * See the License for the specific language governing permissions and
15 * limitations under the License.
16 */
17
26#ifndef _H_GRB_ALGORITHMS_SPARSE_NN_SINGLE_INFERENCE
27#define _H_GRB_ALGORITHMS_SPARSE_NN_SINGLE_INFERENCE
28
29#include <limits>
30#include <graphblas.hpp>
31
32
33namespace grb {
34
35 namespace algorithms {
36
37 namespace internal {
38
47 template<
48 Descriptor descr,
49 bool thresholded, typename ThresholdType,
50 typename IOType, typename WeightType, typename BiasType,
51 class ReluMonoid, class Ring, class MinMonoid
52 >
55 const grb::Vector< IOType > &in,
56 const std::vector< grb::Matrix< WeightType > > &layers,
57 const std::vector< BiasType > &biases,
58 const ThresholdType threshold,
60 const ReluMonoid &relu,
61 const MinMonoid &min,
62 const Ring &ring
63 ) {
64 static_assert( !(descr & descriptors::no_casting) ||
65 (
66 std::is_same< IOType, WeightType >::value &&
67 std::is_same< IOType, BiasType >::value
68 ), "Input containers have different domains even though the no_casting"
69 "descriptor was given"
70 );
71
72 const size_t num_layers = layers.size();
73
74 // run-time checks
75 {
76 const size_t n = grb::size( out );
77 if( num_layers == 0 ) {
78 return ILLEGAL;
79 }
80 if( biases.size() != num_layers ) {
81 return ILLEGAL;
82 }
83 if( grb::size( in ) != grb::nrows( ( layers[ 0 ] ) ) ||
84 grb::size( out ) != grb::ncols( ( layers[ num_layers - 1 ] ) ) ||
85 grb::size( out ) != grb::size( temp )
86 ) {
87 return MISMATCH;
88 }
89 for( size_t i = 1; i < num_layers; ++i ) {
90 if( grb::ncols( ( layers[ i - 1 ] ) ) != grb::nrows( ( layers[ i ] ) ) ) {
91 return MISMATCH;
92 }
93 }
94 for( size_t i = 0; i < num_layers; ++i ) {
95 if( grb::ncols( ( layers[ i ] ) ) != grb::nrows( ( layers[ i ] ) ) ) {
96 return ILLEGAL;
97 }
98 }
99 assert( n == grb::size( in ) );
100 assert( n == grb::size( temp ) );
101 if( grb::capacity( out ) != n ) {
102 return ILLEGAL;
103 }
104 if( grb::capacity( temp ) != n ) {
105 return ILLEGAL;
106 }
107 }
108
109 grb::RC ret = SUCCESS;
110
111 /*
112 // this is a correct implementation that does not unroll the first and the last iterations
113 // we do not use it because it requires setting the input vector to the output vector
114 // which results in copying data for 2*n elements
115
116 ret = grb::set( out, in );
117
118 for( size_t i = 1; ret == SUCCESS && i < num_layers ; ++i ) {
119
120 std::swap( out, temp );
121 ret = ret ? ret : grb::set( out, 0 );
122 ret = ret ? ret : grb::vxm( out, temp, *(layers[ i - 1 ]), ring );
123 ret = ret ? ret : grb::foldl< descriptors::dense >( out, biases[ i ], ring.getAdditiveMonoid() );
124 ret = ret ? ret : grb::foldl< descriptors::dense >( out, 0, relu );
125 if( thresholded ) {
126 ret = ret ? ret : grb::foldl< descriptors::dense >( out, threshold, min );
127 }
128 }
129 */
130
131 ret = grb::set( out, 0 ); assert( ret == SUCCESS );
132
133 ret = ret ? ret : grb::vxm( out, in, ( layers[ 0 ] ), ring );
134 assert( ret == SUCCESS );
135
136 ret = ret ? ret : grb::foldl< descriptors::dense >(
137 out, biases[ 1 ], ring.getAdditiveMonoid()
138 ); assert( ret == SUCCESS );
139
140 for( size_t i = 1; ret == SUCCESS && i < num_layers - 1; ++i ) {
141
142 ret = ret ? ret : grb::foldl< descriptors::dense >( out, 0, relu );
143 assert( ret == SUCCESS );
144
145 if( thresholded ) {
146 ret = ret ? ret : grb::foldl< descriptors::dense >( out, threshold, min );
147 assert( ret == SUCCESS );
148 }
149
150 if( ret == SUCCESS ) {
151 std::swap( out, temp );
152 }
153
154 ret = grb::set( out, 0 );
155 assert( ret == SUCCESS );
156
157 ret = ret ? ret : grb::vxm< descriptors::dense >(
158 out, temp, ( layers[ i ] ), ring
159 ); assert( ret == SUCCESS );
160
161 ret = ret ? ret : grb::foldl< descriptors::dense >(
162 out, biases[ i + 1 ], ring.getAdditiveMonoid()
163 ); assert( ret == SUCCESS );
164 }
165
166 ret = ret ? ret : grb::foldl< descriptors::dense >( out, 0, relu );
167 assert( ret == SUCCESS );
168
169 if( thresholded ) {
170 ret = ret ? ret : grb::foldl< descriptors::dense >( out, threshold, min );
171 assert( ret == SUCCESS );
172 }
173
174 return ret;
175 }
176
177 } // end namespace ``grb::internal''
178
257 template< Descriptor descr = descriptors::no_operation,
258 typename IOType,
259 typename WeightType,
260 typename BiasType,
261 class ReluMonoid = Monoid<
264 >,
265 class Ring = Semiring<
268 >
269 >
272 const grb::Vector< IOType > &in,
273 const std::vector< grb::Matrix< WeightType > > &layers,
274 const std::vector< BiasType > &biases,
276 const ReluMonoid &relu = ReluMonoid(),
277 const Ring &ring = Ring()
278 ) {
279 static_assert( !(descr & descriptors::no_casting) ||
280 (
281 std::is_same< IOType, WeightType >::value &&
282 std::is_same< IOType, BiasType >::value
283 ), "Input containers have different domains even though the no_casting "
284 "descriptor was given" );
285 Monoid<
287 > dummyTresholdMonoid;
289 descr, false, double
290 > (
291 out, in, layers,
292 biases, 0.0,
293 temp,
294 relu, dummyTresholdMonoid, ring
295 );
296 }
297
383 template< Descriptor descr = descriptors::no_operation,
384 typename IOType,
385 typename WeightType,
386 typename BiasType,
387 typename ThresholdType = IOType,
388 class MinMonoid = Monoid<
390 >,
391 class ReluMonoid = Monoid<
394 >,
395 class Ring = Semiring<
398 >
399 >
402 const grb::Vector< IOType > &in,
403 const std::vector< grb::Matrix< WeightType > > &layers,
404 const std::vector< BiasType > &biases,
405 const ThresholdType threshold,
407 const ReluMonoid &relu = ReluMonoid(),
408 const MinMonoid &min = MinMonoid(),
409 const Ring &ring = Ring()
410 ) {
411 static_assert( !(descr & descriptors::no_casting) ||
412 (
413 std::is_same< IOType, WeightType >::value &&
414 std::is_same< IOType, BiasType >::value
415 ), "Input containers have different domains even though the no_casting "
416 "descriptor was given" );
418 descr, true
419 > (
420 out, in, layers,
421 biases, threshold,
422 temp,
423 relu, min, ring
424 );
425 }
426
427 } // namespace algorithms
428
429} // end namespace grb
430
431#endif // end _H_GRB_ALGORITHMS_SPARSE_NN_SINGLE_INFERENCE
432
An ALP/GraphBLAS matrix.
Definition: matrix.hpp:71
A generalised monoid.
Definition: monoid.hpp:54
A generalised semiring.
Definition: semiring.hpp:186
A GraphBLAS vector.
Definition: vector.hpp:64
Standard identity for the minimum operator.
Definition: identities.hpp:101
Standard identity for the maximum operator.
Definition: identities.hpp:124
Standard identity for numerical multiplication.
Definition: identities.hpp:79
Standard identity for numerical addition.
Definition: identities.hpp:57
This operator takes the sum of the two input parameters and writes it to the output variable.
Definition: ops.hpp:177
This operator takes the minimum of the two input parameters and writes the result to the output varia...
Definition: ops.hpp:276
This operator multiplies the two input parameters and writes the result to the output variable.
Definition: ops.hpp:210
This operation is equivalent to grb::operators::min.
Definition: ops.hpp:516
The main header to include in order to use the ALP/GraphBLAS API.
RC vxm(Vector< IOType, backend, Coords > &u, const Vector< InputType3, backend, Coords > &u_mask, const Vector< InputType1, backend, Coords > &v, const Vector< InputType4, backend, Coords > &v_mask, const Matrix< InputType2, backend, RIT, CIT, NIT > &A, const Semiring &semiring=Semiring(), const Phase &phase=EXECUTE, typename std::enable_if< grb::is_semiring< Semiring >::value &&!grb::is_object< InputType1 >::value &&!grb::is_object< InputType2 >::value &&!grb::is_object< InputType3 >::value &&!grb::is_object< InputType4 >::value &&!grb::is_object< IOType >::value, void >::type *=nullptr)
Left-handed in-place doubly-masked sparse matrix times vector multiplication, .
Definition: blas2.hpp:307
RC set(Vector< DataType, backend, Coords > &x, const T val, const Phase &phase=EXECUTE, const typename std::enable_if< !grb::is_object< DataType >::value &&!grb::is_object< T >::value, void >::type *const =nullptr) noexcept
Sets all elements of a vector to the given value.
Definition: io.hpp:857
size_t size(const Vector< DataType, backend, Coords > &x) noexcept
Request the size of a given vector.
Definition: io.hpp:235
size_t ncols(const Matrix< InputType, backend, RIT, CIT, NIT > &A) noexcept
Requests the column size of a given matrix.
Definition: io.hpp:339
size_t capacity(const Vector< InputType, backend, Coords > &x) noexcept
Queries the capacity of the given ALP/GraphBLAS container.
Definition: io.hpp:388
size_t nrows(const Matrix< InputType, backend, RIT, CIT, NIT > &A) noexcept
Requests the row size of a given matrix.
Definition: io.hpp:286
grb::RC sparse_nn_single_inference(grb::Vector< IOType > &out, const grb::Vector< IOType > &in, const std::vector< grb::Matrix< WeightType > > &layers, const std::vector< BiasType > &biases, const ThresholdType threshold, grb::Vector< IOType > &temp, const ReluMonoid &relu=ReluMonoid(), const MinMonoid &min=MinMonoid(), const Ring &ring=Ring())
Performs an inference step of a single data element through a Sparse Neural Network defined by num_la...
Definition: sparse_nn_single_inference.hpp:400
grb::RC sparse_nn_single_inference(grb::Vector< IOType > &out, const grb::Vector< IOType > &in, const std::vector< grb::Matrix< WeightType > > &layers, const std::vector< BiasType > &biases, grb::Vector< IOType > &temp, const ReluMonoid &relu=ReluMonoid(), const Ring &ring=Ring())
Performs an inference step of a single data element through a Sparse Neural Network defined by num_la...
Definition: sparse_nn_single_inference.hpp:270
static constexpr Descriptor no_casting
Disallows the standard casting of input parameters to a compatible domain in case they did not match ...
Definition: descriptors.hpp:196
static constexpr Descriptor no_operation
Indicates no additional pre- or post-processing on any of the GraphBLAS function arguments.
Definition: descriptors.hpp:63
The ALP/GraphBLAS namespace.
Definition: graphblas.hpp:452
RC
Return codes of ALP primitives.
Definition: rc.hpp:47
@ ILLEGAL
A call to a primitive has determined that one of its arguments was illegal as per the specification o...
Definition: rc.hpp:143
@ MISMATCH
One or more of the ALP/GraphBLAS objects passed to the primitive that returned this error have mismat...
Definition: rc.hpp:90
@ SUCCESS
Indicates the primitive has executed successfully.
Definition: rc.hpp:54
unsigned int Descriptor
Descriptors indicate pre- or post-processing for some or all of the arguments to an ALP/GraphBLAS cal...
Definition: descriptors.hpp:54