>
> Yes I apologize for not including code. A program that runs is included 
> below. Principally, I'm not sure whether deal.ii means to create an 
> Epetra_FEVector with LinearMap() equal to true or not. In the source code, 
> the word "linear" is used when constructing Epetra_FEVector object.


The reason I need a Epetra_FEVector for which LinearMap() returns true is 
that Hypre seems to need this. I am using Ifpack_Hypre to access the 
BoomerAMG preconditioner. The problem code from Ifpack_Hypre is included 
below. 

int Ifpack_Hypre::ApplyInverse(const Epetra_MultiVector& X, 
Epetra_MultiVector& Y) const{
  if(IsComputed() == false){
    IFPACK_CHK_ERR(-1);
  }
  // These are hypre requirements
  // hypre needs A, X, and Y to have the same contiguous distribution
  // NOTE: Maps are only considered to be contiguous if they were generated 
using a
  // particular constructor.  Otherwise, LinearMap() will not detect 
whether they are
  // actually contiguous.
  if(!X.Map().LinearMap() || !Y.Map().LinearMap()) {
    std::cerr << "ERROR: X and Y must have contiguous maps.\n";
    IFPACK_CHK_ERR(-1);
  }
  if(!X.Map().PointSameAs(*MySimpleMap_) ||
     !Y.Map().PointSameAs(*MySimpleMap_)) {
    std::cerr << "ERROR: X, Y, and A must have the same distribution.\n";
    IFPACK_CHK_ERR(-1);
  }


 



#include <deal.II/grid/tria.h>
#include <deal.II/dofs/dof_handler.h>
#include <deal.II/grid/grid_generator.h>
#include <deal.II/grid/tria_accessor.h>
#include <deal.II/grid/tria_iterator.h>
#include <deal.II/dofs/dof_accessor.h>
#include <deal.II/fe/fe_q.h>
#include <deal.II/dofs/dof_tools.h>
#include <deal.II/fe/fe_values.h>
#include <deal.II/base/quadrature_lib.h>
#include <deal.II/base/function.h>
#include <deal.II/numerics/vector_tools.h>
#include <deal.II/numerics/matrix_tools.h>
#include <deal.II/lac/vector.h>
#include <deal.II/lac/full_matrix.h>
#include <deal.II/lac/sparse_matrix.h>
#include <deal.II/lac/dynamic_sparsity_pattern.h>
#include <deal.II/lac/precondition.h>
#include <deal.II/numerics/data_out.h>
#include <fstream>
#include <iostream>
#include <vector>
//
///////////////////////////////////////////////
//
#include <deal.II/base/mpi.h>
#include <deal.II/lac/generic_linear_algebra.h>
#include <deal.II/base/utilities.h>
#include <deal.II/base/conditional_ostream.h>
#include <deal.II/base/index_set.h>
#include <deal.II/lac/sparsity_tools.h>
//
#include <deal.II/distributed/tria.h>
#include <deal.II/distributed/grid_refinement.h>
//
#include <deal.II/base/timer.h>
//
#include <deal.II/grid/grid_tools.h>
#include <deal.II/dofs/dof_renumbering.h>
#include <deal.II/lac/trilinos_solver.h>
//
#include <deal.II/lac/trilinos_vector.h>
#include <deal.II/lac/trilinos_vector_base.h>
#include <deal.II/lac/trilinos_sparse_matrix.h>

#include<Ifpack_Hypre.h>
#include <Teuchos_ParameterList.hpp>

#include "HYPRE_IJ_mv.h"
#include "HYPRE_parcsr_ls.h"
#include "krylov.h"
#include "_hypre_parcsr_mv.h"
#include "_hypre_IJ_mv.h"
#include "HYPRE_parcsr_mv.h"
#include "HYPRE.h"

#include "Epetra_MultiVector.h"


namespace LA =  dealii::LinearAlgebraTrilinos;

using namespace dealii;


class AmgElliptic
{
public:
  AmgElliptic ();

  void run ();
  void set_uniform_grid (unsigned int n_cells, double length);


  void setup_system ();
  void assemble_system ();
  void solve ();
  void output_results () const;

  Triangulation<2>     triangulation;
  FE_Q<2>              fe;
  DoFHandler<2>        dof_handler;


  SparsityPattern      sparsity_pattern;

  MPI_Comm     mpi_communicator;

  const unsigned int n_mpi_processes;
  const unsigned int this_mpi_process;

  ConditionalOStream        pcout;
 // TimerOutput               computing_timer;

  LA::MPI::SparseMatrix     system_matrix;
  LA::MPI::Vector           system_rhs;
  LA::MPI::Vector     solution;

};


AmgElliptic::AmgElliptic ()
  :
  mpi_communicator (MPI_COMM_WORLD),
  n_mpi_processes (Utilities::MPI::n_mpi_processes(mpi_communicator)),
  this_mpi_process (Utilities::MPI::this_mpi_process(mpi_communicator)),
  pcout (std::cout,
         (this_mpi_process == 0)),
  fe (1),
  dof_handler (triangulation)
{ set_uniform_grid(10,1.0); }

void AmgElliptic :: set_uniform_grid(unsigned int n_cells, double length){
dealii::Point<2> lower_left_pt(0.0,0.0),upper_right_pt(length,length);
std::vector<unsigned int> repetitions(2);
repetitions[0] = n_cells;
repetitions[1] = n_cells;

dealii::GridGenerator::subdivided_hyper_rectangle(triangulation , 
repetitions , lower_left_pt , upper_right_pt);

}

void AmgElliptic::setup_system ()
{

  GridTools::partition_triangulation (n_mpi_processes, triangulation);

  dof_handler.distribute_dofs (fe);
  DoFRenumbering::subdomain_wise (dof_handler);

  
  
    const std::vector<IndexSet> locally_owned_dofs_per_proc = 
DoFTools::locally_owned_dofs_per_subdomain(dof_handler);
    const IndexSet locally_owned_dofs = 
locally_owned_dofs_per_proc[this_mpi_process];

    DynamicSparsityPattern dsp(dof_handler.n_dofs());
    DoFTools::make_sparsity_pattern (dof_handler, dsp);
    sparsity_pattern.copy_from(dsp);


    system_matrix.reinit (locally_owned_dofs,
                          locally_owned_dofs,
                          dsp,
                          mpi_communicator);
   //
   // Here I check if is_ascending_and_one_to_one(mpi_communicator) will be 
true and when I run
   // the result printed is 1
   //
   std::cout<<"is_ascending_and_one_to_one = "<<
    
locally_owned_dofs.is_ascending_and_one_to_one(mpi_communicator)<<"\n\n";
   //
   //
   //
    solution.reinit(locally_owned_dofs, mpi_communicator);
    system_rhs.reinit(locally_owned_dofs, mpi_communicator);
   //
   // Here I check if LinearMap_ for these vectors is true. When I run, 
this prints 0
   //
    std::cout<<"solution.trilinos_vector().Map().LinearMap() = "<<
    solution.trilinos_vector().Map().LinearMap()<<"\n\n";

    //
    // Here I try something else and make a new vector with a copy 
constructor to see
    // if this will have LinearMap_ set to 1
    //
Epetra_MultiVector e_soln_vector(solution.trilinos_vector());
Epetra_MultiVector e_rhs_vector(system_rhs.trilinos_vector());
//
    std::cout<<"e_soln_vector.trilinos_vector().Map().LinearMap() = "<<
    e_soln_vector.Map().LinearMap()<<"\n";

}

void AmgElliptic::run ()
{
  setup_system ();
}

int main (int argc, char **argv)
{
  Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1);   
    
  deallog.depth_console (2);

  AmgElliptic test_simulation;
  
  test_simulation.run();

  return 0;
}


 

-- 
The deal.II project is located at http://www.dealii.org/
For mailing list/forum options, see 
https://groups.google.com/d/forum/dealii?hl=en
--- 
You received this message because you are subscribed to the Google Groups 
"deal.II User Group" group.
To unsubscribe from this group and stop receiving emails from it, send an email 
to [email protected].
For more options, visit https://groups.google.com/d/optout.

Reply via email to