// ============================================================================
int ML_Epetra::MultiLevelPreconditioner::
VisualizeSmoothers(int NumPreCycles, int NumPostCycles)
{

  bool viz = List_.get("viz: enable",false);
  if (viz) {
    if (IsPreconditionerComputed() == false)
      ML_CHK_ERR(-1); // need an already computed preconditioner

    bool VizPreSmoother = false;
    bool VizPostSmoother = false;

    if (NumPreCycles != 0)
      VizPreSmoother = true;
    if (NumPostCycles != 0)
      VizPostSmoother = true;

    int ierr = Visualize(false, VizPreSmoother, VizPostSmoother,
                 false, NumPreCycles, NumPostCycles, -1);

    ML_CHK_ERR(ierr);
  }
  else
  {
    std::cout << PrintMsg_ << "You need to specify `viz: enable' = true" << std::endl;
    std::cout << PrintMsg_ << "in the parameter list before building the ML" << std::endl;
    std::cout << PrintMsg_ << "preconditioner in order to visualize" << std::endl;
    ML_CHK_ERR(-1);
  }

  return(0);
}
Exemple #2
0
// ================================================ ====== ==== ==== == = 
int ML_Epetra::FaceMatrixFreePreconditioner::NodeAggregate(ML_Aggregate_Struct *&MLAggr,ML_Operator *&P,ML_Operator* TMT_ML,int &NumAggregates){
  /* Pull Teuchos Options */
  string CoarsenType = List_.get("aggregation: type", "Uncoupled");
  double Threshold   = List_.get("aggregation: threshold", 0.0);  
  int    NodesPerAggr = List_.get("aggregation: nodes per aggregate", 
                                  ML_Aggregate_Get_OptimalNumberOfNodesPerAggregate());

  string PrintMsg_ = "FMFP (Level 0): ";

  ML_Aggregate_Create(&MLAggr);
  ML_Aggregate_Set_MaxLevels(MLAggr, 2);
  ML_Aggregate_Set_StartLevel(MLAggr, 0);
  ML_Aggregate_Set_Threshold(MLAggr, Threshold);
  ML_Aggregate_Set_MaxCoarseSize(MLAggr,1);
  MLAggr->cur_level = 0;
  ML_Aggregate_Set_Reuse(MLAggr); 
  MLAggr->keep_agg_information = 1;  
  P = ML_Operator_Create(ml_comm_);
  
  /* Process Teuchos Options */
  if (CoarsenType == "Uncoupled")
    ML_Aggregate_Set_CoarsenScheme_Uncoupled(MLAggr);
  else if (CoarsenType == "Uncoupled-MIS"){
    ML_Aggregate_Set_CoarsenScheme_UncoupledMIS(MLAggr);
  }
  else if (CoarsenType == "METIS"){
    ML_Aggregate_Set_CoarsenScheme_METIS(MLAggr);
    ML_Aggregate_Set_NodesPerAggr(0, MLAggr, 0, NodesPerAggr);
  }/*end if*/
  else {
    if(!Comm_->MyPID()) printf("FMFP: Unsupported (1,1) block aggregation type(%s), resetting to uncoupled-mis\n",CoarsenType.c_str());
    ML_Aggregate_Set_CoarsenScheme_UncoupledMIS(MLAggr);
  }

  /* Aggregate Nodes */
  int printlevel=ML_Get_PrintLevel();
  ML_Set_PrintLevel(10);
  NumAggregates = ML_Aggregate_Coarsen(MLAggr, TMT_ML, &P, ml_comm_);
  ML_Set_PrintLevel(printlevel);

  if (NumAggregates == 0){
    cerr << "Found 0 aggregates, perhaps the problem is too small." << endl;
    ML_CHK_ERR(-2);
  }/*end if*/
  else if(very_verbose_) printf("[%d] FMFP: %d aggregates created invec_leng=%d\n",Comm_->MyPID(),NumAggregates,P->invec_leng);

  int globalAggs;
  Comm_->SumAll(&NumAggregates,&globalAggs,1);
  if( verbose_ && !Comm_->MyPID()) {
    std::cout << PrintMsg_ << "Aggregation threshold = " << Threshold << std::endl;
    std::cout << PrintMsg_ << "Global aggregates = " << globalAggs << std::endl;
    //ML_Aggregate_Print_Complexity(MLAggr);
  }


  if(P==0) {fprintf(stderr,"%s","ERROR: No tentative prolongator found\n");ML_CHK_ERR(-5);}
  return 0;
}  
// ================================================ ====== ==== ==== == =
// Forms the coarse matrix, given the prolongator
int  ML_Epetra::FaceMatrixFreePreconditioner::FormCoarseMatrix()
{
  CoarseMat_ML = ML_Operator_Create(ml_comm_);
  CoarseMat_ML->data_destroy=free;
  ML_Operator *Temp_ML=0;
  ML_Operator *R= ML_Operator_Create(ml_comm_);
  ML_Operator *P= ML_Operator_Create(ml_comm_);

  /* Build ML_Operator version of Prolongator_, Restriction Operator */
  ML_CHK_ERR(ML_Operator_WrapEpetraCrsMatrix(Prolongator_,P,verbose_));
  P->num_rigid=P->num_PDEs=dim;

  //NTS: ML_CHK_ERR won't work on this: it returns 1
  ML_Operator_Transpose_byrow(P, R);

  /* OPTION: Disable the addon */
  const Epetra_CrsMatrix *Op11crs = dynamic_cast<const Epetra_CrsMatrix*>(&*Operator_);
  const Epetra_Operator_With_MatMat *Op11mm = dynamic_cast<const Epetra_Operator_With_MatMat*>(&*Operator_);

  /* Do the A*P  with or without addon*/
  if(Op11crs){
    if(verbose_ && !Comm_->MyPID()) printf("FMFP: Running *without* addon\n");
    ML_Operator *SM_ML = ML_Operator_Create(ml_comm_);
    Temp_ML = ML_Operator_Create(ml_comm_);
    ML_Operator_WrapEpetraCrsMatrix((Epetra_CrsMatrix*)Op11crs,SM_ML,verbose_);
    ML_2matmult(SM_ML,P,Temp_ML,ML_CSR_MATRIX);
    ML_Operator_Destroy(&SM_ML);
  }
  else if(Op11mm){
    if(verbose_ && !Comm_->MyPID()) printf("FMFP: Running with addon\n");
    ML_CHK_ERR(Op11mm->MatrixMatrix_Multiply(*Prolongator_,ml_comm_,&Temp_ML));
  }
  else{
    if(!Comm_->MyPID()) printf("ERROR: FMFP Illegal Operator\n");
    delete R;
    ML_CHK_ERR(-1);
  }

  /* Do R * AP */
  R->num_rigid=R->num_PDEs=dim;
  ML_2matmult_block(R, Temp_ML,CoarseMat_ML,ML_CSR_MATRIX);

  /* Wrap to Epetra-land */
  int nnz=100;
  double time;
  ML_Operator2EpetraCrsMatrix(CoarseMat_ML,CoarseMatrix,nnz,true,time,0,verbose_);
  // NTS: This is a hack to get around the sticking ones on the diagonal issue;

  /* Cleanup */
  ML_Operator_Destroy(&P);
  ML_Operator_Destroy(&R);
  ML_Operator_Destroy(&Temp_ML);
  ML_Operator_Destroy(&CoarseMat_ML);CoarseMat_ML=0;//HAX
  return 0;
}/*end FormCoarseMatrix*/
// ============================================================================
int ML_Epetra::MultiLevelPreconditioner::
VisualizeAggregates()
{
  bool viz = List_.get("viz: enable",false);
  if (viz) {
    if (IsPreconditionerComputed() == false)
      ML_CHK_ERR(-1); // need an already computed preconditioner

    ML_CHK_ERR(Visualize(true, false, false, false, -1, -1, -1));

  }
  return(0);
}
int ML_Amesos_Solve( void *data, double x[], double rhs[] )
{
  Amesos_Handle_Type *Amesos_Handle = (Amesos_Handle_Type *) data;
  if (Amesos_Handle->A_Base == 0) return 0;

  Amesos_BaseSolver *A_Base = (Amesos_BaseSolver *) Amesos_Handle->A_Base ;

  Epetra_Time Time(A_Base->Comm());  

  Epetra_LinearProblem *Amesos_LinearProblem = (Epetra_LinearProblem *)A_Base->GetProblem() ;
  
  const Epetra_BlockMap & map = Amesos_LinearProblem->GetOperator()->OperatorDomainMap() ; 

  Epetra_Vector EV_rhs( View, map, rhs ) ;
  Epetra_Vector EV_lhs( View, map, x ) ;

  Amesos_LinearProblem->SetRHS( &EV_rhs ) ; 
  Amesos_LinearProblem->SetLHS( &EV_lhs ) ;

  A_Base->Solve() ; 

  TimeForSolve__ += Time.ElapsedTime();
  NumSolves__++;

#ifdef ML_AMESOS_DEBUG
  // verify that the residual is actually small (and print the max
  // in the destruction phase)
  Epetra_Vector Ax(map);
  
  (Amesos_LinearProblem->GetMatrix())->Multiply(false,EV_lhs,Ax);
  
  ML_CHK_ERR(Ax.Update(1.0, EV_rhs, -1.0));
  
  double residual;
  ML_CHK_ERR(Ax.Norm2(&residual));
  if( residual > MaxError__ ) MaxError__ = residual;
#endif
  
  return 0;
} //ML_Amesos_Solve()
// ============================================================================
int ML_Epetra::MultiLevelPreconditioner::
VisualizeCycle(int NumCycles)
{

  bool viz = List_.get("viz: enable",false);
  if (viz) {
    if (IsPreconditionerComputed() == false)
      ML_CHK_ERR(-1); // need an already computed preconditioner

    int ierr = Visualize(false, false, false, true,
               -1, -1, NumCycles);

    ML_CHK_ERR(ierr);
  }
  else
  {
    cout << PrintMsg_ << "You need to specify `viz: enable' = true" << endl;
    cout << PrintMsg_ << "in the parameter list before building the ML" << endl;
    cout << PrintMsg_ << "preconditioner in order to visualize" << endl;
    ML_CHK_ERR(-1);
  }
  return(0);
}
// ================================================ ====== ==== ==== == =
// Computes C= <me> * A
int  ML_Epetra::Epetra_Multi_CrsMatrix::MatrixMatrix_Multiply(const Epetra_CrsMatrix & A, ML_Comm *comm, ML_Operator **C) const
{

  int rv=0;
  ML_Comm* temp = global_comm;

  /* Setup for 1st Matmat */
  ML_Operator * MV[2]={0,0},*CV;
  MV[(NumMatrices_-1)%2]= ML_Operator_Create(comm);
  rv=ML_Operator_WrapEpetraCrsMatrix((Epetra_CrsMatrix*)&A,MV[(NumMatrices_-1)%2]);
  ML_CHK_ERR(rv);

  /* Do the matmats */
  for(int i=NumMatrices_-1;rv==0 && i>=0;i--){
    /* Do pre-wraps */
    if(MV[(i+1)%2] && i!=NumMatrices_-1) ML_Operator_Destroy(&MV[(i+1)%2]);
    MV[(i+1)%2]=ML_Operator_Create(comm);
    CV=ML_Operator_Create(comm);
    rv=ML_Operator_WrapEpetraCrsMatrix(CrsMatrices_[i],CV);
    ML_CHK_ERR(rv);

    /* Do matmat */
    ML_2matmult(CV,MV[i%2],MV[(i+1)%2],ML_CSR_MATRIX);

    ML_Operator_Destroy(&CV);
  }/*end for*/
  global_comm = temp;

  /* Final Postwrap */
  *C=MV[1];

  /* Cleanup */
  if(MV[0]) ML_Operator_Destroy(&MV[0]);

  return rv;
}/*end MatrixMatrix_Multiply*/
// ================================================ ====== ==== ==== == =
//! Build the face-to-node prolongator described by Bochev, Siefert, Tuminaro, Xu and Zhu (2007).
int ML_Epetra::FaceMatrixFreePreconditioner::PBuildSparsity(ML_Operator *P, Epetra_CrsMatrix *&Psparse){

  /* Create wrapper to do abs(T) */
  // NTS: Assume D0 has already been reindexed by now.
  ML_Operator* AbsFN_ML = ML_Operator_Create(ml_comm_);
  ML_CHK_ERR(ML_Operator_WrapEpetraCrsMatrix(const_cast<Epetra_CrsMatrix*>(&*FaceNode_Matrix_),AbsFN_ML,verbose_));
  ML_Operator_Set_Getrow(AbsFN_ML,AbsFN_ML->outvec_leng,CSR_getrow_ones);

  /* Form abs(T) * P_n */
  ML_Operator* AbsFNP = ML_Operator_Create(ml_comm_);
  ML_2matmult(AbsFN_ML,P,AbsFNP, ML_CSR_MATRIX);

  /* Wrap P_n into Epetra-land */
  Epetra_CrsMatrix_Wrap_ML_Operator(AbsFNP,*Comm_,*FaceRangeMap_,&Psparse,Copy,0);

  /* Nuke the rows in Psparse */
  if(BCfaces_.size()>0) Apply_BCsToMatrixRows(BCfaces_.get(),BCfaces_.size(),*Psparse);

  // Cleanup
  ML_Operator_Destroy(&AbsFN_ML);
  ML_Operator_Destroy(&AbsFNP);

  return 0;
}
// ================================================ ====== ==== ==== == =
//! Apply the preconditioner to an Epetra_MultiVector X, puts the result in Y
int ML_Epetra::FaceMatrixFreePreconditioner::ApplyInverse(const Epetra_MultiVector& B_, Epetra_MultiVector& X) const{
  const Epetra_MultiVector *B;
  Epetra_MultiVector *Bcopy=0;

  /* Sanity Checks */
  int NumVectors=B_.NumVectors();
  if (!B_.Map().SameAs(*FaceDomainMap_)) ML_CHK_ERR(-1);
  if (NumVectors != X.NumVectors()) ML_CHK_ERR(-1);

  Epetra_MultiVector r_edge(*FaceDomainMap_,NumVectors,false);
  Epetra_MultiVector e_edge(*FaceDomainMap_,NumVectors,false);
  Epetra_MultiVector e_node(*CoarseMap_,NumVectors,false);
  Epetra_MultiVector r_node(*CoarseMap_,NumVectors,false);

  /* Deal with the B==X case */
  if (B_.Pointers()[0] == X.Pointers()[0]){
    Bcopy=new Epetra_MultiVector(B_);
    B=Bcopy;
    X.PutScalar(0.0);
  }
  else B=&B_;


  for(int i=0;i<num_cycles;i++){
    /* Pre-smoothing */
#ifdef HAVE_ML_IFPACK
    if(Smoother_) ML_CHK_ERR(Smoother_->ApplyInverse(*B,X));
#endif

    if(MaxLevels > 0){
      if(i != 0
#ifdef HAVE_ML_IFPACK
         || Smoother_
#endif
         ){
        /* Calculate Residual (r_e = b - (S+M+Addon) * x) */
        ML_CHK_ERR(Operator_->Apply(X,r_edge));
        ML_CHK_ERR(r_edge.Update(1.0,*B,-1.0));

        /* Xfer to coarse grid (r_n = P' * r_e) */
        ML_CHK_ERR(Prolongator_->Multiply(true,r_edge,r_node));
      }
      else{
        /* Xfer to coarse grid (r_n = P' * r_e) */
        ML_CHK_ERR(Prolongator_->Multiply(true,*B,r_node));
      }

      /* AMG on coarse grid  (e_n = (CoarseMatrix)^{-1} r_n) */
      ML_CHK_ERR(CoarsePC->ApplyInverse(r_node,e_node));

      /* Xfer back to fine grid (e_e = P * e_n) */
      ML_CHK_ERR(Prolongator_->Multiply(false,e_node,e_edge));

      /* Add in correction (x = x + e_e)        */
      ML_CHK_ERR(X.Update(1.0,e_edge,1.0));
    }/*end if*/

    /* Post-Smoothing */
#ifdef HAVE_ML_IFPACK
    if(Smoother_) ML_CHK_ERR(Smoother_->ApplyInverse(*B,X));
#endif

  }/*end for*/

  /* Cleanup */
  if(Bcopy) delete Bcopy;

  return 0;
}/*end ApplyInverse*/
// ================================================ ====== ==== ==== == =
//! Build the face-to-node prolongator described by Bochev, Siefert, Tuminaro, Xu and Zhu (2007).
int ML_Epetra::FaceMatrixFreePreconditioner::BuildProlongator()
{

  /* Wrap TMT_Matrix in a ML_Operator */
  ML_Operator* TMT_ML = ML_Operator_Create(ml_comm_);
  ML_Operator_WrapEpetraCrsMatrix(const_cast<Epetra_CrsMatrix*>(&*TMT_Matrix_),TMT_ML);

  /* Nodal Aggregation */
  ML_Aggregate_Struct *MLAggr=0;
  ML_Operator *P=0;
  int NumAggregates;
  int rv=ML_Epetra::RefMaxwell_Aggregate_Nodes(*TMT_Matrix_,List_,ml_comm_,std::string("FMFP (level 0) :"),MLAggr,P,NumAggregates);
  if(rv || !P) {if(!Comm_->MyPID()) printf("ERROR: Building nodal P\n");ML_CHK_ERR(-1);}

  /* Build 1-unknown sparsity of prolongator */
  Epetra_CrsMatrix *Psparse=0;
  PBuildSparsity(P,Psparse);
  if(!Psparse) {if(!Comm_->MyPID()) printf("ERROR: Building Psparse\n");ML_CHK_ERR(-2);}

  /* Build the "nullspace" */
  Epetra_MultiVector *nullspace;
  BuildNullspace(nullspace);
  if(!nullspace) {if(!Comm_->MyPID()) printf("ERROR: Building Nullspace\n");ML_CHK_ERR(-3);}


  /* Build the DomainMap of the new operator*/
  const Epetra_Map & FineColMap = Psparse->ColMap();
  CoarseMap_=new Epetra_Map(-1,NumAggregates*dim,0,*Comm_);

  /* Allocate the Prolongator_ */
  int max_nz_per_row=Psparse->MaxNumEntries();
  Prolongator_=new Epetra_CrsMatrix(Copy,*FaceRangeMap_,0);
  int ne1, *idx1, *idx2;
  idx2=new int [dim*max_nz_per_row];
  double *vals1, *vals2;
  vals2=new double[dim*max_nz_per_row];
  int nonzeros;

  for(int i=0;i<Prolongator_->NumMyRows();i++){
    Psparse->ExtractMyRowView(i,ne1,vals1,idx1);
    nonzeros=0;
    for(int j=0;j<ne1;j++) nonzeros+=ABS(vals1[j])>0;

    for(int j=0;j<ne1;j++){
      for(int k=0;k<dim;k++) {
        idx2[j*dim+k]=FineColMap.GID(idx1[j])*dim+k;
        //FIX: This works only because there's an implicit linear mapping which
        //we're exploiting.
        if(idx2[j*dim+k]==-1) printf("[%d] ERROR: idx1[j]=%d / idx1[j]*dim+k=%d does not have a GID!\n",Comm_->MyPID(),idx1[j],idx1[j]*dim+k);
        if(vals1[j]==0 ) vals2[j*dim+k]=0;
	else vals2[j*dim+k]=(*nullspace)[k][i] / nonzeros;
      }/*end for*/
    }/*end for*/
    Prolongator_->InsertGlobalValues(FaceRangeMap_->GID(i),dim*ne1,vals2,idx2);
  }/*end for*/


  /* FillComplete / OptimizeStorage for Prolongator*/
  Prolongator_->FillComplete(*CoarseMap_,*FaceRangeMap_);
  Prolongator_->OptimizeStorage();

#ifndef NO_OUTPUT
  /* DEBUG: Dump aggregates */
  Epetra_IntVector AGG(View,*NodeDomainMap_,MLAggr->aggr_info[0]);
  IVOUT(AGG,"agg.dat");
  EpetraExt::RowMatrixToMatlabFile("fmfp_psparse.dat",*Psparse);
  EpetraExt::RowMatrixToMatlabFile("fmfp_prolongator.dat",*Prolongator_);
  EpetraExt::VectorToMatrixMarketFile("fmfp_null0.dat",*(*nullspace)(0));
  EpetraExt::VectorToMatrixMarketFile("fmfp_null1.dat",*(*nullspace)(1));
  EpetraExt::VectorToMatrixMarketFile("fmfp_null2.dat",*(*nullspace)(2));

#endif

  /* EXPERIMENTAL: Normalize Prolongator Columns */
  bool normalize_prolongator=List_.get("face matrix free: normalize prolongator",false);
  if(normalize_prolongator){
    Epetra_Vector n_vector(*CoarseMap_,false);
    Prolongator_->InvColSums(n_vector);
    Prolongator_->RightScale(n_vector);
  }/*end if*/

  /* Post-wrapping to convert to ML indexing */
#ifdef HAVE_ML_EPETRAEXT
  Prolongator_ = dynamic_cast<Epetra_CrsMatrix*>(ModifyEpetraMatrixColMap(*Prolongator_,ProlongatorColMapTrans_,"Prolongator",(verbose_&&!Comm_->MyPID())));
#endif

  /* Cleanup */
  ML_qr_fix_Destroy();
  ML_Aggregate_Destroy(&MLAggr);
  ML_Operator_Destroy(&TMT_ML);
  ML_Operator_Destroy(&P);

  delete nullspace;
  delete Psparse;
  delete [] idx2;
  delete [] vals2;
  return 0;
}/*end BuildProlongator_*/
// ================================================ ====== ==== ==== == =
//! Build pi operator described by Bochev, Siefert, Tuminaro, Xu and Zhu (2007).
int ML_Epetra::FaceMatrixFreePreconditioner::BuildNullspace(Epetra_MultiVector *& nullspace){
  int Nf=FaceRangeMap_->NumMyElements();

  /* Pull the coordinates from Teuchos */
  double * xcoord=List_.get("x-coordinates",(double*)0);
  double * ycoord=List_.get("y-coordinates",(double*)0);
  double * zcoord=List_.get("z-coordinates",(double*)0);
  dim=(xcoord!=0) + (ycoord!=0) + (zcoord!=0);

  // Sanity check
  if(dim!=3){if(!Comm_->MyPID()) printf("ERROR: FaceMatrixFreePreconditioner only works in 3D"); ML_CHK_ERR(-1);}

  // Build the (unimported) coordinate multivector
  double **d_coords=new double* [dim];
  d_coords[0]=xcoord; d_coords[1]=ycoord;
  if(dim==3) d_coords[2]=zcoord;
  Epetra_MultiVector n_coords_domain(View,*NodeDomainMap_,d_coords,dim);
  Epetra_MultiVector *n_coords;

  // Import coordinate info
  if(FaceNode_Matrix_->Importer()){
    n_coords=new Epetra_MultiVector(FaceNode_Matrix_->ColMap(),dim);
    n_coords->PutScalar(0.0);
    n_coords->Import(n_coords_domain,*FaceNode_Matrix_->Importer(),Add);
  }
  else n_coords=&n_coords_domain;

  // Sanity HAQ - Only works on Hexes
  if(FaceNode_Matrix_->GlobalMaxNumEntries()!=4)
    {if(!Comm_->MyPID()) printf("ERROR: FaceMatrixFreePreconditioner only works on Hexes"); ML_CHK_ERR(-2);}

  // Allocate vector
  nullspace=new Epetra_MultiVector(*FaceDomainMap_,3);

  // Fill matrix - NTS this will *NOT* do periodic BCs correctly.
  double *a=new double[dim];
  double *b=new double[dim];
  double *c=new double[dim];
  for(int i=0;i<Nf;i++){
    int Ni,*indices;
    double *values;
    FaceNode_Matrix_->ExtractMyRowView(i,Ni,values,indices);
    if(Ni != 4){
      printf("ERROR: Face %d has only %d nodes\n",i,Ni);
      ML_free(a);
      ML_free(b);
      ML_free(c);
      ML_CHK_ERR(-1);
    }

    a[0] = (*n_coords)[0][indices[1]] - (*n_coords)[0][indices[0]];
    a[1] = (*n_coords)[1][indices[1]] - (*n_coords)[1][indices[0]];
    a[2] = (*n_coords)[2][indices[1]] - (*n_coords)[2][indices[0]];

    b[0] = (*n_coords)[0][indices[2]] - (*n_coords)[0][indices[0]];
    b[1] = (*n_coords)[1][indices[2]] - (*n_coords)[1][indices[0]];
    b[2] = (*n_coords)[2][indices[2]] - (*n_coords)[2][indices[0]];

    cross_product(a,b,c);

    // HAQ - Hardwiring for hexes
    // HAQ - Absolute value, presuming all hexes are actually pointed the same way.  This is a HAQ!!!
    (*nullspace)[0][i]=ABS(c[0])/6.0;
    (*nullspace)[1][i]=ABS(c[1])/6.0;
    (*nullspace)[2][i]=ABS(c[2])/6.0;
  }

  /* Cleanup */
  if(FaceNode_Matrix_->Importer()) delete n_coords;
  delete [] a; delete [] b; delete [] c;
  delete [] d_coords;
  return 0;
}
// ================================================ ====== ==== ==== == =
int ML_Epetra::RefMaxwell_Aggregate_Nodes(const Epetra_CrsMatrix & A, Teuchos::ParameterList & List, ML_Comm * ml_comm, std::string PrintMsg,
					  ML_Aggregate_Struct *& MLAggr,ML_Operator *&P, int &NumAggregates){

  /* Output level */
  bool verbose, very_verbose;
  int OutputLevel = List.get("ML output", -47);
  if(OutputLevel == -47) OutputLevel = List.get("output", 1);
  if(OutputLevel>=15) very_verbose=verbose=true;
  if(OutputLevel > 5) {very_verbose=false;verbose=true;}
  else very_verbose=verbose=false;

  /* Wrap A in a ML_Operator */
  ML_Operator* A_ML = ML_Operator_Create(ml_comm);
  ML_Operator_WrapEpetraCrsMatrix(const_cast<Epetra_CrsMatrix*>(&A),A_ML);

 /* Pull Teuchos Options */
  std::string CoarsenType  = List.get("aggregation: type", "Uncoupled");
  double Threshold    = List.get("aggregation: threshold", 0.0);
  int    NodesPerAggr = List.get("aggregation: nodes per aggregate",
                                  ML_Aggregate_Get_OptimalNumberOfNodesPerAggregate());
  bool UseAux         = List.get("aggregation: aux: enable",false);
  double AuxThreshold = List.get("aggregation: aux: threshold",0.0);
  int  MaxAuxLevels   = List.get("aggregation: aux: max levels",10);


  ML_Aggregate_Create(&MLAggr);
  ML_Aggregate_Set_MaxLevels(MLAggr, 2);
  ML_Aggregate_Set_StartLevel(MLAggr, 0);
  ML_Aggregate_Set_Threshold(MLAggr, Threshold);
  ML_Aggregate_Set_MaxCoarseSize(MLAggr,1);
  MLAggr->cur_level = 0;
  ML_Aggregate_Set_Reuse(MLAggr);
  MLAggr->keep_agg_information = 1;
  P = ML_Operator_Create(ml_comm);

  /* Process Teuchos Options */
  if (CoarsenType == "Uncoupled")
    ML_Aggregate_Set_CoarsenScheme_Uncoupled(MLAggr);
  else if (CoarsenType == "Uncoupled-MIS"){
    ML_Aggregate_Set_CoarsenScheme_UncoupledMIS(MLAggr);
  }
  else if (CoarsenType == "METIS"){
    ML_Aggregate_Set_CoarsenScheme_METIS(MLAggr);
    ML_Aggregate_Set_NodesPerAggr(0, MLAggr, 0, NodesPerAggr);
  }/*end if*/
  else {
    if(!A.Comm().MyPID()) printf("%s Unsupported (1,1) block aggregation type(%s), resetting to uncoupled-mis\n",PrintMsg.c_str(),CoarsenType.c_str());
    ML_Aggregate_Set_CoarsenScheme_UncoupledMIS(MLAggr);
  }

  /* Setup Aux Data */
  if(UseAux) {
    A_ML->aux_data->enable=1;
    A_ML->aux_data->threshold=AuxThreshold;
    A_ML->aux_data->max_level=MaxAuxLevels;
    ML_Init_Aux(A_ML,List);
    if(verbose && !A.Comm().MyPID()) {
      printf("%s Using auxiliary matrix\n",PrintMsg.c_str());
      printf("%s aux threshold = %e\n",PrintMsg.c_str(),A_ML->aux_data->threshold);
    }
  }

  /* Aggregate Nodes */
  int printlevel=ML_Get_PrintLevel();
  if(verbose) ML_Set_PrintLevel(10);
  NumAggregates = ML_Aggregate_Coarsen(MLAggr,A_ML, &P, ml_comm);
  if(verbose) ML_Set_PrintLevel(printlevel);

  if (NumAggregates == 0){
    std::cerr << "Found 0 aggregates, perhaps the problem is too small." << std::endl;
    ML_CHK_ERR(-2);
  }/*end if*/
  else if(very_verbose) printf("[%d] %s %d aggregates created invec_leng=%d\n",A.Comm().MyPID(),PrintMsg.c_str(),NumAggregates,P->invec_leng);

  if(verbose){
    int globalAggs=0;
    A.Comm().SumAll(&NumAggregates,&globalAggs,1);
    if(!A.Comm().MyPID()) {
      printf("%s Aggregation threshold = %e\n",PrintMsg.c_str(),Threshold);
      printf("%s Global aggregates     = %d\n",PrintMsg.c_str(),globalAggs);

    }
  }

  /* Cleanup */
  ML_qr_fix_Destroy();
  if(UseAux) ML_Finalize_Aux(A_ML);
  ML_Operator_Destroy(&A_ML);

  return 0;
}
// ============================================================================
int ML_Epetra::MultiLevelPreconditioner::
CreateAuxiliaryMatrixVbr(Epetra_VbrMatrix* &FakeMatrix)
{

  // FakeMatrix has already been created before
  if (FakeMatrix == 0)
    ML_CHK_ERR(-1); // something went wrong

  int NumDimensions = 0;

  double* x_coord = List_.get("x-coordinates", (double *)0);
  if (x_coord != 0) ++NumDimensions;

  // at least x-coordinates must be not null
  if( NumDimensions == 0 ) {
    std::cerr << ErrorMsg_ << "Option `aggregation: use auxiliary matrix' == true" << std::endl
         << ErrorMsg_ << "requires x-, y-, or z-coordinates." << std::endl
         << ErrorMsg_ << "You must specify them using options" << std::endl
         << ErrorMsg_ << "`x-coordinates' (and equivalently for" << std::endl
         << ErrorMsg_ << "y- and z-)." << std::endl;
    ML_CHK_ERR(-2); // wrong parameters
  }

  double* y_coord = List_.get("y-coordinates", (double *)0);
  if (y_coord != 0) ++NumDimensions;

  double* z_coord = List_.get("z-coordinates", (double *)0);
  if (z_coord != 0) ++NumDimensions;

  // small check to avoid strange behavior
  if( z_coord != 0 && y_coord == 0 ) {
    std::cerr << ErrorMsg_ << "Something wrong: `y-coordinates'" << std::endl
         << ErrorMsg_ << "is null, while `z-coordinates' is not null" << std::endl;
    ML_CHK_ERR(-3); // something went wrong
  }

  // usual crap to clutter the output
  if( verbose_ ) {
    std::cout << std::endl;
    std::cout << PrintMsg_ << "*** Using auxiliary matrix to create the aggregates" << std::endl
         << PrintMsg_ << "*** Number of dimensions = " << NumDimensions << std::endl
         << PrintMsg_ << "*** (the version for Epetra_VbrMatrix is currently used)" << std::endl;
    std::cout << std::endl;
  }

  const Epetra_BlockMap& RowMap = FakeMatrix->RowMap();
  const Epetra_BlockMap& ColMap = FakeMatrix->ColMap();
  int NumMyRowElements = RowMap.NumMyElements();
  int NumMyColElements = ColMap.NumMyElements();
  int* MyGlobalRowElements = RowMap.MyGlobalElements();
  int* MyGlobalColElements = ColMap.MyGlobalElements();

  // use point map to exchange coordinates
  Epetra_Map PointRowMap(-1,NumMyRowElements,MyGlobalRowElements,0,Comm());
  Epetra_Map PointColMap(-1,NumMyColElements,MyGlobalColElements,0,Comm());

  Epetra_Vector RowX(PointRowMap);
  Epetra_Vector RowY(PointRowMap);
  Epetra_Vector RowZ(PointRowMap);

  for (int i = 0 ; i < NumMyRowElements ; ++i) {
    RowX[i] = x_coord[i];
    if (NumDimensions > 1) RowY[i] = y_coord[i];
    if (NumDimensions > 2) RowZ[i] = z_coord[i];
  }

  // create vectors containing coordinates for columns
  // (this is useful only if MIS/ParMETIS are used)
  Epetra_Vector ColX(PointColMap);
  Epetra_Vector ColY(PointColMap);
  Epetra_Vector ColZ(PointColMap);

  // get coordinates for non-local nodes (in column map)
  Epetra_Import Importer(PointColMap,PointRowMap);
  ColX.Import(RowX,Importer,Insert);
  if (NumDimensions > 1) ColY.Import(RowY,Importer,Insert);
  if (NumDimensions > 2) ColZ.Import(RowZ,Importer,Insert);

  // room for getrow()
  int MaxNnz = FakeMatrix->MaxNumEntries();
  std::vector<int> colInd(MaxNnz);
  std::vector<double> colVal(MaxNnz);
  std::vector<double> coord_i(3);
  std::vector<double> coord_j(3);

  // change the entries of FakeMatrix so that it corresponds to a discrete
  // Laplacian. Note: This is not exactly the same as in the Crs case.

  FakeMatrix->PutScalar(0.0);

  for (int LocalRow = 0; LocalRow < NumMyRowElements ; ++LocalRow) {

    int RowDim, NumBlockEntries;
    int* BlockIndices;
    Epetra_SerialDenseMatrix** RowValues;

    switch (NumDimensions) {
    case 3:
      coord_i[2] = RowZ[LocalRow];
    case 2:
      coord_i[1] = RowY[LocalRow];
    case 1:
      coord_i[0] = RowX[LocalRow];
    }

    FakeMatrix->ExtractMyBlockRowView(LocalRow,RowDim,NumBlockEntries,
                                      BlockIndices,RowValues);

    // accumulator for zero row-sum
    double total = 0.0;

    for (int j = 0 ; j < NumBlockEntries ; ++j) {

      int LocalCol = BlockIndices[j];

      // insert diagonal later
      if (LocalCol != LocalRow) {

        // get coordinates of this node
        switch (NumDimensions) {
        case 3:
          coord_j[2] = ColZ[LocalCol];
        case 2:
          coord_j[1] = ColY[LocalCol];
        case 1:
          coord_j[0] = ColX[LocalCol];
        }

        // d2 is the square of the distance between node `i' and
        // node `j'
        double d2 = (coord_i[0] - coord_j[0]) * (coord_i[0] - coord_j[0]) +
          (coord_i[1] - coord_j[1]) * (coord_i[1] - coord_j[1]) +
          (coord_i[2] - coord_j[2]) * (coord_i[2] - coord_j[2]);

        if (d2 == 0.0) {
          std::cerr << std::endl;
          std::cerr << ErrorMsg_ << "distance between node " << LocalRow << " and node "
            << LocalCol << std::endl
            << ErrorMsg_ << "is zero. Coordinates of these nodes are" << std::endl
            << ErrorMsg_ << "x_i = " << coord_i[0] << ", x_j = " << coord_j[0] << std::endl
            << ErrorMsg_ << "y_i = " << coord_i[1] << ", y_j = " << coord_j[1] << std::endl
            << ErrorMsg_ << "z_i = " << coord_i[2] << ", z_j = " << coord_j[2] << std::endl
            << ErrorMsg_ << "Now proceeding with distance = 1.0" << std::endl;
          std::cerr << std::endl;
          d2 = 1.0;
        }

        for (int k = 0 ; k < RowValues[j]->M() ; ++k) {
          for (int h = 0 ; h < RowValues[j]->N() ; ++h) {
            if (k == h)
              (*RowValues[j])(k,h) = - 1.0 / d2;
          }
        }

        total += 1.0 /d2;

      }
    }

    // check that the diagonal block exists
    bool ok = false;
    int DiagonalBlock = 0;
    for (int j = 0 ; j < NumBlockEntries ; ++j) {
      if (BlockIndices[j] == LocalRow) {
        DiagonalBlock = j;
        ok = true;
        break;
      }
    }
    assert (ok == true);

    for (int k = 0 ; k < RowValues[DiagonalBlock]->N() ; ++k) {
      (*RowValues[DiagonalBlock])(k,k) = total;
    }
  }

  // tell ML to keep the tentative prolongator
  ML_Aggregate_Set_Reuse(agg_);

  // pray that no bugs will tease us
  return(0);
}
// ================================================ ====== ==== ==== == =
// the tentative null space is in input because the user
// has to remember to allocate and fill it, and then to delete
// it after calling this method.
int ML_Epetra::MultiLevelPreconditioner::
ComputeAdaptivePreconditioner(int TentativeNullSpaceSize,
                              double* TentativeNullSpace)
{

  if ((TentativeNullSpaceSize == 0) || (TentativeNullSpace == 0))
    ML_CHK_ERR(-1);
   
  // ================================== //
  // get parameters from the input list //
  // ================================== //
  
  // maximum number of relaxation sweeps
  int MaxSweeps = List_.get("adaptive: max sweeps", 10);
  // number of std::vector to be added to the tentative null space
  int NumAdaptiveVectors = List_.get("adaptive: num vectors", 1);

  if (verbose_) {
    std::cout << PrintMsg_ << "*** Adaptive Smoother Aggregation setup ***" << std::endl;
    std::cout << PrintMsg_ << "    Maximum relaxation sweeps     = " << MaxSweeps << std::endl;
    std::cout << PrintMsg_ << "    Additional vectors to compute = " << NumAdaptiveVectors << std::endl;
  }

  // ==================================================== //
  // compute the preconditioner, set null space from user //
  // (who will have to delete std::vector TentativeNullSpace)  //
  // ==================================================== //
  
  double* NewNullSpace = 0;
  double* OldNullSpace = TentativeNullSpace;
  int OldNullSpaceSize = TentativeNullSpaceSize;

  // need some work otherwise matvec() with Epetra_Vbr fails.
  // Also, don't differentiate between range and domain here
  // as ML will not work if range != domain
  const Epetra_VbrMatrix* VbrA = NULL;
  VbrA = dynamic_cast<const Epetra_VbrMatrix*>(RowMatrix_);

  Epetra_Vector* LHS = 0;
  Epetra_Vector* RHS = 0;

  if (VbrA != 0) {
    LHS = new Epetra_Vector(VbrA->DomainMap());
    RHS = new Epetra_Vector(VbrA->DomainMap());
  } else {
    LHS = new Epetra_Vector(RowMatrix_->OperatorDomainMap());
    RHS = new Epetra_Vector(RowMatrix_->OperatorDomainMap());
  }

  // destroy what we may already have
  if (IsComputePreconditionerOK_ == true) {
    DestroyPreconditioner();
  }

  // build the preconditioner for the first time
  List_.set("null space: type", "pre-computed");
  List_.set("null space: dimension", OldNullSpaceSize);
  List_.set("null space: vectors", OldNullSpace);
  ComputePreconditioner();

  // ====================== //
  // add one std::vector at time //
  // ====================== //
  
  for (int istep = 0 ; istep < NumAdaptiveVectors ; ++istep) {

    if (verbose_) {
      std::cout << PrintMsg_ << "\tAdaptation step " << istep << std::endl;
      std::cout << PrintMsg_ << "\t---------------" << std::endl;
    }

    // ==================== //
    // look for "bad" modes //
    // ==================== //

    // note: should an error occur, ML_CHK_ERR will return,
    // and LHS and RHS will *not* be delete'd (--> memory leak).
    // Anyway, this means that something wrong happened in the code
    // and should be fixed by the user.

    LHS->Random();
    double Norm2;

    for (int i = 0 ; i < MaxSweeps ; ++i) {
      // RHS = (I - ML^{-1} A) LHS
      ML_CHK_ERR(RowMatrix_->Multiply(false,*LHS,*RHS));
      // FIXME: can do something slightly better here
      ML_CHK_ERR(ApplyInverse(*RHS,*RHS));
      ML_CHK_ERR(LHS->Update(-1.0,*RHS,1.0));
      LHS->Norm2(&Norm2);
      if (verbose_) {
        std::cout << PrintMsg_ << "\titer " << i << ", ||x||_2 = ";
        std::cout << Norm2 << std::endl;
      }
    }

    // scaling vectors
    double NormInf;
    LHS->NormInf(&NormInf);
    LHS->Scale(1.0 / NormInf);

    // ========================================================= //
    // copy tentative and computed null space into NewNullSpace, //
    // which now becomes the standard null space                 //
    // ========================================================= //

    int NewNullSpaceSize = OldNullSpaceSize + 1;
    NewNullSpace = new double[NumMyRows() * NewNullSpaceSize];
    assert (NewNullSpace != 0);
    int itmp = OldNullSpaceSize * NumMyRows();
    for (int i = 0 ; i < itmp ; ++i) {
      NewNullSpace[i] = OldNullSpace[i];
    }

    for (int j = 0 ; j < NumMyRows() ; ++j) {
      NewNullSpace[itmp + j] = (*LHS)[j];
    }

    // =============== //
    // visualize modes //
    // =============== //

    if (List_.get("adaptive: visualize", false)) {

      double* x_coord = List_.get("viz: x-coordinates", (double*)0);
      double* y_coord = List_.get("viz: y-coordinates", (double*)0);
      double* z_coord = List_.get("viz: z-coordinates", (double*)0);
      assert (x_coord != 0);

      std::vector<double> plot_me(NumMyRows()/NumPDEEqns_);
      ML_Aggregate_Viz_Stats info;
      info.Amatrix = &(ml_->Amat[LevelID_[0]]);
      info.x = x_coord;
      info.y = y_coord;
      info.z = z_coord;
      info.Nlocal = NumMyRows() / NumPDEEqns_;
      info.Naggregates = 1;
      ML_Operator_AmalgamateAndDropWeak(&(ml_->Amat[LevelID_[0]]),
                                        NumPDEEqns_, 0.0);

      for (int ieqn = 0 ; ieqn < NumPDEEqns_ ; ++ieqn) {
        for (int j = 0 ; j < NumMyRows() ; j+=NumPDEEqns_) {
          plot_me[j / NumPDEEqns_] = (*LHS)[j + ieqn];
        }
        char FileName[80];
        sprintf(FileName,"nullspace-mode%d-eq%d.xyz", istep, ieqn);
        if (verbose_)
          std::cout << PrintMsg_ << "writing file " << FileName << "..." << std::endl;
        ML_Aggregate_VisualizeXYZ(info,FileName,
                                  ml_->comm,&plot_me[0]);
      }

      ML_Operator_UnAmalgamateAndDropWeak(&(ml_->Amat[LevelID_[0]]),
                                          NumPDEEqns_, 0.0);
    }
    
    // Destroy the old preconditioner
    DestroyPreconditioner();

    // ==================================================== //
    // build the new preconditioner with the new null space //
    // ==================================================== //

    List_.set("null space: type", "pre-computed");
    List_.set("null space: dimension", NewNullSpaceSize);
    List_.set("null space: vectors", NewNullSpace);

    ML_CHK_ERR(ComputePreconditioner());

    if (istep && (istep != NumAdaptiveVectors))
      delete OldNullSpace;

    OldNullSpace = NewNullSpace;
    OldNullSpaceSize = NewNullSpaceSize;

  }

  // keep trace of this pointer, it will be delete'd later
  NullSpaceToFree_ = NewNullSpace;

  delete LHS;
  delete RHS;

  return(0);

}
// ================================================ ====== ==== ==== == =
// Computes the preconditioner
int ML_Epetra::FaceMatrixFreePreconditioner::ComputePreconditioner(const bool CheckFiltering)
{
  Teuchos::ParameterList & ListCoarse=List_.sublist("face matrix free: coarse");

  /* ML Communicator */
  ML_Comm_Create(&ml_comm_);

  /* Parameter List Options */
  int SmootherSweeps = List_.get("smoother: sweeps", 0);
  MaxLevels = List_.get("max levels",10);
  print_hierarchy= List_.get("print hierarchy",false);
  num_cycles  = List_.get("cycle applications",1);

  /* Sanity Checking*/
  int OperatorDomainPoints =  OperatorDomainMap().NumGlobalPoints();
  int OperatorRangePoints =  OperatorRangeMap().NumGlobalPoints();
  if (OperatorDomainPoints != OperatorRangePoints)
    ML_CHK_ERR(-1); // only square matrices

 /* Output Header */
  if(verbose_ && !Comm_->MyPID()) {
    printf("------------------------------------------------------------------------------\n");
    printf("***\n");
    printf("*** ML_Epetra::FaceMatrixFreePreconditioner [%s]\n",Label());
    printf("***\n");
  }

  /* Invert non-zeros on the diagonal */
  if(SmootherSweeps){
    for (int i = 0; i < InvDiagonal_->MyLength(); ++i)
      if ((*InvDiagonal_)[i] != 0.0)
        (*InvDiagonal_)[i] = 1.0 / (*InvDiagonal_)[i];
    double nrm;
    InvDiagonal_->Norm2(&nrm);
    if(verbose_ && !Comm_->MyPID()) printf("Inverse Diagonal Norm = %6.4e\n",nrm);
  }/*end if*/

  /* Do the eigenvalue estimation for Chebyshev */
  if(SmootherSweeps)
    ML_CHK_ERR(SetupSmoother());


  if(MaxLevels > 0) {
    /* Build the prolongator */
    ML_CHK_ERR(BuildProlongator());

    /* DEBUG: Output matrices */
    if(print_hierarchy) EpetraExt::RowMatrixToMatlabFile("prolongator.dat",*Prolongator_);

    /* Form the coarse matrix */
    ML_CHK_ERR(FormCoarseMatrix());

    /* DEBUG: Output matrices */
    if(print_hierarchy) EpetraExt::RowMatrixToMatlabFile("coarsemat.dat",*CoarseMatrix);

    /* Setup Preconditioner on Coarse Matrix */
    CoarsePC = new MultiLevelPreconditioner(*CoarseMatrix,ListCoarse);
    if(!CoarsePC) ML_CHK_ERR(-2);

    /* Clean Up */
    //delete nullspace;
  }/*end if*/

  return 0;
}/*end ComputePreconditioner*/
Exemple #16
0
  Epetra_Operator* ML_Gen_Smoother_Ifpack_Epetra(const Epetra_Operator *A,const Epetra_Vector *InvDiagonal, Teuchos::ParameterList & List,string printMsg,bool verbose){
  /* Variables */
  double lambda_min = 0.0;
  double lambda_max = 0.0;
  Teuchos::ParameterList IFPACKList=List.sublist("smoother: ifpack list");

  /* Parameter-list Options */
  string PreOrPostSmoother = List.get("smoother: pre or post","both");
  string SmooType = List.get("smoother: type", "Chebyshev");
  if(SmooType=="IFPACK") SmooType=List.get("smoother: ifpack type","Chebyshev");
  int Sweeps = List.get("smoother: sweeps", 3);
  int IfpackOverlap = List.get("smoother: ifpack overlap",0);
  double omega = List.get("smoother: damping factor",1.0);
  Ifpack_Chebyshev* SmootherC_=0;
  Ifpack_Preconditioner* SmootherP_=0;

  /* Sanity Check*/
  if(Sweeps==0) return 0;

  /* Early Output*/
  int Nrows=A->OperatorDomainMap().NumGlobalElements();
  int Nnz=-1;
  const Epetra_RowMatrix *Arow=dynamic_cast<const Epetra_RowMatrix*>(A);
  if(Arow) Nnz=Arow->NumGlobalNonzeros();

  if(verbose && !A->Comm().MyPID())
    cout <<printMsg<<"# global rows = "<<Nrows<<" # estim. global nnz = "<<Nnz<<endl;


  /**********************************************/
  /***      Chebyshev (Including Block)       ***/
  /**********************************************/
  if(SmooType=="Chebyshev" || SmooType=="MLS" || SmooType=="IFPACK-Chebyshev" || SmooType=="IFPACK-Block Chebyshev"){
    bool allocated_inv_diagonal=false;
    int MaximumIterations = List.get("eigen-analysis: max iters", 10);
    string EigenType_ = List.get("eigen-analysis: type", "cg");
    double boost = List.get("eigen-analysis: boost for lambda max", 1.0);
    double alpha = List.get("chebyshev: alpha",30.0001);  
    Epetra_Vector *InvDiagonal_=0;

    /* Block Chebyshev stuff if needed */    
    int  MyCheby_nBlocks=     List.get("smoother: Block Chebyshev number of blocks",0);
    int* MyCheby_blockIndices=List.get("smoother: Block Chebyshev block list",(int*)0);
    int* MyCheby_blockStarts= List.get("smoother: Block Chebyshev block starts",(int*)0);
    bool MyCheby_NE=          List.get("smoother: chebyshev solve normal equations",false);

    if(SmooType == "IFPACK-Block Chebyshev" && MyCheby_blockIndices && MyCheby_blockStarts){
      // If we're using Block Chebyshev, it can compute it's own eigenvalue estimate..
      Teuchos::ParameterList PermuteList,BlockList;
      BlockList.set("apply mode","invert");
      PermuteList.set("number of local blocks",MyCheby_nBlocks);
      PermuteList.set("block start index",MyCheby_blockStarts);
      //        if(is_lid) PermuteList.set("block entry lids",Blockids_);
      //NTS: Add LID support
      PermuteList.set("block entry gids",MyCheby_blockIndices);        
      PermuteList.set("blockdiagmatrix: list",BlockList);
      
      IFPACKList.set("chebyshev: use block mode",true);
      IFPACKList.set("chebyshev: block list",PermuteList);
      IFPACKList.set("chebyshev: eigenvalue max iterations",10);
      
      // EXPERIMENTAL: Cheby-NE
      IFPACKList.set("chebyshev: solve normal equations",MyCheby_NE);
    }
    else {
      /* Non-Blocked Chebyshev */
      /* Grab Diagonal & invert if not provided */
      if(InvDiagonal) InvDiagonal_=const_cast<Epetra_Vector *>(InvDiagonal);
      else{
	const Epetra_CrsMatrix* Acrs=dynamic_cast<const Epetra_CrsMatrix*>(A);
	if(!Acrs) return 0;
	allocated_inv_diagonal=true;
	InvDiagonal_ = new Epetra_Vector(Acrs->RowMap());  
	Acrs->ExtractDiagonalCopy(*InvDiagonal_);
	for (int i = 0; i < InvDiagonal_->MyLength(); ++i)
	  if ((*InvDiagonal_)[i] != 0.0)
	    (*InvDiagonal_)[i] = 1.0 / (*InvDiagonal_)[i];   
      }

      /* Do the eigenvalue estimation*/
      if (EigenType_ == "power-method") Ifpack_Chebyshev::PowerMethod(*A,*InvDiagonal_,MaximumIterations,lambda_max);
      else if(EigenType_ == "cg") Ifpack_Chebyshev::CG(*A,*InvDiagonal_,MaximumIterations,lambda_min,lambda_max);
      else ML_CHK_ERR(0); // not recognized
    
      lambda_min=lambda_max / alpha;      

      /* Setup the Smoother's List*/
      IFPACKList.set("chebyshev: min eigenvalue", lambda_min);
      IFPACKList.set("chebyshev: max eigenvalue", boost * lambda_max);
      IFPACKList.set("chebyshev: ratio eigenvalue",alpha);
      IFPACKList.set("chebyshev: operator inv diagonal", InvDiagonal_);
    }

    /* Setup the Smoother's List*/    
    IFPACKList.set("chebyshev: ratio eigenvalue", alpha);
    IFPACKList.set("chebyshev: degree", Sweeps);
    IFPACKList.set("chebyshev: zero starting solution",false);

    // Setup
    SmootherC_= new Ifpack_Chebyshev(A);
    if (SmootherC_ == 0) return 0;
    SmootherC_->SetParameters(IFPACKList);
    SmootherC_->Initialize();
    SmootherC_->Compute();

    // Grab the lambda's if needed
    if(SmooType=="IFPACK-Block Chebyshev"){
      lambda_min=SmootherC_->GetLambdaMin();
      lambda_max=SmootherC_->GetLambdaMax();
    }    

    // Smoother Info Output
    if(verbose && !A->Comm().MyPID()){
      if(SmooType=="IFPACK-Block Chebyshev") {
	cout << printMsg << "MLS/Block-Chebyshev, polynomial order = "
	     <<  Sweeps
	     << ", alpha = " << alpha << endl;
	cout << printMsg << "lambda_min = " << lambda_min
	     << ", lambda_max = " << lambda_max << endl;
      }
      else {
	cout << printMsg << "MLS/Chebyshev, polynomial order = "
	     <<  Sweeps
	     << ", alpha = " << alpha << endl;
	cout << printMsg << "lambda_min = " << lambda_min
             << ", lambda_max = " << boost*lambda_max << endl;
      }
    }


    // Cleanup:  Since Chebyshev will keep it's own copy of the Inverse Diagonal...
    if (allocated_inv_diagonal) delete InvDiagonal_;

    return SmootherC_;
  }
  /**********************************************/
  /***           Point Relaxation             ***/
  /**********************************************/
  else if(SmooType=="Gauss-Seidel" || SmooType=="symmetric Gauss-Seidel" || SmooType=="Jacobi"
	  || SmooType=="point relaxation stand-alone" || SmooType=="point relaxation" ){      
    const Epetra_CrsMatrix* Acrs=dynamic_cast<const Epetra_CrsMatrix*>(A);
    if(!Acrs) return 0;
    string MyIfpackType="point relaxation stand-alone";
    if(IfpackOverlap > 0) MyIfpackType="point relaxation";
    string MyRelaxType="symmetric Gauss-Seidel";
    if(SmooType=="symmetric Gauss-Seidel") MyRelaxType=SmooType;
    else if(SmooType=="Jacobi") MyRelaxType=SmooType;

    IFPACKList.set("relaxation: type", IFPACKList.get("relaxation: type",MyRelaxType));
    IFPACKList.set("relaxation: sweeps", Sweeps);
    IFPACKList.set("relaxation: damping factor", omega);
    IFPACKList.set("relaxation: zero starting solution",false);

    if(verbose && !A->Comm().MyPID()){
      cout << printMsg << IFPACKList.get("relaxation: type",MyRelaxType).c_str()<<" (sweeps="
	   << Sweeps << ",omega=" << omega <<  ")" <<endl;
    }
    	
    Ifpack Factory;
    SmootherP_ = Factory.Create(MyIfpackType,const_cast<Epetra_CrsMatrix*>(Acrs),IfpackOverlap);
    if (SmootherP_ == 0) return 0;
    SmootherP_->SetParameters(IFPACKList);
    SmootherP_->Initialize();
    SmootherP_->Compute();
    return SmootherP_;
  }
  /**********************************************/
  /***           Block Relaxation             ***/
  /**********************************************/
  else if(SmooType=="block Gauss-Seidel" || SmooType=="symmetric block Gauss-Seidel" || SmooType=="block Jacobi" 
	  || SmooType=="block relaxation stand-alone" || SmooType=="block relaxation" ){      
    const Epetra_CrsMatrix* Acrs=dynamic_cast<const Epetra_CrsMatrix*>(A);
    if(!Acrs) return 0;
    string MyIfpackType="block relaxation stand-alone";
    if(IfpackOverlap > 0) MyIfpackType="block relaxation";
    string MyRelaxType="symmetric Gauss-Seidel";
    if(SmooType=="block Gauss-Seidel") MyRelaxType="Gauss-Seidel";
    else if(SmooType=="block Jacobi") MyRelaxType="Jacobi";

    IFPACKList.set("relaxation: type", IFPACKList.get("relaxation: type",MyRelaxType));
    IFPACKList.set("relaxation: sweeps", Sweeps);
    IFPACKList.set("relaxation: damping factor", omega);
    IFPACKList.set("relaxation: zero starting solution",false);
   
    if(verbose && !A->Comm().MyPID()){
      cout << printMsg << "block " << IFPACKList.get("relaxation: type",MyRelaxType).c_str()<<" (sweeps="
	   << Sweeps << ",omega=" << omega << ")" <<endl;
    }

    Ifpack Factory;
    SmootherP_ = Factory.Create(MyIfpackType,const_cast<Epetra_CrsMatrix*>(Acrs),IfpackOverlap);
    if (SmootherP_ == 0) return 0;
    SmootherP_->SetParameters(IFPACKList);
    SmootherP_->Initialize();
    SmootherP_->Compute();
    return SmootherP_;
  }
  /**********************************************/
  /***        Incomplete Factorization        ***/
  /**********************************************/
  else if(SmooType == "ILU" || SmooType == "IC" || SmooType == "ILUT"   ||
	  SmooType == "ICT" || SmooType == "SILU") {
    const Epetra_RowMatrix* Arow=dynamic_cast<const Epetra_RowMatrix*>(A);
    double MyLOF=0.0;
    if(SmooType=="ILUT" || SmooType=="ICT") MyLOF=List.get("smoother: ifpack level-of-fill",1.0);
    else MyLOF=List.get("smoother: ifpack level-of-fill",0.0);

    int MyIfpackOverlap = List.get("smoother: ifpack overlap", 0);
    double MyIfpackRT = List.get("smoother: ifpack relative threshold", 1.0);
    double MyIfpackAT = List.get("smoother: ifpack absolute threshold", 0.0);
    IFPACKList.set("ILU: sweeps",Sweeps);
    
    // Set the fact: LOF options, but only if they're not set already... All this sorcery is because level-of-fill
    // is an int for ILU and a double for ILUT.  Lovely.
    if(SmooType=="ILUT" || SmooType=="ICT"){
	IFPACKList.set("fact: level-of-fill", IFPACKList.get("fact: level-of-fill",MyLOF));
	IFPACKList.set("fact: ilut level-of-fill", IFPACKList.get("fact: ilut level-of-fill",MyLOF));
	IFPACKList.set("fact: ict level-of-fill", IFPACKList.get("fact: ict level-of-fill",MyLOF));
	MyLOF=IFPACKList.get("fact: level-of-fill",MyLOF);
    }
    else{
      IFPACKList.set("fact: level-of-fill", (int) IFPACKList.get("fact: level-of-fill",(int)MyLOF));
      MyLOF=IFPACKList.get("fact: level-of-fill",(int)MyLOF);
    }
    
    IFPACKList.set("fact: relative threshold", MyIfpackRT);
    IFPACKList.set("fact: absolute threshold", MyIfpackAT);

    if(verbose && !A->Comm().MyPID()){
      cout << printMsg << "IFPACK, type=`" << SmooType << "'," << endl
	   << printMsg << PreOrPostSmoother  << ",overlap=" << MyIfpackOverlap << endl;
      cout << printMsg << "level-of-fill=" << MyLOF;
      cout << ",rel. threshold=" << MyIfpackRT
	   << ",abs. threshold=" << MyIfpackAT << endl;
    }

    Ifpack Factory;
    SmootherP_ = Factory.Create(SmooType,const_cast<Epetra_RowMatrix*>(Arow),IfpackOverlap);
    if (SmootherP_ == 0) return 0;
    SmootherP_->SetParameters(IFPACKList);
    SmootherP_->Initialize();
    SmootherP_->Compute();
    return SmootherP_;
  }
  /**********************************************/
  /***                  SORa                  ***/
  /**********************************************/
  else if(SmooType=="SORa"){
    const Epetra_RowMatrix* Arow=dynamic_cast<const Epetra_RowMatrix*>(A);
    if(verbose && !A->Comm().MyPID()){
      cout << printMsg << "IFPACK/SORa("<<IFPACKList.get("sora: alpha",1.5)<<","<<IFPACKList.get("sora: gamma",1.0)<<")"
	   << ", sweeps = " <<IFPACKList.get("sora: sweeps",1)<<endl;
      if(IFPACKList.get("sora: oaz boundaries",false))
	cout << printMsg << "oaz boundary handling enabled"<<endl;
      if(IFPACKList.get("sora: use interproc damping",false))
	cout << printMsg << "interproc damping enabled"<<endl;
      if(IFPACKList.get("sora: use global damping",false))
	cout << printMsg << "global damping enabled"<<endl;
    }
    Ifpack Factory;
    SmootherP_ = Factory.Create(SmooType,const_cast<Epetra_RowMatrix*>(Arow),IfpackOverlap);
    if (SmootherP_ == 0) return 0;
    SmootherP_->SetParameters(IFPACKList);
    SmootherP_->Initialize();
    SmootherP_->Compute();
    return SmootherP_;
  }
  else{
    printf("ML_Gen_Smoother_Ifpack_New: Unknown preconditioner\n");
    ML_CHK_ERR(0);
  }
  return 0;

}/*ML_Gen_Smoother_Ifpack_New*/
// ============================================================================
int ML_Epetra::MultiLevelPreconditioner::
CreateAuxiliaryMatrixCrs(Epetra_FECrsMatrix* &FakeMatrix)
{

  int NumMyRows = RowMatrix_->NumMyRows();

  const Epetra_Map& RowMap = RowMatrix_->RowMatrixRowMap();
  const Epetra_Map& ColMap = RowMatrix_->RowMatrixColMap();
  FakeMatrix = new Epetra_FECrsMatrix(Copy,RowMap,
				      2*RowMatrix_->MaxNumEntries());

  if (FakeMatrix == 0)
    ML_CHK_ERR(-1); // something went wrong

  int NumDimensions = 0;

  double* x_coord = List_.get("x-coordinates", (double *)0);
  if (x_coord != 0) ++NumDimensions;

  // at least x-coordinates must be not null
  if( NumDimensions == 0 ) {
    std::cerr << ErrorMsg_ << "Option `aggregation: use auxiliary matrix' == true" << std::endl
         << ErrorMsg_ << "requires x-, y-, or z-coordinates." << std::endl
         << ErrorMsg_ << "You must specify them using options" << std::endl
         << ErrorMsg_ << "`x-coordinates' (and equivalently for" << std::endl
         << ErrorMsg_ << "y- and z-." << std::endl;
    ML_CHK_ERR(-2); // wrong parameters
  }

  double* y_coord = List_.get("y-coordinates", (double *)0);
  if (y_coord != 0) ++NumDimensions;

  double* z_coord = List_.get("z-coordinates", (double *)0);
  if (z_coord != 0) ++NumDimensions;

  // small check to avoid strange behavior
  if( z_coord != 0 && y_coord == 0 ) {
    std::cerr << ErrorMsg_ << "Something wrong: `y-coordinates'" << std::endl
         << ErrorMsg_ << "is null, while `z-coordinates' is null" << std::endl;
    ML_CHK_ERR(-3); // something went wrong
  }

  double theta = List_.get("aggregation: theta",0.0);

  bool SymmetricPattern = List_.get("aggregation: use symmetric pattern",false);

  // usual crap to clutter the output
  if( verbose_ ) {
    std::cout << std::endl;
    std::cout << PrintMsg_ << "*** Using auxiliary matrix to create the aggregates" << std::endl
         << PrintMsg_ << "*** Number of dimensions = " << NumDimensions << std::endl
         << PrintMsg_ << "*** theta = " << theta;
    if( SymmetricPattern ) std::cout << ", using symmetric pattern" << std::endl;
    else                   std::cout << ", using original pattern" << std::endl;
    std::cout << std::endl;
  }

  // create vectors containing coordinates, replicated for all unknonws
  // FIXME: I don't really need Z in all cases
  // for west
  // I am over-allocating, for large number of equations per node
  // this is not optimal. However, it is a only-once importing
  // of some more data. It should harm too much...
  //
  // The following will work for constant number of equations per
  // node only. The fix should be easy, though.
  Epetra_Vector RowX(RowMap); RowX.PutScalar(0.0);
  Epetra_Vector RowY(RowMap); RowY.PutScalar(0.0);
  Epetra_Vector RowZ(RowMap); RowZ.PutScalar(0.0);

  for (int i = 0 ; i < NumMyRows ; i += NumPDEEqns_) {
    RowX[i] = x_coord[i / NumPDEEqns_];
    if (NumDimensions > 1) RowY[i] = y_coord[i / NumPDEEqns_];
    if (NumDimensions > 2) RowZ[i] = z_coord[i / NumPDEEqns_];
  }

  // create vectors containing coordinates for columns
  // (this is useful only if MIS/ParMETIS are used)
  Epetra_Vector ColX(ColMap); ColX.PutScalar(0.0);
  Epetra_Vector ColY(ColMap); ColY.PutScalar(0.0);
  Epetra_Vector ColZ(ColMap); ColZ.PutScalar(0.0);

  // get coordinates for non-local nodes (in column map)
  Epetra_Import Importer(ColMap,RowMap);
  ColX.Import(RowX,Importer,Insert);
  if (NumDimensions > 1) ColY.Import(RowY,Importer,Insert);
  if (NumDimensions > 2) ColZ.Import(RowZ,Importer,Insert);

  // global row and column numbering
  int* MyGlobalRowElements = RowMap.MyGlobalElements();
  int* MyGlobalColElements = ColMap.MyGlobalElements();

  // room for getrow()
  int MaxNnz = RowMatrix_->MaxNumEntries();
  std::vector<int> colInd(MaxNnz);
  std::vector<double> colVal(MaxNnz);
  std::vector<double> coord_i(3);
  std::vector<double> coord_j(3);

  // =================== //
  // cycle over all rows //
  // =================== //

  for (int i = 0; i < NumMyRows ; i += NumPDEEqns_) {

    int GlobalRow = MyGlobalRowElements[i];

    if( i%NumPDEEqns_ == 0 ) { // do it just once for each block row
      switch( NumDimensions ) {
      case 3:
	coord_i[2] = RowZ[i];
      case 2:
	coord_i[1] = RowY[i];
      case 1:
	coord_i[0] = RowX[i];

      }

      int NumEntries;
      ML_CHK_ERR(RowMatrix_->ExtractMyRowCopy(i,MaxNnz,NumEntries,
                                              &colVal[0],&colInd[0]));

      // NOTE: for VBR matrices, the "real" value that will be used in
      // the subsequent part of the code is only the one for the first
      // equations. For each block, I replace values with the sum of
      // the std::abs of each block entry.
      for (int j = 0 ; j < NumEntries ; j += NumPDEEqns_) {
	colVal[j] = std::fabs(colVal[j]);
	for (int k = 1 ; k < NumPDEEqns_ ; ++k) {
	  colVal[j] += std::fabs(colVal[j+k]);
	}
      }

      // work only on the first equations. Theta will blend the
      // coordinate part with the sub of std::abs of row elements.
      int GlobalCol;
      double total = 0.0;

      for (int j = 0 ; j < NumEntries ; j += NumPDEEqns_) {

        if (colInd[j] >= NumMyRows)
          continue;

	if (colInd[j]%NumPDEEqns_ == 0) {

	  // insert diagonal later
	  if (colInd[j] != i) {

	    // get coordinates of this node
	    switch (NumDimensions) {
	    case 3:
	      coord_j[2] = ColZ[colInd[j]];
	    case 2:
	      coord_j[1] = ColY[colInd[j]];
	    case 1:
	      coord_j[0] = ColX[colInd[j]];
	    }

	    // d2 is the square of the distance between node `i' and
	    // node `j'
	    double d2 = (coord_i[0] - coord_j[0]) * (coord_i[0] - coord_j[0]) +
	      (coord_i[1] - coord_j[1]) * (coord_i[1] - coord_j[1]) +
	      (coord_i[2] - coord_j[2]) * (coord_i[2] - coord_j[2]);

	    if (d2 == 0.0) {
	      std::cerr << std::endl;
	      std::cerr << ErrorMsg_ << "distance between node " << i/NumPDEEqns_ << " and node "
                   << colInd[j]/NumPDEEqns_ << std::endl
                   << ErrorMsg_ << "is zero. Coordinates of these nodes are" << std::endl
	           << ErrorMsg_ << "x_i = " << coord_i[0] << ", x_j = " << coord_j[0] << std::endl
		   << ErrorMsg_ << "y_i = " << coord_i[1] << ", y_j = " << coord_j[1] << std::endl
		   << ErrorMsg_ << "z_i = " << coord_i[2] << ", z_j = " << coord_j[2] << std::endl
		   << ErrorMsg_ << "Now proceeding with distance = 1.0" << std::endl;
	      std::cerr << std::endl;
	      d2 = 1.0;
	    }

	    // blend d2 with the actual values of the matrix
	    // FIXME: am I useful?
	    double val = -(1.0 - theta) * (1.0 / d2) + theta * (colVal[j]);

	    GlobalCol = MyGlobalColElements[colInd[j]];

	    // insert this value on all rows
	    for (int k = 0 ; k < NumPDEEqns_ ; ++k) {
	      int row = GlobalRow + k;
	      int col = GlobalCol + k;

	      if (FakeMatrix->SumIntoGlobalValues(1,&row,1,&col,&val) != 0) {
		ML_CHK_ERR(FakeMatrix->InsertGlobalValues(1,&row,1,&col,&val));
	      }

	    }

	    total -= val;

	    // put (j,i) element as well, only for in-process stuff.
	    // I have some problems with off-processor elements.
            // It is here that I need the FE matrix.
	    if (SymmetricPattern == true && colInd[j] < NumMyRows ) {

	      for( int k=0 ; k<NumPDEEqns_ ; ++k ) {
		int row = GlobalCol+k;
		int col = GlobalRow+k;

		if( FakeMatrix->SumIntoGlobalValues(1,&row,1,&col,&val) != 0 ) {
		  ML_CHK_ERR(FakeMatrix->InsertGlobalValues(1,&row,1,&col,&val));
		}

	      }
	      total -= val;
	    }
	  }
	}
      }

      // create lines with zero-row sum
      for (int k = 0 ; k < NumPDEEqns_ ; ++k) {
	int row = GlobalRow + k;
	if (FakeMatrix->SumIntoGlobalValues(1,&row,1,&row,&total) != 0) {
	  if (FakeMatrix->InsertGlobalValues(1,&row,1,&row,&total) != 0)
	    ML_CHK_ERR(-9); // something went wrong
	}

      }
    }
  }

  if (FakeMatrix->FillComplete())
    ML_CHK_ERR(-5); // something went wrong

  // stick pointer in Amat for level 0 (finest level)
  ml_->Amat[LevelID_[0]].data = (void *)FakeMatrix;

  // tell ML to keep the tentative prolongator
  ML_Aggregate_Set_Reuse(agg_);

  // pray that no bugs will tease us
  return(0);
}