phaseMap.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2020-2022 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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.
15 
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.
20 
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/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "phaseMap.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace functionObjects
34 {
35  defineTypeNameAndDebug(phaseMap, 0);
36  addToRunTimeSelectionTable(functionObject, phaseMap, dictionary);
37 }
38 }
39 
40 
41 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42 
44 (
45  const word& name,
46  const Time& runTime,
47  const dictionary& dict
48 )
49 :
50  fvMeshFunctionObject(name, runTime, dict),
51  phases_
52  (
53  mesh_.lookupObject<phaseSystem>(phaseSystem::propertiesName).phases()
54  )
55 {}
56 
57 
58 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
59 
61 {}
62 
63 
64 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
65 
67 {
68  return true;
69 }
70 
71 
73 {
74  volScalarField phaseMap
75  (
76  IOobject
77  (
78  IOobject::groupName(phases_[0].member(), "map"),
79  mesh_.time().timeName(),
80  mesh_,
83  ),
84  mesh_,
86  );
87 
88  scalar level = 0;
89 
90  forAll(phases_, i)
91  {
92  phaseMap += level*phases_[i];
93  level += 1;
94  }
95 
96  phaseMap.write();
97 
98  return true;
99 }
100 
101 
102 // ************************************************************************* //
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
addToRunTimeSelectionTable(functionObject, Qdot, dictionary)
const dimensionSet dimless
virtual ~phaseMap()
Destructor.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:372
phaseMap(const word &name, const Time &runTime, const dictionary &)
Construct from Time and dictionary.
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
static word timeName(const scalar, const int precision=curPrecision_)
Return time name of given scalar time.
Definition: Time.C:666
static word groupName(Name name, const word &group)
virtual bool write()
Write the force fields.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
defineTypeNameAndDebug(Qdot, 0)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual bool execute()
Calculate the force fields.
phaseSystem::phaseModelList & phases
Definition: createFields.H:12
const fvMesh & mesh_
Reference to the fvMesh.
Namespace for OpenFOAM.