processorCyclicPointPatchField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2013 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 "transformField.H"
28 #include "processorPolyPatch.H"
29 
30 
31 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32 
33 template<class Type>
35 (
36  const pointPatch& p,
38 )
39 :
41  procPatch_(refCast<const processorCyclicPointPatch>(p)),
42  receiveBuf_(0)
43 {}
44 
45 
46 template<class Type>
48 (
49  const pointPatch& p,
51  const dictionary& dict
52 )
53 :
55  procPatch_(refCast<const processorCyclicPointPatch>(p)),
56  receiveBuf_(0)
57 {}
58 
59 
60 template<class Type>
62 (
64  const pointPatch& p,
66  const pointPatchFieldMapper& mapper
67 )
68 :
69  coupledPointPatchField<Type>(ptf, p, iF, mapper),
70  procPatch_(refCast<const processorCyclicPointPatch>(ptf.patch())),
71  receiveBuf_(0)
72 {}
73 
74 
75 template<class Type>
77 (
80 )
81 :
83  procPatch_(refCast<const processorCyclicPointPatch>(ptf.patch())),
84  receiveBuf_(0)
85 {}
86 
87 
88 // * * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * //
89 
90 template<class Type>
92 {}
93 
94 
95 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
96 
97 template<class Type>
99 (
100  const Pstream::commsTypes commsType,
101  Field<Type>& pField
102 ) const
103 {
104  if (Pstream::parRun())
105  {
106  // Get internal field into correct order for opposite side
107  Field<Type> pf
108  (
109  this->patchInternalField
110  (
111  pField,
112  procPatch_.reverseMeshPoints()
113  )
114  );
115 
116  if (commsType == Pstream::nonBlocking)
117  {
118  receiveBuf_.setSize(pf.size());
120  (
121  commsType,
122  procPatch_.neighbProcNo(),
123  reinterpret_cast<char*>(receiveBuf_.begin()),
124  receiveBuf_.byteSize(),
125  procPatch_.tag(),
126  procPatch_.comm()
127  );
128  }
130  (
131  commsType,
132  procPatch_.neighbProcNo(),
133  reinterpret_cast<const char*>(pf.begin()),
134  pf.byteSize(),
135  procPatch_.tag(),
136  procPatch_.comm()
137  );
138  }
139 }
140 
141 
142 template<class Type>
144 (
145  const Pstream::commsTypes commsType,
146  Field<Type>& pField
147 ) const
148 {
149  if (Pstream::parRun())
150  {
151  // If nonblocking data has already been received into receiveBuf_
152  if (commsType != Pstream::nonBlocking)
153  {
154  receiveBuf_.setSize(this->size());
156  (
157  commsType,
158  procPatch_.neighbProcNo(),
159  reinterpret_cast<char*>(receiveBuf_.begin()),
160  receiveBuf_.byteSize(),
161  procPatch_.tag(),
162  procPatch_.comm()
163  );
164  }
165 
166  if (doTransform())
167  {
168  const processorCyclicPolyPatch& ppp =
169  procPatch_.procCyclicPolyPatch();
170  const tensor& forwardT = ppp.forwardT()[0];
171 
172  transform(receiveBuf_, forwardT, receiveBuf_);
173  }
174 
175  // All points are separated
176  this->addToInternalField(pField, receiveBuf_);
177  }
178 }
179 
180 
181 // ************************************************************************* //
processorCyclicPointPatchField(const pointPatch &, const DimensionedField< Type, pointMesh > &)
Construct from patch and internal field.
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.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
A Coupled boundary condition for pointField.
void size(const label)
Override size to be inconsistent with allocated storage.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
runTime write()
virtual void initSwapAddSeparated(const Pstream::commsTypes commsType, Field< Type > &) const
Initialise swap of non-collocated patch point values.
iterator begin()
Return an iterator to begin traversing the UList.
dictionary dict
volScalarField & p
Definition: createFields.H:51
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:465
commsTypes
Types of communications.
Definition: UPstream.H:64
Spatial transformation functions for primitive fields.
Pre-declare SubField and related Field type.
Definition: Field.H:57
virtual void swapAddSeparated(const Pstream::commsTypes commsType, Field< Type > &) const
Complete swap of patch point values and add to local values.
Foam::pointPatchFieldMapper.
bool read(const char *, int32_t &)
Definition: int32IO.C:87
const pointPatch & patch() const
Return patch.
Foam::processorCyclicPointPatchField.
virtual const tensorField & forwardT() const
Return face transformation tensor.