noBlending.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) 2014-2018 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 "noBlending.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace blendingMethods
34 {
35  defineTypeNameAndDebug(noBlending, 0);
36 
38  (
39  blendingMethod,
40  noBlending,
41  dictionary
42  );
43 }
44 }
45 
46 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
47 
49 (
50  const dictionary& dict,
51  const wordList& phaseNames
52 )
53 :
54  blendingMethod(dict),
55  continuousPhase_(dict.lookup("continuousPhase"))
56 {}
57 
58 
59 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
60 
62 {}
63 
64 
65 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
66 
68 (
69  const phaseModel& phase1,
70  const phaseModel& phase2
71 ) const
72 {
73  const fvMesh& mesh(phase1.mesh());
74 
75  return
76  tmp<volScalarField>
77  (
78  new volScalarField
79  (
80  IOobject
81  (
82  "f",
83  mesh.time().timeName(),
84  mesh
85  ),
86  mesh,
88  (
89  "f",
90  dimless,
91  phase2.name() != continuousPhase_
92  )
93  )
94  );
95 }
96 
97 
99 (
100  const phaseModel& phase1,
101  const phaseModel& phase2
102 ) const
103 {
104  const fvMesh& mesh(phase1.mesh());
105 
106  return
107  tmp<volScalarField>
108  (
109  new volScalarField
110  (
111  IOobject
112  (
113  "f",
114  mesh.time().timeName(),
115  mesh
116  ),
117  mesh,
119  (
120  "f",
121  dimless,
122  phase1.name() == continuousPhase_
123  )
124  )
125  );
126 }
127 
128 
129 // ************************************************************************* //
dictionary dict
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
stressControl lookup("compactNormalStress") >> compactNormalStress
dynamicFvMesh & mesh
virtual tmp< volScalarField > f1(const phaseModel &phase1, const phaseModel &phase2) const
Factor for primary phase.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
List< word > wordList
A List of words.
Definition: fileName.H:54
virtual tmp< volScalarField > f2(const phaseModel &phase1, const phaseModel &phase2) const
Factor for secondary phase.
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
noBlending(const dictionary &dict, const wordList &phaseNames)
Construct from a dictionary and two phases.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.