Initial commit for version 2.0.x patch release
[OpenFOAM-2.0.x.git] / applications / utilities / postProcessing / miscellaneous / writeCellCentres / writeCellCentres.C
blob57fcf209b55e8c7422ec65bb2284dc876acef0a6
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2004-2010 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
11     OpenFOAM is free software: you can redistribute it and/or modify it
12     under the terms of the GNU General Public License as published by
13     the Free Software Foundation, either version 3 of the License, or
14     (at your option) any later version.
16     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
19     for more details.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
24 Description
25     Write the three components of the cell centres as volScalarFields so
26     they can be used in postprocessing in thresholding.
28 \*---------------------------------------------------------------------------*/
30 #include "argList.H"
31 #include "timeSelector.H"
32 #include "Time.H"
33 #include "fvMesh.H"
34 #include "vectorIOField.H"
35 #include "volFields.H"
37 using namespace Foam;
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 // Main program:
43 int main(int argc, char *argv[])
45     timeSelector::addOptions();
47 #   include "setRootCase.H"
48 #   include "createTime.H"
50     instantList timeDirs = timeSelector::select0(runTime, args);
52 #   include "createMesh.H"
54     forAll(timeDirs, timeI)
55     {
56         runTime.setTime(timeDirs[timeI], timeI);
58         Info<< "Time = " << runTime.timeName() << endl;
60         // Check for new mesh
61         mesh.readUpdate();
63         volVectorField cc
64         (
65             IOobject
66             (
67                 "cellCentres",
68                 runTime.timeName(),
69                 mesh,
70                 IOobject::NO_READ,
71                 IOobject::AUTO_WRITE
72             ),
73             mesh.C()
74         );
76         // Info<< "Writing cellCentre positions to " << cc.name() << " in "
77         //     << runTime.timeName() << endl;
78         //
79         // cc.write();
81         Info<< "Writing components of cellCentre positions to volScalarFields"
82             << " ccx, ccy, ccz in " <<  runTime.timeName() << endl;
84         for (direction i=0; i<vector::nComponents; i++)
85         {
86             volScalarField cci
87             (
88                 IOobject
89                 (
90                     "cc" + word(vector::componentNames[i]),
91                     runTime.timeName(),
92                     mesh,
93                     IOobject::NO_READ,
94                     IOobject::AUTO_WRITE
95                 ),
96                 mesh.C().component(i)
97             );
99             cci.write();
100         }
101     }
103     Info<< "\nEnd\n" << endl;
105     return 0;
109 // ************************************************************************* //