freestreamPressureFvPatchScalarField.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2018 OpenFOAM Foundation
9  Copyright (C) 2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
30 #include "volFields.H"
32 
33 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34 
37 (
38  const fvPatch& p,
40 )
41 :
42  mixedFvPatchScalarField(p, iF),
43  UName_("U")
44 {}
45 
46 
49 (
50  const fvPatch& p,
52  const dictionary& dict
53 )
54 :
55  mixedFvPatchScalarField(p, iF),
56  UName_(dict.getOrDefault<word>("U", "U"))
57 {
58  freestreamValue() = scalarField("freestreamValue", dict, p.size());
59 
60  if (dict.found("value"))
61  {
63  (
64  scalarField("value", dict, p.size())
65  );
66  }
67  else
68  {
70  }
71 
72  refGrad() = Zero;
73  valueFraction() = 0;
74 }
75 
76 
79 (
81  const fvPatch& p,
83  const fvPatchFieldMapper& mapper
84 )
85 :
86  mixedFvPatchScalarField(ptf, p, iF, mapper),
87  UName_(ptf.UName_)
88 {}
89 
90 
93 (
95 )
96 :
97  mixedFvPatchScalarField(wbppsf),
98  UName_(wbppsf.UName_)
99 {}
100 
101 
104 (
107 )
108 :
109  mixedFvPatchScalarField(wbppsf, iF),
110  UName_(wbppsf.UName_)
111 {}
112 
113 
114 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
115 
117 {
118  if (updated())
119  {
120  return;
121  }
122 
123  const Field<vector>& Up =
124  patch().template lookupPatchField<volVectorField, vector>
125  (
126  UName_
127  );
128 
129  valueFraction() = 0.5 + 0.5*(Up & patch().nf())/mag(Up);
130 
132 }
133 
134 
136 {
138  os.writeEntryIfDifferent<word>("U", "U", UName_);
139  freestreamValue().writeEntry("freestreamValue", os);
140  writeEntry("value", os);
141 }
142 
143 
144 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
145 
146 namespace Foam
147 {
149  (
151  freestreamPressureFvPatchScalarField
152  );
153 }
154 
155 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: volSurfaceMapping.H:47
Foam::fvPatchScalarField
fvPatchField< scalar > fvPatchScalarField
Definition: fvPatchFieldsFwd.H:33
volFields.H
Foam::fvPatchField::write
virtual void write(Ostream &) const
Definition: fvPatchField.C:377
Foam::fvPatchField::operator=
virtual void operator=(const UList< Type > &)
Foam::Ostream::writeEntryIfDifferent
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Definition: Ostream.H:244
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:46
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:63
freestreamPressureFvPatchScalarField.H
Foam::Zero
static constexpr const zero Zero
Definition: zero.H:131
Foam::freestreamPressureFvPatchScalarField
This boundary condition provides a free-stream condition for pressure.
Definition: freestreamPressureFvPatchScalarField.H:92
Foam::freestreamPressureFvPatchScalarField::write
virtual void write(Ostream &) const
Definition: freestreamPressureFvPatchScalarField.C:128
Foam::fvPatchField::operator
friend Ostream & operator(Ostream &, const fvPatchField< Type > &)
Foam::freestreamPressureFvPatchScalarField::freestreamValue
const scalarField & freestreamValue() const
Definition: freestreamPressureFvPatchScalarField.H:172
Foam::freestreamPressureFvPatchScalarField::freestreamPressureFvPatchScalarField
freestreamPressureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Definition: freestreamPressureFvPatchScalarField.C:30
Foam::Field< vector >
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:119
Foam::freestreamPressureFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Definition: freestreamPressureFvPatchScalarField.C:109
os
OBJstream os(runTime.globalPath()/outputName)
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Definition: atmBoundaryLayer.C:26
Foam::fvPatchField< Type >::updateCoeffs
virtual void updateCoeffs()
Definition: fvPatchField.C:314
Foam::foamVersion::patch
const std::string patch
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:41
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:52
Foam::makePatchTypeField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:50