processorCyclicPointPatchField.C
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 
27 #include "processorPolyPatch.H"
28 
29 
30 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31 
32 template<class Type>
34 (
35  const pointPatch& p,
37 )
38 :
40  procPatch_(refCast<const processorCyclicPointPatch>(p)),
41  receiveBuf_(0)
42 {}
43 
44 
45 template<class Type>
47 (
48  const pointPatch& p,
50  const dictionary& dict
51 )
52 :
54  procPatch_(refCast<const processorCyclicPointPatch>(p)),
55  receiveBuf_(0)
56 {}
57 
58 
59 template<class Type>
61 (
63  const pointPatch& p,
65  const pointPatchFieldMapper& mapper
66 )
67 :
68  coupledPointPatchField<Type>(ptf, p, iF, mapper),
69  procPatch_(refCast<const processorCyclicPointPatch>(ptf.patch())),
70  receiveBuf_(0)
71 {}
72 
73 
74 template<class Type>
76 (
79 )
80 :
82  procPatch_(refCast<const processorCyclicPointPatch>(ptf.patch())),
83  receiveBuf_(0)
84 {}
85 
86 
87 // * * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * //
88 
89 template<class Type>
91 {}
92 
93 
94 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
95 
96 template<class Type>
98 (
99  const Pstream::commsTypes commsType,
100  Field<Type>& pField
101 ) const
102 {
103  if (Pstream::parRun())
104  {
105  // Get internal field into correct order for opposite side
106  Field<Type> pf
107  (
108  this->patchInternalField
109  (
110  pField,
111  procPatch_.reverseMeshPoints()
112  )
113  );
114 
115  if (commsType == Pstream::commsTypes::nonBlocking)
116  {
117  receiveBuf_.setSize(pf.size());
119  (
120  commsType,
121  procPatch_.neighbProcNo(),
122  reinterpret_cast<char*>(receiveBuf_.begin()),
123  receiveBuf_.byteSize(),
124  procPatch_.tag(),
125  procPatch_.comm()
126  );
127  }
129  (
130  commsType,
131  procPatch_.neighbProcNo(),
132  reinterpret_cast<const char*>(pf.begin()),
133  pf.byteSize(),
134  procPatch_.tag(),
135  procPatch_.comm()
136  );
137  }
138 }
139 
140 
141 template<class Type>
143 (
144  const Pstream::commsTypes commsType,
145  Field<Type>& pField
146 ) const
147 {
148  if (Pstream::parRun())
149  {
150  // If nonblocking data has already been received into receiveBuf_
151  if (commsType != Pstream::commsTypes::nonBlocking)
152  {
153  receiveBuf_.setSize(this->size());
155  (
156  commsType,
157  procPatch_.neighbProcNo(),
158  reinterpret_cast<char*>(receiveBuf_.begin()),
159  receiveBuf_.byteSize(),
160  procPatch_.tag(),
161  procPatch_.comm()
162  );
163  }
164 
165  procPatch_.procCyclicPolyPatch().transform().transform
166  (
167  receiveBuf_,
168  receiveBuf_
169  );
170 
171  // All points are separated
172  this->addToInternalField(pField, receiveBuf_);
173  }
174 }
175 
176 
177 // ************************************************************************* //
Foam::processorCyclicPointPatchField.
dictionary dict
processorCyclicPointPatchField(const pointPatch &, const DimensionedField< Type, pointMesh > &)
Construct from patch and internal field.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
commsTypes
Types of communications.
Definition: UPstream.H:64
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
Foam::pointPatchFieldMapper.
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
virtual void initSwapAddSeparated(const Pstream::commsTypes commsType, Field< Type > &) const
Initialise swap of non-collocated patch point values.
Pre-declare SubField and related Field type.
Definition: Field.H:56
const pointPatch & patch() const
Return patch.
iterator begin()
Return an iterator to begin traversing the UList.
Definition: UListI.H:216
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
Basic pointPatch represents a set of points from the mesh.
Definition: pointPatch.H:56
std::streamsize byteSize() const
Return the binary size in number of characters of the UList.
Definition: UList.C:100
A Coupled boundary condition for pointField.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
volScalarField & p
virtual void swapAddSeparated(const Pstream::commsTypes commsType, Field< Type > &) const
Complete swap of patch point values and add to local values.