diff --git a/applications/solvers/heatTransfer/solidFoam/Make/files b/applications/solvers/heatTransfer/solidFoam/Make/files new file mode 100644 index 0000000000..5e2a2eeaf0 --- /dev/null +++ b/applications/solvers/heatTransfer/solidFoam/Make/files @@ -0,0 +1,3 @@ +solidFoam.C + +EXE = $(FOAM_APPBIN)/solidFoam diff --git a/applications/solvers/heatTransfer/solidFoam/Make/options b/applications/solvers/heatTransfer/solidFoam/Make/options new file mode 100644 index 0000000000..f8e6b9af40 --- /dev/null +++ b/applications/solvers/heatTransfer/solidFoam/Make/options @@ -0,0 +1,19 @@ +EXE_INC = \ + -DFULLDEBUG -g -O0 \ + -I$(LIB_SRC)/finiteVolume/lnInclude \ + -I$(LIB_SRC)/meshTools/lnInclude \ + -I$(LIB_SRC)/sampling/lnInclude \ + -I$(LIB_SRC)/transportModels/compressible/lnInclude \ + -I$(LIB_SRC)/thermophysicalModels/basic/lnInclude \ + -I$(LIB_SRC)/thermophysicalModels/solidThermo/lnInclude \ + -I$(LIB_SRC)/thermophysicalModels/radiation/lnInclude + +EXE_LIBS = \ + -lfiniteVolume \ + -lfvOptions \ + -lmeshTools \ + -lsampling \ + -lcompressibleTransportModels \ + -lfluidThermophysicalModels \ + -lradiationModels \ + -lspecie diff --git a/applications/solvers/heatTransfer/solidFoam/createFields.H b/applications/solvers/heatTransfer/solidFoam/createFields.H new file mode 100644 index 0000000000..1c6a939db7 --- /dev/null +++ b/applications/solvers/heatTransfer/solidFoam/createFields.H @@ -0,0 +1,102 @@ +Info<< "Reading thermophysical properties\n" << endl; + +autoPtr pThermo(solidThermo::New(mesh)); +solidThermo& thermo = pThermo(); + +tmp trho = thermo.rho(); +const volScalarField& rho = trho(); + +tmp tcp = thermo.Cp(); +const volScalarField& cp = tcp(); + +volScalarField& p = thermo.p(); +volScalarField& h = thermo.he(); + +autoPtr coordinatesPtr; +autoPtr taniAlpha; + +if (!thermo.isotropic()) +{ + Info<< "Adding coordinateSystem\n" << endl; + coordinatesPtr = coordinateSystem::New + ( + mesh, + thermo, + coordinateSystem::typeName_() + ); + + tmp tkappaByCp = thermo.Kappa()/thermo.Cp(); + + taniAlpha.reset + ( + new volSymmTensorField + ( + IOobject + ( + "Anialpha", + runTime.timeName(), + mesh, + IOobject::NO_READ, + IOobject::NO_WRITE + ), + mesh, + dimensionedSymmTensor(tkappaByCp().dimensions(), Zero), + zeroGradientFvPatchSymmTensorField::typeName + ) + ); + volSymmTensorField& aniAlpha = *taniAlpha; + + aniAlpha.primitiveFieldRef() = + coordinatesPtr->transformPrincipal + ( + mesh.cellCentres(), + tkappaByCp() + ); + aniAlpha.correctBoundaryConditions(); +} + + +IOobject betavSolidIO +( + "betavSolid", + runTime.timeName(), + mesh, + IOobject::MUST_READ, + IOobject::AUTO_WRITE +); + +autoPtr betavPtr; +if (betavSolidIO.typeHeaderOk(true)) +{ + betavPtr.reset + ( + new volScalarField + ( + betavSolidIO, + mesh + ) + ); +} +else +{ + betavPtr.reset + ( + new volScalarField + ( + IOobject + ( + "betavSolid", + runTime.timeName(), + mesh, + IOobject::NO_READ, + IOobject::NO_WRITE + ), + mesh, + dimensionedScalar("1", dimless, scalar(1)) + ) + ); +} +const volScalarField& betav = *betavPtr; + +#include "createRadiationModel.H" +#include "createFvOptions.H" diff --git a/applications/solvers/heatTransfer/solidFoam/hEqn.H b/applications/solvers/heatTransfer/solidFoam/hEqn.H new file mode 100644 index 0000000000..587c1028b2 --- /dev/null +++ b/applications/solvers/heatTransfer/solidFoam/hEqn.H @@ -0,0 +1,28 @@ +{ + fvScalarMatrix hEqn + ( + fvm::ddt(betav*rho, h) + - ( + thermo.isotropic() + ? fvm::laplacian(betav*thermo.alpha(), h, "laplacian(alpha,h)") + : fvm::laplacian(betav*taniAlpha(), h, "laplacian(alpha,h)") + ) + == + fvOptions(rho, h) + ); + + hEqn.relax(); + + fvOptions.constrain(hEqn); + + hEqn.solve(); //mesh.solver(h.select(finalIter))); + + fvOptions.correct(h); + + thermo.correct(); + + Info<< "Min/max T:" << min(thermo.T()).value() << ' ' + << max(thermo.T()).value() << endl; + + radiation->correct(); +} diff --git a/applications/solvers/heatTransfer/solidFoam/solidFoam.C b/applications/solvers/heatTransfer/solidFoam/solidFoam.C new file mode 100644 index 0000000000..771eb1fe64 --- /dev/null +++ b/applications/solvers/heatTransfer/solidFoam/solidFoam.C @@ -0,0 +1,112 @@ +/*---------------------------------------------------------------------------*\ + ========= | + \\ / F ield | OpenFOAM: The Open Source CFD Toolbox + \\ / O peration | + \\ / A nd | www.openfoam.com + \\/ M anipulation | +------------------------------------------------------------------------------- + Copyright (C) 2020 OpenCFD Ltd. +------------------------------------------------------------------------------- +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 3 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, see . + +Application + solidFoam + +Group + grpHeatTransferSolvers + +Description + Solver for energy transport and thermodynamics on a solid. + +\*---------------------------------------------------------------------------*/ + +#include "fvCFD.H" +#include "solidThermo.H" +#include "radiationModel.H" +#include "fvOptions.H" +#include "simpleControl.H" +#include "pimpleControl.H" +#include "coordinateSystem.H" + +// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // + +int main(int argc, char *argv[]) +{ + argList::addNote + ( + "Solver for energy transport and thermodynamics on a solid" + ); + + #define NO_CONTROL + #include "postProcess.H" + + #include "addCheckCaseOptions.H" + #include "setRootCaseLists.H" + #include "createTime.H" + #include "createMesh.H" + #include "createFields.H" + + // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // + + Info<< "\nEvolving thermodynamics\n" << endl; + + if (mesh.solutionDict().found("SIMPLE")) + { + simpleControl simple(mesh); + + while (simple.loop()) + { + Info<< "Time = " << runTime.timeName() << nl << endl; + + while (simple.correctNonOrthogonal()) + { + #include "hEqn.H" + } + + runTime.write(); + + runTime.printExecutionTime(Info); + } + } + else + { + pimpleControl pimple(mesh); + + while (runTime.run()) + { + ++runTime; + + Info<< "Time = " << runTime.timeName() << nl << endl; + + while (pimple.correctNonOrthogonal()) + { + #include "hEqn.H" + } + + runTime.write(); + + runTime.printExecutionTime(Info); + } + } + + Info<< "End\n" << endl; + + return 0; +} + + +// ************************************************************************* //