From 8a1dd2eb9adcd210b61a66b3f85937c8373980ad Mon Sep 17 00:00:00 2001 From: Henry Weller Date: Tue, 30 May 2023 14:59:58 +0100 Subject: [PATCH] solver modules: Moved constructing the face velocity/momentum into preSolve() The velocity boundary conditions are corrected before the construction of the face velocity or momentum but for multi-region cases with interacting velocity boundary conditions this is only possible after all the region solver modules have been constructed so it is better to delay the optional construction of the face velocity/momentum until preSolve(). --- applications/modules/VoFSolver/VoFSolver.C | 42 +++++++++---------- .../incompressibleDenseParticleFluid.C | 38 ++++++++--------- .../incompressibleFluid/incompressibleFluid.C | 42 +++++++++---------- .../modules/isothermalFluid/isothermalFluid.C | 42 +++++++++---------- 4 files changed, 82 insertions(+), 82 deletions(-) diff --git a/applications/modules/VoFSolver/VoFSolver.C b/applications/modules/VoFSolver/VoFSolver.C index 97b2b9bce7..099c960340 100644 --- a/applications/modules/VoFSolver/VoFSolver.C +++ b/applications/modules/VoFSolver/VoFSolver.C @@ -146,27 +146,6 @@ Foam::solvers::VoFSolver::VoFSolver { mesh.schemes().setFluxRequired(p_rgh.name()); - if (mesh.dynamic() || MRF.size()) - { - Info<< "Constructing face momentum Uf" << endl; - - // Ensure the U BCs are up-to-date before constructing Uf - U_.correctBoundaryConditions(); - - Uf = new surfaceVectorField - ( - IOobject - ( - "Uf", - runTime.name(), - mesh, - IOobject::READ_IF_PRESENT, - IOobject::AUTO_WRITE - ), - fvc::interpolate(U_) - ); - } - if (LTS) { Info<< "Using LTS" << endl; @@ -221,6 +200,27 @@ void Foam::solvers::VoFSolver::preSolve() // Read the controls readControls(); + if ((mesh.dynamic() || MRF.size()) && !Uf.valid()) + { + Info<< "Constructing face momentum Uf" << endl; + + // Ensure the U BCs are up-to-date before constructing Uf + U_.correctBoundaryConditions(); + + Uf = new surfaceVectorField + ( + IOobject + ( + "Uf", + runTime.name(), + mesh, + IOobject::READ_IF_PRESENT, + IOobject::AUTO_WRITE + ), + fvc::interpolate(U_) + ); + } + if (transient()) { correctCoNum(); diff --git a/applications/modules/incompressibleDenseParticleFluid/incompressibleDenseParticleFluid.C b/applications/modules/incompressibleDenseParticleFluid/incompressibleDenseParticleFluid.C index b225510a34..75b833f809 100644 --- a/applications/modules/incompressibleDenseParticleFluid/incompressibleDenseParticleFluid.C +++ b/applications/modules/incompressibleDenseParticleFluid/incompressibleDenseParticleFluid.C @@ -230,7 +230,25 @@ incompressibleDenseParticleFluid alphacf = fvc::interpolate(alphac); alphaPhic = alphacf*phic; - if (mesh.dynamic()) + correctCoNum(); +} + + +// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * // + +Foam::solvers::incompressibleDenseParticleFluid:: +~incompressibleDenseParticleFluid() +{} + + +// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * // + +void Foam::solvers::incompressibleDenseParticleFluid::preSolve() +{ + // Read the controls + readControls(); + + if (mesh.dynamic() && !Ucf.valid()) { Info<< "Constructing face momentum Ucf" << endl; @@ -251,24 +269,6 @@ incompressibleDenseParticleFluid ); } - correctCoNum(); -} - - -// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * // - -Foam::solvers::incompressibleDenseParticleFluid:: -~incompressibleDenseParticleFluid() -{} - - -// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * // - -void Foam::solvers::incompressibleDenseParticleFluid::preSolve() -{ - // Read the controls - readControls(); - // Store the particle positions if (mesh.topoChanging()) { diff --git a/applications/modules/incompressibleFluid/incompressibleFluid.C b/applications/modules/incompressibleFluid/incompressibleFluid.C index 3aaa9b424f..1aad253f4b 100644 --- a/applications/modules/incompressibleFluid/incompressibleFluid.C +++ b/applications/modules/incompressibleFluid/incompressibleFluid.C @@ -123,27 +123,6 @@ Foam::solvers::incompressibleFluid::incompressibleFluid(fvMesh& mesh) momentumTransport->validate(); - if (mesh.dynamic() || MRF.size()) - { - Info<< "Constructing face momentum Uf" << endl; - - // Ensure the U BCs are up-to-date before constructing Uf - U_.correctBoundaryConditions(); - - Uf = new surfaceVectorField - ( - IOobject - ( - "Uf", - runTime.name(), - mesh, - IOobject::READ_IF_PRESENT, - IOobject::AUTO_WRITE - ), - fvc::interpolate(U) - ); - } - if (transient()) { correctCoNum(); @@ -186,6 +165,27 @@ void Foam::solvers::incompressibleFluid::preSolve() // Read the controls readControls(); + if ((mesh.dynamic() || MRF.size()) && !Uf.valid()) + { + Info<< "Constructing face momentum Uf" << endl; + + // Ensure the U BCs are up-to-date before constructing Uf + U_.correctBoundaryConditions(); + + Uf = new surfaceVectorField + ( + IOobject + ( + "Uf", + runTime.name(), + mesh, + IOobject::READ_IF_PRESENT, + IOobject::AUTO_WRITE + ), + fvc::interpolate(U) + ); + } + fvModels().preUpdateMesh(); if (transient()) diff --git a/applications/modules/isothermalFluid/isothermalFluid.C b/applications/modules/isothermalFluid/isothermalFluid.C index 119a6840d5..04c5ed53fa 100644 --- a/applications/modules/isothermalFluid/isothermalFluid.C +++ b/applications/modules/isothermalFluid/isothermalFluid.C @@ -223,27 +223,6 @@ Foam::solvers::isothermalFluid::isothermalFluid ); } - if (mesh.dynamic() || MRF.size()) - { - Info<< "Constructing face momentum rhoUf" << endl; - - // Ensure the U BCs are up-to-date before constructing Uf - U_.correctBoundaryConditions(); - - rhoUf = new surfaceVectorField - ( - IOobject - ( - "rhoUf", - runTime.name(), - mesh, - IOobject::READ_IF_PRESENT, - IOobject::AUTO_WRITE - ), - fvc::interpolate(rho*U) - ); - } - if (transient()) { correctCoNum(); @@ -292,6 +271,27 @@ void Foam::solvers::isothermalFluid::preSolve() // Read the controls readControls(); + if ((mesh.dynamic() || MRF.size()) && !rhoUf.valid()) + { + Info<< "Constructing face momentum rhoUf" << endl; + + // Ensure the U BCs are up-to-date before constructing Uf + U_.correctBoundaryConditions(); + + rhoUf = new surfaceVectorField + ( + IOobject + ( + "rhoUf", + runTime.name(), + mesh, + IOobject::READ_IF_PRESENT, + IOobject::AUTO_WRITE + ), + fvc::interpolate(rho*U) + ); + } + if (transient()) { correctCoNum();