/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | foam-extend: Open Source CFD \\ / O peration | Version: 4.0 \\ / 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 solidDisplacementFoam Description Transient segregated finite-volume solver of linear-elastic, small-strain deformation of a solid body, with optional thermal diffusion and thermal stresses. Simple linear elasticity structural analysis code. Solves for the displacement vector field D, also generating the stress tensor field sigma. \*---------------------------------------------------------------------------*/ #include "fvCFD.H" #include "Switch.H" // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // int main(int argc, char *argv[]) { # include "setRootCase.H" # include "createTime.H" # include "createMesh.H" # include "readMechanicalProperties.H" # include "readThermalProperties.H" # include "readSolidDisplacementFoamControls.H" # include "createFields.H" // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // Info<< "\nCalculating displacement field\n" << endl; while (runTime.loop()) { Info<< "Iteration: " << runTime.value() << nl << endl; # include "readSolidDisplacementFoamControls.H" int iCorr = 0; scalar initialResidual = 0; do { if (thermalStress) { volScalarField& T = Tptr(); solve ( fvm::ddt(T) == fvm::laplacian(DT, T) ); } { fvVectorMatrix DEqn ( fvm::d2dt2(D) == fvm::laplacian(2*mu + lambda, D, "laplacian(DD,D)") + divSigmaExp ); if (thermalStress) { const volScalarField& T = Tptr(); DEqn += fvc::grad(threeKalpha*T); } //DEqn.setComponentReference(1, 0, vector::X, 0); //DEqn.setComponentReference(1, 0, vector::Z, 0); initialResidual = DEqn.solve().initialResidual(); if (!compactNormalStress) { divSigmaExp = fvc::div(DEqn.flux()); } } { volTensorField gradD = fvc::grad(D); sigmaD = mu*twoSymm(gradD) + (lambda*I)*tr(gradD); if (compactNormalStress) { divSigmaExp = fvc::div ( sigmaD - (2*mu + lambda)*gradD, "div(sigmaD)" ); } else { divSigmaExp += fvc::div(sigmaD); } } } while (initialResidual > convergenceTolerance && ++iCorr < nCorr); # include "calculateStress.H" Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s" << " ClockTime = " << runTime.elapsedClockTime() << " s" << nl << endl; } Info<< "End\n" << endl; return 0; } // ************************************************************************* //