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/coupled/conjugateHeatFoam/conjugateHeatFoam.C
Henrik Rusche c0dc374d98 Merge PISO/PIMPLE solver updates. Author: Vanja Skuric. Merge: Henrik Rusche
Conflicts:
	applications/solvers/compressible/rhoPorousMRFPimpleFoam/rhoPorousMRFPimpleFoam.C
	applications/solvers/compressible/steadyCompressibleFoam/pEqn.H
	applications/solvers/compressible/steadyCompressibleMRFFoam/createFields.H
	applications/solvers/compressible/steadyCompressibleMRFFoam/pEqn.H
	applications/solvers/compressible/steadyCompressibleSRFFoam/createFields.H
	applications/solvers/compressible/steadyCompressibleSRFFoam/pEqn.H
	applications/solvers/coupled/MRFPorousFoam/createFields.H
	applications/solvers/coupled/pUCoupledFoam/createFields.H
	tutorials/compressible/rhoPorousMRFPimpleFoam/mixerVessel2D/system/fvSolution
	tutorials/compressible/steadyCompressibleFoam/transonicBump/system/fvSchemes
	tutorials/compressible/steadyCompressibleSRFFoam/bentBlade/system/fvSchemes
	tutorials/compressible/steadyCompressibleSRFFoam/simpleBlade/system/fvSchemes
	tutorials/coupled/conjugateHeatFoam/conjugateCavity/system/fvSchemes
	tutorials/coupled/conjugateHeatSimpleFoam/conjugateCavity/system/fvSchemes
	tutorials/equationReader/equationReaderDemo/pitzDaily/system/fvSchemes
	tutorials/heatTransfer/buoyantBoussinesqPisoFoam/hotRoom/system/fvSchemes
	tutorials/heatTransfer/buoyantBoussinesqSimpleFoam/hotRoom/system/fvSchemes
	tutorials/heatTransfer/buoyantBoussinesqSimpleFoam/iglooWithFridges/system/fvSchemes
	tutorials/incompressible/pUCoupledFoam/backwardFacingStepLaminar/system/fvSchemes
	tutorials/incompressible/pUCoupledFoam/backwardFacingStepTurbulent/system/fvSchemes
	tutorials/incompressible/pUCoupledFoam/cavity/system/fvSchemes
	tutorials/viscoelastic/viscoelasticFluidFoam/XPP_SE/system/fvSchemes
2016-05-29 16:37:07 +02:00

131 lines
3.8 KiB
C

/*---------------------------------------------------------------------------*\
========= |
\\ / 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
This file is part of foam-extend.
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
Free Software Foundation, either version 3 of the License, or (at your
option) any later version.
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
along with foam-extend. If not, see <http://www.gnu.org/licenses/>.
Application
conjugateHeatFoam
Description
Transient solver for buoyancy-driven turbulent flow of incompressible
Newtonian fluids with conjugate heat transfer, complex heat conduction
and radiation
\*---------------------------------------------------------------------------*/
#include "fvCFD.H"
#include "coupledFvMatrices.H"
#include "regionCouplePolyPatch.H"
#include "radiationModel.H"
#include "thermalModel.H"
#include "singlePhaseTransportModel.H"
#include "RASModel.H"
#include "pisoControl.H"
#include "simpleControl.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
int main(int argc, char *argv[])
{
# include "setRootCase.H"
# include "createTime.H"
# include "createFluidMesh.H"
# include "createSolidMesh.H"
pisoControl piso(mesh);
# include "readGravitationalAcceleration.H"
# include "createFields.H"
# include "createSolidFields.H"
# include "initContinuityErrs.H"
# include "createTimeControls.H"
# include "CourantNo.H"
# include "setInitialDeltaT.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
Info<< "\nStarting time loop\n" << endl;
while (runTime.loop())
{
Info<< "Time = " << runTime.timeName() << nl << endl;
# include "readTimeControls.H"
# include "CourantNo.H"
# include "setDeltaT.H"
// Detach patches
# include "detachPatches.H"
# include "UEqn.H"
p_rgh.storePrevIter();
while (piso.correct())
{
# include "pEqn.H"
}
// Update turbulent quantities
turbulence->correct();
radiation->correct();
// Update thermal conductivity in the fluid
kappaEff = rho*Cp*(turbulence->nu()/Pr + turbulence->nut()/Prt);
// Update thermal conductivity in the solid
solidThermo.correct();
ksolid = solidThermo.k();
rhoCpsolid.oldTime();
rhoCpsolid = solidThermo.rho()*solidThermo.C();
// Coupled patches
# include "attachPatches.H"
kappaEff.correctBoundaryConditions();
ksolid.correctBoundaryConditions();
// Interpolate to the faces and add thermal resistance
surfaceScalarField ksolidf = fvc::interpolate(ksolid);
solidThermo.modifyResistance(ksolidf);
# include "solveEnergy.H"
// Update density according to Boussinesq approximation
rhok = 1.0 - beta*(T - TRef);
rhok.correctBoundaryConditions();
runTime.write();
Info<< "ExecutionTime = "
<< runTime.elapsedCpuTime()
<< " s\n\n" << endl;
}
Info<< "End\n" << endl;
return(0);
}
// ************************************************************************* //