initial commit for version 1.5.x patch release
[OpenFOAM-1.5.x.git] / applications / utilities / errorEstimation / momentScalarError / momentScalarError.C
blob762a18aa46f9b986a00e6595fea211d7cfa0dd0c
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2008 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[])
44 #   include "addTimeOptions.H"
45 #   include "setRootCase.H"
47     Info<< "\nEstimating error in scalar transport equation\n" << endl;
49 #   include "createTime.H"
51     // Get times list
52     instantList Times = runTime.times();
54 #   include "checkTimeOptions.H"
56     runTime.setTime(Times[startTime], startTime);
58 #   include "createMesh.H"
60     Info<< "Reading transportProperties\n" << endl;
62     IOdictionary transportProperties
63     (
64         IOobject
65         (
66             "transportProperties",
67             runTime.constant(),
68             mesh,
69             IOobject::MUST_READ,
70             IOobject::NO_WRITE
71         )
72     );
75     Info<< "Reading diffusivity DT\n" << endl;
77     dimensionedScalar DT
78     (
79         transportProperties.lookup("DT")
80     );
83     for (label i=startTime; i<endTime; i++)
84     {
85         runTime.setTime(Times[i], i);
87         Info<< "Time = " << runTime.timeName() << endl;
89         mesh.readUpdate();
91         IOobject THeader
92         (
93             "T",
94             runTime.timeName(),
95             mesh,
96             IOobject::MUST_READ
97         );
99         IOobject Uheader
100         (
101             "U",
102             runTime.timeName(),
103             mesh,
104             IOobject::MUST_READ
105         );
107         if (THeader.headerOk() && Uheader.headerOk())
108         {
109             Info << "Reading T" << endl;
110             volScalarField T(THeader, mesh);
112             Info << "Reading U" << endl;
113             volVectorField U(Uheader, mesh);
115 #           include "createPhi.H"
117             volVectorField gradT = fvc::grad(T);
119             volScalarField TE = 0.5*sqr(T);
121             volScalarField L
122             (
123                 IOobject
124                 (
125                     "L",
126                     mesh.time().timeName(),
127                     mesh,
128                     IOobject::NO_READ,
129                     IOobject::NO_WRITE
130                 ),
131                 mesh,
132                 dimensionedScalar("one", dimLength, 1.0)
133             );
135             L.internalField() =
136                 mesh.V()/fvc::surfaceSum(mesh.magSf())().internalField();
138             // Divergence of the error in the T squared
139             volScalarField momError
140             (
141                 IOobject
142                 (
143                     "momErrorL" + T.name(),
144                     mesh.time().timeName(),
145                     mesh,
146                     IOobject::NO_READ,
147                     IOobject::NO_WRITE
148                 ),
149                 sqrt
150                 (
151                     2.0*mag
152                     (
153                         (
154                             fv::gaussConvectionScheme<scalar>
155                             (
156                                 mesh,
157                                 phi,
158                                 tmp<surfaceInterpolationScheme<scalar> >
159                                 (
160                                     new linear<scalar>(mesh)
161                                 )
162                             ).fvcDiv(phi, TE)
164                           - DT*
165                             fv::gaussLaplacianScheme<scalar, scalar>(mesh)
166                            .fvcLaplacian
167                             (
168                                 TE
169                             )
170                           + DT*(gradT & gradT)
171                         )*L/(mag(U) + DT/L)
172                     )
173                 )
174             );
176             momError.boundaryField() = 0.0;
177             momError.write();
178         }
179         else
180         {
181             Info<< "    No T or U" << endl;
182         }
184         Info<< endl;
185     }
187     Info<< "End\n" << endl;
189     return(0);
193 // ************************************************************************* //