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 #include "Ifpack_ConfigDefs.h"
00031 #include "Ifpack_OverlappingRowMatrix.h"
00032 #include "Epetra_RowMatrix.h"
00033 #include "Epetra_Map.h"
00034 #include "Epetra_CrsMatrix.h"
00035 #include "Epetra_Comm.h"
00036 #include "Epetra_MultiVector.h"
00037
00038 #ifdef IFPACK_NODE_AWARE_CODE
00039 #include "Epetra_IntVector.h"
00040 #include "Epetra_MpiComm.h"
00041 #include "Teuchos_Hashtable.hpp"
00042
00043 #include "Teuchos_Array.hpp"
00044 #include "EpetraExt_OperatorOut.h"
00045 #include "EpetraExt_BlockMapOut.h"
00046 extern int ML_NODE_ID;
00047 int IFPACK_NODE_ID;
00048 #endif
00049
00050 using namespace Teuchos;
00051
00052
00053
00054 Ifpack_OverlappingRowMatrix::
00055 Ifpack_OverlappingRowMatrix(const RCP<const Epetra_RowMatrix>& Matrix_in,
00056 int OverlapLevel_in) :
00057 Matrix_(Matrix_in),
00058 OverlapLevel_(OverlapLevel_in)
00059 {
00060
00061 if (OverlapLevel_in == 0)
00062 IFPACK_CHK_ERRV(-1);
00063
00064
00065 if (Comm().NumProc() == 1)
00066 IFPACK_CHK_ERRV(-1);
00067
00068 NumMyRowsA_ = A().NumMyRows();
00069
00070
00071 vector<int> ExtElements;
00072
00073 RCP<Epetra_Map> TmpMap;
00074 RCP<Epetra_CrsMatrix> TmpMatrix;
00075 RCP<Epetra_Import> TmpImporter;
00076
00077
00078
00079 const Epetra_Map *RowMap;
00080 const Epetra_Map *ColMap;
00081
00082 for (int overlap = 0 ; overlap < OverlapLevel_in ; ++overlap) {
00083 if (TmpMatrix != Teuchos::null) {
00084 RowMap = &(TmpMatrix->RowMatrixRowMap());
00085 ColMap = &(TmpMatrix->RowMatrixColMap());
00086 }
00087 else {
00088 RowMap = &(A().RowMatrixRowMap());
00089 ColMap = &(A().RowMatrixColMap());
00090 }
00091
00092 int size = ColMap->NumMyElements() - RowMap->NumMyElements();
00093 vector<int> list(size);
00094
00095 int count = 0;
00096
00097
00098 for (int i = 0 ; i < ColMap->NumMyElements() ; ++i) {
00099 int GID = ColMap->GID(i);
00100 if (A().RowMatrixRowMap().LID(GID) == -1) {
00101 vector<int>::iterator pos
00102 = find(ExtElements.begin(),ExtElements.end(),GID);
00103 if (pos == ExtElements.end()) {
00104 ExtElements.push_back(GID);
00105 list[count] = GID;
00106 ++count;
00107 }
00108 }
00109 }
00110
00111 TmpMap = rcp( new Epetra_Map(-1,count, &list[0],0,Comm()) );
00112
00113 TmpMatrix = rcp( new Epetra_CrsMatrix(Copy,*TmpMap,0) );
00114
00115 TmpImporter = rcp( new Epetra_Import(*TmpMap,A().RowMatrixRowMap()) );
00116
00117 TmpMatrix->Import(A(),*TmpImporter,Insert);
00118 TmpMatrix->FillComplete(A().OperatorDomainMap(),*TmpMap);
00119
00120 }
00121
00122
00123
00124 vector<int> list(NumMyRowsA_ + ExtElements.size());
00125 for (int i = 0 ; i < NumMyRowsA_ ; ++i)
00126 list[i] = A().RowMatrixRowMap().GID(i);
00127 for (int i = 0 ; i < (int)ExtElements.size() ; ++i)
00128 list[i + NumMyRowsA_] = ExtElements[i];
00129
00130 Map_ = rcp( new Epetra_Map(-1, NumMyRowsA_ + ExtElements.size(),
00131 &list[0], 0, Comm()) );
00132 # ifdef IFPACK_NODE_AWARE_CODE
00133 colMap_ = &*Map_;
00134 # endif
00135
00136
00137 ExtMap_ = rcp( new Epetra_Map(-1,ExtElements.size(),
00138 &ExtElements[0],0,A().Comm()) );
00139 ExtMatrix_ = rcp( new Epetra_CrsMatrix(Copy,*ExtMap_,*Map_,0) );
00140
00141 ExtImporter_ = rcp( new Epetra_Import(*ExtMap_,A().RowMatrixRowMap()) );
00142 ExtMatrix_->Import(A(),*ExtImporter_,Insert);
00143 ExtMatrix_->FillComplete(A().OperatorDomainMap(),*Map_);
00144
00145 Importer_ = rcp( new Epetra_Import(*Map_,A().RowMatrixRowMap()) );
00146
00147
00148 NumMyRowsB_ = B().NumMyRows();
00149 NumMyRows_ = NumMyRowsA_ + NumMyRowsB_;
00150 NumMyCols_ = NumMyRows_;
00151
00152 NumMyDiagonals_ = A().NumMyDiagonals() + B().NumMyDiagonals();
00153
00154 NumMyNonzeros_ = A().NumMyNonzeros() + B().NumMyNonzeros();
00155 Comm().SumAll(&NumMyNonzeros_,&NumGlobalNonzeros_,1);
00156 MaxNumEntries_ = A().MaxNumEntries();
00157
00158 if (MaxNumEntries_ < B().MaxNumEntries())
00159 MaxNumEntries_ = B().MaxNumEntries();
00160
00161 }
00162
00163 #ifdef IFPACK_NODE_AWARE_CODE
00164
00165
00166
00167 Ifpack_OverlappingRowMatrix::
00168 Ifpack_OverlappingRowMatrix(const RCP<const Epetra_RowMatrix>& Matrix_in,
00169 int OverlapLevel_in, int nodeID) :
00170 Matrix_(Matrix_in),
00171 OverlapLevel_(OverlapLevel_in)
00172 {
00173
00174
00175 #ifdef IFPACK_OVA_TIME_BUILD
00176 double t0 = MPI_Wtime();
00177 double t1, true_t0=t0;
00178 #endif
00179
00180
00181 const int votesMax = INT_MAX;
00182
00183
00184 if (OverlapLevel_in == 0)
00185 IFPACK_CHK_ERRV(-1);
00186
00187
00188 if (Comm().NumProc() == 1)
00189 IFPACK_CHK_ERRV(-1);
00190
00191
00192
00193
00194 # ifdef HAVE_MPI
00195 const Epetra_MpiComm *pComm = dynamic_cast<const Epetra_MpiComm*>( &Comm() );
00196 assert(pComm != NULL);
00197 MPI_Comm nodeMPIComm;
00198 MPI_Comm_split(pComm->Comm(),nodeID,pComm->MyPID(),&nodeMPIComm);
00199 Epetra_MpiComm *nodeComm = new Epetra_MpiComm(nodeMPIComm);
00200 # else
00201 Epetra_SerialComm *nodeComm = dynamic_cast<Epetra_MpiComm*>( &(Matrix_in->RowMatrixRowMap().Comm()) );
00202 # endif
00203
00204 NumMyRowsA_ = A().NumMyRows();
00205
00206 RCP<Epetra_Map> TmpMap;
00207 RCP<Epetra_CrsMatrix> TmpMatrix;
00208 RCP<Epetra_Import> TmpImporter;
00209
00210
00211
00212 const Epetra_Map *RowMap;
00213 const Epetra_Map *ColMap;
00214 const Epetra_Map *DomainMap;
00215
00216
00217
00218
00219
00220 #ifdef IFPACK_OVA_TIME_BUILD
00221 t1 = MPI_Wtime();
00222 fprintf(stderr,"[node %d]: time for initialization %2.3e\n", nodeID, t1-t0);
00223 t0=t1;
00224 #endif
00225
00226
00227 #ifdef IFPACK_OVA_TIME_BUILD
00228 t1 = MPI_Wtime();
00229 fprintf(stderr,"[node %d]: overlap hash table allocs %2.3e\n", nodeID, t1-t0);
00230 t0=t1;
00231 #endif
00232
00233 Teuchos::Hashtable<int,int> colMapTable(3 * A().RowMatrixColMap().NumMyElements() );
00234
00235
00236
00237 Teuchos::Hashtable<int,int> ghostTable(3 * A().RowMatrixColMap().NumMyElements() );
00238
00239
00240
00241
00242 for (int overlap = 0 ; overlap < OverlapLevel_in ; ++overlap)
00243 {
00244
00245
00246
00247 if (TmpMatrix != Teuchos::null) {
00248 RowMap = &(TmpMatrix->RowMatrixRowMap());
00249 ColMap = &(TmpMatrix->RowMatrixColMap());
00250 DomainMap = &(TmpMatrix->OperatorDomainMap());
00251 }
00252 else {
00253 RowMap = &(A().RowMatrixRowMap());
00254 ColMap = &(A().RowMatrixColMap());
00255 DomainMap = &(A().OperatorDomainMap());
00256 }
00257
00258 #ifdef IFPACK_OVA_TIME_BUILD
00259 t1 = MPI_Wtime();
00260 fprintf(stderr,"[node %d]: overlap pointer pulls %2.3e\n", nodeID, t1-t0);
00261 t0=t1;
00262 #endif
00263
00264
00265
00266 Epetra_IntVector colIdList( *ColMap );
00267 Epetra_IntVector rowIdList(*DomainMap);
00268 rowIdList.PutValue(nodeID);
00269 Teuchos::RCP<Epetra_Import> nodeIdImporter = rcp(new Epetra_Import( *ColMap, *DomainMap ));
00270
00271 #ifdef IFPACK_OVA_TIME_BUILD
00272 t1 = MPI_Wtime();
00273 fprintf(stderr,"[node %d]: overlap intvector/importer allocs %2.3e\n", nodeID, t1-t0);
00274 t0=t1;
00275 #endif
00276
00277 colIdList.Import(rowIdList,*nodeIdImporter,Insert);
00278
00279 int size = ColMap->NumMyElements() - RowMap->NumMyElements();
00280 int count = 0;
00281
00282 #ifdef IFPACK_OVA_TIME_BUILD
00283 t1 = MPI_Wtime();
00284 fprintf(stderr,"[node %d]: overlap col/row ID imports %2.3e\n", nodeID, t1-t0);
00285 t0=t1;
00286 #endif
00287
00288
00289
00290
00291 for (int i = 0 ; i < ColMap->NumMyElements() ; ++i) {
00292 int GID = ColMap->GID(i);
00293 int myvotes=0;
00294 if (ghostTable.containsKey(GID)) myvotes = ghostTable.get(GID);
00295 if ( colIdList[i] != nodeID && myvotes < votesMax)
00296 {
00297 int votes;
00298 if (ghostTable.containsKey(GID)) {
00299 votes = ghostTable.get(GID);
00300 votes++;
00301 ghostTable.put(GID,votes);
00302 } else {
00303 ghostTable.put(GID,1);
00304 }
00305 }
00306 }
00307
00308 Teuchos::Array<int> gidsarray,votesarray;
00309 ghostTable.arrayify(gidsarray,votesarray);
00310 int *gids = gidsarray.getRawPtr();
00311 int *votes = votesarray.getRawPtr();
00312
00313
00314
00315
00316
00317
00318
00319
00320
00321
00322 # ifdef HAVE_MPI //FIXME What if we build in serial?! This file will likely not compile.
00323 int *cullList;
00324 int ncull;
00325
00326
00327
00328 #ifdef IFPACK_OVA_TIME_BUILD
00329 t1 = MPI_Wtime();
00330 fprintf(stderr,"[node %d]: overlap before culling %2.3e\n", nodeID, t1-t0);
00331 t0=t1;
00332 #endif
00333
00334
00335
00336 if (nodeComm->MyPID() == 0)
00337 {
00338
00339 MPI_Status status;
00340 int *lengths = new int[nodeComm->NumProc()+1];
00341 lengths[0] = 0;
00342 lengths[1] = ghostTable.size();
00343 for (int i=1; i<nodeComm->NumProc(); i++) {
00344 int leng;
00345 MPI_Recv( &leng, 1, MPI_INT, i, MPI_ANY_TAG, nodeComm->Comm(), &status);
00346 lengths[i+1] = lengths[i] + leng;
00347 }
00348 int total = lengths[nodeComm->NumProc()];
00349
00350 int* ghosts = new int[total];
00351 for (int i=0; i<total; i++) ghosts[i] = -9;
00352 int *round = new int[total];
00353 int *owningpid = new int[total];
00354
00355 for (int i=1; i<nodeComm->NumProc(); i++) {
00356 int count = lengths[i+1] - lengths[i];
00357 MPI_Recv( ghosts+lengths[i], count, MPI_INT, i, MPI_ANY_TAG, nodeComm->Comm(), &status);
00358 MPI_Recv( round+lengths[i], count, MPI_INT, i, MPI_ANY_TAG, nodeComm->Comm(), &status);
00359 }
00360
00361
00362 for (int i=0; i<lengths[1]; i++) {
00363 ghosts[i] = gids[i];
00364 round[i] = votes[i];
00365 owningpid[i] = 0;
00366 }
00367
00368
00369 for (int j=1; j<nodeComm->NumProc(); j++) {
00370 for (int i=lengths[j]; i<lengths[j+1]; i++) {
00371 owningpid[i] = j;
00372 }
00373 }
00374
00375
00376 int* roundpid[2];
00377 roundpid[0] = round; roundpid[1] = owningpid;
00378 Epetra_Util epetraUtil;
00379 epetraUtil.Sort(true,total,ghosts,0,0,2,roundpid);
00380
00381
00382 int *nlosers = new int[nodeComm->NumProc()];
00383 int** losers = new int*[nodeComm->NumProc()];
00384 for (int i=0; i<nodeComm->NumProc(); i++) {
00385 nlosers[i] = 0;
00386 losers[i] = new int[ lengths[i+1]-lengths[i] ];
00387 }
00388
00389
00390
00391
00392
00393
00394
00395 int max=0;
00396
00397
00398 for (int i=1; i<total; i++) {
00399 if (ghosts[i] == ghosts[i-1]) {
00400 int idx = i;
00401 if (round[i] > round[max]) {
00402 idx = max;
00403 max=i;
00404 }
00405 int j = owningpid[idx];
00406 losers[j][nlosers[j]++] = ghosts[idx];
00407 } else {
00408 max=i;
00409 }
00410 }
00411
00412 delete [] round;
00413 delete [] ghosts;
00414 delete [] owningpid;
00415
00416
00417 for (int i=1; i<nodeComm->NumProc(); i++) {
00418 MPI_Send( nlosers+i, 1, MPI_INT, i, 8675, nodeComm->Comm());
00419 MPI_Send( losers[i], nlosers[i], MPI_INT, i, 8675, nodeComm->Comm());
00420 }
00421
00422
00423
00424
00425 ncull = nlosers[0];
00426 cullList = new int[ncull+1];
00427 for (int i=0; i<nlosers[0]; i++)
00428 cullList[i] = losers[0][i];
00429
00430 for (int i=0; i<nodeComm->NumProc(); i++)
00431 delete [] losers[i];
00432
00433 delete [] lengths;
00434 delete [] losers;
00435 delete [] nlosers;
00436
00437 } else {
00438
00439
00440 int hashsize = ghostTable.size();
00441 MPI_Send( &hashsize, 1, MPI_INT, 0, 8675, nodeComm->Comm());
00442 MPI_Send( gids, hashsize, MPI_INT, 0, 8675, nodeComm->Comm());
00443 MPI_Send( votes, hashsize, MPI_INT, 0, 8675, nodeComm->Comm());
00444
00445
00446 MPI_Status status;
00447 MPI_Recv( &ncull, 1, MPI_INT, 0, 8675, nodeComm->Comm(), &status);
00448 cullList = new int[ncull+1];
00449 MPI_Recv( cullList, ncull, MPI_INT, 0, 8675, nodeComm->Comm(), &status);
00450 }
00451
00452
00453
00454
00455
00456 for (int i=0; i<ncull; i++) {
00457 try{ghostTable.remove(cullList[i]);}
00458
00459 catch(...) {
00460 printf("pid %d: In OverlappingRowMatrix ctr, problem removing ghost elt %d from ghostTable\n",
00461 Comm().MyPID(),cullList[i]);
00462 fflush(stdout);
00463 }
00464 }
00465
00466 delete [] cullList;
00467
00468
00469
00470 gidsarray.clear(); votesarray.clear();
00471 ghostTable.arrayify(gidsarray,votesarray);
00472
00473 vector<int> list(size);
00474 count=0;
00475 for (int i=0; i<ghostTable.size(); i++) {
00476
00477 if (votesarray[i] < votesMax) {
00478 list[count] = gidsarray[i];
00479 ghostTable.put(gidsarray[i],votesMax);
00480 count++;
00481 }
00482 }
00483
00484
00485 #ifdef IFPACK_OVA_TIME_BUILD
00486 t1 = MPI_Wtime();
00487 fprintf(stderr,"[node %d]: overlap duplicate removal %2.3e\n", nodeID, t1-t0);
00488 t0=t1;
00489 #endif
00490
00491
00492 # endif //ifdef HAVE_MPI
00493
00494 TmpMap = rcp( new Epetra_Map(-1,count, &list[0],0,Comm()) );
00495
00496 TmpMatrix = rcp( new Epetra_CrsMatrix(Copy,*TmpMap,0) );
00497
00498 TmpImporter = rcp( new Epetra_Import(*TmpMap,A().RowMatrixRowMap()) );
00499
00500 TmpMatrix->Import(A(),*TmpImporter,Insert);
00501 TmpMatrix->FillComplete(A().OperatorDomainMap(),*TmpMap);
00502
00503
00504 #ifdef IFPACK_OVA_TIME_BUILD
00505 t1 = MPI_Wtime();
00506 fprintf(stderr,"[node %d]: overlap TmpMatrix fillcomplete %2.3e\n", nodeID, t1-t0);
00507 t0=t1;
00508 #endif
00509
00510
00511
00512
00513
00514
00515
00516
00517 Epetra_IntVector ov_colIdList( TmpMatrix->ColMap() );
00518 ov_colIdList.PutValue(-1);
00519 Epetra_IntVector ov_rowIdList( TmpMatrix->RowMap() );
00520 ov_rowIdList.PutValue(nodeID);
00521 Teuchos::RCP<Epetra_Import> ov_nodeIdImporter = rcp(new Epetra_Import( TmpMatrix->ColMap(), TmpMatrix->RowMap()));
00522 ov_colIdList.Import(ov_rowIdList,*ov_nodeIdImporter,Insert);
00523
00524 for (int i=0 ; i<ov_colIdList.MyLength(); i++) {
00525 if (ov_colIdList[i] == nodeID) {
00526 int GID = (ov_colIdList.Map()).GID(i);
00527 colMapTable.put(GID,1);
00528 }
00529 }
00530
00531
00532
00533
00534 ov_colIdList.PutValue(-1);
00535 Epetra_IntVector ArowIdList( A().RowMatrixRowMap() );
00536 ArowIdList.PutValue(nodeID);
00537 nodeIdImporter = rcp(new Epetra_Import( TmpMatrix->ColMap(), A().RowMatrixRowMap() ));
00538 ov_colIdList.Import(ArowIdList,*nodeIdImporter,Insert);
00539
00540 for (int i=0 ; i<ov_colIdList.MyLength(); i++) {
00541 if (ov_colIdList[i] == nodeID) {
00542 int GID = (ov_colIdList.Map()).GID(i);
00543 colMapTable.put(GID,1);
00544 }
00545 }
00546
00547
00548 #ifdef IFPACK_OVA_TIME_BUILD
00549 t1 = MPI_Wtime();
00550 fprintf(stderr,"[node %d]: overlap 2 imports to fix up colmap %2.3e\n", nodeID, t1-t0);
00551 t0=t1;
00552 #endif
00553
00554
00555 }
00556
00557
00558
00559
00560
00561
00562 vector<int> ghostElements;
00563
00564 Teuchos::Array<int> gidsarray,votesarray;
00565 ghostTable.arrayify(gidsarray,votesarray);
00566 for (int i=0; i<ghostTable.size(); i++) {
00567 ghostElements.push_back(gidsarray[i]);
00568 }
00569
00570 for (int i = 0 ; i < A().RowMatrixColMap().NumMyElements() ; ++i) {
00571 int GID = A().RowMatrixColMap().GID(i);
00572
00573 if (colMapTable.containsKey(GID)) {
00574 try{colMapTable.remove(GID);}
00575 catch(...) {
00576 printf("pid %d: In OverlappingRowMatrix ctr, problem removing entry %d from colMapTable\n", Comm().MyPID(),GID);
00577 fflush(stdout);
00578 }
00579 }
00580 }
00581
00582
00583 vector<int> colMapElements;
00584
00585 gidsarray.clear(); votesarray.clear();
00586 colMapTable.arrayify(gidsarray,votesarray);
00587 for (int i=0; i<colMapTable.size(); i++)
00588 colMapElements.push_back(gidsarray[i]);
00589
00590
00591
00592
00593
00594
00595
00596
00597
00598
00599
00600
00601 vector<int> rowList(NumMyRowsA_ + ghostElements.size());
00602 for (int i = 0 ; i < NumMyRowsA_ ; ++i)
00603 rowList[i] = A().RowMatrixRowMap().GID(i);
00604 for (int i = 0 ; i < (int)ghostElements.size() ; ++i)
00605 rowList[i + NumMyRowsA_] = ghostElements[i];
00606
00607
00608 Map_ = rcp( new Epetra_Map(-1, NumMyRowsA_ + ghostElements.size(), &rowList[0], 0, Comm()) );
00609
00610
00611
00612
00613
00614
00615
00616 vector<int> colList(A().RowMatrixColMap().NumMyElements() + colMapElements.size());
00617 int nc = A().RowMatrixColMap().NumMyElements();
00618 for (int i = 0 ; i < nc; i++)
00619 colList[i] = A().RowMatrixColMap().GID(i);
00620 for (int i = 0 ; i < (int)colMapElements.size() ; i++)
00621 colList[nc+i] = colMapElements[i];
00622
00623
00624
00625 colMap_ = new Epetra_Map(-1, A().RowMatrixColMap().NumMyElements() + colMapElements.size(), &colList[0], 0, Comm());
00626
00627
00628
00629 #ifdef IFPACK_OVA_TIME_BUILD
00630 t1 = MPI_Wtime();
00631 fprintf(stderr,"[node %d]: time to init B() col/row maps %2.3e\n", nodeID, t1-t0);
00632 t0=t1;
00633 #endif
00634
00635
00636
00637
00638
00639
00640
00641 try {
00642
00643 Epetra_Map* nodeMap = new Epetra_Map(-1,NumMyRowsA_ + ghostElements.size(),&rowList[0],0,*nodeComm);
00644 int numMyElts = colMap_->NumMyElements();
00645 assert(numMyElts!=0);
00646
00647
00648
00649
00650
00651 int* myGlobalElts = new int[numMyElts];
00652 colMap_->MyGlobalElements(myGlobalElts);
00653 int* pidList = new int[numMyElts];
00654 nodeMap->RemoteIDList(numMyElts, myGlobalElts, pidList, 0);
00655
00656
00657 Epetra_Util Util;
00658 int *tt[1];
00659 tt[0] = myGlobalElts;
00660 Util.Sort(true, numMyElts, pidList, 0, (double**)0, 1, tt);
00661
00662
00663
00664 int localStart=0;
00665 while (localStart<numMyElts) {
00666 int currPID = (pidList)[localStart];
00667 int i=localStart;
00668 while (i<numMyElts && (pidList)[i] == currPID) i++;
00669 if (currPID != nodeComm->MyPID())
00670 Util.Sort(true, i-localStart, (myGlobalElts)+localStart, 0, 0, 0, 0);
00671 localStart = i;
00672 }
00673
00674
00675 localStart=0;
00676 while (localStart<numMyElts && (pidList)[localStart] != nodeComm->MyPID()) localStart++;
00677 assert(localStart != numMyElts);
00678 int localEnd=localStart;
00679 while (localEnd<numMyElts && (pidList)[localEnd] == nodeComm->MyPID()) localEnd++;
00680 int* mySortedGlobalElts = new int[numMyElts];
00681 int leng = localEnd - localStart;
00682
00683
00684
00685 int *rowGlobalElts = nodeMap->MyGlobalElements();
00686
00687
00688
00689
00690
00691
00692
00693
00694
00695
00696
00697
00698
00699
00700
00701
00702
00703
00704
00705 for (int i=0; i<leng; i++) {
00706
00707 (mySortedGlobalElts)[i] = rowGlobalElts[i];
00708
00709
00710 }
00711 for (int i=leng; i< localEnd; i++) {
00712 (mySortedGlobalElts)[i] = (myGlobalElts)[i-leng];
00713 }
00714 for (int i=localEnd; i<numMyElts; i++) {
00715 (mySortedGlobalElts)[i] = (myGlobalElts)[i];
00716 }
00717
00718
00719 #ifdef IFPACK_OVA_TIME_BUILD
00720 t1 = MPI_Wtime();
00721 fprintf(stderr,"[node %d]: time to sort col map arrays (cp=1) %2.3e\n", nodeID, t1-t0);
00722 t0=t1;
00723 #endif
00724
00725
00726 int indexBase = colMap_->IndexBase();
00727 if (colMap_) delete colMap_;
00728 colMap_ = new Epetra_Map(-1,numMyElts,mySortedGlobalElts,indexBase,Comm());
00729
00730 delete nodeMap;
00731 delete [] myGlobalElts;
00732 delete [] pidList;
00733 delete [] mySortedGlobalElts;
00734
00735 }
00736 catch(...) {
00737 printf("** * Ifpack_OverlappingRowmatrix ctor: problem creating column map * **\n\n");
00738 }
00739
00740
00741
00742
00743
00744
00745 #ifdef IFPACK_OVA_TIME_BUILD
00746 t1 = MPI_Wtime();
00747 fprintf(stderr,"[node %d]: time to sort col map arrays (cp=2) %2.3e\n", nodeID, t1-t0);
00748 t0=t1;
00749 #endif
00750
00751
00752
00753
00754
00755
00756
00757
00758
00759
00760
00761
00762
00763
00764
00765
00766
00767
00768
00769
00770
00771
00772
00773
00774
00775
00776
00777 ExtMap_ = rcp( new Epetra_Map(-1,ghostElements.size(), &ghostElements[0],0,Comm()) );
00778 ExtMatrix_ = rcp( new Epetra_CrsMatrix(Copy,*ExtMap_,*colMap_,0) );
00779
00780 ExtImporter_ = rcp( new Epetra_Import(*ExtMap_,A().RowMatrixRowMap()) );
00781 ExtMatrix_->Import(A(),*ExtImporter_,Insert);
00782
00783
00784 #ifdef IFPACK_OVA_TIME_BUILD
00785 t1 = MPI_Wtime();
00786 fprintf(stderr,"[node %d]: time to import and setup ExtMap_ %2.3e\n", nodeID, t1-t0);
00787 t0=t1;
00788 #endif
00789
00790
00791
00792
00793
00794
00795
00796
00797
00798
00799 ExtMatrix_->FillComplete( *colMap_ , *Map_ );
00800
00801
00802 #ifdef IFPACK_OVA_TIME_BUILD
00803 t1 = MPI_Wtime();
00804 fprintf(stderr,"[node %d]: time to FillComplete B() %2.3e\n", nodeID, t1-t0);
00805 t0=t1;
00806 #endif
00807
00808
00809
00810
00811
00812 Importer_ = rcp( new Epetra_Import(*Map_,A().RowMatrixRowMap()) );
00813
00814
00815 NumMyRowsB_ = B().NumMyRows();
00816 NumMyRows_ = NumMyRowsA_ + NumMyRowsB_;
00817
00818
00819 NumMyCols_ = B().NumMyCols();
00820
00821
00822
00823 NumMyDiagonals_ = A().NumMyDiagonals() + B().NumMyDiagonals();
00824
00825 NumMyNonzeros_ = A().NumMyNonzeros() + B().NumMyNonzeros();
00826 Comm().SumAll(&NumMyNonzeros_,&NumGlobalNonzeros_,1);
00827 MaxNumEntries_ = A().MaxNumEntries();
00828
00829 if (MaxNumEntries_ < B().MaxNumEntries())
00830 MaxNumEntries_ = B().MaxNumEntries();
00831
00832 # ifdef HAVE_MPI
00833 delete nodeComm;
00834 # endif
00835
00836
00837 #ifdef IFPACK_OVA_TIME_BUILD
00838 t1 = MPI_Wtime();
00839 fprintf(stderr,"[node %d]: time for final calcs %2.3e\n", nodeID, t1-t0);
00840 fprintf(stderr,"[node %d]: total IORM ctor time %2.3e\n", nodeID, t1-true_t0);
00841 #endif
00842
00843
00844
00845 }
00846
00847
00848 Ifpack_OverlappingRowMatrix::~Ifpack_OverlappingRowMatrix() {
00849 delete colMap_;
00850 }
00851
00852 #endif //ifdef IFPACK_NODE_AWARE_CODE
00853
00854
00855
00856 int Ifpack_OverlappingRowMatrix::
00857 NumMyRowEntries(int MyRow, int & NumEntries) const
00858 {
00859 if (MyRow < NumMyRowsA_)
00860 return(A().NumMyRowEntries(MyRow,NumEntries));
00861 else
00862 return(B().NumMyRowEntries(MyRow - NumMyRowsA_, NumEntries));
00863 }
00864
00865 #ifdef IFPACK_NODE_AWARE_CODE
00866
00867 int Ifpack_OverlappingRowMatrix::
00868 ExtractMyRowCopy(int LocRow, int Length, int & NumEntries, double *Values,
00869 int * Indices) const
00870 {
00871 assert(1==0);
00872 int ierr;
00873 const Epetra_Map *Themap;
00874 if (LocRow < NumMyRowsA_) {
00875 ierr = A().ExtractMyRowCopy(LocRow,Length,NumEntries,Values,Indices);
00876 Themap=&A().RowMatrixColMap();
00877 }
00878 else {
00879 ierr = B().ExtractMyRowCopy(LocRow-NumMyRowsA_,Length,NumEntries,Values,Indices);
00880 Themap=&B().RowMatrixColMap();
00881 }
00882
00883 IFPACK_RETURN(ierr);
00884 }
00885
00886
00887 int Ifpack_OverlappingRowMatrix::
00888 ExtractGlobalRowCopy(int GlobRow, int Length, int & NumEntries, double *Values,
00889 int * Indices) const
00890 {
00891 int ierr;
00892 const Epetra_Map *Themap;
00893 int LocRow = A().RowMatrixRowMap().LID(GlobRow);
00894 if (LocRow < NumMyRowsA_ && LocRow != -1) {
00895 ierr = A().ExtractMyRowCopy(LocRow,Length,NumEntries,Values,Indices);
00896 Themap=&A().RowMatrixColMap();
00897 }
00898 else {
00899 LocRow = B().RowMatrixRowMap().LID(GlobRow);
00900 assert(LocRow!=-1);
00901
00902 ierr = B().ExtractMyRowCopy(LocRow,Length,NumEntries,Values,Indices);
00903 Themap=&B().RowMatrixColMap();
00904 }
00905
00906 for (int i=0; i<NumEntries; i++) {
00907 Indices[i]=Themap->GID(Indices[i]);
00908 assert(Indices[i] != -1);
00909 }
00910
00911 IFPACK_RETURN(ierr);
00912 }
00913 #else
00914
00915
00916 int Ifpack_OverlappingRowMatrix::
00917 ExtractMyRowCopy(int MyRow, int Length, int & NumEntries, double *Values,
00918 int * Indices) const
00919 {
00920 int ierr;
00921 if (MyRow < NumMyRowsA_)
00922 ierr = A().ExtractMyRowCopy(MyRow,Length,NumEntries,Values,Indices);
00923 else
00924 ierr = B().ExtractMyRowCopy(MyRow - NumMyRowsA_,Length,NumEntries,
00925 Values,Indices);
00926
00927 IFPACK_RETURN(ierr);
00928 }
00929 #endif
00930
00931
00932 int Ifpack_OverlappingRowMatrix::
00933 ExtractDiagonalCopy(Epetra_Vector & Diagonal) const
00934 {
00935 IFPACK_CHK_ERR(-1);
00936 }
00937
00938
00939
00940 int Ifpack_OverlappingRowMatrix::
00941 Multiply(bool TransA, const Epetra_MultiVector& X, Epetra_MultiVector& Y) const
00942 {
00943 int NumVectors = X.NumVectors();
00944 vector<int> Ind(MaxNumEntries_);
00945 vector<double> Val(MaxNumEntries_);
00946
00947 Y.PutScalar(0.0);
00948
00949
00950 for (int i = 0 ; i < NumMyRowsA_ ; ++i) {
00951 for (int k = 0 ; k < NumVectors ; ++k) {
00952 int Nnz;
00953 IFPACK_CHK_ERR(A().ExtractMyRowCopy(i,MaxNumEntries_,Nnz,
00954 &Val[0], &Ind[0]));
00955 for (int j = 0 ; j < Nnz ; ++j) {
00956 Y[k][i] += Val[j] * X[k][Ind[j]];
00957 }
00958 }
00959 }
00960
00961
00962 for (int i = 0 ; i < NumMyRowsB_ ; ++i) {
00963 for (int k = 0 ; k < NumVectors ; ++k) {
00964 int Nnz;
00965 IFPACK_CHK_ERR(B().ExtractMyRowCopy(i,MaxNumEntries_,Nnz,
00966 &Val[0], &Ind[0]));
00967 for (int j = 0 ; j < Nnz ; ++j) {
00968 Y[k][i + NumMyRowsA_] += Val[j] * X[k][Ind[j]];
00969 }
00970 }
00971 }
00972 return(0);
00973 }
00974
00975
00976 int Ifpack_OverlappingRowMatrix::
00977 Apply(const Epetra_MultiVector& X, Epetra_MultiVector& Y) const
00978 {
00979 IFPACK_CHK_ERR(Multiply(UseTranspose(),X,Y));
00980 return(0);
00981 }
00982
00983
00984 int Ifpack_OverlappingRowMatrix::
00985 ApplyInverse(const Epetra_MultiVector& X, Epetra_MultiVector& Y) const
00986 {
00987 IFPACK_CHK_ERR(-1);
00988 }
00989
00990
00991 #ifndef IFPACK_NODE_AWARE_CODE
00992 Epetra_RowMatrix& Ifpack_OverlappingRowMatrix::B() const
00993 {
00994 return(*ExtMatrix_);
00995 }
00996 #endif
00997
00998 const Epetra_BlockMap& Ifpack_OverlappingRowMatrix::Map() const
00999 {
01000 return(*Map_);
01001 }
01002
01003
01004 int Ifpack_OverlappingRowMatrix::
01005 ImportMultiVector(const Epetra_MultiVector& X, Epetra_MultiVector& OvX,
01006 Epetra_CombineMode CM)
01007 {
01008 OvX.Import(X,*Importer_,CM);
01009 return(0);
01010 }
01011
01012
01013 int Ifpack_OverlappingRowMatrix::
01014 ExportMultiVector(const Epetra_MultiVector& OvX, Epetra_MultiVector& X,
01015 Epetra_CombineMode CM)
01016 {
01017 X.Export(OvX,*Importer_,CM);
01018 return(0);
01019 }
01020