int extract_matrices ( Epetra_CrsMatrix *A, // i/p: A matrix shylu_symbolic *ssym, // symbolic structure shylu_data *data, // numeric structure, TODO: Required ? shylu_config *config, // i/p: library configuration bool insertValues // true implies values will be inserted and fill // complete will be called. false implies values // will be replaced. ) { Teuchos::RCP<Epetra_CrsMatrix> D = ssym->D; Teuchos::RCP<Epetra_CrsMatrix> C = ssym->C; Teuchos::RCP<Epetra_CrsMatrix> R = ssym->R; Teuchos::RCP<Epetra_CrsMatrix> G = ssym->G; Teuchos::RCP<Epetra_CrsGraph> Sg = ssym->Sg; int *DColElems = data->DColElems; int *gvals = data->gvals; double Sdiagfactor = config->Sdiagfactor; int *LeftIndex = new int[data->lmax]; double *LeftValues = new double[data->lmax]; int *RightIndex = new int[data->rmax]; double *RightValues = new double[data->rmax]; int err; int lcnt, rcnt ; int gcid; int gid; int *Ai; double *Ax; int nrows = A->RowMap().NumMyElements(); int *rows = A->RowMap().MyGlobalElements(); for (int i = 0; i < nrows ; i++) { int NumEntries; err = A->ExtractMyRowView(i, NumEntries, Ax, Ai); lcnt = 0; rcnt = 0; // Place the entry in the correct sub matrix, Works only for sym gid = rows[i]; int lcid; for (int j = 0 ; j < NumEntries ; j++) { // O(nnz) ! Careful what you do inside // Row permutation does not matter here gcid = A->GCID(Ai[j]); assert(gcid != -1); //Either in D or R if ((gvals[gid] != 1 && gvals[gcid] == 1) || (gvals[gid] == 1 && A->LRID(gcid) != -1 && gvals[gcid] == 1)) { assert(lcnt < data->lmax); if (insertValues) LeftIndex[lcnt] = gcid; else { //local column id lcid = (gvals[gid] == 1 ? D->LCID(gcid) : R->LCID(gcid)); assert(lcid != -1); LeftIndex[lcnt] = lcid; } LeftValues[lcnt++] = Ax[j]; } else { assert(rcnt < data->rmax); if (insertValues) RightIndex[rcnt] = gcid; else { //local column id lcid = (gvals[gid] == 1 ? C->LCID(gcid) : G->LCID(gcid)); assert(lcid != -1); RightIndex[rcnt] = lcid; } RightValues[rcnt++] = Ax[j]; } } if (gvals[gid] == 1) { // D or C row if (insertValues) { err = D->InsertGlobalValues(gid, lcnt, LeftValues, LeftIndex); assert(err == 0); err = C->InsertGlobalValues(gid, rcnt, RightValues, RightIndex); assert(err == 0); } else { err = D->ReplaceMyValues(D->LRID(gid), lcnt, LeftValues, LeftIndex); assert(err == 0); err = C->ReplaceMyValues(C->LRID(gid), rcnt, RightValues, RightIndex); assert(err == 0); } } else { // R or S row //assert(lcnt > 0); // TODO: Enable this once using narrow sep. if (insertValues) { assert(rcnt > 0); err = R->InsertGlobalValues(gid, lcnt, LeftValues, LeftIndex); assert(err == 0); err = G->InsertGlobalValues(gid, rcnt, RightValues, RightIndex); assert(err == 0); if (config->schurApproxMethod == 1) { err = Sg->InsertGlobalIndices(gid, rcnt, RightIndex); assert(err == 0); } } else { assert(rcnt > 0); err = R->ReplaceMyValues(R->LRID(gid), lcnt, LeftValues, LeftIndex); assert(err == 0); err = G->ReplaceMyValues(G->LRID(gid), rcnt, RightValues, RightIndex); assert(err == 0); } } } if (insertValues) { /* ------------- Create the maps for the DBBD form ------------------ */ Epetra_Map *DRowMap, *SRowMap, *DColMap; Epetra_SerialComm LComm; if (config->sep_type != 1) { DRowMap = new Epetra_Map(-1, data->Dnr, data->DRowElems, 0, A->Comm()); SRowMap = new Epetra_Map(-1, data->Snr, data->SRowElems, 0, A->Comm()); DColMap = new Epetra_Map(-1, data->Dnc, DColElems, 0, A->Comm()); } else { DRowMap = new Epetra_Map(-1, data->Dnr, data->DRowElems, 0, LComm); SRowMap = new Epetra_Map(-1, data->Snr, data->SRowElems, 0, LComm); DColMap = new Epetra_Map(-1, data->Dnc, DColElems, 0, LComm); } D->FillComplete(); //config->dm.print(5, "Done D fillcomplete"); G->FillComplete(); //config->dm.print(5, "Done G fillcomplete"); C->FillComplete(*SRowMap, *DRowMap); //TODO:Won't work if permutation is // unsymmetric SRowMap //config->dm.print(5, "Done C fillcomplete"); R->FillComplete(*DColMap, *SRowMap); //config->dm.print(5, "Done R fillcomplete"); int Sdiag = (int) data->Snr * Sdiagfactor; Sdiag = MIN(Sdiag, data->Snr-1); Sdiag = MAX(Sdiag, 0); // Add the diagonals to Sg for (int i = 0; config->schurApproxMethod == 1 && i < nrows ; i++) { gid = rows[i]; if (gvals[gid] == 1) continue; // not a row in S if (data->Snr == 0) assert(0 == 1); rcnt = 0; //TODO Will be trouble if SNumGlobalCols != Snc //assert(SNumGlobalCols == Snc); //for (int j = MAX(i-Sdiag,0) ; j<MIN(SNumGlobalCols, i+Sdiag); j++) for (int j = MAX(i-Sdiag, 0) ; j < MIN(data->Snr, i+Sdiag); j++) { // find the adjacent columns from the row map of S //assert (j >= 0 && j < Snr); RightIndex[rcnt++] = data->SRowElems[j]; } err = Sg->InsertGlobalIndices(gid, rcnt, RightIndex); assert(err == 0); // Always insert the diagonals, if it is added twice that is fine. err = Sg->InsertGlobalIndices(gid, 1, &gid); assert(err == 0); } if (config->schurApproxMethod == 1) Sg->FillComplete(); delete DRowMap; delete SRowMap; delete DColMap; } #if 0 if (insertValues) { #ifdef TIMING_OUTPUT Teuchos::Time ttime("transpose time"); ttime.start(); #endif bool MakeDataContiguous = true; ssym->transposer = Teuchos::RCP<EpetraExt::RowMatrix_Transpose>(new EpetraExt::RowMatrix_Transpose(MakeDataContiguous)); ssym->DT = Teuchos::rcp( dynamic_cast<Epetra_CrsMatrix *>(&(*ssym->transposer)(*D)), false); #ifdef TIMING_OUTPUT ttime.stop(); cout << "Transpose Time" << ttime.totalElapsedTime() << endl; ttime.reset(); #endif } else { ssym->transposer->fwd(); //ssym->ReIdx_LP->fwd(); // TODO: Needed ? } #endif // A is no longer needed delete[] LeftIndex; delete[] LeftValues; delete[] RightIndex; delete[] RightValues; //cout << msg << "S rows=" << S.NumGlobalRows() << " S cols=" << //S.NumGlobalCols() << "#cols in column map="<< //S.ColMap().NumMyElements() << endl; //cout << msg << "C rows=" << Cptr->NumGlobalRows() << " C cols=" << //Cptr->NumGlobalCols() << "#cols in column map="<< //Cptr->ColMap().NumMyElements() << endl; //cout << msg << "D rows=" << D.NumGlobalRows() << " D cols=" << //D.NumGlobalCols() << "#cols in column map="<< //D.ColMap().NumMyElements() << endl; //cout << msg << "R rows=" << Rptr->NumGlobalRows() << " R cols=" << //Rptr->NumGlobalCols() << "#cols in column map="<< //Rptr->ColMap().NumMyElements() << endl; // ] return 0; }