GNU Octave  3.8.0
A high-level interpreted language, primarily intended for numerical computations, mostly compatible with Matlab
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
zdotc3.f
Go to the documentation of this file.
1 c Copyright (C) 2009-2013 VZLU Prague, a.s., Czech Republic
2 c
3 c Author: Jaroslav Hajek <highegg@gmail.com>
4 c
5 c This file is part of Octave.
6 c
7 c Octave is free software; you can redistribute it and/or modify
8 c it under the terms of the GNU General Public License as published by
9 c the Free Software Foundation; either version 3 of the License, or
10 c (at your option) any later version.
11 c
12 c This program is distributed in the hope that it will be useful,
13 c but WITHOUT ANY WARRANTY; without even the implied warranty of
14 c MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15 c GNU General Public License for more details.
16 c
17 c You should have received a copy of the GNU General Public License
18 c along with this software; see the file COPYING. If not, see
19 c <http://www.gnu.org/licenses/>.
20 c
21  subroutine zdotc3(m,n,k,a,b,c)
22 c purpose: a 3-dimensional dot product.
23 c c = sum (conj (a) .* b, 2), where a and b are 3d arrays.
24 c arguments:
25 c m,n,k (in) the dimensions of a and b
26 c a,b (in) double complex input arrays of size (m,k,n)
27 c c (out) double complex output array, size (m,n)
28  integer m,n,k,i,j,l
29  double complex a(m,k,n),b(m,k,n)
30  double complex c(m,n)
31 
32  double complex zdotc
33  external zdotc
34 
35 c quick return if possible.
36  if (m <= 0 .or. n <= 0) return
37 
38  if (m == 1) then
39 c the column-major case.
40  do j = 1,n
41  c(1,j) = zdotc(k,a(1,1,j),1,b(1,1,j),1)
42  end do
43  else
44 c We prefer performance here, because that's what we generally
45 c do by default in reduction functions. Besides, the accuracy
46 c of xDOT is questionable. Hence, do a cache-aligned nested loop.
47  do j = 1,n
48  do i = 1,m
49  c(i,j) = 0d0
50  end do
51  do l = 1,k
52  do i = 1,m
53  c(i,j) = c(i,j) + conjg(a(i,l,j))*b(i,l,j)
54  end do
55  end do
56  end do
57  end if
58 
59  end subroutine