EpetraExt Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
EpetraExt_BlockVector.cpp
Go to the documentation of this file.
1 //@HEADER
2 // ***********************************************************************
3 //
4 // EpetraExt: Epetra Extended - Linear Algebra Services Package
5 // Copyright (2011) Sandia Corporation
6 //
7 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
8 // the U.S. Government retains certain rights in this software.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
38 //
39 // ***********************************************************************
40 //@HEADER
41 
42 #include "EpetraExt_BlockVector.h"
43 #include "EpetraExt_BlockUtility.h"
44 #include "Epetra_Map.h"
45 #include "Epetra_Comm.h"
46 
47 namespace EpetraExt {
48 
49 //=============================================================================
50 // EpetraExt::BlockVector Constructor
52  const Epetra_BlockMap & BaseMap,
53  const Epetra_BlockMap & GlobalMap)
54  : Epetra_Vector( GlobalMap ),
55  BaseMap_( BaseMap ),
56  Offset_( BlockUtility::CalculateOffset64( BaseMap ) )
57 {
58 }
59 
60 //=============================================================================
61 // EpetraExt::BlockVector Constructor
64  const Epetra_BlockMap & BaseMap,
65  const Epetra_Vector & BlockVec)
66  : Epetra_Vector( CV, BlockVec, 0 ),
67  BaseMap_( BaseMap ),
68  Offset_( BlockUtility::CalculateOffset64( BaseMap ) )
69 {
70 }
71 
72 //==========================================================================
73 // Copy Constructor
75  : Epetra_Vector( dynamic_cast<const Epetra_Vector &>(Source) ),
76  BaseMap_( Source.BaseMap_ ),
77  Offset_( Source.Offset_ )
78 {
79 }
80 
81 //=========================================================================
83 {
84 }
85 
86 //=========================================================================
87 int BlockVector::ExtractBlockValues(Epetra_Vector & BaseVector, long long GlobalBlockRow) const
88 {
89  long long IndexOffset = GlobalBlockRow * Offset_;
90  int localIndex=0;
91 
92  // For each entry in the base vector, translate its global ID
93  // by the IndexOffset and extract the value from this blockVector
94  for (int i=0; i<BaseMap_.NumMyElements(); i++) {
95  localIndex = this->Map().LID((IndexOffset + BaseMap_.GID64(i)));
96  if (localIndex==-1) {
97  std::cout << "Error in BlockVector::GetBlock: " << i << " "
98  << IndexOffset << " " << BaseMap_.GID64(i) << std::endl;
99  return -1;
100  }
101  BaseVector[i] = Values_[localIndex];
102  }
103 
104  return 0;
105 }
106 
107 //=========================================================================
108 int BlockVector::LoadBlockValues(const Epetra_Vector & BaseVector, long long GlobalBlockRow)
109 {
110  long long IndexOffset = GlobalBlockRow * Offset_;
111  int localIndex=0;
112 
113  // For each entry in the base vector, translate its global ID
114  // by the IndexOffset and load into this blockVector
115  for (int i=0; i<BaseMap_.NumMyElements(); i++) {
116  localIndex = this->Map().LID((IndexOffset + BaseMap_.GID64(i)));
117  if (localIndex==-1) {
118  std::cout << "Error in BlockVector::GetBlock: " << i << " "
119  << IndexOffset << " " << BaseMap_.GID64(i) << std::endl;
120  return -1;
121  }
122  (*this)[localIndex] = BaseVector[i];
123  }
124 
125  return 0;
126 }
127 //=========================================================================
128 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
129 int BlockVector::BlockSumIntoGlobalValues(int NumIndices, double* theValues,
130  int* Indices, int GlobalBlockRow)
131 {
132  long long IndexOffset = GlobalBlockRow * Offset_;
133  int localIndex=0;
134 
135  // For each entry in the base vector, translate its global ID
136  // by the IndexOffset and load into this blockVector
137  for (int i=0; i<NumIndices; i++) {
138  localIndex = this->Map().LID((IndexOffset + Indices[i]));
139  if (localIndex==-1) {
140  std::cout << "Error in BlockVector::BlockSumIntoGlobalValues: " << i
141  << " " << IndexOffset << " " << Indices[i] << std::endl;
142  return -1;
143  }
144  (*this)[localIndex] += theValues[i];
145  }
146 
147  return 0;
148 }
149 //=========================================================================
150 int BlockVector::BlockReplaceGlobalValues(int NumIndices, double* theValues,
151  int* Indices, int GlobalBlockRow)
152 {
153  long long IndexOffset = GlobalBlockRow * Offset_;
154  int localIndex=0;
155 
156  // For each entry in the base vector, translate its global ID
157  // by the IndexOffset and load into this blockVector
158  for (int i=0; i<NumIndices; i++) {
159  localIndex = this->Map().LID((IndexOffset + Indices[i]));
160  if (localIndex==-1) {
161  std::cout << "Error in BlockVector::BlockReplaceGlobalValues: " << i
162  << " " << IndexOffset << " " << Indices[i] << std::endl;
163  return -1;
164  }
165  (*this)[localIndex] = theValues[i];
166  }
167 
168  return 0;
169 }
170 #endif
171 //=========================================================================
172 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
173 int BlockVector::BlockSumIntoGlobalValues(int NumIndices, double* theValues,
174  long long* Indices, long long GlobalBlockRow)
175 {
176  long long IndexOffset = GlobalBlockRow * Offset_;
177  int localIndex=0;
178 
179  // For each entry in the base vector, translate its global ID
180  // by the IndexOffset and load into this blockVector
181  for (int i=0; i<NumIndices; i++) {
182  localIndex = this->Map().LID((IndexOffset + Indices[i]));
183  if (localIndex==-1) {
184  std::cout << "Error in BlockVector::BlockSumIntoGlobalValues: " << i
185  << " " << IndexOffset << " " << Indices[i] << std::endl;
186  return -1;
187  }
188  (*this)[localIndex] += theValues[i];
189  }
190 
191  return 0;
192 }
193 //=========================================================================
194 int BlockVector::BlockReplaceGlobalValues(int NumIndices, double* theValues,
195  long long* Indices, long long GlobalBlockRow)
196 {
197  long long IndexOffset = GlobalBlockRow * Offset_;
198  int localIndex=0;
199 
200  // For each entry in the base vector, translate its global ID
201  // by the IndexOffset and load into this blockVector
202  for (int i=0; i<NumIndices; i++) {
203  localIndex = this->Map().LID((IndexOffset + Indices[i]));
204  if (localIndex==-1) {
205  std::cout << "Error in BlockVector::BlockReplaceGlobalValues: " << i
206  << " " << IndexOffset << " " << Indices[i] << std::endl;
207  return -1;
208  }
209  (*this)[localIndex] = theValues[i];
210  }
211 
212  return 0;
213 }
214 #endif
215 //=========================================================================
217 BlockVector::GetBlock(long long GlobalBlockRow) const
218 {
219  long long offset = GlobalBlockRow * BaseMap_.NumMyElements();
220  return Teuchos::rcp(new Epetra_Vector(View, BaseMap_, Values_+offset));
221 }
222 
223 //=========================================================================
225 BlockVector::GetBlock(long long GlobalBlockRow)
226 {
227  long long offset = GlobalBlockRow * BaseMap_.NumMyElements();
228  return Teuchos::rcp(new Epetra_Vector(View, BaseMap_, Values_+offset));
229 }
230 
231 //=========================================================================
232 const Epetra_BlockMap&
234 {
235  return BaseMap_;
236 }
237 
238 
239 } //namespace EpetraExt
int BlockSumIntoGlobalValues(int NumIndices, double *Values, int *Indices, int BlockRow)
Load entries into BlockVector with base vector indices offset by BlockRow.
int BlockReplaceGlobalValues(int NumIndices, double *Values, int *Indices, int BlockRow)
Load entries into BlockVector with base vector indices offset by BlockRow.
virtual ~BlockVector()
Destructor.
int NumMyElements() const
long long GID64(int LID) const
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
virtual const Epetra_BlockMap & Map() const =0
int ExtractBlockValues(Epetra_Vector &BaseVec, long long BlockRow) const
Extract a single block from a Block Vector: block row is global, not a stencil value.
int LID(int GID) const
Epetra_Vector(const Epetra_BlockMap &Map, bool zeroOut=true)
int LoadBlockValues(const Epetra_Vector &BaseVec, long long BlockRow)
Load a single block into a Block Vector: block row is global, not a stencil value.
BlockVector(const Epetra_BlockMap &BaseMap, const Epetra_BlockMap &GlobalMap)
BlockVector constuctor with one block row per processor.
Teuchos::RCP< const Epetra_Vector > GetBlock(long long BlockRow) const
Return Epetra_Vector for given block row.
Epetra_DataAccess
const Epetra_BlockMap & GetBaseMap() const
Return base map.