From 8ac29612de4d6a6c71b7a298309d1b21f9d85359 Mon Sep 17 00:00:00 2001
From: Henry <Henry>
Date: Fri, 9 Jan 2015 22:05:32 +0000
Subject: [PATCH] advectionDiffusionPatchDistMethod: Only run the predictor
 method on first call Improves efficiency for moving-mesh cases

---
 .../advectionDiffusionPatchDistMethod.C              | 12 +++++++++---
 .../advectionDiffusionPatchDistMethod.H              |  3 +++
 2 files changed, 12 insertions(+), 3 deletions(-)

diff --git a/src/finiteVolume/fvMesh/wallDist/patchDistMethods/advectionDiffusion/advectionDiffusionPatchDistMethod.C b/src/finiteVolume/fvMesh/wallDist/patchDistMethods/advectionDiffusion/advectionDiffusionPatchDistMethod.C
index 85100b39ecc..a3797f1fa2e 100644
--- a/src/finiteVolume/fvMesh/wallDist/patchDistMethods/advectionDiffusion/advectionDiffusionPatchDistMethod.C
+++ b/src/finiteVolume/fvMesh/wallDist/patchDistMethods/advectionDiffusion/advectionDiffusionPatchDistMethod.C
@@ -68,7 +68,8 @@ Foam::patchDistMethods::advectionDiffusion::advectionDiffusion
     ),
     epsilon_(coeffs_.lookupOrDefault<scalar>("epsilon", 0.1)),
     tolerance_(coeffs_.lookupOrDefault<scalar>("tolerance", 1e-3)),
-    maxIter_(coeffs_.lookupOrDefault<int>("maxIter", 10))
+    maxIter_(coeffs_.lookupOrDefault<int>("maxIter", 10)),
+    predicted_(false)
 {}
 
 
@@ -86,7 +87,11 @@ bool Foam::patchDistMethods::advectionDiffusion::correct
     volVectorField& n
 )
 {
-    pdmPredictor_->correct(y);
+    if (!predicted_)
+    {
+        pdmPredictor_->correct(y);
+        predicted_ = true;
+    }
 
     volVectorField ny
     (
@@ -97,7 +102,7 @@ bool Foam::patchDistMethods::advectionDiffusion::correct
             mesh_
         ),
         mesh_,
-        dimensionedVector("nWall", dimless, vector::zero),
+        dimensionedVector("ny", dimless, vector::zero),
         patchTypes<vector>(mesh_, patchIDs_)
     );
 
@@ -133,6 +138,7 @@ bool Foam::patchDistMethods::advectionDiffusion::correct
 
         yEqn.relax();
         initialResidual = yEqn.solve().initialResidual();
+
     } while (initialResidual > tolerance_ && ++iter < maxIter_);
 
     // Only calculate n if the field is defined
diff --git a/src/finiteVolume/fvMesh/wallDist/patchDistMethods/advectionDiffusion/advectionDiffusionPatchDistMethod.H b/src/finiteVolume/fvMesh/wallDist/patchDistMethods/advectionDiffusion/advectionDiffusionPatchDistMethod.H
index 590ae65757c..c62f7e79f9f 100644
--- a/src/finiteVolume/fvMesh/wallDist/patchDistMethods/advectionDiffusion/advectionDiffusionPatchDistMethod.H
+++ b/src/finiteVolume/fvMesh/wallDist/patchDistMethods/advectionDiffusion/advectionDiffusionPatchDistMethod.H
@@ -162,6 +162,9 @@ class advectionDiffusion
         //  to correct the distance-to-patch and normal-to-patch fields
         int maxIter_;
 
+        //- Flag to indicate the predictor step has been executed
+        bool predicted_;
+
 
     // Private Member Functions
 
-- 
GitLab