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
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
|
// SPDX-License-Identifier: Apache-2.0
//
// Copyright 2008-2016 Conrad Sanderson (http://conradsanderson.id.au)
// Copyright 2008-2016 National ICT Australia (NICTA)
//
// Licensed under the Apache License, Version 2.0 (the "License");
// you may not use this file except in compliance with the License.
// You may obtain a copy of the License at
// http://www.apache.org/licenses/LICENSE-2.0
//
// Unless required by applicable law or agreed to in writing, software
// distributed under the License is distributed on an "AS IS" BASIS,
// WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
// See the License for the specific language governing permissions and
// limitations under the License.
// ------------------------------------------------------------------------
//! \addtogroup gemm_mixed
//! @{
//! \brief
//! Matrix multplication where the matrices have differing element types.
//! Uses caching for speedup.
//! Matrix 'C' is assumed to have been set to the correct size (ie. taking into account transposes)
template<const bool do_trans_A=false, const bool do_trans_B=false, const bool use_alpha=false, const bool use_beta=false>
class gemm_mixed_large
{
public:
template<typename out_eT, typename in_eT1, typename in_eT2>
arma_hot
inline
static
void
apply
(
Mat<out_eT>& C,
const Mat<in_eT1>& A,
const Mat<in_eT2>& B,
const out_eT alpha = out_eT(1),
const out_eT beta = out_eT(0)
)
{
arma_extra_debug_sigprint();
const uword A_n_rows = A.n_rows;
const uword A_n_cols = A.n_cols;
const uword B_n_rows = B.n_rows;
const uword B_n_cols = B.n_cols;
if( (do_trans_A == false) && (do_trans_B == false) )
{
podarray<in_eT1> tmp(A_n_cols);
in_eT1* A_rowdata = tmp.memptr();
#if defined(ARMA_USE_OPENMP)
const bool use_mp = (B_n_cols >= 2) && (B.n_elem >= 8192) && (mp_thread_limit::in_parallel() == false);
#else
const bool use_mp = false;
#endif
if(use_mp)
{
#if defined(ARMA_USE_OPENMP)
{
const int n_threads = int( (std::min)( uword(mp_thread_limit::get()), uword(B_n_cols) ) );
for(uword row_A=0; row_A < A_n_rows; ++row_A)
{
tmp.copy_row(A, row_A);
#pragma omp parallel for schedule(static) num_threads(n_threads)
for(uword col_B=0; col_B < B_n_cols; ++col_B)
{
const in_eT2* B_coldata = B.colptr(col_B);
out_eT acc = out_eT(0);
for(uword i=0; i < B_n_rows; ++i)
{
acc += upgrade_val<in_eT1,in_eT2>::apply(A_rowdata[i]) * upgrade_val<in_eT1,in_eT2>::apply(B_coldata[i]);
}
if( (use_alpha == false) && (use_beta == false) ) { C.at(row_A,col_B) = acc; }
else if( (use_alpha == true ) && (use_beta == false) ) { C.at(row_A,col_B) = alpha*acc; }
else if( (use_alpha == false) && (use_beta == true ) ) { C.at(row_A,col_B) = acc + beta*C.at(row_A,col_B); }
else if( (use_alpha == true ) && (use_beta == true ) ) { C.at(row_A,col_B) = alpha*acc + beta*C.at(row_A,col_B); }
}
}
}
#endif
}
else
{
for(uword row_A=0; row_A < A_n_rows; ++row_A)
{
tmp.copy_row(A, row_A);
for(uword col_B=0; col_B < B_n_cols; ++col_B)
{
const in_eT2* B_coldata = B.colptr(col_B);
out_eT acc = out_eT(0);
for(uword i=0; i < B_n_rows; ++i)
{
acc += upgrade_val<in_eT1,in_eT2>::apply(A_rowdata[i]) * upgrade_val<in_eT1,in_eT2>::apply(B_coldata[i]);
}
if( (use_alpha == false) && (use_beta == false) ) { C.at(row_A,col_B) = acc; }
else if( (use_alpha == true ) && (use_beta == false) ) { C.at(row_A,col_B) = alpha*acc; }
else if( (use_alpha == false) && (use_beta == true ) ) { C.at(row_A,col_B) = acc + beta*C.at(row_A,col_B); }
else if( (use_alpha == true ) && (use_beta == true ) ) { C.at(row_A,col_B) = alpha*acc + beta*C.at(row_A,col_B); }
}
}
}
}
else
if( (do_trans_A == true) && (do_trans_B == false) )
{
#if defined(ARMA_USE_OPENMP)
const bool use_mp = (B_n_cols >= 2) && (B.n_elem >= 8192) && (mp_thread_limit::in_parallel() == false);
#else
const bool use_mp = false;
#endif
if(use_mp)
{
#if defined(ARMA_USE_OPENMP)
{
const int n_threads = int( (std::min)( uword(mp_thread_limit::get()), uword(B_n_cols) ) );
for(uword col_A=0; col_A < A_n_cols; ++col_A)
{
// col_A is interpreted as row_A when storing the results in matrix C
const in_eT1* A_coldata = A.colptr(col_A);
#pragma omp parallel for schedule(static) num_threads(n_threads)
for(uword col_B=0; col_B < B_n_cols; ++col_B)
{
const in_eT2* B_coldata = B.colptr(col_B);
out_eT acc = out_eT(0);
for(uword i=0; i < B_n_rows; ++i)
{
acc += upgrade_val<in_eT1,in_eT2>::apply(A_coldata[i]) * upgrade_val<in_eT1,in_eT2>::apply(B_coldata[i]);
}
if( (use_alpha == false) && (use_beta == false) ) { C.at(col_A,col_B) = acc; }
else if( (use_alpha == true ) && (use_beta == false) ) { C.at(col_A,col_B) = alpha*acc; }
else if( (use_alpha == false) && (use_beta == true ) ) { C.at(col_A,col_B) = acc + beta*C.at(col_A,col_B); }
else if( (use_alpha == true ) && (use_beta == true ) ) { C.at(col_A,col_B) = alpha*acc + beta*C.at(col_A,col_B); }
}
}
}
#endif
}
else
{
for(uword col_A=0; col_A < A_n_cols; ++col_A)
{
// col_A is interpreted as row_A when storing the results in matrix C
const in_eT1* A_coldata = A.colptr(col_A);
for(uword col_B=0; col_B < B_n_cols; ++col_B)
{
const in_eT2* B_coldata = B.colptr(col_B);
out_eT acc = out_eT(0);
for(uword i=0; i < B_n_rows; ++i)
{
acc += upgrade_val<in_eT1,in_eT2>::apply(A_coldata[i]) * upgrade_val<in_eT1,in_eT2>::apply(B_coldata[i]);
}
if( (use_alpha == false) && (use_beta == false) ) { C.at(col_A,col_B) = acc; }
else if( (use_alpha == true ) && (use_beta == false) ) { C.at(col_A,col_B) = alpha*acc; }
else if( (use_alpha == false) && (use_beta == true ) ) { C.at(col_A,col_B) = acc + beta*C.at(col_A,col_B); }
else if( (use_alpha == true ) && (use_beta == true ) ) { C.at(col_A,col_B) = alpha*acc + beta*C.at(col_A,col_B); }
}
}
}
}
else
if( (do_trans_A == false) && (do_trans_B == true) )
{
Mat<in_eT2> B_tmp;
op_strans::apply_mat_noalias(B_tmp, B);
gemm_mixed_large<false, false, use_alpha, use_beta>::apply(C, A, B_tmp, alpha, beta);
}
else
if( (do_trans_A == true) && (do_trans_B == true) )
{
// mat B_tmp = trans(B);
// dgemm_arma<true, false, use_alpha, use_beta>::apply(C, A, B_tmp, alpha, beta);
// By using the trans(A)*trans(B) = trans(B*A) equivalency,
// transpose operations are not needed
podarray<in_eT2> tmp(B_n_cols);
in_eT2* B_rowdata = tmp.memptr();
for(uword row_B=0; row_B < B_n_rows; ++row_B)
{
tmp.copy_row(B, row_B);
for(uword col_A=0; col_A < A_n_cols; ++col_A)
{
const in_eT1* A_coldata = A.colptr(col_A);
out_eT acc = out_eT(0);
for(uword i=0; i < A_n_rows; ++i)
{
acc += upgrade_val<in_eT1,in_eT2>::apply(B_rowdata[i]) * upgrade_val<in_eT1,in_eT2>::apply(A_coldata[i]);
}
if( (use_alpha == false) && (use_beta == false) ) { C.at(col_A,row_B) = acc; }
else if( (use_alpha == true ) && (use_beta == false) ) { C.at(col_A,row_B) = alpha*acc; }
else if( (use_alpha == false) && (use_beta == true ) ) { C.at(col_A,row_B) = acc + beta*C.at(col_A,row_B); }
else if( (use_alpha == true ) && (use_beta == true ) ) { C.at(col_A,row_B) = alpha*acc + beta*C.at(col_A,row_B); }
}
}
}
}
};
//! \brief
//! Matrix multplication where the matrices have differing element types.
template<const bool do_trans_A=false, const bool do_trans_B=false, const bool use_alpha=false, const bool use_beta=false>
class gemm_mixed
{
public:
//! immediate multiplication of matrices A and B, storing the result in C
template<typename out_eT, typename in_eT1, typename in_eT2>
inline
static
void
apply
(
Mat<out_eT>& C,
const Mat<in_eT1>& A,
const Mat<in_eT2>& B,
const out_eT alpha = out_eT(1),
const out_eT beta = out_eT(0)
)
{
arma_extra_debug_sigprint();
if((is_cx<in_eT1>::yes && do_trans_A) || (is_cx<in_eT2>::yes && do_trans_B))
{
// better-than-nothing handling of hermitian transpose
Mat<in_eT1> tmp_A;
Mat<in_eT2> tmp_B;
const bool predo_trans_A = ( (do_trans_A == true) && (is_cx<in_eT1>::yes) );
const bool predo_trans_B = ( (do_trans_B == true) && (is_cx<in_eT2>::yes) );
if(predo_trans_A) { op_htrans::apply_mat_noalias(tmp_A, A); }
if(predo_trans_B) { op_htrans::apply_mat_noalias(tmp_B, B); }
const Mat<in_eT1>& AA = (predo_trans_A == false) ? A : tmp_A;
const Mat<in_eT2>& BB = (predo_trans_B == false) ? B : tmp_B;
gemm_mixed_large<((predo_trans_A) ? false : do_trans_A), ((predo_trans_B) ? false : do_trans_B), use_alpha, use_beta>::apply(C, AA, BB, alpha, beta);
}
else
{
gemm_mixed_large<do_trans_A, do_trans_B, use_alpha, use_beta>::apply(C, A, B, alpha, beta);
}
}
};
//! @}
|