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/solvers/engine/icoDyMEngineFoam/icoDyMEngineFoam.C
Hrvoje Jasak b2cd738861 Merge branch 'parallelTopo'
Conflicts:
	applications/utilities/mesh/manipulation/moveDyMEngineMesh/Make/options
	applications/utilities/mesh/manipulation/moveDyMEngineMesh/moveDyMEngineMesh.C
	applications/utilities/mesh/manipulation/setsToZones/setsToZones.C
	applications/utilities/parallelProcessing/decomposePar/decomposePar.C
	src/OpenFOAM/meshes/polyMesh/polyPatches/constraint/processor/processorPolyPatch.C
	src/decompositionMethods/decompositionMethods/patchConstrainedDecomp/patchConstrainedDecomp.C
	src/decompositionMethods/metisDecomp/metisDecomp.C
	src/dynamicMesh/topoChangerFvMesh/Make/files
	src/dynamicMesh/topoChangerFvMesh/mixerFvMesh/mixerFvMesh.C
	src/finiteArea/faMesh/faMesh.C
	src/finiteVolume/cfdTools/general/include/checkVolContinuity.H
	src/tetDecompositionFiniteElement/tetPolyMeshCellDecomp/tetPolyPatches/constraint/processor/calcProcessorTetPolyPatchCellDecompPointAddr.C
	src/turbulenceModels/compressible/RAS/derivedFvPatchFields/wallFunctions/epsilonWallFunctions/epsilonWallFunction/epsilonWallFunctionFvPatchScalarField.C
	src/turbulenceModels/compressible/RAS/derivedFvPatchFields/wallFunctions/omegaWallFunctions/omegaWallFunction/omegaWallFunctionFvPatchScalarField.C
	src/turbulenceModels/incompressible/RAS/derivedFvPatchFields/wallFunctions/epsilonWallFunctions/epsilonWallFunction/epsilonWallFunctionFvPatchScalarField.C
	src/turbulenceModels/incompressible/RAS/derivedFvPatchFields/wallFunctions/omegaWallFunctions/omegaWallFunction/omegaWallFunctionFvPatchScalarField.C
	tutorials/incompressible/icoDyMFoam/movingConeTopo/system/decomposeParDict
2012-04-25 10:38:12 +01:00

149 lines
4.2 KiB
C

/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright held by original author
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software; you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by the
Free Software Foundation; either version 2 of the License, or (at your
option) any later version.
OpenFOAM 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
along with OpenFOAM; if not, write to the Free Software Foundation,
Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
Application
icoDyMFoamEngine
Description
Transient solver for incompressible, laminar flow of Newtonian fluids
with dynamic mesh for in-cylinder internal combustion engine simulations
Author
Hrvoje Jasak, Wikki Ltd. All rights reserved
\*---------------------------------------------------------------------------*/
#include "fvCFD.H"
#include "dynamicFvMesh.H"
#include "engineTime.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
int main(int argc, char *argv[])
{
# include "setRootCase.H"
# include "createEngineTime.H"
# include "createDynamicFvMesh.H"
# include "initContinuityErrs.H"
# include "createFields.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
Info<< "\nStarting time loop\n" << endl;
while (runTime.run())
{
# include "readControls.H"
# include "checkTotalVolume.H"
# include "CourantNo.H"
// Make the fluxes absolute
fvc::makeAbsolute(phi, U);
# include "setDeltaT.H"
runTime++;
Info<< "Time = " << runTime.timeName() << nl << endl;
bool meshChanged = mesh.update();
reduce(meshChanged, orOp<bool>());
if (meshChanged)
{
# include "checkTotalVolume.H"
# include "correctPhi.H"
# include "CourantNo.H"
}
// Make the fluxes relative to the mesh motion
fvc::makeRelative(phi, U);
if (checkMeshCourantNo)
{
# include "meshCourantNo.H"
}
# include "UEqn.H"
// --- PISO loop
for (int corr = 0; corr < nCorr; corr++)
{
rUA = 1.0/UEqn.A();
U = rUA*UEqn.H();
phi = (fvc::interpolate(U) & mesh.Sf());
//+ fvc::ddtPhiCorr(rUA, U, phi);
adjustPhi(phi, U, p);
for (int nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++)
{
fvScalarMatrix pEqn
(
fvm::laplacian(rUA, p) == fvc::div(phi)
);
pEqn.setReference(pRefCell, pRefValue);
if (corr == nCorr - 1 && nonOrth == nNonOrthCorr)
{
pEqn.solve(mesh.solutionDict().solver(p.name() + "Final"));
}
else
{
pEqn.solve(mesh.solutionDict().solver(p.name()));
}
if (nonOrth == nNonOrthCorr)
{
phi -= pEqn.flux();
}
}
# include "continuityErrs.H"
// Make the fluxes relative to the mesh motion
fvc::makeRelative(phi, U);
U -= rUA*fvc::grad(p);
U.correctBoundaryConditions();
}
runTime.write();
Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
<< " ClockTime = " << runTime.elapsedClockTime() << " s"
<< nl << endl;
}
Info<< "End\n" << endl;
return(0);
}
// ************************************************************************* //