|
EpetraExt
Development
|
00001 //@HEADER 00002 // *********************************************************************** 00003 // 00004 // EpetraExt: Epetra Extended - Linear Algebra Services Package 00005 // Copyright (2011) Sandia Corporation 00006 // 00007 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation, 00008 // the U.S. Government retains certain rights in this software. 00009 // 00010 // Redistribution and use in source and binary forms, with or without 00011 // modification, are permitted provided that the following conditions are 00012 // met: 00013 // 00014 // 1. Redistributions of source code must retain the above copyright 00015 // notice, this list of conditions and the following disclaimer. 00016 // 00017 // 2. Redistributions in binary form must reproduce the above copyright 00018 // notice, this list of conditions and the following disclaimer in the 00019 // documentation and/or other materials provided with the distribution. 00020 // 00021 // 3. Neither the name of the Corporation nor the names of the 00022 // contributors may be used to endorse or promote products derived from 00023 // this software without specific prior written permission. 00024 // 00025 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY 00026 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE 00027 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR 00028 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE 00029 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, 00030 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, 00031 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR 00032 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF 00033 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING 00034 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS 00035 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE. 00036 // 00037 // Questions? Contact Michael A. Heroux (maherou@sandia.gov) 00038 // 00039 // *********************************************************************** 00040 //@HEADER 00041 00042 #include "EpetraExt_BlockVector.h" 00043 #include "EpetraExt_BlockUtility.h" 00044 #include "Epetra_Map.h" 00045 #include "Epetra_Comm.h" 00046 00047 namespace EpetraExt { 00048 00049 //============================================================================= 00050 // EpetraExt::BlockVector Constructor 00051 BlockVector::BlockVector( 00052 const Epetra_BlockMap & BaseMap, 00053 const Epetra_BlockMap & GlobalMap) 00054 : Epetra_Vector( GlobalMap ), 00055 BaseMap_( BaseMap ), 00056 Offset_( BlockUtility::CalculateOffset( BaseMap ) ) 00057 { 00058 } 00059 00060 //============================================================================= 00061 // EpetraExt::BlockVector Constructor 00062 BlockVector::BlockVector( 00063 Epetra_DataAccess CV, 00064 const Epetra_BlockMap & BaseMap, 00065 const Epetra_Vector & BlockVec) 00066 : Epetra_Vector( CV, BlockVec, 0 ), 00067 BaseMap_( BaseMap ), 00068 Offset_( BlockUtility::CalculateOffset( BaseMap ) ) 00069 { 00070 } 00071 00072 //========================================================================== 00073 // Copy Constructor 00074 BlockVector::BlockVector(const BlockVector& Source) 00075 : Epetra_Vector( dynamic_cast<const Epetra_Vector &>(Source) ), 00076 BaseMap_( Source.BaseMap_ ), 00077 Offset_( Source.Offset_ ) 00078 { 00079 } 00080 00081 //========================================================================= 00082 BlockVector::~BlockVector() 00083 { 00084 } 00085 00086 //========================================================================= 00087 int BlockVector::ExtractBlockValues(Epetra_Vector & BaseVector, int GlobalBlockRow) const 00088 { 00089 int IndexOffset = GlobalBlockRow * Offset_; 00090 int localIndex=0; 00091 00092 // For each entry in the base vector, translate its global ID 00093 // by the IndexOffset and extract the value from this blockVector 00094 for (int i=0; i<BaseMap_.NumMyElements(); i++) { 00095 localIndex = this->Map().LID((IndexOffset + BaseMap_.GID(i))); 00096 if (localIndex==-1) { 00097 cout << "Error in BlockVector::GetBlock: " << i << " " 00098 << IndexOffset << " " << BaseMap_.GID(i) << endl; 00099 return -1; 00100 } 00101 BaseVector[i] = Values_[localIndex]; 00102 } 00103 00104 return 0; 00105 } 00106 00107 //========================================================================= 00108 int BlockVector::LoadBlockValues(const Epetra_Vector & BaseVector, int GlobalBlockRow) 00109 { 00110 int IndexOffset = GlobalBlockRow * Offset_; 00111 int localIndex=0; 00112 00113 // For each entry in the base vector, translate its global ID 00114 // by the IndexOffset and load into this blockVector 00115 for (int i=0; i<BaseMap_.NumMyElements(); i++) { 00116 localIndex = this->Map().LID((IndexOffset + BaseMap_.GID(i))); 00117 if (localIndex==-1) { 00118 cout << "Error in BlockVector::GetBlock: " << i << " " 00119 << IndexOffset << " " << BaseMap_.GID(i) << endl; 00120 return -1; 00121 } 00122 (*this)[localIndex] = BaseVector[i]; 00123 } 00124 00125 return 0; 00126 } 00127 //========================================================================= 00128 int BlockVector::BlockSumIntoGlobalValues(int NumIndices, double* Values, 00129 int* Indices, int GlobalBlockRow) 00130 { 00131 int IndexOffset = GlobalBlockRow * Offset_; 00132 int localIndex=0; 00133 00134 // For each entry in the base vector, translate its global ID 00135 // by the IndexOffset and load into this blockVector 00136 for (int i=0; i<NumIndices; i++) { 00137 localIndex = this->Map().LID((IndexOffset + Indices[i])); 00138 if (localIndex==-1) { 00139 cout << "Error in BlockVector::BlockSumIntoGlobalValues: " << i 00140 << " " << IndexOffset << " " << Indices[i] << endl; 00141 return -1; 00142 } 00143 (*this)[localIndex] += Values[i]; 00144 } 00145 00146 return 0; 00147 } 00148 //========================================================================= 00149 int BlockVector::BlockReplaceGlobalValues(int NumIndices, double* Values, 00150 int* Indices, int GlobalBlockRow) 00151 { 00152 int IndexOffset = GlobalBlockRow * Offset_; 00153 int localIndex=0; 00154 00155 // For each entry in the base vector, translate its global ID 00156 // by the IndexOffset and load into this blockVector 00157 for (int i=0; i<NumIndices; i++) { 00158 localIndex = this->Map().LID((IndexOffset + Indices[i])); 00159 if (localIndex==-1) { 00160 cout << "Error in BlockVector::BlockReplaceGlobalValues: " << i 00161 << " " << IndexOffset << " " << Indices[i] << endl; 00162 return -1; 00163 } 00164 (*this)[localIndex] = Values[i]; 00165 } 00166 00167 return 0; 00168 } 00169 00170 //========================================================================= 00171 Teuchos::RCP<const Epetra_Vector> 00172 BlockVector::GetBlock(int GlobalBlockRow) const 00173 { 00174 int offset = GlobalBlockRow * BaseMap_.NumMyElements(); 00175 return Teuchos::rcp(new Epetra_Vector(View, BaseMap_, Values_+offset)); 00176 } 00177 00178 //========================================================================= 00179 Teuchos::RCP<Epetra_Vector> 00180 BlockVector::GetBlock(int GlobalBlockRow) 00181 { 00182 int offset = GlobalBlockRow * BaseMap_.NumMyElements(); 00183 return Teuchos::rcp(new Epetra_Vector(View, BaseMap_, Values_+offset)); 00184 } 00185 00186 //========================================================================= 00187 const Epetra_BlockMap& 00188 BlockVector::GetBaseMap() const 00189 { 00190 return BaseMap_; 00191 } 00192 00193 00194 } //namespace EpetraExt
1.7.6.1