ENH: autoLayerDriver: better layering information message
[OpenFOAM-2.0.x.git] / applications / solvers / basic / potentialFoam / potentialFoam.C
blobe968fa5c16ac8e513a1fb56a260a28c9880b04e2
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
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
13     the Free Software Foundation, either version 3 of the License, or
14     (at your 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, see <http://www.gnu.org/licenses/>.
24 Application
25     potentialFoam
27 Description
28     Simple potential flow solver which can be used to generate starting fields
29     for full Navier-Stokes codes.
31 \*---------------------------------------------------------------------------*/
33 #include "fvCFD.H"
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38 int main(int argc, char *argv[])
40     argList::addBoolOption("writep", "write the final pressure field");
41     argList::addBoolOption
42     (
43         "initialiseUBCs",
44         "initialise U boundary conditions"
45     );
47     #include "setRootCase.H"
48     #include "createTime.H"
49     #include "createMesh.H"
50     #include "readControls.H"
51     #include "createFields.H"
53     // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
55     Info<< nl << "Calculating potential flow" << endl;
57     // Since solver contains no time loop it would never execute
58     // function objects so do it ourselves.
59     runTime.functionObjects().start();
61     adjustPhi(phi, U, p);
63     for (int nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++)
64     {
65         fvScalarMatrix pEqn
66         (
67             fvm::laplacian
68             (
69                 dimensionedScalar
70                 (
71                     "1",
72                     dimTime/p.dimensions()*dimensionSet(0, 2, -2, 0, 0),
73                     1
74                 ),
75                 p
76             )
77          ==
78             fvc::div(phi)
79         );
81         pEqn.setReference(pRefCell, pRefValue);
82         pEqn.solve();
84         if (nonOrth == nNonOrthCorr)
85         {
86             phi -= pEqn.flux();
87         }
88     }
90     Info<< "continuity error = "
91         << mag(fvc::div(phi))().weightedAverage(mesh.V()).value()
92         << endl;
94     U = fvc::reconstruct(phi);
95     U.correctBoundaryConditions();
97     Info<< "Interpolated U error = "
98         << (sqrt(sum(sqr((fvc::interpolate(U) & mesh.Sf()) - phi)))
99           /sum(mesh.magSf())).value()
100         << endl;
102     // Force the write
103     U.write();
104     phi.write();
106     if (args.optionFound("writep"))
107     {
108         p.write();
109     }
111     runTime.functionObjects().end();
113     Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
114         << "  ClockTime = " << runTime.elapsedClockTime() << " s"
115         << nl << endl;
117     Info<< "End\n" << endl;
119     return 0;
123 // ************************************************************************* //