//EpetraCrsMatrix_To_TpetraCrsMatrix: copies Epetra_CrsMatrix to its analogous Tpetra_CrsMatrix Teuchos::RCP<Tpetra_CrsMatrix> Petra::EpetraCrsMatrix_To_TpetraCrsMatrix(const Epetra_CrsMatrix& epetraCrsMatrix_, const Teuchos::RCP<const Teuchos::Comm<int> >& commT_) { //get row map of Epetra::CrsMatrix & convert to Tpetra::Map auto tpetraRowMap_ = EpetraMap_To_TpetraMap(epetraCrsMatrix_.RowMap(), commT_); //get col map of Epetra::CrsMatrix & convert to Tpetra::Map auto tpetraColMap_ = EpetraMap_To_TpetraMap(epetraCrsMatrix_.ColMap(), commT_); //get CrsGraph of Epetra::CrsMatrix & convert to Tpetra::CrsGraph const Epetra_CrsGraph epetraCrsGraph_ = epetraCrsMatrix_.Graph(); std::size_t maxEntries = epetraCrsGraph_.GlobalMaxNumIndices(); Teuchos::RCP<Tpetra_CrsGraph> tpetraCrsGraph_ = Teuchos::rcp(new Tpetra_CrsGraph(tpetraRowMap_, tpetraColMap_, maxEntries)); for (LO i=0; i<epetraCrsGraph_.NumMyRows(); i++) { LO NumEntries; LO *Indices; epetraCrsGraph_.ExtractMyRowView(i, NumEntries, Indices); tpetraCrsGraph_->insertLocalIndices(i, NumEntries, Indices); } tpetraCrsGraph_->fillComplete(); //convert Epetra::CrsMatrix to Tpetra::CrsMatrix, after creating Tpetra::CrsMatrix based on above Tpetra::CrsGraph Teuchos::RCP<Tpetra_CrsMatrix> tpetraCrsMatrix_ = Teuchos::rcp(new Tpetra_CrsMatrix(tpetraCrsGraph_)); tpetraCrsMatrix_->setAllToScalar(0.0); for (LO i=0; i<epetraCrsMatrix_.NumMyRows(); i++) { LO NumEntries; LO *Indices; ST *Values; epetraCrsMatrix_.ExtractMyRowView(i, NumEntries, Values, Indices); tpetraCrsMatrix_->replaceLocalValues(i, NumEntries, Values, Indices); } tpetraCrsMatrix_->fillComplete(); return tpetraCrsMatrix_; }
/* Debugging utility to check if columns have been Inserted into the * matrix that do not correspond to a row on any processor */ void check_for_rogue_columns( Epetra_CrsMatrix& mat) { // Set up rowVector of 0s and column vector of 1s const Epetra_Map& rowMap = mat.RowMap(); const Epetra_Map& colMap = mat.ColMap(); Epetra_Vector rowVec(rowMap); rowVec.PutScalar(0.0); Epetra_Vector colVec(colMap); colVec.PutScalar(1.0); Epetra_Import importer(colMap, rowMap); // Overwrite colVec 1s with rowVec 0s colVec.Import(rowVec, importer, Insert); // Check that all 1s have been overwritten double nrm=0.0; colVec.Norm1(&nrm); // nrm = number of columns not overwritten by rows // If any rogue columns, exit now (or just get nans later) if (nrm>=1.0) { *out << "ERROR: Column map has " << nrm << " rogue entries that are not associated with any row." << endl; rowMap.Comm().Barrier(); exit(-3); } }
//TpetraCrsMatrix_To_EpetraCrsMatrix: copies Tpetra::CrsMatrix object into its analogous //Epetra_CrsMatrix object void Petra::TpetraCrsMatrix_To_EpetraCrsMatrix(const Teuchos::RCP<const Tpetra_CrsMatrix>& tpetraCrsMatrix_, Epetra_CrsMatrix& epetraCrsMatrix_, const Teuchos::RCP<const Epetra_Comm>& comm_) { //check if row maps of epetraCrsMatrix_ and tpetraCrsMatrix_ are the same const Epetra_BlockMap epetraRowMap_ = epetraCrsMatrix_.RowMap(); Teuchos::RCP<const Tpetra_Map> tpetraRowMap_ = tpetraCrsMatrix_->getRowMap(); Teuchos::RCP<const Epetra_Map> tpetraRowMapE_ = TpetraMap_To_EpetraMap(tpetraRowMap_, comm_); bool isRowSame = tpetraRowMapE_->SameAs(epetraRowMap_); //if epetraCrsMatrix_ and tpetraCrsMatrix_ do not have the same row map, throw an exception if (isRowSame != true) { EpetraExt::BlockMapToMatrixMarketFile("epetraRowMap.mm", epetraRowMap_); EpetraExt::BlockMapToMatrixMarketFile("tpetraRowMapE.mm", *tpetraRowMapE_); } TEUCHOS_TEST_FOR_EXCEPTION((isRowSame != true), std::logic_error, "Error in Petra::TpetraCrsMatrix_To_EpetraCrsMatrix! Arguments Epetra_CrsMatrix and Tpetra::CrsMatrix do not have same row map." << std::endl) ; //check if column maps of epetraCrsMatrix_ and tpetraCrsMatrix_ are the same const Epetra_BlockMap epetraColMap_ = epetraCrsMatrix_.ColMap(); Teuchos::RCP<const Tpetra_Map> tpetraColMap_ = tpetraCrsMatrix_->getColMap(); Teuchos::RCP<const Epetra_Map> tpetraColMapE_ = TpetraMap_To_EpetraMap(tpetraColMap_, comm_); bool isColSame = tpetraColMapE_->SameAs(epetraColMap_); //if epetraCrsMatrix_ and tpetraCrsMatrix_ do not have the same column map, throw an exception TEUCHOS_TEST_FOR_EXCEPTION((isColSame != true), std::logic_error, "Error in Petra::TpetraCrsMatrix_To_EpetraCrsMatrix! Arguments Epetra_CrsMatrix and Tpetra::CrsMatrix do not have same column map." << std::endl) ; epetraCrsMatrix_.PutScalar(0.0); for (LO i = 0; i<tpetraCrsMatrix_->getNodeNumRows(); i++) { LO NumEntries; const LO *Indices; const ST *Values; tpetraCrsMatrix_->getLocalRowView(i, NumEntries, Values, Indices); epetraCrsMatrix_.ReplaceMyValues(i, NumEntries, Values, Indices); } }
/* Computes the approximate Schur complement for the wide separator using guided probing*/ Teuchos::RCP<Epetra_CrsMatrix> computeSchur_GuidedProbing ( shylu_config *config, shylu_symbolic *ssym, // symbolic structure shylu_data *data, // numeric structure Epetra_Map *localDRowMap ) { int i; double relative_thres = config->relative_threshold; Epetra_CrsMatrix *G = ssym->G.getRawPtr(); Epetra_CrsMatrix *R = ssym->R.getRawPtr(); Epetra_LinearProblem *LP = ssym->LP.getRawPtr(); Amesos_BaseSolver *solver = ssym->Solver.getRawPtr(); Ifpack_Preconditioner *ifSolver = ssym->ifSolver.getRawPtr(); Epetra_CrsMatrix *C = ssym->C.getRawPtr(); // Need to create local G (block diagonal portion) , R, C // Get row map of G Epetra_Map CrMap = C->RowMap(); int *c_rows = CrMap.MyGlobalElements(); int *c_cols = (C->ColMap()).MyGlobalElements(); //int c_totalElems = CrMap.NumGlobalElements(); int c_localElems = CrMap.NumMyElements(); int c_localcolElems = (C->ColMap()).NumMyElements(); Epetra_Map GrMap = G->RowMap(); int *g_rows = GrMap.MyGlobalElements(); //int g_totalElems = GrMap.NumGlobalElements(); int g_localElems = GrMap.NumMyElements(); Epetra_Map RrMap = R->RowMap(); int *r_rows = RrMap.MyGlobalElements(); int *r_cols = (R->ColMap()).MyGlobalElements(); //int r_totalElems = RrMap.NumGlobalElements(); int r_localElems = RrMap.NumMyElements(); int r_localcolElems = (R->ColMap()).NumMyElements(); Epetra_SerialComm LComm; Epetra_Map C_localRMap (-1, c_localElems, c_rows, 0, LComm); Epetra_Map C_localCMap (-1, c_localcolElems, c_cols, 0, LComm); Epetra_Map G_localRMap (-1, g_localElems, g_rows, 0, LComm); Epetra_Map R_localRMap (-1, r_localElems, r_rows, 0, LComm); Epetra_Map R_localCMap (-1, r_localcolElems, r_cols, 0, LComm); //cout << "#local rows" << g_localElems << "#non zero local cols" << c_localcolElems << endl; #ifdef DEBUG cout << "DEBUG MODE" << endl; int nrows = C->RowMap().NumMyElements(); assert(nrows == localDRowMap->NumGlobalElements()); int gids[nrows], gids1[nrows]; C_localRMap.MyGlobalElements(gids); localDRowMap->MyGlobalElements(gids1); cout << "Comparing R's domain map with D's row map" << endl; for (int i = 0; i < nrows; i++) { assert(gids[i] == gids1[i]); } #endif int nentries1, gid; // maxentries is the maximum of all three possible matrices as the arrays // are reused between the three int maxentries = max(C->MaxNumEntries(), R->MaxNumEntries()); maxentries = max(maxentries, G->MaxNumEntries()); double *values1 = new double[maxentries]; double *values2 = new double[maxentries]; double *values3 = new double[maxentries]; int *indices1 = new int[maxentries]; int *indices2 = new int[maxentries]; int *indices3 = new int[maxentries]; //cout << "Creating local matrices" << endl; int err; Epetra_CrsMatrix localC(Copy, C_localRMap, C->MaxNumEntries(), false); for (i = 0; i < c_localElems ; i++) { gid = c_rows[i]; err = C->ExtractGlobalRowCopy(gid, maxentries, nentries1, values1, indices1); assert (err == 0); //if (nentries1 > 0) // TODO: Later //{ err = localC.InsertGlobalValues(gid, nentries1, values1, indices1); assert (err == 0); //} } localC.FillComplete(G_localRMap, C_localRMap); //cout << "Created local C matrix" << endl; Epetra_CrsMatrix localR(Copy, R_localRMap, R->MaxNumEntries(), false); for (i = 0; i < r_localElems ; i++) { gid = r_rows[i]; R->ExtractGlobalRowCopy(gid, maxentries, nentries1, values1, indices1); localR.InsertGlobalValues(gid, nentries1, values1, indices1); } localR.FillComplete(*localDRowMap, R_localRMap); //cout << "Created local R matrix" << endl; // Sbar - Approximate Schur complement Teuchos::RCP<Epetra_CrsMatrix> Sbar = Teuchos::rcp(new Epetra_CrsMatrix( Copy, GrMap, g_localElems)); // Include only the block diagonal elements of G in localG Epetra_CrsMatrix localG(Copy, G_localRMap, G->MaxNumEntries(), false); int cnt, scnt; for (i = 0; i < g_localElems ; i++) { gid = g_rows[i]; G->ExtractGlobalRowCopy(gid, maxentries, nentries1, values1, indices1); cnt = 0; scnt = 0; for (int j = 0 ; j < nentries1 ; j++) { if (G->LRID(indices1[j]) != -1) { values2[cnt] = values1[j]; indices2[cnt++] = indices1[j]; } else { // Add it to Sbar immediately values3[scnt] = values1[j]; indices3[scnt++] = indices1[j]; } } localG.InsertGlobalValues(gid, cnt, values2, indices2); Sbar->InsertGlobalValues(gid, scnt, values3, indices3); } localG.FillComplete(); cout << "Created local G matrix" << endl; int nvectors = 16; ShyLU_Probing_Operator probeop(config, ssym, &localG, &localR, LP, solver, ifSolver, &localC, localDRowMap, nvectors); #ifdef DUMP_MATRICES //ostringstream fnamestr; //fnamestr << "localC" << C->Comm().MyPID() << ".mat"; //string Cfname = fnamestr.str(); //EpetraExt::RowMatrixToMatlabFile(Cfname.c_str(), localC); //Epetra_Map defMapg(-1, g_localElems, 0, localG.Comm()); //EpetraExt::ViewTransform<Epetra_CrsMatrix> * ReIdx_MatTransg = //new EpetraExt::CrsMatrix_Reindex( defMapg ); //Epetra_CrsMatrix t2G = (*ReIdx_MatTransg)( localG ); //ReIdx_MatTransg->fwd(); //EpetraExt::RowMatrixToMatlabFile("localG.mat", t2G); #endif //cout << " totalElems in Schur Complement" << totalElems << endl; //cout << myPID << " localElems" << localElems << endl; // **************** Two collectives here ********************* #ifdef TIMING_OUTPUT Teuchos::Time ftime("setup time"); #endif #ifdef TIMING_OUTPUT Teuchos::Time app_time("setup time"); #endif Teuchos::RCP<Epetra_CrsGraph> lSGraph = Teuchos::RCP<Epetra_CrsGraph> ( new Epetra_CrsGraph(Copy, G_localRMap, maxentries)); if (data->num_compute % config->reset_iter == 0) { int nentries; // size > maxentries as there could be fill // TODO: Currently the size of the two arrays can be one, Even if we switch // the loop below the size of the array required is nvectors. Fix it double *values = new double[g_localElems]; int *indices = new int[g_localElems]; double *vecvalues; int dropped = 0; double *maxvalue = new double[nvectors]; #ifdef TIMING_OUTPUT ftime.start(); #endif int findex = g_localElems / nvectors ; int cindex; // int mypid = C->Comm().MyPID(); // unused Epetra_MultiVector probevec(G_localRMap, nvectors); Epetra_MultiVector Scol(G_localRMap, nvectors); for (i = 0 ; i < findex*nvectors ; i+=nvectors) { probevec.PutScalar(0.0); // TODO: Move it out for (int k = 0; k < nvectors; k++) { cindex = k+i; // TODO: Can do better than this, just need to go to the column map // of C, there might be null columns in C // Not much of use for Shasta 2x2 .. Later. probevec.ReplaceGlobalValue(g_rows[cindex], k, 1.0); //if (mypid == 0) //cout << "Changing row to 1.0 " << g_rows[cindex] << endl; } #ifdef TIMING_OUTPUT app_time.start(); #endif probeop.Apply(probevec, Scol); #ifdef TIMING_OUTPUT app_time.stop(); #endif Scol.MaxValue(maxvalue); for (int k = 0; k < nvectors; k++) //TODO:Need to switch these loops { cindex = k+i; vecvalues = Scol[k]; //cout << "MAX" << maxvalue << endl; for (int j = 0 ; j < g_localElems ; j++) { nentries = 0; // inserting one entry in each row for now if (g_rows[cindex] == g_rows[j]) // diagonal entry { values[nentries] = vecvalues[j]; indices[nentries] = g_rows[cindex]; nentries++; err = Sbar->InsertGlobalValues(g_rows[j], nentries, values, indices); assert(err >= 0); err = lSGraph->InsertGlobalIndices(g_rows[j], nentries, indices); assert(err >= 0); } else if (abs(vecvalues[j]/maxvalue[k]) > relative_thres) { values[nentries] = vecvalues[j]; indices[nentries] = g_rows[cindex]; nentries++; err = Sbar->InsertGlobalValues(g_rows[j], nentries, values, indices); assert(err >= 0); err = lSGraph->InsertGlobalIndices(g_rows[j], nentries, indices); assert(err >= 0); } else { if (vecvalues[j] != 0.0) { dropped++; //cout << "vecvalues[j]" << vecvalues[j] << // " max" << maxvalue[k] << endl; } } } } } probeop.ResetTempVectors(1); for ( ; i < g_localElems ; i++) { // TODO: Can move the next two decalarations outside the loop Epetra_MultiVector probevec(G_localRMap, 1); Epetra_MultiVector Scol(G_localRMap, 1); probevec.PutScalar(0.0); // TODO: Can do better than this, just need to go to the column map // of C, there might be null columns in C probevec.ReplaceGlobalValue(g_rows[i], 0, 1.0); #ifdef TIMING_OUTPUT app_time.start(); #endif probeop.Apply(probevec, Scol); #ifdef TIMING_OUTPUT app_time.stop(); #endif vecvalues = Scol[0]; Scol.MaxValue(maxvalue); //cout << "MAX" << maxvalue << endl; for (int j = 0 ; j < g_localElems ; j++) { nentries = 0; // inserting one entry in each row for now if (g_rows[i] == g_rows[j]) // diagonal entry { values[nentries] = vecvalues[j]; indices[nentries] = g_rows[i]; nentries++; err = Sbar->InsertGlobalValues(g_rows[j], nentries, values, indices); assert(err >= 0); err = lSGraph->InsertGlobalIndices(g_rows[j], nentries, indices); assert(err >= 0); } else if (abs(vecvalues[j]/maxvalue[0]) > relative_thres) { values[nentries] = vecvalues[j]; indices[nentries] = g_rows[i]; nentries++; err = Sbar->InsertGlobalValues(g_rows[j], nentries, values, indices); assert(err >= 0); err = lSGraph->InsertGlobalIndices(g_rows[j], nentries, indices); assert(err >= 0); } else { if (vecvalues[j] != 0.0) dropped++; } } } #ifdef TIMING_OUTPUT ftime.stop(); cout << "Time in finding and dropping entries" << ftime.totalElapsedTime() << endl; ftime.reset(); #endif #ifdef TIMING_OUTPUT cout << "Time in Apply of probing" << app_time.totalElapsedTime() << endl; #endif probeop.PrintTimingInfo(); Sbar->FillComplete(); lSGraph->FillComplete(); data->localSbargraph = lSGraph; #ifdef DUMP_MATRICES Epetra_Map defMap2(-1, g_localElems, 0, C->Comm()); EpetraExt::ViewTransform<Epetra_CrsMatrix> * ReIdx_MatTrans2 = new EpetraExt::CrsMatrix_Reindex( defMap2 ); Epetra_CrsMatrix t2S = (*ReIdx_MatTrans2)( *Sbar ); ReIdx_MatTrans2->fwd(); EpetraExt::RowMatrixToMatlabFile("Schur.mat", t2S); #endif cout << "#dropped entries" << dropped << endl; delete[] values; delete[] indices; delete[] maxvalue; } else { if (((data->num_compute-1) % config->reset_iter) == 0) { // We recomputed the Schur complement with dropping for the last // compute. Reset the prober with the new orthogonal vectors for // the Sbar from the previous iteration. Teuchos::ParameterList pList; Teuchos::RCP<Isorropia::Epetra::Prober> gprober = Teuchos::RCP<Isorropia::Epetra::Prober> (new Isorropia::Epetra::Prober( data->localSbargraph.getRawPtr(), pList, false)); gprober->color(); data->guided_prober = gprober; } // Use the prober to probe the probeop for the sparsity pattern // add that to Sbar and call Fill complete int nvectors = data->guided_prober->getNumOrthogonalVectors(); cout << "Number of Orthogonal Vectors for guided probing" << nvectors << endl; probeop.ResetTempVectors(nvectors); Teuchos::RCP<Epetra_CrsMatrix> blockdiag_Sbar = data->guided_prober->probe(probeop); int maxentries = blockdiag_Sbar->GlobalMaxNumEntries(); int *indices = new int[maxentries]; double *values = new double[maxentries]; int numentries; for (int i = 0; i < blockdiag_Sbar->NumGlobalRows() ; i++) { int gid = blockdiag_Sbar->GRID(i); blockdiag_Sbar->ExtractGlobalRowCopy(gid, maxentries, numentries, values, indices); Sbar->InsertGlobalValues(gid, numentries, values, indices); } Sbar->FillComplete(); delete[] indices; delete[] values; } delete[] values1; delete[] indices1; delete[] values2; delete[] indices2; delete[] values3; delete[] indices3; return Sbar; }
// // Diagonal: 0=no change, 1=eliminate entry // from the map for the largest row element in process 0 // 2=add diagonal entries to the matrix, with a zero value // (assume row map contains all diagonal entries). // // ReindexRowMap: // 0=no change, 1= add 2 (still contiguous), 2=non-contiguous // // ReindexColMap // 0=same as RowMap, 1=add 4 - Different From RowMap, but contiguous) // // RangeMap: // 0=no change, 1=serial map, 2=bizarre distribution, 3=replicated map // // DomainMap: // 0=no change, 1=serial map, 2=bizarre distribution, 3=replicated map // RCP<Epetra_CrsMatrix> NewMatNewMap(Epetra_CrsMatrix& In, int Diagonal, int ReindexRowMap, int ReindexColMap, int RangeMapType, int DomainMapType ) { // // If we are making no change, return the original matrix (which has a linear map) // #if 0 std::cout << __FILE__ << "::" << __LINE__ << " " << Diagonal << " " << ReindexRowMap << " " << ReindexColMap << " " << RangeMapType << " " << DomainMapType << " " << std::endl ; #endif if ( Diagonal + ReindexRowMap + ReindexColMap + RangeMapType + DomainMapType == 0 ) { RCP<Epetra_CrsMatrix> ReturnOrig = rcp( &In, false ); return ReturnOrig ; } // // Diagonal==2 is used for a different purpose - // Making sure that the diagonal of the matrix is non-empty. // Note: The diagonal must exist in In.RowMap(). // if ( Diagonal == 2 ) { assert( ReindexRowMap==0 && ReindexColMap == 0 ) ; } int (*RowPermute)(int in) = 0; int (*ColPermute)(int in) = 0; assert( Diagonal >= 0 && Diagonal <= 2 ); assert( ReindexRowMap>=0 && ReindexRowMap<=2 ); assert( ReindexColMap>=0 && ReindexColMap<=1 ); assert( RangeMapType>=0 && RangeMapType<=3 ); assert( DomainMapType>=0 && DomainMapType<=3 ); Epetra_Map DomainMap = In.DomainMap(); Epetra_Map RangeMap = In.RangeMap(); Epetra_Map ColMap = In.ColMap(); Epetra_Map RowMap = In.RowMap(); int NumMyRowElements = RowMap.NumMyElements(); int NumMyColElements = ColMap.NumMyElements(); int NumMyRangeElements = RangeMap.NumMyElements(); int NumMyDomainElements = DomainMap.NumMyElements(); int NumGlobalRowElements = RowMap.NumGlobalElements(); int NumGlobalColElements = ColMap.NumGlobalElements(); int NumGlobalRangeElements = RangeMap.NumGlobalElements(); int NumGlobalDomainElements = DomainMap.NumGlobalElements(); assert( NumGlobalRangeElements == NumGlobalDomainElements ) ; std::vector<int> MyGlobalRowElements( NumMyRowElements ) ; std::vector<int> NumEntriesPerRow( NumMyRowElements ) ; std::vector<int> MyPermutedGlobalRowElements( NumMyRowElements ) ; std::vector<int> MyGlobalColElements( NumMyColElements ) ; std::vector<int> MyPermutedGlobalColElements( NumMyColElements ) ; // Used to create the column map std::vector<int> MyPermutedGlobalColElementTable( NumMyColElements ) ; // To convert local indices to global std::vector<int> MyGlobalRangeElements( NumMyRangeElements ) ; std::vector<int> MyPermutedGlobalRangeElements( NumMyRangeElements ) ; std::vector<int> MyGlobalDomainElements( NumMyDomainElements ) ; std::vector<int> MyPermutedGlobalDomainElements( NumMyDomainElements ) ; RowMap.MyGlobalElements(&MyGlobalRowElements[0]); ColMap.MyGlobalElements(&MyGlobalColElements[0]); RangeMap.MyGlobalElements(&MyGlobalRangeElements[0]); DomainMap.MyGlobalElements(&MyGlobalDomainElements[0]); switch( ReindexRowMap ) { case 0: RowPermute = &NoPermute ; break; case 1: RowPermute = &SmallRowPermute ; break; case 2: RowPermute = BigRowPermute ; break; } switch( ReindexColMap ) { case 0: ColPermute = RowPermute ; break; case 1: ColPermute = &SmallColPermute ; break; } // // Create Serial Range and Domain Maps based on the permuted indexing // int nlocal = 0; if (In.Comm().MyPID()==0) nlocal = NumGlobalRangeElements; std::vector<int> AllIDs( NumGlobalRangeElements ) ; for ( int i = 0; i < NumGlobalRangeElements ; i++ ) AllIDs[i] = (*RowPermute)( i ) ; Epetra_Map SerialRangeMap( -1, nlocal, &AllIDs[0], 0, In.Comm()); std::vector<int> AllIDBs( NumGlobalRangeElements ) ; for ( int i = 0; i < NumGlobalRangeElements ; i++ ) AllIDBs[i] = (*ColPermute)( i ) ; Epetra_Map SerialDomainMap( -1, nlocal, &AllIDBs[0], 0, In.Comm()); // // Create Bizarre Range and Domain Maps based on the permuted indexing // These are nearly serial, having all but one element on process 0 // The goal here is to make sure that we can use Domain and Range maps // that are neither serial, nor distributed in the normal manner. // std::vector<int> AllIDCs( NumGlobalRangeElements ) ; for ( int i = 0; i < NumGlobalRangeElements ; i++ ) AllIDCs[i] = (*ColPermute)( i ) ; if ( In.Comm().NumProc() > 1 ) { if (In.Comm().MyPID()==0) nlocal = NumGlobalRangeElements-1; if (In.Comm().MyPID()==1) { nlocal = 1; AllIDCs[0] = (*ColPermute)( NumGlobalRangeElements - 1 ); } } int iam = In.Comm().MyPID(); Epetra_Map BizarreDomainMap( -1, nlocal, &AllIDCs[0], 0, In.Comm()); std::vector<int> AllIDDs( NumGlobalRangeElements ) ; for ( int i = 0; i < NumGlobalRangeElements ; i++ ) AllIDDs[i] = (*RowPermute)( i ) ; if ( In.Comm().NumProc() > 1 ) { if (In.Comm().MyPID()==0) nlocal = NumGlobalRangeElements-1; if (In.Comm().MyPID()==1) { nlocal = 1; AllIDDs[0] = (*RowPermute)( NumGlobalRangeElements -1 ) ; } } Epetra_Map BizarreRangeMap( -1, nlocal, &AllIDDs[0], 0, In.Comm()); // // Compute the column map // // If Diagonal==1, remove the column corresponding to the last row owned // by process 0. Removing this column from a tridiagonal matrix, leaves // a disconnected, but non-singular matrix. // int NumMyColElementsOut = 0 ; int NumGlobalColElementsOut ; if ( Diagonal == 1 ) NumGlobalColElementsOut = NumGlobalColElements-1; else NumGlobalColElementsOut = NumGlobalColElements; if ( Diagonal == 1 && iam==0 ) { for ( int i=0; i < NumMyColElements ; i++ ) { if ( MyGlobalColElements[i] != MyGlobalRowElements[NumMyRowElements-1] ) { MyPermutedGlobalColElements[NumMyColElementsOut++] = (*ColPermute)( MyGlobalColElements[i] ) ; } } assert( NumMyColElementsOut == NumMyColElements-1 ); } else { for ( int i=0; i < NumMyColElements ; i++ ) MyPermutedGlobalColElements[i] = (*ColPermute)( MyGlobalColElements[i] ) ; NumMyColElementsOut = NumMyColElements ; if ( Diagonal == 2 ) { // For each row, make sure that the column map has this row in it, // if it doesn't, add it to the column map. // Note: MyPermutedGlobalColElements == MyGlobalColElements when // Diagonal==2 because ( Diagonal == 2 ) implies: // ReindexRowMap==0 && ReindexColMap == 0 - see assert above for ( int i=0; i < NumMyRowElements ; i++ ) { bool MissingDiagonal = true; for ( int j=0; j < NumMyColElements; j++ ) { if ( MyGlobalRowElements[i] == MyGlobalColElements[j] ) { MissingDiagonal = false; } } if ( MissingDiagonal ) { MyPermutedGlobalColElements.resize(NumMyColElements+1); MyPermutedGlobalColElements[NumMyColElementsOut] = MyGlobalRowElements[i]; NumMyColElementsOut++; } } In.Comm().SumAll(&NumMyColElementsOut,&NumGlobalColElementsOut,1); } } // // These tables are used both as the permutation tables and to create the maps. // for ( int i=0; i < NumMyColElements ; i++ ) MyPermutedGlobalColElementTable[i] = (*ColPermute)( MyGlobalColElements[i] ) ; for ( int i=0; i < NumMyRowElements ; i++ ) MyPermutedGlobalRowElements[i] = (*RowPermute)( MyGlobalRowElements[i] ) ; for ( int i=0; i < NumMyRangeElements ; i++ ) MyPermutedGlobalRangeElements[i] = (*RowPermute)( MyGlobalRangeElements[i] ) ; for ( int i=0; i < NumMyDomainElements ; i++ ) MyPermutedGlobalDomainElements[i] = (*ColPermute)( MyGlobalDomainElements[i] ) ; RCP<Epetra_Map> PermutedRowMap = rcp( new Epetra_Map( NumGlobalRowElements, NumMyRowElements, &MyPermutedGlobalRowElements[0], 0, In.Comm() ) ); RCP<Epetra_Map> PermutedColMap = rcp( new Epetra_Map( NumGlobalColElementsOut, NumMyColElementsOut, &MyPermutedGlobalColElements[0], 0, In.Comm() ) ); RCP<Epetra_Map> PermutedRangeMap = rcp( new Epetra_Map( NumGlobalRangeElements, NumMyRangeElements, &MyPermutedGlobalRangeElements[0], 0, In.Comm() ) ); RCP<Epetra_Map> PermutedDomainMap = rcp( new Epetra_Map( NumGlobalDomainElements, NumMyDomainElements, &MyPermutedGlobalDomainElements[0], 0, In.Comm() ) ); // // These vectors are filled and then passed to InsertGlobalValues // std::vector<int> ThisRowIndices( In.MaxNumEntries() ); std::vector<double> ThisRowValues( In.MaxNumEntries() ); std::vector<int> PermutedGlobalColIndices( In.MaxNumEntries() ); //std::cout << __FILE__ << "::" <<__LINE__ << std::endl ; RCP<Epetra_CrsMatrix> Out = rcp( new Epetra_CrsMatrix( Copy, *PermutedRowMap, *PermutedColMap, 0 ) ); for (int i=0; i<NumMyRowElements; i++) { int NumIndicesThisRow = 0; assert( In.ExtractMyRowCopy( i, In.MaxNumEntries(), NumIndicesThisRow, &ThisRowValues[0], &ThisRowIndices[0] ) == 0 ) ; for (int j = 0 ; j < NumIndicesThisRow ; j++ ) { PermutedGlobalColIndices[j] = MyPermutedGlobalColElementTable[ ThisRowIndices[j] ] ; } bool MissingDiagonal = false; if ( Diagonal==2 ) { // assert( MyGlobalRowElements[i] == MyPermutedGlobalRowElements[i] ); MissingDiagonal = true; for( int j =0 ; j < NumIndicesThisRow ; j++ ) { if ( PermutedGlobalColIndices[j] == MyPermutedGlobalRowElements[i] ) { MissingDiagonal = false ; } } #if 0 std::cout << __FILE__ << "::" << __LINE__ << " i = " << i << " MyPermutedGlobalRowElements[i] = " << MyPermutedGlobalRowElements[i] << " MissingDiagonal = " << MissingDiagonal << std::endl ; #endif } if ( MissingDiagonal ) { ThisRowValues.resize(NumIndicesThisRow+1) ; ThisRowValues[NumIndicesThisRow] = 0.0; PermutedGlobalColIndices.resize(NumIndicesThisRow+1); PermutedGlobalColIndices[NumIndicesThisRow] = MyPermutedGlobalRowElements[i] ; #if 0 std::cout << __FILE__ << "::" << __LINE__ << " i = " << i << "NumIndicesThisRow = " << NumIndicesThisRow << "ThisRowValues[NumIndicesThisRow = " << ThisRowValues[NumIndicesThisRow] << " PermutedGlobalColIndices[NumIndcesThisRow] = " << PermutedGlobalColIndices[NumIndicesThisRow] << std::endl ; #endif NumIndicesThisRow++ ; } assert( Out->InsertGlobalValues( MyPermutedGlobalRowElements[i], NumIndicesThisRow, &ThisRowValues[0], &PermutedGlobalColIndices[0] ) >= 0 ); } // Epetra_LocalMap ReplicatedMap( NumGlobalRangeElements, 0, In.Comm() ); RCP<Epetra_Map> OutRangeMap ; RCP<Epetra_Map> OutDomainMap ; switch( RangeMapType ) { case 0: OutRangeMap = PermutedRangeMap ; break; case 1: OutRangeMap = rcp(&SerialRangeMap, false); break; case 2: OutRangeMap = rcp(&BizarreRangeMap, false); break; case 3: OutRangeMap = rcp(&ReplicatedMap, false); break; } // switch( DomainMapType ) { switch( DomainMapType ) { case 0: OutDomainMap = PermutedDomainMap ; break; case 1: OutDomainMap = rcp(&SerialDomainMap, false); break; case 2: OutDomainMap = rcp(&BizarreDomainMap, false); break; case 3: OutDomainMap = rcp(&ReplicatedMap, false); break; } #if 0 assert(Out->FillComplete( *PermutedDomainMap, *PermutedRangeMap )==0); #else assert(Out->FillComplete( *OutDomainMap, *OutRangeMap )==0); #endif #if 0 std::cout << __FILE__ << "::" << __LINE__ << std::endl ; Out->Print( std::cout ) ; #endif return Out; }
bool FiniteDifferenceColoringWithUpdate::differenceProbe(const Epetra_Vector& x, Epetra_CrsMatrix& jac,const Epetra_MapColoring& colors){ // Allocate space for perturbation, get column version of x for scaling Epetra_Vector xp(x); Epetra_Vector *xcol; int N=jac.NumMyRows(); if(jac.ColMap().SameAs(x.Map())) xcol=const_cast<Epetra_Vector*>(&x); else{ xcol=new Epetra_Vector(jac.ColMap(),true);//zeros out by default xcol->Import(x,*jac.Importer(),InsertAdd); } // Counters for probing diagnostics double tmp,probing_error_lower_bound=0.0,jc_norm=0.0; // Grab coloring info (being very careful to ignore color 0) int Ncolors=colors.MaxNumColors()+1; int num_c0_global,num_c0_local=colors.NumElementsWithColor(0); colors.Comm().MaxAll(&num_c0_local,&num_c0_global,1); if(num_c0_global>0) Ncolors--; if(Ncolors==0) return false; // Pointers for Matrix Info int entries, *indices; double *values; // NTS: Fix me if ( diffType == Centered ) exit(1); double scaleFactor = 1.0; if ( diffType == Backward ) scaleFactor = -1.0; // Compute RHS at initial solution computeF(x,fo,NOX::Epetra::Interface::Required::FD_Res); /* Probing, vector by vector since computeF does not have a MultiVector interface */ // Assume that anything with Color 0 gets ignored. for(int j=1;j<Ncolors;j++){ xp=x; for(int i=0;i<N;i++){ if(colors[i]==j) xp[i] += scaleFactor*(alpha*abs(x[i])+beta); } computeF(xp, fp, NOX::Epetra::Interface::Required::FD_Res); // Do the subtraction to estimate the Jacobian (w/o including step length) Jc.Update(1.0, fp, -1.0, fo, 0.0); // Relative error in probing if(use_probing_diags){ Jc.Norm2(&tmp); jc_norm+=tmp*tmp; } for(int i=0;i<N;i++){ // Skip for uncolored row/columns, else update entries if(colors[i]==0) continue; jac.ExtractMyRowView(i,entries,values,indices); for(int k=0;k<jac.NumMyEntries(i);k++){ if(colors[indices[k]]==j){ values[k]=Jc[i] / (scaleFactor*(alpha*abs((*xcol)[indices[k]])+beta)); // If probing diagnostics are on, zero out the entries as they are used if(use_probing_diags) Jc[i]=0.0; break;// Only one value per row... } } } if(use_probing_diags){ Jc.Norm2(&tmp); probing_error_lower_bound+=tmp*tmp; } } // If diagnostics are requested, output Frobenius norm lower bound if(use_probing_diags && !x.Comm().MyPID()) printf("Probing Error Lower Bound (Frobenius) abs = %6.4e rel = %6.4e\n",sqrt(probing_error_lower_bound),sqrt(probing_error_lower_bound)/sqrt(jc_norm)); // Cleanup if(!jac.ColMap().SameAs(x.Map())) delete xcol; return true; }
int TPackAndPrepareWithOwningPIDs(const Epetra_CrsMatrix & A, int NumExportIDs, int * ExportLIDs, int & LenExports, char *& Exports, int & SizeOfPacket, int * Sizes, bool & VarSizes, std::vector<int>& pids) { int i, j; VarSizes = true; //enable variable block size data comm int TotalSendLength = 0; int * IntSizes = 0; if( NumExportIDs>0 ) IntSizes = new int[NumExportIDs]; int SizeofIntType = sizeof(int_type); for(i = 0; i < NumExportIDs; ++i) { int NumEntries; A.NumMyRowEntries( ExportLIDs[i], NumEntries ); // Will have NumEntries doubles, 2*NumEntries +2 ints pack them interleaved Sizes[i] = NumEntries; // NTS: We add the owning PID as the SECOND int of the pair for each entry Sizes[i] = NumEntries; // NOTE: Mixing and matching Int Types would be more efficient, BUT what about variable alignment? IntSizes[i] = 1 + (((2*NumEntries+2)*SizeofIntType)/(int)sizeof(double)); TotalSendLength += (Sizes[i]+IntSizes[i]); } double * DoubleExports = 0; SizeOfPacket = (int)sizeof(double); //setup buffer locally for memory management by this object if( TotalSendLength*SizeOfPacket > LenExports ) { if( LenExports > 0 ) delete [] Exports; LenExports = TotalSendLength*SizeOfPacket; DoubleExports = new double[TotalSendLength]; for( int i = 0; i < TotalSendLength; ++i ) DoubleExports[i] = 0.0; Exports = (char *) DoubleExports; } int NumEntries; double * values; double * valptr, * dintptr; // Each segment of Exports will be filled by a packed row of information for each row as follows: // 1st int: GRID of row where GRID is the global row ID for the source matrix // next int: NumEntries, Number of indices in row // next 2*NumEntries: The actual indices and owning [1] PID each for the row in (GID,PID) pairs with the GID first. // [1] Owning is defined in the sense of "Who owns the GID in the DomainMap," aka, who sends the GID in the importer const Epetra_Map & rowMap = A.RowMap(); const Epetra_Map & colMap = A.ColMap(); if( NumExportIDs > 0 ) { int_type * Indices; int_type FromRow; int_type * intptr; int maxNumEntries = A.MaxNumEntries(); std::vector<int> MyIndices(maxNumEntries); dintptr = (double *) Exports; valptr = dintptr + IntSizes[0]; intptr = (int_type *) dintptr; for (i=0; i<NumExportIDs; i++) { FromRow = (int_type) rowMap.GID64(ExportLIDs[i]); intptr[0] = FromRow; values = valptr; Indices = intptr + 2; EPETRA_CHK_ERR(A.ExtractMyRowCopy(ExportLIDs[i], maxNumEntries, NumEntries, values, &MyIndices[0])); for (j=0; j<NumEntries; j++) { Indices[2*j] = (int_type) colMap.GID64(MyIndices[j]); // convert to GIDs Indices[2*j+1] = pids[MyIndices[j]]; // PID owning the entry. } intptr[1] = NumEntries; // Load second slot of segment if( i < (NumExportIDs-1) ) { dintptr += (IntSizes[i]+Sizes[i]); valptr = dintptr + IntSizes[i+1]; intptr = (int_type *) dintptr; } } for(i = 0; i < NumExportIDs; ++i ) Sizes[i] += IntSizes[i]; } if( IntSizes ) delete [] IntSizes; return(0); }