initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / applications / utilities / errorEstimation / momentScalarError / momentScalarError.C
blobf1767e75b01e6f2aa025bf39be577c686f2effa4
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 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 the
13     Free Software Foundation; either version 2 of the License, or (at your
14     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, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 Application
26     momentScalarError
28 Description
29     Estimates the error in the solution for a scalar transport equation in the
30     standard form
32 \*---------------------------------------------------------------------------*/
34 #include "fvCFD.H"
35 #include "linear.H"
36 #include "gaussConvectionScheme.H"
37 #include "gaussLaplacianScheme.H"
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 int main(int argc, char *argv[])
43     timeSelector::addOptions();
45 #   include "setRootCase.H"
46 #   include "createTime.H"
48     instantList timeDirs = timeSelector::select0(runTime, args);
50 #   include "createMesh.H"
52     Info<< "\nEstimating error in scalar transport equation\n"
53         << "Reading transportProperties\n" << endl;
55     IOdictionary transportProperties
56     (
57         IOobject
58         (
59             "transportProperties",
60             runTime.constant(),
61             mesh,
62             IOobject::MUST_READ,
63             IOobject::NO_WRITE
64         )
65     );
68     Info<< "Reading diffusivity DT\n" << endl;
70     dimensionedScalar DT
71     (
72         transportProperties.lookup("DT")
73     );
76     forAll(timeDirs, timeI)
77     {
78         runTime.setTime(timeDirs[timeI], timeI);
80         Info<< "Time = " << runTime.timeName() << endl;
82         mesh.readUpdate();
84         IOobject THeader
85         (
86             "T",
87             runTime.timeName(),
88             mesh,
89             IOobject::MUST_READ
90         );
92         IOobject Uheader
93         (
94             "U",
95             runTime.timeName(),
96             mesh,
97             IOobject::MUST_READ
98         );
100         if (THeader.headerOk() && Uheader.headerOk())
101         {
102             Info << "Reading T" << endl;
103             volScalarField T(THeader, mesh);
105             Info << "Reading U" << endl;
106             volVectorField U(Uheader, mesh);
108 #           include "createPhi.H"
110             volVectorField gradT = fvc::grad(T);
112             volScalarField TE = 0.5*sqr(T);
114             volScalarField L
115             (
116                 IOobject
117                 (
118                     "L",
119                     mesh.time().timeName(),
120                     mesh,
121                     IOobject::NO_READ,
122                     IOobject::NO_WRITE
123                 ),
124                 mesh,
125                 dimensionedScalar("one", dimLength, 1.0)
126             );
128             L.internalField() =
129                 mesh.V()/fvc::surfaceSum(mesh.magSf())().internalField();
131             // Divergence of the error in the T squared
132             volScalarField momError
133             (
134                 IOobject
135                 (
136                     "momErrorL" + T.name(),
137                     mesh.time().timeName(),
138                     mesh,
139                     IOobject::NO_READ,
140                     IOobject::NO_WRITE
141                 ),
142                 sqrt
143                 (
144                     2.0*mag
145                     (
146                         (
147                             fv::gaussConvectionScheme<scalar>
148                             (
149                                 mesh,
150                                 phi,
151                                 tmp<surfaceInterpolationScheme<scalar> >
152                                 (
153                                     new linear<scalar>(mesh)
154                                 )
155                             ).fvcDiv(phi, TE)
157                           - DT*
158                             fv::gaussLaplacianScheme<scalar, scalar>(mesh)
159                            .fvcLaplacian
160                             (
161                                 TE
162                             )
163                           + DT*(gradT & gradT)
164                         )*L/(mag(U) + DT/L)
165                     )
166                 )
167             );
169             momError.boundaryField() = 0.0;
170             momError.write();
171         }
172         else
173         {
174             Info<< "    No T or U" << endl;
175         }
177         Info<< endl;
178     }
180     Info<< "End\n" << endl;
182     return 0;
186 // ************************************************************************* //