remoteI.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) 2022-2024 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 #include "remote.H"
27 
28 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
29 
31 :
32  proci(-1),
33  elementi(-1)
34 {}
35 
36 
37 inline Foam::remote::remote(const label p, const label e)
38 :
39  proci(p),
40  elementi(e)
41 {}
42 
43 
45 {
46  is >> *this;
47 }
48 
49 
50 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
51 
52 inline Foam::remote Foam::remote::firstProcOp::operator()
53 (
54  const remote& a,
55  const remote& b
56 ) const
57 {
58  if (b.proci != -1 && (a.proci == -1 || a.proci < b.proci))
59  {
60  return b;
61  }
62  else
63  {
64  return a;
65  }
66 }
67 
68 
69 inline void Foam::remote::firstProcEqOp::operator()
70 (
71  remote& a,
72  const remote& b
73 ) const
74 {
75  if (b.proci != -1 && (a.proci == -1 || a.proci < b.proci))
76  {
77  a = b;
78  }
79 }
80 
81 
82 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
83 
84 inline bool Foam::operator==(const remote& a, const remote& b)
85 {
86  return a.proci == b.proci && a.elementi == b.elementi;
87 }
88 
89 
90 inline bool Foam::operator!=(const remote& a, const remote& b)
91 {
92  return !(a == b);
93 }
94 
95 
96 // * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * * //
97 
99 {
100  return os << p.proci << token::SPACE << p.elementi;
101 }
102 
103 
105 {
106  return is >> p.proci >> p.elementi;
107 }
108 
109 
110 // * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
111 
112 namespace Foam
113 {
114 
115 template<class Cmpt>
116 class Tensor;
117 
118 inline remote transform(const Tensor<scalar>&, const remote& p)
119 {
120  return p;
121 }
122 
123 }
124 
125 
126 // ************************************************************************* //
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:60
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
Struct for keeping processor, element (cell, face, point) index.
Definition: remote.H:57
remote()
Construct null.
Definition: remoteI.H:30
label elementi
Element index.
Definition: remote.H:70
label proci
Processor index.
Definition: remote.H:67
volScalarField & b
Definition: createFields.H:25
Namespace for OpenFOAM.
bool operator!=(const particle &, const particle &)
Definition: particle.C:1210
const doubleScalar e
Definition: doubleScalar.H:106
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
Istream & operator>>(Istream &, pistonPointEdgeData &)
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:504
Ostream & operator<<(Ostream &os, const fvConstraints &constraints)
volScalarField & p