initial commit for version 1.5.x patch release
[OpenFOAM-1.5.x.git] / applications / solvers / basic / potentialFoam / potentialFoam.C
blob68cc358327767429846c01eac667cba7fb14fb53
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2008 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
11     OpenFOAM is free software; you can redistribute it and/or modify it
12     under the terms of the GNU General Public License as published by the
13     Free Software Foundation; either version 2 of the License, or (at your
14     option) any later version.
16     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
19     for more details.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM; if not, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 Application
26     potentialFoam
28 Description
29     Simple potential flow solver which can be used to generate starting fields
30     for full Navier-Stokes codes.
32 \*---------------------------------------------------------------------------*/
34 #include "fvCFD.H"
37 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39 int main(int argc, char *argv[])
42     argList::validOptions.insert("writep", "");
44 #   include "setRootCase.H"
46 #   include "createTime.H"
47 #   include "createMesh.H"
48 #   include "createFields.H"
49 #   include "readSIMPLEControls.H"
51 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
53     Info<< nl << "Calculating potential flow" << endl;
55     adjustPhi(phi, U, p);
57     for (int nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++)
58     {
59         fvScalarMatrix pEqn
60         (
61             fvm::laplacian
62             (
63                 dimensionedScalar
64                 (
65                     "1",
66                     dimTime/p.dimensions()*dimensionSet(0, 2, -2, 0, 0),
67                     1
68                 ),
69                 p
70             )
71          ==
72             fvc::div(phi)
73         );
75         pEqn.setReference(pRefCell, pRefValue);
76         pEqn.solve();
78         if (nonOrth == nNonOrthCorr)
79         {
80             phi -= pEqn.flux();
81         }
82     }
84     Info<< "continuity error = "
85         << mag(fvc::div(phi))().weightedAverage(mesh.V()).value()
86         << endl;
88     U = fvc::reconstruct(phi);
89     U.correctBoundaryConditions();
91     Info<< "Interpolated U error = "
92         << (sqrt(sum(sqr((fvc::interpolate(U) & mesh.Sf()) - phi)))
93           /sum(mesh.magSf())).value()
94         << endl;
96     // Force the write
97     U.write();
98     phi.write();
100     if (args.options().found("writep"))
101     {
102         p.write();
103     }
105     Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
106         << "  ClockTime = " << runTime.elapsedClockTime() << " s"
107         << nl << endl;
109     Info<< "End\n" << endl;
111     return(0);
115 // ************************************************************************* //