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

Problems with time derivative

Register Blogs Community New Posts Updated Threads Search

Like Tree8Likes
  • 1 Post By ahmmedshakil
  • 2 Post By alexeym
  • 1 Post By demichie
  • 1 Post By alexeym
  • 1 Post By Franko
  • 2 Post By alexeym

Reply
 
LinkBack Thread Tools Search this Thread Display Modes
Old   January 25, 2015, 18:28
Default Problems with time derivative
  #1
New Member
 
Join Date: Sep 2014
Posts: 15
Rep Power: 11
Franko is on a distinguished road
Hello foamers!
I have two equations (assume p and T are dimensionless volScalarFields)

1) \frac{\partial p}{\partial t} + p = 0;

2) \frac{\partial p T}{\partial t} + T = 0;

I already know T_i^n and p_i^n, i = 0,1,...,I

The first equation is easy:

Code:
fvScalarMatrix pEqn
(   
fvm::ddt(p) + p
);
pEqn.solve();
As you know, this is equivalent to

\frac{p_i^{n+1} - p_i^n}{\tau} + p_i^n = 0

Now I know p_i^{n+1} for all i=0,1,...,I


If I write the second equation using finite differences, I will get

\frac{p_i^{n+1}T_i^{n+1} - p_i^nT_i^n}{\tau} + T_i^n = 0

But what is the analogue of this equation in OpenFOAM?

If I try to use
Code:
fvm::ddt(p,T) + T
I will get the wrong equation:

\frac{p_i^{n+1}T_i^{n+1} - p_i^{n+1}T_i^n}{\tau} + T_i^n = 0

So how to use OpenFOAM correctly here?
Please help!
Franko is offline   Reply With Quote

Old   January 26, 2015, 02:39
Default
  #2
Senior Member
 
Alexey Matveichev
Join Date: Aug 2011
Location: Nancy, France
Posts: 1,930
Rep Power: 38
alexeym has a spectacular aura aboutalexeym has a spectacular aura about
Send a message via Skype™ to alexeym
Hi,

If I am not mistaken

\frac{\partial{pT}}{\partial{t}} = p\frac{\partial{T}}{\partial{t}} + T\frac{\partial{p}}{\partial{t}}

and then using first equation:

\frac{\partial{pT}}{\partial{t}} = p\frac{\partial{T}}{\partial{t}} - pT

and

p\frac{\partial{T}}{\partial{t}}

is

Code:
fvm::ddt(p, T)
alexeym is offline   Reply With Quote

Old   January 26, 2015, 03:18
Default
  #3
Member
 
Mattia de\' Michieli Vitturi
Join Date: Mar 2009
Posts: 50
Rep Power: 17
demichie is on a distinguished road
Hi,
are you sure that about this?

Quote:
Originally Posted by alexeym View Post
p\frac{\partial{T}}{\partial{t}}

is

Code:
fvm::ddt(p, T)

I've always tought that ddt(p,T) discretizes the derivative in time of pT. Otherwise, if it was p times the temporal derivative of T, why not using p*ddt(t)?
Also looking at Table 2.2 of the programming guide it seems that ddt(p,T) should be the discretization of the time derivative of (p*T).

Can you clarify it?

Thank you
Mattia
demichie is offline   Reply With Quote

Old   January 26, 2015, 03:55
Default
  #4
Senior Member
 
Mohammad Shakil Ahmmed
Join Date: Oct 2012
Location: AUS
Posts: 137
Rep Power: 14
ahmmedshakil is on a distinguished road
Quote:
Originally Posted by demichie View Post
Hi,
are you sure that about this?




I've always tought that ddt(p,T) discretizes the derivative in time of pT. Otherwise, if it was p times the temporal derivative of T, why not using p*ddt(t)?
Also looking at Table 2.2 of the programming guide it seems that ddt(p,T) should be the discretization of the time derivative of (p*T).

Can you clarify it?

Thank you
Mattia
Yes, I think Mattia de\' Michieli Vitturi is right. To the best knowledge, it's discretised as time derivative of (p*T). To be confirmed, I guess, you can check the book "An Introduction to Computational Fluid Dynamics: The Finite Volume Method"-by H.K. Versteeg and W. Malalasekera.

Cheers
shakil
emjay likes this.
ahmmedshakil is offline   Reply With Quote

Old   January 26, 2015, 04:13
Default
  #5
New Member
 
Join Date: Sep 2014
Posts: 15
Rep Power: 11
Franko is on a distinguished road
So what is the correct formula:
1)fvm::ddt(p,T)=\frac{p_i^{n+1}T_i^{n+1} - p_i^nT_i^n}{\tau}

or 2) fvm::ddt(p,T)=\frac{p_i^{n+1}T_i^{n+1} - p_i^{n+1}T_i^n}{\tau}?

If you think that the first formula is right, are you sure that OpenFOAM "remember" the old p (p_i^n, i = 1,2,...,I)?
Franko is offline   Reply With Quote

Old   January 26, 2015, 04:18
Default
  #6
Senior Member
 
Alexey Matveichev
Join Date: Aug 2011
Location: Nancy, France
Posts: 1,930
Rep Power: 38
alexeym has a spectacular aura aboutalexeym has a spectacular aura about
Send a message via Skype™ to alexeym
Well,

Here fvc::ddt method in EulerDdtScheme:

Code:
template<class Type>
tmp<GeometricField<Type, fvPatchField, volMesh> >
EulerDdtScheme<Type>::fvcDdt
(
    const dimensionedScalar& rho,
    const GeometricField<Type, fvPatchField, volMesh>& vf
)
{
    dimensionedScalar rDeltaT = 1.0/mesh().time().deltaT();

    IOobject ddtIOobject
    (
        "ddt("+rho.name()+','+vf.name()+')',
        mesh().time().timeName(),
        mesh()
    );
    if (mesh().moving())
    {
        ...
    }
    else
    {
        return tmp<GeometricField<Type, fvPatchField, volMesh> >
        (
            new GeometricField<Type, fvPatchField, volMesh>
            (
                ddtIOobject,
                rDeltaT*rho*(vf - vf.oldTime())
            )
        );
    }
}
While fvm::ddt

Code:
template<class Type>
tmp<fvMatrix<Type> >
EulerDdtScheme<Type>::fvmDdt
(
    const volScalarField& rho,
    const GeometricField<Type, fvPatchField, volMesh>& vf
)
{
    ...
    scalar rDeltaT = 1.0/mesh().time().deltaTValue();

    fvm.diag() = rDeltaT*rho.internalField()*mesh().Vsc();
    if (mesh().moving())
    {
        ...
    }
    else
    {
        fvm.source() = rDeltaT
            *rho.oldTime().internalField()
            *vf.oldTime().internalField()*mesh().Vsc();
    }
    ...
}
demichie and Franko like this.
alexeym is offline   Reply With Quote

Old   January 26, 2015, 04:18
Default
  #7
Member
 
Mattia de\' Michieli Vitturi
Join Date: Mar 2009
Posts: 50
Rep Power: 17
demichie is on a distinguished road
The correct formula is the first one. Take a look at the programmers guide (page 37 for version 2.3.1) and you will find everything. Eq. 2.21 is exactly what you are looking for.

Mattia
Franko likes this.
demichie is offline   Reply With Quote

Old   January 27, 2015, 18:35
Default
  #8
New Member
 
Join Date: Sep 2014
Posts: 15
Rep Power: 11
Franko is on a distinguished road
I am not sure that I understand fvc::ddt correctly...

I think that in the beginning of the time loop (p and T are dimensionless volScalarFields) we have formula below:

fvc::ddt(p,T) = \frac{p_i^n*(T_i^{n+1} - T_i^n)}{\tau}

But after solving this equation:

Code:
fvScalarMatrix pEqn
(  
fvm::ddt(p) + p
);
pEqn.solve();
fvc::ddt(p,T) = \frac{p_i^{n+1}*(T_i^{n+1} - T_i^n)}{\tau}

So fvc::ddt(rho, phi) = \frac{rho*(phi^{n+1} - phi^n)}{\tau} ?
I mean when we use fvc::ddt(rho, phi), rho always behaves like a scalar, right?
Please correct me, if I am wrong.
Franko is offline   Reply With Quote

Old   January 28, 2015, 03:26
Default
  #9
Senior Member
 
Alexey Matveichev
Join Date: Aug 2011
Location: Nancy, France
Posts: 1,930
Rep Power: 38
alexeym has a spectacular aura aboutalexeym has a spectacular aura about
Send a message via Skype™ to alexeym
It behaves like scalar field i.e. it can vary in space.
Franko likes this.
alexeym is offline   Reply With Quote

Old   January 28, 2015, 04:10
Default
  #10
New Member
 
Join Date: Sep 2014
Posts: 15
Rep Power: 11
Franko is on a distinguished road
alexeym, yes, sure.

And if we have only one argument (volScalarField phi) we will get fvc::ddt(phi) = \frac{phi^{n+1} - phi^n}{\tau}

Now it is clear.
Franko is offline   Reply With Quote

Old   January 28, 2015, 16:44
Default new problem
  #11
New Member
 
Join Date: Sep 2014
Posts: 15
Rep Power: 11
Franko is on a distinguished road
Hi, I have a new big problem with temporal derivative!

Again, A and B are dimensionless volScalarFields
I also have
Code:
dimensionedScalar c ("c", dimensionSet(0, 0, 1, 0, 0, 0 ,0), 1.0);
In the beginning, I know A_i^n and B_i^n, i = 0,1,...,I

I solved the first equation:

Code:
fvScalarMatrix AEqn
        (
        c*fvm::ddt(A) + A
        );
AEqn.solve();
As you know, this is equivalent to

\frac{A_i^{n+1} - A_i^n}{\tau} + A_i^n = 0

Now I know A_i^{n+1}, A_i^n and B_i^n, i=0,1,...,I

After that I tried to solve the second equation:

Code:
fvScalarMatrix BEqn
        (
            c*fvm::ddt(A) +c*fvm::ddt(B) + B
        );
BEqn.solve();
It compiles, but then

Code:
--> FOAM FATAL ERROR : incompatible fields
for operation [A] + [B]
It seems that OpenFOAM don't understand that I have already calculated field A_i^{n+1} in the first equation!
How can I tell him that?

Thank you very much for any hint ...
Annier likes this.
Franko is offline   Reply With Quote

Old   January 29, 2015, 01:45
Default
  #12
Senior Member
 
Alexey Matveichev
Join Date: Aug 2011
Location: Nancy, France
Posts: 1,930
Rep Power: 38
alexeym has a spectacular aura aboutalexeym has a spectacular aura about
Send a message via Skype™ to alexeym
Hi,

fvm:: methods usually operate on the same field, so OpenFOAM is not quite happy about your fvm::ddt(A) + fvm::ddt(B). If you need time derivative of A in equation for B, you can use it as a source term with fvc::ddt(A). So your code should be:

Code:
fvm::ddt(B) + B/c + fvc::ddt(A)
(if c is scalar != 0, why not divide equations by c?)
Annier and Franko like this.
alexeym is offline   Reply With Quote

Old   January 29, 2015, 02:53
Default
  #13
New Member
 
Join Date: Sep 2014
Posts: 15
Rep Power: 11
Franko is on a distinguished road
alexeym, exactly, I need to use fvc::ddt(A) in my equation...
Thanks a lot!
Franko is offline   Reply With Quote

Reply


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
Elastic or Plastic Deformations in OpenFOAM NickolasPl OpenFOAM Pre-Processing 36 September 23, 2023 08:22
AMI speed performance danny123 OpenFOAM 21 October 24, 2020 04:13
simpleFoam error - "Floating point exception" mbcx4jc2 OpenFOAM Running, Solving & CFD 12 August 4, 2015 02:20
Moving mesh Niklas Wikstrom (Wikstrom) OpenFOAM Running, Solving & CFD 122 June 15, 2014 06:20
Micro Scale Pore, icoFoam gooya_kabir OpenFOAM Running, Solving & CFD 2 November 2, 2013 13:58


All times are GMT -4. The time now is 18:38.