/*! Create an Epetra_MpiComm from Fortran */
CT_Epetra_MpiComm_ID_t Epetra_MpiComm_Fortran_Create ( int fcomm )
{
    MPI_Fint mfcomm = (MPI_Fint) fcomm;
    MPI_Comm ccomm = MPI_Comm_f2c(mfcomm);

    /* duplicate the communicator so that we won't get any cross-talk
     * from the application */
    MPI_Comm dupcomm;
    int ret = MPI_Comm_dup(ccomm, &dupcomm);
    if (ret != MPI_SUCCESS)
        throw CTrilinos::CTrilinosMiscException("Error on MPI_Comm_dup");

    return Epetra_MpiComm_Create(dupcomm);
}
int main(int argc, char *argv[])
{
  int ierr = 0, i, forierr = 0;
#ifdef HAVE_MPI

  CT_Epetra_MpiComm_ID_Flex_t Comm;

  // Initialize MPI

  MPI_Init(&argc,&argv);
  int rank; // My process ID

  MPI_Comm_rank(MPI_COMM_WORLD, &rank);
  Comm.Epetra_MpiComm = Epetra_MpiComm_Create( MPI_COMM_WORLD );

#else

  int rank = 0;
  CT_Epetra_SerialComm_ID_Flex_t Comm;
  Comm.Epetra_SerialComm = Epetra_SerialComm_Create();

#endif

  bool verbose = false;

  // Check if we should print results to standard out
  if (argc>1) if (argv[1][0]=='-' && argv[1][1]=='v') verbose = true;

  int verbose_int = verbose ? 1 : 0;
  Epetra_Comm_Broadcast_Int(Comm.Epetra_Comm, &verbose_int, 1, 0);
  verbose = verbose_int==1 ? true : false;


  //  char tmp;
  //  if (rank==0) cout << "Press any key to continue..."<< endl;
  //  if (rank==0) cin >> tmp;
  //  Epetra_Comm_Barrier(Comm.Epetra_Comm);

//  Epetra_Comm_SetTracebackMode(Comm.Epetra_Comm, 0); // This should shut down any error traceback reporting
  int MyPID = Epetra_Comm_MyPID(Comm.Epetra_Comm);
  int NumProc = Epetra_Comm_NumProc(Comm.Epetra_Comm);

  if(verbose && MyPID==0)
    cout << Epetra_Version() << endl << endl;

  if (verbose) cout << "Processor "<<MyPID<<" of "<< NumProc
		    << " is alive."<<endl;

  // Redefine verbose to only print on PE 0
  if(verbose && rank!=0) 
		verbose = false;

  int NumMyEquations = 10000;
  int NumGlobalEquations = (NumMyEquations * NumProc) + EPETRA_MIN(NumProc,3);
  if(MyPID < 3) 
    NumMyEquations++;

  // Construct a Map that puts approximately the same Number of equations on each processor

  CT_Epetra_Map_ID_Flex_t Map;
  Map.Epetra_Map = Epetra_Map_Create_Linear(NumGlobalEquations, NumMyEquations, 0, Comm.Epetra_Comm);
  
  // Get update list and number of local equations from newly created Map
  vector<int> MyGlobalElements(Epetra_BlockMap_NumMyElements(Map.Epetra_BlockMap));
  Epetra_BlockMap_MyGlobalElements_Fill(Map.Epetra_BlockMap, &MyGlobalElements[0]);

  // Create an integer vector NumNz that is used to build the Petra Matrix.
  // NumNz[i] is the Number of OFF-DIAGONAL term for the ith global equation on this processor

  vector<int> NumNz(NumMyEquations);

  // We are building a tridiagonal matrix where each row has (-1 2 -1)
  // So we need 2 off-diagonal terms (except for the first and last equation)

  for(i = 0; i < NumMyEquations; i++)
    if((MyGlobalElements[i] == 0) || (MyGlobalElements[i] == NumGlobalEquations - 1))
      NumNz[i] = 1;
    else
      NumNz[i] = 2;

  // Create a Epetra_Matrix

  CT_Epetra_CrsMatrix_ID_Flex_t A;
  A.Epetra_CrsMatrix = Epetra_CrsMatrix_Create_VarPerRow(
      CT_Epetra_DataAccess_E_Copy, Map.Epetra_Map, &NumNz[0], FALSE);
  EPETRA_TEST_ERR(Epetra_CrsMatrix_IndicesAreGlobal(A.Epetra_CrsMatrix),ierr);
  EPETRA_TEST_ERR(Epetra_CrsMatrix_IndicesAreLocal(A.Epetra_CrsMatrix),ierr);
  
  // Add  rows one-at-a-time
  // Need some vectors to help
  // Off diagonal Values will always be -1


  vector<double> Values(2);
  Values[0] = -1.0; 
	Values[1] = -1.0;
	vector<int> Indices(2);
  double two = 2.0;
  int NumEntries;

  forierr = 0;
  for(i = 0; i < NumMyEquations; i++) {
    if(MyGlobalElements[i] == 0) {
			Indices[0] = 1;
			NumEntries = 1;
		}
    else if (MyGlobalElements[i] == NumGlobalEquations-1) {
			Indices[0] = NumGlobalEquations-2;
			NumEntries = 1;
		}
    else {
			Indices[0] = MyGlobalElements[i]-1;
			Indices[1] = MyGlobalElements[i]+1;
			NumEntries = 2;
		}
		forierr += !(Epetra_CrsMatrix_InsertGlobalValues(A.Epetra_CrsMatrix,
                             MyGlobalElements[i], NumEntries, &Values[0], &Indices[0])==0);
		forierr += !(Epetra_CrsMatrix_InsertGlobalValues(A.Epetra_CrsMatrix,
                             MyGlobalElements[i], 1, &two, &MyGlobalElements[i])>0); // Put in the diagonal entry
  }
  EPETRA_TEST_ERR(forierr,ierr);

  // Finish up
  Epetra_CrsMatrix_FillComplete(A.Epetra_CrsMatrix, TRUE);
  Epetra_CrsMatrix_OptimizeStorage(A.Epetra_CrsMatrix);

  CT_Epetra_JadMatrix_ID_Flex_t JadA, JadA1, JadA2;
  JadA.Epetra_JadMatrix = Epetra_JadMatrix_Create(A.Epetra_RowMatrix);
  JadA1.Epetra_JadMatrix = Epetra_JadMatrix_Create(A.Epetra_RowMatrix);
  JadA2.Epetra_JadMatrix = Epetra_JadMatrix_Create(A.Epetra_RowMatrix);

  // Create vectors for Power method

  CT_Epetra_Vector_ID_Flex_t q, z, resid;
  q.Epetra_Vector = Epetra_Vector_Create(Map.Epetra_BlockMap, TRUE);
  z.Epetra_Vector = Epetra_Vector_Create(Map.Epetra_BlockMap, TRUE);
  Epetra_MultiVector_Random(z.Epetra_MultiVector);
  resid.Epetra_Vector = Epetra_Vector_Create(Map.Epetra_BlockMap, TRUE);

  CT_Epetra_Flops_ID_t flopcounter = Epetra_Flops_Create();
  Epetra_CompObject_SetFlopCounter(A.Epetra_CompObject, flopcounter);

  Epetra_CompObject_SetFlopCounter_Matching(q.Epetra_CompObject, A.Epetra_CompObject);
  Epetra_CompObject_SetFlopCounter_Matching(z.Epetra_CompObject, A.Epetra_CompObject);
  Epetra_CompObject_SetFlopCounter_Matching(resid.Epetra_CompObject, A.Epetra_CompObject);

  Epetra_CompObject_SetFlopCounter_Matching(JadA.Epetra_CompObject, A.Epetra_CompObject);
  Epetra_CompObject_SetFlopCounter_Matching(JadA1.Epetra_CompObject, A.Epetra_CompObject);
  Epetra_CompObject_SetFlopCounter_Matching(JadA2.Epetra_CompObject, A.Epetra_CompObject);
  

  if (verbose) cout << "=======================================" << endl
		    << "Testing Jad using CrsMatrix as input..." << endl
		    << "=======================================" << endl;

  Epetra_CompObject_ResetFlops(A.Epetra_CompObject);
  powerMethodTests(A.Epetra_RowMatrix, JadA.Epetra_RowMatrix, Map, q, z, resid, verbose);

  // Increase diagonal dominance

  if (verbose) cout << "\n\nIncreasing the magnitude of first diagonal term and solving again\n\n"
		    << endl;

  
  if (Epetra_CrsMatrix_MyGlobalRow(A.Epetra_CrsMatrix, 0)) {
    int numvals = Epetra_CrsMatrix_NumGlobalEntries(A.Epetra_CrsMatrix, 0);
    vector<double> Rowvals(numvals);
    vector<int> Rowinds(numvals);
    Epetra_CrsMatrix_ExtractGlobalRowCopy_WithIndices(A.Epetra_CrsMatrix, 0,
        numvals, &numvals, &Rowvals[0], &Rowinds[0]); // Get A[0,0]

    for (i=0; i<numvals; i++) if (Rowinds[i] == 0) Rowvals[i] *= 10.0;
    
    Epetra_CrsMatrix_ReplaceGlobalValues(A.Epetra_CrsMatrix, 0, numvals, &Rowvals[0], &Rowinds[0]);
  }
  Epetra_JadMatrix_UpdateValues(JadA.Epetra_JadMatrix, A.Epetra_RowMatrix, FALSE);
  Epetra_CompObject_ResetFlops(A.Epetra_CompObject);
  powerMethodTests(A.Epetra_RowMatrix, JadA.Epetra_RowMatrix, Map, q, z, resid, verbose);

  if (verbose) cout << "================================================================" << endl
		          << "Testing Jad using Jad matrix as input matrix for construction..." << endl
		          << "================================================================" << endl;
  Epetra_CompObject_ResetFlops(JadA1.Epetra_CompObject);
  powerMethodTests(JadA1.Epetra_RowMatrix, JadA2.Epetra_RowMatrix, Map, q, z, resid, verbose);

#ifdef HAVE_MPI
  MPI_Finalize() ;
#endif

return ierr ;
}
int main(int argc, char *argv[])
{
  CT_Epetra_Map_ID_Flex_t Map;
  CT_Epetra_CrsMatrix_ID_Flex_t A;
  CT_Epetra_Vector_ID_Flex_t x, b;
  CT_Epetra_LinearProblem_ID_t problem;
  CT_AztecOO_ID_t solver;

  int NumMyElements, NumGlobalElements, i, GlobalRow, RowLess1, RowPlus1;
  double negOne, posTwo;

#ifdef HAVE_MPI
  /* Initialize MPI */
  CT_Epetra_MpiComm_ID_Flex_t Comm;
  MPI_Init(&argc, &argv);
  Comm.Epetra_MpiComm = Epetra_MpiComm_Create(MPI_COMM_WORLD);
#else
  CT_Epetra_SerialComm_ID_Flex_t Comm;
  Comm.Epetra_SerialComm = Epetra_SerialComm_Create();
#endif

  NumMyElements = 1000;

  /* Construct a Map that puts same number of equations on each processor */
  Map.Epetra_Map = Epetra_Map_Create_Linear(-1, NumMyElements, 0, Comm.Epetra_Comm);

  NumGlobalElements = Epetra_BlockMap_NumGlobalElements(Map.Epetra_BlockMap);

  /* Create a Epetra_Matrix */
  A.Epetra_CrsMatrix = Epetra_CrsMatrix_Create(CT_Epetra_DataAccess_E_Copy, Map.Epetra_Map, 3, FALSE);
  
  /* Add rows one-at-a-time */
  negOne = -1.0;
  posTwo = 2.0;
  for (i=0; i<NumMyElements; i++) {
    GlobalRow = Epetra_CrsMatrix_GRID(A.Epetra_CrsMatrix, i);
    RowLess1 = GlobalRow - 1;
    RowPlus1 = GlobalRow + 1;

    if (RowLess1 != -1)
      Epetra_CrsMatrix_InsertGlobalValues(A.Epetra_CrsMatrix, GlobalRow, 1, &negOne, &RowLess1);
    if (RowPlus1 != NumGlobalElements)
      Epetra_CrsMatrix_InsertGlobalValues(A.Epetra_CrsMatrix, GlobalRow, 1, &negOne, &RowPlus1);
    Epetra_CrsMatrix_InsertGlobalValues(A.Epetra_CrsMatrix, GlobalRow, 1, &posTwo, &GlobalRow);
  }
  
  /* Finish up */
  Epetra_CrsMatrix_FillComplete(A.Epetra_CrsMatrix, TRUE);

  /* Create x and b vectors */
  x.Epetra_Vector = Epetra_Vector_Create(Map.Epetra_BlockMap, TRUE);
  b.Epetra_Vector = Epetra_Vector_Create(Map.Epetra_BlockMap, TRUE);
  Epetra_MultiVector_Random(b.Epetra_MultiVector);

  /* Create Linear Problem */
  problem = Epetra_LinearProblem_Create_FromMatrix(A.Epetra_RowMatrix, x.Epetra_MultiVector, b.Epetra_MultiVector);

  /* Create AztecOO instance */
  solver = AztecOO_Create_FromLinearProblem(problem);

  AztecOO_SetAztecOption(solver, AZ_precond, AZ_Jacobi);
  AztecOO_Iterate_Current(solver, 1000, 1.0E-8);

  printf("Solver performed %d iterations.\n", AztecOO_NumIters(solver));
  printf("Norm of true residual = %g\n", AztecOO_TrueResidual(solver));

#ifdef HAVE_MPI
  MPI_Finalize() ;
#endif
  return 0;
}