FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

solidDisplacementFoam.C

Go to the documentation of this file.
00001 /*---------------------------------------------------------------------------*\
00002   =========                 |
00003   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
00004    \\    /   O peration     |
00005     \\  /    A nd           | Copyright (C) 1991-2010 OpenCFD Ltd.
00006      \\/     M anipulation  |
00007 -------------------------------------------------------------------------------
00008 License
00009     This file is part of OpenFOAM.
00010 
00011     OpenFOAM is free software: you can redistribute it and/or modify it
00012     under the terms of the GNU General Public License as published by
00013     the Free Software Foundation, either version 3 of the License, or
00014     (at your option) any later version.
00015 
00016     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
00017     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00019     for more details.
00020 
00021     You should have received a copy of the GNU General Public License
00022     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
00023 
00024 Application
00025     solidDisplacementFoam
00026 
00027 Description
00028     Transient segregated finite-volume solver of linear-elastic,
00029     small-strain deformation of a solid body, with optional thermal
00030     diffusion and thermal stresses.
00031 
00032     Simple linear elasticity structural analysis code.
00033     Solves for the displacement vector field D, also generating the
00034     stress tensor field sigma.
00035 
00036 Usage
00037     - solidDisplacementFoam [OPTION]
00038 
00039     @param -case <dir> \n
00040     Specify the case directory
00041 
00042     @param -parallel \n
00043     Run the case in parallel
00044 
00045     @param -help \n
00046     Display short usage message
00047 
00048     @param -doc \n
00049     Display Doxygen documentation page
00050 
00051     @param -srcDoc \n
00052     Display source code
00053 
00054 \*---------------------------------------------------------------------------*/
00055 
00056 #include <finiteVolume/fvCFD.H>
00057 #include <OpenFOAM/Switch.H>
00058 
00059 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00060 
00061 int main(int argc, char *argv[])
00062 {
00063     #include <OpenFOAM/setRootCase.H>
00064 
00065     #include <OpenFOAM/createTime.H>
00066     #include <OpenFOAM/createMesh.H>
00067     #include "readMechanicalProperties.H"
00068     #include "readThermalProperties.H"
00069     #include "readSolidDisplacementFoamControls.H"
00070     #include "createFields.H"
00071 
00072     // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00073 
00074     Info<< "\nCalculating displacement field\n" << endl;
00075 
00076     while (runTime.loop())
00077     {
00078         Info<< "Iteration: " << runTime.value() << nl << endl;
00079 
00080         #include "readSolidDisplacementFoamControls.H"
00081 
00082         int iCorr = 0;
00083         scalar initialResidual = 0;
00084 
00085         do
00086         {
00087             if (thermalStress)
00088             {
00089                 volScalarField& T = Tptr();
00090                 solve
00091                 (
00092                     fvm::ddt(T) == fvm::laplacian(DT, T)
00093                 );
00094             }
00095 
00096             {
00097                 fvVectorMatrix DEqn
00098                 (
00099                     fvm::d2dt2(D)
00100                  ==
00101                     fvm::laplacian(2*mu + lambda, D, "laplacian(DD,D)")
00102                   + divSigmaExp
00103                 );
00104 
00105                 if (thermalStress)
00106                 {
00107                     const volScalarField& T = Tptr();
00108                     DEqn += fvc::grad(threeKalpha*T);
00109                 }
00110 
00111                 //DEqn.setComponentReference(1, 0, vector::X, 0);
00112                 //DEqn.setComponentReference(1, 0, vector::Z, 0);
00113 
00114                 initialResidual = DEqn.solve().initialResidual();
00115 
00116                 if (!compactNormalStress)
00117                 {
00118                     divSigmaExp = fvc::div(DEqn.flux());
00119                 }
00120             }
00121 
00122             {
00123                 volTensorField gradD = fvc::grad(D);
00124                 sigmaD = mu*twoSymm(gradD) + (lambda*I)*tr(gradD);
00125 
00126                 if (compactNormalStress)
00127                 {
00128                     divSigmaExp = fvc::div
00129                     (
00130                         sigmaD - (2*mu + lambda)*gradD,
00131                         "div(sigmaD)"
00132                     );
00133                 }
00134                 else
00135                 {
00136                     divSigmaExp += fvc::div(sigmaD);
00137                 }
00138             }
00139 
00140         } while (initialResidual > convergenceTolerance && ++iCorr < nCorr);
00141 
00142         #include "calculateStress.H"
00143 
00144         Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
00145             << "  ClockTime = " << runTime.elapsedClockTime() << " s"
00146             << nl << endl;
00147     }
00148 
00149     Info<< "End\n" << endl;
00150 
00151     return 0;
00152 }
00153 
00154 
00155 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines