Forward compatibility: flex
[foam-extend-3.2.git] / applications / utilities / errorEstimation / icoErrorEstimate / icoErrorEstimate.C
blobb9cd522bd7ea4f41d5684ec33de28fff9f863f55
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | foam-extend: Open Source CFD
4    \\    /   O peration     | Version:     3.2
5     \\  /    A nd           | Web:         http://www.foam-extend.org
6      \\/     M anipulation  | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
8 License
9     This file is part of foam-extend.
11     foam-extend 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 3 of the License, or (at your
14     option) any later version.
16     foam-extend is distributed in the hope that it will be useful, but
17     WITHOUT ANY WARRANTY; without even the implied warranty of
18     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19     General Public License for more details.
21     You should have received a copy of the GNU General Public License
22     along with foam-extend.  If not, see <http://www.gnu.org/licenses/>.
24 Application
25     icoErrorEstimate
27 Description
28     Estimates error for the incompressible laminar CFD application icoFoam.
30 \*---------------------------------------------------------------------------*/
32 #include "fvCFD.H"
33 #include "errorEstimate.H"
34 #include "resError.H"
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38 int main(int argc, char *argv[])
40     timeSelector::addOptions();
42 #   include "setRootCase.H"
43 #   include "createTime.H"
45     instantList timeDirs = timeSelector::select0(runTime, args);
47 #   include "createMesh.H"
49     Info<< "\nEstimating error in the incompressible momentum equation\n"
50         << "Reading transportProperties\n" << endl;
52     IOdictionary transportProperties
53     (
54         IOobject
55         (
56             "transportProperties",
57             runTime.constant(),
58             mesh,
59             IOobject::MUST_READ,
60             IOobject::NO_WRITE
61         )
62     );
64     dimensionedScalar nu
65     (
66         transportProperties.lookup("nu")
67     );
69     forAll(timeDirs, timeI)
70     {
71         runTime.setTime(timeDirs[timeI], timeI);
73         Info<< "Time = " << runTime.timeName() << endl;
75         mesh.readUpdate();
77         IOobject pHeader
78         (
79             "p",
80             runTime.timeName(),
81             mesh,
82             IOobject::MUST_READ
83         );
85         IOobject Uheader
86         (
87             "U",
88             runTime.timeName(),
89             mesh,
90             IOobject::MUST_READ
91         );
93         if (pHeader.headerOk() && Uheader.headerOk())
94         {
95             Info << "Reading p" << endl;
96             volScalarField p(pHeader, mesh);
98             Info << "Reading U" << endl;
99             volVectorField U(Uheader, mesh);
101 #           include "createPhi.H"
103             errorEstimate<vector> ee
104             (
105                 resError::div(phi, U)
106               - resError::laplacian(nu, U)
107              ==
108                -fvc::grad(p)
109             );
111             volVectorField e = ee.error();
112             e.write();
113             mag(e)().write();
114         }
115         else
116         {
117             Info<< "    No p or U" << endl;
118         }
120         Info<< endl;
121     }
123     Info<< "End\n" << endl;
125     return 0;
129 // ************************************************************************* //