CFD Online Logo CFD Online URL
www.cfd-online.com
[Sponsors]
Home > Forums > Software User Forums > OpenFOAM > OpenFOAM Programming & Development

Study of the EEqn.H in rhoPimpleDyMFoam.

Register Blogs Community New Posts Updated Threads Search

Like Tree6Likes
  • 5 Post By Horacio Aguerre
  • 1 Post By Horacio Aguerre

 
 
LinkBack Thread Tools Search this Thread Display Modes
Prev Previous Post   Next Post Next
Old   May 27, 2013, 08:03
Default Study of the EEqn.H in rhoPimpleDyMFoam.
  #1
New Member
 
Join Date: May 2013
Location: Santa Fé, Argentina
Posts: 5
Rep Power: 13
Horacio Aguerre is on a distinguished road
Hi,

I'm working with the solver "rhoPimpleDyMFoam" from the OpenFOAM 2-2-x version. In specific i was trying to replicate the polytropic process in an adiabatic compression.
The case i studied is a pistone-cylinder sistem where air is compressed by the pistone patch which is a moving wall. The pistone velocity is adjusted to 0,1(m/s) to aproximate a reversible process.

When i choose energy "he" as internal energy i could reproduce the polytropic process in an adiabatic compression.

However when i choose energy "he" as entalphy I couldn't achieve the theoretical expression pV^{1,4}=constant for air. Instead the system looks as it would be non adiabatic or like the conversion between work and internal energy was not made properly.

Therefore i started to study the EEqn.H.

I compare the theoretical balance of total energy in terms of entalphy with what is written in the EEqn.H.

Balance of a tensorial property in an arbritary moving domain.

\dfrac{d}{dt} \int_{V}\rho \phi dV + \int_{S}\rho \phi\left(\mathbf{v}-\mathbf{v_{b}}\right) \cdot \mathbf{ds}=-\int_{S}\rho \mathbf{Q_{\phi}} \cdot \mathbf{ds}+\int_{V}S_{\phi}dV

Balance of mechanical energy.

\phi \equiv \frac{\mathbf{v} \cdot \mathbf{v}}{2} = K

\dfrac{d}{dt} \int_{V}\rho K dV + \int_{S}\rho K \left(\mathbf{v}-\mathbf{v_{b}}\right) \cdot \mathbf{ds} =-\int_{S}\left(p\mathbf{v}+ \boldsymbol{\tau} \cdot \mathbf{v}\right) \cdot \mathbf{ds} + \int_{V} \left[p \left(\boldsymbol{\nabla} \cdot \mathbf{v}\right)+\boldsymbol{\tau}:\boldsymbol{\nabla v}\right] dV+ \int_{V}\rho\left(\mathbf{v} \cdot \mathbf{g}\right) dV


Balance of internal energy.

\phi \equiv \widehat{u}


\dfrac{d}{dt} \int_{V}\rho \widehat{u} dV + \int_{S}\rho \widehat{u} \left(\mathbf{v}-\mathbf{v_{b}}\right) \cdot \mathbf{ds}=-\int_{S} \mathbf{q} \cdot \mathbf{ds}-\int_{V} \left[p \left(\boldsymbol{\nabla} \cdot \mathbf{v}\right)+\boldsymbol{\tau}:\boldsymbol{\nabla v}\right] dV


The balance of internal energy can be expressed in terms of entalphy.

\widehat{h}=\widehat{u}+\frac{p}{\rho}\implies \widehat{u}=\widehat{h}-\frac{p}{\rho}

Replacing in the balance of internal energy.

\dfrac{d}{dt} \int_{V}\rho \left(\widehat{h}-\frac{p}{\rho}\right) dV + \int_{S}\rho \left(\widehat{h}-\frac{p}{\rho}\right) \left(\mathbf{v}-\mathbf{v_{b}}\right) \cdot \mathbf{ds}=-\int_{S} \mathbf{q} \cdot \mathbf{ds}-\int_{V} \left[p \left(\boldsymbol{\nabla} \cdot \mathbf{v}\right)+\boldsymbol{\tau}:\boldsymbol{\nabla v}\right] dV

\dfrac{d}{dt} \int_{V}\rho \widehat{h} dV+ \int_{S}\rho \widehat{h} \left(\mathbf{v}-\mathbf{v_{b}}\right) \cdot \mathbf{ds}-\dfrac{d}{dt} \int_{V}p dV - \int_{S} p \left(\mathbf{v}-\mathbf{v_{b}}\right) \cdot \mathbf{ds}=-\int_{S} \mathbf{q} \cdot \mathbf{ds}-\int_{V} \left[p \left(\boldsymbol{\nabla} \cdot \mathbf{v}\right)+\boldsymbol{\tau}:\boldsymbol{\nabla v}\right] dV


Balance of total energy.

The sum of the balances of mechanical and internal energy results in the balance of total energy.

\dfrac{d}{dt} \int_{V}\rho \left(\widehat{h}+ K \right) dV + \int_{S}\rho \left(\widehat{h}+ K\right) \left(\mathbf{v}-\mathbf{v_{b}}\right) \cdot \mathbf{ds} -\dfrac{d}{dt} \int_{V}p dV - \int_{S} p \left(\mathbf{v}-\mathbf{v_{b}}\right) \cdot \mathbf{ds}=-\int_{S} \mathbf{q} \cdot \mathbf{ds}-\int_{S}\left(p\mathbf{v}+ \boldsymbol{\tau} \cdot \mathbf{v}\right) \cdot \mathbf{ds} + \int_{V}\rho\left(\mathbf{v} \cdot \mathbf{g}\right) dV

If the terms of work of viscous forces and the gravity are neglected, the balance is:

\dfrac{d}{dt} \int_{V}\rho \left(\widehat{h}+ K \right) dV + \int_{S}\rho \left(\widehat{h}+ K\right) \left(\mathbf{v}-\mathbf{v_{b}}\right) \cdot \mathbf{ds} -\dfrac{d}{dt} \int_{V}p dV + \int_{S} p \mathbf{v_{b}} \cdot \mathbf{ds}=
-\int_{S} \mathbf{q} \cdot \mathbf{ds}

rhoPimpleDyMFoam: EEqn.H

The balance of total energy in the EEqn.H file with he considered as entalphy \widehat{h} is:
Code:
fvm::ddt(rho, he) + fvm::div(phi, he) + fvc::ddt(rho, K) + fvc::div(phi, K)
- dpdt - fvm::laplacian(turbulence ->alphaEff(), he) == fvOptions(rho, he)
Where, \phi=\rho_{f} \left(\mathbf{v}-\mathbf{v_{b}}\right)_{f} \cdot \mathbf{S_{f}}, is relative to the mesh.

If fvOptions(rho, he) are considered as zero, the term \int_{S} p \mathbf{v_{b}} \cdot \mathbf{ds} from the final balance is no considered in EEqn.H.

Then i try adding in the EEqn.H the therm \int_{S} p \mathbf{v_{b}} \cdot \mathbf{ds}.

Code:
fvm::ddt(rho, he) + fvm::div(phi, he) + fvc::ddt(rho, K) + fvc::div(phi, K)
- dpdt + fvc::div(fvc::meshPhi(U),p) - fvm::laplacian(turbulence ->alphaEff(), he) == fvOptions(rho, he)
The result was correct. As in the first case (he=internal energy) the adiabatic compression was reproduced with the relation pV^{1,4}=constant.

The question is,

should the term \int_{S} p \mathbf{v_{b}}  \cdot \mathbf{ds} be introduced in any way in fvOptions(rho, he) or it must be explicit in the EEqn.H?

Thanks,

Horacio
mturcios777, fredo490, mgg and 2 others like this.
Horacio Aguerre is offline   Reply With Quote

 

Tags
rhopimpledymfoam


Posting Rules
You may not post new threads
You may not post replies
You may not post attachments
You may not edit your posts

BB code is On
Smilies are On
[IMG] code is On
HTML code is Off
Trackbacks are Off
Pingbacks are On
Refbacks are On


Similar Threads
Thread Thread Starter Forum Replies Last Post
grid study, important question hamid1 FLUENT 1 August 4, 2013 00:14
[ANSYS Meshing] grid study, important questions hamid1 ANSYS Meshing & Geometry 2 February 10, 2012 13:28
Loop outlet data to inlet of new study pggow FLUENT 0 December 21, 2011 14:31
dpm study of solid-liquid study. bhanwar Main CFD Forum 0 September 29, 2011 08:21
Study Mario Main CFD Forum 0 August 24, 2009 09:59


All times are GMT -4. The time now is 21:06.