/*----------------------------------------------------------------------*
 |  Constructor (public)                                     m.gee 12/04|
 |  IMPORTANT:                                                          |
 |  No matter on which level we are here, the vector xfine is ALWAYS    |
 |  a fine grid vector here!                                            |
 *----------------------------------------------------------------------*/
ML_NOX::ML_Nox_MatrixfreeLevel::ML_Nox_MatrixfreeLevel(int level, int nlevel, int plevel, 
                                 ML* ml, ML_Aggregate* ag, Epetra_CrsMatrix** P,
                                 ML_NOX::Ml_Nox_Fineinterface& interface, const Epetra_Comm& comm,
                                 const Epetra_Vector& xfine, double fd_alpha, double fd_beta,
                                 bool fd_centered, bool isDiagonalOnly, int bsize) 
: fineinterface_(interface),
comm_(comm)                                                          
{
   level_          = level;
   nlevel_         = nlevel;
   ml_printlevel_  = plevel;
   ml_             = ml;
   ag_             = ag;
   fd_alpha_       = fd_alpha;
   fd_beta_        = fd_beta;
   fd_centered_    = fd_centered;
   isDiagonalOnly_ = isDiagonalOnly;
   A_              = 0;
   coarseinterface_= 0;
   bsize_          = bsize;
   

   // we need the graph of the operator on this level. On the fine grid we can just ask the
   // fineinterface for it, on the coarser levels it has to be extracted from the ML-hierachy
   if (level_==0)
   {
      // the Epetra_CrsGraph-copy-constructor shares data with the original one.
      // We want a really deep copy here so we cannot use it
      // graph_ will be given to the FiniteDifferencing class and will be destroyed by it
      graph_ = ML_NOX::deepcopy_graph(interface.getGraph());
   }
   else
   {
      // Note that ML has no understanding of global indices, so it makes up new GIDs
      // (This also holds for the Prolongators P)
      Epetra_CrsMatrix* tmpMat  = 0;
      int               maxnnz  = 0;
      double            cputime = 0.0;
      ML_Operator2EpetraCrsMatrix(&(ml_->Amat[level_]), tmpMat, maxnnz, false, cputime);
      // copy the graph
      double t0 = GetClock();
      graph_ = ML_NOX::deepcopy_graph(&(tmpMat->Graph()));
      // delete the copy of the Epetra_CrsMatrix
      if (tmpMat) delete tmpMat; tmpMat = 0;
      double t1 = GetClock();
      if (ml_printlevel_ > 0 && 0 == comm_.MyPID())
         cout << "matrixfreeML (level " << level_ << "): extraction/copy of Graph in " << cputime+t1-t0 << " sec\n"
              << "                        max-nonzeros in Graph: " << maxnnz << "\n";
   }
   
   // create this levels coarse interface
   coarseinterface_ = new ML_NOX::Nox_CoarseProblem_Interface(fineinterface_,level_,ml_printlevel_,
                                                              P,&(graph_->RowMap()),nlevel_);

   // restrict the xfine-vector to this level 
   Epetra_Vector* xthis = coarseinterface_->restrict_fine_to_this(xfine);
   if (!xthis)
   {
     cout << "**ERR**: ML_Epetra::ML_Nox_MatrixfreeLevel::ML_Nox_MatrixfreeLevel:\n"
          << "**ERR**: ML_Epetra::Nox_CoarseProblem_Interface::restrict_fine_to_this returned NULL on level " 
          << level_ << "\n"
          << "**ERR**: file/line: " << __FILE__ << "/" << __LINE__ << "\n"; throw -1;
   }

   Epetra_Vector* xc = new Epetra_Vector(graph_->RowMap(),false);
   // FIXME: after intesive testing, this test might be obsolet
#if 0
   bool samemap = xc->Map().PointSameAs(xthis->Map());
   if (samemap)
   {
#endif
      xc->Update(1.0,*xthis,0.0);
#if 0
   }
   else
   {
      cout << "**WRN** Maps are not equal in\n"
           << "**WRN** file/line: " << __FILE__ << "/" << __LINE__ << "\n";
      // import the xthis vector in the Map that ML produced for graph_
      Epetra_Import* importer = new Epetra_Import(graph_->RowMap(),xthis->Map());  
      int ierr = xc->Import(*xthis,*importer,Insert);
      if (ierr)
      {
         cout << "**ERR**: ML_Epetra::ML_Nox_MatrixfreeLevel::ML_Nox_MatrixfreeLevel:\n"
              << "**ERR**: export from xthis to xc returned err=" << ierr <<"\n"
              << "**ERR**: file/line: " << __FILE__ << "/" << __LINE__ << "\n"; throw -1;
      }
      if (importer) delete importer; importer = 0;
   }
#endif   
   if (xthis) delete xthis; xthis = 0;

   // create the coloring of the graph
   if (ml_printlevel_>0 && comm_.MyPID()==0) 
   {
      cout << "matrixfreeML (level " << level_ << "): Entering Coloring on level " << level_ << "\n";
      fflush(stdout);
   }
   double t0 = GetClock();
   colorMap_ = ML_NOX::ML_Nox_collapsedcoloring(graph_,bsize_,isDiagonalOnly,ml_printlevel_);
   if (!colorMap_) colorMap_ = ML_NOX::ML_Nox_standardcoloring(graph_,isDiagonalOnly);
   colorMapIndex_ = new EpetraExt::CrsGraph_MapColoringIndex(*colorMap_);
   colorcolumns_  = &(*colorMapIndex_)(*graph_);
   double t1 = GetClock();
   if (ml_printlevel_>0 && comm_.MyPID()==0)
   {
      cout << "matrixfreeML (level " << level_ << "): Proc " << comm_.MyPID() <<" Coloring time is " << (t1-t0) << " sec\n";
      fflush(stdout);
   }

   // construct the FiniteDifferenceColoring-Matrix
   if (ml_printlevel_>0 && comm_.MyPID()==0)
   {
      cout << "matrixfreeML (level " << level_ << "): Entering Construction FD-Operator on level " << level_ << "\n";
      fflush(stdout);
   }

   t0 = GetClock();

#if 1 // FD-operator with coloring
   FD_ = new NOX::EpetraNew::FiniteDifferenceColoring(*coarseinterface_,
                                                      *xc,
                                                      *graph_,
                                                      *colorMap_,
                                                      *colorcolumns_,
                                                      true,
                                                      isDiagonalOnly_,
                                                      fd_beta_,fd_alpha_);
#else // FD-operator without coloring
   FD_ = new NOX::EpetraNew::FiniteDifference(*coarseinterface_,
                                              *xc,
                                              *graph_,
                                              fd_beta_,fd_alpha_);
#endif
   // set differencing method
   if (fd_centered_)
      FD_->setDifferenceMethod(NOX::EpetraNew::FiniteDifferenceColoring::Centered);

   bool err = FD_->computeJacobian(*xc); 
   if (err==false)
   {
     cout << "**ERR**: ML_NOX::ML_Nox_MatrixfreeLevel::ML_Nox_MatrixfreeLevel:\n"
          << "**ERR**: NOX::Epetra::FiniteDifferenceColoring returned an error on level " 
          << level_ << "\n"
          << "**ERR**: file/line: " << __FILE__ << "/" << __LINE__ << "\n"; throw -1;
   }

   // print number of calls to the coarse interface
   if (ml_printlevel_>0 && comm_.MyPID()==0)
      cout << "matrixfreeML (level " << level_ << "): Calls to coarse-computeF in FD-Operator: "
           << coarseinterface_->numcallscomputeF() << "\n"; 

   t1 = GetClock();
   if (ml_printlevel_>0 && comm_.MyPID()==0)
   {
      cout << "matrixfreeML (level " << level_ << "): Proc " << comm_.MyPID() <<" colored Finite Differencing time is " << (t1-t0) << " sec\n";
      fflush(stdout);
   }

   // get ref to computed Epetra_CrsMatrix  
   A_ = dynamic_cast<Epetra_CrsMatrix*>(&(FD_->getUnderlyingMatrix()));

   // set counter for number of calls to the coarseinterface_->computeF back to zero
   coarseinterface_->resetnumcallscomputeF();   
   
   // tidy up 
   if (xc) delete xc; xc = 0;
   
   return;
}
/*----------------------------------------------------------------------*
 |  recreate this level (public)                             m.gee 01/05|
 |  this function assumes, that the graph of the fine level problem has |
 |  not changed since call to the constructor and therefore             |
 | the graph and it's coloring do not have to be recomputed             |  
 |  IMPORTANT:                                                          |
 |  No matter on which level we are here, the vector xfine is ALWAYS    |
 |  a fine grid vector here!                                            |
 *----------------------------------------------------------------------*/
bool ML_NOX::ML_Nox_MatrixfreeLevel::recreateLevel(int level, int nlevel, int plevel, 
                                                   ML* ml, ML_Aggregate* ag, Epetra_CrsMatrix** P,
                                                   ML_NOX::Ml_Nox_Fineinterface& interface, 
                                                   const Epetra_Comm& comm, const Epetra_Vector& xfine) 
{
   // make some tests
   if (level != level_)
   {
     cout << "**ERR**: ML_Epetra::ML_Nox_MatrixfreeLevel::recreateLevel:\n"
          << "**ERR**: level_ " << level_ << " not equal level " << level << "\n"
          << "**ERR**: file/line: " << __FILE__ << "/" << __LINE__ << "\n"; throw -1;
   }
   if (nlevel != nlevel_)
   {
     cout << "**ERR**: ML_Epetra::ML_Nox_MatrixfreeLevel::recreateLevel:\n"
          << "**ERR**: nlevel_ " << nlevel_ << " not equal nlevel " << nlevel << "\n" 
          << "**ERR**: file/line: " << __FILE__ << "/" << __LINE__ << "\n"; throw -1;
   }
   // printlevel might have changed
   ml_printlevel_ = plevel;
   ml_            = ml;
   ag_            = ag;
   destroyP();  // safer to use the new Ps
   setP(NULL);

   // we need the graph of the operator on this level. On the fine grid we can just ask the
   // fineinterface for it, on the coarser levels it has to be extracted from the ML-hierachy
   bool same;
   if (level_==0)
   {
      const Epetra_CrsGraph* graph = interface.getGraph();
      // check whether the old graph matches the new one
      same = compare_graphs(graph,graph_);
      destroyFD(); // we are here to recompute the FD-operator (this destroys graph_)
      graph_ = ML_NOX::deepcopy_graph(graph);
   }
   else
   {
      // Note that ML has no understanding of global indices, so it makes up new GIDs
      // (This also holds for the Prolongators P)
      Epetra_CrsMatrix* tmpMat  = 0;
      int               maxnnz  = 0;
      double            cputime = 0.0;
      ML_Operator2EpetraCrsMatrix(&(ml_->Amat[level_]), tmpMat, maxnnz, false, cputime);
      // get a view from the graph
      const Epetra_CrsGraph& graph = tmpMat->Graph();
      // compare the graph to the existing one
      same = compare_graphs(&graph,graph_);
      destroyFD(); // we are here to recompute the FD-operator (this destroys graph_)
      double t0 = GetClock();
      graph_ = ML_NOX::deepcopy_graph(&graph);
      // delete the copy of the Epetra_CrsMatrix
      if (tmpMat) delete tmpMat; tmpMat = 0;
      double t1 = GetClock();
      if (ml_printlevel_ > 0 && 0 == comm_.MyPID())
         cout << "matrixfreeML (level " << level_ << "): extraction/copy of Graph in " << cputime+t1-t0 << " sec\n"
              << "                        max-nonzeros in Graph: " << maxnnz << "\n";
   }
   
   // recreate this levels coarse interface
   if (same)
      coarseinterface_->recreate(ml_printlevel_,P,&(graph_->RowMap()));
   else
   {
      delete coarseinterface_;
      coarseinterface_ = new ML_NOX::Nox_CoarseProblem_Interface(fineinterface_,level_,ml_printlevel_,
                                                                 P,&(graph_->RowMap()),nlevel_);
   }
   
   // restrict the xfine-vector to this level 
   Epetra_Vector* xthis = coarseinterface_->restrict_fine_to_this(xfine);
   if (!xthis)
   {
     cout << "**ERR**: ML_Epetra::ML_Nox_MatrixfreeLevel::ML_Nox_MatrixfreeLevel:\n"
          << "**ERR**: ML_Epetra::Nox_CoarseProblem_Interface::restrict_fine_to_this returned NULL on level " 
          << level_ << "\n"
          << "**ERR**: file/line: " << __FILE__ << "/" << __LINE__ << "\n"; throw -1;
   }

   Epetra_Vector* xc = new Epetra_Vector(graph_->RowMap(),false);
   // FIXME: after intesive testing, this test might be obsolet
#if 0
   bool samemap = xc->Map().PointSameAs(xthis->Map());
   if (samemap)
   {
#endif
      xc->Update(1.0,*xthis,0.0);
#if 0
   }
   else
   {
      cout << "**WRN** Maps are not equal in\n"
           << "**WRN** file/line: " << __FILE__ << "/" << __LINE__ << "\n";
      // import the xthis vector in the Map that ML produced for graph_
      Epetra_Import* importer = new Epetra_Import(graph_->RowMap(),xthis->Map());  
      int ierr = xc->Import(*xthis,*importer,Insert);
      if (ierr)
      {
         cout << "**ERR**: ML_Epetra::ML_Nox_MatrixfreeLevel::ML_Nox_MatrixfreeLevel:\n"
              << "**ERR**: export from xthis to xc returned err=" << ierr <<"\n"
              << "**ERR**: file/line: " << __FILE__ << "/" << __LINE__ << "\n"; throw -1;
      }
      if (importer) delete importer; importer = 0;
   }
#endif   
   if (xthis) delete xthis; xthis = 0;

   // create the coloring of the graph
   if (ml_printlevel_>0 && comm_.MyPID()==0)
   {
      cout << "matrixfreeML (level " << level_ << "): Entering Recoloring on level " << level_ << "\n";
      fflush(stdout);
   }
   double t0 = GetClock();
   if (!same) // te graph has obviously changed, so we need to recolor
   {
      if (colorMap_) delete colorMap_; colorMap_ = 0;
      if (colorMapIndex_) delete colorMapIndex_; colorMapIndex_ = 0;
      if (colorcolumns_) delete colorcolumns_; colorcolumns_ = 0;

      colorMap_ = ML_NOX::ML_Nox_collapsedcoloring(graph_,bsize_,isDiagonalOnly_,ml_printlevel_);
      if (!colorMap_) colorMap_ = ML_NOX::ML_Nox_standardcoloring(graph_,isDiagonalOnly_);
      colorMapIndex_ = new EpetraExt::CrsGraph_MapColoringIndex(*colorMap_);
      colorcolumns_  = &(*colorMapIndex_)(*graph_);
    }
    else if (ml_printlevel_>0 && comm_.MyPID()==0)
      cout << "matrixfreeML (level " << level_ << "): Reusing existing Coloring on level " << level_ << "\n";
    double t1 = GetClock();
    if (ml_printlevel_>5)
    {
      cout << "matrixfreeML (level " << level_ << "): Proc " << comm_.MyPID() <<" (Re)Coloring time is " << (t1-t0) << " sec\n";
      fflush(stdout);
    }

#if 0
   // print the colorMap_
   if (comm_.MyPID()==0)
   cout << "colorMap_\n";
   cout << *colorMap_;
   for (int i=0; i<colorcolumns_->size(); i++)
   {
      if (comm_.MyPID()==0)
         cout << "the " << i << " th colorcolumn_ - vector\n";
      cout << colorcolumns_->operator[](i);
   }
#endif
   
   // construct the FiniteDifferenceColoring-Matrix
   if (ml_printlevel_>0 && comm_.MyPID()==0)
   {
      cout << "matrixfreeML (level " << level_ << "): Entering Construction FD-Operator on level " << level_ << "\n";
      fflush(stdout);
   }

   t0 = GetClock();
#if 1 // FD-operator with coloring (see the #if 1 in ml_nox_matrixfreelevel.H as well!)
   FD_ = new NOX::EpetraNew::FiniteDifferenceColoring(*coarseinterface_,
                                                      *xc,
                                                      *graph_,
                                                      *colorMap_,
                                                      *colorcolumns_,
                                                      true,
                                                      isDiagonalOnly_,
                                                      fd_beta_,fd_alpha_);
#else // FD-operator without coloring
   FD_ = new NOX::EpetraNew::FiniteDifference(*coarseinterface_,
                                              *xc,
                                              *graph_,
                                              fd_beta_,fd_alpha_);
#endif

   // set differencing method
   if (fd_centered_)
   {
      FD_->setDifferenceMethod(NOX::EpetraNew::FiniteDifferenceColoring::Centered);
   }

   bool err = FD_->computeJacobian(*xc); 
   if (err==false)
   {
     cout << "**ERR**: ML_Epetra::ML_Nox_MatrixfreeLevel::ML_Nox_MatrixfreeLevel:\n"
          << "**ERR**: NOX::Epetra::FiniteDifferenceColoring returned an error on level " 
          << level_ << "\n"
          << "**ERR**: file/line: " << __FILE__ << "/" << __LINE__ << "\n"; throw -1;
   }

   t1 = GetClock();
   if (ml_printlevel_>5)
      cout << "matrixfreeML (level " << level_ << "): Proc " << comm_.MyPID() 
           <<" Finite Differencing operator constr. in " << (t1-t0) << " sec\n";

   // get ref to computed Epetra_CrsMatrix  
   A_ = dynamic_cast<Epetra_CrsMatrix*>(&(FD_->getUnderlyingMatrix()));
   
   // print number of calls to the coarse interface
   if (ml_printlevel_>5 && comm_.MyPID()==0)
      cout << "matrixfreeML (level " << level_ << "): Calls to coarse-computeF in FD-Operator: "
           << coarseinterface_->numcallscomputeF() << "\n"; 
   
   // set counter for number of calls to the coarseinterface_->computeF back to zero
   coarseinterface_->resetnumcallscomputeF();   
   
   // tidy up 
   if (xc)       delete xc;       xc = 0;
   
   return true;
}
Beispiel #3
0
int 
main (int argc, char *argv[])
{
  using Teuchos::ArrayRCP;
  using Teuchos::ArrayView;
  using Teuchos::Comm;
  using Teuchos::CommandLineProcessor;
  using Teuchos::FancyOStream;
  using Teuchos::getFancyOStream;
  using Teuchos::OSTab;
  using Teuchos::ptr;
  using Teuchos::RCP;
  using Teuchos::rcp;
  using Teuchos::rcpFromRef;
  using std::cout;
  using std::endl;

  bool success = true; // May be changed by tests

  Teuchos::oblackholestream blackHole;
  //Teuchos::GlobalMPISession (&argc, &argv, &blackHole);
  MPI_Init (&argc, &argv);

  //
  // Construct communicators, and verify that we are on 4 processors.
  //

  // Construct a Teuchos Comm object.
  RCP<const Comm<int> > teuchosComm = Teuchos::DefaultComm<int>::getComm();
  const int numProcs = teuchosComm->getSize();
  const int pid = teuchosComm->getRank();
  RCP<FancyOStream> pOut = 
    getFancyOStream (rcpFromRef ((pid == 0) ? std::cout : blackHole));
  FancyOStream& out = *pOut;
  // Verify that we are on four processors (which manifests the bug).
  if (teuchosComm->getSize() != 4) {
    out << "This test must be run on four processors.  Exiting ..." << endl;
    return EXIT_FAILURE;
  }

  // We also need an Epetra Comm, so that we can compare Tpetra and
  // Epetra results.
  Epetra_MpiComm epetraComm (MPI_COMM_WORLD);

  //
  // Default values of command-line options.
  //
  bool verbose = false;
  bool printEpetra = false;
  bool printTpetra = false;
  CommandLineProcessor cmdp (false,true);
  //
  // Set command-line options.
  //
  cmdp.setOption ("verbose", "quiet", &verbose, "Print verbose output.");
  // Epetra and Tpetra output will ask the Maps and Import objects to
  // print themselves in distributed, maximally verbose fashion.  It's
  // best to turn on either Epetra or Tpetra, but not both.  Then you
  // can compare their output side by side.
  cmdp.setOption ("printEpetra", "dontPrintEpetra", &printEpetra, 
		  "Print Epetra output (in verbose mode only).");
  cmdp.setOption ("printTpetra", "dontPrintTpetra", &printTpetra, 
		  "Print Tpetra output (in verbose mode only).");
  // Parse command-line options.
  if (cmdp.parse (argc,argv) != CommandLineProcessor::PARSE_SUCCESSFUL) {
    out << "End Result: TEST FAILED" << endl;
    MPI_Finalize ();
    return EXIT_FAILURE;
  }

  if (verbose) {
    out << "Running test on " << numProcs << " process" 
	<< (numProcs != 1 ? "es" : "") << "." << endl;
  }

  // The maps for this problem are derived from a 3D structured mesh.
  // In this example, the dimensions are 4x4x2 and there are 2
  // processors assigned to the first dimension and 2 processors
  // assigned to the second dimension, with no parallel decomposition
  // along the third dimension.  The "owned" arrays represent the
  // one-to-one map, with each array representing a 2x2x2 slice.  If
  // DIMENSIONS == 2, then only the first 4 values will be used,
  // representing a 2x2(x1) slice.
  int owned0[8] = { 0, 1, 4, 5,16,17,20,21};
  int owned1[8] = { 2, 3, 6, 7,18,19,22,23};
  int owned2[8] = { 8, 9,12,13,24,25,28,29};
  int owned3[8] = {10,11,14,15,26,27,30,31};

  // The "overlap" arrays represent the map with communication
  // elements, with each array representing a 3x3x2 slice.  If
  // DIMENSIONS == 2, then only the first 9 values will be used,
  // representing a 3x3(x1) slice.
  int overlap0[18] = {0,1,2,4, 5, 6, 8, 9,10,16,17,18,20,21,22,24,25,26};
  int overlap1[18] = {1,2,3,5, 6, 7, 9,10,11,17,18,19,21,22,23,25,26,27};
  int overlap2[18] = {4,5,6,8, 9,10,12,13,14,20,21,22,24,25,26,28,29,30};
  int overlap3[18] = {5,6,7,9,10,11,13,14,15,21,22,23,25,26,27,29,30,31};

  // Construct the owned and overlap maps for both Epetra and Tpetra.
  int* owned;
  int* overlap;
  if (pid == 0) {
    owned   = owned0;
    overlap = overlap0;
  }
  else if (pid == 1) {
    owned   = owned1;
    overlap = overlap1;
  }
  else if (pid == 2) {
    owned   = owned2;
    overlap = overlap2;
  }
  else {
    owned   = owned3;
    overlap = overlap3;
  }

#if DIMENSIONS == 2
  int ownedSize   = 4;
  int overlapSize = 9;
#elif DIMENSIONS == 3
  int ownedSize   =  8;
  int overlapSize = 18;
#endif

  // Create the two Epetra Maps.  Source for the Import is the owned
  // map; target for the Import is the overlap map.
  Epetra_Map epetraOwnedMap (  -1, ownedSize,   owned,   0, epetraComm);
  Epetra_Map epetraOverlapMap (-1, overlapSize, overlap, 0, epetraComm);

  if (verbose && printEpetra) {
    // Have the Epetra_Map objects describe themselves.
    //
    // Epetra_BlockMap::Print() takes an std::ostream&, and expects
    // all MPI processes to be able to write to it.  (The method
    // handles its own synchronization.)
    out << "Epetra owned map:" << endl;
    epetraOwnedMap.Print (std::cout);
    out << "Epetra overlap map:" << endl;
    epetraOverlapMap.Print (std::cout);
  }

  // Create the two Tpetra Maps.  The "invalid" global element count
  // input tells Tpetra::Map to compute the global number of elements
  // itself.
  const int invalid = Teuchos::OrdinalTraits<int>::invalid();
  RCP<Tpetra::Map<int> > tpetraOwnedMap = 
    rcp (new Tpetra::Map<int> (invalid, ArrayView<int> (owned, ownedSize), 
			       0, teuchosComm));
  tpetraOwnedMap->setObjectLabel ("Owned Map");
  RCP<Tpetra::Map<int> > tpetraOverlapMap =
    rcp (new Tpetra::Map<int> (invalid, ArrayView<int> (overlap, overlapSize),
			       0, teuchosComm));
  tpetraOverlapMap->setObjectLabel ("Overlap Map");

  // In verbose mode, have the Tpetra::Map objects describe themselves.
  if (verbose && printTpetra) {
    Teuchos::EVerbosityLevel verb = Teuchos::VERB_EXTREME;

    // Tpetra::Map::describe() takes a FancyOStream, but expects all
    // MPI processes to be able to write to it.  (The method handles
    // its own synchronization.)
    RCP<FancyOStream> globalOut = getFancyOStream (rcpFromRef (std::cout));
    out << "Tpetra owned map:" << endl;
    {
      OSTab tab (globalOut);
      tpetraOwnedMap->describe (*globalOut, verb);
    }
    out << "Tpetra overlap map:" << endl;
    {
      OSTab tab (globalOut);
      tpetraOverlapMap->describe (*globalOut, verb);
    }
  }

  // Use the owned and overlap maps to construct an importer for both
  // Epetra and Tpetra.
  Epetra_Import       epetraImporter (epetraOverlapMap, epetraOwnedMap  );
  Tpetra::Import<int> tpetraImporter (tpetraOwnedMap  , tpetraOverlapMap);

  // In verbose mode, have the Epetra_Import object describe itself.
  if (verbose && printEpetra) {
    out << "Epetra importer:" << endl;
    // The importer's Print() method takes an std::ostream& and plans
    // to write to it on all MPI processes (handling synchronization
    // itself).
    epetraImporter.Print (std::cout);
    out << endl;
  }

  // In verbose mode, have the Tpetra::Import object describe itself.
  if (verbose && printTpetra) {
    out << "Tpetra importer:" << endl;
    // The importer doesn't implement Teuchos::Describable.  It wants
    // std::cout and plans to write to it on all MPI processes (with
    // its own synchronization).
    tpetraImporter.print (std::cout);
    out << endl;
  }

  // Construct owned and overlap vectors for both Epetra and Tpetra.
  Epetra_Vector epetraOwnedVector   (epetraOwnedMap  );
  Epetra_Vector epetraOverlapVector (epetraOverlapMap);
  Tpetra::Vector<double,int> tpetraOwnedVector   (tpetraOwnedMap  );
  Tpetra::Vector<double,int> tpetraOverlapVector (tpetraOverlapMap);

  // The test is as follows: initialize the owned and overlap vectors
  // with global IDs in the owned regions.  Initialize the overlap
  // vectors to equal -1 in the overlap regions.  Then perform a
  // communication from the owned vectors to the overlap vectors.  The
  // resulting overlap vectors should have global IDs everywhere and
  // all of the -1 values should be overwritten.

  // Initialize.  We cannot assign directly to the Tpetra Vectors;
  // instead, we extract nonconst views and assign to those.  The
  // results aren't guaranteed to be committed to the vector unless
  // the views are released (by assigning Teuchos::null to them).
  epetraOverlapVector.PutScalar(-1);
  tpetraOverlapVector.putScalar(-1);
  ArrayRCP<double> tpetraOwnedArray   = tpetraOwnedVector.getDataNonConst(0);
  ArrayRCP<double> tpetraOverlapArray = tpetraOverlapVector.getDataNonConst(0);
  for (int owned_lid = 0; 
       owned_lid < tpetraOwnedMap->getNodeElementList().size(); 
       ++owned_lid) {
    int gid         = tpetraOwnedMap->getGlobalElement(owned_lid);
    int overlap_lid = tpetraOverlapMap->getLocalElement(gid);
    epetraOwnedVector[owned_lid]     = gid;
    epetraOverlapVector[overlap_lid] = gid;
    tpetraOwnedArray[owned_lid]      = gid;
    tpetraOverlapArray[overlap_lid]  = gid;
  }
  // Make sure that the changes to the Tpetra Vector were committed,
  // by releasing the nonconst views.
  tpetraOwnedArray = Teuchos::null;
  tpetraOverlapArray = Teuchos::null;

  // Test the Epetra and Tpetra Import.
  if (verbose) {
    out << "Testing Import from owned Map to overlap Map:" << endl << endl;
  }
  epetraOverlapVector.Import(  epetraOwnedVector, epetraImporter, Insert);
  tpetraOverlapVector.doImport(tpetraOwnedVector, tpetraImporter, 
			       Tpetra::INSERT);
  // Check the Import results.
  success = countFailures (teuchosComm, epetraOwnedMap, epetraOwnedVector, 
			   epetraOverlapMap, epetraOverlapVector, 
			   tpetraOwnedMap, tpetraOwnedVector, 
			   tpetraOverlapMap, tpetraOverlapVector, verbose);

  const bool testOtherDirections = false;
  if (testOtherDirections) {
    //
    // Reinitialize the Tpetra vectors and test whether Export works.
    //
    tpetraOverlapVector.putScalar(-1);
    tpetraOwnedArray   = tpetraOwnedVector.getDataNonConst(0);
    tpetraOverlapArray = tpetraOverlapVector.getDataNonConst(0);
    for (int owned_lid = 0; 
	 owned_lid < tpetraOwnedMap->getNodeElementList().size(); 
	 ++owned_lid) 
      {
	int gid         = tpetraOwnedMap->getGlobalElement(owned_lid);
	int overlap_lid = tpetraOverlapMap->getLocalElement(gid);
	tpetraOwnedArray[owned_lid]      = gid;
	tpetraOverlapArray[overlap_lid]  = gid;
      }
    // Make sure that the changes to the Tpetra Vector were committed,
    // by releasing the nonconst views.
    tpetraOwnedArray = Teuchos::null;
    tpetraOverlapArray = Teuchos::null;

    // Make a Tpetra Export object, and test the export.
    Tpetra::Export<int> tpetraExporter1 (tpetraOwnedMap, tpetraOverlapMap);
    if (verbose) {
      out << "Testing Export from owned Map to overlap Map:" << endl << endl;
    }
    tpetraOverlapVector.doExport (tpetraOwnedVector, tpetraExporter1, 
				  Tpetra::INSERT);

    // Check the Export results.
    success = countFailures (teuchosComm, epetraOwnedMap, epetraOwnedVector, 
			     epetraOverlapMap, epetraOverlapVector, 
			     tpetraOwnedMap, tpetraOwnedVector, 
			     tpetraOverlapMap, tpetraOverlapVector, verbose);
    //
    // Reinitialize the Tpetra vectors and see what Import in the
    // other direction does.
    //
    tpetraOverlapVector.putScalar(-1);
    tpetraOwnedArray   = tpetraOwnedVector.getDataNonConst(0);
    tpetraOverlapArray = tpetraOverlapVector.getDataNonConst(0);
    for (int owned_lid = 0; 
	 owned_lid < tpetraOwnedMap->getNodeElementList().size(); 
	 ++owned_lid) 
      {
	int gid         = tpetraOwnedMap->getGlobalElement(owned_lid);
	int overlap_lid = tpetraOverlapMap->getLocalElement(gid);
	tpetraOwnedArray[owned_lid]      = gid;
	tpetraOverlapArray[overlap_lid]  = gid;
      }
    // Make sure that the changes to the Tpetra Vector were committed,
    // by releasing the nonconst views.
    tpetraOwnedArray = Teuchos::null;
    tpetraOverlapArray = Teuchos::null;

    if (verbose) {
      out << "Testing Import from overlap Map to owned Map:" << endl << endl;
    }
    Tpetra::Import<int> tpetraImporter2 (tpetraOverlapMap, tpetraOwnedMap);
    tpetraOwnedVector.doImport (tpetraOverlapVector, tpetraImporter2, 
				Tpetra::INSERT);
    // Check the Import results.
    success = countFailures (teuchosComm, epetraOwnedMap, epetraOwnedVector, 
			     epetraOverlapMap, epetraOverlapVector, 
			     tpetraOwnedMap, tpetraOwnedVector, 
			     tpetraOverlapMap, tpetraOverlapVector, verbose);
  } // if testOtherDirections

  out << "End Result: TEST " << (success ? "PASSED" : "FAILED") << endl;
  MPI_Finalize ();
  return success ? EXIT_SUCCESS : EXIT_FAILURE;
}
bool FiniteDifferenceColoringWithUpdate::differenceProbe(const Epetra_Vector& x, Epetra_CrsMatrix& jac,const Epetra_MapColoring& colors){

  // Allocate space for perturbation, get column version of x for scaling
  Epetra_Vector xp(x);
  Epetra_Vector *xcol;
  int N=jac.NumMyRows();

  if(jac.ColMap().SameAs(x.Map()))
     xcol=const_cast<Epetra_Vector*>(&x);
  else{
    xcol=new Epetra_Vector(jac.ColMap(),true);//zeros out by default
    xcol->Import(x,*jac.Importer(),InsertAdd);
  }

  // Counters for probing diagnostics
  double tmp,probing_error_lower_bound=0.0,jc_norm=0.0;

  // Grab coloring info (being very careful to ignore color 0)
  int Ncolors=colors.MaxNumColors()+1;
  int num_c0_global,num_c0_local=colors.NumElementsWithColor(0);
  colors.Comm().MaxAll(&num_c0_local,&num_c0_global,1);
  if(num_c0_global>0) Ncolors--;

  if(Ncolors==0) return false;

  // Pointers for Matrix Info
  int entries, *indices;
  double *values;

  // NTS: Fix me
  if ( diffType == Centered ) exit(1);

  double scaleFactor = 1.0;
  if ( diffType == Backward )
    scaleFactor = -1.0;

  // Compute RHS at initial solution
  computeF(x,fo,NOX::Epetra::Interface::Required::FD_Res);

  /* Probing, vector by vector since computeF does not have a MultiVector interface */
  // Assume that anything with Color 0 gets ignored.
  for(int j=1;j<Ncolors;j++){
    xp=x;
    for(int i=0;i<N;i++){
      if(colors[i]==j)
    xp[i] += scaleFactor*(alpha*abs(x[i])+beta);
    }

    computeF(xp, fp, NOX::Epetra::Interface::Required::FD_Res);

    // Do the subtraction to estimate the Jacobian (w/o including step length)
    Jc.Update(1.0, fp, -1.0, fo, 0.0);

    // Relative error in probing
     if(use_probing_diags){
       Jc.Norm2(&tmp);
       jc_norm+=tmp*tmp;
     }

    for(int i=0;i<N;i++){
      // Skip for uncolored row/columns, else update entries
      if(colors[i]==0) continue;

      jac.ExtractMyRowView(i,entries,values,indices);
      for(int k=0;k<jac.NumMyEntries(i);k++){
    if(colors[indices[k]]==j){
      values[k]=Jc[i] / (scaleFactor*(alpha*abs((*xcol)[indices[k]])+beta));
      // If probing diagnostics are on, zero out the entries as they are used
      if(use_probing_diags) Jc[i]=0.0;
      break;// Only one value per row...
    }
      }
    }
    if(use_probing_diags){
      Jc.Norm2(&tmp);
      probing_error_lower_bound+=tmp*tmp;
    }
  }

  // If diagnostics are requested, output Frobenius norm lower bound
  if(use_probing_diags && !x.Comm().MyPID()) printf("Probing Error Lower Bound (Frobenius) abs = %6.4e rel = %6.4e\n",sqrt(probing_error_lower_bound),sqrt(probing_error_lower_bound)/sqrt(jc_norm));

  // Cleanup
  if(!jac.ColMap().SameAs(x.Map()))
    delete xcol;

  return true;
}
Beispiel #5
0
// ====================================================================== 
int Amesos_Lapack::GEEV(Epetra_Vector& Er, Epetra_Vector& Ei)
{
  if (IsSymbolicFactorizationOK_ == false)
    AMESOS_CHK_ERR(SymbolicFactorization());

  if (MyPID_ == 0)
    AMESOS_CHK_ERR(DenseMatrix_.Shape(static_cast<int>(NumGlobalRows64()),static_cast<int>(NumGlobalRows64())));

  AMESOS_CHK_ERR(DistributedToSerial());
  AMESOS_CHK_ERR(SerialToDense());

  Teuchos::RCP<Epetra_Vector> LocalEr;
  Teuchos::RCP<Epetra_Vector> LocalEi;

  if (NumProcs_ == 1)
  {
    LocalEr = Teuchos::rcp(&Er, false);
    LocalEi = Teuchos::rcp(&Ei, false);
  }
  else
  {
    LocalEr = Teuchos::rcp(new Epetra_Vector(*SerialMap_));
    LocalEi = Teuchos::rcp(new Epetra_Vector(*SerialMap_));
  }

  if (MyPID_ == 0) 
  {
    int n = static_cast<int>(NumGlobalRows64());
    char jobvl = 'N'; /* V/N to calculate/not calculate left eigenvectors
                         of matrix H.*/
    char jobvr = 'N'; /* As above, but for right eigenvectors. */
    int info = 0;
    int ldvl = n;
    int ldvr = n;

    Er.PutScalar(0.0);
    Ei.PutScalar(0.0);

    Epetra_LAPACK LAPACK;

    std::vector<double> work(1);
    int lwork = -1;

    LAPACK.GEEV(jobvl, jobvr, n, DenseMatrix_.A(), n, 
                LocalEr->Values(), LocalEi->Values(), NULL,
                ldvl, NULL, 
                ldvr, &work[0], 
                lwork, &info);

    lwork = (int)work[0];
    work.resize(lwork);
    LAPACK.GEEV(jobvl, jobvr, n, DenseMatrix_.A(), n, 
                LocalEr->Values(), LocalEi->Values(), NULL,
                ldvl, NULL, 
                ldvr, &work[0], 
                lwork, &info);

    if (info)
      AMESOS_CHK_ERR(info);
  }

  if (NumProcs_ != 1)
  {
    // I am not really sure that exporting the results make sense... 
    // It is just to be coherent with the other parts of the code.
    Er.Import(*LocalEr, Epetra_Import(Er.Map(), SerialMap()), Insert);
    Ei.Import(*LocalEi, Epetra_Import(Ei.Map(), SerialMap()), Insert);
  }

  return(0);
}