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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
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 
29 #include "transformField.H"
30 #include "processorPolyPatch.H"
31 
32 
33 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34 
35 template<class Type>
37 (
38  const pointPatch& p,
40 )
41 :
43  procPatch_(refCast<const processorCyclicPointPatch>(p)),
44  receiveBuf_(0)
45 {}
46 
47 
48 template<class Type>
50 (
51  const pointPatch& p,
53  const dictionary& dict
54 )
55 :
57  procPatch_(refCast<const processorCyclicPointPatch>(p, dict)),
58  receiveBuf_(0)
59 {}
60 
61 
62 template<class Type>
64 (
66  const pointPatch& p,
68  const pointPatchFieldMapper& mapper
69 )
70 :
71  coupledPointPatchField<Type>(ptf, p, iF, mapper),
72  procPatch_(refCast<const processorCyclicPointPatch>(ptf.patch())),
73  receiveBuf_(0)
74 {}
75 
76 
77 template<class Type>
79 (
82 )
83 :
85  procPatch_(refCast<const processorCyclicPointPatch>(ptf.patch())),
86  receiveBuf_(0)
87 {}
88 
89 
90 // * * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * //
91 
92 template<class Type>
94 {}
95 
96 
97 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
98 
99 template<class Type>
101 (
102  const Pstream::commsTypes commsType,
103  Field<Type>& pField
104 ) const
105 {
106  if (Pstream::parRun())
107  {
108  // Get internal field into correct order for opposite side
109  Field<Type> pf
110  (
111  this->patchInternalField
112  (
113  pField,
114  procPatch_.reverseMeshPoints()
115  )
116  );
117 
118  if (commsType == Pstream::commsTypes::nonBlocking)
119  {
120  receiveBuf_.setSize(pf.size());
122  (
123  commsType,
124  procPatch_.neighbProcNo(),
125  reinterpret_cast<char*>(receiveBuf_.begin()),
126  receiveBuf_.byteSize(),
127  procPatch_.tag(),
128  procPatch_.comm()
129  );
130  }
132  (
133  commsType,
134  procPatch_.neighbProcNo(),
135  reinterpret_cast<const char*>(pf.begin()),
136  pf.byteSize(),
137  procPatch_.tag(),
138  procPatch_.comm()
139  );
140  }
141 }
142 
143 
144 template<class Type>
146 (
147  const Pstream::commsTypes commsType,
148  Field<Type>& pField
149 ) const
150 {
151  if (Pstream::parRun())
152  {
153  // If nonblocking data has already been received into receiveBuf_
154  if (commsType != Pstream::commsTypes::nonBlocking)
155  {
156  receiveBuf_.setSize(this->size());
158  (
159  commsType,
160  procPatch_.neighbProcNo(),
161  reinterpret_cast<char*>(receiveBuf_.begin()),
162  receiveBuf_.byteSize(),
163  procPatch_.tag(),
164  procPatch_.comm()
165  );
166  }
167 
168  if (doTransform())
169  {
170  const processorCyclicPolyPatch& ppp =
171  procPatch_.procCyclicPolyPatch();
172  const tensor& forwardT = ppp.forwardT()[0];
173 
174  transform(receiveBuf_, forwardT, receiveBuf_);
175  }
176 
177  // All points are separated
178  this->addToInternalField(pField, receiveBuf_);
179  }
180 }
181 
182 
183 // ************************************************************************* //
Foam::Tensor< scalar >
Foam::processorCyclicPointPatchField::swapAddSeparated
virtual void swapAddSeparated(const Pstream::commsTypes commsType, Field< Type > &) const
Complete swap of patch point values and add to local values.
Definition: processorCyclicPointPatchField.C:146
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::coupledPointPatchField
A Coupled boundary condition for pointField.
Definition: coupledPointPatchField.H:52
Foam::processorCyclicPolyPatch
Neighbour processor patch.
Definition: processorCyclicPolyPatch.H:52
Foam::processorCyclicPointPatchField::initSwapAddSeparated
virtual void initSwapAddSeparated(const Pstream::commsTypes commsType, Field< Type > &) const
Initialise swap of non-collocated patch point values.
Definition: processorCyclicPointPatchField.C:101
Foam::processorCyclicPointPatchField
Foam::processorCyclicPointPatchField.
Definition: processorCyclicPointPatchField.H:52
Foam::pointPatch
Basic pointPatch represents a set of points from the mesh.
Definition: pointPatch.H:58
Foam::transform
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:519
transformField.H
Spatial transformation functions for primitive fields.
Foam::pointPatchField::patch
const pointPatch & patch() const
Return patch.
Definition: pointPatchField.H:268
Foam::pointPatchFieldMapper
Foam::pointPatchFieldMapper.
Definition: pointPatchFieldMapper.H:48
processorCyclicPointPatchField.H
Foam::Field
Generic templated field type.
Definition: Field.H:63
Foam::processorCyclicPointPatchField::~processorCyclicPointPatchField
virtual ~processorCyclicPointPatchField()
Destructor.
Definition: processorCyclicPointPatchField.C:93
Foam::blockMeshTools::read
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
Definition: blockMeshTools.C:33
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
processorPolyPatch.H
Foam::UPstream::commsTypes
commsTypes
Types of communications.
Definition: UPstream.H:66
Foam::processorCyclicPointPatchField::processorCyclicPointPatchField
processorCyclicPointPatchField(const pointPatch &, const DimensionedField< Type, pointMesh > &)
Construct from patch and internal field.
Definition: processorCyclicPointPatchField.C:37
Foam::vtk::write
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
Definition: foamVtkOutputTemplates.C:35
Foam::processorCyclicPolyPatch::forwardT
virtual const tensorField & forwardT() const
Return face transformation tensor.
Definition: processorCyclicPolyPatch.H:368
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54