Merge branch 'upstream/OpenFOAM' into master
[freefoam.git] / applications / utilities / errorEstimation / momentScalarError / momentScalarError.C
blob27efb5a594013d7a1cfa57a3e1722e0f5f4165e3
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 Usage
34     - momentScalarError [OPTIONS]
36     @param -noZero \n
37     Ignore timestep 0.
39     @param -constant \n
40     Include the constant directory.
42     @param -time \<time\>\n
43     Apply only to specific time.
45     @param -latestTime \n
46     Only apply to latest time step.
48     @param -case \<dir\>\n
49     Case directory.
51     @param -parallel \n
52     Run in parallel.
54     @param -help \n
55     Display help message.
57     @param -doc \n
58     Display Doxygen API documentation page for this application.
60     @param -srcDoc \n
61     Display Doxygen source documentation page for this application.
63 \*---------------------------------------------------------------------------*/
65 #include <finiteVolume/fvCFD.H>
66 #include <finiteVolume/linear.H>
67 #include <finiteVolume/gaussConvectionScheme.H>
68 #include <finiteVolume/gaussLaplacianScheme.H>
70 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
72 int main(int argc, char *argv[])
74     timeSelector::addOptions();
76 #   include <OpenFOAM/setRootCase.H>
77 #   include <OpenFOAM/createTime.H>
79     instantList timeDirs = timeSelector::select0(runTime, args);
81 #   include <OpenFOAM/createMesh.H>
83     Info<< "\nEstimating error in scalar transport equation\n"
84         << "Reading transportProperties\n" << endl;
86     IOdictionary transportProperties
87     (
88         IOobject
89         (
90             "transportProperties",
91             runTime.constant(),
92             mesh,
93             IOobject::MUST_READ,
94             IOobject::NO_WRITE
95         )
96     );
99     Info<< "Reading diffusivity DT\n" << endl;
101     dimensionedScalar DT
102     (
103         transportProperties.lookup("DT")
104     );
107     forAll(timeDirs, timeI)
108     {
109         runTime.setTime(timeDirs[timeI], timeI);
111         Info<< "Time = " << runTime.timeName() << endl;
113         mesh.readUpdate();
115         IOobject THeader
116         (
117             "T",
118             runTime.timeName(),
119             mesh,
120             IOobject::MUST_READ
121         );
123         IOobject Uheader
124         (
125             "U",
126             runTime.timeName(),
127             mesh,
128             IOobject::MUST_READ
129         );
131         if (THeader.headerOk() && Uheader.headerOk())
132         {
133             Info << "Reading T" << endl;
134             volScalarField T(THeader, mesh);
136             Info << "Reading U" << endl;
137             volVectorField U(Uheader, mesh);
139 #           include <finiteVolume/createPhi.H>
141             volVectorField gradT = fvc::grad(T);
143             volScalarField TE = 0.5*sqr(T);
145             volScalarField L
146             (
147                 IOobject
148                 (
149                     "L",
150                     mesh.time().timeName(),
151                     mesh,
152                     IOobject::NO_READ,
153                     IOobject::NO_WRITE
154                 ),
155                 mesh,
156                 dimensionedScalar("one", dimLength, 1.0)
157             );
159             L.internalField() =
160                 mesh.V()/fvc::surfaceSum(mesh.magSf())().internalField();
162             // Divergence of the error in the T squared
163             volScalarField momError
164             (
165                 IOobject
166                 (
167                     "momErrorL" + T.name(),
168                     mesh.time().timeName(),
169                     mesh,
170                     IOobject::NO_READ,
171                     IOobject::NO_WRITE
172                 ),
173                 sqrt
174                 (
175                     2.0*mag
176                     (
177                         (
178                             fv::gaussConvectionScheme<scalar>
179                             (
180                                 mesh,
181                                 phi,
182                                 tmp<surfaceInterpolationScheme<scalar> >
183                                 (
184                                     new linear<scalar>(mesh)
185                                 )
186                             ).fvcDiv(phi, TE)
188                           - DT*
189                             fv::gaussLaplacianScheme<scalar, scalar>(mesh)
190                            .fvcLaplacian
191                             (
192                                 TE
193                             )
194                           + DT*(gradT & gradT)
195                         )*L/(mag(U) + DT/L)
196                     )
197                 )
198             );
200             momError.boundaryField() = 0.0;
201             momError.write();
202         }
203         else
204         {
205             Info<< "    No T or U" << endl;
206         }
208         Info<< endl;
209     }
211     Info<< "End\n" << endl;
213     return 0;
217 // ************************ vim: set sw=4 sts=4 et: ************************ //