From e2294e539bedf050c0642fc49a6fb2af9ed238cd Mon Sep 17 00:00:00 2001
From: andy <a.heather@opencfd.co.uk>
Date: Wed, 20 Oct 2010 13:28:07 +0100
Subject: [PATCH] ENH: Updated thermo trackData to better handle Cp

---
 .../Templates/ReactingCloud/ReactingCloud.C   | 10 +--------
 .../ReactingMultiphaseCloud.C                 | 10 +--------
 .../Templates/ThermoCloud/ThermoCloud.C       | 10 +--------
 .../ReactingMultiphaseParcel.H                |  1 -
 .../ReactingMultiphaseParcelI.H               |  3 +--
 .../Templates/ReactingParcel/ReactingParcel.H |  1 -
 .../ReactingParcel/ReactingParcelI.H          |  3 +--
 .../Templates/ThermoParcel/ThermoParcel.H     |  6 ++++-
 .../Templates/ThermoParcel/ThermoParcelI.H    | 22 +++++++++++++++++--
 9 files changed, 30 insertions(+), 36 deletions(-)

diff --git a/src/lagrangian/intermediate/clouds/Templates/ReactingCloud/ReactingCloud.C b/src/lagrangian/intermediate/clouds/Templates/ReactingCloud/ReactingCloud.C
index cf3fec1732d..dfec32bc4bc 100644
--- a/src/lagrangian/intermediate/clouds/Templates/ReactingCloud/ReactingCloud.C
+++ b/src/lagrangian/intermediate/clouds/Templates/ReactingCloud/ReactingCloud.C
@@ -67,15 +67,7 @@ void Foam::ReactingCloud<ParcelType>::preEvolve()
 template<class ParcelType>
 void Foam::ReactingCloud<ParcelType>::evolveCloud()
 {
-    const volScalarField Cp = this->thermo().thermo().Cp();
-
-    autoPtr<interpolation<scalar> > CpInterp = interpolation<scalar>::New
-    (
-        this->solution().interpolationSchemes(),
-        Cp
-    );
-
-    typename ParcelType::trackData td(*this, CpInterp());
+    typename ParcelType::trackData td(*this);
 
     label preInjectionSize = this->size();
 
diff --git a/src/lagrangian/intermediate/clouds/Templates/ReactingMultiphaseCloud/ReactingMultiphaseCloud.C b/src/lagrangian/intermediate/clouds/Templates/ReactingMultiphaseCloud/ReactingMultiphaseCloud.C
index 9dc75d51389..b383f38f10a 100644
--- a/src/lagrangian/intermediate/clouds/Templates/ReactingMultiphaseCloud/ReactingMultiphaseCloud.C
+++ b/src/lagrangian/intermediate/clouds/Templates/ReactingMultiphaseCloud/ReactingMultiphaseCloud.C
@@ -40,15 +40,7 @@ void Foam::ReactingMultiphaseCloud<ParcelType>::preEvolve()
 template<class ParcelType>
 void Foam::ReactingMultiphaseCloud<ParcelType>::evolveCloud()
 {
-    const volScalarField Cp = this->thermo().thermo().Cp();
-
-    autoPtr<interpolation<scalar> > CpInterp = interpolation<scalar>::New
-    (
-        this->solution().interpolationSchemes(),
-        Cp
-    );
-
-    typename ParcelType::trackData td(*this, CpInterp());
+    typename ParcelType::trackData td(*this);
 
     label preInjectionSize = this->size();
 
diff --git a/src/lagrangian/intermediate/clouds/Templates/ThermoCloud/ThermoCloud.C b/src/lagrangian/intermediate/clouds/Templates/ThermoCloud/ThermoCloud.C
index 5296583c2b0..17469774af3 100644
--- a/src/lagrangian/intermediate/clouds/Templates/ThermoCloud/ThermoCloud.C
+++ b/src/lagrangian/intermediate/clouds/Templates/ThermoCloud/ThermoCloud.C
@@ -41,15 +41,7 @@ void Foam::ThermoCloud<ParcelType>::preEvolve()
 template<class ParcelType>
 void Foam::ThermoCloud<ParcelType>::evolveCloud()
 {
-    const volScalarField Cp = thermo_.thermo().Cp();
-
-    autoPtr<interpolation<scalar> > CpInterp = interpolation<scalar>::New
-    (
-        this->solution().interpolationSchemes(),
-        Cp
-    );
-
-    typename ParcelType::trackData td(*this, CpInterp());
+    typename ParcelType::trackData td(*this);
 
     label preInjectionSize = this->size();
 
diff --git a/src/lagrangian/intermediate/parcels/Templates/ReactingMultiphaseParcel/ReactingMultiphaseParcel.H b/src/lagrangian/intermediate/parcels/Templates/ReactingMultiphaseParcel/ReactingMultiphaseParcel.H
index b4b9ff213a2..97cfb935dc4 100644
--- a/src/lagrangian/intermediate/parcels/Templates/ReactingMultiphaseParcel/ReactingMultiphaseParcel.H
+++ b/src/lagrangian/intermediate/parcels/Templates/ReactingMultiphaseParcel/ReactingMultiphaseParcel.H
@@ -127,7 +127,6 @@ public:
             inline trackData
             (
                 ReactingMultiphaseCloud<ParcelType>& cloud,
-                const interpolation<scalar>& CpInterp,
                 typename ReactingParcel<ParcelType>::trackData::trackPart
                     part = ReactingParcel<ParcelType>::trackData::tpLinearTrack
              );
diff --git a/src/lagrangian/intermediate/parcels/Templates/ReactingMultiphaseParcel/ReactingMultiphaseParcelI.H b/src/lagrangian/intermediate/parcels/Templates/ReactingMultiphaseParcel/ReactingMultiphaseParcelI.H
index 454ec13f0ef..4d67ca43d46 100644
--- a/src/lagrangian/intermediate/parcels/Templates/ReactingMultiphaseParcel/ReactingMultiphaseParcelI.H
+++ b/src/lagrangian/intermediate/parcels/Templates/ReactingMultiphaseParcel/ReactingMultiphaseParcelI.H
@@ -54,11 +54,10 @@ template<class ParcelType>
 inline Foam::ReactingMultiphaseParcel<ParcelType>::trackData::trackData
 (
     ReactingMultiphaseCloud<ParcelType>& cloud,
-    const interpolation<scalar>& CpInterp,
     typename ReactingParcel<ParcelType>::trackData::trackPart part
 )
 :
-    ReactingParcel<ParcelType>::trackData(cloud, CpInterp, part),
+    ReactingParcel<ParcelType>::trackData(cloud, part),
     cloud_(cloud),
     constProps_(cloud.constProps())
 {}
diff --git a/src/lagrangian/intermediate/parcels/Templates/ReactingParcel/ReactingParcel.H b/src/lagrangian/intermediate/parcels/Templates/ReactingParcel/ReactingParcel.H
index 0d6cc119ab5..d4423bc8c14 100644
--- a/src/lagrangian/intermediate/parcels/Templates/ReactingParcel/ReactingParcel.H
+++ b/src/lagrangian/intermediate/parcels/Templates/ReactingParcel/ReactingParcel.H
@@ -139,7 +139,6 @@ public:
             inline trackData
             (
                 ReactingCloud<ParcelType>& cloud,
-                const interpolation<scalar>& CpInterp,
                 typename ThermoParcel<ParcelType>::trackData::trackPart
                     part = ThermoParcel<ParcelType>::trackData::tpLinearTrack
             );
diff --git a/src/lagrangian/intermediate/parcels/Templates/ReactingParcel/ReactingParcelI.H b/src/lagrangian/intermediate/parcels/Templates/ReactingParcel/ReactingParcelI.H
index 285a1af0dad..4206975a70c 100644
--- a/src/lagrangian/intermediate/parcels/Templates/ReactingParcel/ReactingParcelI.H
+++ b/src/lagrangian/intermediate/parcels/Templates/ReactingParcel/ReactingParcelI.H
@@ -43,11 +43,10 @@ template<class ParcelType>
 inline Foam::ReactingParcel<ParcelType>::trackData::trackData
 (
     ReactingCloud<ParcelType>& cloud,
-    const interpolation<scalar>& CpInterp,
     typename ThermoParcel<ParcelType>::trackData::trackPart part
 )
 :
-    ThermoParcel<ParcelType>::trackData(cloud, CpInterp, part),
+    ThermoParcel<ParcelType>::trackData(cloud, part),
     cloud_(cloud),
     constProps_(cloud.constProps()),
     pInterp_
diff --git a/src/lagrangian/intermediate/parcels/Templates/ThermoParcel/ThermoParcel.H b/src/lagrangian/intermediate/parcels/Templates/ThermoParcel/ThermoParcel.H
index 38de73df1bc..2924dc3eb0a 100644
--- a/src/lagrangian/intermediate/parcels/Templates/ThermoParcel/ThermoParcel.H
+++ b/src/lagrangian/intermediate/parcels/Templates/ThermoParcel/ThermoParcel.H
@@ -141,6 +141,11 @@ public:
             //- Particle constant properties
             const constantProperties& constProps_;
 
+            //- Local copy of specific heat field
+            //  Cp not stored on acrrier thermo, but returned as tmp<...>
+            const volScalarField Cp_;
+
+
             // Interpolators for continuous phase fields
 
                 //- Temperature field interpolator
@@ -158,7 +163,6 @@ public:
             inline trackData
             (
                 ThermoCloud<ParcelType>& cloud,
-                const interpolation<scalar>& CpInterp,
                 typename KinematicParcel<ParcelType>::trackData::trackPart
                     part = KinematicParcel<ParcelType>::trackData::tpLinearTrack
             );
diff --git a/src/lagrangian/intermediate/parcels/Templates/ThermoParcel/ThermoParcelI.H b/src/lagrangian/intermediate/parcels/Templates/ThermoParcel/ThermoParcelI.H
index a7b9b9d68b9..69ba4e61094 100644
--- a/src/lagrangian/intermediate/parcels/Templates/ThermoParcel/ThermoParcelI.H
+++ b/src/lagrangian/intermediate/parcels/Templates/ThermoParcel/ThermoParcelI.H
@@ -45,13 +45,24 @@ template<class ParcelType>
 inline Foam::ThermoParcel<ParcelType>::trackData::trackData
 (
     ThermoCloud<ParcelType>& cloud,
-    const interpolation<scalar>& CpInterp,
     typename KinematicParcel<ParcelType>::trackData::trackPart part
 )
 :
     KinematicParcel<ParcelType>::trackData(cloud, part),
     cloud_(cloud),
     constProps_(cloud.constProps()),
+    Cp_
+    (
+        IOobject
+        (
+            "Cp",
+            cloud.db().time().timeName(),
+            cloud.db(),
+            IOobject::NO_READ,
+            IOobject::NO_WRITE
+        ),
+        cloud.thermo().thermo().Cp()
+    ),
     TInterp_
     (
         interpolation<scalar>::New
@@ -60,7 +71,14 @@ inline Foam::ThermoParcel<ParcelType>::trackData::trackData
             cloud.T()
         )
     ),
-    CpInterp_(CpInterp)
+    CpInterp_
+    (
+        interpolation<scalar>::New
+        (
+            cloud.solution().interpolationSchemes(),
+            Cp_
+        )
+    )
 {}
 
 
-- 
GitLab