From 8cd5eb72d2b4dc0d3ebe20a2344f472097ff5568 Mon Sep 17 00:00:00 2001
From: Bernd Flemisch <bernd@iws.uni-stuttgart.de>
Date: Thu, 25 Oct 2012 12:35:03 +0000
Subject: [PATCH] Bring the joy of parallelism also to the common decoupled
 stuff. Fix the pressure matrix assembly in parallel if faces are visited only
 once. Can be tested with test/decoupled/2p/test_parallel in dumux-devel.
 Reviewed by Markus.

git-svn-id: svn://svn.iws.uni-stuttgart.de/DUMUX/dumux/trunk@9434 2fb0f335-1f38-0410-981e-8018bf24f1b0
---
 dumux/decoupled/common/fv/fvpressure.hh  | 18 ++++++++++++------
 dumux/decoupled/common/fv/fvtransport.hh | 13 +++++++++++++
 2 files changed, 25 insertions(+), 6 deletions(-)

diff --git a/dumux/decoupled/common/fv/fvpressure.hh b/dumux/decoupled/common/fv/fvpressure.hh
index 1b05b20502..49674e43dc 100644
--- a/dumux/decoupled/common/fv/fvpressure.hh
+++ b/dumux/decoupled/common/fv/fvpressure.hh
@@ -393,11 +393,15 @@ void FVPressure<TypeTag>::assemble(bool first)
 
                     int globalIdxJ = problem_.variables().index(*elementNeighbor);
 
-                    //check for hanging nodes
-                    //take a hanging node never from the element with smaller level!
+                    // check for hanging nodes
+                    // take a hanging node never from the element with smaller level!
                     bool haveSameLevel = (eIt->level() == elementNeighbor->level());
-                    //calculate only from one side, but add matrix entries for both sides
-                    if (GET_PROP_VALUE(TypeTag, VisitFacesOnlyOnce) && (globalIdxI > globalIdxJ) && haveSameLevel)
+                    // calculate only from one side, but add matrix entries for both sides
+                    // the last condition is needed to properly assemble in the presence 
+                    // of ghost elements
+                    if (GET_PROP_VALUE(TypeTag, VisitFacesOnlyOnce) 
+                        && (globalIdxI > globalIdxJ) && haveSameLevel
+                        && elementNeighbor->partitionType() == Dune::InteriorEntity)
                         continue;
 
                     //check for hanging nodes
@@ -410,10 +414,12 @@ void FVPressure<TypeTag>::assemble(bool first)
                     // set diagonal entry
                     A_[globalIdxI][globalIdxI] += entries[matrix];
 
-                        // set off-diagonal entry
+                    // set off-diagonal entry
                     A_[globalIdxI][globalIdxJ] -= entries[matrix];
 
-                    if (GET_PROP_VALUE(TypeTag, VisitFacesOnlyOnce))
+                    // The second condition is needed to not spoil the ghost element entries
+                    if (GET_PROP_VALUE(TypeTag, VisitFacesOnlyOnce) 
+                        && elementNeighbor->partitionType() == Dune::InteriorEntity)
                     {
                         f_[globalIdxJ] += entries[rhs];
                         A_[globalIdxJ][globalIdxJ] += entries[matrix];
diff --git a/dumux/decoupled/common/fv/fvtransport.hh b/dumux/decoupled/common/fv/fvtransport.hh
index 5a260d61fb..8123012022 100644
--- a/dumux/decoupled/common/fv/fvtransport.hh
+++ b/dumux/decoupled/common/fv/fvtransport.hh
@@ -22,6 +22,7 @@
 #include <dune/grid/common/gridenums.hh>
 #include <dumux/decoupled/common/transportproperties.hh>
 #include <dumux/decoupled/common/decoupledproperties.hh>
+#include <dumux/linear/vectorexchange.hh>
 
 /**
  * @file
@@ -279,6 +280,18 @@ void FVTransport<TypeTag>::update(const Scalar t, Scalar& dt, TransportSolutionT
         //store update
         cellDataI.setUpdate(updateVec[globalIdxI]);
     } // end grid traversal
+    
+#if HAVE_MPI        
+    // communicate updated values
+    typedef typename GET_PROP(TypeTag, SolutionTypes) SolutionTypes;
+    typedef typename SolutionTypes::ElementMapper ElementMapper;
+    typedef VectorExchange<ElementMapper, Dune::BlockVector<Dune::FieldVector<Scalar, 1> > > DataHandle;
+    DataHandle dataHandle(problem_.variables().elementMapper(), updateVec);
+    problem_.gridView().template communicate<DataHandle>(dataHandle, 
+                                                         Dune::InteriorBorder_All_Interface, 
+                                                         Dune::ForwardCommunication);
+    dt = problem_.gridView().comm().min(dt);
+#endif
 }
 
 }
-- 
GitLab