// ====================================================================== 
Operator GetJacobiIterationOperator(const Operator& Amat, double Damping)
{

  struct ML_AGG_Matrix_Context* widget = new struct ML_AGG_Matrix_Context;
  widget->near_bdry = 0;
  widget->aggr_info = 0;
  widget->drop_tol  = 0.0;

  widget->Amat = Amat.GetML_Operator();
  widget->omega = Damping;

  ML_Operator* tmp_ML = ML_Operator_Create(GetML_Comm());
  ML_Operator_Set_ApplyFuncData(tmp_ML, widget->Amat->invec_leng,
                                widget->Amat->outvec_leng, widget,
                                widget->Amat->matvec->Nrows, NULL, 0);

  tmp_ML->data_destroy = widget_destroy;

  ML_Operator_Set_Getrow(tmp_ML, widget->Amat->getrow->Nrows, 
                         ML_AGG_JacobiSmoother_Getrows);

  // Creates a new copy of pre_comm, so that the old pre_comm
  // can be destroyed without worry
  ML_CommInfoOP_Clone(&(tmp_ML->getrow->pre_comm),
                      widget->Amat->getrow->pre_comm);

  Operator tmp(Amat.GetDomainSpace(), Amat.GetRangeSpace(), tmp_ML, true,
               Amat.GetRCPOperatorBox());

  return(tmp);
}
// ====================================================================== 
Operator GetIdentity(const Space& DomainSpace, const Space& RangeSpace)
{
  ML_Operator* ML_eye = ML_Operator_Create(GetML_Comm());
  int size = DomainSpace.GetNumMyElements();
  ML_Operator_Set_ApplyFuncData(ML_eye, size, size,
            NULL, size, eye_matvec, 0);
  ML_Operator_Set_Getrow(ML_eye, size, eye_getrows);
  Operator eye(DomainSpace,DomainSpace,ML_eye,true);
  return(eye);
}
// ====================================================================== 
Operator GetTranspose(const Operator& A, const bool byrow = true) 
{
  ML_Operator* ML_transp;
  ML_transp = ML_Operator_Create(GetML_Comm());
  if (byrow)
    ML_Operator_Transpose_byrow(A.GetML_Operator(),ML_transp);
  else
    ML_Operator_Transpose(A.GetML_Operator(),ML_transp);

  Operator transp(A.GetRangeSpace(),A.GetDomainSpace(), ML_transp,true);
  return(transp);
}
// ====================================================================== 
Operator GetPtent1D(const MultiVector& D, const int offset = 0)
{
  if (D.GetNumVectors() != 1)
    ML_THROW("D.GetNumVectors() != 1", -1);

  int size = D.GetMyLength();
  if (size == 0)
    ML_THROW("empty diagonal vector in input", -1);

  double* diag = new double[size];
  for (int i = 0 ; i < size ; ++i)
    diag[i] = D(i);

  // creates the ML operator and store the diag pointer,
  // as well as the function pointers
  ML_Operator* MLDiag = ML_Operator_Create(GetML_Comm());

  int invec_leng = size / 3 + size % 3;
  int outvec_leng = size;

  MLDiag->invec_leng = invec_leng;
  MLDiag->outvec_leng = outvec_leng;
  MLDiag->data = (void*)diag;
  MLDiag->data_destroy = Ptent1D_destroy;
  MLDiag->matvec->func_ptr = Ptent1D_matvec;

  MLDiag->matvec->ML_id = ML_NONEMPTY;
  MLDiag->matvec->Nrows = outvec_leng;
  MLDiag->from_an_ml_operator = 0;

  MLDiag->getrow->func_ptr = Ptent1D_getrows;

  MLDiag->getrow->ML_id = ML_NONEMPTY;
  MLDiag->getrow->Nrows = outvec_leng;

  // creates the domain space
  vector<int> MyGlobalElements(invec_leng);
  for (int i = 0 ; i < invec_leng ; ++i) 
    MyGlobalElements[i] = D.GetVectorSpace()(i * 3) / 3;
  Space DomainSpace(invec_leng, -1, &MyGlobalElements[0]);
  Space RangeSpace = D.GetVectorSpace();

  // creates the MLAPI wrapper
  Operator Diag(DomainSpace,RangeSpace,MLDiag,true);
  return(Diag);
}
Esempio n. 5
0
// ====================================================================== 
double MaxEigPowerMethod(const Operator& Op, const bool DiagonalScaling) 
{

  ML_Krylov *kdata;
  double MaxEigen;

    kdata = ML_Krylov_Create(GetML_Comm());
    if (DiagonalScaling == false)
      kdata->ML_dont_scale_by_diag = ML_TRUE;
    else
      kdata->ML_dont_scale_by_diag = ML_FALSE;
    ML_Krylov_Set_PrintFreq(kdata, 0);
    ML_Krylov_Set_ComputeNonSymEigenvalues(kdata);
    ML_Krylov_Set_Amatrix(kdata, Op.GetML_Operator());
    ML_Krylov_Solve(kdata, Op.GetML_Operator()->outvec_leng, NULL, NULL);
    MaxEigen = ML_Krylov_Get_MaxEigenvalue(kdata);
  ML_Krylov_Destroy(&kdata);

  return(MaxEigen);
}
// ====================================================================== 
Operator GetDiagonal(const MultiVector& D)
{
  if (D.GetNumVectors() != 1)
    ML_THROW("D.GetNumVectors() != 1", -1);

  int size = D.GetMyLength();
  if (size == 0)
    ML_THROW("empty diagonal vector in input", -1);

  double* diag = new double[size];
  for (int i = 0 ; i < size ; ++i) 
    diag[i] = D(i);

  // creates the ML operator and store the diag pointer,
  // as well as the function pointers
  ML_Operator* MLDiag = ML_Operator_Create(GetML_Comm());

  MLDiag->invec_leng = size;
  MLDiag->outvec_leng = size;
  MLDiag->data = (void*)diag;
  MLDiag->matvec->func_ptr = diag_matvec;

  MLDiag->matvec->ML_id = ML_NONEMPTY;
  MLDiag->matvec->Nrows = size;
  MLDiag->from_an_ml_operator = 0;
  MLDiag->data_destroy = diag_destroy;

  MLDiag->getrow->func_ptr = diag_getrows;

  MLDiag->getrow->ML_id = ML_NONEMPTY;
  MLDiag->getrow->Nrows = size;

  // creates the MLAPI wrapper
  Operator Diag(D.GetVectorSpace(),D.GetVectorSpace(),MLDiag,true);
  return(Diag);
}
Esempio n. 7
0
// ======================================================================
void GetPtent(const Operator& A, Teuchos::ParameterList& List,
              const MultiVector& ThisNS,
              Operator& Ptent, MultiVector& NextNS)
{
  std::string CoarsenType     = List.get("aggregation: type", "Uncoupled");
  /* old version
  int    NodesPerAggr    = List.get("aggregation: per aggregate", 64);
  */
  double Threshold       = List.get("aggregation: threshold", 0.0);
  int    NumPDEEquations = List.get("PDE equations", 1);

  ML_Aggregate* agg_object;
  ML_Aggregate_Create(&agg_object);
  ML_Aggregate_Set_MaxLevels(agg_object,2);
  ML_Aggregate_Set_StartLevel(agg_object,0);
  ML_Aggregate_Set_Threshold(agg_object,Threshold);
  //agg_object->curr_threshold = 0.0;

  ML_Operator* ML_Ptent = 0;
  ML_Ptent = ML_Operator_Create(GetML_Comm());

  if (ThisNS.GetNumVectors() == 0)
    ML_THROW("zero-dimension null space", -1);

  int size = ThisNS.GetMyLength();

  double* null_vect = 0;
  ML_memory_alloc((void **)(&null_vect), sizeof(double) * size * ThisNS.GetNumVectors(), "ns");

  int incr = 1;
  for (int v = 0 ; v < ThisNS.GetNumVectors() ; ++v)
    DCOPY_F77(&size, (double*)ThisNS.GetValues(v), &incr,
              null_vect + v * ThisNS.GetMyLength(), &incr);


  ML_Aggregate_Set_NullSpace(agg_object, NumPDEEquations,
                             ThisNS.GetNumVectors(), null_vect,
                             ThisNS.GetMyLength());

  if (CoarsenType == "Uncoupled")
    agg_object->coarsen_scheme = ML_AGGR_UNCOUPLED;
  else if (CoarsenType == "Uncoupled-MIS")
    agg_object->coarsen_scheme = ML_AGGR_HYBRIDUM;
  else if (CoarsenType == "MIS") {
   /* needed for MIS, otherwise it sets the number of equations to
    * the null space dimension */
    agg_object->max_levels  = -7;
    agg_object->coarsen_scheme = ML_AGGR_MIS;
  }
  else if (CoarsenType == "METIS")
    agg_object->coarsen_scheme = ML_AGGR_METIS;
  else {
    ML_THROW("Requested aggregation scheme (" + CoarsenType +
             ") not recognized", -1);
  }

  int NextSize = ML_Aggregate_Coarsen(agg_object, A.GetML_Operator(),
                                      &ML_Ptent, GetML_Comm());

  /* This is the old version
  int NextSize;

  if (CoarsenType == "Uncoupled") {
    NextSize = ML_Aggregate_CoarsenUncoupled(agg_object, A.GetML_Operator(),
  }
  else if (CoarsenType == "MIS") {
    NextSize = ML_Aggregate_CoarsenMIS(agg_object, A.GetML_Operator(),
                                       &ML_Ptent, GetML_Comm());
  }
  else if (CoarsenType == "METIS") {
    ML ml_object;
    ml_object.ML_num_levels = 1; // crap for line below
    ML_Aggregate_Set_NodesPerAggr(&ml_object,agg_object,0,NodesPerAggr);
    NextSize = ML_Aggregate_CoarsenMETIS(agg_object, A.GetML_Operator(),
                                         &ML_Ptent, GetML_Comm());
  }
  else {
    ML_THROW("Requested aggregation scheme (" + CoarsenType +
             ") not recognized", -1);
  }
  */

  ML_Operator_ChangeToSinglePrecision(ML_Ptent);

  int NumMyElements = NextSize;
  Space CoarseSpace(-1,NumMyElements);
  Ptent.Reshape(CoarseSpace,A.GetRangeSpace(),ML_Ptent,true);

  assert (NextSize * ThisNS.GetNumVectors() != 0);

  NextNS.Reshape(CoarseSpace, ThisNS.GetNumVectors());

  size = NextNS.GetMyLength();
  for (int v = 0 ; v < NextNS.GetNumVectors() ; ++v)
    DCOPY_F77(&size, agg_object->nullspace_vect + v * size, &incr,
              NextNS.GetValues(v), &incr);

  ML_Aggregate_Destroy(&agg_object);
  ML_memory_free((void**)(&null_vect));
}
Esempio n. 8
0
// ======================================================================
int GetAggregates(Epetra_RowMatrix& A, Teuchos::ParameterList& List,
                  double* thisns, Epetra_IntVector& aggrinfo)
{
  if (!A.RowMatrixRowMap().SameAs(aggrinfo.Map()))
    ML_THROW("map of aggrinfo must match row map of operator", -1);

  std::string CoarsenType     = List.get("aggregation: type", "Uncoupled");
  double Threshold       = List.get("aggregation: threshold", 0.0);
  int    NumPDEEquations = List.get("PDE equations", 1);
  int    nsdim           = List.get("null space: dimension",-1);
  if (nsdim==-1) ML_THROW("dimension of nullspace not set", -1);
  int size = A.RowMatrixRowMap().NumMyElements();

  ML_Aggregate* agg_object;
  ML_Aggregate_Create(&agg_object);
  ML_Aggregate_KeepInfo(agg_object,1);
  ML_Aggregate_Set_MaxLevels(agg_object,2);
  ML_Aggregate_Set_StartLevel(agg_object,0);
  ML_Aggregate_Set_Threshold(agg_object,Threshold);
  //agg_object->curr_threshold = 0.0;

  ML_Operator* ML_Ptent = 0;
  ML_Ptent = ML_Operator_Create(GetML_Comm());

  if (!thisns)
    ML_THROW("nullspace is NULL", -1);

  ML_Aggregate_Set_NullSpace(agg_object, NumPDEEquations,
                             nsdim, thisns,size);

  if (CoarsenType == "Uncoupled")
    agg_object->coarsen_scheme = ML_AGGR_UNCOUPLED;
  else if (CoarsenType == "Uncoupled-MIS")
    agg_object->coarsen_scheme = ML_AGGR_HYBRIDUM;
  else if (CoarsenType == "MIS") {
   /* needed for MIS, otherwise it sets the number of equations to
    * the null space dimension */
    agg_object->max_levels  = -7;
    agg_object->coarsen_scheme = ML_AGGR_MIS;
  }
  else if (CoarsenType == "METIS")
    agg_object->coarsen_scheme = ML_AGGR_METIS;
  else {
    ML_THROW("Requested aggregation scheme (" + CoarsenType +
             ") not recognized", -1);
  }

  ML_Operator* ML_A = ML_Operator_Create(GetML_Comm());
  ML_Operator_WrapEpetraMatrix(&A,ML_A);

  int NextSize = ML_Aggregate_Coarsen(agg_object, ML_A,
                                      &ML_Ptent, GetML_Comm());

  int* aggrmap = NULL;
  ML_Aggregate_Get_AggrMap(agg_object,0,&aggrmap);
  if (!aggrmap) ML_THROW("aggr_info not available", -1);

#if 0 // debugging
  fflush(stdout);
  for (int proc=0; proc<A.GetRowMatrix()->Comm().NumProc(); ++proc)
  {
    if (A.GetRowMatrix()->Comm().MyPID()==proc)
    {
      std::cout << "Proc " << proc << ":" << std::endl;
      std::cout << "aggrcount " << aggrcount << std::endl;
      std::cout << "NextSize " << NextSize << std::endl;
      for (int i=0; i<size; ++i)
        std::cout << "aggrmap[" << i << "] = " << aggrmap[i] << std::endl;
      fflush(stdout);
    }
    A.GetRowMatrix()->Comm().Barrier();
  }
#endif

  assert (NextSize * nsdim != 0);
  for (int i=0; i<size; ++i) aggrinfo[i] = aggrmap[i];

  ML_Aggregate_Destroy(&agg_object);

  return (NextSize/nsdim);
}