diff --git a/include/cantera/oneD/StFlow.h b/include/cantera/oneD/StFlow.h index cd995e5e3c3..9ec91b00826 100644 --- a/include/cantera/oneD/StFlow.h +++ b/include/cantera/oneD/StFlow.h @@ -9,6 +9,8 @@ #include "Domain1D.h" #include "cantera/base/Array.h" #include "cantera/thermo/IdealGasPhase.h" +#include "cantera/thermo/RedlichKwongMFTP.h" +#include "cantera/thermo/PengRobinson.h" #include "cantera/kinetics/Kinetics.h" namespace Cantera @@ -71,7 +73,7 @@ class StFlow : public Domain1D * Set the thermo manager. Note that the flow equations assume * the ideal gas equation. */ - void setThermo(IdealGasPhase& th) { + void setThermo(ThermoPhase& th) { m_thermo = &th; } @@ -393,6 +395,10 @@ class StFlow : public Domain1D //! Update the diffusive mass fluxes. virtual void updateDiffFluxes(const doublereal* x, size_t j0, size_t j1); + //Get the gradient of speies specific molar enthalpies + virtual vector_fp grad_hk(const doublereal* x, size_t j); + virtual void updateMolarEnthalpies(const doublereal* x, size_t j); + //--------------------------------------------------------- // member data //--------------------------------------------------------- @@ -423,7 +429,7 @@ class StFlow : public Domain1D size_t m_nsp; - IdealGasPhase* m_thermo; + ThermoPhase* m_thermo; Kinetics* m_kin; Transport* m_trans; @@ -471,6 +477,11 @@ class StFlow : public Domain1D //! Temperature at the point used to fix the flame location double m_tfixed; + //Pointers to save molar enthalpies + vector_fp m_hk_current; + vector_fp m_hk_right; + vector_fp m_hk_left; + private: vector_fp m_ybar; }; diff --git a/samples/cxx/flamespeed/flamespeed.cpp b/samples/cxx/flamespeed/flamespeed.cpp index 9f1be1c0b0e..5cd66fbad26 100644 --- a/samples/cxx/flamespeed/flamespeed.cpp +++ b/samples/cxx/flamespeed/flamespeed.cpp @@ -16,7 +16,7 @@ #include "cantera/oneD/Sim1D.h" #include "cantera/oneD/Boundary1D.h" #include "cantera/oneD/StFlow.h" -#include "cantera/thermo/IdealGasPhase.h" +#include "cantera/thermo/ThermoPhase.h" #include "cantera/transport.h" #include "cantera/base/Solution.h" #include "cantera/base/stringUtils.h" diff --git a/src/oneD/StFlow.cpp b/src/oneD/StFlow.cpp index e0635dc302c..158a645e1b4 100644 --- a/src/oneD/StFlow.cpp +++ b/src/oneD/StFlow.cpp @@ -364,6 +364,7 @@ void StFlow::evalResidual(double* x, double* rsd, int* diag, } } + vector_fp h_molar(m_nsp); for (size_t j = jmin; j <= jmax; j++) { //---------------------------------------------- // left boundary @@ -451,22 +452,22 @@ void StFlow::evalResidual(double* x, double* rsd, int* diag, // - sum_k(J_k c_p_k / M_k) dT/dz //----------------------------------------------- if (m_do_energy[j]) { + // Update vectors m_hk_current, m_hk_left and m_hk_right + updateMolarEnthalpies(x, j); setGas(x,j); // heat release term - const vector_fp& h_RT = m_thermo->enthalpy_RT_ref(); - const vector_fp& cp_R = m_thermo->cp_R_ref(); + vector_fp dHk_dz = grad_hk(x, j); + m_thermo->getPartialMolarEnthalpies(&h_molar[0]); + double sum = 0.0; double sum2 = 0.0; for (size_t k = 0; k < m_nsp; k++) { double flxk = 0.5*(m_flux(k,j-1) + m_flux(k,j)); - sum += wdot(k,j)*h_RT[k]; - sum2 += flxk*cp_R[k]/m_wt[k]; + sum += wdot(k,j)*h_molar[k]; + sum2 += flxk * dHk_dz[k] / m_wt[k]; } - sum *= GasConstant * T(x,j); double dtdzj = dTdz(x,j); - sum2 *= GasConstant * dtdzj; - rsd[index(c_offset_T, j)] = - m_cp[j]*rho_u(x,j)*dtdzj - divHeatFlux(x,j) - sum - sum2; rsd[index(c_offset_T, j)] /= (m_rho[j]*m_cp[j]); @@ -898,4 +899,41 @@ void StFlow::evalContinuity(size_t j, double* x, double* rsd, int* diag, double } } +vector_fp StFlow::grad_hk(const doublereal* x, size_t j) +{ + vector_fp dhk_dz(m_nsp, 0.0); + + for(size_t k = 0; k < m_nsp; k++) + { + if (u(x, j) > 0.0) + { + dhk_dz[k] = (m_hk_current[k] - m_hk_left[k]) / m_dz[j - 1]; + } + else + { + dhk_dz[k] = (m_hk_right[k] - m_hk_current[k]) / m_dz[j]; + } + } + return dhk_dz; +} + +void StFlow::updateMolarEnthalpies(const doublereal* x, size_t j) +{ + m_hk_left.resize(m_nsp,0.0); + m_hk_current.resize(m_nsp,0.0); + m_hk_right.resize(m_nsp,0.0); + + // Node (j-1) + setGas(x, j-1); + m_thermo->getPartialMolarEnthalpies(&m_hk_left[0]); + + // Current Node j + setGas(x, j); + m_thermo->getPartialMolarEnthalpies(&m_hk_current[0]); + + // Node (j+1) + setGas(x, j + 1); + m_thermo->getPartialMolarEnthalpies(&m_hk_right[0]); +} + } // namespace