forces.C 9.42 KB
Newer Older
1
2
3
4
/*---------------------------------------------------------------------------*\
  =========                 |
  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
   \\    /   O peration     |
5
    \\  /    A nd           | Copyright (C) 1991-2008 OpenCFD Ltd.
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
     \\/     M anipulation  |
-------------------------------------------------------------------------------
License
    This file is part of OpenFOAM.

    OpenFOAM is free software; you can redistribute it and/or modify it
    under the terms of the GNU General Public License as published by the
    Free Software Foundation; either version 2 of the License, or (at your
    option) any later version.

    OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
    ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
    FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
    for more details.

    You should have received a copy of the GNU General Public License
    along with OpenFOAM; if not, write to the Free Software Foundation,
    Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA

\*---------------------------------------------------------------------------*/

#include "forces.H"
#include "volFields.H"
#include "dictionary.H"
#include "Time.H"

#include "incompressible/singlePhaseTransportModel/singlePhaseTransportModel.H"
33
34
#include "incompressible/RASModel/RASModel.H"
#include "incompressible/LESModel/LESModel.H"
35
#include "basicThermo.H"
36
37
#include "compressible/RASModel/RASModel.H"
#include "compressible/LESModel/LESModel.H"
38
39
40
41
42
43
44
45
46
47
48
49
50


// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //

namespace Foam
{
    defineTypeNameAndDebug(forces, 0);
}

// * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //

Foam::tmp<Foam::volSymmTensorField> Foam::forces::devRhoReff() const
{
51
    if (obr_.foundObject<compressible::RASModel>("RASProperties"))
52
    {
53
54
        const compressible::RASModel& ras
            = obr_.lookupObject<compressible::RASModel>("RASProperties");
55
56
57

        return ras.devRhoReff();
    }
58
    else if (obr_.foundObject<incompressible::RASModel>("RASProperties"))
59
    {
60
61
        const incompressible::RASModel& ras
            = obr_.lookupObject<incompressible::RASModel>("RASProperties");
62
63
64

        return rhoRef_*ras.devReff();
    }
65
    else if (obr_.foundObject<compressible::LESModel>("LESProperties"))
66
    {
67
68
        const compressible::LESModel& les =
        obr_.lookupObject<compressible::LESModel>("LESProperties");
69
70
71

        return les.devRhoBeff();
    }
72
    else if (obr_.foundObject<incompressible::LESModel>("LESProperties"))
73
    {
74
75
        const incompressible::LESModel& les
            = obr_.lookupObject<incompressible::LESModel>("LESProperties");
76
77
78
79
80
81
82
83

        return rhoRef_*les.devBeff();
    }
    else if (obr_.foundObject<basicThermo>("thermophysicalProperties"))
    {
        const basicThermo& thermo =
             obr_.lookupObject<basicThermo>("thermophysicalProperties");

Andrew Heather's avatar
Andrew Heather committed
84
        const volVectorField& U = obr_.lookupObject<volVectorField>(UName_);
85
86
87
88
89
90
91
92
93
94
95
96

        return -thermo.mu()*dev(twoSymm(fvc::grad(U)));
    }
    else if
    (
        obr_.foundObject<singlePhaseTransportModel>("transportProperties")
    )
    {
        const singlePhaseTransportModel& laminarT =
            obr_.lookupObject<singlePhaseTransportModel>
            ("transportProperties");

Andrew Heather's avatar
Andrew Heather committed
97
        const volVectorField& U = obr_.lookupObject<volVectorField>(UName_);
98
99
100
101
102
103
104
105
106
107

        return -rhoRef_*laminarT.nu()*dev(twoSymm(fvc::grad(U)));
    }
    else if (obr_.foundObject<dictionary>("transportProperties"))
    {
        const dictionary& transportProperties =
             obr_.lookupObject<dictionary>("transportProperties");

        dimensionedScalar nu(transportProperties.lookup("nu"));

Andrew Heather's avatar
Andrew Heather committed
108
        const volVectorField& U = obr_.lookupObject<volVectorField>(UName_);
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151

        return -rhoRef_*nu*dev(twoSymm(fvc::grad(U)));
    }
    else
    {
        FatalErrorIn("forces::devRhoReff()")
            << "No valid model for viscous stress calculation."
            << exit(FatalError);

        return volSymmTensorField::null();
    }
}


Foam::scalar Foam::forces::rho(const volScalarField& p) const
{
    if (p.dimensions() == dimPressure)
    {
        return 1.0;
    }
    else
    {
        return rhoRef_;
    }
}


// * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //

Foam::forces::forces
(
    const word& name,
    const objectRegistry& obr,
    const dictionary& dict,
    const bool loadFromFiles
)
:
    name_(name),
    obr_(obr),
    active_(true),
    log_(false),
    patchSet_(),
    pName_(""),
Andrew Heather's avatar
Andrew Heather committed
152
    UName_(""),
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
    rhoRef_(0),
    CofR_(vector::zero),
    forcesFilePtr_(NULL)
{
    // Check if the available mesh is an fvMesh otherise deactivate
    if (!isA<fvMesh>(obr_))
    {
        active_ = false;
        WarningIn
        (
            "forces::forces(const objectRegistry& obr, const dictionary& dict)"
        )   << "No fvMesh available, deactivating."
            << endl;
    }

    read(dict);
}


// * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //

Foam::forces::~forces()
{}


// * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //

void Foam::forces::read(const dictionary& dict)
{
    if (active_)
    {
        log_ = dict.lookupOrDefault<Switch>("log", false);

        const fvMesh& mesh = refCast<const fvMesh>(obr_);

        patchSet_ =
            mesh.boundaryMesh().patchSet(wordList(dict.lookup("patches")));

        // Optional entries U and p
        pName_ = dict.lookupOrDefault<word>("pName", "p");
Andrew Heather's avatar
Andrew Heather committed
193
        UName_ = dict.lookupOrDefault<word>("UName", "U");
194

Andrew Heather's avatar
Andrew Heather committed
195
        // Check whether UName and pName exists, if not deactivate forces
196
197
        if
        (
Andrew Heather's avatar
Andrew Heather committed
198
            !obr_.foundObject<volVectorField>(UName_)
199
200
201
202
203
         || !obr_.foundObject<volScalarField>(pName_)
        )
        {
            active_ = false;
            WarningIn("void forces::read(const dictionary& dict)")
Andrew Heather's avatar
Andrew Heather committed
204
                << "Could not find " << UName_ << " or "
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
                << pName_ << " in database." << nl
                << "    De-activating forces."
                << endl;
        }

        // Reference density needed for incompressible calculations
        rhoRef_ = readScalar(dict.lookup("rhoInf"));

        // Centre of rotation for moment calculations
        CofR_ = dict.lookup("CofR");
    }
}


void Foam::forces::makeFile()
{
    // Create the forces file if not already created
    if (!forcesFilePtr_.valid())
    {
        if (debug)
        {
            Info<< "Creating forces file." << endl;
        }

        // File update
        if (Pstream::master())
        {
            fileName forcesDir;
            if (Pstream::parRun())
            {
                // Put in undecomposed case (Note: gives problems for
                // distributed data running)
                forcesDir =
                obr_.time().path()/".."/name_/obr_.time().timeName();
            }
            else
            {
                forcesDir = obr_.time().path()/name_/obr_.time().timeName();
            }

            // Create directory if does not exist.
            mkDir(forcesDir);

            // Open new file at start up
            forcesFilePtr_.reset(new OFstream(forcesDir/(type() + ".dat")));

            // Add headers to output data
            writeFileHeader();
        }
    }
}


void Foam::forces::writeFileHeader()
{
    if (forcesFilePtr_.valid())
    {
        forcesFilePtr_()
            << "# Time" << tab
            << "forces(pressure, viscous) moment(pressure, viscous)"
            << endl;
    }
}


Andrew Heather's avatar
Andrew Heather committed
270
271
272
273
274
void Foam::forces::execute()
{
    // Do nothing - only valid on write
}

275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
void Foam::forces::write()
{
    if (active_)
    {
        // Create the forces file if not already created
        makeFile();

        forcesMoments fm = calcForcesMoment();

        if (Pstream::master())
        {
            forcesFilePtr_() << obr_.time().value() << tab << fm << endl;

            if (log_)
            {
                Info<< "forces output:" << nl
                    << "    forces(pressure, viscous)" << fm.first() << nl
                    << "    moment(pressure, viscous)" << fm.second() << nl
                    << endl;
            }
        }
    }
}


Foam::forces::forcesMoments Foam::forces::calcForcesMoment() const
{
Andrew Heather's avatar
Andrew Heather committed
302
    const volVectorField& U = obr_.lookupObject<volVectorField>(UName_);
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
    const volScalarField& p = obr_.lookupObject<volScalarField>(pName_);

    const fvMesh& mesh = U.mesh();

    forcesMoments fm
    (
        pressureViscous(vector::zero, vector::zero),
        pressureViscous(vector::zero, vector::zero)
    );

    const surfaceVectorField::GeometricBoundaryField& Sfb =
        mesh.Sf().boundaryField();

    tmp<volSymmTensorField> tdevRhoReff = devRhoReff();
    const volSymmTensorField::GeometricBoundaryField& devRhoReffb
        = tdevRhoReff().boundaryField();

    forAllConstIter(labelHashSet, patchSet_, iter)
    {
        label patchi = iter.key();

        vectorField Md = mesh.C().boundaryField()[patchi] - CofR_;

        vectorField pf =
            mesh.Sf().boundaryField()[patchi]*p.boundaryField()[patchi];

        fm.first().first() += rho(p)*sum(pf);
        fm.second().first() += rho(p)*sum(Md ^ pf);

        vectorField vf = Sfb[patchi] & devRhoReffb[patchi];

        fm.first().second() += sum(vf);
        fm.second().second() += sum(Md ^ vf);
    }

    reduce(fm, sumOp());

    return fm;
}


// ************************************************************************* //