Loading...
Searching...
No Matches
WallCollisionRecordI.H
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-2015 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
29
30template<class Type>
32(
33 const vector& pRel,
34 scalar radius
35)
36{
37 scalar magpRel_= mag(pRel_);
38
39 scalar magpRel = mag(pRel);
40
41 // Using the new data as the acceptance criterion
42 scalar cosAcceptanceAngle = magpRel/radius;
43
44 if (cosAcceptanceAngle > errorCosAngle)
45 {
46 Info<< "pRel_ " << pRel_ << " " << magpRel_ << nl
47 << "pRel " << pRel << " " << magpRel << nl
48 << "unit vector dot product " << (pRel & pRel_)/(magpRel_*magpRel)
49 << nl << "cosAcceptanceAngle " << cosAcceptanceAngle
50 << endl;
51
53 << "Problem with matching WallCollisionRecord." << nl
54 << "The given radius, " << radius << ", is smaller than distance "
55 << "to the relative position of the WallInteractionSite, "
56 << magpRel << nl
57 << abort(FatalError);
58 }
59
60 // Are the test and recorded pRel (relative position vectors)
61 // aligned to within the calculated tolerance?
62 bool matched = (pRel & pRel_)/(magpRel_*magpRel) > cosAcceptanceAngle;
63
64 if (matched)
65 {
66 pRel_ = pRel;
67 }
69 return matched;
70}
71
72
73template<class Type>
74inline const Foam::vector&
77 return pRel_;
78}
79
80
81template<class Type>
82inline const Type&
84{
85 return data_;
86}
87
88
89template<class Type>
91{
92 return data_;
93}
94
95
96template<class Type>
98{
99 return accessed_;
100}
101
102
103template<class Type>
105{
106 accessed_ = true;
107}
108
109
110template<class Type>
112{
113 accessed_ = false;
114}
115
116
117// * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * * //
118
119template<class Type>
120inline bool Foam::operator==
121(
122 const WallCollisionRecord<Type>& a,
123 const WallCollisionRecord<Type>& b
124)
125{
126 return
127 (
128 a.accessed_ == b.accessed_
129 && a.pRel_ == b.pRel_
130 && a.data_ == b.data_
131 );
132}
133
134
135template<class Type>
136inline bool Foam::operator!=
137(
138 const WallCollisionRecord<Type>& a,
139 const WallCollisionRecord<Type>& b
140)
141{
142 return !(a == b);
143}
144
145
146// ************************************************************************* //
Record of a collision between the particle holding the record and a wall face at the position relativ...
const vector & pRel() const
Return the pRel data.
bool match(const vector &pRel, scalar radius)
void setAccessed()
Set the accessed property of the record to accessed.
const Type & collisionData() const
Return access to the collision data.
bool accessed() const
Return the accessed status of the record.
void setUnaccessed()
Set the accessed property of the record to unaccessed.
static const scalar errorCosAngle
Tolerance for detecting seriously erroneous wall matches.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
messageStream Info
Information stream (stdout output on master, null elsewhere).
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
errorManip< error > abort(error &err)
Definition errorManip.H:139
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
Vector< scalar > vector
Definition vector.H:57
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
volScalarField & b