annotate libgfortran/m4/matmul_internal.m4 @ 111:04ced10e8804

gcc 7
author kono
date Fri, 27 Oct 2017 22:46:09 +0900
parents
children 84e7813d76e9
Ignore whitespace changes - Everywhere: Within whitespace: At end of lines:
rev   line source
111
kono
parents:
diff changeset
1 `void
kono
parents:
diff changeset
2 'matmul_name` ('rtype` * const restrict retarray,
kono
parents:
diff changeset
3 'rtype` * const restrict a, 'rtype` * const restrict b, int try_blas,
kono
parents:
diff changeset
4 int blas_limit, blas_call gemm)
kono
parents:
diff changeset
5 {
kono
parents:
diff changeset
6 const 'rtype_name` * restrict abase;
kono
parents:
diff changeset
7 const 'rtype_name` * restrict bbase;
kono
parents:
diff changeset
8 'rtype_name` * restrict dest;
kono
parents:
diff changeset
9
kono
parents:
diff changeset
10 index_type rxstride, rystride, axstride, aystride, bxstride, bystride;
kono
parents:
diff changeset
11 index_type x, y, n, count, xcount, ycount;
kono
parents:
diff changeset
12
kono
parents:
diff changeset
13 assert (GFC_DESCRIPTOR_RANK (a) == 2
kono
parents:
diff changeset
14 || GFC_DESCRIPTOR_RANK (b) == 2);
kono
parents:
diff changeset
15
kono
parents:
diff changeset
16 /* C[xcount,ycount] = A[xcount, count] * B[count,ycount]
kono
parents:
diff changeset
17
kono
parents:
diff changeset
18 Either A or B (but not both) can be rank 1:
kono
parents:
diff changeset
19
kono
parents:
diff changeset
20 o One-dimensional argument A is implicitly treated as a row matrix
kono
parents:
diff changeset
21 dimensioned [1,count], so xcount=1.
kono
parents:
diff changeset
22
kono
parents:
diff changeset
23 o One-dimensional argument B is implicitly treated as a column matrix
kono
parents:
diff changeset
24 dimensioned [count, 1], so ycount=1.
kono
parents:
diff changeset
25 */
kono
parents:
diff changeset
26
kono
parents:
diff changeset
27 if (retarray->base_addr == NULL)
kono
parents:
diff changeset
28 {
kono
parents:
diff changeset
29 if (GFC_DESCRIPTOR_RANK (a) == 1)
kono
parents:
diff changeset
30 {
kono
parents:
diff changeset
31 GFC_DIMENSION_SET(retarray->dim[0], 0,
kono
parents:
diff changeset
32 GFC_DESCRIPTOR_EXTENT(b,1) - 1, 1);
kono
parents:
diff changeset
33 }
kono
parents:
diff changeset
34 else if (GFC_DESCRIPTOR_RANK (b) == 1)
kono
parents:
diff changeset
35 {
kono
parents:
diff changeset
36 GFC_DIMENSION_SET(retarray->dim[0], 0,
kono
parents:
diff changeset
37 GFC_DESCRIPTOR_EXTENT(a,0) - 1, 1);
kono
parents:
diff changeset
38 }
kono
parents:
diff changeset
39 else
kono
parents:
diff changeset
40 {
kono
parents:
diff changeset
41 GFC_DIMENSION_SET(retarray->dim[0], 0,
kono
parents:
diff changeset
42 GFC_DESCRIPTOR_EXTENT(a,0) - 1, 1);
kono
parents:
diff changeset
43
kono
parents:
diff changeset
44 GFC_DIMENSION_SET(retarray->dim[1], 0,
kono
parents:
diff changeset
45 GFC_DESCRIPTOR_EXTENT(b,1) - 1,
kono
parents:
diff changeset
46 GFC_DESCRIPTOR_EXTENT(retarray,0));
kono
parents:
diff changeset
47 }
kono
parents:
diff changeset
48
kono
parents:
diff changeset
49 retarray->base_addr
kono
parents:
diff changeset
50 = xmallocarray (size0 ((array_t *) retarray), sizeof ('rtype_name`));
kono
parents:
diff changeset
51 retarray->offset = 0;
kono
parents:
diff changeset
52 }
kono
parents:
diff changeset
53 else if (unlikely (compile_options.bounds_check))
kono
parents:
diff changeset
54 {
kono
parents:
diff changeset
55 index_type ret_extent, arg_extent;
kono
parents:
diff changeset
56
kono
parents:
diff changeset
57 if (GFC_DESCRIPTOR_RANK (a) == 1)
kono
parents:
diff changeset
58 {
kono
parents:
diff changeset
59 arg_extent = GFC_DESCRIPTOR_EXTENT(b,1);
kono
parents:
diff changeset
60 ret_extent = GFC_DESCRIPTOR_EXTENT(retarray,0);
kono
parents:
diff changeset
61 if (arg_extent != ret_extent)
kono
parents:
diff changeset
62 runtime_error ("Incorrect extent in return array in"
kono
parents:
diff changeset
63 " MATMUL intrinsic: is %ld, should be %ld",
kono
parents:
diff changeset
64 (long int) ret_extent, (long int) arg_extent);
kono
parents:
diff changeset
65 }
kono
parents:
diff changeset
66 else if (GFC_DESCRIPTOR_RANK (b) == 1)
kono
parents:
diff changeset
67 {
kono
parents:
diff changeset
68 arg_extent = GFC_DESCRIPTOR_EXTENT(a,0);
kono
parents:
diff changeset
69 ret_extent = GFC_DESCRIPTOR_EXTENT(retarray,0);
kono
parents:
diff changeset
70 if (arg_extent != ret_extent)
kono
parents:
diff changeset
71 runtime_error ("Incorrect extent in return array in"
kono
parents:
diff changeset
72 " MATMUL intrinsic: is %ld, should be %ld",
kono
parents:
diff changeset
73 (long int) ret_extent, (long int) arg_extent);
kono
parents:
diff changeset
74 }
kono
parents:
diff changeset
75 else
kono
parents:
diff changeset
76 {
kono
parents:
diff changeset
77 arg_extent = GFC_DESCRIPTOR_EXTENT(a,0);
kono
parents:
diff changeset
78 ret_extent = GFC_DESCRIPTOR_EXTENT(retarray,0);
kono
parents:
diff changeset
79 if (arg_extent != ret_extent)
kono
parents:
diff changeset
80 runtime_error ("Incorrect extent in return array in"
kono
parents:
diff changeset
81 " MATMUL intrinsic for dimension 1:"
kono
parents:
diff changeset
82 " is %ld, should be %ld",
kono
parents:
diff changeset
83 (long int) ret_extent, (long int) arg_extent);
kono
parents:
diff changeset
84
kono
parents:
diff changeset
85 arg_extent = GFC_DESCRIPTOR_EXTENT(b,1);
kono
parents:
diff changeset
86 ret_extent = GFC_DESCRIPTOR_EXTENT(retarray,1);
kono
parents:
diff changeset
87 if (arg_extent != ret_extent)
kono
parents:
diff changeset
88 runtime_error ("Incorrect extent in return array in"
kono
parents:
diff changeset
89 " MATMUL intrinsic for dimension 2:"
kono
parents:
diff changeset
90 " is %ld, should be %ld",
kono
parents:
diff changeset
91 (long int) ret_extent, (long int) arg_extent);
kono
parents:
diff changeset
92 }
kono
parents:
diff changeset
93 }
kono
parents:
diff changeset
94 '
kono
parents:
diff changeset
95 sinclude(`matmul_asm_'rtype_code`.m4')dnl
kono
parents:
diff changeset
96 `
kono
parents:
diff changeset
97 if (GFC_DESCRIPTOR_RANK (retarray) == 1)
kono
parents:
diff changeset
98 {
kono
parents:
diff changeset
99 /* One-dimensional result may be addressed in the code below
kono
parents:
diff changeset
100 either as a row or a column matrix. We want both cases to
kono
parents:
diff changeset
101 work. */
kono
parents:
diff changeset
102 rxstride = rystride = GFC_DESCRIPTOR_STRIDE(retarray,0);
kono
parents:
diff changeset
103 }
kono
parents:
diff changeset
104 else
kono
parents:
diff changeset
105 {
kono
parents:
diff changeset
106 rxstride = GFC_DESCRIPTOR_STRIDE(retarray,0);
kono
parents:
diff changeset
107 rystride = GFC_DESCRIPTOR_STRIDE(retarray,1);
kono
parents:
diff changeset
108 }
kono
parents:
diff changeset
109
kono
parents:
diff changeset
110
kono
parents:
diff changeset
111 if (GFC_DESCRIPTOR_RANK (a) == 1)
kono
parents:
diff changeset
112 {
kono
parents:
diff changeset
113 /* Treat it as a a row matrix A[1,count]. */
kono
parents:
diff changeset
114 axstride = GFC_DESCRIPTOR_STRIDE(a,0);
kono
parents:
diff changeset
115 aystride = 1;
kono
parents:
diff changeset
116
kono
parents:
diff changeset
117 xcount = 1;
kono
parents:
diff changeset
118 count = GFC_DESCRIPTOR_EXTENT(a,0);
kono
parents:
diff changeset
119 }
kono
parents:
diff changeset
120 else
kono
parents:
diff changeset
121 {
kono
parents:
diff changeset
122 axstride = GFC_DESCRIPTOR_STRIDE(a,0);
kono
parents:
diff changeset
123 aystride = GFC_DESCRIPTOR_STRIDE(a,1);
kono
parents:
diff changeset
124
kono
parents:
diff changeset
125 count = GFC_DESCRIPTOR_EXTENT(a,1);
kono
parents:
diff changeset
126 xcount = GFC_DESCRIPTOR_EXTENT(a,0);
kono
parents:
diff changeset
127 }
kono
parents:
diff changeset
128
kono
parents:
diff changeset
129 if (count != GFC_DESCRIPTOR_EXTENT(b,0))
kono
parents:
diff changeset
130 {
kono
parents:
diff changeset
131 if (count > 0 || GFC_DESCRIPTOR_EXTENT(b,0) > 0)
kono
parents:
diff changeset
132 runtime_error ("dimension of array B incorrect in MATMUL intrinsic");
kono
parents:
diff changeset
133 }
kono
parents:
diff changeset
134
kono
parents:
diff changeset
135 if (GFC_DESCRIPTOR_RANK (b) == 1)
kono
parents:
diff changeset
136 {
kono
parents:
diff changeset
137 /* Treat it as a column matrix B[count,1] */
kono
parents:
diff changeset
138 bxstride = GFC_DESCRIPTOR_STRIDE(b,0);
kono
parents:
diff changeset
139
kono
parents:
diff changeset
140 /* bystride should never be used for 1-dimensional b.
kono
parents:
diff changeset
141 The value is only used for calculation of the
kono
parents:
diff changeset
142 memory by the buffer. */
kono
parents:
diff changeset
143 bystride = 256;
kono
parents:
diff changeset
144 ycount = 1;
kono
parents:
diff changeset
145 }
kono
parents:
diff changeset
146 else
kono
parents:
diff changeset
147 {
kono
parents:
diff changeset
148 bxstride = GFC_DESCRIPTOR_STRIDE(b,0);
kono
parents:
diff changeset
149 bystride = GFC_DESCRIPTOR_STRIDE(b,1);
kono
parents:
diff changeset
150 ycount = GFC_DESCRIPTOR_EXTENT(b,1);
kono
parents:
diff changeset
151 }
kono
parents:
diff changeset
152
kono
parents:
diff changeset
153 abase = a->base_addr;
kono
parents:
diff changeset
154 bbase = b->base_addr;
kono
parents:
diff changeset
155 dest = retarray->base_addr;
kono
parents:
diff changeset
156
kono
parents:
diff changeset
157 /* Now that everything is set up, we perform the multiplication
kono
parents:
diff changeset
158 itself. */
kono
parents:
diff changeset
159
kono
parents:
diff changeset
160 #define POW3(x) (((float) (x)) * ((float) (x)) * ((float) (x)))
kono
parents:
diff changeset
161 #define min(a,b) ((a) <= (b) ? (a) : (b))
kono
parents:
diff changeset
162 #define max(a,b) ((a) >= (b) ? (a) : (b))
kono
parents:
diff changeset
163
kono
parents:
diff changeset
164 if (try_blas && rxstride == 1 && (axstride == 1 || aystride == 1)
kono
parents:
diff changeset
165 && (bxstride == 1 || bystride == 1)
kono
parents:
diff changeset
166 && (((float) xcount) * ((float) ycount) * ((float) count)
kono
parents:
diff changeset
167 > POW3(blas_limit)))
kono
parents:
diff changeset
168 {
kono
parents:
diff changeset
169 const int m = xcount, n = ycount, k = count, ldc = rystride;
kono
parents:
diff changeset
170 const 'rtype_name` one = 1, zero = 0;
kono
parents:
diff changeset
171 const int lda = (axstride == 1) ? aystride : axstride,
kono
parents:
diff changeset
172 ldb = (bxstride == 1) ? bystride : bxstride;
kono
parents:
diff changeset
173
kono
parents:
diff changeset
174 if (lda > 0 && ldb > 0 && ldc > 0 && m > 1 && n > 1 && k > 1)
kono
parents:
diff changeset
175 {
kono
parents:
diff changeset
176 assert (gemm != NULL);
kono
parents:
diff changeset
177 gemm (axstride == 1 ? "N" : "T", bxstride == 1 ? "N" : "T", &m,
kono
parents:
diff changeset
178 &n, &k, &one, abase, &lda, bbase, &ldb, &zero, dest,
kono
parents:
diff changeset
179 &ldc, 1, 1);
kono
parents:
diff changeset
180 return;
kono
parents:
diff changeset
181 }
kono
parents:
diff changeset
182 }
kono
parents:
diff changeset
183
kono
parents:
diff changeset
184 if (rxstride == 1 && axstride == 1 && bxstride == 1)
kono
parents:
diff changeset
185 {
kono
parents:
diff changeset
186 /* This block of code implements a tuned matmul, derived from
kono
parents:
diff changeset
187 Superscalar GEMM-based level 3 BLAS, Beta version 0.1
kono
parents:
diff changeset
188
kono
parents:
diff changeset
189 Bo Kagstrom and Per Ling
kono
parents:
diff changeset
190 Department of Computing Science
kono
parents:
diff changeset
191 Umea University
kono
parents:
diff changeset
192 S-901 87 Umea, Sweden
kono
parents:
diff changeset
193
kono
parents:
diff changeset
194 from netlib.org, translated to C, and modified for matmul.m4. */
kono
parents:
diff changeset
195
kono
parents:
diff changeset
196 const 'rtype_name` *a, *b;
kono
parents:
diff changeset
197 'rtype_name` *c;
kono
parents:
diff changeset
198 const index_type m = xcount, n = ycount, k = count;
kono
parents:
diff changeset
199
kono
parents:
diff changeset
200 /* System generated locals */
kono
parents:
diff changeset
201 index_type a_dim1, a_offset, b_dim1, b_offset, c_dim1, c_offset,
kono
parents:
diff changeset
202 i1, i2, i3, i4, i5, i6;
kono
parents:
diff changeset
203
kono
parents:
diff changeset
204 /* Local variables */
kono
parents:
diff changeset
205 'rtype_name` f11, f12, f21, f22, f31, f32, f41, f42,
kono
parents:
diff changeset
206 f13, f14, f23, f24, f33, f34, f43, f44;
kono
parents:
diff changeset
207 index_type i, j, l, ii, jj, ll;
kono
parents:
diff changeset
208 index_type isec, jsec, lsec, uisec, ujsec, ulsec;
kono
parents:
diff changeset
209 'rtype_name` *t1;
kono
parents:
diff changeset
210
kono
parents:
diff changeset
211 a = abase;
kono
parents:
diff changeset
212 b = bbase;
kono
parents:
diff changeset
213 c = retarray->base_addr;
kono
parents:
diff changeset
214
kono
parents:
diff changeset
215 /* Parameter adjustments */
kono
parents:
diff changeset
216 c_dim1 = rystride;
kono
parents:
diff changeset
217 c_offset = 1 + c_dim1;
kono
parents:
diff changeset
218 c -= c_offset;
kono
parents:
diff changeset
219 a_dim1 = aystride;
kono
parents:
diff changeset
220 a_offset = 1 + a_dim1;
kono
parents:
diff changeset
221 a -= a_offset;
kono
parents:
diff changeset
222 b_dim1 = bystride;
kono
parents:
diff changeset
223 b_offset = 1 + b_dim1;
kono
parents:
diff changeset
224 b -= b_offset;
kono
parents:
diff changeset
225
kono
parents:
diff changeset
226 /* Empty c first. */
kono
parents:
diff changeset
227 for (j=1; j<=n; j++)
kono
parents:
diff changeset
228 for (i=1; i<=m; i++)
kono
parents:
diff changeset
229 c[i + j * c_dim1] = ('rtype_name`)0;
kono
parents:
diff changeset
230
kono
parents:
diff changeset
231 /* Early exit if possible */
kono
parents:
diff changeset
232 if (m == 0 || n == 0 || k == 0)
kono
parents:
diff changeset
233 return;
kono
parents:
diff changeset
234
kono
parents:
diff changeset
235 /* Adjust size of t1 to what is needed. */
kono
parents:
diff changeset
236 index_type t1_dim;
kono
parents:
diff changeset
237 t1_dim = (a_dim1-1) * 256 + b_dim1;
kono
parents:
diff changeset
238 if (t1_dim > 65536)
kono
parents:
diff changeset
239 t1_dim = 65536;
kono
parents:
diff changeset
240
kono
parents:
diff changeset
241 t1 = malloc (t1_dim * sizeof('rtype_name`));
kono
parents:
diff changeset
242
kono
parents:
diff changeset
243 /* Start turning the crank. */
kono
parents:
diff changeset
244 i1 = n;
kono
parents:
diff changeset
245 for (jj = 1; jj <= i1; jj += 512)
kono
parents:
diff changeset
246 {
kono
parents:
diff changeset
247 /* Computing MIN */
kono
parents:
diff changeset
248 i2 = 512;
kono
parents:
diff changeset
249 i3 = n - jj + 1;
kono
parents:
diff changeset
250 jsec = min(i2,i3);
kono
parents:
diff changeset
251 ujsec = jsec - jsec % 4;
kono
parents:
diff changeset
252 i2 = k;
kono
parents:
diff changeset
253 for (ll = 1; ll <= i2; ll += 256)
kono
parents:
diff changeset
254 {
kono
parents:
diff changeset
255 /* Computing MIN */
kono
parents:
diff changeset
256 i3 = 256;
kono
parents:
diff changeset
257 i4 = k - ll + 1;
kono
parents:
diff changeset
258 lsec = min(i3,i4);
kono
parents:
diff changeset
259 ulsec = lsec - lsec % 2;
kono
parents:
diff changeset
260
kono
parents:
diff changeset
261 i3 = m;
kono
parents:
diff changeset
262 for (ii = 1; ii <= i3; ii += 256)
kono
parents:
diff changeset
263 {
kono
parents:
diff changeset
264 /* Computing MIN */
kono
parents:
diff changeset
265 i4 = 256;
kono
parents:
diff changeset
266 i5 = m - ii + 1;
kono
parents:
diff changeset
267 isec = min(i4,i5);
kono
parents:
diff changeset
268 uisec = isec - isec % 2;
kono
parents:
diff changeset
269 i4 = ll + ulsec - 1;
kono
parents:
diff changeset
270 for (l = ll; l <= i4; l += 2)
kono
parents:
diff changeset
271 {
kono
parents:
diff changeset
272 i5 = ii + uisec - 1;
kono
parents:
diff changeset
273 for (i = ii; i <= i5; i += 2)
kono
parents:
diff changeset
274 {
kono
parents:
diff changeset
275 t1[l - ll + 1 + ((i - ii + 1) << 8) - 257] =
kono
parents:
diff changeset
276 a[i + l * a_dim1];
kono
parents:
diff changeset
277 t1[l - ll + 2 + ((i - ii + 1) << 8) - 257] =
kono
parents:
diff changeset
278 a[i + (l + 1) * a_dim1];
kono
parents:
diff changeset
279 t1[l - ll + 1 + ((i - ii + 2) << 8) - 257] =
kono
parents:
diff changeset
280 a[i + 1 + l * a_dim1];
kono
parents:
diff changeset
281 t1[l - ll + 2 + ((i - ii + 2) << 8) - 257] =
kono
parents:
diff changeset
282 a[i + 1 + (l + 1) * a_dim1];
kono
parents:
diff changeset
283 }
kono
parents:
diff changeset
284 if (uisec < isec)
kono
parents:
diff changeset
285 {
kono
parents:
diff changeset
286 t1[l - ll + 1 + (isec << 8) - 257] =
kono
parents:
diff changeset
287 a[ii + isec - 1 + l * a_dim1];
kono
parents:
diff changeset
288 t1[l - ll + 2 + (isec << 8) - 257] =
kono
parents:
diff changeset
289 a[ii + isec - 1 + (l + 1) * a_dim1];
kono
parents:
diff changeset
290 }
kono
parents:
diff changeset
291 }
kono
parents:
diff changeset
292 if (ulsec < lsec)
kono
parents:
diff changeset
293 {
kono
parents:
diff changeset
294 i4 = ii + isec - 1;
kono
parents:
diff changeset
295 for (i = ii; i<= i4; ++i)
kono
parents:
diff changeset
296 {
kono
parents:
diff changeset
297 t1[lsec + ((i - ii + 1) << 8) - 257] =
kono
parents:
diff changeset
298 a[i + (ll + lsec - 1) * a_dim1];
kono
parents:
diff changeset
299 }
kono
parents:
diff changeset
300 }
kono
parents:
diff changeset
301
kono
parents:
diff changeset
302 uisec = isec - isec % 4;
kono
parents:
diff changeset
303 i4 = jj + ujsec - 1;
kono
parents:
diff changeset
304 for (j = jj; j <= i4; j += 4)
kono
parents:
diff changeset
305 {
kono
parents:
diff changeset
306 i5 = ii + uisec - 1;
kono
parents:
diff changeset
307 for (i = ii; i <= i5; i += 4)
kono
parents:
diff changeset
308 {
kono
parents:
diff changeset
309 f11 = c[i + j * c_dim1];
kono
parents:
diff changeset
310 f21 = c[i + 1 + j * c_dim1];
kono
parents:
diff changeset
311 f12 = c[i + (j + 1) * c_dim1];
kono
parents:
diff changeset
312 f22 = c[i + 1 + (j + 1) * c_dim1];
kono
parents:
diff changeset
313 f13 = c[i + (j + 2) * c_dim1];
kono
parents:
diff changeset
314 f23 = c[i + 1 + (j + 2) * c_dim1];
kono
parents:
diff changeset
315 f14 = c[i + (j + 3) * c_dim1];
kono
parents:
diff changeset
316 f24 = c[i + 1 + (j + 3) * c_dim1];
kono
parents:
diff changeset
317 f31 = c[i + 2 + j * c_dim1];
kono
parents:
diff changeset
318 f41 = c[i + 3 + j * c_dim1];
kono
parents:
diff changeset
319 f32 = c[i + 2 + (j + 1) * c_dim1];
kono
parents:
diff changeset
320 f42 = c[i + 3 + (j + 1) * c_dim1];
kono
parents:
diff changeset
321 f33 = c[i + 2 + (j + 2) * c_dim1];
kono
parents:
diff changeset
322 f43 = c[i + 3 + (j + 2) * c_dim1];
kono
parents:
diff changeset
323 f34 = c[i + 2 + (j + 3) * c_dim1];
kono
parents:
diff changeset
324 f44 = c[i + 3 + (j + 3) * c_dim1];
kono
parents:
diff changeset
325 i6 = ll + lsec - 1;
kono
parents:
diff changeset
326 for (l = ll; l <= i6; ++l)
kono
parents:
diff changeset
327 {
kono
parents:
diff changeset
328 f11 += t1[l - ll + 1 + ((i - ii + 1) << 8) - 257]
kono
parents:
diff changeset
329 * b[l + j * b_dim1];
kono
parents:
diff changeset
330 f21 += t1[l - ll + 1 + ((i - ii + 2) << 8) - 257]
kono
parents:
diff changeset
331 * b[l + j * b_dim1];
kono
parents:
diff changeset
332 f12 += t1[l - ll + 1 + ((i - ii + 1) << 8) - 257]
kono
parents:
diff changeset
333 * b[l + (j + 1) * b_dim1];
kono
parents:
diff changeset
334 f22 += t1[l - ll + 1 + ((i - ii + 2) << 8) - 257]
kono
parents:
diff changeset
335 * b[l + (j + 1) * b_dim1];
kono
parents:
diff changeset
336 f13 += t1[l - ll + 1 + ((i - ii + 1) << 8) - 257]
kono
parents:
diff changeset
337 * b[l + (j + 2) * b_dim1];
kono
parents:
diff changeset
338 f23 += t1[l - ll + 1 + ((i - ii + 2) << 8) - 257]
kono
parents:
diff changeset
339 * b[l + (j + 2) * b_dim1];
kono
parents:
diff changeset
340 f14 += t1[l - ll + 1 + ((i - ii + 1) << 8) - 257]
kono
parents:
diff changeset
341 * b[l + (j + 3) * b_dim1];
kono
parents:
diff changeset
342 f24 += t1[l - ll + 1 + ((i - ii + 2) << 8) - 257]
kono
parents:
diff changeset
343 * b[l + (j + 3) * b_dim1];
kono
parents:
diff changeset
344 f31 += t1[l - ll + 1 + ((i - ii + 3) << 8) - 257]
kono
parents:
diff changeset
345 * b[l + j * b_dim1];
kono
parents:
diff changeset
346 f41 += t1[l - ll + 1 + ((i - ii + 4) << 8) - 257]
kono
parents:
diff changeset
347 * b[l + j * b_dim1];
kono
parents:
diff changeset
348 f32 += t1[l - ll + 1 + ((i - ii + 3) << 8) - 257]
kono
parents:
diff changeset
349 * b[l + (j + 1) * b_dim1];
kono
parents:
diff changeset
350 f42 += t1[l - ll + 1 + ((i - ii + 4) << 8) - 257]
kono
parents:
diff changeset
351 * b[l + (j + 1) * b_dim1];
kono
parents:
diff changeset
352 f33 += t1[l - ll + 1 + ((i - ii + 3) << 8) - 257]
kono
parents:
diff changeset
353 * b[l + (j + 2) * b_dim1];
kono
parents:
diff changeset
354 f43 += t1[l - ll + 1 + ((i - ii + 4) << 8) - 257]
kono
parents:
diff changeset
355 * b[l + (j + 2) * b_dim1];
kono
parents:
diff changeset
356 f34 += t1[l - ll + 1 + ((i - ii + 3) << 8) - 257]
kono
parents:
diff changeset
357 * b[l + (j + 3) * b_dim1];
kono
parents:
diff changeset
358 f44 += t1[l - ll + 1 + ((i - ii + 4) << 8) - 257]
kono
parents:
diff changeset
359 * b[l + (j + 3) * b_dim1];
kono
parents:
diff changeset
360 }
kono
parents:
diff changeset
361 c[i + j * c_dim1] = f11;
kono
parents:
diff changeset
362 c[i + 1 + j * c_dim1] = f21;
kono
parents:
diff changeset
363 c[i + (j + 1) * c_dim1] = f12;
kono
parents:
diff changeset
364 c[i + 1 + (j + 1) * c_dim1] = f22;
kono
parents:
diff changeset
365 c[i + (j + 2) * c_dim1] = f13;
kono
parents:
diff changeset
366 c[i + 1 + (j + 2) * c_dim1] = f23;
kono
parents:
diff changeset
367 c[i + (j + 3) * c_dim1] = f14;
kono
parents:
diff changeset
368 c[i + 1 + (j + 3) * c_dim1] = f24;
kono
parents:
diff changeset
369 c[i + 2 + j * c_dim1] = f31;
kono
parents:
diff changeset
370 c[i + 3 + j * c_dim1] = f41;
kono
parents:
diff changeset
371 c[i + 2 + (j + 1) * c_dim1] = f32;
kono
parents:
diff changeset
372 c[i + 3 + (j + 1) * c_dim1] = f42;
kono
parents:
diff changeset
373 c[i + 2 + (j + 2) * c_dim1] = f33;
kono
parents:
diff changeset
374 c[i + 3 + (j + 2) * c_dim1] = f43;
kono
parents:
diff changeset
375 c[i + 2 + (j + 3) * c_dim1] = f34;
kono
parents:
diff changeset
376 c[i + 3 + (j + 3) * c_dim1] = f44;
kono
parents:
diff changeset
377 }
kono
parents:
diff changeset
378 if (uisec < isec)
kono
parents:
diff changeset
379 {
kono
parents:
diff changeset
380 i5 = ii + isec - 1;
kono
parents:
diff changeset
381 for (i = ii + uisec; i <= i5; ++i)
kono
parents:
diff changeset
382 {
kono
parents:
diff changeset
383 f11 = c[i + j * c_dim1];
kono
parents:
diff changeset
384 f12 = c[i + (j + 1) * c_dim1];
kono
parents:
diff changeset
385 f13 = c[i + (j + 2) * c_dim1];
kono
parents:
diff changeset
386 f14 = c[i + (j + 3) * c_dim1];
kono
parents:
diff changeset
387 i6 = ll + lsec - 1;
kono
parents:
diff changeset
388 for (l = ll; l <= i6; ++l)
kono
parents:
diff changeset
389 {
kono
parents:
diff changeset
390 f11 += t1[l - ll + 1 + ((i - ii + 1) << 8) -
kono
parents:
diff changeset
391 257] * b[l + j * b_dim1];
kono
parents:
diff changeset
392 f12 += t1[l - ll + 1 + ((i - ii + 1) << 8) -
kono
parents:
diff changeset
393 257] * b[l + (j + 1) * b_dim1];
kono
parents:
diff changeset
394 f13 += t1[l - ll + 1 + ((i - ii + 1) << 8) -
kono
parents:
diff changeset
395 257] * b[l + (j + 2) * b_dim1];
kono
parents:
diff changeset
396 f14 += t1[l - ll + 1 + ((i - ii + 1) << 8) -
kono
parents:
diff changeset
397 257] * b[l + (j + 3) * b_dim1];
kono
parents:
diff changeset
398 }
kono
parents:
diff changeset
399 c[i + j * c_dim1] = f11;
kono
parents:
diff changeset
400 c[i + (j + 1) * c_dim1] = f12;
kono
parents:
diff changeset
401 c[i + (j + 2) * c_dim1] = f13;
kono
parents:
diff changeset
402 c[i + (j + 3) * c_dim1] = f14;
kono
parents:
diff changeset
403 }
kono
parents:
diff changeset
404 }
kono
parents:
diff changeset
405 }
kono
parents:
diff changeset
406 if (ujsec < jsec)
kono
parents:
diff changeset
407 {
kono
parents:
diff changeset
408 i4 = jj + jsec - 1;
kono
parents:
diff changeset
409 for (j = jj + ujsec; j <= i4; ++j)
kono
parents:
diff changeset
410 {
kono
parents:
diff changeset
411 i5 = ii + uisec - 1;
kono
parents:
diff changeset
412 for (i = ii; i <= i5; i += 4)
kono
parents:
diff changeset
413 {
kono
parents:
diff changeset
414 f11 = c[i + j * c_dim1];
kono
parents:
diff changeset
415 f21 = c[i + 1 + j * c_dim1];
kono
parents:
diff changeset
416 f31 = c[i + 2 + j * c_dim1];
kono
parents:
diff changeset
417 f41 = c[i + 3 + j * c_dim1];
kono
parents:
diff changeset
418 i6 = ll + lsec - 1;
kono
parents:
diff changeset
419 for (l = ll; l <= i6; ++l)
kono
parents:
diff changeset
420 {
kono
parents:
diff changeset
421 f11 += t1[l - ll + 1 + ((i - ii + 1) << 8) -
kono
parents:
diff changeset
422 257] * b[l + j * b_dim1];
kono
parents:
diff changeset
423 f21 += t1[l - ll + 1 + ((i - ii + 2) << 8) -
kono
parents:
diff changeset
424 257] * b[l + j * b_dim1];
kono
parents:
diff changeset
425 f31 += t1[l - ll + 1 + ((i - ii + 3) << 8) -
kono
parents:
diff changeset
426 257] * b[l + j * b_dim1];
kono
parents:
diff changeset
427 f41 += t1[l - ll + 1 + ((i - ii + 4) << 8) -
kono
parents:
diff changeset
428 257] * b[l + j * b_dim1];
kono
parents:
diff changeset
429 }
kono
parents:
diff changeset
430 c[i + j * c_dim1] = f11;
kono
parents:
diff changeset
431 c[i + 1 + j * c_dim1] = f21;
kono
parents:
diff changeset
432 c[i + 2 + j * c_dim1] = f31;
kono
parents:
diff changeset
433 c[i + 3 + j * c_dim1] = f41;
kono
parents:
diff changeset
434 }
kono
parents:
diff changeset
435 i5 = ii + isec - 1;
kono
parents:
diff changeset
436 for (i = ii + uisec; i <= i5; ++i)
kono
parents:
diff changeset
437 {
kono
parents:
diff changeset
438 f11 = c[i + j * c_dim1];
kono
parents:
diff changeset
439 i6 = ll + lsec - 1;
kono
parents:
diff changeset
440 for (l = ll; l <= i6; ++l)
kono
parents:
diff changeset
441 {
kono
parents:
diff changeset
442 f11 += t1[l - ll + 1 + ((i - ii + 1) << 8) -
kono
parents:
diff changeset
443 257] * b[l + j * b_dim1];
kono
parents:
diff changeset
444 }
kono
parents:
diff changeset
445 c[i + j * c_dim1] = f11;
kono
parents:
diff changeset
446 }
kono
parents:
diff changeset
447 }
kono
parents:
diff changeset
448 }
kono
parents:
diff changeset
449 }
kono
parents:
diff changeset
450 }
kono
parents:
diff changeset
451 }
kono
parents:
diff changeset
452 free(t1);
kono
parents:
diff changeset
453 return;
kono
parents:
diff changeset
454 }
kono
parents:
diff changeset
455 else if (rxstride == 1 && aystride == 1 && bxstride == 1)
kono
parents:
diff changeset
456 {
kono
parents:
diff changeset
457 if (GFC_DESCRIPTOR_RANK (a) != 1)
kono
parents:
diff changeset
458 {
kono
parents:
diff changeset
459 const 'rtype_name` *restrict abase_x;
kono
parents:
diff changeset
460 const 'rtype_name` *restrict bbase_y;
kono
parents:
diff changeset
461 'rtype_name` *restrict dest_y;
kono
parents:
diff changeset
462 'rtype_name` s;
kono
parents:
diff changeset
463
kono
parents:
diff changeset
464 for (y = 0; y < ycount; y++)
kono
parents:
diff changeset
465 {
kono
parents:
diff changeset
466 bbase_y = &bbase[y*bystride];
kono
parents:
diff changeset
467 dest_y = &dest[y*rystride];
kono
parents:
diff changeset
468 for (x = 0; x < xcount; x++)
kono
parents:
diff changeset
469 {
kono
parents:
diff changeset
470 abase_x = &abase[x*axstride];
kono
parents:
diff changeset
471 s = ('rtype_name`) 0;
kono
parents:
diff changeset
472 for (n = 0; n < count; n++)
kono
parents:
diff changeset
473 s += abase_x[n] * bbase_y[n];
kono
parents:
diff changeset
474 dest_y[x] = s;
kono
parents:
diff changeset
475 }
kono
parents:
diff changeset
476 }
kono
parents:
diff changeset
477 }
kono
parents:
diff changeset
478 else
kono
parents:
diff changeset
479 {
kono
parents:
diff changeset
480 const 'rtype_name` *restrict bbase_y;
kono
parents:
diff changeset
481 'rtype_name` s;
kono
parents:
diff changeset
482
kono
parents:
diff changeset
483 for (y = 0; y < ycount; y++)
kono
parents:
diff changeset
484 {
kono
parents:
diff changeset
485 bbase_y = &bbase[y*bystride];
kono
parents:
diff changeset
486 s = ('rtype_name`) 0;
kono
parents:
diff changeset
487 for (n = 0; n < count; n++)
kono
parents:
diff changeset
488 s += abase[n*axstride] * bbase_y[n];
kono
parents:
diff changeset
489 dest[y*rystride] = s;
kono
parents:
diff changeset
490 }
kono
parents:
diff changeset
491 }
kono
parents:
diff changeset
492 }
kono
parents:
diff changeset
493 else if (axstride < aystride)
kono
parents:
diff changeset
494 {
kono
parents:
diff changeset
495 for (y = 0; y < ycount; y++)
kono
parents:
diff changeset
496 for (x = 0; x < xcount; x++)
kono
parents:
diff changeset
497 dest[x*rxstride + y*rystride] = ('rtype_name`)0;
kono
parents:
diff changeset
498
kono
parents:
diff changeset
499 for (y = 0; y < ycount; y++)
kono
parents:
diff changeset
500 for (n = 0; n < count; n++)
kono
parents:
diff changeset
501 for (x = 0; x < xcount; x++)
kono
parents:
diff changeset
502 /* dest[x,y] += a[x,n] * b[n,y] */
kono
parents:
diff changeset
503 dest[x*rxstride + y*rystride] +=
kono
parents:
diff changeset
504 abase[x*axstride + n*aystride] *
kono
parents:
diff changeset
505 bbase[n*bxstride + y*bystride];
kono
parents:
diff changeset
506 }
kono
parents:
diff changeset
507 else if (GFC_DESCRIPTOR_RANK (a) == 1)
kono
parents:
diff changeset
508 {
kono
parents:
diff changeset
509 const 'rtype_name` *restrict bbase_y;
kono
parents:
diff changeset
510 'rtype_name` s;
kono
parents:
diff changeset
511
kono
parents:
diff changeset
512 for (y = 0; y < ycount; y++)
kono
parents:
diff changeset
513 {
kono
parents:
diff changeset
514 bbase_y = &bbase[y*bystride];
kono
parents:
diff changeset
515 s = ('rtype_name`) 0;
kono
parents:
diff changeset
516 for (n = 0; n < count; n++)
kono
parents:
diff changeset
517 s += abase[n*axstride] * bbase_y[n*bxstride];
kono
parents:
diff changeset
518 dest[y*rxstride] = s;
kono
parents:
diff changeset
519 }
kono
parents:
diff changeset
520 }
kono
parents:
diff changeset
521 else
kono
parents:
diff changeset
522 {
kono
parents:
diff changeset
523 const 'rtype_name` *restrict abase_x;
kono
parents:
diff changeset
524 const 'rtype_name` *restrict bbase_y;
kono
parents:
diff changeset
525 'rtype_name` *restrict dest_y;
kono
parents:
diff changeset
526 'rtype_name` s;
kono
parents:
diff changeset
527
kono
parents:
diff changeset
528 for (y = 0; y < ycount; y++)
kono
parents:
diff changeset
529 {
kono
parents:
diff changeset
530 bbase_y = &bbase[y*bystride];
kono
parents:
diff changeset
531 dest_y = &dest[y*rystride];
kono
parents:
diff changeset
532 for (x = 0; x < xcount; x++)
kono
parents:
diff changeset
533 {
kono
parents:
diff changeset
534 abase_x = &abase[x*axstride];
kono
parents:
diff changeset
535 s = ('rtype_name`) 0;
kono
parents:
diff changeset
536 for (n = 0; n < count; n++)
kono
parents:
diff changeset
537 s += abase_x[n*aystride] * bbase_y[n*bxstride];
kono
parents:
diff changeset
538 dest_y[x*rxstride] = s;
kono
parents:
diff changeset
539 }
kono
parents:
diff changeset
540 }
kono
parents:
diff changeset
541 }
kono
parents:
diff changeset
542 }
kono
parents:
diff changeset
543 #undef POW3
kono
parents:
diff changeset
544 #undef min
kono
parents:
diff changeset
545 #undef max
kono
parents:
diff changeset
546 '