1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | foam-extend: Open Source CFD
4 \\ / O peration | Version: 3.2
5 \\ / A nd | Web: http://www.foam-extend.org
6 \\/ M anipulation | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
9 This file is part of foam-extend.
11 foam-extend 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 3 of the License, or (at your
14 option) any later version.
16 foam-extend is distributed in the hope that it will be useful, but
17 WITHOUT ANY WARRANTY; without even the implied warranty of
18 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 General Public License for more details.
21 You should have received a copy of the GNU General Public License
22 along with foam-extend. If not, see <http://www.gnu.org/licenses/>.
28 Transient solver for potential flow with dynamic mesh.
31 Hrvoje Jasak, Wikki Ltd. All rights reserved.
33 \*---------------------------------------------------------------------------*/
36 #include "dynamicFvMesh.H"
37 #include "pisoControl.H"
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 int main(int argc, char *argv[])
43 argList::validOptions.insert("resetU", "");
44 argList::validOptions.insert("writep", "");
46 # include "setRootCase.H"
47 # include "createTime.H"
48 # include "createDynamicFvMesh.H"
50 pisoControl piso(mesh);
52 # include "createFields.H"
53 # include "initTotalVolume.H"
55 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
57 Info<< "\nStarting time loop\n" << endl;
59 while (runTime.loop())
61 # include "checkTotalVolume.H"
63 Info<< "Time = " << runTime.timeName() << nl << endl;
65 bool meshChanged = mesh.update();
66 reduce(meshChanged, orOp<bool>());
68 p.internalField() = 0;
70 if (args.optionFound("resetU"))
72 U.internalField() = vector::zero;
75 # include "volContinuity.H"
76 # include "meshCourantNo.H"
78 // Solve potential flow equations
81 while (piso.correctNonOrthogonal())
92 dimTime/p.dimensions()*dimensionSet(0, 2, -2, 0, 0),
101 pEqn.setReference(pRefCell, pRefValue);
104 if (piso.finalNonOrthogonalIter())
114 Info<< "continuity error = "
115 << mag(fvc::div(phi))().weightedAverage(mesh.V()).value()
118 U = fvc::reconstruct(phi);
119 U.correctBoundaryConditions();
121 Info<< "Interpolated U error = "
122 << (sqrt(sum(sqr((fvc::interpolate(U) & mesh.Sf()) - phi)))
123 /sum(mesh.magSf())).value()
126 // Calculate velocity magnitude
128 volScalarField magU = mag(U);
130 Info<< "mag(U): max: " << gMax(magU.internalField())
131 << " min: " << gMin(magU.internalField()) << endl;
134 if (args.optionFound("writep"))
136 // Find reference patch
140 // Go through all velocity patches and find the one that fixes
141 // velocity to the largest value
143 forAll (U.boundaryField(), patchI)
145 const fvPatchVectorField& Upatch = U.boundaryField()[patchI];
147 if (Upatch.fixesValue())
149 // Calculate mean velocity
150 scalar u = sum(mag(Upatch));
151 label patchSize = Upatch.size();
153 reduce(u, sumOp<scalar>());
154 reduce(patchSize, sumOp<label>());
158 scalar curMag = u/patchSize;
160 if (curMag > maxMagU)
172 // Calculate reference pressure
173 const fvPatchVectorField& Upatch = U.boundaryField()[refPatch];
174 const fvPatchScalarField& pPatch = p.boundaryField()[refPatch];
176 scalar patchE = sum(mag(pPatch + 0.5*magSqr(Upatch)));
177 label patchSize = Upatch.size();
179 reduce(patchE, sumOp<scalar>());
180 reduce(patchSize, sumOp<label>());
182 scalar e = patchE/patchSize;
184 Info<< "Using reference patch " << refPatch
185 << " with mag(U) = " << maxMagU
186 << " p + 0.5*U^2 = " << e << endl;
188 p.internalField() = e - 0.5*magSqr(U.internalField());
189 p.correctBoundaryConditions();
195 Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
196 << " ClockTime = " << runTime.elapsedClockTime() << " s"
200 Info<< "End\n" << endl;
206 // ************************************************************************* //