mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
Adds overset discretisation to selected physics: - diffusion : overLaplacianDyMFoam - incompressible steady : overSimpleFoam - incompressible transient : overPimpleDyMFoam - compressible transient: overRhoPimpleDyMFoam - two-phase VOF: overInterDyMFoam The overset method chosen is a parallel, fully implicit implementation whereby the interpolation (from donor to acceptor) is inserted as an adapted discretisation on the donor cells, such that the resulting matrix can be solved using the standard linear solvers. Above solvers come with a set of tutorials, showing how to create and set-up simple simulations from scratch.
126 lines
3.5 KiB
C
126 lines
3.5 KiB
C
/*---------------------------------------------------------------------------*\
|
|
========= |
|
|
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
|
\\ / O peration |
|
|
\\ / A nd | Copyright (C) 2011-2015 OpenFOAM Foundation
|
|
\\/ M anipulation | Copyright (C) 2016-2017 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 <http://www.gnu.org/licenses/>.
|
|
|
|
Application
|
|
laplacianFoam
|
|
|
|
Group
|
|
grpBasicSolvers
|
|
|
|
Description
|
|
Laplace equation solver for a scalar quantity.
|
|
|
|
\heading Solver details
|
|
The solver is applicable to, e.g. for thermal diffusion in a solid. The
|
|
equation is given by:
|
|
|
|
\f[
|
|
\ddt{T} = \div \left( D_T \grad T \right)
|
|
\f]
|
|
|
|
Where:
|
|
\vartable
|
|
T | Scalar field which is solved for, e.g. temperature
|
|
D_T | Diffusion coefficient
|
|
\endvartable
|
|
|
|
\heading Required fields
|
|
\plaintable
|
|
T | Scalar field which is solved for, e.g. temperature
|
|
\endplaintable
|
|
|
|
\*---------------------------------------------------------------------------*/
|
|
|
|
#include "fvCFD.H"
|
|
#include "simpleControl.H"
|
|
#include "dynamicFvMesh.H"
|
|
#include "dynamicOversetFvMesh.H"
|
|
|
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
|
|
|
int main(int argc, char *argv[])
|
|
{
|
|
#include "setRootCase.H"
|
|
|
|
#include "createTime.H"
|
|
#include "createNamedDynamicFvMesh.H"
|
|
|
|
simpleControl simple(mesh);
|
|
|
|
#include "createFields.H"
|
|
|
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
|
|
|
Info<< "\nCorrecting boundary conditions on " << T.name() << nl << endl;
|
|
|
|
runTime++;
|
|
|
|
Info<< "Time = " << runTime.timeName() << nl << endl;
|
|
|
|
Info<< "Reading : ExecutionTime = " << runTime.elapsedCpuTime() << " s"
|
|
<< " ClockTime = " << runTime.elapsedClockTime() << " s"
|
|
<< nl << endl;
|
|
|
|
mesh.update();
|
|
|
|
Info<< "Overset calculation : ExecutionTime = "
|
|
<< runTime.elapsedCpuTime() << " s"
|
|
<< " ClockTime = " << runTime.elapsedClockTime() << " s"
|
|
<< nl << endl;
|
|
|
|
|
|
if (false)
|
|
{
|
|
// Test correctBoundaryConditions
|
|
|
|
// Change the internalField
|
|
component(T.ref(), mesh.C(), 0);
|
|
component(T.ref(), mesh.C(), 1);
|
|
// Interpolate + halo swap
|
|
T.correctBoundaryConditions();
|
|
// Check halo swap
|
|
dynamicOversetFvMesh::checkCoupledBC(T);
|
|
}
|
|
if (true)
|
|
{
|
|
// Test solving
|
|
fvScalarMatrix TEqn(fvm::laplacian(DT, T));
|
|
TEqn.solve();
|
|
}
|
|
|
|
runTime.write();
|
|
|
|
|
|
|
|
Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
|
|
<< " ClockTime = " << runTime.elapsedClockTime() << " s"
|
|
<< nl << endl;
|
|
|
|
Info<< "End\n" << endl;
|
|
|
|
return 0;
|
|
}
|
|
|
|
|
|
// ************************************************************************* //
|