basicThermo.C 12 KB
Newer Older
1
2
3
4
/*---------------------------------------------------------------------------*\
  =========                 |
  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
   \\    /   O peration     |
5
    \\  /    A nd           | Copyright (C) 2011-2017 OpenFOAM Foundation
6
     \\/     M anipulation  | Copyright (C) 2017 OpenCFD Ltd.
7
8
9
10
-------------------------------------------------------------------------------
License
    This file is part of OpenFOAM.

11
12
13
14
    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 3 of the License, or
    (at your option) any later version.
15
16
17
18
19
20
21

    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
22
    along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
23
24
25
26

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

#include "basicThermo.H"
27
28
29
30
31
32
33
34
#include "zeroGradientFvPatchFields.H"
#include "fixedEnergyFvPatchScalarField.H"
#include "gradientEnergyFvPatchScalarField.H"
#include "mixedEnergyFvPatchScalarField.H"
#include "fixedJumpFvPatchFields.H"
#include "fixedJumpAMIFvPatchFields.H"
#include "energyJumpFvPatchScalarField.H"
#include "energyJumpAMIFvPatchScalarField.H"
35

36

37
/* * * * * * * * * * * * * * * private static data * * * * * * * * * * * * * */
38
39
40

namespace Foam
{
41
    defineTypeNameAndDebug(basicThermo, 0);
42
    defineRunTimeSelectionTable(basicThermo, fvMesh);
43
}
44

45
const Foam::word Foam::basicThermo::dictName("thermophysicalProperties");
46

47
48
49
50
51

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

Foam::wordList Foam::basicThermo::heBoundaryBaseTypes()
{
52
    const volScalarField::Boundary& tbf =
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
        this->T_.boundaryField();

    wordList hbt(tbf.size(), word::null);

    forAll(tbf, patchi)
    {
        if (isA<fixedJumpFvPatchScalarField>(tbf[patchi]))
        {
            const fixedJumpFvPatchScalarField& pf =
                dynamic_cast<const fixedJumpFvPatchScalarField&>(tbf[patchi]);

            hbt[patchi] = pf.interfaceFieldType();
        }
        else if (isA<fixedJumpAMIFvPatchScalarField>(tbf[patchi]))
        {
            const fixedJumpAMIFvPatchScalarField& pf =
                dynamic_cast<const fixedJumpAMIFvPatchScalarField&>
                (
                    tbf[patchi]
                );

            hbt[patchi] = pf.interfaceFieldType();
        }
    }

    return hbt;
}


Foam::wordList Foam::basicThermo::heBoundaryTypes()
{
84
    const volScalarField::Boundary& tbf =
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
        this->T_.boundaryField();

    wordList hbt = tbf.types();

    forAll(tbf, patchi)
    {
        if (isA<fixedValueFvPatchScalarField>(tbf[patchi]))
        {
            hbt[patchi] = fixedEnergyFvPatchScalarField::typeName;
        }
        else if
        (
            isA<zeroGradientFvPatchScalarField>(tbf[patchi])
         || isA<fixedGradientFvPatchScalarField>(tbf[patchi])
        )
        {
            hbt[patchi] = gradientEnergyFvPatchScalarField::typeName;
        }
        else if (isA<mixedFvPatchScalarField>(tbf[patchi]))
        {
            hbt[patchi] = mixedEnergyFvPatchScalarField::typeName;
        }
        else if (isA<fixedJumpFvPatchScalarField>(tbf[patchi]))
        {
            hbt[patchi] = energyJumpFvPatchScalarField::typeName;
        }
        else if (isA<fixedJumpAMIFvPatchScalarField>(tbf[patchi]))
        {
            hbt[patchi] = energyJumpAMIFvPatchScalarField::typeName;
        }
        else if (tbf[patchi].type() == "energyRegionCoupledFvPatchScalarField")
        {
            hbt[patchi] = "energyRegionCoupledFvPatchScalarField";
        }
    }

    return hbt;
}


125
126
// * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //

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
152
153
154
155
156
157
158
159
160
161
Foam::volScalarField& Foam::basicThermo::lookupOrConstruct
(
    const fvMesh& mesh,
    const char* name
) const
{
    if (!mesh.objectRegistry::foundObject<volScalarField>(name))
    {
        volScalarField* fPtr
        (
            new volScalarField
            (
                IOobject
                (
                    name,
                    mesh.time().timeName(),
                    mesh,
                    IOobject::MUST_READ,
                    IOobject::AUTO_WRITE
                ),
                mesh
            )
        );

        // Transfer ownership of this object to the objectRegistry
        fPtr->store(fPtr);
    }

    return const_cast<volScalarField&>
    (
        mesh.objectRegistry::lookupObject<volScalarField>(name)
    );
}


162
163
164
165
166
167
168
169
170
void Foam::basicThermo::lookupAndCheckout(const char* name) const
{
    if (db().foundObject<volScalarField>(name))
    {
         db().checkOut(*db()[name]);
    }
}


171
172
173
174
175
Foam::basicThermo::basicThermo
(
    const fvMesh& mesh,
    const word& phaseName
)
176
:
177
178
179
180
    IOdictionary
    (
        IOobject
        (
181
            phasePropertyName(dictName, phaseName),
182
183
184
185
186
187
188
            mesh.time().constant(),
            mesh,
            IOobject::MUST_READ_IF_MODIFIED,
            IOobject::NO_WRITE
        )
    ),

189
190
191
    phaseName_(phaseName),

    p_(lookupOrConstruct(mesh, "p")),
192
193
194
195
196

    T_
    (
        IOobject
        (
197
            phasePropertyName("T"),
198
199
200
201
202
203
204
205
206
207
208
209
            mesh.time().timeName(),
            mesh,
            IOobject::MUST_READ,
            IOobject::AUTO_WRITE
        ),
        mesh
    ),

    alpha_
    (
        IOobject
        (
210
            phasePropertyName("thermo:alpha"),
211
212
213
214
215
216
217
            mesh.time().timeName(),
            mesh,
            IOobject::NO_READ,
            IOobject::NO_WRITE
        ),
        mesh,
        dimensionSet(1, -1, -1, 0, 0)
218
    ),
219

220
221
    dpdt_(lookupOrDefault<Switch>("dpdt", true))
{}
222
223


224
225
226
Foam::basicThermo::basicThermo
(
    const fvMesh& mesh,
227
228
    const dictionary& dict,
    const word& phaseName
229
)
230
:
231
232
233
234
    IOdictionary
    (
        IOobject
        (
235
            phasePropertyName(dictName, phaseName),
236
237
238
239
240
241
242
243
            mesh.time().constant(),
            mesh,
            IOobject::NO_READ,
            IOobject::NO_WRITE
        ),
        dict
    ),

244
245
246
    phaseName_(phaseName),

    p_(lookupOrConstruct(mesh, "p")),
247
248
249
250
251

    T_
    (
        IOobject
        (
252
            phasePropertyName("T"),
253
254
255
256
257
258
259
260
261
262
263
264
            mesh.time().timeName(),
            mesh,
            IOobject::MUST_READ,
            IOobject::AUTO_WRITE
        ),
        mesh
    ),

    alpha_
    (
        IOobject
        (
265
            phasePropertyName("thermo:alpha"),
266
267
268
269
270
271
272
273
            mesh.time().timeName(),
            mesh,
            IOobject::NO_READ,
            IOobject::NO_WRITE
        ),
        mesh,
        dimensionSet(1, -1, -1, 0, 0)
    )
274
275
276
{}


277
278
// * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //

279
280
Foam::autoPtr<Foam::basicThermo> Foam::basicThermo::New
(
281
282
    const fvMesh& mesh,
    const word& phaseName
283
284
)
{
285
    return New<basicThermo>(mesh, phaseName);
286
287
288
}


289
290
// * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //

291
Foam::basicThermo::~basicThermo()
292
293
294
{
    lookupAndCheckout("p");
}
295
296


297
298
// * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //

299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
const Foam::basicThermo& Foam::basicThermo::lookupThermo
(
    const fvPatchScalarField& pf
)
{
    if (pf.db().foundObject<basicThermo>(dictName))
    {
        return pf.db().lookupObject<basicThermo>(dictName);
    }
    else
    {
        HashTable<const basicThermo*> thermos =
            pf.db().lookupClass<basicThermo>();

        for
        (
            HashTable<const basicThermo*>::iterator iter = thermos.begin();
            iter != thermos.end();
            ++iter
        )
        {
            if
            (
322
323
                &(iter()->he().internalField())
              == &(pf.internalField())
324
325
326
327
328
329
330
331
332
333
334
            )
            {
                return *iter();
            }
        }
    }

    return pf.db().lookupObject<basicThermo>(dictName);
}


335
336
void Foam::basicThermo::validate
(
337
    const string& app,
338
339
340
    const word& a
) const
{
341
    if (!(he().name() == phasePropertyName(a)))
342
    {
343
        FatalErrorInFunction
344
            << "Supported energy type is " << phasePropertyName(a)
345
346
347
348
349
350
351
            << ", thermodynamics package provides " << he().name()
            << exit(FatalError);
    }
}

void Foam::basicThermo::validate
(
352
    const string& app,
353
354
355
356
    const word& a,
    const word& b
) const
{
357
358
359
360
361
362
363
    if
    (
       !(
            he().name() == phasePropertyName(a)
         || he().name() == phasePropertyName(b)
        )
    )
364
    {
365
        FatalErrorInFunction
366
367
            << "Supported energy types are " << phasePropertyName(a)
            << " and " << phasePropertyName(b)
368
369
370
371
372
373
374
            << ", thermodynamics package provides " << he().name()
            << exit(FatalError);
    }
}

void Foam::basicThermo::validate
(
375
    const string& app,
376
377
378
379
380
381
382
383
    const word& a,
    const word& b,
    const word& c
) const
{
    if
    (
       !(
384
385
386
            he().name() == phasePropertyName(a)
         || he().name() == phasePropertyName(b)
         || he().name() == phasePropertyName(c)
387
388
389
        )
    )
    {
390
        FatalErrorInFunction
391
392
393
            << "Supported energy types are " << phasePropertyName(a)
            << ", " << phasePropertyName(b)
            << " and " << phasePropertyName(c)
394
395
396
397
398
399
400
            << ", thermodynamics package provides " << he().name()
            << exit(FatalError);
    }
}

void Foam::basicThermo::validate
(
401
    const string& app,
402
403
404
405
406
407
408
409
410
    const word& a,
    const word& b,
    const word& c,
    const word& d
) const
{
    if
    (
       !(
411
412
413
414
            he().name() == phasePropertyName(a)
         || he().name() == phasePropertyName(b)
         || he().name() == phasePropertyName(c)
         || he().name() == phasePropertyName(d)
415
416
417
        )
    )
    {
418
        FatalErrorInFunction
419
420
421
422
            << "Supported energy types are " << phasePropertyName(a)
            << ", " << phasePropertyName(b)
            << ", " << phasePropertyName(c)
            << " and " << phasePropertyName(d)
423
424
425
426
427
428
            << ", thermodynamics package provides " << he().name()
            << exit(FatalError);
    }
}


429
430
431
432
433
434
435
436
Foam::wordList Foam::basicThermo::splitThermoName
(
    const word& thermoName,
    const int nCmpt
)
{
    wordList cmpts(nCmpt);

437
    string::size_type beg=0, end=0, endb=0, endc=0;
438
439
440
441
    int i = 0;

    while
    (
442
443
        (endb = thermoName.find('<', beg)) != string::npos
     || (endc = thermoName.find(',', beg)) != string::npos
444
445
    )
    {
446
447
448
449
450
451
452
453
454
455
456
457
458
        if (endb == string::npos)
        {
            end = endc;
        }
        else if ((endc = thermoName.find(',', beg)) != string::npos)
        {
            end = min(endb, endc);
        }
        else
        {
            end = endb;
        }

459
460
461
462
        if (beg < end)
        {
            cmpts[i] = thermoName.substr(beg, end-beg);
            cmpts[i++].replaceAll(">","");
463
464
465
466
467
468
469

            // If the number of number of components in the name
            // is greater than nCmpt return an empty list
            if (i == nCmpt)
            {
                return wordList::null();
            }
470
471
472
473
        }
        beg = end + 1;
    }

474
475
476
477
478
479
480
    // If the number of number of components in the name is not equal to nCmpt
    // return an empty list
    if (i + 1 != nCmpt)
    {
        return wordList::null();
    }

481
482
483
    if (beg < thermoName.size())
    {
        cmpts[i] = thermoName.substr(beg, string::npos);
484
        cmpts[i].replaceAll(">","");
485
486
487
488
489
490
    }

    return cmpts;
}


491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
Foam::volScalarField& Foam::basicThermo::p()
{
    return p_;
}


const Foam::volScalarField& Foam::basicThermo::p() const
{
    return p_;
}


const Foam::volScalarField& Foam::basicThermo::T() const
{
    return T_;
}


509
510
511
512
513
514
Foam::volScalarField& Foam::basicThermo::T()
{
    return T_;
}


515
516
517
518
519
520
const Foam::volScalarField& Foam::basicThermo::alpha() const
{
    return alpha_;
}


521
522
523
524
525
526
const Foam::scalarField& Foam::basicThermo::alpha(const label patchi) const
{
    return alpha_.boundaryField()[patchi];
}


527
528
529
530
531
532
bool Foam::basicThermo::read()
{
    return regIOobject::read();
}


533
// ************************************************************************* //