WallCollisionRecordI.H
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) 2011-2020 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 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
27 
28 template<class Type>
30 (
31  const vector& pRel,
32  scalar radius
33 )
34 {
35  scalar magpRel_= mag(pRel_);
36 
37  scalar magpRel = mag(pRel);
38 
39  // Using the new data as the acceptance criterion
40  scalar cosAcceptanceAngle = magpRel/radius;
41 
42  if (cosAcceptanceAngle > errorCosAngle)
43  {
44  Info<< "pRel_ " << pRel_ << " " << magpRel_ << nl
45  << "pRel " << pRel << " " << magpRel << nl
46  << "unit vector dot product " << (pRel & pRel_)/(magpRel_*magpRel)
47  << nl << "cosAcceptanceAngle " << cosAcceptanceAngle
48  << endl;
49 
51  << "Problem with matching WallCollisionRecord." << nl
52  << "The given radius, " << radius << ", is smaller than distance "
53  << "to the relative position of the WallInteractionSite, "
54  << magpRel << nl
55  << abort(FatalError);
56  }
57 
58  // Are the test and recorded pRel (relative position vectors)
59  // aligned to within the calculated tolerance?
60  bool matched = (pRel & pRel_)/(magpRel_*magpRel) > cosAcceptanceAngle;
61 
62  if (matched)
63  {
64  pRel_ = pRel;
65  }
66 
67  return matched;
68 }
69 
70 
71 template<class Type>
72 inline const Foam::vector&
74 {
75  return pRel_;
76 }
77 
78 
79 template<class Type>
80 inline const Type&
82 {
83  return data_;
84 }
85 
86 
87 template<class Type>
89 {
90  return data_;
91 }
92 
93 
94 template<class Type>
96 {
97  return accessed_;
98 }
99 
100 
101 template<class Type>
103 {
104  accessed_ = true;
105 }
106 
107 
108 template<class Type>
110 {
111  accessed_ = false;
112 }
113 
114 
115 // * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * * //
116 
117 template<class Type>
118 inline bool Foam::operator==
119 (
120  const WallCollisionRecord<Type>& a,
122 )
123 {
124  return
125  (
126  a.accessed_ == b.accessed_
127  && a.pRel_ == b.pRel_
128  && a.data_ == b.data_
129  );
130 }
131 
132 
133 template<class Type>
134 inline bool Foam::operator!=
135 (
136  const WallCollisionRecord<Type>& a,
138 )
139 {
140  return !(a == b);
141 }
142 
143 
144 // ************************************************************************* //
void setUnaccessed()
Set the accessed property of the record to unaccessed.
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
bool accessed() const
Return the accessed status of the record.
const vector & pRel() const
Return the pRel data.
const dimensionedScalar b
Wien displacement law constant: default SI units: [m K].
Definition: createFields.H:27
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
errorManip< error > abort(error &err)
Definition: errorManip.H:131
static const char nl
Definition: Ostream.H:260
Record of a collision between the particle holding the record and a wall face at the position relativ...
bool match(const vector &pRel, scalar radius)
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
const Type & collisionData() const
Return access to the collision data.
void setAccessed()
Set the accessed property of the record to accessed.