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/blockCoupledScalarTransportFoam/blockCoupledScalarTransportFoam.C

157 lines
4.8 KiB
C
Raw Normal View History

2010-09-29 18:57:03 +00:00
/*---------------------------------------------------------------------------*\
========= |
2013-12-11 16:09:41 +00:00
\\ / F ield | foam-extend: Open Source CFD
2010-09-29 18:57:03 +00:00
\\ / O peration |
2013-12-11 16:09:41 +00:00
\\ / A nd | For copyright notice see file Copyright
2010-09-29 18:57:03 +00:00
\\/ M anipulation |
-------------------------------------------------------------------------------
License
2013-12-11 16:09:41 +00:00
This file is part of foam-extend.
2010-09-29 18:57:03 +00:00
2013-12-11 16:09:41 +00:00
foam-extend is free software: you can redistribute it and/or modify it
2010-09-29 18:57:03 +00:00
under the terms of the GNU General Public License as published by the
2013-12-11 16:09:41 +00:00
Free Software Foundation, either version 3 of the License, or (at your
2010-09-29 18:57:03 +00:00
option) any later version.
2013-12-11 16:09:41 +00:00
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.
2010-09-29 18:57:03 +00:00
You should have received a copy of the GNU General Public License
2013-12-11 16:09:41 +00:00
along with foam-extend. If not, see <http://www.gnu.org/licenses/>.
2010-09-29 18:57:03 +00:00
Application
blockCoupledScalarTransportFoam
Description
Solves two coupled transport equations in a block-coupled manner
1) transport equation for a passive scalar
2) diffusion only
This resembles heat exchanging flow through a porous medium
\*---------------------------------------------------------------------------*/
#include "fvCFD.H"
#include "fieldTypes.H"
#include "Time.H"
#include "fvMesh.H"
#include "blockLduSolvers.H"
#include "VectorNFieldTypes.H"
#include "volVectorNFields.H"
#include "blockVectorNMatrices.H"
2010-09-29 18:57:03 +00:00
#include "blockMatrixTools.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
int main(int argc, char *argv[])
{
# include "setRootCase.H"
# include "createTime.H"
# include "createMesh.H"
# include "createFields.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
Info<< "\nCalculating scalar transport\n" << endl;
# include "CourantNo.H"
for (runTime++; !runTime.end(); runTime++)
{
Info<< "Time = " << runTime.timeName() << nl << endl;
# include "readSIMPLEControls.H"
for (int nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++)
{
fvScalarMatrix TEqn
(
fvm::div(phi, T)
- fvm::laplacian(DT, T)
==
alpha*Ts
- fvm::Sp(alpha, T)
);
TEqn.relax();
fvScalarMatrix TsEqn
(
- fvm::laplacian(DTs, Ts)
==
alpha*T
- fvm::Sp(alpha, Ts)
);
TsEqn.relax();
// Prepare block system
BlockLduMatrix<vector2> blockM(mesh);
2011-08-24 10:35:47 +00:00
//- Transfer the coupled interface list for processor/cyclic/etc.
// boundaries
blockM.interfaces() = blockT.boundaryField().blockInterfaces();
2010-09-29 18:57:03 +00:00
// Grab block diagonal and set it to zero
Field<tensor2>& d = blockM.diag().asSquare();
d = tensor2::zero;
// Grab linear off-diagonal and set it to zero
Field<vector2>& l = blockM.lower().asLinear();
2011-08-24 10:35:47 +00:00
Field<vector2>& u = blockM.upper().asLinear();
2010-09-29 18:57:03 +00:00
u = vector2::zero;
l = vector2::zero;
vector2Field& blockX = blockT.internalField();
2010-09-29 18:57:03 +00:00
vector2Field blockB(mesh.nCells(), vector2::zero);
//- Inset equations into block Matrix
blockMatrixTools::insertEquation(0, TEqn, blockM, blockX, blockB);
blockMatrixTools::insertEquation(1, TsEqn, blockM, blockX, blockB);
2011-08-24 10:35:47 +00:00
//- Add off-diagonal terms and remove from block source
2010-09-29 18:57:03 +00:00
forAll(d, i)
{
2011-09-22 14:07:21 +00:00
d[i](0, 1) = -alpha.value()*mesh.V()[i];
d[i](1, 0) = -alpha.value()*mesh.V()[i];
2010-09-29 18:57:03 +00:00
blockB[i][0] -= alpha.value()*blockX[i][1]*mesh.V()[i];
blockB[i][1] -= alpha.value()*blockX[i][0]*mesh.V()[i];
}
//- Block coupled solver call
2010-09-29 18:57:03 +00:00
BlockSolverPerformance<vector2> solverPerf =
BlockLduSolver<vector2>::New
(
2011-08-24 10:35:47 +00:00
blockT.name(),
2010-09-29 18:57:03 +00:00
blockM,
2011-08-24 10:35:47 +00:00
mesh.solutionDict().solver(blockT.name())
2010-09-29 18:57:03 +00:00
)->solve(blockX, blockB);
solverPerf.print();
// Retrieve solution
blockMatrixTools::blockRetrieve(0, T.internalField(), blockX);
blockMatrixTools::blockRetrieve(1, Ts.internalField(), blockX);
T.correctBoundaryConditions();
Ts.correctBoundaryConditions();
}
runTime.write();
}
Info<< "End\n" << endl;
return(0);
}
// ************************************************************************* //