/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | foam-extend: Open Source CFD \\ / O peration | Version: 3.2 \\ / A nd | Web: http://www.foam-extend.org \\/ M anipulation | For copyright notice see file Copyright ------------------------------------------------------------------------------- License This file is part of foam-extend. foam-extend 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. 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. You should have received a copy of the GNU General Public License along with foam-extend. If not, see . Application stressedFoam Description Transient/steady-state segregated finite-volume solver for a solid body. Solves for the displacement vector field U, also generating the stress tensor field sigma. \*---------------------------------------------------------------------------*/ #include "fvCFD.H" #include "rheologyModel.H" #include "componentReferenceList.H" // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // int main(int argc, char *argv[]) { # include "setRootCase.H" # include "createTime.H" # include "createMesh.H" # include "createFields.H" // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // Info<< "\nCalculating displacement field\n" << endl; volTensorField gradU = fvc::grad(U); volScalarField rho = rheology.rho(); while (runTime.loop()) { Info<< "Iteration: " << runTime.timeName() << nl << endl; # include "readStressedFoamControls.H" volScalarField mu = rheology.mu(); volScalarField lambda = rheology.lambda(); int iCorr = 0; scalar initialResidual = 0; do { fvVectorMatrix UEqn ( fvm::d2dt2(rho, U) == fvm::laplacian(2*mu + lambda, U, "laplacian(DU,U)") + fvc::div ( mu*gradU.T() + lambda*(I*tr(gradU)) - (mu + lambda)*gradU, "div(sigma)" ) ); # include "setComponentReference.H" initialResidual = UEqn.solve().initialResidual(); gradU = fvc::grad(U); # include "calculateSigma.H" rheology.correct(); rho = rheology.rho(); mu = rheology.mu(); lambda = rheology.lambda(); } while (initialResidual > convergenceTolerance && ++iCorr < nCorr); # include "calculateStress.H" Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s\n\n" << endl; } Info<< "End\n" << endl; return(0); } // ************************************************************************* //