OpenFOAM logo
The Open Source CFD Toolbox
  Source Guide OpenCFD Solutions Contact OpenFOAM

potentialFoam.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-2009 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 the
00013     Free Software Foundation; either version 2 of the License, or (at your
00014     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, write to the Free Software Foundation,
00023     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
00024 
00025 Application
00026     potentialFoam
00027 
00028 Description
00029     Simple potential flow solver which can be used to generate starting fields
00030     for full Navier-Stokes codes.
00031 
00032 \*---------------------------------------------------------------------------*/
00033 
00034 #include "fvCFD.H"
00035 
00036 
00037 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00038 
00039 int main(int argc, char *argv[])
00040 {
00041 
00042     argList::validOptions.insert("writep", "");
00043 
00044 #   include "setRootCase.H"
00045 
00046 #   include "createTime.H"
00047 #   include "createMesh.H"
00048 #   include "createFields.H"
00049 #   include "readSIMPLEControls.H"
00050 
00051 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00052 
00053     Info<< nl << "Calculating potential flow" << endl;
00054 
00055     adjustPhi(phi, U, p);
00056 
00057     for (int nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++)
00058     {
00059         fvScalarMatrix pEqn
00060         (
00061             fvm::laplacian
00062             (
00063                 dimensionedScalar
00064                 (
00065                     "1",
00066                     dimTime/p.dimensions()*dimensionSet(0, 2, -2, 0, 0),
00067                     1
00068                 ),
00069                 p
00070             )
00071          ==
00072             fvc::div(phi)
00073         );
00074 
00075         pEqn.setReference(pRefCell, pRefValue);
00076         pEqn.solve();
00077 
00078         if (nonOrth == nNonOrthCorr)
00079         {
00080             phi -= pEqn.flux();
00081         }
00082     }
00083 
00084     Info<< "continuity error = "
00085         << mag(fvc::div(phi))().weightedAverage(mesh.V()).value()
00086         << endl;
00087 
00088     U = fvc::reconstruct(phi);
00089     U.correctBoundaryConditions();
00090 
00091     Info<< "Interpolated U error = "
00092         << (sqrt(sum(sqr((fvc::interpolate(U) & mesh.Sf()) - phi)))
00093           /sum(mesh.magSf())).value()
00094         << endl;
00095 
00096     // Force the write
00097     U.write();
00098     phi.write();
00099 
00100     if (args.optionFound("writep"))
00101     {
00102         p.write();
00103     }
00104 
00105     Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
00106         << "  ClockTime = " << runTime.elapsedClockTime() << " s"
00107         << nl << endl;
00108 
00109     Info<< "End\n" << endl;
00110 
00111     return 0;
00112 }
00113 
00114 
00115 // ************************************************************************* //
Copyright © 2000-2009 OpenCFD Ltd