]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Deprecate the serial TrilinosWrappers and PETScWrappers classes instead of the constr...
authorBruno Turcksin <bruno.turcksin@gmail.com>
Fri, 1 May 2015 17:43:32 +0000 (12:43 -0500)
committerBruno Turcksin <bruno.turcksin@gmail.com>
Tue, 5 May 2015 19:17:29 +0000 (14:17 -0500)
examples/step-17/step-17.cc
examples/step-18/step-18.cc
examples/step-31/step-31.cc
examples/step-32/step-32.cc
examples/step-36/step-36.cc
examples/step-41/step-41.cc
examples/step-43/step-43.cc
include/deal.II/lac/petsc_block_vector.h
include/deal.II/lac/petsc_vector.h
include/deal.II/lac/trilinos_block_vector.h
include/deal.II/lac/trilinos_vector.h

index 469b00bca7a6ad965f32dc913b010abbfe6be9c1..5298d9dd8a5ee0551c8c3718370e39f1bcee526a 100644 (file)
@@ -65,7 +65,6 @@
 // simply map to sequential, local vectors and matrices if there is only a
 // single process, i.e. if you are running on only one machine, and without
 // MPI support):
-#include <deal.II/lac/petsc_vector.h>
 #include <deal.II/lac/petsc_parallel_vector.h>
 #include <deal.II/lac/petsc_parallel_sparse_matrix.h>
 // Then we also need interfaces for solvers and preconditioners that PETSc
@@ -143,14 +142,9 @@ namespace Step17
 
     // In step-8, this would have been the place where we would have declared
     // the member variables for the sparsity pattern, the system matrix, right
-    // hand, and solution vector. We change these declarations to use
-    // parallel PETSc objects instead (note that the fact that we use the
-    // parallel versions is denoted the fact that we use the classes from the
-    // <code>PETScWrappers::MPI</code> namespace; sequential versions of these
-    // classes are in the <code>PETScWrappers</code> namespace, i.e. without
-    // the <code>MPI</code> part, be aware that these classes are deprecated).
-    // Note also that we do not use a separate sparsity pattern, since PETSc
-    // manages that as part of its matrix data structures.
+    // hand, and solution vector. We change these declarations to use parallel
+    // PETSc objects instead. Note that we do not use a separate sparsity
+    // pattern, since PETSc manages that as part of its matrix data structures.
     PETScWrappers::MPI::SparseMatrix system_matrix;
 
     PETScWrappers::MPI::Vector       solution;
index ea1906eac41a4684862b6f43f97eb993b9943811..38ae6af245822a1f28feee86df5a2341ad5bd5f2 100644 (file)
@@ -30,7 +30,6 @@
 #include <deal.II/lac/vector.h>
 #include <deal.II/lac/full_matrix.h>
 #include <deal.II/lac/dynamic_sparsity_pattern.h>
-#include <deal.II/lac/petsc_vector.h>
 #include <deal.II/lac/petsc_parallel_vector.h>
 #include <deal.II/lac/petsc_parallel_sparse_matrix.h>
 #include <deal.II/lac/petsc_solver.h>
index 653d8c31d036f303d3132c32644afe08a2e8fb27..be44f33ef529928db6b5912a84aa364ec078c1a1 100644 (file)
@@ -393,6 +393,15 @@ namespace Step31
 
 
 
+    // When using a TrilinosWrappers::MPI::Vector or a
+    // TrilinosWrappers::MPI::BlockVector, the Vector is initialized using an
+    // IndexSet. IndexSet is used not only to resize the
+    // TrilinosWrappers::MPI::Vector but it also associates an index in the
+    // TrilinosWrappers::MPI::Vector with a degree of freedom (see step-40 for
+    // a more detailed explanation). The function complete_index_set() creates
+    // an IndexSet where every valid index is part of the set. Note that this
+    // program can only be run sequentially and will throw an exception if used
+    // in parallel.
     template <class PreconditionerA, class PreconditionerMp>
     BlockSchurPreconditioner<PreconditionerA, PreconditionerMp>::
     BlockSchurPreconditioner(const TrilinosWrappers::BlockSparseMatrix  &S,
@@ -402,19 +411,9 @@ namespace Step31
       :
       stokes_matrix           (&S),
       m_inverse               (&Mpinv),
-      a_preconditioner        (Apreconditioner)
-    {
-      // When using a TrilinosWrappers::MPI::Vector or a
-      // TrilinosWrappers::MPI::BlockVector, the Vector is initialized using an
-      // IndexSet. IndexSet is used not only to resize the
-      // TrilinosWrappers::MPI::Vector but it also associates an index in the
-      // TrilinosWrappers::MPI::Vector with a degree of freedom (see step-40 for
-      // a more detailed explanation). This assocation is done by the add_range()
-      // function.
-      IndexSet tmp_index_set(stokes_matrix->block(1,1).m());
-      tmp_index_set.add_range(0,stokes_matrix->block(1,1).m());
-      tmp.reinit(tmp_index_set, MPI_COMM_WORLD);
-    }
+      a_preconditioner        (Apreconditioner),
+      tmp                     (complete_index_set(stokes_matrix->block(1,1).m()))
+    {}
 
 
     // Next is the <code>vmult</code> function. We implement the action of
@@ -513,7 +512,7 @@ namespace Step31
     DoFHandler<dim>                     stokes_dof_handler;
     ConstraintMatrix                    stokes_constraints;
 
-    std::vector<IndexSet>               stokes_block_sizes;
+    std::vector<IndexSet>               stokes_partitioning;
     TrilinosWrappers::BlockSparseMatrix stokes_matrix;
     TrilinosWrappers::BlockSparseMatrix stokes_preconditioner_matrix;
 
@@ -963,12 +962,9 @@ namespace Step31
     // Trilinos matrices store the sparsity pattern internally, there is no
     // need to keep the sparsity pattern around after the initialization of
     // the matrix.
-    stokes_block_sizes.clear();
-    stokes_block_sizes.resize (2);
-    stokes_block_sizes[0].set_size(n_u);
-    stokes_block_sizes[1].set_size(n_p);
-    stokes_block_sizes[0].add_range(0,n_u);
-    stokes_block_sizes[1].add_range(0,n_p);
+    stokes_partitioning.resize (2);
+    stokes_partitioning[0] = complete_index_set (n_u);
+    stokes_partitioning[1] = complete_index_set (n_p);
     {
       stokes_matrix.clear ();
 
@@ -1054,11 +1050,10 @@ namespace Step31
     // and $\mathbf u^{n-2}$, as well as for the temperatures $T^{n}$,
     // $T^{n-1}$ and $T^{n-2}$ (required for time stepping) and all the system
     // right hand sides to their correct sizes and block structure:
-    IndexSet temperature_partitioning (n_T);
-    temperature_partitioning.add_range(0,n_T);
-    stokes_solution.reinit (stokes_block_sizes, MPI_COMM_WORLD);
-    old_stokes_solution.reinit (stokes_block_sizes, MPI_COMM_WORLD);
-    stokes_rhs.reinit (stokes_block_sizes, MPI_COMM_WORLD);
+    IndexSet temperature_partitioning = complete_index_set (n_T);
+    stokes_solution.reinit (stokes_partitioning, MPI_COMM_WORLD);
+    old_stokes_solution.reinit (stokes_partitioning, MPI_COMM_WORLD);
+    stokes_rhs.reinit (stokes_partitioning, MPI_COMM_WORLD);
 
     temperature_solution.reinit (temperature_partitioning, MPI_COMM_WORLD);
     old_temperature_solution.reinit (temperature_partitioning, MPI_COMM_WORLD);
@@ -2219,10 +2214,8 @@ int main (int argc, char *argv[])
                                                            numbers::invalid_unsigned_int);
 
       // This program can only be run in serial. Otherwise, throw an exception.
-      int size;
-      MPI_Comm_size(MPI_COMM_WORLD,&size);
-      AssertThrow(size==1, ExcMessage("This program can only be run in serial,"
-                                      " use mpirun -np 1 ./step-31"));
+      AssertThrow(Utilities::MPI::n_mpi_processes(MPI_COMM_WORLD)==1,
+                  ExcMessage("This program can only be run in serial, use mpirun -np 1 ./step-31"));
 
       BoussinesqFlowProblem<2> flow_problem;
       flow_problem.run ();
index dab38039d2c7c752f5f124e702be17407dced990..baea0d2bc09baf181fe89aae82ce5a39232848a8 100644 (file)
@@ -861,13 +861,8 @@ namespace Step32
     // - In a bit of naming confusion, you will notice below that some of the
     // variables from namespace TrilinosWrappers are taken from namespace
     // TrilinosWrappers::MPI (such as the right hand side vectors) whereas
-    // others are not (such as the various matrices). For the matrices, we
-    // happen to use the same class names for %parallel and sequential data
-    // structures, i.e., all matrices will actually be considered %parallel
-    // below. On the other hand, for vectors, only those from namespace
-    // TrilinosWrappers::MPI are actually distributed (be aware that
-    // TrilinosWrappers::Vector and TrilinosWrappers::BlockVector are
-    // deprecated). In particular, we will frequently have to query velocities
+    // others are not (such as the various matrices). This is due to legacy
+    // reasons. We will frequently have to query velocities
     // and temperatures at arbitrary quadrature points; consequently, rather
     // than importing ghost information of a vector whenever we need access
     // to degrees of freedom that are relevant locally but owned by another
index 6b5a539508e45af8765173db90dbf4d0a001118f..438ae79a62d732370545d0fe15059c1e4f2171ec 100644 (file)
@@ -49,7 +49,7 @@
 
 // PETSc appears here because SLEPc depends on this library:
 #include <deal.II/lac/petsc_sparse_matrix.h>
-#include <deal.II/lac/petsc_vector.h>
+#include <deal.II/lac/petsc_parallel_vector.h>
 
 // And then we need to actually import the interfaces for solvers that SLEPc
 // provides:
@@ -201,14 +201,15 @@ namespace Step36
     // is initialized using an IndexSet. IndexSet is used not only to resize the
     // PETScWrappers::MPI::Vector but it also associates an index in the
     // PETScWrappers::MPI::Vector with a degree of freedom (see step-40 for a
-    // more detailed explanation). This assocation is done by the add_range()
-    // function:
-    IndexSet eigenfunction_index_set(dof_handler.n_dofs ());
-    eigenfunction_index_set.add_range(0, dof_handler.n_dofs ());
+    // a more detailed explanation). The function complete_index_set() creates
+    // an IndexSet where every valid index is part of the set. Note that this
+    // program can only be run sequentially and will throw an exception if used
+    // in parallel.
+    IndexSet eigenfunction_partitioning = complete_index_set(dof_handler.n_dofs ());
     eigenfunctions
     .resize (parameters.get_integer ("Number of eigenvalues/eigenfunctions"));
     for (unsigned int i=0; i<eigenfunctions.size (); ++i)
-      eigenfunctions[i].reinit (eigenfunction_index_set, MPI_COMM_WORLD);
+      eigenfunctions[i].reinit (eigenfunction_partitioning, MPI_COMM_WORLD);
 
     eigenvalues.resize (eigenfunctions.size ());
   }
@@ -477,11 +478,10 @@ int main (int argc, char **argv)
 
       Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1);
 
+
       // This program can only be run in serial. Otherwise, throw an exception.
-      int size;
-      MPI_Comm_size(MPI_COMM_WORLD,&size);
-      AssertThrow(size==1, ExcMessage("This program can only be run in serial,"
-                                      " use mpirun -np 1 ./step-36"));
+      AssertThrow(Utilities::MPI::n_mpi_processes(MPI_COMM_WORLD)==1,
+                  ExcMessage("This program can only be run in serial, use mpirun -np 1 ./step-36"));
 
       {
         deallog.depth_console (0);
index 279b4c9f0bb48d5516ab5ca38113bd2b8b7c3e7b..14ba2f7cff92501cec92e1ad2a23d9acea7243a5 100644 (file)
@@ -94,7 +94,6 @@ namespace Step41
     DoFHandler<dim>      dof_handler;
     ConstraintMatrix     constraints;
     IndexSet             active_set;
-    IndexSet             solution_index_set;
 
     TrilinosWrappers::SparseMatrix system_matrix;
     TrilinosWrappers::SparseMatrix complete_system_matrix;
@@ -256,12 +255,11 @@ namespace Step41
     system_matrix.reinit (dsp);
     complete_system_matrix.reinit (dsp);
 
-    solution_index_set.set_size(dof_handler.n_dofs());
-    solution_index_set.add_range(0, dof_handler.n_dofs());
-    solution.reinit (solution_index_set, MPI_COMM_WORLD);
-    system_rhs.reinit (solution_index_set, MPI_COMM_WORLD);
-    complete_system_rhs.reinit (solution_index_set, MPI_COMM_WORLD);
-    contact_force.reinit (solution_index_set, MPI_COMM_WORLD);
+    IndexSet solution_partitioning = complete_index_set(dof_handler.n_dofs());
+    solution.reinit (solution_partitioning, MPI_COMM_WORLD);
+    system_rhs.reinit (solution_partitioning, MPI_COMM_WORLD);
+    complete_system_rhs.reinit (solution_partitioning, MPI_COMM_WORLD);
+    contact_force.reinit (solution_partitioning, MPI_COMM_WORLD);
 
     // The only other thing to do here is to compute the factors in the $B$
     // matrix which is used to scale the residual. As discussed in the
@@ -271,7 +269,7 @@ namespace Step41
     TrilinosWrappers::SparseMatrix mass_matrix;
     mass_matrix.reinit (dsp);
     assemble_mass_matrix_diagonal (mass_matrix);
-    diagonal_of_mass_matrix.reinit (solution_index_set);
+    diagonal_of_mass_matrix.reinit (solution_partitioning);
     for (unsigned int j=0; j<solution.size (); j++)
       diagonal_of_mass_matrix (j) = mass_matrix.diag_element (j);
   }
@@ -435,7 +433,7 @@ namespace Step41
 
     const double penalty_parameter = 100.0;
 
-    TrilinosWrappers::MPI::Vector lambda (solution_index_set);
+    TrilinosWrappers::MPI::Vector lambda (complete_index_set(dof_handler.n_dofs()));
     complete_system_matrix.residual (lambda,
                                      solution, complete_system_rhs);
     contact_force.ratio (lambda, diagonal_of_mass_matrix);
@@ -675,10 +673,8 @@ int main (int argc, char *argv[])
                                                            numbers::invalid_unsigned_int);
 
       // This program can only be run in serial. Otherwise, throw an exception.
-      int size;
-      MPI_Comm_size(MPI_COMM_WORLD,&size);
-      AssertThrow(size==1, ExcMessage("This program can only be run in serial,"
-                                      " use mpirun -np 1 ./step-41"));
+      AssertThrow(Utilities::MPI::n_mpi_processes(MPI_COMM_WORLD)==1,
+                  ExcMessage("This program can only be run in serial, use mpirun -np 1 ./step-41"));
 
       ObstacleProblem<2> obstacle_problem;
       obstacle_problem.run ();
index dbea7124ebb65af8654e4df8861316af3ca85eee..d17e8db85b0091e7e5bdf8813580d979be2df95c 100644 (file)
@@ -445,12 +445,9 @@ namespace Step43
       :
       darcy_matrix            (&S),
       m_inverse               (&Mpinv),
-      a_preconditioner        (Apreconditioner)
-    {
-      IndexSet tmp_index_set(darcy_matrix->block(1,1).m());
-      tmp_index_set.add_range(0,darcy_matrix->block(1,1).m());
-      tmp.reinit(tmp_index_set, MPI_COMM_WORLD);
-    }
+      a_preconditioner        (Apreconditioner),
+      tmp                     (complete_index_set(darcy_matrix->block(1,1).m()))
+    {}
 
 
     template <class PreconditionerA, class PreconditionerMp>
@@ -550,7 +547,6 @@ namespace Step43
     FESystem<dim>                        darcy_fe;
     DoFHandler<dim>                      darcy_dof_handler;
     ConstraintMatrix                     darcy_constraints;
-    std::vector<IndexSet>                darcy_index_set;
 
     ConstraintMatrix                     darcy_preconditioner_constraints;
 
@@ -568,7 +564,6 @@ namespace Step43
     FE_Q<dim>                            saturation_fe;
     DoFHandler<dim>                      saturation_dof_handler;
     ConstraintMatrix                     saturation_constraints;
-    IndexSet                             saturation_index_set;
 
     TrilinosWrappers::SparseMatrix       saturation_matrix;
 
@@ -808,35 +803,30 @@ namespace Step43
       saturation_matrix.reinit (dsp);
     }
 
-    darcy_index_set.clear();
-    darcy_index_set.resize(2);
-    darcy_index_set[0].set_size(n_u);
-    darcy_index_set[1].set_size(n_p);
-    darcy_index_set[0].add_range(0,n_u);
-    darcy_index_set[1].add_range(0,n_p);
-    darcy_solution.reinit (darcy_index_set, MPI_COMM_WORLD);
+    std::vector<IndexSet> darcy_partitioning(2);
+    darcy_partitioning[0] = complete_index_set (n_u);
+    darcy_partitioning[1] = complete_index_set (n_p);
+    darcy_solution.reinit (darcy_partitioning, MPI_COMM_WORLD);
     darcy_solution.collect_sizes ();
 
-    last_computed_darcy_solution.reinit (darcy_index_set, MPI_COMM_WORLD);
+    last_computed_darcy_solution.reinit (darcy_partitioning, MPI_COMM_WORLD);
     last_computed_darcy_solution.collect_sizes ();
 
-    second_last_computed_darcy_solution.reinit (darcy_index_set, MPI_COMM_WORLD);
+    second_last_computed_darcy_solution.reinit (darcy_partitioning, MPI_COMM_WORLD);
     second_last_computed_darcy_solution.collect_sizes ();
 
-    darcy_rhs.reinit (darcy_index_set, MPI_COMM_WORLD);
+    darcy_rhs.reinit (darcy_partitioning, MPI_COMM_WORLD);
     darcy_rhs.collect_sizes ();
 
-    saturation_index_set.clear();
-    saturation_index_set.set_size(n_s);
-    saturation_index_set.add_range(0,n_s);
-    saturation_solution.reinit (saturation_index_set, MPI_COMM_WORLD);
-    old_saturation_solution.reinit (saturation_index_set, MPI_COMM_WORLD);
-    old_old_saturation_solution.reinit (saturation_index_set, MPI_COMM_WORLD);
+    IndexSet saturation_partitioning = complete_index_set(n_s);
+    saturation_solution.reinit (saturation_partitioning, MPI_COMM_WORLD);
+    old_saturation_solution.reinit (saturation_partitioning, MPI_COMM_WORLD);
+    old_old_saturation_solution.reinit (saturation_partitioning, MPI_COMM_WORLD);
 
-    saturation_matching_last_computed_darcy_solution.reinit (saturation_index_set,
+    saturation_matching_last_computed_darcy_solution.reinit (saturation_partitioning,
                                                              MPI_COMM_WORLD);
 
-    saturation_rhs.reinit (saturation_index_set, MPI_COMM_WORLD);
+    saturation_rhs.reinit (saturation_partitioning, MPI_COMM_WORLD);
   }
 
 
@@ -2262,10 +2252,8 @@ int main (int argc, char *argv[])
                                                            numbers::invalid_unsigned_int);
 
       // This program can only be run in serial. Otherwise, throw an exception.
-      int size;
-      MPI_Comm_size(MPI_COMM_WORLD,&size);
-      AssertThrow(size==1, ExcMessage("This program can only be run in serial,"
-                                      " use mpirun -np 1 ./step-43"));
+      AssertThrow(Utilities::MPI::n_mpi_processes(MPI_COMM_WORLD)==1,
+                  ExcMessage("This program can only be run in serial, use mpirun -np 1 ./step-43"));
 
       TwoPhaseFlowProblem<2> two_phase_flow_problem(1);
       two_phase_flow_problem.run ();
index 38080427dd547ea7e31b677574a1957e880bfc47..084218813f31cccfc98132e031ad16f6179acb6e 100644 (file)
@@ -43,7 +43,7 @@ namespace PETScWrappers
    * interface, this class handles the actual allocation of vectors and
    * provides functions that are specific to the underlying vector type.
    *
-   * This class is deprecated use PETScWrappers::MPI::BlockVector.
+   * This class is deprecated, use PETScWrappers::MPI::BlockVector.
    *
    * @ingroup Vectors
    *
@@ -87,13 +87,13 @@ namespace PETScWrappers
      * of different sizes.
      */
     explicit BlockVector (const unsigned int num_blocks = 0,
-                          const size_type    block_size = 0) DEAL_II_DEPRECATED;
+                          const size_type    block_size = 0);
 
     /**
      * Copy-Constructor. Dimension set to that of V, all components are copied
      * from V
      */
-    BlockVector (const BlockVector  &V) DEAL_II_DEPRECATED;
+    BlockVector (const BlockVector  &V);
 
     /**
      * Copy-constructor: copy the values from a PETSc wrapper parallel block
@@ -105,13 +105,13 @@ namespace PETScWrappers
      * It is not sufficient if only one processor tries to copy the elements
      * from the other processors over to its own process space.
      */
-    explicit BlockVector (const MPI::BlockVector &v) DEAL_II_DEPRECATED;
+    explicit BlockVector (const MPI::BlockVector &v);
 
     /**
      * Constructor. Set the number of blocks to <tt>n.size()</tt> and
      * initialize each block with <tt>n[i]</tt> zero elements.
      */
-    BlockVector (const std::vector<size_type> &n) DEAL_II_DEPRECATED;
+    BlockVector (const std::vector<size_type> &n);
 
     /**
      * Constructor. Set the number of blocks to <tt>n.size()</tt>. Initialize
@@ -124,7 +124,7 @@ namespace PETScWrappers
     template <typename InputIterator>
     BlockVector (const std::vector<size_type> &n,
                  const InputIterator           first,
-                 const InputIterator           end) DEAL_II_DEPRECATED;
+                 const InputIterator           end);
 
     /**
      * Destructor. Clears memory
@@ -244,7 +244,7 @@ namespace PETScWrappers
      */
     DeclException0 (ExcIteratorRangeDoesNotMatchVectorSize);
     ///@}
-  };
+  } DEAL_II_DEPRECATED;
 
   /*@}*/
 
index 5e953c28108c21bd06c627046e1d84fe3aa6ff7e..c923ffc22e8ff3b8ea8870c9a7ace37704a55994 100644 (file)
@@ -45,7 +45,7 @@ namespace PETScWrappers
    * virtual functions). Only the functions creating a vector of specific type
    * differ, and are implemented in this particular class.
    *
-   * This class is deprecated use PETScWrappers::MPI::Vector instead.
+   * This class is deprecated, use PETScWrappers::MPI::Vector instead.
    *
    * @ingroup Vectors
    * @author Wolfgang Bangerth, 2004
@@ -73,7 +73,7 @@ namespace PETScWrappers
     /**
      * Default constructor. Initialize the vector as empty.
      */
-    Vector () DEAL_II_DEPRECATED;
+    Vector ();
 
     /**
      * Constructor. Set dimension to @p n and initialize all elements with
@@ -85,14 +85,14 @@ namespace PETScWrappers
      * <tt>v=Vector@<number@>(0);</tt>, i.e. the vector is replaced by one of
      * length zero.
      */
-    explicit Vector (const size_type n) DEAL_II_DEPRECATED;
+    explicit Vector (const size_type n);
 
     /**
      * Copy-constructor from deal.II vectors. Sets the dimension to that of
      * the given vector, and copies all elements.
      */
     template <typename Number>
-    explicit Vector (const dealii::Vector<Number> &v) DEAL_II_DEPRECATED;
+    explicit Vector (const dealii::Vector<Number> &v);
 
     /**
      * Construct it from an existing PETSc Vector of type Vec. Note: this does
@@ -100,12 +100,12 @@ namespace PETScWrappers
      * the vector is not used twice at the same time or destroyed while in
      * use. This class does not destroy the PETSc object. Handle with care!
      */
-    explicit Vector (const Vec &v) DEAL_II_DEPRECATED;
+    explicit Vector (const Vec &v);
 
     /**
      * Copy-constructor the values from a PETSc wrapper vector class.
      */
-    Vector (const Vector &v) DEAL_II_DEPRECATED;
+    Vector (const Vector &v);
 
     /**
      * Copy-constructor: copy the values from a PETSc wrapper parallel vector
@@ -116,7 +116,7 @@ namespace PETScWrappers
      * It is not sufficient if only one processor tries to copy the elements
      * from the other processors over to its own process space.
      */
-    explicit Vector (const MPI::Vector &v) DEAL_II_DEPRECATED;
+    explicit Vector (const MPI::Vector &v);
 
     /**
      * Copy the given vector. Resize the present vector if necessary.
@@ -183,7 +183,7 @@ namespace PETScWrappers
      * vector. @p n denotes the total size of the vector to be created.
      */
     void create_vector (const size_type n);
-  };
+  } DEAL_II_DEPRECATED;
 
   /*@}*/
 
index f9eb32422f855ab553d31467e11854e7831c014d..0d388d7c25338d1b1710b85878c66e38143a756b 100644 (file)
@@ -59,7 +59,7 @@ namespace TrilinosWrappers
    * block vector class do only work in case the program is run on only one
    * processor, since the Trilinos matrices are inherently parallel.
    *
-   * This class is deprecated use TrilinosWrappers::MPI::BlockVector instead.
+   * This class is deprecated, use TrilinosWrappers::MPI::BlockVector instead.
    *
    * @ingroup Vectors
    * @ingroup TrilinosWrappers @see
@@ -94,14 +94,14 @@ namespace TrilinosWrappers
     /**
      * Default constructor. Generate an empty vector without any blocks.
      */
-    BlockVector () DEAL_II_DEPRECATED;
+    BlockVector ();
 
     /**
      * Constructor. Generate a block vector with as many blocks as there are
      * entries in Input_Maps.  For this non-distributed vector, the %parallel
      * partitioning is not used, just the global size of the partitioner.
      */
-    explicit BlockVector (const std::vector<Epetra_Map> &partitioner) DEAL_II_DEPRECATED;
+    explicit BlockVector (const std::vector<Epetra_Map> &partitioner);
 
     /**
      * Constructor. Generate a block vector with as many blocks as there are
@@ -109,26 +109,26 @@ namespace TrilinosWrappers
      * partitioning is not used, just the global size of the partitioner.
      */
     explicit BlockVector (const std::vector<IndexSet> &partitioner,
-                          const MPI_Comm              &communicator = MPI_COMM_WORLD) DEAL_II_DEPRECATED;
+                          const MPI_Comm              &communicator = MPI_COMM_WORLD);
 
     /**
      * Copy-Constructor. Set all the properties of the non-%parallel vector to
      * those of the given %parallel vector and import the elements.
      */
-    BlockVector (const MPI::BlockVector &V) DEAL_II_DEPRECATED;
+    BlockVector (const MPI::BlockVector &V);
 
     /**
      * Copy-Constructor. Set all the properties of the vector to those of the
      * given input vector and copy the elements.
      */
-    BlockVector (const BlockVector  &V) DEAL_II_DEPRECATED;
+    BlockVector (const BlockVector  &V);
 
     /**
      * Creates a block vector consisting of <tt>num_blocks</tt> components,
      * but there is no content in the individual components and the user has
      * to fill appropriate data using a reinit of the blocks.
      */
-    explicit BlockVector (const size_type num_blocks) DEAL_II_DEPRECATED;
+    explicit BlockVector (const size_type num_blocks);
 
     /**
      * Constructor. Set the number of blocks to <tt>n.size()</tt> and
@@ -136,7 +136,7 @@ namespace TrilinosWrappers
      *
      * References BlockVector.reinit().
      */
-    explicit BlockVector (const std::vector<size_type> &N) DEAL_II_DEPRECATED;
+    explicit BlockVector (const std::vector<size_type> &N);
 
     /**
      * Constructor. Set the number of blocks to <tt>n.size()</tt>. Initialize
@@ -149,7 +149,7 @@ namespace TrilinosWrappers
     template <typename InputIterator>
     BlockVector (const std::vector<size_type> &n,
                  const InputIterator           first,
-                 const InputIterator           end) DEAL_II_DEPRECATED;
+                 const InputIterator           end);
 
     /**
      * Destructor. Clears memory
@@ -303,7 +303,7 @@ namespace TrilinosWrappers
                     << "local_size = global_size is a necessary condition, but"
                     << arg1 << " != " << arg2 << " was given!");
 
-  };
+  } DEAL_II_DEPRECATED;
 
 
 
index 241a334d84e045fa4a7eb8adf8518e3af61c651d..447065b813a97ce1240435c2b98425d5d96b05d7 100644 (file)
@@ -710,7 +710,7 @@ namespace TrilinosWrappers
    * in order to be able to access all elements in the vector or to apply
    * certain deal.II functions.
    *
-   * This class is deprecated use TrilinosWrappers::MPI::Vector instead.
+   * This class is deprecated, use TrilinosWrappers::MPI::Vector instead.
    *
    * @ingroup TrilinosWrappers
    * @ingroup Vectors
@@ -741,12 +741,12 @@ namespace TrilinosWrappers
      * function <tt>reinit()</tt> will have to give the vector the correct
      * size.
      */
-    Vector () DEAL_II_DEPRECATED;
+    Vector ();
 
     /**
      * This constructor takes as input the number of elements in the vector.
      */
-    explicit Vector (const size_type n) DEAL_II_DEPRECATED;
+    explicit Vector (const size_type n);
 
     /**
      * This constructor takes as input the number of elements in the vector.
@@ -757,7 +757,7 @@ namespace TrilinosWrappers
      * ignored, the only thing that matters is the size of the index space
      * described by this argument.
      */
-    explicit Vector (const Epetra_Map &partitioning) DEAL_II_DEPRECATED;
+    explicit Vector (const Epetra_Map &partitioning);
 
     /**
      * This constructor takes as input the number of elements in the vector.
@@ -769,20 +769,20 @@ namespace TrilinosWrappers
      * size of the index space described by this argument.
      */
     explicit Vector (const IndexSet &partitioning,
-                     const MPI_Comm &communicator = MPI_COMM_WORLD) DEAL_II_DEPRECATED;
+                     const MPI_Comm &communicator = MPI_COMM_WORLD);
 
     /**
      * This constructor takes a (possibly parallel) Trilinos Vector and
      * generates a localized version of the whole content on each processor.
      */
-    explicit Vector (const VectorBase &V) DEAL_II_DEPRECATED;
+    explicit Vector (const VectorBase &V);
 
     /**
      * Copy-constructor from deal.II vectors. Sets the dimension to that of
      * the given vector, and copies all elements.
      */
     template <typename Number>
-    explicit Vector (const dealii::Vector<Number> &v) DEAL_II_DEPRECATED;
+    explicit Vector (const dealii::Vector<Number> &v);
 
     /**
      * Reinit function that resizes the vector to the size specified by
@@ -870,7 +870,7 @@ namespace TrilinosWrappers
      * thus an empty function.
      */
     void update_ghost_values () const;
-  };
+  } DEAL_II_DEPRECATED;
 
 
 

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.