|
EpetraExt
Development
|
00001 /* 00002 //@HEADER 00003 // *********************************************************************** 00004 // 00005 // EpetraExt: Epetra Extended - Linear Algebra Services Package 00006 // Copyright (2011) Sandia Corporation 00007 // 00008 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation, 00009 // the U.S. Government retains certain rights in this software. 00010 // 00011 // Redistribution and use in source and binary forms, with or without 00012 // modification, are permitted provided that the following conditions are 00013 // met: 00014 // 00015 // 1. Redistributions of source code must retain the above copyright 00016 // notice, this list of conditions and the following disclaimer. 00017 // 00018 // 2. Redistributions in binary form must reproduce the above copyright 00019 // notice, this list of conditions and the following disclaimer in the 00020 // documentation and/or other materials provided with the distribution. 00021 // 00022 // 3. Neither the name of the Corporation nor the names of the 00023 // contributors may be used to endorse or promote products derived from 00024 // this software without specific prior written permission. 00025 // 00026 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY 00027 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE 00028 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR 00029 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE 00030 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, 00031 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, 00032 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR 00033 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF 00034 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING 00035 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS 00036 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE. 00037 // 00038 // Questions? Contact Michael A. Heroux (maherou@sandia.gov) 00039 // 00040 // *********************************************************************** 00041 //@HEADER 00042 */ 00043 00044 #include "EpetraExt_PointToBlockDiagPermute.h" 00045 #include "EpetraExt_BlockDiagMatrix.h" 00046 #include "Epetra_MultiVector.h" 00047 #include "Epetra_Import.h" 00048 #include "Epetra_Export.h" 00049 #include "Epetra_Comm.h" 00050 00051 #include <stdio.h> 00052 #include <fstream> 00053 00054 #define MAX(x,y) ((x)>(y)?(x):(y)) 00055 00056 //========================================================================= 00057 EpetraExt_PointToBlockDiagPermute::EpetraExt_PointToBlockDiagPermute(const Epetra_CrsMatrix& MAT) 00058 :Epetra_DistObject(MAT.RowMap()), 00059 Matrix_(&MAT), 00060 PurelyLocalMode_(true), 00061 ContiguousBlockMode_(false), 00062 ContiguousBlockSize_(0), 00063 NumBlocks_(0), 00064 Blockstart_(0), 00065 Blockids_(0), 00066 BDMap_(0), 00067 CompatibleMap_(0), 00068 BDMat_(0), 00069 Importer_(0), 00070 Exporter_(0), 00071 ImportVector_(0), 00072 ExportVector_(0) 00073 { 00074 00075 } 00076 00077 //========================================================================= 00078 // Destructor 00079 EpetraExt_PointToBlockDiagPermute::~EpetraExt_PointToBlockDiagPermute() 00080 { 00081 if(BDMap_) delete BDMap_; 00082 if(CompatibleMap_) delete CompatibleMap_; 00083 if(BDMat_) delete BDMat_; 00084 if(Importer_) delete Importer_; 00085 if(Exporter_) delete Exporter_; 00086 if(ImportVector_) delete ImportVector_; 00087 if(ExportVector_) delete ExportVector_; 00088 } 00089 00090 00091 00092 //========================================================================= 00093 // Set list 00094 int EpetraExt_PointToBlockDiagPermute::SetParameters(Teuchos::ParameterList & List){ 00095 List_=List; 00096 00097 // Check for contiguous blocking first 00098 ContiguousBlockSize_=List_.get("contiguous block size",0); 00099 if(ContiguousBlockSize_!=0){ 00100 ContiguousBlockMode_=true; 00101 PurelyLocalMode_=false; 00102 } 00103 00104 // Local vs. global ids & mode 00105 NumBlocks_=List_.get("number of local blocks",0); 00106 Blockstart_=List_.get("block start index",(int*)0); 00107 Blockids_=List_.get("block entry lids",(int*)0); 00108 if(Blockids_) 00109 PurelyLocalMode_=true; 00110 else{ 00111 Blockids_=List_.get("block entry gids",(int*)0); 00112 PurelyLocalMode_=false; 00113 } 00114 00115 // Sanity checks 00116 if(ContiguousBlockMode_){ 00117 // Can't use contiguous at the same time as the other modes 00118 if(NumBlocks_ || Blockstart_ || Blockids_) EPETRA_CHK_ERR(-4); 00119 } 00120 else { 00121 if(NumBlocks_ <= 0) EPETRA_CHK_ERR(-1); 00122 if(!Blockstart_) EPETRA_CHK_ERR(-2); 00123 if(!Blockids_) EPETRA_CHK_ERR(-3); 00124 } 00125 00126 return 0; 00127 } 00128 00129 00130 //========================================================================= 00131 // Extracts the block-diagonal, builds maps, etc. 00132 int EpetraExt_PointToBlockDiagPermute::Compute(){ 00133 int rv=ExtractBlockDiagonal(); 00134 return rv; 00135 } 00136 00137 //========================================================================= 00138 // Returns the result of a Epetra_Operator applied to a Epetra_MultiVector X in Y. 00139 int EpetraExt_PointToBlockDiagPermute::Apply(const Epetra_MultiVector& X, Epetra_MultiVector& Y) const{ 00140 return -1; 00141 00142 } 00143 00144 //========================================================================= 00145 // Returns the result of a Epetra_Operator inverse applied to an Epetra_MultiVector X in Y. 00146 int EpetraExt_PointToBlockDiagPermute::ApplyInverse(const Epetra_MultiVector& X, Epetra_MultiVector& Y) const{ 00147 // Stuff borrowed from Epetra_CrsMatrix 00148 int NumVectors = X.NumVectors(); 00149 if (NumVectors!=Y.NumVectors()) { 00150 EPETRA_CHK_ERR(-2); // Need same number of vectors in each MV 00151 } 00152 00153 const Epetra_MultiVector *Xp=&X; 00154 Epetra_MultiVector *Yp=&Y; 00155 00156 // Allocate temp workspace if X==Y and there are no imports or exports 00157 Epetra_MultiVector * Xcopy = 0; 00158 if (&X==&Y && Importer_==0 && Exporter_==0) { 00159 Xcopy = new Epetra_MultiVector(X); 00160 Xp=Xcopy; 00161 } 00162 00163 UpdateImportVector(NumVectors); // Make sure Import and Export Vectors are compatible 00164 UpdateExportVector(NumVectors); 00165 00166 // If we have a non-trivial importer, we must import elements that are permuted or are on other processors 00167 if (Importer_){ 00168 EPETRA_CHK_ERR(ImportVector_->Import(X, *Importer_, Insert)); 00169 Xp=ImportVector_; 00170 } 00171 00172 // If we have a non-trivial exporter, we must export elements that are permuted or belong to other processors 00173 if (Exporter_) { 00174 Yp=ExportVector_; 00175 } 00176 00177 // Do the matvec 00178 BDMat_->ApplyInverse(*Xp,*Yp); 00179 00180 // Export if needed 00181 if (Exporter_) { 00182 Y.PutScalar(0.0); // Make sure target is zero 00183 Y.Export(*ExportVector_, *Exporter_, Add); // Fill Y with Values from export vector 00184 } 00185 00186 // Cleanup 00187 if(Xcopy) { 00188 delete Xcopy; 00189 EPETRA_CHK_ERR(1); // Return positive code to alert the user about needing extra copy of X 00190 return 1; 00191 } 00192 00193 return 0; 00194 } 00195 00196 //========================================================================= 00197 // Print method 00198 void EpetraExt_PointToBlockDiagPermute::Print(ostream& os) const{ 00199 if(Importer_) cout<<*Importer_<<endl; 00200 if(Exporter_) cout<<*Exporter_<<endl; 00201 if(BDMat_) cout<<*BDMat_<<endl; 00202 } 00203 00204 //========================================================================= 00205 // Pulls the block diagonal of the matrix and then builds the BDMat_ 00206 int EpetraExt_PointToBlockDiagPermute::ExtractBlockDiagonal(){ 00207 int i,j; 00208 std::vector<int> ExtRows; 00209 int* LocalColIDS=0; 00210 int ExtSize=0,ExtCols=0,MainCols=0; 00211 int Nrows=Matrix_->NumMyRows(); 00212 int *l2b,*block_offset,*l2blockid; 00213 const Epetra_Map &RowMap=Matrix_->RowMap(); 00214 int index,col,row_in_block,col_in_block,length,*colind; 00215 double *values; 00216 00217 bool verbose=(bool)(List_.get("output",0) > 0); 00218 00219 // Contiguous Setup 00220 SetupContiguousMode(); 00221 00222 // Compute block size lists 00223 int *bsize=new int[NumBlocks_]; 00224 for(i=0;i<NumBlocks_;i++) 00225 bsize[i]=Blockstart_[i+1]-Blockstart_[i]; 00226 00227 // Use the ScanSum function to compute a prefix sum of the number of points 00228 int MyMinGID; 00229 Matrix_->Comm().ScanSum(&NumBlocks_,&MyMinGID, 1); 00230 MyMinGID-=NumBlocks_; 00231 int *MyBlockGIDs=new int[NumBlocks_]; 00232 for(i=0;i<NumBlocks_;i++) 00233 MyBlockGIDs[i]=MyMinGID+i; 00234 00235 BDMap_=new Epetra_BlockMap(-1,NumBlocks_,MyBlockGIDs,bsize,0,Matrix_->Comm()); 00236 00237 // Allocate the Epetra_DistBlockMatrix 00238 BDMat_=new EpetraExt_BlockDiagMatrix(*BDMap_,true); 00239 00240 // First check to see if we can switch back to PurelyLocalMode_ if we're not 00241 // in it 00242 if(!PurelyLocalMode_){ 00243 // Find the non-local row IDs 00244 ExtSize=MAX(0,Blockstart_[NumBlocks_]-RowMap.NumMyElements());// estimate 00245 ExtRows.reserve(ExtSize); 00246 00247 for(i=0;i<Blockstart_[NumBlocks_];i++){ 00248 if(RowMap.LID(Blockids_[i])==-1) 00249 ExtRows.push_back(Blockids_[i]); 00250 } 00251 00252 // Switch to PurelyLocalMode_ if we never need GIDs 00253 int size_sum; 00254 ExtSize=ExtRows.size(); 00255 RowMap.Comm().SumAll(&ExtSize,&size_sum,1); 00256 if(size_sum==0){ 00257 if(verbose && !Matrix_->Comm().MyPID()) printf("EpetraExt_PointToBlockDiagPermute: Switching to purely local mode\n"); 00258 PurelyLocalMode_=true; 00259 for(i=0;i<Blockstart_[NumBlocks_];i++){ 00260 Blockids_[i]=RowMap.LID(Blockids_[i]); 00261 } 00262 } 00263 } 00264 00265 if(PurelyLocalMode_){ 00266 /*******************************************************/ 00267 // Allocations 00268 l2b=new int[Nrows]; 00269 block_offset=new int[Nrows]; 00270 00271 // Build the local-id-to-block-id list, block offset list 00272 for(i=0;i<NumBlocks_;i++) { 00273 for(j=Blockstart_[i];j<Blockstart_[i+1];j++){ 00274 block_offset[Blockids_[j]]=j-Blockstart_[i]; 00275 l2b[Blockids_[j]]=i; 00276 } 00277 } 00278 00279 // Copy the data to the EpetraExt_BlockDiagMatrix 00280 for(i=0;i<Nrows;i++) { 00281 int block_num = l2b[i]; 00282 if(block_num>=0 && block_num<NumBlocks_) { 00283 row_in_block=block_offset[i]; 00284 Matrix_->ExtractMyRowView(i,length,values,colind); 00285 int Nnz=0; 00286 for (j = 0; j < length; j++) { 00287 col = colind[j]; 00288 if(col < Nrows && l2b[col]==block_num){ 00289 Nnz++; 00290 col_in_block = block_offset[col]; 00291 index = col_in_block * bsize[block_num] + row_in_block; 00292 (*BDMat_)[block_num][index]=values[j]; 00293 } 00294 } 00295 /* Handle the case of a zero row. */ 00296 /* By just putting a 1 on the diagonal. */ 00297 if (Nnz == 0) { 00298 index = row_in_block * bsize[block_num] + row_in_block; 00299 (*BDMat_)[block_num][index]=1.0; 00300 } 00301 } 00302 } 00303 00304 // Build the compatible map for import/export 00305 l2blockid=new int[Nrows]; 00306 for(i=0;i<Nrows;i++) 00307 l2blockid[Blockstart_[l2b[i]]+block_offset[i]]=Matrix_->RowMap().GID(i); 00308 00309 // Build the Compatible Map, Import/Export Objects 00310 CompatibleMap_=new Epetra_Map(-1,Nrows,l2blockid,0,Matrix_->Comm()); 00311 00312 } 00313 else{ 00314 /*******************************************************/ 00315 // Do the import to grab matrix entries 00316 // Allocate temporaries for import 00317 int* ExtRowsPtr = ExtRows.size()>0 ? &ExtRows[0] : NULL; 00318 Epetra_Map TmpMap(-1,ExtSize, ExtRowsPtr,0,Matrix_->Comm());; 00319 Epetra_CrsMatrix TmpMatrix(Copy,TmpMap,0); 00320 Epetra_Import TmpImporter(TmpMap,RowMap); 00321 00322 TmpMatrix.Import(*Matrix_,TmpImporter,Insert); 00323 TmpMatrix.FillComplete(); 00324 00325 ExtCols=TmpMatrix.NumMyCols(); 00326 MainCols=Matrix_->NumMyCols(); 00327 00328 00329 // Build the column reidex - main matrix 00330 LocalColIDS=new int[MainCols+ExtCols]; 00331 for(i=0;i<MainCols;i++){ 00332 int GID=Matrix_->GCID(i); 00333 int MainLID=RowMap.LID(GID); 00334 if(MainLID!=-1) LocalColIDS[i]=MainLID; 00335 else{ 00336 int ExtLID=TmpMatrix.LRID(GID); 00337 if(ExtLID!=-1) LocalColIDS[i]=Nrows+ExtLID; 00338 else LocalColIDS[i]=-1; 00339 } 00340 } 00341 00342 // Build the column reidex - ext matrix 00343 for(i=0;i<ExtCols;i++){ 00344 int GID=TmpMatrix.GCID(i); 00345 int MainLID=RowMap.LID(GID); 00346 if(MainLID!=-1) LocalColIDS[MainCols+i]=MainLID; 00347 else{ 00348 int ExtLID=TmpMatrix.LRID(GID); 00349 if(ExtLID!=-1) LocalColIDS[MainCols+i]=Nrows+ExtLID; 00350 else LocalColIDS[MainCols+i]=-1; 00351 } 00352 } 00353 00354 // Allocations 00355 l2b=new int[Nrows+ExtSize]; 00356 block_offset=new int[Nrows+ExtSize]; 00357 00358 // Build l2b/block_offset with the expanded local index 00359 //NTS: Uses the same ordering of operation as the above routine, which is why 00360 //it works. 00361 for(i=0;i<Nrows+ExtSize;i++) block_offset[i]=l2b[i]=-1; 00362 int ext_idx=0; 00363 for(i=0;i<NumBlocks_;i++) { 00364 for(j=Blockstart_[i];j<Blockstart_[i+1];j++){ 00365 int LID=RowMap.LID(Blockids_[j]); 00366 if(LID==-1) {LID=Nrows+ext_idx;ext_idx++;} 00367 block_offset[LID]=j-Blockstart_[i]; 00368 l2b[LID]=i; 00369 } 00370 } 00371 00372 // Copy the data to the EpetraExt_BlockDiagMatrix from Matrix_ 00373 for(i=0;i<Nrows;i++) { 00374 int block_num = l2b[i]; 00375 if(block_num>=0 && block_num<NumBlocks_) { 00376 row_in_block=block_offset[i]; 00377 Matrix_->ExtractMyRowView(i,length,values,colind); 00378 int Nnz=0; 00379 for (j = 0; j < length; j++) { 00380 col = LocalColIDS[colind[j]]; 00381 if(col!=-1 && l2b[col]==block_num){ 00382 Nnz++; 00383 col_in_block = block_offset[col]; 00384 index = col_in_block * bsize[block_num] + row_in_block; 00385 (*BDMat_)[block_num][index]=values[j]; 00386 } 00387 } 00388 /* Handle the case of a zero row. */ 00389 /* By just putting a 1 on the diagonal. */ 00390 if (Nnz == 0) { 00391 index = row_in_block * bsize[block_num] + row_in_block; 00392 (*BDMat_)[block_num][index]=values[j]; 00393 } 00394 } 00395 } 00396 00397 // Copy the data to the EpetraExt_BlockDiagMatrix from TmpMatrix 00398 for(i=0;i<ExtSize;i++) { 00399 int block_num = l2b[Nrows+i]; 00400 if(block_num>=0 && block_num<NumBlocks_) { 00401 row_in_block=block_offset[Nrows+i]; 00402 TmpMatrix.ExtractMyRowView(i,length,values,colind); 00403 int Nnz=0; 00404 for (j = 0; j < length; j++) { 00405 col = LocalColIDS[MainCols+colind[j]]; 00406 if(col!=-1 && l2b[col]==block_num){ 00407 Nnz++; 00408 col_in_block = block_offset[col]; 00409 index = col_in_block * bsize[block_num] + row_in_block; 00410 (*BDMat_)[block_num][index]=values[j]; 00411 } 00412 } 00413 /* Handle the case of a zero row. */ 00414 /* By just putting a 1 on the diagonal. */ 00415 if (Nnz == 0) { 00416 index = row_in_block * bsize[block_num] + row_in_block; 00417 (*BDMat_)[block_num][index]=1.0; 00418 } 00419 } 00420 } 00421 00422 // Build the compatible map for import/export 00423 l2blockid=new int[Blockstart_[NumBlocks_]]; 00424 for(i=0;i<Nrows;i++){ 00425 int bkid=l2b[i]; 00426 if(bkid>-1) l2blockid[Blockstart_[bkid]+block_offset[i]]=RowMap.GID(i); 00427 } 00428 // NTS: This is easier - if we imported it, we know we need it 00429 for(i=0;i<ExtSize;i++) 00430 l2blockid[Blockstart_[l2b[Nrows+i]]+block_offset[Nrows+i]]=TmpMatrix.GRID(i); 00431 00432 // Build the Compatible Map, Import/Export Objects 00433 CompatibleMap_=new Epetra_Map(-1,Blockstart_[NumBlocks_],l2blockid,0,Matrix_->Comm()); 00434 00435 }//end else 00436 00437 // Set BDMat_'s parameter list and compute! 00438 Teuchos::ParameterList dummy,inList; 00439 inList=List_.get("blockdiagmatrix: list",dummy); 00440 BDMat_->SetParameters(inList); 00441 BDMat_->Compute(); 00442 00443 // Build importer/exporter 00444 if(!CompatibleMap_->SameAs(Matrix_->DomainMap())) 00445 Importer_ = new Epetra_Import(*CompatibleMap_,Matrix_->DomainMap()); 00446 if(!CompatibleMap_->SameAs(Matrix_->RangeMap())) 00447 Exporter_ = new Epetra_Export(*CompatibleMap_,Matrix_->RangeMap()); 00448 00449 // Cleanup 00450 delete [] LocalColIDS; 00451 delete [] block_offset; 00452 delete [] l2b; 00453 delete [] l2blockid; 00454 delete [] bsize; 00455 delete [] MyBlockGIDs; 00456 00457 // Contiguous Cleanup 00458 CleanupContiguousMode(); 00459 00460 return 0; 00461 } 00462 //======================================================================================================= 00463 int EpetraExt_PointToBlockDiagPermute::SetupContiguousMode(){ 00464 if(!ContiguousBlockMode_) return 0; 00465 // NTS: In case of processor-crossing blocks, the lowest PID always gets the block; 00466 const Epetra_Map &RowMap=Matrix_->RowMap(); 00467 00468 int MinMyGID=RowMap.MinMyGID(); 00469 int MaxMyGID=RowMap.MaxMyGID(); 00470 int Base=Matrix_->IndexBase(); 00471 00472 // Find the GID that begins my first block 00473 int MyFirstBlockGID=ContiguousBlockSize_*(int)ceil(((double)(MinMyGID - Base))/ContiguousBlockSize_)+Base; 00474 NumBlocks_=(int)ceil((double)((MaxMyGID-MyFirstBlockGID+1.0)) / ContiguousBlockSize_); 00475 00476 // Allocate memory 00477 Blockstart_=new int[NumBlocks_+1]; 00478 Blockids_=new int[NumBlocks_*ContiguousBlockSize_]; 00479 Blockstart_[NumBlocks_]=NumBlocks_*ContiguousBlockSize_; 00480 00481 // Fill the arrays 00482 for(int i=0,ct=0;i<NumBlocks_;i++){ 00483 Blockstart_[i]=ct; 00484 for(int j=0;j<ContiguousBlockSize_;j++,ct++){ 00485 Blockids_[ct]=MyFirstBlockGID+ct; 00486 } 00487 } 00488 00489 return 0; 00490 } 00491 //======================================================================================================= 00492 int EpetraExt_PointToBlockDiagPermute::CleanupContiguousMode(){ 00493 if(!ContiguousBlockMode_) return 0; 00494 NumBlocks_=0; 00495 if(Blockstart_) {delete [] Blockstart_; Blockstart_=0;} 00496 if(Blockids_) {delete [] Blockids_; Blockids_=0;} 00497 return 0; 00498 } 00499 00500 00501 //======================================================================================================= 00502 // Creates an Epetra_CrsMatrix from the BlockDiagMatrix. This is generally only useful if you want to do a matrix-matrix multiply. 00503 Epetra_FECrsMatrix * EpetraExt_PointToBlockDiagPermute::CreateFECrsMatrix(){ 00504 Epetra_FECrsMatrix * NewMat=new Epetra_FECrsMatrix(Copy,Matrix_->RowMap(),0); 00505 00506 const Epetra_BlockMap &BlockMap=BDMat_->BlockMap(); 00507 const Epetra_BlockMap &DataMap=BDMat_->DataMap(); 00508 const int *vlist=DataMap.FirstPointInElementList(); 00509 const int *xlist=BlockMap.FirstPointInElementList(); 00510 const int *blocksize=BlockMap.ElementSizeList(); 00511 const double *values=BDMat_->Values(); 00512 int NumBlocks=BDMat_->NumMyBlocks(); 00513 00514 // Maximum size vector for stashing GIDs 00515 std::vector<int> GIDs; 00516 GIDs.resize(BlockMap.MaxMyElementSize()); 00517 00518 00519 for(int i=0;i<NumBlocks;i++){ 00520 int Nb=blocksize[i]; 00521 int vidx0=vlist[i]; 00522 int xidx0=xlist[i]; 00523 // Get global indices 00524 for(int j=0;j<Nb;j++) 00525 GIDs[j]=CompatibleMap_->GID(xidx0+j); 00526 00527 // Remember: We're using column-major storage for LAPACK's benefit 00528 int ierr=NewMat->InsertGlobalValues(Nb,&GIDs[0],&values[vidx0],Epetra_FECrsMatrix::COLUMN_MAJOR); 00529 if(ierr < 0) throw "CreateFECrsMatrix: ERROR in InsertGlobalValues"; 00530 } 00531 NewMat->GlobalAssemble(); 00532 return NewMat; 00533 } 00534 00535 00536 //======================================================================================================= 00537 void EpetraExt_PointToBlockDiagPermute::UpdateImportVector(int NumVectors) const { 00538 if(Importer_ != 0) { 00539 if(ImportVector_ != 0) { 00540 if(ImportVector_->NumVectors() != NumVectors) { 00541 delete ImportVector_; 00542 ImportVector_= 0; 00543 } 00544 } 00545 if(ImportVector_ == 0) 00546 ImportVector_ = new Epetra_MultiVector(*CompatibleMap_,NumVectors); // Create import vector if needed 00547 } 00548 return; 00549 } 00550 //======================================================================================================= 00551 void EpetraExt_PointToBlockDiagPermute::UpdateExportVector(int NumVectors) const { 00552 if(Exporter_ != 0) { 00553 if(ExportVector_ != 0) { 00554 if(ExportVector_->NumVectors() != NumVectors) { 00555 delete ExportVector_; 00556 ExportVector_= 0; 00557 } 00558 } 00559 if(ExportVector_ == 0) 00560 ExportVector_ = new Epetra_MultiVector(*CompatibleMap_,NumVectors); // Create Export vector if needed 00561 } 00562 return; 00563 } 00564 00565 00566 00567 00568 //========================================================================= 00569 int EpetraExt_PointToBlockDiagPermute::Import(const Epetra_SrcDistObject& A, const Epetra_Import& Importer, Epetra_CombineMode CombineMode, const Epetra_OffsetIndex *Indexor){ 00570 return -1; 00571 } 00572 00573 //========================================================================= 00574 int EpetraExt_PointToBlockDiagPermute::Import(const Epetra_SrcDistObject& A, const Epetra_Export& Exporter, Epetra_CombineMode CombineMode, const Epetra_OffsetIndex *Indexor){ 00575 return -1; 00576 } 00577 00578 //========================================================================= 00579 int EpetraExt_PointToBlockDiagPermute::Export(const Epetra_SrcDistObject& A, const Epetra_Import & Importer, Epetra_CombineMode CombineMode, const Epetra_OffsetIndex *Indexor){ 00580 return -1; 00581 } 00582 00583 //========================================================================= 00584 int EpetraExt_PointToBlockDiagPermute::Export(const Epetra_SrcDistObject& A, const Epetra_Export& Exporter, Epetra_CombineMode CombineMode, const Epetra_OffsetIndex *Indexor){ 00585 return -1; 00586 } 00587 00588 //========================================================================= 00589 // Allows the source and target (\e this) objects to be compared for compatibility, return nonzero if not. 00590 int EpetraExt_PointToBlockDiagPermute::CheckSizes(const Epetra_SrcDistObject& Source){ 00591 return -1; 00592 } 00593 00594 //========================================================================= 00595 // Perform ID copies and permutations that are on processor. 00596 int EpetraExt_PointToBlockDiagPermute::CopyAndPermute(const Epetra_SrcDistObject& Source, 00597 int NumSameIDs, 00598 int NumPermuteIDs, 00599 int * PermuteToLIDs, 00600 int * PermuteFromLIDs, 00601 const Epetra_OffsetIndex * Indexor){ 00602 return -1; 00603 } 00604 00605 //========================================================================= 00606 // Perform any packing or preparation required for call to DoTransfer(). 00607 int EpetraExt_PointToBlockDiagPermute::PackAndPrepare(const Epetra_SrcDistObject& Source, 00608 int NumExportIDs, 00609 int* ExportLIDs, 00610 int& LenExports, 00611 char*& Exports, 00612 int& SizeOfPacket, 00613 int* Sizes, 00614 bool & VarSizes, 00615 Epetra_Distributor& Distor){ 00616 return -1; 00617 } 00618 00619 //========================================================================= 00620 // Perform any unpacking and combining after call to DoTransfer(). 00621 int EpetraExt_PointToBlockDiagPermute::UnpackAndCombine(const Epetra_SrcDistObject& Source, 00622 int NumImportIDs, 00623 int* ImportLIDs, 00624 int LenImports, 00625 char* Imports, 00626 int& SizeOfPacket, 00627 Epetra_Distributor& Distor, 00628 Epetra_CombineMode CombineMode, 00629 const Epetra_OffsetIndex * Indexor){ 00630 return -1; 00631 } 00632 00633
1.7.6.1