2010-05-12 13:27:55 +00:00
|
|
|
/*---------------------------------------------------------------------------*\
|
|
|
|
========= |
|
2013-12-11 16:09:41 +00:00
|
|
|
\\ / F ield | foam-extend: Open Source CFD
|
2015-05-17 13:32:07 +00:00
|
|
|
\\ / O peration | Version: 3.2
|
|
|
|
\\ / A nd | Web: http://www.foam-extend.org
|
|
|
|
\\/ M anipulation | For copyright notice see file Copyright
|
2010-05-12 13:27:55 +00:00
|
|
|
-------------------------------------------------------------------------------
|
|
|
|
License
|
2013-12-11 16:09:41 +00:00
|
|
|
This file is part of foam-extend.
|
2010-05-12 13:27:55 +00:00
|
|
|
|
2013-12-11 16:09:41 +00:00
|
|
|
foam-extend is free software: you can redistribute it and/or modify it
|
2010-05-12 13:27:55 +00:00
|
|
|
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
|
2010-05-12 13:27:55 +00:00
|
|
|
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.
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
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/>.
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
Application
|
|
|
|
icoMomentError
|
|
|
|
|
|
|
|
Description
|
|
|
|
Estimates error for the incompressible laminar CFD application icoFoam.
|
|
|
|
|
|
|
|
\*---------------------------------------------------------------------------*/
|
|
|
|
|
|
|
|
#include "fvCFD.H"
|
|
|
|
#include "linear.H"
|
|
|
|
#include "gaussConvectionScheme.H"
|
|
|
|
#include "gaussLaplacianScheme.H"
|
|
|
|
|
|
|
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
|
|
|
|
|
|
|
int main(int argc, char *argv[])
|
|
|
|
{
|
2010-09-21 14:32:04 +00:00
|
|
|
timeSelector::addOptions();
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
# include "setRootCase.H"
|
|
|
|
# include "createTime.H"
|
|
|
|
|
2010-09-21 14:32:04 +00:00
|
|
|
instantList timeDirs = timeSelector::select0(runTime, args);
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
# include "createMesh.H"
|
|
|
|
|
2010-09-21 14:32:04 +00:00
|
|
|
Info<< "\nEstimating error in the incompressible momentum equation\n"
|
|
|
|
<< "Reading transportProperties\n" << endl;
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
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-05-12 13:27:55 +00:00
|
|
|
{
|
2010-09-21 14:32:04 +00:00
|
|
|
runTime.setTime(timeDirs[timeI], timeI);
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
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"
|
|
|
|
|
|
|
|
volScalarField ek = 0.5*magSqr(U);
|
|
|
|
volTensorField gradU = fvc::grad(U);
|
|
|
|
|
|
|
|
// Divergence of the error in U squared
|
|
|
|
|
|
|
|
volScalarField L
|
|
|
|
(
|
|
|
|
IOobject
|
|
|
|
(
|
|
|
|
"L",
|
|
|
|
mesh.time().timeName(),
|
|
|
|
mesh,
|
|
|
|
IOobject::NO_READ,
|
|
|
|
IOobject::NO_WRITE
|
|
|
|
),
|
|
|
|
mesh,
|
|
|
|
dimensionedScalar("one", dimLength, 1.0)
|
|
|
|
);
|
|
|
|
|
|
|
|
L.internalField() =
|
|
|
|
mesh.V()/fvc::surfaceSum(mesh.magSf())().internalField();
|
|
|
|
|
|
|
|
// Warning: 4th row of this equation specially modified
|
|
|
|
// for the momentum equation. The "real" formulation would
|
|
|
|
// have diffusivity*(gradV && gradV)
|
|
|
|
volScalarField momError
|
|
|
|
(
|
|
|
|
IOobject
|
|
|
|
(
|
|
|
|
"momErrorL" + U.name(),
|
|
|
|
mesh.time().timeName(),
|
|
|
|
mesh,
|
|
|
|
IOobject::NO_READ,
|
|
|
|
IOobject::NO_WRITE
|
|
|
|
),
|
|
|
|
sqrt
|
|
|
|
(
|
|
|
|
2.0*mag
|
|
|
|
(
|
|
|
|
(
|
|
|
|
fv::gaussConvectionScheme<scalar>
|
|
|
|
(
|
|
|
|
mesh,
|
|
|
|
phi,
|
|
|
|
tmp<surfaceInterpolationScheme<scalar> >
|
|
|
|
(
|
|
|
|
new linear<scalar>(mesh)
|
|
|
|
)
|
|
|
|
).fvcDiv(phi, ek)
|
|
|
|
|
|
|
|
- nu*
|
|
|
|
fv::gaussLaplacianScheme<scalar, scalar>(mesh)
|
|
|
|
.fvcLaplacian
|
|
|
|
(
|
|
|
|
ek
|
|
|
|
)
|
|
|
|
- (U & fvc::grad(p))
|
|
|
|
// + nu*(gradU && gradU)
|
|
|
|
+ 0.5*nu*
|
|
|
|
(
|
|
|
|
gradU && (gradU + gradU.T())
|
|
|
|
)
|
|
|
|
)*L/(mag(U) + nu/L)
|
|
|
|
)
|
|
|
|
)
|
|
|
|
);
|
|
|
|
|
|
|
|
momError.boundaryField() = 0.0;
|
|
|
|
momError.write();
|
|
|
|
}
|
|
|
|
else
|
|
|
|
{
|
|
|
|
Info<< " No p or U" << endl;
|
|
|
|
}
|
|
|
|
|
|
|
|
Info<< endl;
|
|
|
|
}
|
|
|
|
|
|
|
|
Info<< "End\n" << endl;
|
|
|
|
|
2010-09-21 14:32:04 +00:00
|
|
|
return 0;
|
2010-05-12 13:27:55 +00:00
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
// ************************************************************************* //
|