This repository has been archived on 2023-11-20. You can view files and clone it, but cannot push or open issues or pull requests.
foam-extend4.1-coherent-io/applications/utilities/errorEstimation/icoErrorEstimate/icoErrorEstimate.C

130 lines
3.3 KiB
C
Raw Normal View History

/*---------------------------------------------------------------------------*\
========= |
2013-12-11 16:09:41 +00:00
\\ / F ield | foam-extend: Open Source CFD
\\ / O peration | Version: 3.2
\\ / A nd | Web: http://www.foam-extend.org
\\/ M anipulation | For copyright notice see file Copyright
-------------------------------------------------------------------------------
License
2013-12-11 16:09:41 +00:00
This file is part of foam-extend.
2013-12-11 16:09:41 +00:00
foam-extend is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by the
2013-12-11 16:09:41 +00:00
Free Software Foundation, either version 3 of the License, or (at your
option) any later version.
2013-12-11 16:09:41 +00:00
foam-extend is distributed in the hope that it will be useful, but
WITHOUT ANY WARRANTY; without even the implied warranty of
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
General Public License for more details.
You should have received a copy of the GNU General Public License
2013-12-11 16:09:41 +00:00
along with foam-extend. If not, see <http://www.gnu.org/licenses/>.
Application
icoErrorEstimate
Description
Estimates error for the incompressible laminar CFD application icoFoam.
\*---------------------------------------------------------------------------*/
#include "fvCFD.H"
#include "errorEstimate.H"
#include "resError.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
int main(int argc, char *argv[])
{
2010-09-21 14:32:04 +00:00
timeSelector::addOptions();
# include "setRootCase.H"
# include "createTime.H"
2010-09-21 14:32:04 +00:00
instantList timeDirs = timeSelector::select0(runTime, args);
# include "createMesh.H"
2010-09-21 14:32:04 +00:00
Info<< "\nEstimating error in the incompressible momentum equation\n"
<< "Reading transportProperties\n" << endl;
IOdictionary transportProperties
(
IOobject
(
"transportProperties",
runTime.constant(),
mesh,
IOobject::MUST_READ,
IOobject::NO_WRITE
)
);
dimensionedScalar nu
(
transportProperties.lookup("nu")
);
2010-09-21 14:32:04 +00:00
forAll(timeDirs, timeI)
{
2010-09-21 14:32:04 +00:00
runTime.setTime(timeDirs[timeI], timeI);
Info<< "Time = " << runTime.timeName() << endl;
mesh.readUpdate();
IOobject pHeader
(
"p",
runTime.timeName(),
mesh,
IOobject::MUST_READ
);
IOobject Uheader
(
"U",
runTime.timeName(),
mesh,
IOobject::MUST_READ
);
if (pHeader.headerOk() && Uheader.headerOk())
{
Info << "Reading p" << endl;
volScalarField p(pHeader, mesh);
Info << "Reading U" << endl;
volVectorField U(Uheader, mesh);
# include "createPhi.H"
errorEstimate<vector> ee
(
resError::div(phi, U)
- resError::laplacian(nu, U)
==
-fvc::grad(p)
);
volVectorField e = ee.error();
e.write();
mag(e)().write();
}
else
{
Info<< " No p or U" << endl;
}
Info<< endl;
}
Info<< "End\n" << endl;
2010-09-21 14:32:04 +00:00
return 0;
}
// ************************************************************************* //