/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | \\ / A nd | Copyright (C) 1991-2007 OpenCFD Ltd. \\/ 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 2 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, write to the Free Software Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA Description Sample field data with a choice of interpolation schemes, sampling options and write formats. interpolationScheme : choice of cell : use cell-centre value only; constant over cells cellPoint : use cell-centre and vertex values cellPointFace : use cell-centre, vertex and face values. sample: choice of uniform evenly distributed points on line face one point per face intersection midPoint one point per cell, inbetween two face intersections midPointAndFace combination of face and midPoint curve specified points, not nessecary on line, uses tracking cloud specified points, uses findCell writeFormat : choice of xmgr jplot gnuplot raw \*---------------------------------------------------------------------------*/ #include "Pstream.H" #include "argList.H" #include "OSspecific.H" #include "Cloud.H" #include "passiveParticle.H" #include "meshSearch.H" #include "interpolation.H" #include "volPointInterpolation.H" #include "writer.H" #include "sampleSet.H" #include "volFieldSampler.H" #include "dictionaryEntry.H" #include "combineSampleSets.H" #include "combineSampleValues.H" using namespace Foam; template void writeSampleFile ( const coordSet& masterSampleSet, const PtrList >& masterFields, const label setI, const fileName& timeDir, const word& writeFormat ) { if (masterFields.size() > 0) { wordList valueSetNames(masterFields.size()); List*> valueSets(masterFields.size()); forAll(masterFields, fieldI) { valueSetNames[fieldI] = masterFields[fieldI].name(); valueSets[fieldI] = &masterFields[fieldI][setI]; } autoPtr > formatter ( writer::New(writeFormat) ); fileName fName ( timeDir/formatter().getFileName(masterSampleSet, valueSetNames) ); Info<< "Writing fields to " << fName << endl; formatter().write ( masterSampleSet, valueSetNames, valueSets, OFstream(fName)() ); } } // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // // Main program: int main(int argc, char *argv[]) { # include "addTimeOptions.H" # include "setRootCase.H" # include "createTime.H" // Get times list instantList Times = runTime.times(); // set startTime and endTime depending on -time and -latestTime options # include "checkTimeOptions.H" runTime.setTime(Times[startTime], startTime); # include "createMesh.H" // // Hack: initialize Cloud to initialize the processor table so from // now on we can use cloud on single processors only. // Cloud dummyCloud(mesh, IDLList()); // Read control dictionary IOdictionary sampleDict ( IOobject ( "sampleDict", runTime.system(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE ) ); const word interpolationScheme(sampleDict.lookup("interpolationScheme")); const wordList fieldNames = sampleDict.lookup("fields"); // // Construct writers // word writeFormat(sampleDict.lookup("writeFormat")); // // Construct interpolation dictionary (same interpolation for all fields) // dictionary interpolationSchemes; forAll(fieldNames, fieldI) { interpolationSchemes.add ( fieldNames[fieldI], interpolationScheme ); } // Set up interpolation autoPtr pMeshPtr(new pointMesh(mesh)); autoPtr pInterpPtr ( new volPointInterpolation(mesh, pMeshPtr()) ); // Set up mesh searching meshSearch searchEngine(mesh, true); fileName samplePath; if (Pstream::master()) { if (Pstream::parRun()) { samplePath = runTime.path()/".."/"samples"; } else { samplePath = runTime.path()/"samples"; } if (exists(samplePath)) { Info<< "Deleting samples/ directory" << endl << endl; rmDir(samplePath); } } fileName oldPointsDir("constant"); for (label i=startTime; i > sampledScalarFields ( fieldNames.size() ); PtrList > sampledVectorFields ( fieldNames.size() ); PtrList > sampledSphericalTensorFields ( fieldNames.size() ); PtrList > sampledSymmTensorFields ( fieldNames.size() ); PtrList > sampledTensorFields ( fieldNames.size() ); // // Do actual interpolation // label nScalarFields = 0; label nVectorFields = 0; label nSphericalTensorFields = 0; label nSymmTensorFields = 0; label nTensorFields = 0; forAll(fieldNames, fieldI) { const word& fieldName = fieldNames[fieldI]; IOobject fieldHeader ( fieldName, runTime.timeName(), mesh, IOobject::MUST_READ ); // Determine number of processor actually having this field label fieldFound = (fieldHeader.headerOk() ? 1 : 0); reduce(fieldFound, sumOp