DenseLinAlgPack_PermVecMat.cpp

Go to the documentation of this file.
00001 // @HEADER
00002 // ***********************************************************************
00003 // 
00004 // Moocho: Multi-functional Object-Oriented arCHitecture for Optimization
00005 //                  Copyright (2003) Sandia Corporation
00006 // 
00007 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
00008 // license for use of this work by or on behalf of the U.S. Government.
00009 // 
00010 // This library is free software; you can redistribute it and/or modify
00011 // it under the terms of the GNU Lesser General Public License as
00012 // published by the Free Software Foundation; either version 2.1 of the
00013 // License, or (at your option) any later version.
00014 //  
00015 // This library is distributed in the hope that it will be useful, but
00016 // WITHOUT ANY WARRANTY; without even the implied warranty of
00017 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00018 // Lesser General Public License for more details.
00019 //  
00020 // You should have received a copy of the GNU Lesser General Public
00021 // License along with this library; if not, write to the Free Software
00022 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307
00023 // USA
00024 // Questions? Contact Roscoe A. Bartlett (rabartl@sandia.gov) 
00025 // 
00026 // ***********************************************************************
00027 // @HEADER
00028 
00029 #include "DenseLinAlgPack_PermVecMat.hpp"
00030 #include "DenseLinAlgPack_DMatrixClass.hpp"
00031 #include "DenseLinAlgPack_IVector.hpp"
00032 
00033 #ifdef TEUCHOS_DEBUG   // Debug only!
00034 bool DenseLinAlgPack::PermVecMat_print = false;
00035 #include <iostream>
00036 #include "DenseLinAlgPack_PermOut.hpp"
00037 #include "DenseLinAlgPack_DVectorOut.hpp"
00038 #endif
00039 
00040 // Local assert function
00041 namespace {
00042 inline void i_assert_perm_size(size_t size1, size_t size2)
00043 {
00044 #ifdef LINALGPACK_CHECK_RHS_SIZES
00045   if(size1 != size2)
00046     throw std::length_error("The size of the permutation vector is not correct");
00047 #endif
00048 }
00049 
00050 } // end namespace
00051 
00052 void DenseLinAlgPack::identity_perm(IVector* perm) {
00053   if(!perm->size())
00054     throw std::length_error("DenseLinAlgPack::identity_perm(): perm must be sized");
00055   IVector::iterator itr_perm    = perm->begin();
00056   for(size_type i = 1; i <= perm->size(); ++i)
00057     *itr_perm++ = i;
00058 }
00059 
00060 void DenseLinAlgPack::inv_perm(const IVector& perm, IVector* inv_perm) {
00061   inv_perm->resize(perm.size());
00062   for(size_type i = 1; i <= perm.size(); ++i)
00063     (*inv_perm)(perm(i)) = i;
00064 }
00065 
00066 void DenseLinAlgPack::perm_ele(const IVector& perm, DVectorSlice* vs)
00067 {
00068   i_assert_perm_size(vs->dim(),perm.size());
00069   DVector tmp_v(vs->dim());
00070   DVector::iterator   v_itr   = tmp_v.begin(),
00071               v_itr_end = tmp_v.end();
00072   IVector::const_iterator perm_itr = perm.begin();
00073   // Copy the elements in the permed order into the temp vector
00074   for(; v_itr != v_itr_end; ++v_itr, ++perm_itr)
00075     *v_itr = (*vs)(*perm_itr);
00076     
00077   // Copy the permed vector back
00078   (*vs) = tmp_v();
00079 }
00080 
00081 void DenseLinAlgPack::perm_ele(const DVectorSlice& x, const IVector& perm, DVectorSlice* y)
00082 {
00083   i_assert_perm_size(x.dim(),perm.size());
00084   i_assert_perm_size(y->dim(),perm.size());
00085 
00086   IVector::const_iterator
00087     perm_itr  = perm.begin();
00088   DVectorSlice::iterator
00089     y_itr   = y->begin(),
00090     y_end   = y->end();
00091   while(y_itr != y_end)
00092     *y_itr++ = x(*perm_itr++);
00093 }
00094 
00095 void DenseLinAlgPack::inv_perm_ele(const DVectorSlice& y, const IVector& perm, DVectorSlice* x)
00096 {
00097   i_assert_perm_size(y.dim(),perm.size());
00098   i_assert_perm_size(x->dim(),perm.size());
00099 #ifdef TEUCHOS_DEBUG
00100   if( PermVecMat_print ) {
00101     std::cerr
00102       << "enter inv_perm_ele(y,perm,x):\n"
00103       << "before:\n"
00104       << "y =\n" << y
00105       << "perm =\n" << perm
00106       << "x =\n" << *x;
00107   }
00108 #endif  
00109   DVectorSlice::const_iterator
00110     y_itr   = y.begin(),
00111     y_end   = y.end();
00112   IVector::const_iterator
00113     perm_itr  = perm.begin();
00114   while(y_itr != y_end)
00115     (*x)(*perm_itr++) = *y_itr++;
00116 #ifdef TEUCHOS_DEBUG
00117   if( PermVecMat_print ) {
00118     std::cerr
00119       << "inv_perm_ele(y,perm,x):\n"
00120       << "after:\n"
00121       << "x =\n" << *x
00122       << "exit inv_perm_ele(...) ...\n";
00123   }
00124 #endif  
00125 }
00126 
00127 void DenseLinAlgPack::perm_rows(const IVector& row_perm, DMatrixSlice* gms)
00128 {
00129   i_assert_perm_size(gms->rows(),row_perm.size());
00130   DMatrix tmp_gm(gms->rows(),gms->cols());
00131   DMatrixSlice::size_type rows = gms->rows(), i;
00132   // Copy the rows in the correct order into the temp matrix.
00133   for(i = 1; i <= rows; ++i)
00134     tmp_gm.row(i) = gms->row(row_perm(i));
00135   // Copy the permed matrix back
00136   (*gms) = tmp_gm();
00137 }
00138 
00139 void DenseLinAlgPack::perm_cols(const IVector& col_perm, DMatrixSlice* gms)
00140 {
00141   i_assert_perm_size(gms->cols(),col_perm.size());
00142   DMatrix tmp_gm(gms->rows(),gms->cols());
00143   DMatrixSlice::size_type cols = gms->cols(), i;
00144   // Copy the columns in the correct order into the temp matrix.
00145   for(i = 1; i <= cols; ++i)
00146     tmp_gm.col(i) = gms->col(col_perm(i));
00147   // Copy the permed matrix back
00148   (*gms) = tmp_gm();
00149 }
00150 
00151 void DenseLinAlgPack::perm_rows_cols(const IVector& row_perm, const IVector& col_perm
00152   , DMatrixSlice* gms)
00153 {
00154   i_assert_perm_size(gms->rows(),row_perm.size());
00155   i_assert_perm_size(gms->cols(),col_perm.size());
00156   DMatrix tmp_gm(gms->rows(),gms->cols());
00157   DMatrixSlice::size_type rows = gms->rows(), cols = gms->cols(), i;
00158   // Copy the rows in the correct order into the temp matrix.
00159   for(i = 1; i <= rows; ++i)
00160     tmp_gm.row(i) = gms->row(row_perm(i));
00161   // Copy the columns in the correct order back into matrix.
00162   for(i = 1; i <= cols; ++i)
00163     gms->col(i) = tmp_gm.col(col_perm(i));
00164 }

Generated on Tue Oct 20 12:51:46 2009 for MOOCHO (Single Doxygen Collection) by doxygen 1.4.7