FreeFOAM The Cross-Platform CFD Toolkit
inverseDistanceDiffusivity.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 1991-2010 OpenCFD Ltd.
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 
28 #include <meshTools/patchWave.H>
29 #include <OpenFOAM/HashSet.H>
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37  defineTypeNameAndDebug(inverseDistanceDiffusivity, 0);
38 
40  (
41  motionDiffusivity,
42  inverseDistanceDiffusivity,
43  Istream
44  );
45 }
46 
47 
48 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
49 
50 Foam::inverseDistanceDiffusivity::inverseDistanceDiffusivity
51 (
52  const fvMotionSolver& mSolver,
53  Istream& mdData
54 )
55 :
56  uniformDiffusivity(mSolver, mdData),
57  patchNames_(mdData)
58 {
59  correct();
60 }
61 
62 
63 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
64 
66 {}
67 
68 
69 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
70 
72 {
73  const polyMesh& mesh = mSolver().mesh();
74 
75  labelHashSet patchSet(mesh.boundaryMesh().patchSet(patchNames_));
76 
77  if (patchSet.size())
78  {
79  return tmp<scalarField>
80  (
81  new scalarField(patchWave(mesh, patchSet, false).distance())
82  );
83  }
84  else
85  {
86  return tmp<scalarField>(new scalarField(mesh.nCells(), 1.0));
87  }
88 }
89 
90 
92 {
93  const fvMesh& mesh = mSolver().mesh();
94 
96  (
97  IOobject
98  (
99  "y",
100  mesh.time().timeName(),
101  mesh
102  ),
103  mesh,
104  dimless,
105  zeroGradientFvPatchScalarField::typeName
106  );
107  y_.internalField() = y();
109 
110  faceDiffusivity_ = 1.0/fvc::interpolate(y_);
111 }
112 
113 
114 // ************************ vim: set sw=4 sts=4 et: ************************ //