// @(#)root/matrix:$Name: $:$Id: TVectorT.cxx,v 1.4 2006/01/26 16:31:01 brun Exp $
// Authors: Fons Rademakers, Eddy Offermann Nov 2003
/*************************************************************************
* Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. *
* All rights reserved. *
* *
* For the licensing terms see $ROOTSYS/LICENSE. *
* For the list of contributors see $ROOTSYS/README/CREDITS. *
*************************************************************************/
//////////////////////////////////////////////////////////////////////////
// //
// TVectorT //
// //
// Template class of Vectors in the linear algebra package //
// //
// Unless otherwise specified, vector indices always start with 0, //
// spanning up to the specified limit-1. //
// //
// For (n) vectors where n <= kSizeMax (5 currently) storage space is //
// available on the stack, thus avoiding expensive allocation/ //
// deallocation of heap space . However, this introduces of course //
// kSizeMax overhead for each vector object . If this is an issue //
// recompile with a new appropriate value (>=0) for kSizeMax //
// //
// Another way to assign and store vector data is through Use //
// see for instance stressLinear.cxx file . //
// //
// Note that Constructors/assignments exists for all different matrix //
// views //
// //
// For usage examples see $ROOTSYS/test/stressLinear.cxx //
// //
//////////////////////////////////////////////////////////////////////////
#include "TVectorT.h"
#ifndef R__ALPHA
templateClassImp(TVectorT)
#endif
//______________________________________________________________________________
template<class Element>
void TVectorT<Element>::Delete_m(Int_t size,Element *&m)
{
if (m) {
if (size > kSizeMax)
delete [] m;
m = 0;
}
}
//______________________________________________________________________________
template<class Element>
Element* TVectorT<Element>::New_m(Int_t size)
{
if (size == 0) return 0;
else {
if ( size <= kSizeMax )
return fDataStack;
else {
Element *heap = new Element[size];
return heap;
}
}
}
//______________________________________________________________________________
template<class Element>
Int_t TVectorT<Element>::Memcpy_m(Element *newp,const Element *oldp,Int_t copySize,
Int_t newSize,Int_t oldSize)
{
if (copySize == 0 || oldp == newp) return 0;
else {
if ( newSize <= kSizeMax && oldSize <= kSizeMax ) {
// both pointers are inside fDataStack, be careful with copy direction !
if (newp > oldp) {
for (Int_t i = copySize-1; i >= 0; i--)
newp[i] = oldp[i];
} else {
for (Int_t i = 0; i < copySize; i++)
newp[i] = oldp[i];
}
}
else {
memcpy(newp,oldp,copySize*sizeof(Element));
}
}
return 0;
}
//______________________________________________________________________________
template<class Element>
void TVectorT<Element>::Allocate(Int_t nrows,Int_t row_lwb,Int_t init)
{
// Allocate new vector. Arguments are number of rows and row
// lowerbound (0 default).
Assert(nrows >= 0);
MakeValid();
fNrows = nrows;
fRowLwb = row_lwb;
fIsOwner = kTRUE;
fElements = New_m(fNrows);
if (init)
memset(fElements,0,fNrows*sizeof(Element));
}
//______________________________________________________________________________
template<class Element>
TVectorT<Element>::TVectorT(Int_t n)
{
Allocate(n,0,1);
}
//______________________________________________________________________________
template<class Element>
TVectorT<Element>::TVectorT(Int_t lwb,Int_t upb)
{
Allocate(upb-lwb+1,lwb,1);
}
//______________________________________________________________________________
template<class Element>
TVectorT<Element>::TVectorT(Int_t n,const Element *elements)
{
Allocate(n,0);
SetElements(elements);
}
//______________________________________________________________________________
template<class Element>
TVectorT<Element>::TVectorT(Int_t lwb,Int_t upb,const Element *elements)
{
Allocate(upb-lwb+1,lwb);
SetElements(elements);
}
//______________________________________________________________________________
template<class Element>
TVectorT<Element>::TVectorT(const TVectorT &another) : TObject(another)
{
Assert(another.IsValid());
Allocate(another.GetUpb()-another.GetLwb()+1,another.GetLwb());
*this = another;
}
//______________________________________________________________________________
template<class Element>
TVectorT<Element>::TVectorT(const TMatrixTRow_const<Element> &mr) : TObject()
{
const TMatrixTBase<Element> *mt = mr.GetMatrix();
Assert(mt->IsValid());
Allocate(mt->GetColUpb()-mt->GetColLwb()+1,mt->GetColLwb());
*this = mr;
}
//______________________________________________________________________________
template<class Element>
TVectorT<Element>::TVectorT(const TMatrixTColumn_const<Element> &mc) : TObject()
{
const TMatrixTBase<Element> *mt = mc.GetMatrix();
Assert(mt->IsValid());
Allocate(mt->GetRowUpb()-mt->GetRowLwb()+1,mt->GetRowLwb());
*this = mc;
}
//______________________________________________________________________________
template<class Element>
TVectorT<Element>::TVectorT(const TMatrixTDiag_const<Element> &md) : TObject()
{
const TMatrixTBase<Element> *mt = md.GetMatrix();
Assert(mt->IsValid());
Allocate(TMath::Min(mt->GetNrows(),mt->GetNcols()));
*this = md;
}
//______________________________________________________________________________
template<class Element>
TVectorT<Element>::TVectorT(Int_t lwb,Int_t upb,Element va_(iv1), ...)
{
// Make a vector and assign initial values. Argument list should contain
// Element values to assign to vector elements. The list must be
// terminated by the string "END". Example:
// TVectorT foo(1,3,0.0,1.0,1.5,"END");
const Int_t no_rows = upb-lwb+1;
Assert(no_rows);
Allocate(no_rows,lwb);
va_list args;
va_start(args,va_(iv1)); // Init 'args' to the beginning of
// the variable length list of args
(*this)(lwb) = iv1;
for (Int_t i = lwb+1; i <= upb; i++)
(*this)(i) = (Element)va_arg(args,Double_t);
if (strcmp((char *)va_arg(args,char *),"END"))
Error("TVectorT(Int_t, Int_t, ...)", "argument list must be terminated by \"END\"");
va_end(args);
}
//______________________________________________________________________________
template<class Element>
TVectorT<Element> &TVectorT<Element>::ResizeTo(Int_t lwb,Int_t upb)
{
// Resize the vector to [lwb:upb] .
// New dynamic elemenst are created, the overlapping part of the old ones are
// copied to the new structures, then the old elements are deleleted.
Assert(IsValid());
if (!fIsOwner) {
Error("ResizeTo(lwb,upb)","Not owner of data array,cannot resize");
Invalidate();
return *this;
}
const Int_t new_nrows = upb-lwb+1;
if (fNrows > 0) {
if (fNrows == new_nrows && fRowLwb == lwb)
return *this;
else if (new_nrows == 0) {
Clear();
return *this;
}
Element *elements_old = GetMatrixArray();
const Int_t nrows_old = fNrows;
const Int_t rowLwb_old = fRowLwb;
Allocate(new_nrows,lwb);
Assert(IsValid());
if (fNrows > kSizeMax || nrows_old > kSizeMax)
memset(GetMatrixArray(),0,fNrows*sizeof(Element));
else if (fNrows > nrows_old)
memset(GetMatrixArray()+nrows_old,0,(fNrows-nrows_old)*sizeof(Element));
// Copy overlap
const Int_t rowLwb_copy = TMath::Max(fRowLwb,rowLwb_old);
const Int_t rowUpb_copy = TMath::Min(fRowLwb+fNrows-1,rowLwb_old+nrows_old-1);
const Int_t nrows_copy = rowUpb_copy-rowLwb_copy+1;
const Int_t nelems_new = fNrows;
Element *elements_new = GetMatrixArray();
if (nrows_copy > 0) {
const Int_t rowOldOff = rowLwb_copy-rowLwb_old;
const Int_t rowNewOff = rowLwb_copy-fRowLwb;
Memcpy_m(elements_new+rowNewOff,elements_old+rowOldOff,nrows_copy,nelems_new,nrows_old);
}
Delete_m(nrows_old,elements_old);
} else {
Allocate(upb-lwb+1,lwb,1);
}
return *this;
}
//______________________________________________________________________________
template<class Element>
TVectorT<Element> &TVectorT<Element>::Use(Int_t lwb,Int_t upb,Element *data)
{
Assert(upb >= lwb);
Clear();
fNrows = upb-lwb+1;
fRowLwb = lwb;
fElements = data;
fIsOwner = kFALSE;
return *this;
}
//______________________________________________________________________________
template<class Element>
TVectorT<Element> &TVectorT<Element>::GetSub(Int_t row_lwb,Int_t row_upb,TVectorT<Element> &target,Option_t *option) const
{
// Get subvector [row_lwb..row_upb]; The indexing range of the
// returned vector depends on the argument option:
//
// option == "S" : return [0..row_upb-row_lwb+1] (default)
// else : return [row_lwb..row_upb]
Assert(IsValid());
if (row_lwb < fRowLwb || row_lwb > fRowLwb+fNrows-1) {
Error("GetSub","row_lwb out of bounds");
target.Invalidate();
return target;
}
if (row_upb < fRowLwb || row_upb > fRowLwb+fNrows-1) {
Error("GetSub","row_upb out of bounds");
target.Invalidate();
return target;
}
if (row_upb < row_lwb) {
Error("GetSub","row_upb < row_lwb");
target.Invalidate();
return target;
}
TString opt(option);
opt.ToUpper();
const Int_t shift = (opt.Contains("S")) ? 1 : 0;
Int_t row_lwb_sub;
Int_t row_upb_sub;
if (shift) {
row_lwb_sub = 0;
row_upb_sub = row_upb-row_lwb;
} else {
row_lwb_sub = row_lwb;
row_upb_sub = row_upb;
}
target.ResizeTo(row_lwb_sub,row_upb_sub);
const Int_t nrows_sub = row_upb_sub-row_lwb_sub+1;
const Element *ap = this->GetMatrixArray()+(row_lwb-fRowLwb);
Element *bp = target.GetMatrixArray();
for (Int_t irow = 0; irow < nrows_sub; irow++)
*bp++ = *ap++;
return target;
}
//______________________________________________________________________________
template<class Element>
TVectorT<Element> &TVectorT<Element>::SetSub(Int_t row_lwb,const TVectorT<Element> &source)
{
// Insert vector source starting at [row_lwb], thereby overwriting the part
// [row_lwb..row_lwb+nrows_source];
Assert(IsValid());
Assert(source.IsValid());
if (row_lwb < fRowLwb && row_lwb > fRowLwb+fNrows-1) {
Error("SetSub","row_lwb outof bounds");
Invalidate();
return *this;
}
const Int_t nRows_source = source.GetNrows();
if (row_lwb+nRows_source > fRowLwb+fNrows) {
Error("SetSub","source vector too large");
Invalidate();
return *this;
}
const Element *bp = source.GetMatrixArray();
Element *ap = this->GetMatrixArray()+(row_lwb-fRowLwb);
for (Int_t irow = 0; irow < nRows_source; irow++)
*ap++ = *bp++;
return *this;
}
//______________________________________________________________________________
template<class Element>
TVectorT<Element> &TVectorT<Element>::Zero()
{
Assert(IsValid());
memset(this->GetMatrixArray(),0,fNrows*sizeof(Element));
return *this;
}
//______________________________________________________________________________
template<class Element>
TVectorT<Element> &TVectorT<Element>::Abs()
{
// Take an absolute value of a vector, i.e. apply Abs() to each element.
Assert(IsValid());
Element *ep = this->GetMatrixArray();
const Element * const fp = ep+fNrows;
while (ep < fp) {
*ep = TMath::Abs(*ep);
ep++;
}
return *this;
}
//______________________________________________________________________________
template<class Element>
TVectorT<Element> &TVectorT<Element>::Sqr()
{
// Square each element of the vector.
Assert(IsValid());
Element *ep = this->GetMatrixArray();
const Element * const fp = ep+fNrows;
while (ep < fp) {
*ep = (*ep) * (*ep);
ep++;
}
return *this;
}
//______________________________________________________________________________
template<class Element>
TVectorT<Element> &TVectorT<Element>::Sqrt()
{
// Take square root of all elements.
Assert(IsValid());
Element *ep = this->GetMatrixArray();
const Element * const fp = ep+fNrows;
while (ep < fp) {
Assert(*ep >= 0);
*ep = TMath::Sqrt(*ep);
ep++;
}
return *this;
}
//______________________________________________________________________________
template<class Element>
TVectorT<Element> &TVectorT<Element>::Invert()
{
// v[i] = 1/v[i]
Assert(IsValid());
Element *ep = this->GetMatrixArray();
const Element * const fp = ep+fNrows;
while (ep < fp) {
Assert(*ep != 0.0);
*ep = 1./ *ep;
ep++;
}
return *this;
}
//______________________________________________________________________________
template<class Element>
TVectorT<Element> &TVectorT<Element>::SelectNonZeros(const TVectorT<Element> &select)
{
if (!AreCompatible(*this,select)) {
Error("SelectNonZeros(const TVectorT<Element> &","vector's not compatible");
Invalidate();
return *this;
}
const Element *sp = select.GetMatrixArray();
Element *ep = this->GetMatrixArray();
const Element * const fp = ep+fNrows;
while (ep < fp) {
if (*sp == 0.0)
*ep = 0.0;
sp++; ep++;
}
return *this;
}
//______________________________________________________________________________
template<class Element>
Element TVectorT<Element>::Norm1() const
{
// Compute the 1-norm of the vector SUM{ |v[i]| }.
Assert(IsValid());
Element norm = 0;
const Element *ep = this->GetMatrixArray();
const Element * const fp = ep+fNrows;
while (ep < fp)
norm += TMath::Abs(*ep++);
return norm;
}
//______________________________________________________________________________
template<class Element>
Element TVectorT<Element>::Norm2Sqr() const
{
// Compute the square of the 2-norm SUM{ v[i]^2 }.
Assert(IsValid());
Element norm = 0;
const Element *ep = this->GetMatrixArray();
const Element * const fp = ep+fNrows;
while (ep < fp) {
norm += (*ep) * (*ep);
ep++;
}
return norm;
}
//______________________________________________________________________________
template<class Element>
Element TVectorT<Element>::NormInf() const
{
// Compute the infinity-norm of the vector MAX{ |v[i]| }.
Assert(IsValid());
Element norm = 0;
const Element *ep = this->GetMatrixArray();
const Element * const fp = ep+fNrows;
while (ep < fp)
norm = TMath::Max(norm,TMath::Abs(*ep++));
return norm;
}
//______________________________________________________________________________
template<class Element>
Int_t TVectorT<Element>::NonZeros() const
{
// Compute the number of elements != 0.0
Assert(IsValid());
Int_t nr_nonzeros = 0;
const Element *ep = this->GetMatrixArray();
const Element * const fp = ep+fNrows;
while (ep < fp)
if (*ep++) nr_nonzeros++;
return nr_nonzeros;
}
//______________________________________________________________________________
template<class Element>
Element TVectorT<Element>::Sum() const
{
// Compute sum of elements
Assert(IsValid());
Element sum = 0.0;
const Element *ep = this->GetMatrixArray();
const Element * const fp = ep+fNrows;
while (ep < fp)
sum += *ep++;
return sum;
}
//______________________________________________________________________________
template<class Element>
Element TVectorT<Element>::Min() const
{
// return minimum vector element value
Assert(IsValid());
const Int_t index = TMath::LocMin(fNrows,fElements);
return fElements[index];
}
//______________________________________________________________________________
template<class Element>
Element TVectorT<Element>::Max() const
{
// return maximum vector element value
Assert(IsValid());
const Int_t index = TMath::LocMax(fNrows,fElements);
return fElements[index];
}
//______________________________________________________________________________
template<class Element>
TVectorT<Element> &TVectorT<Element>::operator=(const TVectorT<Element> &source)
{
// Notice that this assignment does NOT change the ownership :
// if the storage space was adopted, source is copied to
// this space .
if (!AreCompatible(*this,source)) {
Error("operator=(const TVectorT<Element> &)","vectors not compatible");
Invalidate();
return *this;
}
if (this != &source) {
TObject::operator=(source);
memcpy(fElements,source.GetMatrixArray(),fNrows*sizeof(Element));
}
return *this;
}
//______________________________________________________________________________
template<class Element>
TVectorT<Element> &TVectorT<Element>::operator=(const TMatrixTRow_const<Element> &mr)
{
// Assign a matrix row to a vector.
Assert(IsValid());
const TMatrixTBase<Element> *mt = mr.GetMatrix();
Assert(mt->IsValid());
if (mt->GetColLwb() != fRowLwb || mt->GetNcols() != fNrows) {
Error("operator=(const TMatrixTRow_const &)","vector and row not compatible");
Invalidate();
return *this;
}
const Int_t inc = mr.GetInc();
const Element *rp = mr.GetPtr(); // Row ptr
Element *ep = this->GetMatrixArray(); // Vector ptr
const Element * const fp = ep+fNrows;
while (ep < fp) {
*ep++ = *rp;
rp += inc;
}
Assert(rp == mr.GetPtr()+mt->GetNcols());
return *this;
}
//______________________________________________________________________________
template<class Element>
TVectorT<Element> &TVectorT<Element>::operator=(const TMatrixTColumn_const<Element> &mc)
{
// Assign a matrix column to a vector.
Assert(IsValid());
const TMatrixTBase<Element> *mt = mc.GetMatrix();
Assert(mt->IsValid());
if (mt->GetRowLwb() != fRowLwb || mt->GetNrows() != fNrows) {
Error("operator=(const TMatrixTColumn_const &)","vector and column not compatible");
Invalidate();
return *this;
}
const Int_t inc = mc.GetInc();
const Element *cp = mc.GetPtr(); // Column ptr
Element *ep = this->GetMatrixArray(); // Vector ptr
const Element * const fp = ep+fNrows;
while (ep < fp) {
*ep++ = *cp;
cp += inc;
}
Assert(cp == mc.GetPtr()+mt->GetNoElements());
return *this;
}
//______________________________________________________________________________
template<class Element>
TVectorT<Element> &TVectorT<Element>::operator=(const TMatrixTDiag_const<Element> &md)
{
// Assign the matrix diagonal to a vector.
Assert(IsValid());
const TMatrixTBase<Element> *mt = md.GetMatrix();
Assert(mt->IsValid());
if (md.GetNdiags() != fNrows) {
Error("operator=(const TMatrixTDiag_const &)","vector and matrix-diagonal not compatible");
Invalidate();
return *this;
}
const Int_t inc = md.GetInc();
const Element *dp = md.GetPtr(); // Diag ptr
Element *ep = this->GetMatrixArray(); // Vector ptr
const Element * const fp = ep+fNrows;
while (ep < fp) {
*ep++ = *dp;
dp += inc;
}
Assert(dp < md.GetPtr()+mt->GetNoElements()+inc);
return *this;
}
//______________________________________________________________________________
template<class Element>
TVectorT<Element> &TVectorT<Element>::operator=(const TMatrixTSparseRow_const<Element> &mr)
{
// Assign a sparse matrix row to a vector. The matrix row is implicitly transposed
// to allow the assignment in the strict sense.
Assert(IsValid());
const TMatrixTBase<Element> *mt = mr.GetMatrix();
Assert(mt->IsValid());
if (mt->GetColLwb() != fRowLwb || mt->GetNcols() != fNrows) {
Error("operator=(const TMatrixTSparseRow_const &)","vector and row not compatible");
Invalidate();
return *this;
}
const Int_t nIndex = mr.GetNindex();
const Element * const prData = mr.GetDataPtr(); // Row Data ptr
const Int_t * const prCol = mr.GetColPtr(); // Col ptr
Element * const pvData = this->GetMatrixArray(); // Vector ptr
memset(pvData,0,fNrows*sizeof(Element));
for (Int_t index = 0; index < nIndex; index++) {
const Int_t icol = prCol[index];
pvData[icol] = prData[index];
}
return *this;
}
//______________________________________________________________________________
template<class Element>
TVectorT<Element> &TVectorT<Element>::operator=(const TMatrixTSparseDiag_const<Element> &md)
{
// Assign a sparse matrix diagonal to a vector.
Assert(IsValid());
const TMatrixTBase<Element> *mt = md.GetMatrix();
Assert(mt->IsValid());
if (md.GetNdiags() != fNrows) {
Error("operator=(const TMatrixTSparseDiag_const &)","vector and matrix-diagonal not compatible");
Invalidate();
return *this;
}
Element * const pvData = this->GetMatrixArray();
for (Int_t idiag = 0; idiag < fNrows; idiag++)
pvData[idiag] = md(idiag);
return *this;
}
//______________________________________________________________________________
template<class Element>
TVectorT<Element> &TVectorT<Element>::operator=(Element val)
{
// Assign val to every element of the vector.
Assert(IsValid());
Element *ep = this->GetMatrixArray();
const Element * const fp = ep+fNrows;
while (ep < fp)
*ep++ = val;
return *this;
}
//______________________________________________________________________________
template<class Element>
TVectorT<Element> &TVectorT<Element>::operator+=(Element val)
{
// Add val to every element of the vector.
Assert(IsValid());
Element *ep = this->GetMatrixArray();
const Element * const fp = ep+fNrows;
while (ep < fp)
*ep++ += val;
return *this;
}
//______________________________________________________________________________
template<class Element>
TVectorT<Element> &TVectorT<Element>::operator-=(Element val)
{
// Subtract val from every element of the vector.
Assert(IsValid());
Element *ep = this->GetMatrixArray();
const Element * const fp = ep+fNrows;
while (ep < fp)
*ep++ -= val;
return *this;
}
//______________________________________________________________________________
template<class Element>
TVectorT<Element> &TVectorT<Element>::operator*=(Element val)
{
// Multiply every element of the vector with val.
Assert(IsValid());
Element *ep = this->GetMatrixArray();
const Element * const fp = ep+fNrows;
while (ep < fp)
*ep++ *= val;
return *this;
}
//______________________________________________________________________________
template<class Element>
TVectorT<Element> &TVectorT<Element>::operator+=(const TVectorT<Element> &source)
{
// Add vector source
if (!AreCompatible(*this,source)) {
Error("operator+=(const TVectorT<Element> &)","vector's not compatible");
Invalidate();
return *this;
}
const Element *sp = source.GetMatrixArray();
Element *tp = this->GetMatrixArray();
const Element * const tp_last = tp+fNrows;
while (tp < tp_last)
*tp++ += *sp++;
return *this;
}
//______________________________________________________________________________
template<class Element>
TVectorT<Element> &TVectorT<Element>::operator-=(const TVectorT<Element> &source)
{
// Subtract vector source
if (!AreCompatible(*this,source)) {
Error("operator-=(const TVectorT<Element> &)","vector's not compatible");
Invalidate();
return *this;
}
const Element *sp = source.GetMatrixArray();
Element *tp = this->GetMatrixArray();
const Element * const tp_last = tp+fNrows;
while (tp < tp_last)
*tp++ -= *sp++;
return *this;
}
//______________________________________________________________________________
template<class Element>
TVectorT<Element> &TVectorT<Element>::operator*=(const TMatrixT<Element> &a)
{
// "Inplace" multiplication target = A*target. A needn't be a square one
// If target has to be resized, it should own the storage: fIsOwner = kTRUE
Assert(IsValid());
Assert(a.IsValid());
if (a.GetNcols() != fNrows || a.GetColLwb() != fRowLwb) {
Error("operator*=(const TMatrixT &)","vector and matrix incompatible");
Invalidate();
return *this;
}
const Bool_t doResize = (fNrows != a.GetNrows() || fRowLwb != a.GetRowLwb());
if (doResize && !fIsOwner) {
Error("operator*=(const TMatrixT &)","vector has to be resized but not owner");
Invalidate();
return *this;
}
Element work[kWorkMax];
Bool_t isAllocated = kFALSE;
Element *elements_old = work;
const Int_t nrows_old = fNrows;
if (nrows_old > kWorkMax) {
isAllocated = kTRUE;
elements_old = new Element[nrows_old];
}
memcpy(elements_old,fElements,nrows_old*sizeof(Element));
if (doResize) {
const Int_t rowlwb_new = a.GetRowLwb();
const Int_t nrows_new = a.GetNrows();
ResizeTo(rowlwb_new,rowlwb_new+nrows_new-1);
}
memset(fElements,0,fNrows*sizeof(Element));
const Element *mp = a.GetMatrixArray(); // Matrix row ptr
Element *tp = this->GetMatrixArray(); // Target vector ptr
#ifdef CBLAS
if (typeid(Element) == typeid(Double_t))
cblas_dgemv(CblasRowMajor,CblasNoTrans,a.GetNrows(),a.GetNcols(),1.0,mp,
a.GetNcols(),elements_old,1,0.0,tp,1);
else if (typeid(Element) != typeid(Float_t))
cblas_sgemv(CblasRowMajor,CblasNoTrans,a.GetNrows(),a.GetNcols(),1.0,mp,
a.GetNcols(),elements_old,1,0.0,tp,1);
else
Error("operator*=","type %s not implemented in BLAS library",typeid(Element));
#else
const Element * const tp_last = tp+fNrows;
while (tp < tp_last) {
Element sum = 0;
for (const Element *sp = elements_old; sp < elements_old+nrows_old; )
sum += *sp++ * *mp++;
*tp++ = sum;
}
Assert(mp == a.GetMatrixArray()+a.GetNoElements());
#endif
if (isAllocated)
delete [] elements_old;
return *this;
}
//______________________________________________________________________________
template<class Element>
TVectorT<Element> &TVectorT<Element>::operator*=(const TMatrixTSparse<Element> &a)
{
// "Inplace" multiplication target = A*target. A needn't be a square one
// If target has to be resized, it should own the storage: fIsOwner = kTRUE
Assert(IsValid());
Assert(a.IsValid());
if (a.GetNcols() != fNrows || a.GetColLwb() != fRowLwb) {
Error("operator*=(const TMatrixTSparse &)","vector and matrix incompatible");
Invalidate();
return *this;
}
const Bool_t doResize = (fNrows != a.GetNrows() || fRowLwb != a.GetRowLwb());
if (doResize && !fIsOwner) {
Error("operator*=(const TMatrixTSparse &)","vector has to be resized but not owner");
Invalidate();
return *this;
}
Element work[kWorkMax];
Bool_t isAllocated = kFALSE;
Element *elements_old = work;
const Int_t nrows_old = fNrows;
if (nrows_old > kWorkMax) {
isAllocated = kTRUE;
elements_old = new Element[nrows_old];
}
memcpy(elements_old,fElements,nrows_old*sizeof(Element));
if (doResize) {
const Int_t rowlwb_new = a.GetRowLwb();
const Int_t nrows_new = a.GetNrows();
ResizeTo(rowlwb_new,rowlwb_new+nrows_new-1);
}
memset(fElements,0,fNrows*sizeof(Element));
const Int_t * const pRowIndex = a.GetRowIndexArray();
const Int_t * const pColIndex = a.GetColIndexArray();
const Element * const mp = a.GetMatrixArray(); // Matrix row ptr
const Element * const sp = elements_old;
Element * tp = this->GetMatrixArray(); // Target vector ptr
for (Int_t irow = 0; irow < fNrows; irow++) {
const Int_t sIndex = pRowIndex[irow];
const Int_t eIndex = pRowIndex[irow+1];
Element sum = 0.0;
for (Int_t index = sIndex; index < eIndex; index++) {
const Int_t icol = pColIndex[index];
sum += mp[index]*sp[icol];
}
tp[irow] = sum;
}
if (isAllocated)
delete [] elements_old;
return *this;
}
//______________________________________________________________________________
template<class Element>
TVectorT<Element> &TVectorT<Element>::operator*=(const TMatrixTSym<Element> &a)
{
// "Inplace" multiplication target = A*target. A is symmetric .
// vector size will not change
Assert(IsValid());
Assert(a.IsValid());
if (a.GetNcols() != fNrows || a.GetColLwb() != fRowLwb) {
Error("operator*=(const TMatrixTSym &)","vector and matrix incompatible");
Invalidate();
return *this;
}
Element work[kWorkMax];
Bool_t isAllocated = kFALSE;
Element *elements_old = work;
const Int_t nrows_old = fNrows;
if (nrows_old > kWorkMax) {
isAllocated = kTRUE;
elements_old = new Element[nrows_old];
}
memcpy(elements_old,fElements,nrows_old*sizeof(Element));
memset(fElements,0,fNrows*sizeof(Element));
const Element *mp = a.GetMatrixArray(); // Matrix row ptr
Element *tp = this->GetMatrixArray(); // Target vector ptr
#ifdef CBLAS
if (typeid(Element) == typeid(Double_t))
cblas_dsymv(CblasRowMajor,CblasUpper,fNrows,1.0,mp,
fNrows,elements_old,1,0.0,tp,1);
else if (typeid(Element) != typeid(Float_t))
cblas_ssymv(CblasRowMajor,CblasUpper,fNrows,1.0,mp,
fNrows,elements_old,1,0.0,tp,1);
else
Error("operator*=","type %s not implemented in BLAS library",typeid(Element));
#else
const Element * const tp_last = tp+fNrows;
while (tp < tp_last) {
Element sum = 0;
for (const Element *sp = elements_old; sp < elements_old+nrows_old; )
sum += *sp++ * *mp++;
*tp++ = sum;
}
Assert(mp == a.GetMatrixArray()+a.GetNoElements());
#endif
if (isAllocated)
delete [] elements_old;
return *this;
}
//______________________________________________________________________________
template<class Element>
Bool_t TVectorT<Element>::operator==(Element val) const
{
// Are all vector elements equal to val?
Assert(IsValid());
const Element *ep = this->GetMatrixArray();
const Element * const fp = ep+fNrows;
while (ep < fp)
if (!(*ep++ == val))
return kFALSE;
return kTRUE;
}
//______________________________________________________________________________
template<class Element>
Bool_t TVectorT<Element>::operator!=(Element val) const
{
// Are all vector elements not equal to val?
Assert(IsValid());
const Element *ep = this->GetMatrixArray();
const Element * const fp = ep+fNrows;
while (ep < fp)
if (!(*ep++ != val))
return kFALSE;
return kTRUE;
}
//______________________________________________________________________________
template<class Element>
Bool_t TVectorT<Element>::operator<(Element val) const
{
// Are all vector elements < val?
Assert(IsValid());
const Element *ep = this->GetMatrixArray();
const Element * const fp = ep+fNrows;
while (ep < fp)
if (!(*ep++ < val))
return kFALSE;
return kTRUE;
}
//______________________________________________________________________________
template<class Element>
Bool_t TVectorT<Element>::operator<=(Element val) const
{
// Are all vector elements <= val?
Assert(IsValid());
const Element *ep = this->GetMatrixArray();
const Element * const fp = ep+fNrows;
while (ep < fp)
if (!(*ep++ <= val))
return kFALSE;
return kTRUE;
}
//______________________________________________________________________________
template<class Element>
Bool_t TVectorT<Element>::operator>(Element val) const
{
// Are all vector elements > val?
Assert(IsValid());
const Element *ep = this->GetMatrixArray();
const Element * const fp = ep+fNrows;
while (ep < fp)
if (!(*ep++ > val))
return kFALSE;
return kTRUE;
}
//______________________________________________________________________________
template<class Element>
Bool_t TVectorT<Element>::operator>=(Element val) const
{
// Are all vector elements >= val?
Assert(IsValid());
const Element *ep = this->GetMatrixArray();
const Element * const fp = ep+fNrows;
while (ep < fp)
if (!(*ep++ >= val))
return kFALSE;
return kTRUE;
}
//______________________________________________________________________________
template<class Element>
Bool_t TVectorT<Element>::MatchesNonZeroPattern(const TVectorT<Element> &select)
{
if (!AreCompatible(*this,select)) {
Error("MatchesNonZeroPattern(const TVectorT&)","vector's not compatible");
return kFALSE;
}
const Element *sp = select.GetMatrixArray();
const Element *ep = this->GetMatrixArray();
const Element * const fp = ep+fNrows;
while (ep < fp) {
if (*sp == 0.0 && *ep != 0.0)
return kFALSE;
sp++; ep++;
}
return kTRUE;
}
//______________________________________________________________________________
template<class Element>
Bool_t TVectorT<Element>::SomePositive(const TVectorT<Element> &select)
{
if (!AreCompatible(*this,select)) {
Error("SomePositive(const TVectorT&)","vector's not compatible");
return kFALSE;
}
const Element *sp = select.GetMatrixArray();
const Element *ep = this->GetMatrixArray();
const Element * const fp = ep+fNrows;
while (ep < fp) {
if (*sp != 0.0 && *ep <= 0.0)
return kFALSE;
sp++; ep++;
}
return kTRUE;
}
//______________________________________________________________________________
template<class Element>
void TVectorT<Element>::AddSomeConstant(Element val,const TVectorT<Element> &select)
{
if (!AreCompatible(*this,select))
Error("AddSomeConstant(Element,const TVectorT&)(const TVectorT&)","vector's not compatible");
const Element *sp = select.GetMatrixArray();
Element *ep = this->GetMatrixArray();
const Element * const fp = ep+fNrows;
while (ep < fp) {
if (*sp)
*ep += val;
sp++; ep++;
}
}
extern Double_t Drand(Double_t &ix);
//______________________________________________________________________________
template<class Element>
void TVectorT<Element>::Randomize(Element alpha,Element beta,Double_t &seed)
{
// randomize vector elements value
Assert(IsValid());
const Element scale = beta-alpha;
const Element shift = alpha/scale;
Element * ep = GetMatrixArray();
const Element * const fp = ep+fNrows;
while (ep < fp)
*ep++ = scale*(Drand(seed)+shift);
}
//______________________________________________________________________________
template<class Element>
TVectorT<Element> &TVectorT<Element>::Apply(const TElementActionT<Element> &action)
{
// Apply action to each element of the vector.
Assert(IsValid());
for (Element *ep = fElements; ep < fElements+fNrows; ep++)
action.Operation(*ep);
return *this;
}
//______________________________________________________________________________
template<class Element>
TVectorT<Element> &TVectorT<Element>::Apply(const TElementPosActionT<Element> &action)
{
// Apply action to each element of the vector. In action the location
// of the current element is known.
Assert(IsValid());
Element *ep = fElements;
for (action.fI = fRowLwb; action.fI < fRowLwb+fNrows; action.fI++)
action.Operation(*ep++);
Assert(ep == fElements+fNrows);
return *this;
}
//______________________________________________________________________________
template<class Element>
void TVectorT<Element>::Draw(Option_t *option)
{
// Draw this vector using an intermediate histogram
// The histogram is named "TVectorT" by default and no title
//create the hist utility manager (a plugin)
TVirtualUtilHist *util = (TVirtualUtilHist*)gROOT->GetListOfSpecials()->FindObject("R__TVirtualUtilHist");
if (!util) {
TPluginHandler *h;
if ((h = gROOT->GetPluginManager()->FindHandler("TVirtualUtilHist"))) {
if (h->LoadPlugin() == -1)
return;
h->ExecPlugin(0);
util = (TVirtualUtilHist*)gROOT->GetListOfSpecials()->FindObject("R__TVirtualUtilHist");
}
}
util->PaintVector(*this,option);
}
//______________________________________________________________________________
template<class Element>
void TVectorT<Element>::Print(Option_t *flag) const
{
// Print the vector as a list of elements.
Assert(IsValid());
printf("\nVector (%d) %s is as follows",fNrows,flag);
printf("\n\n | %6d |", 1);
printf("\n%s\n", "------------------");
for (Int_t i = 0; i < fNrows; i++) {
printf("%4d |",i+fRowLwb);
//printf("%11.4g \n",(*this)(i+fRowLwb));
printf("%g \n",(*this)(i+fRowLwb));
}
printf("\n");
}
//______________________________________________________________________________
template<class Element>
Bool_t operator==(const TVectorT<Element> &v1,const TVectorT<Element> &v2)
{
// Check to see if two vectors are identical.
if (!AreCompatible(v1,v2)) return kFALSE;
return (memcmp(v1.GetMatrixArray(),v2.GetMatrixArray(),v1.GetNrows()*sizeof(Element)) == 0);
}
//______________________________________________________________________________
template<class Element>
Element operator*(const TVectorT<Element> &v1,const TVectorT<Element> &v2)
{
// Compute the scalar product.
if (!AreCompatible(v1,v2)) {
Error("operator*(const TVectorT<Element> &,const TVectorT<Element> &)","vector's are incompatible");
return 0.0;
}
const Element *v1p = v1.GetMatrixArray();
const Element *v2p = v2.GetMatrixArray();
Element sum = 0.0;
const Element * const fv1p = v1p+v1.GetNrows();
while (v1p < fv1p)
sum += *v1p++ * *v2p++;
return sum;
}
//______________________________________________________________________________
template<class Element>
TVectorT<Element> operator+(const TVectorT<Element> &source1,const TVectorT<Element> &source2)
{
TVectorT<Element> target = source1;
target += source2;
return target;
}
//______________________________________________________________________________
template<class Element>
TVectorT<Element> operator-(const TVectorT<Element> &source1,const TVectorT<Element> &source2)
{
TVectorT<Element> target = source1;
target -= source2;
return target;
}
//______________________________________________________________________________
template<class Element>
TVectorT<Element> operator*(const TMatrixT<Element> &a,const TVectorT<Element> &source)
{
Assert(a.IsValid());
TVectorT<Element> target(a.GetRowLwb(),a.GetRowUpb());
return Add(target,Element(1.0),a,source);
}
//______________________________________________________________________________
template<class Element>
TVectorT<Element> operator*(const TMatrixTSym<Element> &a,const TVectorT<Element> &source)
{
Assert(a.IsValid());
TVectorT<Element> target(a.GetRowLwb(),a.GetRowUpb());
return Add(target,Element(1.0),a,source);
}
//______________________________________________________________________________
template<class Element>
TVectorT<Element> operator*(const TMatrixTSparse<Element> &a,const TVectorT<Element> &source)
{
Assert(a.IsValid());
TVectorT<Element> target(a.GetRowLwb(),a.GetRowUpb());
return Add(target,Element(1.0),a,source);
}
//______________________________________________________________________________
template<class Element>
TVectorT<Element> operator*(Element val,const TVectorT<Element> &source)
{
TVectorT<Element> target = source;
target *= val;
return target;
}
//______________________________________________________________________________
template<class Element>
TVectorT<Element> &Add(TVectorT<Element> &target,Element scalar,const TVectorT<Element> &source)
{
// Modify addition: target += scalar * source.
if (!AreCompatible(target,source)) {
Error("Add(TVectorT<Element> &,Element,const TVectorT<Element> &)","vector's are incompatible");
target.Invalidate();
return target;
}
const Element * sp = source.GetMatrixArray();
Element * tp = target.GetMatrixArray();
const Element * const ftp = tp+target.GetNrows();
if (scalar == 1.0 ) {
while ( tp < ftp )
*tp++ += *sp++;
} else if (scalar == -1.0) {
while ( tp < ftp )
*tp++ -= *sp++;
} else {
while ( tp < ftp )
*tp++ += scalar * *sp++;
}
return target;
}
//______________________________________________________________________________
template<class Element>
TVectorT<Element> &Add(TVectorT<Element> &target,Element scalar,
const TMatrixT<Element> &a,const TVectorT<Element> &source)
{
// Modify addition: target += scalar * A * source.
Assert(target.IsValid());
Assert(a.IsValid());
if (a.GetNrows() != target.GetNrows() || a.GetRowLwb() != target.GetLwb()) {
Error("Add(TVectorT &,const TMatrixT &,const TVectorT &)","target vector and matrix are incompatible");
target.Invalidate();
return target;
}
Assert(source.IsValid());
if (a.GetNcols() != source.GetNrows() || a.GetColLwb() != source.GetLwb()) {
Error("Add(TVectorT &,const TMatrixT &,const TVectorT &)","source vector and matrix are incompatible");
target.Invalidate();
return target;
}
const Element * const sp = source.GetMatrixArray(); // sources vector ptr
const Element * mp = a.GetMatrixArray(); // Matrix row ptr
Element * tp = target.GetMatrixArray(); // Target vector ptr
#ifdef CBLAS
if (typeid(Element) == typeid(Double_t))
cblas_dsymv(CblasRowMajor,CblasUpper,fNrows,scalar,mp,
fNrows,sp,1,0.0,tp,1);
else if (typeid(Element) != typeid(Float_t))
cblas_ssymv(CblasRowMajor,CblasUpper,fNrows,scalar,mp,
fNrows,sp,1,0.0,tp,1);
else
Error("operator*=","type %s not implemented in BLAS library",typeid(Element));
#else
const Element * const sp_last = sp+source.GetNrows();
const Element * const tp_last = tp+target.GetNrows();
if (scalar == 1.0) {
while (tp < tp_last) {
const Element *sp1 = sp;
Element sum = 0;
while (sp1 < sp_last)
sum += *sp1++ * *mp++;
*tp++ += sum;
}
} else if (scalar == -1.0) {
while (tp < tp_last) {
const Element *sp1 = sp;
Element sum = 0;
while (sp1 < sp_last)
sum += *sp1++ * *mp++;
*tp++ -= sum;
}
} else {
while (tp < tp_last) {
const Element *sp1 = sp;
Element sum = 0;
while (sp1 < sp_last)
sum += *sp1++ * *mp++;
*tp++ += scalar * sum;
}
}
Assert(mp == a.GetMatrixArray()+a.GetNoElements());
#endif
return target;
}
//______________________________________________________________________________
template<class Element>
TVectorT<Element> &Add(TVectorT<Element> &target,Element scalar,
const TMatrixTSym<Element> &a,const TVectorT<Element> &source)
{
// Modify addition: target += A * source.
Assert(target.IsValid());
Assert(source.IsValid());
Assert(a.IsValid());
if (a.GetNrows() != target.GetNrows() || a.GetRowLwb() != target.GetLwb()) {
Error("Add(TVectorT &,const TMatrixT &,const TVectorT &)","target vector and matrix are incompatible");
target.Invalidate();
return target;
}
const Element * const sp = source.GetMatrixArray(); // sources vector ptr
const Element * mp = a.GetMatrixArray(); // Matrix row ptr
Element * tp = target.GetMatrixArray(); // Target vector ptr
#ifdef CBLAS
if (typeid(Element) == typeid(Double_t))
cblas_dsymv(CblasRowMajor,CblasUpper,fNrows,1.0,mp,
fNrows,sp,1,0.0,tp,1);
else if (typeid(Element) != typeid(Float_t))
cblas_ssymv(CblasRowMajor,CblasUpper,fNrows,1.0,mp,
fNrows,sp,1,0.0,tp,1);
else
Error("operator*=","type %s not implemented in BLAS library",typeid(Element));
#else
const Element * const sp_last = sp+source.GetNrows();
const Element * const tp_last = tp+target.GetNrows();
if (scalar == 1.0) {
while (tp < tp_last) {
const Element *sp1 = sp;
Element sum = 0;
while (sp1 < sp_last)
sum += *sp1++ * *mp++;
*tp++ += sum;
}
} else if (scalar == -1.0) {
while (tp < tp_last) {
const Element *sp1 = sp;
Element sum = 0;
while (sp1 < sp_last)
sum += *sp1++ * *mp++;
*tp++ -= sum;
}
} else {
while (tp < tp_last) {
const Element *sp1 = sp;
Element sum = 0;
while (sp1 < sp_last)
sum += *sp1++ * *mp++;
*tp++ += scalar * sum;
}
}
Assert(mp == a.GetMatrixArray()+a.GetNoElements());
#endif
return target;
}
//______________________________________________________________________________
template<class Element>
TVectorT<Element> &Add(TVectorT<Element> &target,Element scalar,
const TMatrixTSparse<Element> &a,const TVectorT<Element> &source)
{
// Modify addition: target += A * source.
Assert(target.IsValid());
Assert(a.IsValid());
if (a.GetNrows() != target.GetNrows() || a.GetRowLwb() != target.GetLwb()) {
Error("Add(TVectorT &,const TMatrixT &,const TVectorT &)","target vector and matrix are incompatible");
target.Invalidate();
return target;
}
Assert(source.IsValid());
if (a.GetNcols() != source.GetNrows() || a.GetColLwb() != source.GetLwb()) {
Error("Add(TVectorT &,const TMatrixT &,const TVectorT &)","source vector and matrix are incompatible");
target.Invalidate();
return target;
}
const Int_t * const pRowIndex = a.GetRowIndexArray();
const Int_t * const pColIndex = a.GetColIndexArray();
const Element * const mp = a.GetMatrixArray(); // Matrix row ptr
const Element * const sp = source.GetMatrixArray(); // Source vector ptr
Element * tp = target.GetMatrixArray(); // Target vector ptr
if (scalar == 1.0) {
for (Int_t irow = 0; irow < a.GetNrows(); irow++) {
const Int_t sIndex = pRowIndex[irow];
const Int_t eIndex = pRowIndex[irow+1];
Element sum = 0.0;
for (Int_t index = sIndex; index < eIndex; index++) {
const Int_t icol = pColIndex[index];
sum += mp[index]*sp[icol];
}
tp[irow] += sum;
}
} else if (scalar == -1.0) {
for (Int_t irow = 0; irow < a.GetNrows(); irow++) {
const Int_t sIndex = pRowIndex[irow];
const Int_t eIndex = pRowIndex[irow+1];
Element sum = 0.0;
for (Int_t index = sIndex; index < eIndex; index++) {
const Int_t icol = pColIndex[index];
sum += mp[index]*sp[icol];
}
tp[irow] -= sum;
}
} else {
for (Int_t irow = 0; irow < a.GetNrows(); irow++) {
const Int_t sIndex = pRowIndex[irow];
const Int_t eIndex = pRowIndex[irow+1];
Element sum = 0.0;
for (Int_t index = sIndex; index < eIndex; index++) {
const Int_t icol = pColIndex[index];
sum += mp[index]*sp[icol];
}
tp[irow] += scalar * sum;
}
}
return target;
}
//______________________________________________________________________________
template<class Element>
TVectorT<Element> &AddElemMult(TVectorT<Element> &target,Element scalar,
const TVectorT<Element> &source1,const TVectorT<Element> &source2)
{
// Modify addition: target += scalar * ElementMult(source1,source2) .
if (!(AreCompatible(target,source1) && AreCompatible(target,source1))) {
Error("AddElemMult(TVectorT<Element> &,Element,const TVectorT<Element> &,const TVectorT<Element> &)",
"vector's are incompatible");
target.Invalidate();
return target;
}
const Element * sp1 = source1.GetMatrixArray();
const Element * sp2 = source2.GetMatrixArray();
Element * tp = target.GetMatrixArray();
const Element * const ftp = tp+target.GetNrows();
if (scalar == 1.0 ) {
while ( tp < ftp )
*tp++ += *sp1++ * *sp2++;
} else if (scalar == -1.0) {
while ( tp < ftp )
*tp++ -= *sp1++ * *sp2++;
} else {
while ( tp < ftp )
*tp++ += scalar * *sp1++ * *sp2++;
}
return target;
}
//______________________________________________________________________________
template<class Element>
TVectorT<Element> &AddElemMult(TVectorT<Element> &target,Element scalar,
const TVectorT<Element> &source1,const TVectorT<Element> &source2,const TVectorT<Element> &select)
{
// Modify addition: target += scalar * ElementMult(source1,source2) only for those elements
// where select[i] != 0.0
if (!( AreCompatible(target,source1) && AreCompatible(target,source1) &&
AreCompatible(target,select) )) {
Error("AddElemMult(TVectorT<Element> &,Element,const TVectorT<Element> &,const TVectorT<Element> &,onst TVectorT<Element> &)",
"vector's are incompatible");
target.Invalidate();
return target;
}
const Element * sp1 = source1.GetMatrixArray();
const Element * sp2 = source2.GetMatrixArray();
const Element * mp = select.GetMatrixArray();
Element * tp = target.GetMatrixArray();
const Element * const ftp = tp+target.GetNrows();
if (scalar == 1.0 ) {
while ( tp < ftp ) {
if (*mp) *tp += *sp1 * *sp2;
mp++; tp++; sp1++; sp2++;
}
} else if (scalar == -1.0) {
while ( tp < ftp ) {
if (*mp) *tp -= *sp1 * *sp2;
mp++; tp++; sp1++; sp2++;
}
} else {
while ( tp < ftp ) {
if (*mp) *tp += scalar * *sp1 * *sp2;
mp++; tp++; sp1++; sp2++;
}
}
return target;
}
//______________________________________________________________________________
template<class Element>
TVectorT<Element> &AddElemDiv(TVectorT<Element> &target,Element scalar,
const TVectorT<Element> &source1,const TVectorT<Element> &source2)
{
// Modify addition: target += scalar * ElementMult(source1,source2) .
if (!(AreCompatible(target,source1) && AreCompatible(target,source1))) {
Error("AddElemDiv(TVectorT<Element> &,Element,const TVectorT<Element> &,const TVectorT<Element> &)",
"vector's are incompatible");
target.Invalidate();
return target;
}
const Element * sp1 = source1.GetMatrixArray();
const Element * sp2 = source2.GetMatrixArray();
Element * tp = target.GetMatrixArray();
const Element * const ftp = tp+target.GetNrows();
if (scalar == 1.0 ) {
while ( tp < ftp )
*tp++ += *sp1++ / *sp2++;
} else if (scalar == -1.0) {
while ( tp < ftp )
*tp++ -= *sp1++ / *sp2++;
} else {
while ( tp < ftp )
*tp++ += scalar * *sp1++ / *sp2++;
}
return target;
}
//______________________________________________________________________________
template<class Element>
TVectorT<Element> &AddElemDiv(TVectorT<Element> &target,Element scalar,
const TVectorT<Element> &source1,const TVectorT<Element> &source2,const TVectorT<Element> &select)
{
// Modify addition: target += scalar * ElementMult(source1,source2) only for those elements
// where select[i] != 0.0
if (!( AreCompatible(target,source1) && AreCompatible(target,source1) &&
AreCompatible(target,select) )) {
Error("AddElemDiv(TVectorT<Element> &,Element,const TVectorT<Element> &,const TVectorT<Element> &,onst TVectorT<Element> &)",
"vector's are incompatible");
target.Invalidate();
return target;
}
const Element * sp1 = source1.GetMatrixArray();
const Element * sp2 = source2.GetMatrixArray();
const Element * mp = select.GetMatrixArray();
Element * tp = target.GetMatrixArray();
const Element * const ftp = tp+target.GetNrows();
if (scalar == 1.0 ) {
while ( tp < ftp ) {
if (*mp) *tp += *sp1 / *sp2;
mp++; tp++; sp1++; sp2++;
}
} else if (scalar == -1.0) {
while ( tp < ftp ) {
if (*mp) *tp -= *sp1 / *sp2;
mp++; tp++; sp1++; sp2++;
}
} else {
while ( tp < ftp ) {
if (*mp) *tp += scalar * *sp1 / *sp2;
mp++; tp++; sp1++; sp2++;
}
}
return target;
}
//______________________________________________________________________________
template<class Element>
TVectorT<Element> &ElementMult(TVectorT<Element> &target,const TVectorT<Element> &source)
{
// Multiply target by the source, element-by-element.
if (!AreCompatible(target,source)) {
Error("ElementMult(TVectorT<Element> &,const TVectorT<Element> &)","vector's are incompatible");
target.Invalidate();
return target;
}
const Element * sp = source.GetMatrixArray();
Element * tp = target.GetMatrixArray();
const Element * const ftp = tp+target.GetNrows();
while ( tp < ftp )
*tp++ *= *sp++;
return target;
}
//______________________________________________________________________________
template<class Element>
TVectorT<Element> &ElementMult(TVectorT<Element> &target,const TVectorT<Element> &source,const TVectorT<Element> &select)
{
// Multiply target by the source, element-by-element only where select[i] != 0.0
if (!(AreCompatible(target,source) && AreCompatible(target,select))) {
Error("ElementMult(TVectorT<Element> &,const TVectorT<Element> &,const TVectorT<Element> &)","vector's are incompatible");
target.Invalidate();
return target;
}
const Element * sp = source.GetMatrixArray();
const Element * mp = select.GetMatrixArray();
Element * tp = target.GetMatrixArray();
const Element * const ftp = tp+target.GetNrows();
while ( tp < ftp ) {
if (*mp) *tp *= *sp;
mp++; tp++; sp++;
}
return target;
}
//______________________________________________________________________________
template<class Element>
TVectorT<Element> &ElementDiv(TVectorT<Element> &target,const TVectorT<Element> &source)
{
// Divide target by the source, element-by-element.
if (!AreCompatible(target,source)) {
Error("ElementDiv(TVectorT<Element> &,const TVectorT<Element> &)","vector's are incompatible");
target.Invalidate();
return target;
}
const Element * sp = source.GetMatrixArray();
Element * tp = target.GetMatrixArray();
const Element * const ftp = tp+target.GetNrows();
while ( tp < ftp )
*tp++ /= *sp++;
return target;
}
//______________________________________________________________________________
template<class Element>
TVectorT<Element> &ElementDiv(TVectorT<Element> &target,const TVectorT<Element> &source,const TVectorT<Element> &select)
{
// Divide target by the source, element-by-element only where select[i] != 0.0
if (!AreCompatible(target,source)) {
Error("ElementDiv(TVectorT<Element> &,const TVectorT<Element> &,const TVectorT<Element> &)","vector's are incompatible");
target.Invalidate();
return target;
}
const Element * sp = source.GetMatrixArray();
const Element * mp = select.GetMatrixArray();
Element * tp = target.GetMatrixArray();
const Element * const ftp = tp+target.GetNrows();
while ( tp < ftp ) {
if (*mp) *tp /= *sp;
mp++; tp++; sp++;
}
return target;
}
//______________________________________________________________________________
template<class Element1,class Element2>
Bool_t AreCompatible(const TVectorT<Element1> &v1,const TVectorT<Element2> &v2,Int_t verbose)
{
if (!v1.IsValid()) {
if (verbose)
::Error("AreCompatible", "vector 1 not valid");
return kFALSE;
}
if (!v2.IsValid()) {
if (verbose)
::Error("AreCompatible", "vector 2 not valid");
return kFALSE;
}
if (v1.GetNrows() != v2.GetNrows() || v1.GetLwb() != v2.GetLwb()) {
if (verbose)
::Error("AreCompatible", "matrices 1 and 2 not compatible");
return kFALSE;
}
return kTRUE;
}
//______________________________________________________________________________
template<class Element>
void Compare(const TVectorT<Element> &v1,const TVectorT<Element> &v2)
{
// Compare two vectors and print out the result of the comparison.
if (!AreCompatible(v1,v2)) {
Error("Compare(const TVectorT<Element> &,const TVectorT<Element> &)","vectors are incompatible");
return;
}
printf("\n\nComparison of two TVectorTs:\n");
Element norm1 = 0; // Norm of the Matrices
Element norm2 = 0; // Norm of the Matrices
Element ndiff = 0; // Norm of the difference
Int_t imax = 0; // For the elements that differ most
Element difmax = -1;
const Element *mp1 = v1.GetMatrixArray(); // Vector element pointers
const Element *mp2 = v2.GetMatrixArray();
for (Int_t i = 0; i < v1.GetNrows(); i++) {
const Element mv1 = *mp1++;
const Element mv2 = *mp2++;
const Element diff = TMath::Abs(mv1-mv2);
if (diff > difmax) {
difmax = diff;
imax = i;
}
norm1 += TMath::Abs(mv1);
norm2 += TMath::Abs(mv2);
ndiff += TMath::Abs(diff);
}
imax += v1.GetLwb();
printf("\nMaximal discrepancy \t\t%g",difmax);
printf("\n occured at the point\t\t(%d)",imax);
const Element mv1 = v1(imax);
const Element mv2 = v2(imax);
printf("\n Vector 1 element is \t\t%g",mv1);
printf("\n Vector 2 element is \t\t%g",mv2);
printf("\n Absolute error v2[i]-v1[i]\t\t%g",mv2-mv1);
printf("\n Relative error\t\t\t\t%g\n",
(mv2-mv1)/TMath::Max(TMath::Abs(mv2+mv1)/2,(Element)1e-7));
printf("\n||Vector 1|| \t\t\t%g",norm1);
printf("\n||Vector 2|| \t\t\t%g",norm2);
printf("\n||Vector1-Vector2||\t\t\t\t%g",ndiff);
printf("\n||Vector1-Vector2||/sqrt(||Vector1|| ||Vector2||)\t%g\n\n",
ndiff/TMath::Max(TMath::Sqrt(norm1*norm2),1e-7));
}
//______________________________________________________________________________
template<class Element>
Bool_t VerifyVectorValue(const TVectorT<Element> &v,Element val,
Int_t verbose,Element maxDevAllow)
{
// Validate that all elements of vector have value val within maxDevAllow .
Int_t imax = 0;
Element maxDevObs = 0;
if (TMath::Abs(maxDevAllow) <= 0.0)
maxDevAllow = std::numeric_limits<Element>::epsilon();
for (Int_t i = v.GetLwb(); i <= v.GetUpb(); i++) {
const Element dev = TMath::Abs(v(i)-val);
if (dev > maxDevObs) {
imax = i;
maxDevObs = dev;
}
}
if (maxDevObs == 0)
return kTRUE;
if (verbose) {
printf("Largest dev for (%d); dev = |%g - %g| = %g\n",imax,v(imax),val,maxDevObs);
if(maxDevObs > maxDevAllow)
Error("VerifyVectorValue","Deviation > %g\n",maxDevAllow);
}
if(maxDevObs > maxDevAllow)
return kFALSE;
return kTRUE;
}
//______________________________________________________________________________
template<class Element>
Bool_t VerifyVectorIdentity(const TVectorT<Element> &v1,const TVectorT<Element> &v2,
Int_t verbose, Element maxDevAllow)
{
// Verify that elements of the two vectors are equal within maxDevAllow .
Int_t imax = 0;
Element maxDevObs = 0;
if (!AreCompatible(v1,v2))
return kFALSE;
if (TMath::Abs(maxDevAllow) <= 0.0)
maxDevAllow = std::numeric_limits<Element>::epsilon();
for (Int_t i = v1.GetLwb(); i <= v1.GetUpb(); i++) {
const Element dev = TMath::Abs(v1(i)-v2(i));
if (dev > maxDevObs) {
imax = i;
maxDevObs = dev;
}
}
if (maxDevObs == 0)
return kTRUE;
if (verbose) {
printf("Largest dev for (%d); dev = |%g - %g| = %g\n",imax,v1(imax),v2(imax),maxDevObs);
if(maxDevObs > maxDevAllow)
Error("VerifyVectorIdentity","Deviation > %g\n",maxDevAllow);
}
if(maxDevObs > maxDevAllow) {
return kFALSE;
}
return kTRUE;
}
//______________________________________________________________________________
template<class Element>
void TVectorT<Element>::Streamer(TBuffer &R__b)
{
// Stream an object of class TVectorT.
if (R__b.IsReading()) {
UInt_t R__s, R__c;
Version_t R__v = R__b.ReadVersion(&R__s,&R__c);
if (R__v > 1) {
Clear();
TVectorT<Element>::Class()->ReadBuffer(R__b,this,R__v,R__s,R__c);
} else { //====process old versions before automatic schema evolution
TObject::Streamer(R__b);
R__b >> fRowLwb;
fNrows = R__b.ReadArray(fElements);
R__b.CheckByteCount(R__s, R__c, TVectorT<Element>::IsA());
}
if (fNrows > 0 && fNrows <= kSizeMax) {
memcpy(fDataStack,fElements,fNrows*sizeof(Element));
delete [] fElements;
fElements = fDataStack;
} else if (fNrows < 0)
Invalidate();
if (R__v < 3)
MakeValid();
} else {
TVectorT<Element>::Class()->WriteBuffer(R__b,this);
}
}
#ifndef ROOT_TMatrixFfwd
#include "TMatrixFfwd.h"
#endif
#ifndef ROOT_TMatrixFSymfwd
#include "TMatrixFSymfwd.h"
#endif
#ifndef ROOT_TMatrixFSparsefwd
#include "TMatrixFSparsefwd.h"
#endif
template class TVectorT<Float_t>;
template Bool_t operator== <Float_t>(const TVectorF &source1,const TVectorF &source2);
template TVectorF operator+ <Float_t>(const TVectorF &source1,const TVectorF &source2);
template TVectorF operator- <Float_t>(const TVectorF &source1,const TVectorF &source2);
template Float_t operator* <Float_t>(const TVectorF &source1,const TVectorF &source2);
template TVectorF operator* <Float_t>(const TMatrixF &a, const TVectorF &source);
template TVectorF operator* <Float_t>(const TMatrixFSym &a, const TVectorF &source);
template TVectorF operator* <Float_t>(const TMatrixFSparse &a, const TVectorF &source);
template TVectorF operator* <Float_t>( Float_t val, const TVectorF &source);
template TVectorF &Add <Float_t>( TVectorF &target, Float_t scalar, const TVectorF &source);
template TVectorF &Add <Float_t>( TVectorF &target, Float_t scalar, const TMatrixF &a,
const TVectorF &source);
template TVectorF &Add <Float_t>( TVectorF &target, Float_t scalar, const TMatrixFSym &a,
const TVectorF &source);
template TVectorF &Add <Float_t>( TVectorF &target, Float_t scalar, const TMatrixFSparse &a,
const TVectorF &source);
template TVectorF &AddElemMult <Float_t>( TVectorF &target, Float_t scalar, const TVectorF &source1,
const TVectorF &source2);
template TVectorF &AddElemMult <Float_t>( TVectorF &target, Float_t scalar, const TVectorF &source1,
const TVectorF &source2,const TVectorF &select);
template TVectorF &AddElemDiv <Float_t>( TVectorF &target, Float_t scalar, const TVectorF &source1,
const TVectorF &source2);
template TVectorF &AddElemDiv <Float_t>( TVectorF &target, Float_t scalar, const TVectorF &source1,
const TVectorF &source2,const TVectorF &select);
template TVectorF &ElementMult <Float_t>( TVectorF &target, const TVectorF &source);
template TVectorF &ElementMult <Float_t>( TVectorF &target, const TVectorF &source, const TVectorF &select);
template TVectorF &ElementDiv <Float_t>( TVectorF &target, const TVectorF &source);
template TVectorF &ElementDiv <Float_t>( TVectorF &target, const TVectorF &source, const TVectorF &select);
template Bool_t AreCompatible <Float_t,Float_t>
(const TVectorF &source1,const TVectorF &source2, Int_t verbose);
template Bool_t AreCompatible <Float_t,Double_t>
(const TVectorF &source1,const TVectorD &source2, Int_t verbose);
template void Compare <Float_t>(const TVectorF &source1,const TVectorF &source2);
template Bool_t VerifyVectorValue <Float_t>(const TVectorF &m, Float_t val, Int_t verbose,
Float_t maxDevAllow);
template Bool_t VerifyVectorIdentity<Float_t>(const TVectorF &m1, const TVectorF &m2, Int_t verbose,
Float_t maxDevAllow);
#ifndef ROOT_TMatrixDfwd
#include "TMatrixDfwd.h"
#endif
#ifndef ROOT_TMatrixDSymfwd
#include "TMatrixDSymfwd.h"
#endif
#ifndef ROOT_TMatrixDSparsefwd
#include "TMatrixDSparsefwd.h"
#endif
template class TVectorT<Double_t>;
template Bool_t operator== <Double_t>(const TVectorD &source1,const TVectorD &source2);
template TVectorD operator+ <Double_t>(const TVectorD &source1,const TVectorD &source2);
template TVectorD operator- <Double_t>(const TVectorD &source1,const TVectorD &source2);
template Double_t operator* <Double_t>(const TVectorD &source1,const TVectorD &source2);
template TVectorD operator* <Double_t>(const TMatrixD &a, const TVectorD &source);
template TVectorD operator* <Double_t>(const TMatrixDSym &a, const TVectorD &source);
template TVectorD operator* <Double_t>(const TMatrixDSparse &a, const TVectorD &source);
template TVectorD operator* <Double_t>( Double_t val, const TVectorD &source);
template TVectorD &Add <Double_t>( TVectorD &target, Double_t scalar, const TVectorD &source);
template TVectorD &Add <Double_t>( TVectorD &target, Double_t scalar, const TMatrixD &a,
const TVectorD &source);
template TVectorD &Add <Double_t>( TVectorD &target, Double_t scalar, const TMatrixDSym &a
, const TVectorD &source);
template TVectorD &Add <Double_t>( TVectorD &target, Double_t scalar, const TMatrixDSparse &a
, const TVectorD &source);
template TVectorD &AddElemMult <Double_t>( TVectorD &target, Double_t scalar, const TVectorD &source1,
const TVectorD &source2);
template TVectorD &AddElemMult <Double_t>( TVectorD &target, Double_t scalar, const TVectorD &source1,
const TVectorD &source2,const TVectorD &select);
template TVectorD &AddElemDiv <Double_t>( TVectorD &target, Double_t scalar, const TVectorD &source1,
const TVectorD &source2);
template TVectorD &AddElemDiv <Double_t>( TVectorD &target, Double_t scalar, const TVectorD &source1,
const TVectorD &source2,const TVectorD &select);
template TVectorD &ElementMult <Double_t>( TVectorD &target, const TVectorD &source);
template TVectorD &ElementMult <Double_t>( TVectorD &target, const TVectorD &source, const TVectorD &select);
template TVectorD &ElementDiv <Double_t>( TVectorD &target, const TVectorD &source);
template TVectorD &ElementDiv <Double_t>( TVectorD &target, const TVectorD &source, const TVectorD &select);
template Bool_t AreCompatible<Double_t,Double_t>
(const TVectorD &source1,const TVectorD &source2, Int_t verbose);
template Bool_t AreCompatible<Double_t,Float_t>
(const TVectorD &source1,const TVectorF &source2, Int_t verbose);
template void Compare <Double_t>(const TVectorD &source1,const TVectorD &source2);
template Bool_t VerifyVectorValue <Double_t>(const TVectorD &m, Double_t val, Int_t verbose,
Double_t maxDevAllow);
template Bool_t VerifyVectorIdentity<Double_t>(const TVectorD &m1, const TVectorD &m2, Int_t verbose,
Double_t maxDevAllow);
ROOT page - Class index - Class Hierarchy - Top of the page
This page has been automatically generated. If you have any comments or suggestions about the page layout send a mail to ROOT support, or contact the developers with any questions or problems regarding ROOT.