/*----------------------------------------------------------------------* | m.gee 11/05| *----------------------------------------------------------------------*/ Epetra_CrsMatrix* ML_NOX::StripZeros(Epetra_CrsMatrix& in, double eps) { Epetra_CrsMatrix* out = new Epetra_CrsMatrix(Copy,in.RowMap(),200); for (int lrow=0; lrow<in.NumMyRows(); ++lrow) { int grow = in.GRID(lrow); if (grow<0) { cout << "ERROR: grow<0 \n"; exit(0); } int numentries; int* lindices; double* values; int err = in.ExtractMyRowView(lrow,numentries,values,lindices); if (err) { cout << "ExtractMyRowView returned " << err << endl; exit(0);} for (int j=0; j<numentries; ++j) { int lcol = lindices[j]; int gcol = in.GCID(lcol); if (gcol<0) { cout << "ERROR: gcol<0 \n"; exit(0); } if (abs(values[j])<eps && gcol != grow) continue; int err = out->InsertGlobalValues(grow,1,&values[j],&gcol); if (err) { cout << "InsertGlobalValues returned " << err << endl; exit(0);} } } out->FillComplete(); out->OptimizeStorage(); return out; }
Epetra_CrsMatrix* TridiagMatrix(const Epetra_Map* Map, const double a, const double b, const double c) { Epetra_CrsMatrix* Matrix = new Epetra_CrsMatrix(Copy, *Map, 3); int NumGlobalElements = Map->NumGlobalElements(); int NumMyElements = Map->NumMyElements(); int* MyGlobalElements = Map->MyGlobalElements(); std::vector<double> Values(2); std::vector<int> Indices(2); int NumEntries; for (int i = 0; i < NumMyElements; ++i) { if (MyGlobalElements[i] == 0) { // off-diagonal for first row Indices[0] = 1; NumEntries = 1; Values[0] = c; } else if (MyGlobalElements[i] == NumGlobalElements - 1) { // off-diagonal for last row Indices[0] = NumGlobalElements - 2; NumEntries = 1; Values[0] = b; } else { // off-diagonal for internal row Indices[0] = MyGlobalElements[i] - 1; Values[1] = b; Indices[1] = MyGlobalElements[i] + 1; Values[0] = c; NumEntries = 2; } Matrix->InsertGlobalValues(MyGlobalElements[i], NumEntries, &Values[0], &Indices[0]); // Put in the diagonal entry Values[0] = a; Matrix->InsertGlobalValues(MyGlobalElements[i], 1, &Values[0], MyGlobalElements + i); } // Finish up, trasforming the matrix entries into local numbering, // to optimize data transfert during matrix-vector products Matrix->FillComplete(); Matrix->OptimizeStorage(); return(Matrix); }
/*----------------------------------------------------------------------* | Create the spd system (public) mwgee 12/05| | Note that this is collective for ALL procs | *----------------------------------------------------------------------*/ Epetra_CrsMatrix* MOERTEL::Manager::MakeSPDProblem() { // time this process Epetra_Time time(Comm()); time.ResetStartTime(); // check whether all interfaces are complete and integrated std::map<int,Teuchos::RCP<MOERTEL::Interface> >::iterator curr; for (curr=interface_.begin(); curr != interface_.end(); ++curr) { if (curr->second->IsComplete() == false) { cout << "***ERR*** MOERTEL::Manager::MakeSPDProblem:\n" << "***ERR*** interface " << curr->second->Id() << " is not Complete()\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; return NULL; } if (curr->second->IsIntegrated() == false) { cout << "***ERR*** MOERTEL::Manager::MakeSPDProblem:\n" << "***ERR*** interface " << curr->second->Id() << " is not integrated yet\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; return NULL; } } // check whether we have a problemmap_ if (problemmap_==Teuchos::null) { cout << "***ERR*** MOERTEL::Manager::MakeSPDProblem:\n" << "***ERR*** No problemmap_ set\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; return NULL; } // check whether we have a constraintsmap_ if (constraintsmap_==Teuchos::null) { cout << "***ERR*** MOERTEL::Manager::MakeSPDProblem:\n" << "***ERR*** onstraintsmap is NULL\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; return NULL; } // check for saddlemap_ if (saddlemap_==Teuchos::null) { cout << "***ERR*** MOERTEL::Manager::MakeSPDProblem:\n" << "***ERR*** saddlemap_==NULL\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; return NULL; } // check for inputmatrix if (inputmatrix_==Teuchos::null) { cout << "***ERR*** MOERTEL::Manager::MakeSPDProblem:\n" << "***ERR*** No inputmatrix set\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; return NULL; } // check whether we have M and D matrices if (D_==Teuchos::null || M_==Teuchos::null) { cout << "***ERR*** MOERTEL::Manager::MakeSPDProblem:\n" << "***ERR*** Matrix M or D is NULL\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; return NULL; } // we need a map from lagrange multiplier dofs to primal dofs on the same node std::vector<MOERTEL::Node*> nodes(0); std::map<int,int> lm_to_dof; for (curr=interface_.begin(); curr!=interface_.end(); ++curr) { Teuchos::RCP<MOERTEL::Interface> inter = curr->second; inter->GetNodeView(nodes); for (int i=0; i<(int)nodes.size(); ++i) { if (!nodes[i]->Nlmdof()) continue; const int* dof = nodes[i]->Dof(); const int* lmdof = nodes[i]->LMDof(); for (int j=0; j<nodes[i]->Nlmdof(); ++j) { //cout << "j " << j << " maps lmdof " << lmdof[j] << " to dof " << dof[j] << endl; lm_to_dof[lmdof[j]] = dof[j]; } } } lm_to_dof_ = Teuchos::rcp(new std::map<int,int>(lm_to_dof)); // this is a very useful map for the moertel_ml_preconditioner /* _ _ | | | Arr Arn Mr | S = | | | Anr Ann D | | | MrT D 0 | |_ _ | _ _ | | | Arr Arn | A = | | | Anr Ann | |_ _| 1) Ann is square and we need it's Range/DomainMap annmap _ _ WT = |_ 0 Dinv _| 2) Build WT (has rowmap/rangemap annmap and domainmap problemmap_) _ _ | | | Mr | B = | | | D | |_ _| 3) Build B (has rowmap/rangemap problemmap_ and domainmap annmap) 4) Build I, the identity matrix with maps problemmap_,problemmap_); */ int err=0; //-------------------------------------------------------------------------- // 1) create the rangemap of Ann std::vector<int> myanngids(problemmap_->NumMyElements()); int count=0; std::map<int,int>::iterator intintcurr; for (intintcurr=lm_to_dof.begin(); intintcurr!=lm_to_dof.end(); ++intintcurr) { if (problemmap_->MyGID(intintcurr->second)==false) continue; if ((int)myanngids.size()<=count) myanngids.resize(myanngids.size()+50); myanngids[count] = intintcurr->second; ++count; } myanngids.resize(count); int numglobalelements; Comm().SumAll(&count,&numglobalelements,1); Epetra_Map* annmap = new Epetra_Map(numglobalelements,count,&myanngids[0],0,Comm()); annmap_ = Teuchos::rcp(annmap); myanngids.clear(); //-------------------------------------------------------------------------- // 2) create WT Epetra_CrsMatrix* WT = new Epetra_CrsMatrix(Copy,*annmap,1,false); for (intintcurr=lm_to_dof.begin(); intintcurr!=lm_to_dof.end(); ++intintcurr) { int lmdof = intintcurr->first; int dof = intintcurr->second; if (D_->MyGRID(lmdof)==false) continue; int lmlrid = D_->LRID(lmdof); int numentries; int* indices; double* values; err = D_->ExtractMyRowView(lmlrid,numentries,values,indices); if (err) cout << "D_->ExtractMyRowView returned err=" << err << endl; bool foundit = false; for (int j=0; j<numentries; ++j) { int gcid = D_->GCID(indices[j]); if (gcid<0) cout << "Cannot find gcid for indices[j]\n"; //cout << "Proc " << Comm().MyPID() << " lmdof " << lmdof << " dof " << dof << " gcid " << gcid << " val " << values[j] << endl; if (gcid==dof) { double val = 1./values[j]; err = WT->InsertGlobalValues(dof,1,&val,&dof); if (err<0) cout << "WT->InsertGlobalValues returned err=" << err << endl; foundit = true; break; } } if (!foundit) { cout << "***ERR*** MOERTEL::Manager::MakeSPDProblem:\n" << "***ERR*** Cannot compute inverse of D_\n" << "***ERR*** file/line: " << __FILE__ << "/" << __LINE__ << "\n"; cout << "lmdof " << lmdof << " dof " << dof << endl; return NULL; } } WT->FillComplete(*problemmap_,*annmap); WT_ = Teuchos::rcp(WT); //-------------------------------------------------------------------------- // 3) create B // create a temporary matrix with rowmap of the Ann block Epetra_CrsMatrix* tmp = new Epetra_CrsMatrix(Copy,*annmap,120); std::vector<int> newindices(100); for (intintcurr=lm_to_dof.begin(); intintcurr!=lm_to_dof.end(); ++intintcurr) { int lmdof = intintcurr->first; int dof = intintcurr->second; if (D_->MyGRID(lmdof)==false) continue; int lmlrid = D_->LRID(lmdof); if (lmlrid<0) cout << "Cannot find lmlrid for lmdof\n"; int numentries; int* indices; double* values; // extract and add values from D err = D_->ExtractMyRowView(lmlrid,numentries,values,indices); if (err) cout << "D_->ExtractMyRowView returned err=" << err << endl; if (numentries>(int)newindices.size()) newindices.resize(numentries); for (int j=0; j<numentries; ++j) { newindices[j] = D_->GCID(indices[j]); if (newindices[j]<0) cout << "Cannot find gcid for indices[j]\n"; } //cout << "Inserting from D in row " << dof << " cols/val "; //for (int j=0; j<numentries; ++j) cout << newindices[j] << "/" << values[j] << " "; //cout << endl; err = tmp->InsertGlobalValues(dof,numentries,values,&newindices[0]); if (err) cout << "tmp->InsertGlobalValues returned err=" << err << endl; // extract and add values from M err = M_->ExtractMyRowView(lmlrid,numentries,values,indices); if (err) cout << "M_->ExtractMyRowView returned err=" << err << endl; if (numentries>(int)newindices.size()) newindices.resize(numentries); for (int j=0; j<numentries; ++j) { newindices[j] = M_->GCID(indices[j]); if (newindices[j]<0) cout << "Cannot find gcid for indices[j]\n"; } //cout << "Inserting from M in row " << dof << " cols/val "; //for (int j=0; j<numentries; ++j) cout << newindices[j] << "/" << values[j] << " "; //cout << endl; err = tmp->InsertGlobalValues(dof,numentries,values,&newindices[0]); if (err) cout << "tmp->InsertGlobalValues returned err=" << err << endl; } tmp->FillComplete(*(problemmap_.get()),*annmap); // B is transposed of tmp EpetraExt::RowMatrix_Transpose* trans = new EpetraExt::RowMatrix_Transpose(false); Epetra_CrsMatrix* B = &(dynamic_cast<Epetra_CrsMatrix&>(((*trans)(const_cast<Epetra_CrsMatrix&>(*tmp))))); delete tmp; tmp = NULL; B_ = Teuchos::rcp(new Epetra_CrsMatrix(*B)); newindices.clear(); //-------------------------------------------------------------------------- // 4) create I Epetra_CrsMatrix* I = MOERTEL::PaddedMatrix(*problemmap_,1.0,1); I->FillComplete(*problemmap_,*problemmap_); I_ = Teuchos::rcp(I); //-------------------------------------------------------------------------- // 5) Build BWT = B * WT Epetra_CrsMatrix* BWT = MOERTEL::MatMatMult(*B,false,*WT,false,OutLevel()); //-------------------------------------------------------------------------- // 6) create CBT = (I-BWT)*A*W*BT Epetra_CrsMatrix* spdrhs = new Epetra_CrsMatrix(Copy,*problemmap_,10,false); MOERTEL::MatrixMatrixAdd(*BWT,false,-1.0,*spdrhs,0.0); MOERTEL::MatrixMatrixAdd(*I,false,1.0,*spdrhs,1.0); spdrhs->FillComplete(); spdrhs->OptimizeStorage(); spdrhs_ = Teuchos::rcp(spdrhs); Epetra_CrsMatrix* tmp1 = MOERTEL::MatMatMult(*inputmatrix_,false,*WT,true,OutLevel()); Epetra_CrsMatrix* tmp2 = MOERTEL::MatMatMult(*tmp1,false,*B,true,OutLevel()); delete tmp1; Epetra_CrsMatrix* CBT = MOERTEL::MatMatMult(*spdrhs,false,*tmp2,false,OutLevel()); delete tmp2; //-------------------------------------------------------------------------- // 7) create spdmatrix_ = A - CBT - (CBT)^T spdmatrix_ = Teuchos::rcp(new Epetra_CrsMatrix(Copy,*problemmap_,10,false)); MOERTEL::MatrixMatrixAdd(*inputmatrix_,false,1.0,*spdmatrix_,0.0); MOERTEL::MatrixMatrixAdd(*CBT,false,-1.0,*spdmatrix_,1.0); MOERTEL::MatrixMatrixAdd(*CBT,true,-1.0,*spdmatrix_,1.0); delete CBT; CBT = NULL; spdmatrix_->FillComplete(); spdmatrix_->OptimizeStorage(); //-------------------------------------------------------------------------- // tidy up lm_to_dof.clear(); delete trans; B = NULL; // time this process double t = time.ElapsedTime(); if (OutLevel()>5 && Comm().MyPID()==0) cout << "MOERTEL (Proc 0): Construct spd system in " << t << " sec\n"; return spdmatrix_.get(); }
int main(int argc, char *argv[]) { #ifdef HAVE_MPI MPI_Init(&argc, &argv); Epetra_MpiComm Comm(MPI_COMM_WORLD); #else Epetra_SerialComm Comm; #endif // initialize the random number generator int ml_one = 1; ML_srandom1(&ml_one); // ===================== // // create linear problem // // ===================== // Epetra_CrsMatrix *BadMatrix; int * i_blockids; int numblocks; EpetraExt::MatlabFileToCrsMatrix("samplemat.dat",Comm,BadMatrix); const Epetra_Map *Map=&BadMatrix->RowMap(); // Read in the block GLOBAL ids Epetra_Vector* d_blockids; int rv=EpetraExt::MatrixMarketFileToVector("blockids.dat",*Map,d_blockids); i_blockids=new int[d_blockids->MyLength()]; numblocks=-1; for(int i=0;i<d_blockids->MyLength();i++){ i_blockids[i]=(int)(*d_blockids)[i]-1; numblocks=MAX(numblocks,i_blockids[i]); } numblocks++; BadMatrix->FillComplete(); BadMatrix->OptimizeStorage(); int N=BadMatrix->RowMatrixRowMap().NumMyElements(); // Create the trivial blockID list int * trivial_blockids=new int[N]; for(int i=0;i<N;i++) trivial_blockids[i]=i; Epetra_Vector LHS(*Map); Epetra_Vector RHS(*Map); Epetra_LinearProblem BadProblem(BadMatrix, &LHS, &RHS); Teuchos::ParameterList MLList; double TotalErrorResidual = 0.0, TotalErrorExactSol = 0.0; char mystring[80]; // ====================== // // ML Cheby // ====================== // if (Comm.MyPID() == 0) PrintLine(); ML_Epetra::SetDefaults("SA",MLList); MLList.set("smoother: type","Chebyshev"); MLList.set("coarse: type","Amesos-KLU"); MLList.set("max levels",2); MLList.set("aggregation: threshold",.02); MLList.set("ML output",10); MLList.set("smoother: polynomial order",2); strcpy(mystring,"Cheby"); TestMultiLevelPreconditioner(mystring, MLList, BadProblem, TotalErrorResidual, TotalErrorExactSol); // These tests only work in serial due to how the block id's are written. if(Comm.NumProc()==1){ // ====================== // // ML Block Cheby (Trivial) // ====================== // if (Comm.MyPID() == 0) PrintLine(); ML_Epetra::SetDefaults("SA",MLList); MLList.set("smoother: type","Block Chebyshev"); MLList.set("smoother: Block Chebyshev number of blocks",N); MLList.set("smoother: Block Chebyshev block list",trivial_blockids); MLList.set("coarse: type","Amesos-KLU"); MLList.set("max levels",2); MLList.set("ML output",10); MLList.set("smoother: polynomial order",2); strcpy(mystring,"ML Block Cheby (Trivial)"); TestMultiLevelPreconditioner(mystring, MLList, BadProblem, TotalErrorResidual, TotalErrorExactSol); // ====================== // // ML Block Cheby (Smart) // ====================== // if (Comm.MyPID() == 0) PrintLine(); ML_Epetra::SetDefaults("SA",MLList); MLList.set("smoother: type","Block Chebyshev"); MLList.set("smoother: Block Chebyshev number of blocks",numblocks); MLList.set("smoother: Block Chebyshev block list",i_blockids); MLList.set("coarse: type","Amesos-KLU"); MLList.set("max levels",2); MLList.set("ML output",10); MLList.set("smoother: polynomial order",2); strcpy(mystring,"ML Block Cheby (Smart)"); TestMultiLevelPreconditioner(mystring, MLList, BadProblem, TotalErrorResidual, TotalErrorExactSol); } #if defined(HAVE_ML_IFPACK) // ====================== // // IFPACK Cheby // ====================== // if (Comm.MyPID() == 0) PrintLine(); ML_Epetra::SetDefaults("SA",MLList); MLList.set("smoother: type","IFPACK-Chebyshev"); MLList.set("coarse: type","Amesos-KLU"); MLList.set("max levels",2); MLList.set("ML output",10); MLList.set("smoother: polynomial order",2); strcpy(mystring,"IFPACK Cheby"); TestMultiLevelPreconditioner(mystring, MLList, BadProblem, TotalErrorResidual, TotalErrorExactSol); // ====================== // // IFPACK Block Cheby (Trivial) // ====================== // int NumBlocks=Map->NumMyElements(); int *BlockStarts=new int[NumBlocks+1]; int *Blockids=new int [NumBlocks]; for(int i=0;i<NumBlocks;i++){ BlockStarts[i]=i; Blockids[i]=Map->GID(i); } BlockStarts[NumBlocks]=NumBlocks; if (Comm.MyPID() == 0) PrintLine(); ML_Epetra::SetDefaults("SA",MLList); MLList.set("smoother: type","IFPACK-Block Chebyshev"); MLList.set("smoother: Block Chebyshev number of blocks",NumBlocks); MLList.set("smoother: Block Chebyshev block starts",&BlockStarts[0]); MLList.set("smoother: Block Chebyshev block list",&Blockids[0]); MLList.set("coarse: type","Amesos-KLU"); MLList.set("max levels",2); MLList.set("ML output",10); MLList.set("smoother: polynomial order",2); strcpy(mystring,"IFPACK Block Cheby (Trivial)"); TestMultiLevelPreconditioner(mystring, MLList, BadProblem, TotalErrorResidual, TotalErrorExactSol); delete [] BlockStarts; delete [] Blockids; // ====================== // // IFPACK Block Cheby (Smart) // ====================== // // Figure out how many blocks we actually have and build a map... Epetra_Map* IfpackMap; int g_NumBlocks=-1,g_MaxSize=-1; if(Comm.MyPID() == 0){ const int lineLength = 1025; char line[lineLength]; FILE * f=fopen("localids_in_blocks.dat","r"); assert(f!=0); // Next, strip off header lines (which start with "%") do { if(fgets(line, lineLength,f)==0) return(-4); } while (line[0] == '%'); // Grab the number we actually care about sscanf(line, "%d %d", &g_NumBlocks, &g_MaxSize); fclose(f); } Comm.Broadcast(&g_NumBlocks,1,0); Comm.Broadcast(&g_MaxSize,1,0); Epetra_Map BlockMap(g_NumBlocks,0,Comm); Epetra_MultiVector *blockids_disk=0; rv=EpetraExt::MatrixMarketFileToMultiVector("localids_in_blocks.dat",BlockMap,blockids_disk); // Put all the block info into the right place NumBlocks=BlockMap.NumMyElements(); BlockStarts=new int[NumBlocks+1]; Blockids= new int[g_MaxSize*NumBlocks]; // NTS: Blockids_ is overallocated because I don't want to write a counting loop int i,cidx; for(i=0,cidx=0;i<NumBlocks;i++){ BlockStarts[i]=cidx; Blockids[cidx]=(int)(*blockids_disk)[0][i];cidx++; if((*blockids_disk)[1][i] > 1e-2){ Blockids[cidx]=(int)(*blockids_disk)[1][i];cidx++; } } BlockStarts[NumBlocks]=cidx; if (Comm.MyPID() == 0) PrintLine(); ML_Epetra::SetDefaults("SA",MLList); MLList.set("smoother: type","IFPACK-Block Chebyshev"); MLList.set("smoother: Block Chebyshev number of blocks",NumBlocks); MLList.set("smoother: Block Chebyshev block starts",&BlockStarts[0]); MLList.set("smoother: Block Chebyshev block list",&Blockids[0]); MLList.set("coarse: type","Amesos-KLU"); MLList.set("max levels",2); MLList.set("ML output",10); MLList.set("smoother: polynomial order",2); strcpy(mystring,"IFPACK Block Cheby (Smart)"); TestMultiLevelPreconditioner(mystring, MLList, BadProblem, TotalErrorResidual, TotalErrorExactSol); delete blockids_disk; delete [] BlockStarts; delete [] Blockids; #endif // ===================== // // print out total error // // ===================== // if (Comm.MyPID() == 0) { cout << endl; cout << "......Total error for residual = " << TotalErrorResidual << endl; cout << "......Total error for exact solution = " << TotalErrorExactSol << endl; cout << endl; } delete [] i_blockids; delete [] trivial_blockids; delete BadMatrix; delete d_blockids; if (TotalErrorResidual > 1e-5) { cerr << "Error: `BlockCheby.exe' failed!" << endl; exit(EXIT_FAILURE); } #ifdef HAVE_MPI MPI_Finalize(); #endif if (Comm.MyPID() == 0) cerr << "`BlockCheby.exe' passed!" << endl; return (EXIT_SUCCESS); }