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

Read temperature dependent thermophysical properties from a file - boundaries false

Register Blogs Community New Posts Updated Threads Search

Like Tree13Likes
  • 3 Post By AnjaMiehe
  • 10 Post By AnjaMiehe

Reply
 
LinkBack Thread Tools Search this Thread Display Modes
Old   June 26, 2012, 11:38
Default [Solved] Read temperature dependent properties from a file and interpolate
  #1
Member
 
Anja Miehe
Join Date: Dec 2009
Location: Freiberg / Germany
Posts: 48
Rep Power: 16
AnjaMiehe is on a distinguished road
Hello everyone,

I am trying to implement the following:
- the thermophysical property "DT" shall be read from a data file and interpolated according to the temperature in the solution domain.

Reading some posts here on cfd-online.com I got quite far, only the bounary values of the "DT" field are not calculated.

#######
The Solver - based on icoFoam with TEqn according to wiki, but "DT" as volScalarField
#######

The icoFileFoam.C file
Code:
/*---------------------------------------------------------------------------*\
  =========                 |
  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
   \\    /   O peration     |
    \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
     \\/     M anipulation  |
-------------------------------------------------------------------------------
License
    This file is part of OpenFOAM.

    OpenFOAM is free software: you can redistribute it and/or modify it
    under the terms of the GNU General Public License as published by
    the Free Software Foundation, either version 3 of the License, or
    (at your option) any later version.

    OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
    ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
    FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
    for more details.

    You should have received a copy of the GNU General Public License
    along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.

Application
    icoFoam

Description
    Transient solver for incompressible, laminar flow of Newtonian fluids.

\*---------------------------------------------------------------------------*/

#include "fvCFD.H"
#include "IFstream.H"
#include "graph.H"
#include "interpolateXY.H"

// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //

int main(int argc, char *argv[])
{
    #include "setRootCase.H"

    #include "createTime.H"
    #include "createMesh.H"
    #include "createFields.H"
    #include "initContinuityErrs.H"
    #include "interpolateProperties.H"

    // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //

    Info<< "\nStarting time loop\n" << endl;

    while (runTime.loop())
    {
        Info<< "Time = " << runTime.timeName() << nl << endl;

        #include "readPISOControls.H"
        #include "CourantNo.H"

        fvVectorMatrix UEqn
        (
            fvm::ddt(U)
          + fvm::div(phi, U)
          - fvm::laplacian(nu, U)
        );

        solve(UEqn == -fvc::grad(p));

        // --- PISO loop

        for (int corr=0; corr<nCorr; corr++)
        {
            volScalarField rAU(1.0/UEqn.A());

            U = rAU*UEqn.H();
            phi = (fvc::interpolate(U) & mesh.Sf())
                + fvc::ddtPhiCorr(rAU, U, phi);

            adjustPhi(phi, U, p);

            for (int nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++)
            {
                fvScalarMatrix pEqn
                (
                    fvm::laplacian(rAU, p) == fvc::div(phi)
                );

                pEqn.setReference(pRefCell, pRefValue);
                pEqn.solve();

                if (nonOrth == nNonOrthCorr)
                {
                    phi -= pEqn.flux();
                }
            }

            #include "continuityErrs.H"

            U -= rAU*fvc::grad(p);
            U.correctBoundaryConditions();
        }

        #include "interpolateProperties.H"

        fvScalarMatrix TEqn
        (
            fvm::ddt(T)
            + fvm::div(phi, T)
            - fvm::laplacian(DT, T)
        );

        TEqn.solve();

        runTime.write();

        Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
            << "  ClockTime = " << runTime.elapsedClockTime() << " s"
            << nl << endl;
    }

    Info<< "End\n" << endl;

    return 0;
}


// ************************************************************************* //
The createFields.H
Code:
    Info<< "Reading transportProperties\n" << endl;

    IOdictionary transportProperties
    (
        IOobject
        (
            "transportProperties",
            runTime.constant(),
            mesh,
            IOobject::MUST_READ_IF_MODIFIED,
            IOobject::NO_WRITE
        )
    );

    dimensionedScalar nu
    (
        transportProperties.lookup("nu")
    );

    Info<< "Reading field p\n" << endl;
    volScalarField p
    (
        IOobject
        (
            "p",
            runTime.timeName(),
            mesh,
            IOobject::MUST_READ,
            IOobject::AUTO_WRITE
        ),
        mesh
    );


    Info<< "Reading field U\n" << endl;
    volVectorField U
    (
        IOobject
        (
            "U",
            runTime.timeName(),
            mesh,
            IOobject::MUST_READ,
            IOobject::AUTO_WRITE
        ),
        mesh
    );

    Info<< "Reading field T\n" <<endl;
    volScalarField T
    (
        IOobject
        (
             "T",
             runTime.timeName(),
             mesh,
             IOobject::MUST_READ,
             IOobject::AUTO_WRITE
         ),
         mesh
    );

    volScalarField DT
    (
        IOobject
        (
             "DT",
             runTime.timeName(),
             mesh,
             IOobject::NO_READ,
             IOobject::AUTO_WRITE
         ),
         mesh,
         dimensionedScalar ("DT",dimensionSet (0,2,-1,0,0,0,0), 1e-3) // this is just for initializing 
    );


#   include "createPhi.H"


    label pRefCell = 0;
    scalar pRefValue = 0.0;
    setRefCell(p, mesh.solutionDict().subDict("PISO"), pRefCell, pRefValue);
the interpolateProperties.H
Code:
    Info << "Reading DT property file and interpolate" << endl; 
    IFstream file_DT(runTime.path()/runTime.constant()/"DT"); 
    
    //Read file
    graph DTGraph
    (
        "DT_data_file",
        "T_data",
        "DT_data",
        file_DT
    ); 
    
    //Interpolation
    DT.field() = interpolateXY
    (
        T.field(),
        DTGraph.x(),
        DTGraph.y()
    );
The complete code is attached as zip.

######
The test case
######
It is conduction only, a block 1 m long, 0.1 m high and 0.01 m in width with 100x10x1 cells. Every boundary is a wall (or empty), no velocity and the temperature is set to 273.15 K except for the right wall being 373.15 K.
I used funkySetFields -time 0 -field T -keepPatches -expression "pos().x*50.0+273.15" to see the effect of interpolated "DT" field. As coded, the file "DT" is in the folder "constant".

The test case is attached as zip, too.

All this works fine. The "DT" field is written to the time directory as demanded and interpolated. Only, the boundary values of the "DT" stay 1e-3 as given in the createFields.H for initialisation (that's also why I chose this odd value). Therefore, the temperature does not build up to a linear slope from left to right, from 273.15 K to 373.15 K as it should. The written-out "DT" file gives sensible values in the internal field, but only "calculated 1e-3" for every boundary.

This means, the boundary values of "DT" are not updated. Can anyone tell me, what line of coding is missing to interpolated these according to temperature, too ? Or maybe, anyone knows that the error is somewhere else?

Thanks in advance
Anja
Attached Files
File Type: zip icoFileFoam.zip (3.2 KB, 222 views)
File Type: zip icoFileFoam_testcase.zip (46.0 KB, 248 views)

Last edited by AnjaMiehe; June 27, 2012 at 08:14.
AnjaMiehe is offline   Reply With Quote

Old   June 27, 2012, 08:06
Default
  #2
Member
 
Anja Miehe
Join Date: Dec 2009
Location: Freiberg / Germany
Posts: 48
Rep Power: 16
AnjaMiehe is on a distinguished road
I got it.
The interpolateProperties.H has to look as follows to account for the boundary values as well:
Code:
    IFstream file_DT(runTime.path()/runTime.constant()/"DT");
    //Read file
    graph DTGraph
    (
        "DT_data_file",
        "T_data",
        "DT_data",
        file_DT
    ); 
    
    //Interpolation internal field
    DT.field() = interpolateXY
    (
        T.field(),
        DTGraph.x(),
        DTGraph.y()
    );

    //Interpolation boundary field
    forAll(DT.boundaryField(), patchi)
    {
        DT.boundaryField()[patchi] = interpolateXY
        (
            T.boundaryField()[patchi],
            DTGraph.x(),
            DTGraph.y()
        );
    }
Now it works, have fun using it.
Anja
AnjaMiehe is offline   Reply With Quote

Old   October 1, 2013, 01:52
Default
  #3
New Member
 
eric
Join Date: Nov 2010
Location: Vancouver, Canada
Posts: 16
Rep Power: 15
armyou is on a distinguished road
Thank you for sharing. Anja
armyou is offline   Reply With Quote

Old   November 10, 2016, 18:55
Default Errors in parallel
  #4
New Member
 
Paul
Join Date: May 2012
Posts: 23
Rep Power: 13
pmdelgado2 is on a distinguished road
When I run this interpolateProperties.H file in parallel, I get the following error (OF2.3.1):

error in IOstream "/path/to/case/directory/processor3/constant/DT" for operation operator >>(Istream&,List<T>&) from function IOstream::fatalCheck(const char*) .

Has anyone else tried this code in parallel? Is there any way I can modify this code so that it runs in parallel?

Last edited by pmdelgado2; November 14, 2016 at 13:18.
pmdelgado2 is offline   Reply With Quote

Old   November 14, 2016, 13:17
Default SOLVED - Parallel Read temperature dependent thermophysical properties from a file
  #5
New Member
 
Paul
Join Date: May 2012
Posts: 23
Rep Power: 13
pmdelgado2 is on a distinguished road
I figured out what the problem was with doing this in parallel.

Apparently, when you run decomposePar, the solver will automatically look for the tabular files in the subfolders

/path/to/case/directory/processor0/constant/
/path/to/case/directory//processor1/constant/
etc..,

instead of just in the directory

/path/to/case/directory/constant/

To resolve this issue, you need to copy the tabular file "DT" from

/path/to/case/directory/constant/DT

to

/path/to/case/directory/processor0/constant/DT
/path/to/case/directory/processor1/constant/DT
etc...

It may be helpful to do this in a loop within your script.
pmdelgado2 is offline   Reply With Quote

Old   February 27, 2018, 07:51
Default
  #6
Member
 
Anil Kunwar
Join Date: Jun 2013
Posts: 64
Rep Power: 11
Annier is an unknown quantity at this point
Thanks for this discussion and it has proved quite helpful for me to define temperature dependent properties in OpenFOAM.
Annier is offline   Reply With Quote

Old   November 14, 2023, 02:42
Default
  #7
New Member
 
Zhao Zhongping
Join Date: Nov 2023
Posts: 1
Rep Power: 0
GooseEast is on a distinguished road
I recurrent it in OpenFOAM-v2306 and get the following error:
interpolateProperties.H:29:9: error: passing ‘const Foam::fvPatchField<double>’ as ‘this’ argument discards qualifiers [-fpermissive] 29 | ); | ^
which the error part is://Interpolation boundary field
forAll(DT.boundaryField(), patchi)
{
DT.boundaryField()[patchi] = interpolateXY
(
T.boundaryField()[patchi],
DTGraph.x(),
DTGraph.y()
);
}
Is there any way I can modify this part code so that it runs in OpenFOAM-v2306?
GooseEast 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
[Other] mesh airfoil NACA0012 anand_30 OpenFOAM Meshing & Mesh Conversion 13 March 7, 2022 17:22
wmake compiling new solver mksca OpenFOAM Programming & Development 14 June 22, 2018 06:29
OpenFOAM Install Script ljsh OpenFOAM Installation 82 October 12, 2009 11:47
Results saving in CFD hawk Main CFD Forum 16 July 21, 2005 20:51


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