View Issue Details

IDProjectCategoryView StatusLast Update
0001861OpenFOAMBugpublic2015-10-05 11:35
Reporterwyldckat Assigned Tohenry  
PrioritylowSeverityfeatureReproducibilityN/A
Status resolvedResolutionfixed 
Product Versiondev 
Summary0001861: applications/utilities/postProcessing: Adding "-region" option to wallGradU, temporalInterpolate and pPrime2
DescriptionIt's a bit hard to fully determine which utilities can properly use the "-region" option, but the attached files provide the ones that seem that can make proper use of this option with minimum modifications:

 - applications/utilities/postProcessing/wall/wallGradU/wallGradU.C

 - applications/utilities/postProcessing/miscellaneous/temporalInterpolate/temporalInterpolate.C

 - applications/utilities/postProcessing/scalarField/pPrime2/pPrime2.C
Steps To ReproduceI've tested the 3 with the tutorial "heatTransfer/chtMultiRegionFoam/multiRegionHeater" for the region "bottomWater", with the following additional lines in "system/controlDict":

    libs
    (
    "libcompressibleTransportModels.so"
    "libfluidThermophysicalModels.so"
    "libsolidThermo.so"
    "libspecie.so"
    "libturbulenceModels.so"
    "libcompressibleTurbulenceModels.so"
    "libmeshTools.so"
    "libfiniteVolume.so"
    "libradiationModels.so"
    "libfvOptions.so"
    "libregionModels.so"
    "libsampling.so"
    );

    functions
    {
        fieldAverage1
        {
            type fieldAverage;
            functionObjectLibs ( "libfieldFunctionObjects.so" );
            enabled true;
            region bottomWater;
            outputControl outputTime;

            fields
            (
                U
                {
                    mean on;
                    prime2Mean on;
                    base time;
                }

                p
                {
                    mean on;
                    prime2Mean on;
                    base time;
                }
            );
        }
    }
Additional InformationThis is a follow up to Issue #1227.

The utilities stressComponents and createTurbulenceFields can also make good use of the option "-region", but for that they need to "-compressible" option as well (or inherit the mechanism from yPlus), which will be addressed in another report in the near future.
TagsNo tags attached.

Relationships

related to 0001227 resolvedhenry "ptot" utility does not support -region 

Activities

wyldckat

2015-10-03 17:50

updater  

wallGradU.C (3,339 bytes)   
/*---------------------------------------------------------------------------*\
  =========                 |
  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
   \\    /   O peration     |
    \\  /    A nd           | Copyright (C) 2011-2013 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
    wallGradU

Description
    Calculates and writes the gradient of U at the wall.

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

#include "fvCFD.H"
#include "wallFvPatch.H"

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

int main(int argc, char *argv[])
{
    timeSelector::addOptions();
    #include "addRegionOption.H"

    #include "setRootCase.H"
    #include "createTime.H"
    
    instantList timeDirs = timeSelector::select0(runTime, args);

    #include "createNamedMesh.H"

    forAll(timeDirs, timeI)
    {
        runTime.setTime(timeDirs[timeI], timeI);
        Info<< "Time = " << runTime.timeName() << endl;

        IOobject Uheader
        (
            "U",
            runTime.timeName(),
            mesh,
            IOobject::MUST_READ
        );

        // Check U exists
        if (Uheader.headerOk())
        {
            mesh.readUpdate();

            Info<< "    Reading U" << endl;
            volVectorField U(Uheader, mesh);

            Info<< "    Calculating wallGradU" << endl;

            volVectorField wallGradU
            (
                IOobject
                (
                    "wallGradU",
                    runTime.timeName(),
                    mesh,
                    IOobject::NO_READ,
                    IOobject::AUTO_WRITE
                ),
                mesh,
                dimensionedVector
                (
                    "wallGradU",
                    U.dimensions()/dimLength,
                    vector::zero
                )
            );

            const fvPatchList& patches = mesh.boundary();

            forAll(wallGradU.boundaryField(), patchi)
            {
                const fvPatch& currPatch = patches[patchi];

                if (isA<wallFvPatch>(currPatch))
                {
                    wallGradU.boundaryField()[patchi] =
                        -U.boundaryField()[patchi].snGrad();
                }
            }

            wallGradU.write();
        }
        else
        {
            Info<< "    No U" << endl;
        }
    }

    Info<< "End" << endl;

    return 0;
}


// ************************************************************************* //
wallGradU.C (3,339 bytes)   

wyldckat

2015-10-03 17:50

updater  

temporalInterpolate.C (8,660 bytes)   
/*---------------------------------------------------------------------------*\
  =========                 |
  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
   \\    /   O peration     |
    \\  /    A nd           | Copyright (C) 2011-2013 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/>.

Description
    Interpolate fields between time-steps e.g. for animation.

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

#include "argList.H"
#include "timeSelector.H"

#include "fvMesh.H"
#include "Time.H"
#include "volMesh.H"
#include "surfaceMesh.H"
#include "volFields.H"
#include "surfaceFields.H"
#include "pointFields.H"
#include "ReadFields.H"
#include "interpolationWeights.H"
#include "uniformInterpolate.H"

using namespace Foam;

class fieldInterpolator
{
    Time& runTime_;
    const fvMesh& mesh_;
    const IOobjectList& objects_;
    const HashSet<word>& selectedFields_;
    instant ti_;
    instant ti1_;
    const interpolationWeights& interpolator_;
    const wordList& timeNames_;
    int divisions_;

public:

    fieldInterpolator
    (
        Time& runTime,
        const fvMesh& mesh,
        const IOobjectList& objects,
        const HashSet<word>& selectedFields,
        const instant& ti,
        const instant& ti1,
        const interpolationWeights& interpolator,
        const wordList& timeNames,
        int divisions
    )
    :
        runTime_(runTime),
        mesh_(mesh),
        objects_(objects),
        selectedFields_(selectedFields),
        ti_(ti),
        ti1_(ti1),
        interpolator_(interpolator),
        timeNames_(timeNames),
        divisions_(divisions)
    {}

    template<class GeoFieldType>
    void interpolate();
};


template<class GeoFieldType>
void fieldInterpolator::interpolate()
{
    const word& fieldClassName = GeoFieldType::typeName;

    IOobjectList fields = objects_.lookupClass(fieldClassName);

    if (fields.size())
    {
        Info<< "    " << fieldClassName << "s:";

        forAllConstIter(IOobjectList, fields, fieldIter)
        {
            if
            (
                selectedFields_.empty()
             || selectedFields_.found(fieldIter()->name())
            )
            {
                Info<< " " << fieldIter()->name() << '(';

                scalar deltaT = (ti1_.value() - ti_.value())/(divisions_ + 1);

                for (int j=0; j<divisions_; j++)
                {
                    instant timej = instant(ti_.value() + (j + 1)*deltaT);

                    runTime_.setTime(instant(timej.name()), 0);

                    Info<< timej.name();

                    if (j < divisions_-1)
                    {
                        Info<< " ";
                    }

                    // Calculate times to read and weights
                    labelList indices;
                    scalarField weights;
                    interpolator_.valueWeights
                    (
                        runTime_.value(),
                        indices,
                        weights
                    );

                    const wordList selectedTimeNames
                    (
                        UIndirectList<word>(timeNames_, indices)()
                    );

                    //Info<< "For time " << runTime_.value()
                    //    << " need times " << selectedTimeNames
                    //    << " need weights " << weights << endl;


                    // Read on the objectRegistry all the required fields
                    ReadFields<GeoFieldType>
                    (
                        fieldIter()->name(),
                        mesh_,
                        selectedTimeNames
                    );

                    GeoFieldType fieldj
                    (
                        uniformInterpolate<GeoFieldType>
                        (
                            IOobject
                            (
                                fieldIter()->name(),
                                runTime_.timeName(),
                                fieldIter()->db(),
                                IOobject::NO_READ,
                                IOobject::NO_WRITE,
                                false
                            ),
                            fieldIter()->name(),
                            selectedTimeNames,
                            weights
                        )
                    );

                    fieldj.write();
                }

                Info<< ')';
            }
        }

        Info<< endl;
    }
}


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

int main(int argc, char *argv[])
{
    timeSelector::addOptions();
    #include "addRegionOption.H"
    argList::addOption
    (
        "fields",
        "list",
        "specify a list of fields to be interpolated. Eg, '(U T p)' - "
        "regular expressions not currently supported"
    );
    argList::addOption
    (
        "divisions",
        "integer",
        "specify number of temporal sub-divisions to create (default = 1)."
    );
    argList::addOption
    (
        "interpolationType",
        "word",
        "specify type of interpolation (linear or spline)"
    );

    #include "setRootCase.H"
    #include "createTime.H"
    runTime.functionObjects().off();

    HashSet<word> selectedFields;
    if (args.optionFound("fields"))
    {
        args.optionLookup("fields")() >> selectedFields;
    }
    if (selectedFields.size())
    {
        Info<< "Interpolating fields " << selectedFields << nl << endl;
    }
    else
    {
        Info<< "Interpolating all fields" << nl << endl;
    }


    int divisions = 1;
    if (args.optionFound("divisions"))
    {
        args.optionLookup("divisions")() >> divisions;
    }
    Info<< "Using " << divisions << " per time interval" << nl << endl;


    const word interpolationType = args.optionLookupOrDefault<word>
    (
        "interpolationType",
        "linear"
    );
    Info<< "Using interpolation " << interpolationType << nl << endl;


    instantList timeDirs = timeSelector::select0(runTime, args);

    scalarField timeVals(timeDirs.size());
    wordList timeNames(timeDirs.size());
    forAll(timeDirs, i)
    {
        timeVals[i] = timeDirs[i].value();
        timeNames[i] = timeDirs[i].name();
    }
    autoPtr<interpolationWeights> interpolatorPtr
    (
        interpolationWeights::New
        (
            interpolationType,
            timeVals
        )
    );


    #include "createNamedMesh.H"

    Info<< "Interpolating fields for times:" << endl;

    for (label timei = 0; timei < timeDirs.size() - 1; timei++)
    {
        runTime.setTime(timeDirs[timei], timei);

        // Read objects in time directory
        IOobjectList objects(mesh, runTime.timeName());

        fieldInterpolator interpolator
        (
            runTime,
            mesh,
            objects,
            selectedFields,
            timeDirs[timei],
            timeDirs[timei+1],
            interpolatorPtr(),
            timeNames,
            divisions
        );

        // Interpolate vol fields
        interpolator.interpolate<volScalarField>();
        interpolator.interpolate<volVectorField>();
        interpolator.interpolate<volSphericalTensorField>();
        interpolator.interpolate<volSymmTensorField>();
        interpolator.interpolate<volTensorField>();

        // Interpolate surface fields
        interpolator.interpolate<surfaceScalarField>();
        interpolator.interpolate<surfaceVectorField>();
        interpolator.interpolate<surfaceSphericalTensorField>();
        interpolator.interpolate<surfaceSymmTensorField>();
        interpolator.interpolate<surfaceTensorField>();
    }

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

    return 0;
}


// ************************************************************************* //
temporalInterpolate.C (8,660 bytes)   

wyldckat

2015-10-03 17:50

updater  

pPrime2.C (2,957 bytes)   
/*---------------------------------------------------------------------------*\
  =========                 |
  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
   \\    /   O peration     |
    \\  /    A nd           | Copyright (C) 2011-2015 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
    pPrime2

Description
    Calculates and writes the scalar field of pPrime2 (sqr(p - pMean)) at
    each time

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

#include "fvCFD.H"


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

int main(int argc, char *argv[])
{
    timeSelector::addOptions();
    #include "addRegionOption.H"

    #include "setRootCase.H"
    #include "createTime.H"

    instantList timeDirs = timeSelector::select0(runTime, args);

    #include "createNamedMesh.H"

    runTime.setTime(timeDirs.last(), timeDirs.size()-1);

    volScalarField pMean
    (
        IOobject
        (
            "pMean",
            runTime.timeName(),
            mesh,
            IOobject::MUST_READ
        ),
        mesh
    );

    forAll(timeDirs, timeI)
    {
        runTime.setTime(timeDirs[timeI], timeI);

        Info<< "Time = " << runTime.timeName() << endl;

        IOobject pheader
        (
            "p",
            runTime.timeName(),
            mesh,
            IOobject::MUST_READ
        );

        // Check p exists
        if (pheader.headerOk())
        {
            mesh.readUpdate();

            Info<< "    Reading p" << endl;
            volScalarField p(pheader, mesh);

            Info<< "    Calculating pPrime2" << endl;
            volScalarField pPrime2
            (
                IOobject
                (
                    "pPrime2",
                    runTime.timeName(),
                    mesh,
                    IOobject::NO_READ
                ),
                sqr(p - pMean)
            );
            pPrime2.write();
        }
        else
        {
            Info<< "    No p" << endl;
        }

        Info<< endl;
    }

    return 0;
}


// ************************************************************************* //
pPrime2.C (2,957 bytes)   

henry

2015-10-05 11:34

manager   ~0005385

Thanks Bruno.

Resolved by commit 445a626805d7f8e98be08f8257a5f18c657a419e

Issue History

Date Modified Username Field Change
2015-10-03 17:50 wyldckat New Issue
2015-10-03 17:50 wyldckat Status new => assigned
2015-10-03 17:50 wyldckat Assigned To => henry
2015-10-03 17:50 wyldckat File Added: wallGradU.C
2015-10-03 17:50 wyldckat File Added: temporalInterpolate.C
2015-10-03 17:50 wyldckat File Added: pPrime2.C
2015-10-03 17:51 wyldckat Relationship added related to 0001227
2015-10-05 11:34 henry Note Added: 0005385
2015-10-05 11:34 henry Status assigned => resolved
2015-10-05 11:34 henry Resolution open => fixed