PairCollisionRecordI.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-2017 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 label queryOrigProcOfOther,
34 label queryOrigIdOfOther
35) const
36{
37 return
38 (
39 queryOrigProcOfOther == origProcOfOther()
40 && queryOrigIdOfOther == origIdOfOther()
41 );
42}
43
44
45template<class Type>
47{
48 return mag(origProcOfOther_) - 1;
49}
50
51
52template<class Type>
54{
55 return origIdOfOther_;
56}
57
58
59template<class Type>
60inline const Type&
62{
63 return data_;
64}
65
66
67template<class Type>
69{
70 return data_;
71}
72
73
74template<class Type>
76{
77 return pos0(origProcOfOther_);
78}
79
80
81template<class Type>
83{
84 origProcOfOther_ = origProcOfOther() + 1;
85}
86
87
88template<class Type>
90{
91 origProcOfOther_ = -(origProcOfOther() + 1);
92}
93
94
95// * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * * //
96
97template<class Type>
98inline bool Foam::operator==
99(
102)
103{
104 return
105 (
106 a.origProcOfOther_ == b.origProcOfOther_
107 && a.origIdOfOther_ == b.origIdOfOther_
108 && a.data_ == b.data_
109 );
110}
111
112
113template<class Type>
114inline bool Foam::operator!=
115(
118)
119{
120 return !(a == b);
121}
122
123
124// ************************************************************************* //
Record of a collision between the particle holding the record and the particle with the stored id.
bool match(label queryOrigProcOfOther, label queryOrigIdOfOther) const
label origIdOfOther() const
Return the origIdOfOther data.
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.
label origProcOfOther() const
Return the origProcOfOther data.
dimensionedScalar pos0(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
volScalarField & b
Definition: createFields.H:27