Files
OpenFOAM-12/applications/solvers/foamMultiRun/foamMultiRun.C
Henry Weller 4d63b39e3e foamMultiRun: Added automatic region prefixing to the Info statements in the log
e.g. for the rivuletBox case the output for a time-step now looks like:

film  Courant Number mean: 0.0003701330848 max: 0.1862204919
panel Diffusion Number mean: 0.007352456305 max: 0.1276468109
box   Courant Number mean: 0.006324172752 max: 0.09030825997
      deltaT = 0.001550908752
      Time = 0.08294s

film  diagonal:  Solving for alpha, Initial residual = 0, Final residual = 0, No Iterations 0
film  diagonal:  Solving for alpha, Initial residual = 0, Final residual = 0, No Iterations 0
box   diagonal:  Solving for rho, Initial residual = 0, Final residual = 0, No Iterations 0
film  DILUPBiCGStab:  Solving for Ux, Initial residual = 0.009869417958, Final residual = 2.132619614e-11, No Iterations 2
film  DILUPBiCGStab:  Solving for Uy, Initial residual = 0.0002799662756, Final residual = 6.101011285e-12, No Iterations 1
film  DILUPBiCGStab:  Solving for Uz, Initial residual = 1, Final residual = 1.854120599e-12, No Iterations 2
box   DILUPBiCGStab:  Solving for Ux, Initial residual = 0.004071057403, Final residual = 4.79249226e-07, No Iterations 1
box   DILUPBiCGStab:  Solving for Uy, Initial residual = 0.006370817152, Final residual = 9.606673696e-07, No Iterations 1
box   DILUPBiCGStab:  Solving for Uz, Initial residual = 0.0158299327, Final residual = 2.104129791e-06, No Iterations 1
film  DILUPBiCGStab:  Solving for e, Initial residual = 0.0002888908396, Final residual = 2.301587523e-11, No Iterations 1
panel GAMG:  Solving for e, Initial residual = 0.00878508958, Final residual = 7.807579738e-12, No Iterations 1
box   DILUPBiCGStab:  Solving for h, Initial residual = 0.004403989559, Final residual = 1.334113552e-06, No Iterations 1
film  DILUPBiCGStab:  Solving for alpha, Initial residual = 0.0002760406755, Final residual = 2.267583256e-14, No Iterations 1
film  time step continuity errors : sum local = 9.01334987e-12, global = 2.296671859e-13, cumulative = 1.907846466e-08
box   GAMG:  Solving for p_rgh, Initial residual = 0.002842335602, Final residual = 1.036572819e-05, No Iterations 4
box   diagonal:  Solving for rho, Initial residual = 0, Final residual = 0, No Iterations 0
box   time step continuity errors : sum local = 4.538744531e-07, global = 1.922637799e-08, cumulative = -6.612579497e-09
box   GAMG:  Solving for p_rgh, Initial residual = 1.283128787e-05, Final residual = 7.063185653e-07, No Iterations 2
box   diagonal:  Solving for rho, Initial residual = 0, Final residual = 0, No Iterations 0
box   time step continuity errors : sum local = 3.069629869e-08, global = 3.780547824e-10, cumulative = -6.234524715e-09
      ExecutionTime = 19.382601 s  ClockTime = 20 s

film  Courant Number mean: 0.0003684434169 max: 0.1840342756
panel Diffusion Number mean: 0.007352456305 max: 0.1276468109
box   Courant Number mean: 0.006292704463 max: 0.09016861809
      deltaT = 0.001550908752
      Time = 0.0844909s

where each line printed by each region solver is prefixed by the region name.
Global messages for the time-step and time are just prefixed with spaces to
align them with the region output.
2023-03-08 10:59:13 +00:00

187 lines
5.1 KiB
C++

/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2022-2023 OpenFOAM Foundation
\\/ 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 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
foamMultiRun
Description
Loads and executes an OpenFOAM solver modules for each region of a
multiregion simulation e.g. for conjugate heat transfer.
The region solvers are specified in the \c regionSolvers dictionary entry in
\c controlDict, containing a list of pairs of region and solver names,
e.g. for a two region case with one fluid region named
liquid and one solid region named tubeWall:
\verbatim
regionSolvers
{
liquid fluid;
tubeWall solid;
}
\endverbatim
The \c regionSolvers entry is a dictionary to support name substitutions to
simplify the specification of a single solver type for a set of
regions, e.g.
\verbatim
fluidSolver fluid;
solidSolver solid;
regionSolvers
{
tube1 $fluidSolver;
tubeWall1 solid;
tube2 $fluidSolver;
tubeWall2 solid;
tube3 $fluidSolver;
tubeWall3 solid;
}
\endverbatim
Uses the flexible PIMPLE (PISO-SIMPLE) solution for time-resolved and
pseudo-transient and steady simulations.
Usage
\b foamMultiRun [OPTION]
- \par -libs '(\"lib1.so\" ... \"libN.so\")'
Specify the additional libraries loaded
Example usage:
- To update and run a \c chtMultiRegion case add the following entries to
the controlDict:
\verbatim
application foamMultiRun;
regionSolvers
{
fluid fluid;
solid solid;
}
\endverbatim
then execute \c foamMultiRun
\*---------------------------------------------------------------------------*/
#include "argList.H"
#include "regionSolvers.H"
#include "pimpleMultiRegionControl.H"
#include "setDeltaT.H"
using namespace Foam;
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
int main(int argc, char *argv[])
{
#include "setRootCase.H"
#include "createTime.H"
// Create the region meshes and solvers
regionSolvers solvers(runTime);
// Create the outer PIMPLE loop and control structure
pimpleMultiRegionControl pimple(runTime, solvers);
// Set the initial time-step
setDeltaT(runTime, solvers);
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
Info<< "\nStarting time loop\n" << endl;
while (pimple.run(runTime))
{
forAll(solvers, i)
{
solvers[i].preSolve();
}
solvers.setGlobalPrefix();
// Adjust the time-step according to the solver maxDeltaT
adjustDeltaT(runTime, solvers);
runTime++;
Info<< "Time = " << runTime.userTimeName() << nl << endl;
// Multi-region PIMPLE corrector loop
while (pimple.loop())
{
forAll(solvers, i)
{
solvers[i].moveMesh();
}
forAll(solvers, i)
{
solvers[i].prePredictor();
}
forAll(solvers, i)
{
solvers[i].momentumPredictor();
}
while (pimple.correctEnergy())
{
forAll(solvers, i)
{
solvers[i].thermophysicalPredictor();
}
}
forAll(solvers, i)
{
solvers[i].pressureCorrector();
}
forAll(solvers, i)
{
solvers[i].postCorrector();
}
}
forAll(solvers, i)
{
solvers[i].postSolve();
}
solvers.setGlobalPrefix();
runTime.write();
Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
<< " ClockTime = " << runTime.elapsedClockTime() << " s"
<< nl << endl;
}
Info<< "End\n" << endl;
return 0;
}
// ************************************************************************* //