1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright held by original author
7 -------------------------------------------------------------------------------
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
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
26 Hrvoje Jasak, Wikki Ltd. All rights reserved
28 \*---------------------------------------------------------------------------*/
30 #include "surfaceFields.H"
31 #include "rasVofForceAndTorqueFunctionObject.H"
32 #include "addToRunTimeSelectionTable.H"
34 #include "incompressible/incompressibleTwoPhaseMixture/twoPhaseMixture.H"
35 #include "incompressible/RAS/RASModel/RASModel.H"
36 #include "incompressible/LES/LESModel/LESModel.H"
38 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42 defineTypeNameAndDebug(rasVofForceAndTorqueFunctionObject, 0);
44 addToRunTimeSelectionTable
47 rasVofForceAndTorqueFunctionObject,
53 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
55 Foam::tmp<Foam::volSymmTensorField>
56 Foam::rasVofForceAndTorqueFunctionObject::devRhoReff() const
59 time_.lookupObject<fvMesh>(regionName_);
61 if (mesh.foundObject<incompressible::RASModel>("RASProperties"))
63 const incompressible::RASModel& ras =
64 mesh.lookupObject<incompressible::RASModel>("RASProperties");
68 else if (mesh.foundObject<incompressible::LESModel>("LESProperties"))
70 const incompressible::LESModel& les =
71 mesh.lookupObject<incompressible::LESModel>("LESProperties");
75 else if (mesh.foundObject<twoPhaseMixture>("transportProperties"))
77 const twoPhaseMixture& twoPhaseProperties =
78 mesh.lookupObject<twoPhaseMixture>("transportProperties");
80 const volVectorField& U = mesh.lookupObject<volVectorField>("U");
82 return twoPhaseProperties.nu()*dev(twoSymm(fvc::grad(U)));
86 FatalErrorIn("floatingBody::devRhoReff()")
87 << "No valid model for viscous stress calculation."
90 return volSymmTensorField::null();
95 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
97 Foam::rasVofForceAndTorqueFunctionObject::
98 rasVofForceAndTorqueFunctionObject
102 const dictionary& dict
105 functionObject(name),
108 regionName_(polyMesh::defaultRegion),
109 patchNames_(dict.lookup("patches")),
110 origin_(dict.lookup("origin")),
111 of_(time_.path()/word(dict.lookup("file")))
113 if (dict.found("region"))
115 dict.lookup("region") >> regionName_;
118 Info<< "Creating rasVofForceAndTorqueFunctionObject for patches "
119 << patchNames_ << endl;
123 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
125 bool Foam::rasVofForceAndTorqueFunctionObject::start()
131 bool Foam::rasVofForceAndTorqueFunctionObject::execute()
134 time_.lookupObject<fvMesh>(regionName_);
136 const volScalarField& alpha = mesh.lookupObject<volScalarField>("alpha1");
137 const volScalarField& p = mesh.lookupObject<volScalarField>("p");
138 const volVectorField& U = mesh.lookupObject<volVectorField>("U");
140 const surfaceVectorField::GeometricBoundaryField& Sfb =
141 mesh.Sf().boundaryField();
143 tmp<volSymmTensorField> tdevRhoReff = devRhoReff();
144 const volSymmTensorField::GeometricBoundaryField& devRhoReffb
145 = tdevRhoReff().boundaryField();
147 vectorField pressureForces(patchNames_.size(), vector::zero);
148 vectorField pressureMoments(patchNames_.size(), vector::zero);
150 vectorField viscousForces(patchNames_.size(), vector::zero);
151 vectorField viscousMoments(patchNames_.size(), vector::zero);
153 forAll (patchNames_, patchI)
155 const label curPatch =
156 mesh.boundaryMesh().findPatchID(patchNames_[patchI]);
162 // Pressure component
163 pressureForces[patchI] =
166 alpha.boundaryField()[curPatch]*
167 p.boundaryField()[curPatch]*Sfb[curPatch]
170 pressureMoments[patchI] =
173 (mesh.Cf().boundaryField()[curPatch] - origin_)
175 alpha.boundaryField()[curPatch]*
176 p.boundaryField()[curPatch]*Sfb[curPatch]
181 viscousForces[patchI] =
184 alpha.boundaryField()[curPatch]*
185 (Sfb[curPatch] & devRhoReffb[curPatch])
188 viscousMoments[patchI] =
191 (mesh.Cf().boundaryField()[curPatch] - origin_)
193 alpha.boundaryField()[curPatch]*
194 (Sfb[curPatch] & devRhoReffb[curPatch])
202 "bool rasVofForceAndTorqueFunctionObject::execute()"
203 ) << "Patch named " << patchNames_[patchI] << " not found. "
204 << "Available patches are: "
205 << mesh.boundaryMesh().names()
206 << abort(FatalError);
210 vectorField totalForces = pressureForces + viscousForces;
211 vectorField totalMoments = pressureMoments + viscousMoments;
213 Info<< "Pressure forces = " << pressureForces << nl
214 << "Pressure moments = " << pressureMoments << nl
215 << "Viscous forces = " << viscousForces << nl
216 << "Viscous moments = " << viscousMoments << nl
217 << "Total forces = " << totalForces << nl
218 << "Total moments = " << totalMoments << endl;
220 of_ << time_.value() << tab
221 << sum(totalForces).x() << tab
222 << sum(totalForces).y() << tab
223 << sum(totalForces).z() << tab
224 << sum(totalMoments).x() << tab
225 << sum(totalMoments).y() << tab
226 << sum(totalMoments).z() << endl;
228 if (patchNames_.size() > 1)
230 Info<< "sum pressure forces = " << sum(pressureForces) << nl
231 << "sum pressure moments = " << sum(pressureMoments) << nl
232 << "sum viscous forces = " << sum(viscousForces) << nl
233 << "sum viscous moments = " << sum(viscousMoments) << nl
234 << "sum total forces = " << sum(totalForces) << nl
235 << "sum total moments = " << sum(totalMoments) << nl << endl;
246 bool Foam::rasVofForceAndTorqueFunctionObject::read(const dictionary& dict)
248 patchNames_ = wordList(dict.lookup("patches"));
249 origin_ = point(dict.lookup("origin"));
254 // ************************************************************************* //