view liboctave/floatLU.cc @ 14327:4d917a6a858b stable

doc: Use Octave coding conventions in @example blocks of docstrings. * accumarray.m, accumdim.m, bar.m, base2dec.m, bincoeff.m, bitcmp.m, bitset.m, celldisp.m, chop.m, clabel.m, cloglog.m, colon.m, compass.m, computer.m, contour3.m, contourc.m, corr.m, cstrcat.m, ctime.m, cylinder.m, date.m, dec2base.m, demo.m, dir.m, dlmwrite.m, expm.m, ezcontourf.m, ezcontour.m, ezmeshc.m, ezmesh.m, ezplot.m, ezsurfc.m, ezsurf.m, feather.m, findobj.m, flipdim.m, fplot.m, genvarname.m, getfield.m, hankel.m, hilb.m, hist.m, idivide.m, index.m, int2str.m, interp1.m, is_leap_year.m, ismember.m, isocolors.m, isonormals.m, isosurface.m, kurtosis.m, legendre.m, linkprop.m, logit.m, logm.m, __makeinfo__.m, __marching_cube__.m, median.m, mkoctfile.m, moment.m, mpoles.m, orderfields.m, pcg.m, pcr.m, plot3.m, plotmatrix.m, polyaffine.m, polygcd.m, poly.m, polyout.m, print.m, qp.m, quadgk.m, qzhess.m, randi.m, rat.m, refreshdata.m, residue.m, rose.m, rot90.m, saveas.m, saveobj.m, shiftdim.m, skewness.m, spaugment.m, spdiags.m, sqp.m, stem.m, str2num.m, strcat.m, strjust.m, strread.m, strsplit.m, structfun.m, subplot.m, subsindex.m, substruct.m, surfl.m, surfnorm.m, svds.m, uimenu.m, union.m, voronoi.m, warning_ids.m, wblpdf.m: Use Octave coding conventions in @example blocks of docstrings.
author Rik <octave@nomad.inbox5.com>
date Sat, 04 Feb 2012 22:12:50 -0800
parents 72c96de7a403
children
line wrap: on
line source

/*

Copyright (C) 1994-2012 John W. Eaton
Copyright (C) 2009 VZLU Prague, a.s.

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 "floatLU.h"
#include "f77-fcn.h"
#include "lo-error.h"
#include "oct-locbuf.h"
#include "fColVector.h"

// Instantiate the base LU class for the types we need.

#include <base-lu.h>
#include <base-lu.cc>

template class base_lu <FloatMatrix>;

// Define the constructor for this particular derivation.

extern "C"
{
  F77_RET_T
  F77_FUNC (sgetrf, SGETRF) (const octave_idx_type&, const octave_idx_type&,
                             float*, const octave_idx_type&, octave_idx_type*,
                             octave_idx_type&);

#ifdef HAVE_QRUPDATE_LUU
  F77_RET_T
  F77_FUNC (slu1up, SLU1UP) (const octave_idx_type&, const octave_idx_type&,
                             float *, const octave_idx_type&,
                             float *, const octave_idx_type&,
                             float *, float *);

  F77_RET_T
  F77_FUNC (slup1up, SLUP1UP) (const octave_idx_type&, const octave_idx_type&,
                               float *, const octave_idx_type&,
                               float *, const octave_idx_type&,
                               octave_idx_type *, const float *,
                               const float *, float *);
#endif
}

FloatLU::FloatLU (const FloatMatrix& a)
{
  octave_idx_type a_nr = a.rows ();
  octave_idx_type a_nc = a.cols ();
  octave_idx_type mn = (a_nr < a_nc ? a_nr : a_nc);

  ipvt.resize (dim_vector (mn, 1));
  octave_idx_type *pipvt = ipvt.fortran_vec ();

  a_fact = a;
  float *tmp_data = a_fact.fortran_vec ();

  octave_idx_type info = 0;

  F77_XFCN (sgetrf, SGETRF, (a_nr, a_nc, tmp_data, a_nr, pipvt, info));

  for (octave_idx_type i = 0; i < mn; i++)
    pipvt[i] -= 1;
}

#ifdef HAVE_QRUPDATE_LUU

void FloatLU::update (const FloatColumnVector& u, const FloatColumnVector& v)
{
  if (packed ())
    unpack ();

  FloatMatrix& l = l_fact;
  FloatMatrix& r = a_fact;

  octave_idx_type m = l.rows ();
  octave_idx_type n = r.columns ();
  octave_idx_type k = l.columns ();

  if (u.length () == m && v.length () == n)
    {
      FloatColumnVector utmp = u, vtmp = v;
      F77_XFCN (slu1up, SLU1UP, (m, n, l.fortran_vec (), m, r.fortran_vec (), k,
                                 utmp.fortran_vec (), vtmp.fortran_vec ()));
    }
  else
    (*current_liboctave_error_handler) ("luupdate: dimensions mismatch");
}

void FloatLU::update (const FloatMatrix& u, const FloatMatrix& v)
{
  if (packed ())
    unpack ();

  FloatMatrix& l = l_fact;
  FloatMatrix& r = a_fact;

  octave_idx_type m = l.rows ();
  octave_idx_type n = r.columns ();
  octave_idx_type k = l.columns ();

  if (u.rows () == m && v.rows () == n && u.cols () == v.cols ())
    {
      for (volatile octave_idx_type i = 0; i < u.cols (); i++)
        {
          FloatColumnVector utmp = u.column (i), vtmp = v.column (i);
          F77_XFCN (slu1up, SLU1UP, (m, n, l.fortran_vec (), m, r.fortran_vec (), k,
                                     utmp.fortran_vec (), vtmp.fortran_vec ()));
        }
    }
  else
    (*current_liboctave_error_handler) ("luupdate: dimensions mismatch");
}

void FloatLU::update_piv (const FloatColumnVector& u, const FloatColumnVector& v)
{
  if (packed ())
    unpack ();

  FloatMatrix& l = l_fact;
  FloatMatrix& r = a_fact;

  octave_idx_type m = l.rows ();
  octave_idx_type n = r.columns ();
  octave_idx_type k = l.columns ();

  if (u.length () == m && v.length () == n)
    {
      FloatColumnVector utmp = u, vtmp = v;
      OCTAVE_LOCAL_BUFFER (float, w, m);
      for (octave_idx_type i = 0; i < m; i++) ipvt(i) += 1; // increment
      F77_XFCN (slup1up, SLUP1UP, (m, n, l.fortran_vec (), m, r.fortran_vec (), k,
                                   ipvt.fortran_vec (), utmp.data (), vtmp.data (), w));
      for (octave_idx_type i = 0; i < m; i++) ipvt(i) -= 1; // decrement
    }
  else
    (*current_liboctave_error_handler) ("luupdate: dimensions mismatch");
}

void FloatLU::update_piv (const FloatMatrix& u, const FloatMatrix& v)
{
  if (packed ())
    unpack ();

  FloatMatrix& l = l_fact;
  FloatMatrix& r = a_fact;

  octave_idx_type m = l.rows ();
  octave_idx_type n = r.columns ();
  octave_idx_type k = l.columns ();

  if (u.rows () == m && v.rows () == n && u.cols () == v.cols ())
    {
      OCTAVE_LOCAL_BUFFER (float, w, m);
      for (octave_idx_type i = 0; i < m; i++) ipvt(i) += 1; // increment
      for (volatile octave_idx_type i = 0; i < u.cols (); i++)
        {
          FloatColumnVector utmp = u.column (i), vtmp = v.column (i);
          F77_XFCN (slup1up, SLUP1UP, (m, n, l.fortran_vec (), m, r.fortran_vec (), k,
                                       ipvt.fortran_vec (), utmp.data (), vtmp.data (), w));
        }
      for (octave_idx_type i = 0; i < m; i++) ipvt(i) -= 1; // decrement
    }
  else
    (*current_liboctave_error_handler) ("luupdate: dimensions mismatch");
}

#else

void FloatLU::update (const FloatColumnVector&, const FloatColumnVector&)
{
  (*current_liboctave_error_handler) ("luupdate: not available in this version");
}

void FloatLU::update (const FloatMatrix&, const FloatMatrix&)
{
  (*current_liboctave_error_handler) ("luupdate: not available in this version");
}

void FloatLU::update_piv (const FloatColumnVector&, const FloatColumnVector&)
{
  (*current_liboctave_error_handler) ("luupdate: not available in this version");
}

void FloatLU::update_piv (const FloatMatrix&, const FloatMatrix&)
{
  (*current_liboctave_error_handler) ("luupdate: not available in this version");
}

#endif