2015-05-11 10:41:43 +00:00
|
|
|
/*---------------------------------------------------------------------------*\
|
|
|
|
========= |
|
2015-05-17 13:32:07 +00:00
|
|
|
\\ / F ield | foam-extend: Open Source CFD
|
2018-05-29 07:35:20 +00:00
|
|
|
\\ / O peration | Version: 4.1
|
2015-05-17 13:32:07 +00:00
|
|
|
\\ / A nd | Web: http://www.foam-extend.org
|
|
|
|
\\/ M anipulation | For copyright notice see file Copyright
|
2015-05-11 10:41:43 +00:00
|
|
|
-------------------------------------------------------------------------------
|
|
|
|
License
|
2015-05-17 13:32:07 +00:00
|
|
|
This file is part of foam-extend.
|
2015-05-11 10:41:43 +00:00
|
|
|
|
2015-05-17 13:32:07 +00:00
|
|
|
foam-extend is free software: you can redistribute it and/or modify it
|
2015-05-11 10:41:43 +00:00
|
|
|
under the terms of the GNU General Public License as published by the
|
2015-05-17 13:32:07 +00:00
|
|
|
Free Software Foundation, either version 3 of the License, or (at your
|
2015-05-11 10:41:43 +00:00
|
|
|
option) any later version.
|
|
|
|
|
2015-05-17 13:32:07 +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.
|
2015-05-11 10:41:43 +00:00
|
|
|
|
|
|
|
You should have received a copy of the GNU General Public License
|
2015-05-17 13:32:07 +00:00
|
|
|
along with foam-extend. If not, see <http://www.gnu.org/licenses/>.
|
2015-05-11 10:41:43 +00:00
|
|
|
|
|
|
|
Application
|
|
|
|
writeIbMasks
|
|
|
|
|
|
|
|
Description
|
|
|
|
Calculate and write immersed boundary masks
|
|
|
|
|
|
|
|
\*---------------------------------------------------------------------------*/
|
|
|
|
|
2017-12-06 21:19:44 +00:00
|
|
|
#include "calc.H"
|
|
|
|
#include "fvc.H"
|
|
|
|
#include "fvMatrices.H"
|
2015-05-11 10:41:43 +00:00
|
|
|
#include "immersedBoundaryFvPatch.H"
|
|
|
|
|
|
|
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
|
|
|
|
2017-12-06 21:19:44 +00:00
|
|
|
void Foam::calc(const argList& args, const Time& runTime, const fvMesh& mesh)
|
2015-05-11 10:41:43 +00:00
|
|
|
{
|
2017-12-06 21:19:44 +00:00
|
|
|
Info<< nl << "Calculating gamma" << endl;
|
|
|
|
volScalarField gamma
|
|
|
|
(
|
|
|
|
IOobject
|
|
|
|
(
|
|
|
|
"gamma",
|
|
|
|
runTime.timeName(),
|
|
|
|
mesh,
|
|
|
|
IOobject::NO_READ,
|
|
|
|
IOobject::AUTO_WRITE
|
|
|
|
),
|
|
|
|
mesh,
|
2017-12-30 09:32:32 +00:00
|
|
|
dimensionedScalar("one", dimless, 1)
|
2017-12-06 21:19:44 +00:00
|
|
|
);
|
|
|
|
gamma.internalField() = mesh.V()/mesh.cellVolumes();
|
2017-12-30 09:32:32 +00:00
|
|
|
gamma.write();
|
|
|
|
|
|
|
|
// Report minimal live cell volume
|
|
|
|
scalar minLiveGamma = GREAT;
|
|
|
|
label minLiveCell = -1;
|
|
|
|
const scalarField& gammaIn = gamma.internalField();
|
|
|
|
|
|
|
|
forAll (mesh.boundary(), patchI)
|
|
|
|
{
|
|
|
|
if (isA<immersedBoundaryFvPatch>(mesh.boundary()[patchI]))
|
|
|
|
{
|
|
|
|
const immersedBoundaryFvPatch& ibPatch =
|
|
|
|
refCast<const immersedBoundaryFvPatch>
|
|
|
|
(
|
|
|
|
mesh.boundary()[patchI]
|
|
|
|
);
|
|
|
|
|
|
|
|
const labelList& ibCells = ibPatch.ibPolyPatch().ibCells();
|
|
|
|
|
|
|
|
forAll (ibCells, dcI)
|
|
|
|
{
|
|
|
|
if (gammaIn[ibCells[dcI]] < minLiveGamma)
|
|
|
|
{
|
|
|
|
minLiveGamma = gammaIn[ibCells[dcI]];
|
|
|
|
minLiveCell = ibCells[dcI];
|
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
Info<< "Min live cell " << minLiveCell
|
|
|
|
<< " gamma = " << minLiveGamma
|
|
|
|
<< endl;
|
2017-12-06 21:19:44 +00:00
|
|
|
|
|
|
|
Info<< nl << "Calculating sGamma" << endl;
|
|
|
|
surfaceScalarField sGamma
|
|
|
|
(
|
|
|
|
IOobject
|
|
|
|
(
|
|
|
|
"sGamma",
|
|
|
|
runTime.timeName(),
|
|
|
|
mesh,
|
|
|
|
IOobject::NO_READ,
|
|
|
|
IOobject::AUTO_WRITE
|
|
|
|
),
|
|
|
|
mesh,
|
|
|
|
dimensionedScalar("one", dimless, 0)
|
|
|
|
);
|
|
|
|
|
|
|
|
const surfaceScalarField& magSf = mesh.magSf();
|
|
|
|
const scalarField magFaceAreas = mag(mesh.faceAreas());
|
2017-12-30 09:32:32 +00:00
|
|
|
|
2017-12-06 21:19:44 +00:00
|
|
|
sGamma.internalField() =
|
|
|
|
magSf.internalField()/
|
|
|
|
scalarField::subField(magFaceAreas, mesh.nInternalFaces());
|
|
|
|
|
|
|
|
forAll (mesh.boundary(), patchI)
|
|
|
|
{
|
|
|
|
if (!isA<immersedBoundaryFvPatch>(mesh.boundary()[patchI]))
|
|
|
|
{
|
|
|
|
sGamma.boundaryField()[patchI] =
|
|
|
|
magSf.boundaryField()[patchI]/
|
|
|
|
mesh.boundary()[patchI].patchSlice(magFaceAreas);
|
|
|
|
|
|
|
|
gamma.boundaryField()[patchI] =
|
|
|
|
sGamma.boundaryField()[patchI];
|
|
|
|
}
|
|
|
|
}
|
2017-12-30 09:32:32 +00:00
|
|
|
|
2017-12-06 21:19:44 +00:00
|
|
|
sGamma.write();
|
|
|
|
|
|
|
|
// Check consistency of face area vectors
|
|
|
|
|
|
|
|
Info<< nl << "Calculating divSf" << endl;
|
|
|
|
volVectorField divSf
|
|
|
|
(
|
|
|
|
"divSf",
|
|
|
|
fvc::div(mesh.Sf())
|
|
|
|
);
|
|
|
|
divSf.write();
|
|
|
|
|
|
|
|
// Check divergence of face area vectors
|
|
|
|
scalarField magDivSf = mag(divSf)().internalField();
|
|
|
|
|
|
|
|
Info<< "Face areas divergence (min, max, average): "
|
|
|
|
<< "(" << min(magDivSf) << " " << max(magDivSf)
|
|
|
|
<< " " << average(magDivSf) << ")"
|
|
|
|
<< endl;
|
|
|
|
|
|
|
|
if (max(magDivSf) > 1e-9)
|
|
|
|
{
|
|
|
|
WarningIn("writeIbMasks")
|
|
|
|
<< "Possible problem with immersed boundary face area vectors: "
|
|
|
|
<< max(magDivSf)
|
|
|
|
<< endl;
|
|
|
|
}
|
2017-12-30 09:32:32 +00:00
|
|
|
|
|
|
|
Info<< endl;
|
2015-05-11 10:41:43 +00:00
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
// ************************************************************************* //
|