Mercurial > jwe > octave
view liboctave/array/dRowVector.cc @ 21202:f7121e111991
maint: indent #ifdef blocks in liboctave and src directories.
* Array-C.cc, Array-b.cc, Array-ch.cc, Array-d.cc, Array-f.cc, Array-fC.cc,
Array-i.cc, Array-idx-vec.cc, Array-s.cc, Array-str.cc, Array-util.cc,
Array-voidp.cc, Array.cc, CColVector.cc, CDiagMatrix.cc, CMatrix.cc,
CNDArray.cc, CRowVector.cc, CSparse.cc, CSparse.h, DiagArray2.cc, MArray-C.cc,
MArray-d.cc, MArray-f.cc, MArray-fC.cc, MArray-i.cc, MArray-s.cc, MArray.cc,
MDiagArray2.cc, MSparse-C.cc, MSparse-d.cc, MSparse.h, MatrixType.cc,
PermMatrix.cc, Range.cc, Sparse-C.cc, Sparse-b.cc, Sparse-d.cc, Sparse.cc,
boolMatrix.cc, boolNDArray.cc, boolSparse.cc, chMatrix.cc, chNDArray.cc,
dColVector.cc, dDiagMatrix.cc, dMatrix.cc, dNDArray.cc, dRowVector.cc,
dSparse.cc, dSparse.h, dim-vector.cc, fCColVector.cc, fCDiagMatrix.cc,
fCMatrix.cc, fCNDArray.cc, fCRowVector.cc, fColVector.cc, fDiagMatrix.cc,
fMatrix.cc, fNDArray.cc, fRowVector.cc, idx-vector.cc, int16NDArray.cc,
int32NDArray.cc, int64NDArray.cc, int8NDArray.cc, intNDArray.cc,
uint16NDArray.cc, uint32NDArray.cc, uint64NDArray.cc, uint8NDArray.cc,
blaswrap.c, cquit.c, f77-extern.cc, f77-fcn.c, f77-fcn.h, lo-error.c, quit.cc,
quit.h, CmplxAEPBAL.cc, CmplxCHOL.cc, CmplxGEPBAL.cc, CmplxHESS.cc, CmplxLU.cc,
CmplxQR.cc, CmplxQRP.cc, CmplxSCHUR.cc, CmplxSVD.cc, CollocWt.cc, DASPK.cc,
DASRT.cc, DASSL.cc, EIG.cc, LSODE.cc, ODES.cc, Quad.cc, base-lu.cc, base-qr.cc,
dbleAEPBAL.cc, dbleCHOL.cc, dbleGEPBAL.cc, dbleHESS.cc, dbleLU.cc, dbleQR.cc,
dbleQRP.cc, dbleSCHUR.cc, dbleSVD.cc, eigs-base.cc, fCmplxAEPBAL.cc,
fCmplxCHOL.cc, fCmplxGEPBAL.cc, fCmplxHESS.cc, fCmplxLU.cc, fCmplxQR.cc,
fCmplxQRP.cc, fCmplxSCHUR.cc, fCmplxSVD.cc, fEIG.cc, floatAEPBAL.cc,
floatCHOL.cc, floatGEPBAL.cc, floatHESS.cc, floatLU.cc, floatQR.cc,
floatQRP.cc, floatSCHUR.cc, floatSVD.cc, lo-mappers.cc, lo-specfun.cc,
oct-convn.cc, oct-fftw.cc, oct-fftw.h, oct-norm.cc, oct-rand.cc,
oct-spparms.cc, randgamma.c, randmtzig.c, randpoisson.c, sparse-chol.cc,
sparse-dmsolve.cc, sparse-lu.cc, sparse-qr.cc, mx-defs.h, dir-ops.cc,
file-ops.cc, file-stat.cc, lo-sysdep.cc, mach-info.cc, oct-env.cc,
oct-group.cc, oct-openmp.h, oct-passwd.cc, oct-syscalls.cc, oct-time.cc,
oct-uname.cc, pathlen.h, sysdir.h, syswait.h, cmd-edit.cc, cmd-hist.cc,
data-conv.cc, f2c-main.c, glob-match.cc, lo-array-errwarn.cc,
lo-array-gripes.cc, lo-cutils.c, lo-cutils.h, lo-ieee.cc, lo-math.h,
lo-regexp.cc, lo-utils.cc, oct-base64.cc, oct-glob.cc, oct-inttypes.cc,
oct-inttypes.h, oct-locbuf.cc, oct-mutex.cc, oct-refcount.h, oct-rl-edit.c,
oct-rl-hist.c, oct-shlib.cc, oct-sort.cc, pathsearch.cc, singleton-cleanup.cc,
sparse-sort.cc, sparse-util.cc, statdefs.h, str-vec.cc, unwind-prot.cc,
url-transfer.cc, display-available.h, main-cli.cc, main-gui.cc, main.in.cc,
mkoctfile.in.cc, octave-config.in.cc, shared-fcns.h:
indent #ifdef blocks in liboctave and src directories.
author | Rik <rik@octave.org> |
---|---|
date | Sat, 06 Feb 2016 06:40:13 -0800 |
parents | 7cac4e7458f2 |
children | 40de9f8f23a6 |
line wrap: on
line source
// RowVector manipulations. /* Copyright (C) 1994-2015 John W. Eaton This file is part of Octave. Octave is free software; you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation; either version 3 of the License, or (at your option) any later version. Octave is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details. You should have received a copy of the GNU General Public License along with Octave; see the file COPYING. If not, see <http://www.gnu.org/licenses/>. */ #ifdef HAVE_CONFIG_H # include <config.h> #endif #include <iostream> #include "Array-util.h" #include "f77-fcn.h" #include "functor.h" #include "lo-error.h" #include "mx-base.h" #include "mx-inlines.cc" #include "oct-cmplx.h" // Fortran functions we call. extern "C" { F77_RET_T F77_FUNC (dgemv, DGEMV) (F77_CONST_CHAR_ARG_DECL, const octave_idx_type&, const octave_idx_type&, const double&, const double*, const octave_idx_type&, const double*, const octave_idx_type&, const double&, double*, const octave_idx_type& F77_CHAR_ARG_LEN_DECL); F77_RET_T F77_FUNC (xddot, XDDOT) (const octave_idx_type&, const double*, const octave_idx_type&, const double*, const octave_idx_type&, double&); } // Row Vector class. bool RowVector::operator == (const RowVector& a) const { octave_idx_type len = numel (); if (len != a.numel ()) return 0; return mx_inline_equal (len, data (), a.data ()); } bool RowVector::operator != (const RowVector& a) const { return !(*this == a); } RowVector& RowVector::insert (const RowVector& a, octave_idx_type c) { octave_idx_type a_len = a.numel (); if (c < 0 || c + a_len > numel ()) (*current_liboctave_error_handler) ("range error for insert"); if (a_len > 0) { make_unique (); for (octave_idx_type i = 0; i < a_len; i++) xelem (c+i) = a.elem (i); } return *this; } RowVector& RowVector::fill (double val) { octave_idx_type len = numel (); if (len > 0) { make_unique (); for (octave_idx_type i = 0; i < len; i++) xelem (i) = val; } return *this; } RowVector& RowVector::fill (double val, octave_idx_type c1, octave_idx_type c2) { octave_idx_type len = numel (); if (c1 < 0 || c2 < 0 || c1 >= len || c2 >= len) (*current_liboctave_error_handler) ("range error for fill"); if (c1 > c2) { std::swap (c1, c2); } if (c2 >= c1) { make_unique (); for (octave_idx_type i = c1; i <= c2; i++) xelem (i) = val; } return *this; } RowVector RowVector::append (const RowVector& a) const { octave_idx_type len = numel (); octave_idx_type nc_insert = len; RowVector retval (len + a.numel ()); retval.insert (*this, 0); retval.insert (a, nc_insert); return retval; } ColumnVector RowVector::transpose (void) const { return MArray<double>::transpose (); } RowVector real (const ComplexRowVector& a) { return do_mx_unary_op<double, Complex> (a, mx_inline_real); } RowVector imag (const ComplexRowVector& a) { return do_mx_unary_op<double, Complex> (a, mx_inline_imag); } RowVector RowVector::extract (octave_idx_type c1, octave_idx_type c2) const { if (c1 > c2) { std::swap (c1, c2); } octave_idx_type new_c = c2 - c1 + 1; RowVector result (new_c); for (octave_idx_type i = 0; i < new_c; i++) result.xelem (i) = elem (c1+i); return result; } RowVector RowVector::extract_n (octave_idx_type r1, octave_idx_type n) const { RowVector result (n); for (octave_idx_type i = 0; i < n; i++) result.xelem (i) = elem (r1+i); return result; } // row vector by matrix -> row vector RowVector operator * (const RowVector& v, const Matrix& a) { RowVector retval; octave_idx_type len = v.numel (); octave_idx_type a_nr = a.rows (); octave_idx_type a_nc = a.cols (); if (a_nr != len) err_nonconformant ("operator *", 1, len, a_nr, a_nc); if (len == 0) retval.resize (a_nc, 0.0); else { // Transpose A to form A'*x == (x'*A)' octave_idx_type ld = a_nr; retval.resize (a_nc); double *y = retval.fortran_vec (); F77_XFCN (dgemv, DGEMV, (F77_CONST_CHAR_ARG2 ("T", 1), a_nr, a_nc, 1.0, a.data (), ld, v.data (), 1, 0.0, y, 1 F77_CHAR_ARG_LEN (1))); } return retval; } // other operations double RowVector::min (void) const { octave_idx_type len = numel (); if (len == 0) return 0; double res = elem (0); for (octave_idx_type i = 1; i < len; i++) if (elem (i) < res) res = elem (i); return res; } double RowVector::max (void) const { octave_idx_type len = numel (); if (len == 0) return 0; double res = elem (0); for (octave_idx_type i = 1; i < len; i++) if (elem (i) > res) res = elem (i); return res; } std::ostream& operator << (std::ostream& os, const RowVector& a) { // int field_width = os.precision () + 7; for (octave_idx_type i = 0; i < a.numel (); i++) os << " " /* setw (field_width) */ << a.elem (i); return os; } std::istream& operator >> (std::istream& is, RowVector& a) { octave_idx_type len = a.numel (); if (len > 0) { double tmp; for (octave_idx_type i = 0; i < len; i++) { is >> tmp; if (is) a.elem (i) = tmp; else break; } } return is; } // other operations RowVector linspace (double x1, double x2, octave_idx_type n) { NoAlias<RowVector> retval; if (n < 1) return retval; else retval.clear (n); retval(0) = x1; double delta = (x2 - x1) / (n - 1); for (octave_idx_type i = 1; i < n-1; i++) retval(i) = x1 + i*delta; retval(n-1) = x2; return retval; } // row vector by column vector -> scalar double operator * (const RowVector& v, const ColumnVector& a) { double retval = 0.0; octave_idx_type len = v.numel (); octave_idx_type a_len = a.numel (); if (len != a_len) err_nonconformant ("operator *", len, a_len); if (len != 0) F77_FUNC (xddot, XDDOT) (len, v.data (), 1, a.data (), 1, retval); return retval; } Complex operator * (const RowVector& v, const ComplexColumnVector& a) { ComplexRowVector tmp (v); return tmp * a; }