00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043 #include "Ifpack_ConfigDefs.h"
00044 #include "Ifpack_OverlappingRowMatrix.h"
00045 #include "Epetra_RowMatrix.h"
00046 #include "Epetra_Map.h"
00047 #include "Epetra_CrsMatrix.h"
00048 #include "Epetra_Comm.h"
00049 #include "Epetra_MultiVector.h"
00050
00051 #ifdef HAVE_IFPACK_PARALLEL_SUBDOMAIN_SOLVERS
00052 #include "Epetra_IntVector.h"
00053 #include "Epetra_MpiComm.h"
00054 #include "Teuchos_Hashtable.hpp"
00055 #include "Teuchos_Array.hpp"
00056 #include "EpetraExt_OperatorOut.h"
00057 #include "EpetraExt_BlockMapOut.h"
00058 #else
00059 # ifdef IFPACK_NODE_AWARE_CODE
00060 # include "Epetra_IntVector.h"
00061 # include "Epetra_MpiComm.h"
00062 # include "Teuchos_Hashtable.hpp"
00063 # include "Teuchos_Array.hpp"
00064 # include "EpetraExt_OperatorOut.h"
00065 # include "EpetraExt_BlockMapOut.h"
00066 extern int ML_NODE_ID;
00067 int IFPACK_NODE_ID;
00068 # endif
00069 #endif
00070
00071 using namespace Teuchos;
00072
00073
00074
00075 template<typename int_type>
00076 void Ifpack_OverlappingRowMatrix::BuildMap(int OverlapLevel_in)
00077 {
00078 RCP<Epetra_Map> TmpMap;
00079 RCP<Epetra_CrsMatrix> TmpMatrix;
00080 RCP<Epetra_Import> TmpImporter;
00081
00082
00083
00084 const Epetra_Map *RowMap;
00085 const Epetra_Map *ColMap;
00086
00087 std::vector<int_type> ExtElements;
00088
00089 for (int overlap = 0 ; overlap < OverlapLevel_in ; ++overlap) {
00090 if (TmpMatrix != Teuchos::null) {
00091 RowMap = &(TmpMatrix->RowMatrixRowMap());
00092 ColMap = &(TmpMatrix->RowMatrixColMap());
00093 }
00094 else {
00095 RowMap = &(A().RowMatrixRowMap());
00096 ColMap = &(A().RowMatrixColMap());
00097 }
00098
00099 int size = ColMap->NumMyElements() - RowMap->NumMyElements();
00100 std::vector<int_type> list(size);
00101
00102 int count = 0;
00103
00104
00105 for (int i = 0 ; i < ColMap->NumMyElements() ; ++i) {
00106 int_type GID = (int_type) ColMap->GID64(i);
00107 if (A().RowMatrixRowMap().LID(GID) == -1) {
00108 typename std::vector<int_type>::iterator pos
00109 = find(ExtElements.begin(),ExtElements.end(),GID);
00110 if (pos == ExtElements.end()) {
00111 ExtElements.push_back(GID);
00112 list[count] = GID;
00113 ++count;
00114 }
00115 }
00116 }
00117
00118 const int_type *listptr = NULL;
00119 if ( ! list.empty() ) listptr = &list[0];
00120 TmpMap = rcp( new Epetra_Map(-1,count, listptr,0,Comm()) );
00121
00122 TmpMatrix = rcp( new Epetra_CrsMatrix(Copy,*TmpMap,0) );
00123
00124 TmpImporter = rcp( new Epetra_Import(*TmpMap,A().RowMatrixRowMap()) );
00125
00126 TmpMatrix->Import(A(),*TmpImporter,Insert);
00127 TmpMatrix->FillComplete(A().OperatorDomainMap(),*TmpMap);
00128
00129 }
00130
00131
00132
00133 std::vector<int_type> list(NumMyRowsA_ + ExtElements.size());
00134 for (int i = 0 ; i < NumMyRowsA_ ; ++i)
00135 list[i] = (int_type) A().RowMatrixRowMap().GID64(i);
00136 for (int i = 0 ; i < (int)ExtElements.size() ; ++i)
00137 list[i + NumMyRowsA_] = ExtElements[i];
00138
00139 const int_type *listptr = NULL;
00140 if ( ! list.empty() ) listptr = &list[0];
00141 {
00142 Map_ = rcp( new Epetra_Map((int_type) -1, NumMyRowsA_ + ExtElements.size(),
00143 listptr, 0, Comm()) );
00144 }
00145 #ifdef HAVE_IFPACK_PARALLEL_SUBDOMAIN_SOLVERS
00146 colMap_ = &*Map_;
00147 #else
00148 # ifdef IFPACK_NODE_AWARE_CODE
00149 colMap_ = &*Map_;
00150 # endif
00151 #endif
00152
00153
00154 {
00155 const int_type * extelsptr = NULL;
00156 if ( ! ExtElements.empty() ) extelsptr = &ExtElements[0];
00157 ExtMap_ = rcp( new Epetra_Map((int_type) -1,ExtElements.size(),
00158 extelsptr,0,A().Comm()) );
00159 }
00160 }
00161
00162
00163
00164 Ifpack_OverlappingRowMatrix::
00165 Ifpack_OverlappingRowMatrix(const RCP<const Epetra_RowMatrix>& Matrix_in,
00166 int OverlapLevel_in) :
00167 Matrix_(Matrix_in),
00168 OverlapLevel_(OverlapLevel_in)
00169 {
00170
00171 if (OverlapLevel_in == 0)
00172 IFPACK_CHK_ERRV(-1);
00173
00174
00175 if (Comm().NumProc() == 1)
00176 IFPACK_CHK_ERRV(-1);
00177
00178 NumMyRowsA_ = A().NumMyRows();
00179
00180
00181
00182 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
00183 if(A().RowMatrixRowMap().GlobalIndicesInt()) {
00184 BuildMap<int>(OverlapLevel_in);
00185 }
00186 else
00187 #endif
00188 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
00189 if(A().RowMatrixRowMap().GlobalIndicesLongLong()) {
00190 BuildMap<long long>(OverlapLevel_in);
00191 }
00192 else
00193 #endif
00194 throw "Ifpack_OverlappingRowMatrix::Ifpack_OverlappingRowMatrix: GlobalIndices type unknown";
00195
00196 ExtMatrix_ = rcp( new Epetra_CrsMatrix(Copy,*ExtMap_,*Map_,0) );
00197
00198 ExtImporter_ = rcp( new Epetra_Import(*ExtMap_,A().RowMatrixRowMap()) );
00199 ExtMatrix_->Import(A(),*ExtImporter_,Insert);
00200 ExtMatrix_->FillComplete(A().OperatorDomainMap(),*Map_);
00201
00202 Importer_ = rcp( new Epetra_Import(*Map_,A().RowMatrixRowMap()) );
00203
00204
00205 NumMyRowsB_ = B().NumMyRows();
00206 NumMyRows_ = NumMyRowsA_ + NumMyRowsB_;
00207 NumMyCols_ = NumMyRows_;
00208
00209 NumMyDiagonals_ = A().NumMyDiagonals() + B().NumMyDiagonals();
00210
00211 NumMyNonzeros_ = A().NumMyNonzeros() + B().NumMyNonzeros();
00212 long long NumMyNonzeros_tmp = NumMyNonzeros_;
00213 Comm().SumAll(&NumMyNonzeros_tmp,&NumGlobalNonzeros_,1);
00214 MaxNumEntries_ = A().MaxNumEntries();
00215
00216 if (MaxNumEntries_ < B().MaxNumEntries())
00217 MaxNumEntries_ = B().MaxNumEntries();
00218
00219 }
00220
00221 #ifdef HAVE_IFPACK_PARALLEL_SUBDOMAIN_SOLVERS
00222
00223
00224
00225 Ifpack_OverlappingRowMatrix::
00226 Ifpack_OverlappingRowMatrix(const RCP<const Epetra_RowMatrix>& Matrix_in,
00227 int OverlapLevel_in, int subdomainID) :
00228 Matrix_(Matrix_in),
00229 OverlapLevel_(OverlapLevel_in)
00230 {
00231
00232
00233 #ifdef IFPACK_OVA_TIME_BUILD
00234 double t0 = MPI_Wtime();
00235 double t1, true_t0=t0;
00236 #endif
00237
00238
00239 const int votesMax = INT_MAX;
00240
00241
00242 if (OverlapLevel_in == 0)
00243 IFPACK_CHK_ERRV(-1);
00244
00245
00246 if (Comm().NumProc() == 1)
00247 IFPACK_CHK_ERRV(-1);
00248
00249
00250
00251
00252 # ifdef HAVE_MPI
00253 const Epetra_MpiComm *pComm = dynamic_cast<const Epetra_MpiComm*>( &Comm() );
00254 assert(pComm != NULL);
00255 MPI_Comm nodeMPIComm;
00256 MPI_Comm_split(pComm->Comm(),subdomainID,pComm->MyPID(),&nodeMPIComm);
00257 Epetra_MpiComm *nodeComm = new Epetra_MpiComm(nodeMPIComm);
00258 # else
00259 Epetra_SerialComm *nodeComm = dynamic_cast<Epetra_MpiComm*>( &(Matrix_in->RowMatrixRowMap().Comm()) );
00260 # endif
00261
00262 NumMyRowsA_ = A().NumMyRows();
00263
00264 RCP<Epetra_Map> TmpMap;
00265 RCP<Epetra_CrsMatrix> TmpMatrix;
00266 RCP<Epetra_Import> TmpImporter;
00267
00268
00269
00270 const Epetra_Map *RowMap;
00271 const Epetra_Map *ColMap;
00272 const Epetra_Map *DomainMap;
00273
00274
00275
00276
00277
00278 #ifdef IFPACK_OVA_TIME_BUILD
00279 t1 = MPI_Wtime();
00280 fprintf(stderr,"[node %d]: time for initialization %2.3e\n", subdomainID, t1-t0);
00281 t0=t1;
00282 #endif
00283
00284
00285 #ifdef IFPACK_OVA_TIME_BUILD
00286 t1 = MPI_Wtime();
00287 fprintf(stderr,"[node %d]: overlap hash table allocs %2.3e\n", subdomainID, t1-t0);
00288 t0=t1;
00289 #endif
00290
00291 Teuchos::Hashtable<int,int> colMapTable(3 * A().RowMatrixColMap().NumMyElements() );
00292
00293
00294
00295 Teuchos::Hashtable<int,int> ghostTable(3 * A().RowMatrixColMap().NumMyElements() );
00296
00297
00298
00299
00300 for (int overlap = 0 ; overlap < OverlapLevel_in ; ++overlap)
00301 {
00302
00303
00304
00305 if (TmpMatrix != Teuchos::null) {
00306 RowMap = &(TmpMatrix->RowMatrixRowMap());
00307 ColMap = &(TmpMatrix->RowMatrixColMap());
00308 DomainMap = &(TmpMatrix->OperatorDomainMap());
00309 }
00310 else {
00311 RowMap = &(A().RowMatrixRowMap());
00312 ColMap = &(A().RowMatrixColMap());
00313 DomainMap = &(A().OperatorDomainMap());
00314 }
00315
00316 #ifdef IFPACK_OVA_TIME_BUILD
00317 t1 = MPI_Wtime();
00318 fprintf(stderr,"[node %d]: overlap pointer pulls %2.3e\n", subdomainID, t1-t0);
00319 t0=t1;
00320 #endif
00321
00322
00323
00324 Epetra_IntVector colIdList( *ColMap );
00325 Epetra_IntVector rowIdList(*DomainMap);
00326 rowIdList.PutValue(subdomainID);
00327 Teuchos::RCP<Epetra_Import> nodeIdImporter = rcp(new Epetra_Import( *ColMap, *DomainMap ));
00328
00329 #ifdef IFPACK_OVA_TIME_BUILD
00330 t1 = MPI_Wtime();
00331 fprintf(stderr,"[node %d]: overlap intvector/importer allocs %2.3e\n", subdomainID, t1-t0);
00332 t0=t1;
00333 #endif
00334
00335 colIdList.Import(rowIdList,*nodeIdImporter,Insert);
00336
00337 int size = ColMap->NumMyElements() - RowMap->NumMyElements();
00338 int count = 0;
00339
00340 #ifdef IFPACK_OVA_TIME_BUILD
00341 t1 = MPI_Wtime();
00342 fprintf(stderr,"[node %d]: overlap col/row ID imports %2.3e\n", subdomainID, t1-t0);
00343 t0=t1;
00344 #endif
00345
00346
00347
00348
00349 for (int i = 0 ; i < ColMap->NumMyElements() ; ++i) {
00350 int GID = ColMap->GID(i);
00351 int myvotes=0;
00352 if (ghostTable.containsKey(GID)) myvotes = ghostTable.get(GID);
00353 if ( colIdList[i] != subdomainID && myvotes < votesMax)
00354 {
00355 int votes;
00356 if (ghostTable.containsKey(GID)) {
00357 votes = ghostTable.get(GID);
00358 votes++;
00359 ghostTable.put(GID,votes);
00360 } else {
00361 ghostTable.put(GID,1);
00362 }
00363 }
00364 }
00365
00366 Teuchos::Array<int> gidsarray,votesarray;
00367 ghostTable.arrayify(gidsarray,votesarray);
00368 int *gids = gidsarray.getRawPtr();
00369 int *votes = votesarray.getRawPtr();
00370
00371
00372
00373
00374
00375
00376
00377
00378
00379
00380 # ifdef HAVE_MPI //FIXME What if we build in serial?! This file will likely not compile.
00381 int *cullList;
00382 int ncull;
00383
00384
00385
00386 #ifdef IFPACK_OVA_TIME_BUILD
00387 t1 = MPI_Wtime();
00388 fprintf(stderr,"[node %d]: overlap before culling %2.3e\n", subdomainID, t1-t0);
00389 t0=t1;
00390 #endif
00391
00392
00393
00394 if (nodeComm->MyPID() == 0)
00395 {
00396
00397 MPI_Status status;
00398 int *lengths = new int[nodeComm->NumProc()+1];
00399 lengths[0] = 0;
00400 lengths[1] = ghostTable.size();
00401 for (int i=1; i<nodeComm->NumProc(); i++) {
00402 int leng;
00403 MPI_Recv( &leng, 1, MPI_INT, i, MPI_ANY_TAG, nodeComm->Comm(), &status);
00404 lengths[i+1] = lengths[i] + leng;
00405 }
00406 int total = lengths[nodeComm->NumProc()];
00407
00408 int* ghosts = new int[total];
00409 for (int i=0; i<total; i++) ghosts[i] = -9;
00410 int *round = new int[total];
00411 int *owningpid = new int[total];
00412
00413 for (int i=1; i<nodeComm->NumProc(); i++) {
00414 int count = lengths[i+1] - lengths[i];
00415 MPI_Recv( ghosts+lengths[i], count, MPI_INT, i, MPI_ANY_TAG, nodeComm->Comm(), &status);
00416 MPI_Recv( round+lengths[i], count, MPI_INT, i, MPI_ANY_TAG, nodeComm->Comm(), &status);
00417 }
00418
00419
00420 for (int i=0; i<lengths[1]; i++) {
00421 ghosts[i] = gids[i];
00422 round[i] = votes[i];
00423 owningpid[i] = 0;
00424 }
00425
00426
00427 for (int j=1; j<nodeComm->NumProc(); j++) {
00428 for (int i=lengths[j]; i<lengths[j+1]; i++) {
00429 owningpid[i] = j;
00430 }
00431 }
00432
00433
00434 int* roundpid[2];
00435 roundpid[0] = round; roundpid[1] = owningpid;
00436 Epetra_Util epetraUtil;
00437 epetraUtil.Sort(true,total,ghosts,0,0,2,roundpid);
00438
00439
00440 int *nlosers = new int[nodeComm->NumProc()];
00441 int** losers = new int*[nodeComm->NumProc()];
00442 for (int i=0; i<nodeComm->NumProc(); i++) {
00443 nlosers[i] = 0;
00444 losers[i] = new int[ lengths[i+1]-lengths[i] ];
00445 }
00446
00447
00448
00449
00450
00451
00452
00453 int max=0;
00454
00455
00456 for (int i=1; i<total; i++) {
00457 if (ghosts[i] == ghosts[i-1]) {
00458 int idx = i;
00459 if (round[i] > round[max]) {
00460 idx = max;
00461 max=i;
00462 }
00463 int j = owningpid[idx];
00464 losers[j][nlosers[j]++] = ghosts[idx];
00465 } else {
00466 max=i;
00467 }
00468 }
00469
00470 delete [] round;
00471 delete [] ghosts;
00472 delete [] owningpid;
00473
00474
00475 for (int i=1; i<nodeComm->NumProc(); i++) {
00476 MPI_Send( nlosers+i, 1, MPI_INT, i, 8675, nodeComm->Comm());
00477 MPI_Send( losers[i], nlosers[i], MPI_INT, i, 8675, nodeComm->Comm());
00478 }
00479
00480
00481
00482
00483 ncull = nlosers[0];
00484 cullList = new int[ncull+1];
00485 for (int i=0; i<nlosers[0]; i++)
00486 cullList[i] = losers[0][i];
00487
00488 for (int i=0; i<nodeComm->NumProc(); i++)
00489 delete [] losers[i];
00490
00491 delete [] lengths;
00492 delete [] losers;
00493 delete [] nlosers;
00494
00495 } else {
00496
00497
00498 int hashsize = ghostTable.size();
00499 MPI_Send( &hashsize, 1, MPI_INT, 0, 8675, nodeComm->Comm());
00500 MPI_Send( gids, hashsize, MPI_INT, 0, 8675, nodeComm->Comm());
00501 MPI_Send( votes, hashsize, MPI_INT, 0, 8675, nodeComm->Comm());
00502
00503
00504 MPI_Status status;
00505 MPI_Recv( &ncull, 1, MPI_INT, 0, 8675, nodeComm->Comm(), &status);
00506 cullList = new int[ncull+1];
00507 MPI_Recv( cullList, ncull, MPI_INT, 0, 8675, nodeComm->Comm(), &status);
00508 }
00509
00510
00511
00512
00513
00514 for (int i=0; i<ncull; i++) {
00515 try{ghostTable.remove(cullList[i]);}
00516
00517 catch(...) {
00518 printf("pid %d: In OverlappingRowMatrix ctr, problem removing ghost elt %d from ghostTable\n",
00519 Comm().MyPID(),cullList[i]);
00520 fflush(stdout);
00521 }
00522 }
00523
00524 delete [] cullList;
00525
00526
00527
00528 gidsarray.clear(); votesarray.clear();
00529 ghostTable.arrayify(gidsarray,votesarray);
00530
00531 std::vector<int> list(size);
00532 count=0;
00533 for (int i=0; i<ghostTable.size(); i++) {
00534
00535 if (votesarray[i] < votesMax) {
00536 list[count] = gidsarray[i];
00537 ghostTable.put(gidsarray[i],votesMax);
00538 count++;
00539 }
00540 }
00541
00542
00543 #ifdef IFPACK_OVA_TIME_BUILD
00544 t1 = MPI_Wtime();
00545 fprintf(stderr,"[node %d]: overlap duplicate removal %2.3e\n", subdomainID, t1-t0);
00546 t0=t1;
00547 #endif
00548
00549
00550 # endif //ifdef HAVE_MPI
00551
00552 TmpMap = rcp( new Epetra_Map(-1,count, &list[0],0,Comm()) );
00553
00554 TmpMatrix = rcp( new Epetra_CrsMatrix(Copy,*TmpMap,0) );
00555
00556 TmpImporter = rcp( new Epetra_Import(*TmpMap,A().RowMatrixRowMap()) );
00557
00558 TmpMatrix->Import(A(),*TmpImporter,Insert);
00559 TmpMatrix->FillComplete(A().OperatorDomainMap(),*TmpMap);
00560
00561
00562 #ifdef IFPACK_OVA_TIME_BUILD
00563 t1 = MPI_Wtime();
00564 fprintf(stderr,"[node %d]: overlap TmpMatrix fillcomplete %2.3e\n", subdomainID, t1-t0);
00565 t0=t1;
00566 #endif
00567
00568
00569
00570
00571
00572
00573
00574
00575 Epetra_IntVector ov_colIdList( TmpMatrix->ColMap() );
00576 ov_colIdList.PutValue(-1);
00577 Epetra_IntVector ov_rowIdList( TmpMatrix->RowMap() );
00578 ov_rowIdList.PutValue(subdomainID);
00579 Teuchos::RCP<Epetra_Import> ov_nodeIdImporter = rcp(new Epetra_Import( TmpMatrix->ColMap(), TmpMatrix->RowMap()));
00580 ov_colIdList.Import(ov_rowIdList,*ov_nodeIdImporter,Insert);
00581
00582 for (int i=0 ; i<ov_colIdList.MyLength(); i++) {
00583 if (ov_colIdList[i] == subdomainID) {
00584 int GID = (ov_colIdList.Map()).GID(i);
00585 colMapTable.put(GID,1);
00586 }
00587 }
00588
00589
00590
00591
00592 ov_colIdList.PutValue(-1);
00593 Epetra_IntVector ArowIdList( A().RowMatrixRowMap() );
00594 ArowIdList.PutValue(subdomainID);
00595 nodeIdImporter = rcp(new Epetra_Import( TmpMatrix->ColMap(), A().RowMatrixRowMap() ));
00596 ov_colIdList.Import(ArowIdList,*nodeIdImporter,Insert);
00597
00598 for (int i=0 ; i<ov_colIdList.MyLength(); i++) {
00599 if (ov_colIdList[i] == subdomainID) {
00600 int GID = (ov_colIdList.Map()).GID(i);
00601 colMapTable.put(GID,1);
00602 }
00603 }
00604
00605
00606 #ifdef IFPACK_OVA_TIME_BUILD
00607 t1 = MPI_Wtime();
00608 fprintf(stderr,"[node %d]: overlap 2 imports to fix up colmap %2.3e\n", subdomainID, t1-t0);
00609 t0=t1;
00610 #endif
00611
00612
00613 }
00614
00615
00616
00617
00618
00619
00620 std::vector<int> ghostElements;
00621
00622 Teuchos::Array<int> gidsarray,votesarray;
00623 ghostTable.arrayify(gidsarray,votesarray);
00624 for (int i=0; i<ghostTable.size(); i++) {
00625 ghostElements.push_back(gidsarray[i]);
00626 }
00627
00628 for (int i = 0 ; i < A().RowMatrixColMap().NumMyElements() ; ++i) {
00629 int GID = A().RowMatrixColMap().GID(i);
00630
00631 if (colMapTable.containsKey(GID)) {
00632 try{colMapTable.remove(GID);}
00633 catch(...) {
00634 printf("pid %d: In OverlappingRowMatrix ctr, problem removing entry %d from colMapTable\n", Comm().MyPID(),GID);
00635 fflush(stdout);
00636 }
00637 }
00638 }
00639
00640
00641 std::vector<int> colMapElements;
00642
00643 gidsarray.clear(); votesarray.clear();
00644 colMapTable.arrayify(gidsarray,votesarray);
00645 for (int i=0; i<colMapTable.size(); i++)
00646 colMapElements.push_back(gidsarray[i]);
00647
00648
00649
00650
00651
00652
00653
00654
00655
00656
00657
00658
00659 std::vector<int> rowList(NumMyRowsA_ + ghostElements.size());
00660 for (int i = 0 ; i < NumMyRowsA_ ; ++i)
00661 rowList[i] = A().RowMatrixRowMap().GID(i);
00662 for (int i = 0 ; i < (int)ghostElements.size() ; ++i)
00663 rowList[i + NumMyRowsA_] = ghostElements[i];
00664
00665
00666 Map_ = rcp( new Epetra_Map(-1, NumMyRowsA_ + ghostElements.size(), &rowList[0], 0, Comm()) );
00667
00668
00669
00670
00671
00672
00673
00674 std::vector<int> colList(A().RowMatrixColMap().NumMyElements() + colMapElements.size());
00675 int nc = A().RowMatrixColMap().NumMyElements();
00676 for (int i = 0 ; i < nc; i++)
00677 colList[i] = A().RowMatrixColMap().GID(i);
00678 for (int i = 0 ; i < (int)colMapElements.size() ; i++)
00679 colList[nc+i] = colMapElements[i];
00680
00681
00682
00683 colMap_ = new Epetra_Map(-1, A().RowMatrixColMap().NumMyElements() + colMapElements.size(), &colList[0], 0, Comm());
00684
00685
00686
00687 #ifdef IFPACK_OVA_TIME_BUILD
00688 t1 = MPI_Wtime();
00689 fprintf(stderr,"[node %d]: time to init B() col/row maps %2.3e\n", subdomainID, t1-t0);
00690 t0=t1;
00691 #endif
00692
00693
00694
00695
00696
00697
00698
00699 try {
00700
00701 Epetra_Map* nodeMap = new Epetra_Map(-1,NumMyRowsA_ + ghostElements.size(),&rowList[0],0,*nodeComm);
00702 int numMyElts = colMap_->NumMyElements();
00703 assert(numMyElts!=0);
00704
00705
00706
00707
00708
00709 int* myGlobalElts = new int[numMyElts];
00710 colMap_->MyGlobalElements(myGlobalElts);
00711 int* pidList = new int[numMyElts];
00712 nodeMap->RemoteIDList(numMyElts, myGlobalElts, pidList, 0);
00713
00714
00715 Epetra_Util Util;
00716 int *tt[1];
00717 tt[0] = myGlobalElts;
00718 Util.Sort(true, numMyElts, pidList, 0, (double**)0, 1, tt);
00719
00720
00721
00722 int localStart=0;
00723 while (localStart<numMyElts) {
00724 int currPID = (pidList)[localStart];
00725 int i=localStart;
00726 while (i<numMyElts && (pidList)[i] == currPID) i++;
00727 if (currPID != nodeComm->MyPID())
00728 Util.Sort(true, i-localStart, (myGlobalElts)+localStart, 0, 0, 0, 0);
00729 localStart = i;
00730 }
00731
00732
00733 localStart=0;
00734 while (localStart<numMyElts && (pidList)[localStart] != nodeComm->MyPID()) localStart++;
00735 assert(localStart != numMyElts);
00736 int localEnd=localStart;
00737 while (localEnd<numMyElts && (pidList)[localEnd] == nodeComm->MyPID()) localEnd++;
00738 int* mySortedGlobalElts = new int[numMyElts];
00739 int leng = localEnd - localStart;
00740
00741
00742
00743 int *rowGlobalElts = nodeMap->MyGlobalElements();
00744
00745
00746
00747
00748
00749
00750
00751
00752
00753
00754
00755
00756
00757
00758
00759
00760
00761
00762
00763 for (int i=0; i<leng; i++) {
00764
00765 (mySortedGlobalElts)[i] = rowGlobalElts[i];
00766
00767
00768 }
00769 for (int i=leng; i< localEnd; i++) {
00770 (mySortedGlobalElts)[i] = (myGlobalElts)[i-leng];
00771 }
00772 for (int i=localEnd; i<numMyElts; i++) {
00773 (mySortedGlobalElts)[i] = (myGlobalElts)[i];
00774 }
00775
00776
00777 #ifdef IFPACK_OVA_TIME_BUILD
00778 t1 = MPI_Wtime();
00779 fprintf(stderr,"[node %d]: time to sort col map arrays (cp=1) %2.3e\n", subdomainID, t1-t0);
00780 t0=t1;
00781 #endif
00782
00783
00784 int indexBase = colMap_->IndexBase();
00785 if (colMap_) delete colMap_;
00786 colMap_ = new Epetra_Map(-1,numMyElts,mySortedGlobalElts,indexBase,Comm());
00787
00788 delete nodeMap;
00789 delete [] myGlobalElts;
00790 delete [] pidList;
00791 delete [] mySortedGlobalElts;
00792
00793 }
00794 catch(...) {
00795 printf("** * Ifpack_OverlappingRowmatrix ctor: problem creating column map * **\n\n");
00796 }
00797
00798
00799
00800
00801
00802
00803 #ifdef IFPACK_OVA_TIME_BUILD
00804 t1 = MPI_Wtime();
00805 fprintf(stderr,"[node %d]: time to sort col map arrays (cp=2) %2.3e\n", subdomainID, t1-t0);
00806 t0=t1;
00807 #endif
00808
00809
00810
00811
00812
00813
00814
00815
00816
00817
00818
00819
00820
00821
00822
00823
00824
00825
00826
00827
00828
00829
00830
00831
00832
00833
00834
00835 ExtMap_ = rcp( new Epetra_Map(-1,ghostElements.size(), &ghostElements[0],0,Comm()) );
00836 ExtMatrix_ = rcp( new Epetra_CrsMatrix(Copy,*ExtMap_,*colMap_,0) );
00837
00838 ExtImporter_ = rcp( new Epetra_Import(*ExtMap_,A().RowMatrixRowMap()) );
00839 ExtMatrix_->Import(A(),*ExtImporter_,Insert);
00840
00841
00842 #ifdef IFPACK_OVA_TIME_BUILD
00843 t1 = MPI_Wtime();
00844 fprintf(stderr,"[node %d]: time to import and setup ExtMap_ %2.3e\n", subdomainID, t1-t0);
00845 t0=t1;
00846 #endif
00847
00848
00849
00850
00851
00852
00853
00854
00855
00856
00857 ExtMatrix_->FillComplete( *colMap_ , *Map_ );
00858
00859
00860 #ifdef IFPACK_OVA_TIME_BUILD
00861 t1 = MPI_Wtime();
00862 fprintf(stderr,"[node %d]: time to FillComplete B() %2.3e\n", subdomainID, t1-t0);
00863 t0=t1;
00864 #endif
00865
00866
00867
00868
00869
00870 Importer_ = rcp( new Epetra_Import(*Map_,A().RowMatrixRowMap()) );
00871
00872
00873 NumMyRowsB_ = B().NumMyRows();
00874 NumMyRows_ = NumMyRowsA_ + NumMyRowsB_;
00875
00876
00877 NumMyCols_ = B().NumMyCols();
00878
00879
00880
00881 NumMyDiagonals_ = A().NumMyDiagonals() + B().NumMyDiagonals();
00882
00883 NumMyNonzeros_ = A().NumMyNonzeros() + B().NumMyNonzeros();
00884 long long NumMyNonZerosTemp = NumMyNonzeros_;
00885 Comm().SumAll(&NumMyNonZerosTemp,&NumGlobalNonzeros_,1);
00886 MaxNumEntries_ = A().MaxNumEntries();
00887
00888 if (MaxNumEntries_ < B().MaxNumEntries())
00889 MaxNumEntries_ = B().MaxNumEntries();
00890
00891 # ifdef HAVE_MPI
00892 delete nodeComm;
00893 # endif
00894
00895
00896 #ifdef IFPACK_OVA_TIME_BUILD
00897 t1 = MPI_Wtime();
00898 fprintf(stderr,"[node %d]: time for final calcs %2.3e\n", subdomainID, t1-t0);
00899 fprintf(stderr,"[node %d]: total IORM ctor time %2.3e\n", subdomainID, t1-true_t0);
00900 #endif
00901
00902
00903
00904 }
00905
00906 #else
00907 # ifdef IFPACK_NODE_AWARE_CODE
00908
00909
00910
00911 Ifpack_OverlappingRowMatrix::
00912 Ifpack_OverlappingRowMatrix(const RCP<const Epetra_RowMatrix>& Matrix_in,
00913 int OverlapLevel_in, int nodeID) :
00914 Matrix_(Matrix_in),
00915 OverlapLevel_(OverlapLevel_in)
00916 {
00917
00918
00919 #ifdef IFPACK_OVA_TIME_BUILD
00920 double t0 = MPI_Wtime();
00921 double t1, true_t0=t0;
00922 #endif
00923
00924
00925 const int votesMax = INT_MAX;
00926
00927
00928 if (OverlapLevel_in == 0)
00929 IFPACK_CHK_ERRV(-1);
00930
00931
00932 if (Comm().NumProc() == 1)
00933 IFPACK_CHK_ERRV(-1);
00934
00935
00936
00937
00938 # ifdef HAVE_MPI
00939 const Epetra_MpiComm *pComm = dynamic_cast<const Epetra_MpiComm*>( &Comm() );
00940 assert(pComm != NULL);
00941 MPI_Comm nodeMPIComm;
00942 MPI_Comm_split(pComm->Comm(),nodeID,pComm->MyPID(),&nodeMPIComm);
00943 Epetra_MpiComm *nodeComm = new Epetra_MpiComm(nodeMPIComm);
00944 # else
00945 Epetra_SerialComm *nodeComm = dynamic_cast<Epetra_MpiComm*>( &(Matrix_in->RowMatrixRowMap().Comm()) );
00946 # endif
00947
00948 NumMyRowsA_ = A().NumMyRows();
00949
00950 RCP<Epetra_Map> TmpMap;
00951 RCP<Epetra_CrsMatrix> TmpMatrix;
00952 RCP<Epetra_Import> TmpImporter;
00953
00954
00955
00956 const Epetra_Map *RowMap;
00957 const Epetra_Map *ColMap;
00958 const Epetra_Map *DomainMap;
00959
00960
00961
00962
00963
00964 #ifdef IFPACK_OVA_TIME_BUILD
00965 t1 = MPI_Wtime();
00966 fprintf(stderr,"[node %d]: time for initialization %2.3e\n", nodeID, t1-t0);
00967 t0=t1;
00968 #endif
00969
00970
00971 #ifdef IFPACK_OVA_TIME_BUILD
00972 t1 = MPI_Wtime();
00973 fprintf(stderr,"[node %d]: overlap hash table allocs %2.3e\n", nodeID, t1-t0);
00974 t0=t1;
00975 #endif
00976
00977 Teuchos::Hashtable<int,int> colMapTable(3 * A().RowMatrixColMap().NumMyElements() );
00978
00979
00980
00981 Teuchos::Hashtable<int,int> ghostTable(3 * A().RowMatrixColMap().NumMyElements() );
00982
00983
00984
00985
00986 for (int overlap = 0 ; overlap < OverlapLevel_in ; ++overlap)
00987 {
00988
00989
00990
00991 if (TmpMatrix != Teuchos::null) {
00992 RowMap = &(TmpMatrix->RowMatrixRowMap());
00993 ColMap = &(TmpMatrix->RowMatrixColMap());
00994 DomainMap = &(TmpMatrix->OperatorDomainMap());
00995 }
00996 else {
00997 RowMap = &(A().RowMatrixRowMap());
00998 ColMap = &(A().RowMatrixColMap());
00999 DomainMap = &(A().OperatorDomainMap());
01000 }
01001
01002 #ifdef IFPACK_OVA_TIME_BUILD
01003 t1 = MPI_Wtime();
01004 fprintf(stderr,"[node %d]: overlap pointer pulls %2.3e\n", nodeID, t1-t0);
01005 t0=t1;
01006 #endif
01007
01008
01009
01010 Epetra_IntVector colIdList( *ColMap );
01011 Epetra_IntVector rowIdList(*DomainMap);
01012 rowIdList.PutValue(nodeID);
01013 Teuchos::RCP<Epetra_Import> nodeIdImporter = rcp(new Epetra_Import( *ColMap, *DomainMap ));
01014
01015 #ifdef IFPACK_OVA_TIME_BUILD
01016 t1 = MPI_Wtime();
01017 fprintf(stderr,"[node %d]: overlap intvector/importer allocs %2.3e\n", nodeID, t1-t0);
01018 t0=t1;
01019 #endif
01020
01021 colIdList.Import(rowIdList,*nodeIdImporter,Insert);
01022
01023 int size = ColMap->NumMyElements() - RowMap->NumMyElements();
01024 int count = 0;
01025
01026 #ifdef IFPACK_OVA_TIME_BUILD
01027 t1 = MPI_Wtime();
01028 fprintf(stderr,"[node %d]: overlap col/row ID imports %2.3e\n", nodeID, t1-t0);
01029 t0=t1;
01030 #endif
01031
01032
01033
01034
01035 for (int i = 0 ; i < ColMap->NumMyElements() ; ++i) {
01036 int GID = ColMap->GID(i);
01037 int myvotes=0;
01038 if (ghostTable.containsKey(GID)) myvotes = ghostTable.get(GID);
01039 if ( colIdList[i] != nodeID && myvotes < votesMax)
01040 {
01041 int votes;
01042 if (ghostTable.containsKey(GID)) {
01043 votes = ghostTable.get(GID);
01044 votes++;
01045 ghostTable.put(GID,votes);
01046 } else {
01047 ghostTable.put(GID,1);
01048 }
01049 }
01050 }
01051
01052 Teuchos::Array<int> gidsarray,votesarray;
01053 ghostTable.arrayify(gidsarray,votesarray);
01054 int *gids = gidsarray.getRawPtr();
01055 int *votes = votesarray.getRawPtr();
01056
01057
01058
01059
01060
01061
01062
01063
01064
01065
01066 # ifdef HAVE_MPI //FIXME What if we build in serial?! This file will likely not compile.
01067 int *cullList;
01068 int ncull;
01069
01070
01071
01072 #ifdef IFPACK_OVA_TIME_BUILD
01073 t1 = MPI_Wtime();
01074 fprintf(stderr,"[node %d]: overlap before culling %2.3e\n", nodeID, t1-t0);
01075 t0=t1;
01076 #endif
01077
01078
01079
01080 if (nodeComm->MyPID() == 0)
01081 {
01082
01083 MPI_Status status;
01084 int *lengths = new int[nodeComm->NumProc()+1];
01085 lengths[0] = 0;
01086 lengths[1] = ghostTable.size();
01087 for (int i=1; i<nodeComm->NumProc(); i++) {
01088 int leng;
01089 MPI_Recv( &leng, 1, MPI_INT, i, MPI_ANY_TAG, nodeComm->Comm(), &status);
01090 lengths[i+1] = lengths[i] + leng;
01091 }
01092 int total = lengths[nodeComm->NumProc()];
01093
01094 int* ghosts = new int[total];
01095 for (int i=0; i<total; i++) ghosts[i] = -9;
01096 int *round = new int[total];
01097 int *owningpid = new int[total];
01098
01099 for (int i=1; i<nodeComm->NumProc(); i++) {
01100 int count = lengths[i+1] - lengths[i];
01101 MPI_Recv( ghosts+lengths[i], count, MPI_INT, i, MPI_ANY_TAG, nodeComm->Comm(), &status);
01102 MPI_Recv( round+lengths[i], count, MPI_INT, i, MPI_ANY_TAG, nodeComm->Comm(), &status);
01103 }
01104
01105
01106 for (int i=0; i<lengths[1]; i++) {
01107 ghosts[i] = gids[i];
01108 round[i] = votes[i];
01109 owningpid[i] = 0;
01110 }
01111
01112
01113 for (int j=1; j<nodeComm->NumProc(); j++) {
01114 for (int i=lengths[j]; i<lengths[j+1]; i++) {
01115 owningpid[i] = j;
01116 }
01117 }
01118
01119
01120 int* roundpid[2];
01121 roundpid[0] = round; roundpid[1] = owningpid;
01122 Epetra_Util epetraUtil;
01123 epetraUtil.Sort(true,total,ghosts,0,0,2,roundpid);
01124
01125
01126 int *nlosers = new int[nodeComm->NumProc()];
01127 int** losers = new int*[nodeComm->NumProc()];
01128 for (int i=0; i<nodeComm->NumProc(); i++) {
01129 nlosers[i] = 0;
01130 losers[i] = new int[ lengths[i+1]-lengths[i] ];
01131 }
01132
01133
01134
01135
01136
01137
01138
01139 int max=0;
01140
01141
01142 for (int i=1; i<total; i++) {
01143 if (ghosts[i] == ghosts[i-1]) {
01144 int idx = i;
01145 if (round[i] > round[max]) {
01146 idx = max;
01147 max=i;
01148 }
01149 int j = owningpid[idx];
01150 losers[j][nlosers[j]++] = ghosts[idx];
01151 } else {
01152 max=i;
01153 }
01154 }
01155
01156 delete [] round;
01157 delete [] ghosts;
01158 delete [] owningpid;
01159
01160
01161 for (int i=1; i<nodeComm->NumProc(); i++) {
01162 MPI_Send( nlosers+i, 1, MPI_INT, i, 8675, nodeComm->Comm());
01163 MPI_Send( losers[i], nlosers[i], MPI_INT, i, 8675, nodeComm->Comm());
01164 }
01165
01166
01167
01168
01169 ncull = nlosers[0];
01170 cullList = new int[ncull+1];
01171 for (int i=0; i<nlosers[0]; i++)
01172 cullList[i] = losers[0][i];
01173
01174 for (int i=0; i<nodeComm->NumProc(); i++)
01175 delete [] losers[i];
01176
01177 delete [] lengths;
01178 delete [] losers;
01179 delete [] nlosers;
01180
01181 } else {
01182
01183
01184 int hashsize = ghostTable.size();
01185 MPI_Send( &hashsize, 1, MPI_INT, 0, 8675, nodeComm->Comm());
01186 MPI_Send( gids, hashsize, MPI_INT, 0, 8675, nodeComm->Comm());
01187 MPI_Send( votes, hashsize, MPI_INT, 0, 8675, nodeComm->Comm());
01188
01189
01190 MPI_Status status;
01191 MPI_Recv( &ncull, 1, MPI_INT, 0, 8675, nodeComm->Comm(), &status);
01192 cullList = new int[ncull+1];
01193 MPI_Recv( cullList, ncull, MPI_INT, 0, 8675, nodeComm->Comm(), &status);
01194 }
01195
01196
01197
01198
01199
01200 for (int i=0; i<ncull; i++) {
01201 try{ghostTable.remove(cullList[i]);}
01202
01203 catch(...) {
01204 printf("pid %d: In OverlappingRowMatrix ctr, problem removing ghost elt %d from ghostTable\n",
01205 Comm().MyPID(),cullList[i]);
01206 fflush(stdout);
01207 }
01208 }
01209
01210 delete [] cullList;
01211
01212
01213
01214 gidsarray.clear(); votesarray.clear();
01215 ghostTable.arrayify(gidsarray,votesarray);
01216
01217 std::vector<int> list(size);
01218 count=0;
01219 for (int i=0; i<ghostTable.size(); i++) {
01220
01221 if (votesarray[i] < votesMax) {
01222 list[count] = gidsarray[i];
01223 ghostTable.put(gidsarray[i],votesMax);
01224 count++;
01225 }
01226 }
01227
01228
01229 #ifdef IFPACK_OVA_TIME_BUILD
01230 t1 = MPI_Wtime();
01231 fprintf(stderr,"[node %d]: overlap duplicate removal %2.3e\n", nodeID, t1-t0);
01232 t0=t1;
01233 #endif
01234
01235
01236 # endif //ifdef HAVE_MPI
01237
01238 TmpMap = rcp( new Epetra_Map(-1,count, &list[0],0,Comm()) );
01239
01240 TmpMatrix = rcp( new Epetra_CrsMatrix(Copy,*TmpMap,0) );
01241
01242 TmpImporter = rcp( new Epetra_Import(*TmpMap,A().RowMatrixRowMap()) );
01243
01244 TmpMatrix->Import(A(),*TmpImporter,Insert);
01245 TmpMatrix->FillComplete(A().OperatorDomainMap(),*TmpMap);
01246
01247
01248 #ifdef IFPACK_OVA_TIME_BUILD
01249 t1 = MPI_Wtime();
01250 fprintf(stderr,"[node %d]: overlap TmpMatrix fillcomplete %2.3e\n", nodeID, t1-t0);
01251 t0=t1;
01252 #endif
01253
01254
01255
01256
01257
01258
01259
01260
01261 Epetra_IntVector ov_colIdList( TmpMatrix->ColMap() );
01262 ov_colIdList.PutValue(-1);
01263 Epetra_IntVector ov_rowIdList( TmpMatrix->RowMap() );
01264 ov_rowIdList.PutValue(nodeID);
01265 Teuchos::RCP<Epetra_Import> ov_nodeIdImporter = rcp(new Epetra_Import( TmpMatrix->ColMap(), TmpMatrix->RowMap()));
01266 ov_colIdList.Import(ov_rowIdList,*ov_nodeIdImporter,Insert);
01267
01268 for (int i=0 ; i<ov_colIdList.MyLength(); i++) {
01269 if (ov_colIdList[i] == nodeID) {
01270 int GID = (ov_colIdList.Map()).GID(i);
01271 colMapTable.put(GID,1);
01272 }
01273 }
01274
01275
01276
01277
01278 ov_colIdList.PutValue(-1);
01279 Epetra_IntVector ArowIdList( A().RowMatrixRowMap() );
01280 ArowIdList.PutValue(nodeID);
01281 nodeIdImporter = rcp(new Epetra_Import( TmpMatrix->ColMap(), A().RowMatrixRowMap() ));
01282 ov_colIdList.Import(ArowIdList,*nodeIdImporter,Insert);
01283
01284 for (int i=0 ; i<ov_colIdList.MyLength(); i++) {
01285 if (ov_colIdList[i] == nodeID) {
01286 int GID = (ov_colIdList.Map()).GID(i);
01287 colMapTable.put(GID,1);
01288 }
01289 }
01290
01291
01292 #ifdef IFPACK_OVA_TIME_BUILD
01293 t1 = MPI_Wtime();
01294 fprintf(stderr,"[node %d]: overlap 2 imports to fix up colmap %2.3e\n", nodeID, t1-t0);
01295 t0=t1;
01296 #endif
01297
01298
01299 }
01300
01301
01302
01303
01304
01305
01306 std::vector<int> ghostElements;
01307
01308 Teuchos::Array<int> gidsarray,votesarray;
01309 ghostTable.arrayify(gidsarray,votesarray);
01310 for (int i=0; i<ghostTable.size(); i++) {
01311 ghostElements.push_back(gidsarray[i]);
01312 }
01313
01314 for (int i = 0 ; i < A().RowMatrixColMap().NumMyElements() ; ++i) {
01315 int GID = A().RowMatrixColMap().GID(i);
01316
01317 if (colMapTable.containsKey(GID)) {
01318 try{colMapTable.remove(GID);}
01319 catch(...) {
01320 printf("pid %d: In OverlappingRowMatrix ctr, problem removing entry %d from colMapTable\n", Comm().MyPID(),GID);
01321 fflush(stdout);
01322 }
01323 }
01324 }
01325
01326
01327 std::vector<int> colMapElements;
01328
01329 gidsarray.clear(); votesarray.clear();
01330 colMapTable.arrayify(gidsarray,votesarray);
01331 for (int i=0; i<colMapTable.size(); i++)
01332 colMapElements.push_back(gidsarray[i]);
01333
01334
01335
01336
01337
01338
01339
01340
01341
01342
01343
01344
01345 std::vector<int> rowList(NumMyRowsA_ + ghostElements.size());
01346 for (int i = 0 ; i < NumMyRowsA_ ; ++i)
01347 rowList[i] = A().RowMatrixRowMap().GID(i);
01348 for (int i = 0 ; i < (int)ghostElements.size() ; ++i)
01349 rowList[i + NumMyRowsA_] = ghostElements[i];
01350
01351
01352 Map_ = rcp( new Epetra_Map(-1, NumMyRowsA_ + ghostElements.size(), &rowList[0], 0, Comm()) );
01353
01354
01355
01356
01357
01358
01359
01360 std::vector<int> colList(A().RowMatrixColMap().NumMyElements() + colMapElements.size());
01361 int nc = A().RowMatrixColMap().NumMyElements();
01362 for (int i = 0 ; i < nc; i++)
01363 colList[i] = A().RowMatrixColMap().GID(i);
01364 for (int i = 0 ; i < (int)colMapElements.size() ; i++)
01365 colList[nc+i] = colMapElements[i];
01366
01367
01368
01369 colMap_ = new Epetra_Map(-1, A().RowMatrixColMap().NumMyElements() + colMapElements.size(), &colList[0], 0, Comm());
01370
01371
01372
01373 #ifdef IFPACK_OVA_TIME_BUILD
01374 t1 = MPI_Wtime();
01375 fprintf(stderr,"[node %d]: time to init B() col/row maps %2.3e\n", nodeID, t1-t0);
01376 t0=t1;
01377 #endif
01378
01379
01380
01381
01382
01383
01384
01385 try {
01386
01387 Epetra_Map* nodeMap = new Epetra_Map(-1,NumMyRowsA_ + ghostElements.size(),&rowList[0],0,*nodeComm);
01388 int numMyElts = colMap_->NumMyElements();
01389 assert(numMyElts!=0);
01390
01391
01392
01393
01394
01395 int* myGlobalElts = new int[numMyElts];
01396 colMap_->MyGlobalElements(myGlobalElts);
01397 int* pidList = new int[numMyElts];
01398 nodeMap->RemoteIDList(numMyElts, myGlobalElts, pidList, 0);
01399
01400
01401 Epetra_Util Util;
01402 int *tt[1];
01403 tt[0] = myGlobalElts;
01404 Util.Sort(true, numMyElts, pidList, 0, (double**)0, 1, tt);
01405
01406
01407
01408 int localStart=0;
01409 while (localStart<numMyElts) {
01410 int currPID = (pidList)[localStart];
01411 int i=localStart;
01412 while (i<numMyElts && (pidList)[i] == currPID) i++;
01413 if (currPID != nodeComm->MyPID())
01414 Util.Sort(true, i-localStart, (myGlobalElts)+localStart, 0, 0, 0, 0);
01415 localStart = i;
01416 }
01417
01418
01419 localStart=0;
01420 while (localStart<numMyElts && (pidList)[localStart] != nodeComm->MyPID()) localStart++;
01421 assert(localStart != numMyElts);
01422 int localEnd=localStart;
01423 while (localEnd<numMyElts && (pidList)[localEnd] == nodeComm->MyPID()) localEnd++;
01424 int* mySortedGlobalElts = new int[numMyElts];
01425 int leng = localEnd - localStart;
01426
01427
01428
01429 int *rowGlobalElts = nodeMap->MyGlobalElements();
01430
01431
01432
01433
01434
01435
01436
01437
01438
01439
01440
01441
01442
01443
01444
01445
01446
01447
01448
01449 for (int i=0; i<leng; i++) {
01450
01451 (mySortedGlobalElts)[i] = rowGlobalElts[i];
01452
01453
01454 }
01455 for (int i=leng; i< localEnd; i++) {
01456 (mySortedGlobalElts)[i] = (myGlobalElts)[i-leng];
01457 }
01458 for (int i=localEnd; i<numMyElts; i++) {
01459 (mySortedGlobalElts)[i] = (myGlobalElts)[i];
01460 }
01461
01462
01463 #ifdef IFPACK_OVA_TIME_BUILD
01464 t1 = MPI_Wtime();
01465 fprintf(stderr,"[node %d]: time to sort col map arrays (cp=1) %2.3e\n", nodeID, t1-t0);
01466 t0=t1;
01467 #endif
01468
01469
01470 int indexBase = colMap_->IndexBase();
01471 if (colMap_) delete colMap_;
01472 colMap_ = new Epetra_Map(-1,numMyElts,mySortedGlobalElts,indexBase,Comm());
01473
01474 delete nodeMap;
01475 delete [] myGlobalElts;
01476 delete [] pidList;
01477 delete [] mySortedGlobalElts;
01478
01479 }
01480 catch(...) {
01481 printf("** * Ifpack_OverlappingRowmatrix ctor: problem creating column map * **\n\n");
01482 }
01483
01484
01485
01486
01487
01488
01489 #ifdef IFPACK_OVA_TIME_BUILD
01490 t1 = MPI_Wtime();
01491 fprintf(stderr,"[node %d]: time to sort col map arrays (cp=2) %2.3e\n", nodeID, t1-t0);
01492 t0=t1;
01493 #endif
01494
01495
01496
01497
01498
01499
01500
01501
01502
01503
01504
01505
01506
01507
01508
01509
01510
01511
01512
01513
01514
01515
01516
01517
01518
01519
01520
01521 ExtMap_ = rcp( new Epetra_Map(-1,ghostElements.size(), &ghostElements[0],0,Comm()) );
01522 ExtMatrix_ = rcp( new Epetra_CrsMatrix(Copy,*ExtMap_,*colMap_,0) );
01523
01524 ExtImporter_ = rcp( new Epetra_Import(*ExtMap_,A().RowMatrixRowMap()) );
01525 ExtMatrix_->Import(A(),*ExtImporter_,Insert);
01526
01527
01528 #ifdef IFPACK_OVA_TIME_BUILD
01529 t1 = MPI_Wtime();
01530 fprintf(stderr,"[node %d]: time to import and setup ExtMap_ %2.3e\n", nodeID, t1-t0);
01531 t0=t1;
01532 #endif
01533
01534
01535
01536
01537
01538
01539
01540
01541
01542
01543 ExtMatrix_->FillComplete( *colMap_ , *Map_ );
01544
01545
01546 #ifdef IFPACK_OVA_TIME_BUILD
01547 t1 = MPI_Wtime();
01548 fprintf(stderr,"[node %d]: time to FillComplete B() %2.3e\n", nodeID, t1-t0);
01549 t0=t1;
01550 #endif
01551
01552
01553
01554
01555
01556 Importer_ = rcp( new Epetra_Import(*Map_,A().RowMatrixRowMap()) );
01557
01558
01559 NumMyRowsB_ = B().NumMyRows();
01560 NumMyRows_ = NumMyRowsA_ + NumMyRowsB_;
01561
01562
01563 NumMyCols_ = B().NumMyCols();
01564
01565
01566
01567 NumMyDiagonals_ = A().NumMyDiagonals() + B().NumMyDiagonals();
01568
01569 NumMyNonzeros_ = A().NumMyNonzeros() + B().NumMyNonzeros();
01570 Comm().SumAll(&NumMyNonzeros_,&NumGlobalNonzeros_,1);
01571 MaxNumEntries_ = A().MaxNumEntries();
01572
01573 if (MaxNumEntries_ < B().MaxNumEntries())
01574 MaxNumEntries_ = B().MaxNumEntries();
01575
01576 # ifdef HAVE_MPI
01577 delete nodeComm;
01578 # endif
01579
01580
01581 #ifdef IFPACK_OVA_TIME_BUILD
01582 t1 = MPI_Wtime();
01583 fprintf(stderr,"[node %d]: time for final calcs %2.3e\n", nodeID, t1-t0);
01584 fprintf(stderr,"[node %d]: total IORM ctor time %2.3e\n", nodeID, t1-true_t0);
01585 #endif
01586
01587
01588
01589 }
01590
01591
01592 Ifpack_OverlappingRowMatrix::~Ifpack_OverlappingRowMatrix() {
01593 delete colMap_;
01594 }
01595
01596 # endif //ifdef IFPACK_NODE_AWARE_CODE
01597 #endif // ifdef HAVE_IFPACK_PARALLEL_SUBDOMAIN_SOLVERS
01598
01599
01600
01601 int Ifpack_OverlappingRowMatrix::
01602 NumMyRowEntries(int MyRow, int & NumEntries) const
01603 {
01604 if (MyRow < NumMyRowsA_)
01605 return(A().NumMyRowEntries(MyRow,NumEntries));
01606 else
01607 return(B().NumMyRowEntries(MyRow - NumMyRowsA_, NumEntries));
01608 }
01609
01610 #ifdef HAVE_IFPACK_PARALLEL_SUBDOMAIN_SOLVERS
01611
01612 int Ifpack_OverlappingRowMatrix::
01613 ExtractMyRowCopy(int LocRow, int Length, int & NumEntries, double *Values,
01614 int * Indices) const
01615 {
01616
01617 int ierr;
01618 const Epetra_Map *Themap;
01619 if (LocRow < NumMyRowsA_) {
01620 ierr = A().ExtractMyRowCopy(LocRow,Length,NumEntries,Values,Indices);
01621 Themap=&A().RowMatrixColMap();
01622 }
01623 else {
01624 ierr = B().ExtractMyRowCopy(LocRow-NumMyRowsA_,Length,NumEntries,Values,Indices);
01625 Themap=&B().RowMatrixColMap();
01626 }
01627
01628 IFPACK_RETURN(ierr);
01629 }
01630
01631
01632 int Ifpack_OverlappingRowMatrix::
01633 ExtractGlobalRowCopy(int GlobRow, int Length, int & NumEntries, double *Values,
01634 int * Indices) const
01635 {
01636 int ierr;
01637 const Epetra_Map *Themap;
01638 int LocRow = A().RowMatrixRowMap().LID(GlobRow);
01639 if (LocRow < NumMyRowsA_ && LocRow != -1) {
01640 ierr = A().ExtractMyRowCopy(LocRow,Length,NumEntries,Values,Indices);
01641 Themap=&A().RowMatrixColMap();
01642 }
01643 else {
01644 LocRow = B().RowMatrixRowMap().LID(GlobRow);
01645 assert(LocRow!=-1);
01646
01647 ierr = B().ExtractMyRowCopy(LocRow,Length,NumEntries,Values,Indices);
01648 Themap=&B().RowMatrixColMap();
01649 }
01650
01651 for (int i=0; i<NumEntries; i++) {
01652 Indices[i]=Themap->GID(Indices[i]);
01653 assert(Indices[i] != -1);
01654 }
01655
01656 IFPACK_RETURN(ierr);
01657 }
01658 #else
01659 # ifdef IFPACK_NODE_AWARE_CODE
01660
01661 int Ifpack_OverlappingRowMatrix::
01662 ExtractMyRowCopy(int LocRow, int Length, int & NumEntries, double *Values,
01663 int * Indices) const
01664 {
01665 assert(1==0);
01666 int ierr;
01667 const Epetra_Map *Themap;
01668 if (LocRow < NumMyRowsA_) {
01669 ierr = A().ExtractMyRowCopy(LocRow,Length,NumEntries,Values,Indices);
01670 Themap=&A().RowMatrixColMap();
01671 }
01672 else {
01673 ierr = B().ExtractMyRowCopy(LocRow-NumMyRowsA_,Length,NumEntries,Values,Indices);
01674 Themap=&B().RowMatrixColMap();
01675 }
01676
01677 IFPACK_RETURN(ierr);
01678 }
01679
01680
01681 int Ifpack_OverlappingRowMatrix::
01682 ExtractGlobalRowCopy(int GlobRow, int Length, int & NumEntries, double *Values,
01683 int * Indices) const
01684 {
01685 int ierr;
01686 const Epetra_Map *Themap;
01687 int LocRow = A().RowMatrixRowMap().LID(GlobRow);
01688 if (LocRow < NumMyRowsA_ && LocRow != -1) {
01689 ierr = A().ExtractMyRowCopy(LocRow,Length,NumEntries,Values,Indices);
01690 Themap=&A().RowMatrixColMap();
01691 }
01692 else {
01693 LocRow = B().RowMatrixRowMap().LID(GlobRow);
01694 assert(LocRow!=-1);
01695
01696 ierr = B().ExtractMyRowCopy(LocRow,Length,NumEntries,Values,Indices);
01697 Themap=&B().RowMatrixColMap();
01698 }
01699
01700 for (int i=0; i<NumEntries; i++) {
01701 Indices[i]=Themap->GID(Indices[i]);
01702 assert(Indices[i] != -1);
01703 }
01704
01705 IFPACK_RETURN(ierr);
01706 }
01707 # else
01708
01709
01710 int Ifpack_OverlappingRowMatrix::
01711 ExtractMyRowCopy(int MyRow, int Length, int & NumEntries, double *Values,
01712 int * Indices) const
01713 {
01714 int ierr;
01715 if (MyRow < NumMyRowsA_)
01716 ierr = A().ExtractMyRowCopy(MyRow,Length,NumEntries,Values,Indices);
01717 else
01718 ierr = B().ExtractMyRowCopy(MyRow - NumMyRowsA_,Length,NumEntries,
01719 Values,Indices);
01720
01721 IFPACK_RETURN(ierr);
01722 }
01723 # endif // ifdef IFPACK_NODE_AWARE_CODE
01724 #endif // ifdef HAVE_IFPACK_PARALLEL_SUBDOMAIN_SOLVERS
01725
01726
01727 int Ifpack_OverlappingRowMatrix::
01728 ExtractDiagonalCopy(Epetra_Vector & Diagonal) const
01729 {
01730 IFPACK_CHK_ERR(-1);
01731 }
01732
01733
01734
01735 int Ifpack_OverlappingRowMatrix::
01736 Multiply(bool TransA, const Epetra_MultiVector& X, Epetra_MultiVector& Y) const
01737 {
01738 int NumVectors = X.NumVectors();
01739 std::vector<int> Ind(MaxNumEntries_);
01740 std::vector<double> Val(MaxNumEntries_);
01741
01742 Y.PutScalar(0.0);
01743
01744
01745 for (int i = 0 ; i < NumMyRowsA_ ; ++i) {
01746 for (int k = 0 ; k < NumVectors ; ++k) {
01747 int Nnz;
01748 IFPACK_CHK_ERR(A().ExtractMyRowCopy(i,MaxNumEntries_,Nnz,
01749 &Val[0], &Ind[0]));
01750 for (int j = 0 ; j < Nnz ; ++j) {
01751 Y[k][i] += Val[j] * X[k][Ind[j]];
01752 }
01753 }
01754 }
01755
01756
01757 for (int i = 0 ; i < NumMyRowsB_ ; ++i) {
01758 for (int k = 0 ; k < NumVectors ; ++k) {
01759 int Nnz;
01760 IFPACK_CHK_ERR(B().ExtractMyRowCopy(i,MaxNumEntries_,Nnz,
01761 &Val[0], &Ind[0]));
01762 for (int j = 0 ; j < Nnz ; ++j) {
01763 Y[k][i + NumMyRowsA_] += Val[j] * X[k][Ind[j]];
01764 }
01765 }
01766 }
01767 return(0);
01768 }
01769
01770
01771 int Ifpack_OverlappingRowMatrix::
01772 Apply(const Epetra_MultiVector& X, Epetra_MultiVector& Y) const
01773 {
01774 IFPACK_CHK_ERR(Multiply(UseTranspose(),X,Y));
01775 return(0);
01776 }
01777
01778
01779 int Ifpack_OverlappingRowMatrix::
01780 ApplyInverse(const Epetra_MultiVector& X, Epetra_MultiVector& Y) const
01781 {
01782 IFPACK_CHK_ERR(-1);
01783 }
01784
01785
01786 #ifndef HAVE_IFPACK_PARALLEL_SUBDOMAIN_SOLVERS
01787 # ifndef IFPACK_NODE_AWARE_CODE
01788 Epetra_RowMatrix& Ifpack_OverlappingRowMatrix::B() const
01789 {
01790 return(*ExtMatrix_);
01791 }
01792 # endif
01793 #endif
01794
01795 const Epetra_BlockMap& Ifpack_OverlappingRowMatrix::Map() const
01796 {
01797 return(*Map_);
01798 }
01799
01800
01801 int Ifpack_OverlappingRowMatrix::
01802 ImportMultiVector(const Epetra_MultiVector& X, Epetra_MultiVector& OvX,
01803 Epetra_CombineMode CM)
01804 {
01805 OvX.Import(X,*Importer_,CM);
01806 return(0);
01807 }
01808
01809
01810 int Ifpack_OverlappingRowMatrix::
01811 ExportMultiVector(const Epetra_MultiVector& OvX, Epetra_MultiVector& X,
01812 Epetra_CombineMode CM)
01813 {
01814 X.Export(OvX,*Importer_,CM);
01815 return(0);
01816 }
01817