Initial commit for version 2.0.x patch release
[OpenFOAM-2.0.x.git] / applications / utilities / postProcessing / stressField / stressComponents / stressComponents.C
blob73cdd13690d02c4257c83d2d0f4d2c9785550217
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 Application
25     stressComponents
27 Description
28     Calculates and writes the scalar fields of the six components of the stress
29     tensor sigma for each time.
31 \*---------------------------------------------------------------------------*/
33 #include "fvCFD.H"
34 #include "incompressible/singlePhaseTransportModel/singlePhaseTransportModel.H"
35 #include "zeroGradientFvPatchFields.H"
37 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39 int main(int argc, char *argv[])
41     timeSelector::addOptions();
43 #   include "setRootCase.H"
44 #   include "createTime.H"
46     instantList timeDirs = timeSelector::select0(runTime, args);
48 #   include "createMesh.H"
50     forAll(timeDirs, timeI)
51     {
52         runTime.setTime(timeDirs[timeI], timeI);
54         Info<< "Time = " << runTime.timeName() << endl;
56         IOobject Uheader
57         (
58             "U",
59             runTime.timeName(),
60             mesh,
61             IOobject::MUST_READ
62         );
64         // Check U exists
65         if (Uheader.headerOk())
66         {
67             mesh.readUpdate();
69             Info<< "    Reading U" << endl;
70             volVectorField U(Uheader, mesh);
72 #           include "createPhi.H"
74             singlePhaseTransportModel laminarTransport(U, phi);
76             volSymmTensorField sigma
77             (
78                 IOobject
79                 (
80                     "sigma",
81                     runTime.timeName(),
82                     mesh,
83                     IOobject::NO_READ,
84                     IOobject::AUTO_WRITE
85                 ),
86                 laminarTransport.nu()*2*dev(symm(fvc::grad(U)))
87             );
90             volScalarField sigmaxx
91             (
92                 IOobject
93                 (
94                     "sigmaxx",
95                     runTime.timeName(),
96                     mesh,
97                     IOobject::NO_READ
98                 ),
99                 sigma.component(symmTensor::XX)
100             );
101             sigmaxx.write();
103             volScalarField sigmayy
104             (
105                 IOobject
106                 (
107                     "sigmayy",
108                     runTime.timeName(),
109                     mesh,
110                     IOobject::NO_READ
111                 ),
112                 sigma.component(symmTensor::YY)
113             );
114             sigmayy.write();
116             volScalarField sigmazz
117             (
118                 IOobject
119                 (
120                     "sigmazz",
121                     runTime.timeName(),
122                     mesh,
123                     IOobject::NO_READ
124                 ),
125                 sigma.component(symmTensor::ZZ)
126             );
127             sigmazz.write();
129             volScalarField sigmaxy
130             (
131                 IOobject
132                 (
133                     "sigmaxy",
134                     runTime.timeName(),
135                     mesh,
136                     IOobject::NO_READ
137                 ),
138                 sigma.component(symmTensor::XY)
139             );
140             sigmaxy.write();
142             volScalarField sigmaxz
143             (
144                 IOobject
145                 (
146                     "sigmaxz",
147                     runTime.timeName(),
148                     mesh,
149                     IOobject::NO_READ
150                 ),
151                 sigma.component(symmTensor::XZ)
152             );
153             sigmaxz.write();
155             volScalarField sigmayz
156             (
157                 IOobject
158                 (
159                     "sigmayz",
160                     runTime.timeName(),
161                     mesh,
162                     IOobject::NO_READ
163                 ),
164                 sigma.component(symmTensor::YZ)
165             );
166             sigmayz.write();
168             volVectorField Ub
169             (
170                 IOobject
171                 (
172                     "Ub",
173                     runTime.timeName(),
174                     mesh,
175                     IOobject::NO_READ
176                 ),
177                 U,
178                 zeroGradientFvPatchVectorField::typeName
179             );
180             Ub.correctBoundaryConditions();
181             Ub.write();
183             volScalarField sigmaUn
184             (
185                 IOobject
186                 (
187                     "sigmaUn",
188                     runTime.timeName(),
189                     mesh,
190                     IOobject::NO_READ
191                 ),
192                 0.0*sigma.component(symmTensor::YZ)
193             );
195             forAll(sigmaUn.boundaryField(), patchI)
196             {
197                 sigmaUn.boundaryField()[patchI] =
198                 (
199                     mesh.boundary()[patchI].nf()
200                   & sigma.boundaryField()[patchI]
201                 )().component(vector::X);
202             }
204             sigmaUn.write();
205         }
206         else
207         {
208             Info<< "    No U" << endl;
209         }
211         Info<< endl;
212     }
214     Info<< "End" << endl;
216     return 0;
220 // ************************************************************************* //