Loading...
Searching...
No Matches
PatchCollisionDensity.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) 2018 OpenFOAM Foundation
9 Copyright (C) 2020-2023 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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 "Pstream.H"
31#include "ListOps.H"
32#include "ListListOps.H"
33
34// * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
35
36template<class CloudType>
38{
39 const scalarField z(this->owner().mesh().nCells(), Zero);
40
41 {
43 (
44 this->owner().mesh().newIOobject
45 (
47 (
48 this->owner().name(),
49 "collisionDensity"
50 )
51 ),
52 this->owner().mesh(),
54 z,
55 collisionDensity_
56 ).write();
57 }
58
59 {
61 (
62 this->owner().mesh().newIOobject
63 (
65 (
66 this->owner().name(),
67 "collisionDensityRate"
68 )
69 ),
70 this->owner().mesh(),
72 z,
73 (collisionDensity_ - collisionDensity0_)
74 /(this->owner().mesh().time().value() - time0_)
75 ).write();
76 }
77
78 collisionDensity0_ == collisionDensity_;
79 time0_ = this->owner().mesh().time().value();
80}
81
82
83// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
84
85template<class CloudType>
87(
88 const dictionary& dict,
89 CloudType& owner,
90 const word& modelName
91)
92:
93 CloudFunctionObject<CloudType>(dict, owner, modelName, typeName),
94 minSpeed_(dict.getOrDefault<scalar>("minSpeed", -1)),
95 collisionDensity_
96 (
97 this->owner().mesh().boundary(),
98 volScalarField::Internal::null(),
99 fvPatchFieldBase::calculatedType()
100 ),
101 collisionDensity0_
102 (
103 this->owner().mesh().boundary(),
104 volScalarField::Internal::null(),
105 fvPatchFieldBase::calculatedType()
106 ),
107 time0_(this->owner().mesh().time().value())
108{
109 collisionDensity_ == 0;
110 collisionDensity0_ == 0;
111
113 (
114 this->owner().mesh().newIOobject
115 (
116 IOobject::scopedName(this->owner().name(), "collisionDensity"),
118 )
119 );
120
121 if (io.typeHeaderOk<volScalarField>())
122 {
123 const volScalarField collisionDensity(io, this->owner().mesh());
124 collisionDensity_ == collisionDensity.boundaryField();
125 collisionDensity0_ == collisionDensity.boundaryField();
126 }
127}
128
129
130template<class CloudType>
132(
133 const PatchCollisionDensity<CloudType>& ppm
134)
135:
136 CloudFunctionObject<CloudType>(ppm),
137 minSpeed_(ppm.minSpeed_),
138 collisionDensity_(ppm.collisionDensity_),
139 collisionDensity0_(ppm.collisionDensity0_),
140 time0_(ppm.time0_)
141{}
142
143
144// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
145
146template<class CloudType>
148(
149 const parcelType& p,
150 const polyPatch& pp,
151 const typename parcelType::trackingData& td
152)
153{
154 const label patchi = pp.index();
155 const label patchFacei = p.face() - pp.start();
156
157 vector nw, Up;
158 this->owner().patchData(p, pp, nw, Up);
159
160 const scalar speed = (p.U() - Up) & nw;
161 if (speed > minSpeed_)
162 {
163 collisionDensity_[patchi][patchFacei] +=
164 1/this->owner().mesh().magSf().boundaryField()[patchi][patchFacei];
165 }
166
167 return true;
168}
169
170
171// ************************************************************************* //
Various functions to operate on Lists.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Templated cloud function object base class.
CloudFunctionObject(CloudType &owner)
Construct null from owner.
const CloudType & owner() const
Return const access to the owner cloud.
const fvMesh & mesh() const
Return reference to the mesh.
Definition DSMCCloudI.H:37
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
@ MUST_READ
Reading required.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition IOobjectI.H:50
PatchCollisionDensity(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
void write()
Write post-processing info.
virtual bool postPatch(const parcelType &p, const polyPatch &pp, const typename parcelType::trackingData &td)
Post-patch hook.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
Template invariant parts for fvPatchField.
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
const dictionary & dict() const
Return const access to the cloud dictionary.
const word & modelName() const
Return const access to the name of the sub-model.
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
faceListList boundary
dynamicFvMesh & mesh
const auto & io
auto & name
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
DSMCCloud< dsmcParcel > CloudType
const dimensionSet dimless
Dimensionless.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const dimensionSet dimArea(sqr(dimLength))
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
Vector< scalar > vector
Definition vector.H:57
dictionary dict