GNU Octave  9.1.0
A high-level interpreted language, primarily intended for numerical computations, mostly compatible with Matlab
oct-convn.cc
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////
2 //
3 // Copyright (C) 2010-2024 The Octave Project Developers
4 //
5 // See the file COPYRIGHT.md in the top-level directory of this
6 // distribution or <https://octave.org/copyright/>.
7 //
8 // This file is part of Octave.
9 //
10 // Octave is free software: you can redistribute it and/or modify it
11 // under the terms of the GNU General Public License as published by
12 // the Free Software Foundation, either version 3 of the License, or
13 // (at your option) any later version.
14 //
15 // Octave is distributed in the hope that it will be useful, but
16 // WITHOUT ANY WARRANTY; without even the implied warranty of
17 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
18 // GNU General Public License for more details.
19 //
20 // You should have received a copy of the GNU General Public License
21 // along with Octave; see the file COPYING. If not, see
22 // <https://www.gnu.org/licenses/>.
23 //
24 ////////////////////////////////////////////////////////////////////////
25 
26 #if defined (HAVE_CONFIG_H)
27 # include "config.h"
28 #endif
29 
30 #include <algorithm>
31 
32 #include "Array.h"
33 #include "CColVector.h"
34 #include "CMatrix.h"
35 #include "CNDArray.h"
36 #include "CRowVector.h"
37 #include "MArray.h"
38 #include "dColVector.h"
39 #include "dMatrix.h"
40 #include "dNDArray.h"
41 #include "dRowVector.h"
42 #include "f77-fcn.h"
43 #include "fCColVector.h"
44 #include "fCMatrix.h"
45 #include "fCNDArray.h"
46 #include "fCRowVector.h"
47 #include "fColVector.h"
48 #include "fMatrix.h"
49 #include "fNDArray.h"
50 #include "fRowVector.h"
51 #include "oct-convn.h"
52 
54 
55 // 2d convolution with a matrix kernel.
56 template <typename T, typename R>
57 static void
58 convolve_2d (const T *a, F77_INT ma, F77_INT na,
59  const R *b, F77_INT mb, F77_INT nb,
60  T *c, bool inner);
61 
62 // Forward instances to our Fortran implementations.
63 #define FORWARD_IMPL(T_CXX, R_CXX, T, R, T_CAST, T_CONST_CAST, \
64  R_CONST_CAST, f, F) \
65  extern "C" \
66  F77_RET_T \
67  F77_FUNC (f##conv2o, F##CONV2O) (const F77_INT&, const F77_INT&, \
68  const T*, const F77_INT&, \
69  const F77_INT&, const R*, T *); \
70  \
71  extern "C" \
72  F77_RET_T \
73  F77_FUNC (f##conv2i, F##CONV2I) (const F77_INT&, const F77_INT&, \
74  const T*, const F77_INT&, \
75  const F77_INT&, const R*, T *); \
76  \
77  template <> void \
78  convolve_2d<T_CXX, R_CXX> (const T_CXX *a, F77_INT ma, F77_INT na, \
79  const R_CXX *b, F77_INT mb, F77_INT nb, \
80  T_CXX *c, bool inner) \
81  { \
82  if (inner) \
83  F77_XFCN (f##conv2i, F##CONV2I, (ma, na, T_CONST_CAST (a), \
84  mb, nb, R_CONST_CAST (b), \
85  T_CAST (c))); \
86  else \
87  F77_XFCN (f##conv2o, F##CONV2O, (ma, na, T_CONST_CAST (a), \
88  mb, nb, R_CONST_CAST (b), \
89  T_CAST (c))); \
90  }
91 
92 FORWARD_IMPL (double, double, F77_DBLE, F77_DBLE,,,, d, D)
93 FORWARD_IMPL (float, float, F77_REAL, F77_REAL,,,, s, S)
94 
95 FORWARD_IMPL (std::complex<double>, std::complex<double>,
98 FORWARD_IMPL (std::complex<float>, std::complex<float>,
101 
102 FORWARD_IMPL (std::complex<double>, double,
104  F77_CONST_DBLE_CMPLX_ARG,, zd, ZD)
105 FORWARD_IMPL (std::complex<float>, float, F77_CMPLX, F77_REAL, F77_CMPLX_ARG,
106  F77_CONST_CMPLX_ARG,, cs, CS)
107 
108 template <typename T, typename R>
109 void convolve_nd (const T *a, const dim_vector& ad, const dim_vector& acd,
110  const R *b, const dim_vector& bd, const dim_vector& bcd,
111  T *c, const dim_vector& ccd, int nd, bool inner)
112 {
113  if (nd == 2)
114  {
115  F77_INT ad0 = to_f77_int (ad(0));
116  F77_INT ad1 = to_f77_int (ad(1));
117 
118  F77_INT bd0 = to_f77_int (bd(0));
119  F77_INT bd1 = to_f77_int (bd(1));
120 
121  convolve_2d<T, R> (a, ad0, ad1, b, bd0, bd1, c, inner);
122  }
123  else
124  {
125  octave_idx_type ma = acd(nd-2);
126  octave_idx_type na = ad(nd-1);
127  octave_idx_type mb = bcd(nd-2);
128  octave_idx_type nb = bd(nd-1);
129  octave_idx_type ldc = ccd(nd-2);
130 
131  if (inner)
132  {
133  for (octave_idx_type ja = 0; ja < na - nb + 1; ja++)
134  for (octave_idx_type jb = 0; jb < nb; jb++)
135  convolve_nd<T, R> (a + ma*(ja+jb), ad, acd,
136  b + mb*(nb-jb-1), bd, bcd,
137  c + ldc*ja, ccd, nd-1, inner);
138  }
139  else
140  {
141  for (octave_idx_type ja = 0; ja < na; ja++)
142  for (octave_idx_type jb = 0; jb < nb; jb++)
143  convolve_nd<T, R> (a + ma*ja, ad, acd, b + mb*jb, bd, bcd,
144  c + ldc*(ja+jb), ccd, nd-1, inner);
145  }
146  }
147 }
148 
149 // Arbitrary convolutor.
150 // The 2nd array is assumed to be the smaller one.
151 template <typename T, typename R>
152 static MArray<T>
153 convolve (const MArray<T>& a, const MArray<R>& b,
154  convn_type ct)
155 {
156  if (a.isempty () || b.isempty ())
157  return MArray<T> ();
158 
159  int nd = std::max (a.ndims (), b.ndims ());
160  const dim_vector adims = a.dims ().redim (nd);
161  const dim_vector bdims = b.dims ().redim (nd);
162  dim_vector cdims = dim_vector::alloc (nd);
163 
164  for (int i = 0; i < nd; i++)
165  {
166  if (ct == convn_valid)
167  cdims(i) = std::max (adims(i) - bdims(i) + 1,
168  static_cast<octave_idx_type> (0));
169  else
170  cdims(i) = std::max (adims(i) + bdims(i) - 1,
171  static_cast<octave_idx_type> (0));
172  }
173 
174  MArray<T> c (cdims, T ());
175 
176  // "valid" shape can sometimes result in empty matrices which must avoid
177  // calling Fortran code which does not expect this (bug #52067)
178  if (c.isempty ())
179  return c;
180 
181  convolve_nd<T, R> (a.data (), adims, adims.cumulative (),
182  b.data (), bdims, bdims.cumulative (),
183  c.fortran_vec (), cdims.cumulative (),
184  nd, ct == convn_valid);
185 
186  if (ct == convn_same)
187  {
188  // Pick the relevant part.
189  Array<idx_vector> sidx (dim_vector (nd, 1));
190 
191  for (int i = 0; i < nd; i++)
192  sidx(i) = idx_vector::make_range (bdims(i)/2, 1, adims(i));
193  c = c.index (sidx);
194  }
195 
196  return c;
197 }
198 
199 #define CONV_DEFS(TPREF, RPREF) \
200  TPREF ## NDArray \
201  convn (const TPREF ## NDArray& a, const RPREF ## NDArray& b, \
202  convn_type ct) \
203  { \
204  return convolve (a, b, ct); \
205  } \
206  TPREF ## Matrix \
207  convn (const TPREF ## Matrix& a, const RPREF ## Matrix& b, \
208  convn_type ct) \
209  { \
210  return convolve (a, b, ct); \
211  } \
212  TPREF ## Matrix \
213  convn (const TPREF ## Matrix& a, const RPREF ## ColumnVector& c, \
214  const RPREF ## RowVector& r, convn_type ct) \
215  { \
216  return convolve (a, c * r, ct); \
217  }
218 
222 CONV_DEFS (Float, Float)
223 CONV_DEFS (FloatComplex, Float)
225 
226 OCTAVE_END_NAMESPACE(octave)
#define C(a, b)
Definition: Faddeeva.cc:259
charNDArray max(char d, const charNDArray &m)
Definition: chNDArray.cc:230
N Dimensional Array with copy-on-write semantics.
Definition: Array.h:130
int ndims() const
Size of the specified dimension.
Definition: Array.h:671
const T * data() const
Size of the specified dimension.
Definition: Array.h:663
bool isempty() const
Size of the specified dimension.
Definition: Array.h:651
const dim_vector & dims() const
Return a const-reference so that dims ()(i) works efficiently.
Definition: Array.h:503
Template for N-dimensional array classes with like-type math operators.
Definition: MArray.h:63
Vector representing the dimensions (size) of an Array.
Definition: dim-vector.h:94
static dim_vector alloc(int n)
Definition: dim-vector.h:202
dim_vector cumulative() const
Return cumulative dimensions.
Definition: dim-vector.h:488
dim_vector redim(int n) const
Force certain dimensionality, preserving numel ().
Definition: dim-vector.cc:226
static idx_vector make_range(octave_idx_type start, octave_idx_type step, octave_idx_type len)
Definition: idx-vector.h:452
OCTAVE_BEGIN_NAMESPACE(octave) static octave_value daspk_fcn
#define F77_CONST_CMPLX_ARG(x)
Definition: f77-fcn.h:313
#define F77_DBLE_CMPLX_ARG(x)
Definition: f77-fcn.h:316
float F77_REAL
Definition: f77-fcn.h:303
#define F77_CMPLX_ARG(x)
Definition: f77-fcn.h:310
double F77_DBLE
Definition: f77-fcn.h:302
double _Complex F77_DBLE_CMPLX
Definition: f77-fcn.h:304
octave_f77_int_type F77_INT
Definition: f77-fcn.h:306
#define F77_CONST_DBLE_CMPLX_ARG(x)
Definition: f77-fcn.h:319
float _Complex F77_CMPLX
Definition: f77-fcn.h:305
F77_RET_T const F77_INT const F77_INT const F77_INT F77_DBLE const F77_INT F77_DBLE const F77_INT F77_DBLE const F77_INT F77_DBLE * Z
F77_RET_T const F77_DBLE const F77_DBLE F77_DBLE * d
std::complex< double > Complex
Definition: oct-cmplx.h:33
std::complex< float > FloatComplex
Definition: oct-cmplx.h:34
#define CONV_DEFS(TPREF, RPREF)
Definition: oct-convn.cc:197
#define FORWARD_IMPL(T_CXX, R_CXX, T, R, T_CAST, T_CONST_CAST, R_CONST_CAST, f, F)
Definition: oct-convn.cc:63
R void convolve_nd(const T *a, const dim_vector &ad, const dim_vector &acd, const R *b, const dim_vector &bd, const dim_vector &bcd, T *c, const dim_vector &ccd, int nd, bool inner)
Definition: oct-convn.cc:107
convn_type
Definition: oct-convn.h:52
@ convn_same
Definition: oct-convn.h:54
@ convn_valid
Definition: oct-convn.h:55