]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Partial merge from the parallel_mg branch: Introduce Tvmult for Trilinos preconditioners.
authorWolfgang Bangerth <bangerth@math.tamu.edu>
Mon, 5 Nov 2012 00:52:50 +0000 (00:52 +0000)
committerWolfgang Bangerth <bangerth@math.tamu.edu>
Mon, 5 Nov 2012 00:52:50 +0000 (00:52 +0000)
git-svn-id: https://svn.dealii.org/trunk@27379 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/doc/news/changes.h
deal.II/include/deal.II/lac/trilinos_precondition.h

index 0a6fcfc76bc20e6d405d294af4d9105f287d027b..db2c26e5a215106e688e914c1e610afe17ad4dfb 100644 (file)
@@ -154,6 +154,12 @@ DoFHandler, in particular removal of specializations.
 <h3>Specific improvements</h3>
 
 <ol>
+<li> New: The TrilinosWrappers::PreconditionBase class now has
+a function TrilinosWrappers::PreconditionBase::Tvmult that
+allows applying the transpose preconditioner.
+<br>
+(Guido Kanschat, 2012/11/04)
+
 <li> New: The parallel::distributed::Triangulation::n_global_levels()
 function returns the maximal refinement level over all involved
 processors.
index aed66aaded1ddbe417682e048802ad220746368c..025817f8bec7fade23172dfc0e47beaec7397074 100644 (file)
@@ -113,6 +113,12 @@ namespace TrilinosWrappers
                                         * Apply the preconditioner.
                                         */
       void vmult (VectorBase       &dst,
+                  const VectorBase &src) const;
+
+                                       /**
+                                        * Apply transpose the preconditioner.
+                                        */
+      void Tvmult (VectorBase       &dst,
                   const VectorBase &src) const;
 
                                        /**
@@ -123,6 +129,16 @@ namespace TrilinosWrappers
                                         * class.
                                         */
       void vmult (dealii::Vector<double>       &dst,
+                  const dealii::Vector<double> &src) const;
+
+                                       /**
+                                        * Apply the transpose preconditioner on
+                                        * deal.II data structures
+                                        * instead of the ones provided
+                                        * in the Trilinos wrapper
+                                        * class.
+                                        */
+      void Tvmult (dealii::Vector<double>       &dst,
                   const dealii::Vector<double> &src) const;
 
                                        /**
@@ -132,6 +148,15 @@ namespace TrilinosWrappers
                                         * wrapper class.
                                         */
       void vmult (dealii::parallel::distributed::Vector<double>       &dst,
+                  const dealii::parallel::distributed::Vector<double> &src) const;
+
+                                       /**
+                                        * Apply the transpose preconditioner on deal.II
+                                        * parallel data structures instead of
+                                        * the ones provided in the Trilinos
+                                        * wrapper class.
+                                        */
+      void Tvmult (dealii::parallel::distributed::Vector<double>       &dst,
                   const dealii::parallel::distributed::Vector<double> &src) const;
 
                                        /**
@@ -1488,6 +1513,23 @@ namespace TrilinosWrappers
     AssertThrow (ierr == 0, ExcTrilinosError(ierr));
   }
 
+  inline
+  void
+  PreconditionBase::Tvmult (VectorBase       &dst,
+                           const VectorBase &src) const
+  {
+    Assert (dst.vector_partitioner().SameAs(preconditioner->OperatorRangeMap()),
+            ExcNonMatchingMaps("dst"));
+    Assert (src.vector_partitioner().SameAs(preconditioner->OperatorDomainMap()),
+            ExcNonMatchingMaps("src"));
+    
+    preconditioner->SetUseTranspose(true);
+    const int ierr = preconditioner->ApplyInverse (src.trilinos_vector(),
+                                                   dst.trilinos_vector());
+    AssertThrow (ierr == 0, ExcTrilinosError(ierr));
+    preconditioner->SetUseTranspose(false);
+  }
+
 
                                         // For the implementation of
                                         // the <code>vmult</code>
@@ -1531,6 +1573,26 @@ namespace TrilinosWrappers
   }
 
 
+  inline
+  void PreconditionBase::Tvmult (dealii::Vector<double>       &dst,
+                                const dealii::Vector<double> &src) const
+  {
+    AssertDimension (static_cast<int>(dst.size()),
+                     preconditioner->OperatorDomainMap().NumMyElements());
+    AssertDimension (static_cast<int>(src.size()),
+                     preconditioner->OperatorRangeMap().NumMyElements());
+    Epetra_Vector tril_dst (View, preconditioner->OperatorDomainMap(),
+                            dst.begin());
+    Epetra_Vector tril_src (View, preconditioner->OperatorRangeMap(),
+                            const_cast<double*>(src.begin()));
+
+    preconditioner->SetUseTranspose(true);
+    const int ierr = preconditioner->ApplyInverse (tril_src, tril_dst);
+    AssertThrow (ierr == 0, ExcTrilinosError(ierr));
+    preconditioner->SetUseTranspose(false);
+  }
+
+
 
   inline
   void
@@ -1550,6 +1612,26 @@ namespace TrilinosWrappers
     AssertThrow (ierr == 0, ExcTrilinosError(ierr));
   }
 
+  inline
+  void
+  PreconditionBase::Tvmult (parallel::distributed::Vector<double>       &dst,
+                           const parallel::distributed::Vector<double> &src) const
+  {
+    AssertDimension (static_cast<int>(dst.local_size()),
+                     preconditioner->OperatorDomainMap().NumMyElements());
+    AssertDimension (static_cast<int>(src.local_size()),
+                     preconditioner->OperatorRangeMap().NumMyElements());
+    Epetra_Vector tril_dst (View, preconditioner->OperatorDomainMap(),
+                            dst.begin());
+    Epetra_Vector tril_src (View, preconditioner->OperatorRangeMap(),
+                            const_cast<double*>(src.begin()));
+
+    preconditioner->SetUseTranspose(true);
+    const int ierr = preconditioner->ApplyInverse (tril_src, tril_dst);
+    AssertThrow (ierr == 0, ExcTrilinosError(ierr));
+    preconditioner->SetUseTranspose(false);
+  }
+
 #endif
 
 }

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.