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  Copyright (C) 2020-2023 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
30 #include "transformField.H"
32 
33 
34 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35 
36 template<class Type>
38 (
39  const pointPatch& p,
41 )
42 :
44  procPatch_(refCast<const processorCyclicPointPatch>(p))
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 {}
59 
60 
61 template<class Type>
63 (
65  const pointPatch& p,
67  const pointPatchFieldMapper& mapper
68 )
69 :
70  coupledPointPatchField<Type>(ptf, p, iF, mapper),
71  procPatch_(refCast<const processorCyclicPointPatch>(ptf.patch()))
72 {}
73 
74 
75 template<class Type>
77 (
80 )
81 :
82  coupledPointPatchField<Type>(ptf, iF),
83  procPatch_(refCast<const processorCyclicPointPatch>(ptf.patch()))
84 {}
85 
86 
87 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
88 
89 template<class Type>
91 (
92  const Pstream::commsTypes commsType,
93  Field<Type>& pField
94 ) const
95 {
96  if (UPstream::parRun())
97  {
98  // Get internal field into correct order for opposite side. Note use
99  // of member data buffer since using non-blocking. Could be optimised
100  // out if not using non-blocking...
101  this->patchInternalField
102  (
103  pField,
104  procPatch_.reverseMeshPoints(),
105  sendBuf_
106  );
107 
108  if (commsType == Pstream::commsTypes::nonBlocking)
109  {
110  recvBuf_.resize_nocopy(sendBuf_.size());
111 
113  (
114  commsType,
115  procPatch_.neighbProcNo(),
116  recvBuf_.data_bytes(),
117  recvBuf_.size_bytes(),
118  procPatch_.tag(),
119  procPatch_.comm()
120  );
121  }
123  (
124  commsType,
125  procPatch_.neighbProcNo(),
126  sendBuf_.cdata_bytes(),
127  sendBuf_.size_bytes(),
128  procPatch_.tag(),
129  procPatch_.comm()
130  );
131  }
132 }
133 
134 
135 template<class Type>
137 (
138  const Pstream::commsTypes commsType,
139  Field<Type>& pField
140 ) const
141 {
142  if (UPstream::parRun())
143  {
144  // If nonblocking, data is already in receive buffer
145 
146  if (commsType != Pstream::commsTypes::nonBlocking)
147  {
148  recvBuf_.resize_nocopy(this->size());
150  (
151  commsType,
152  procPatch_.neighbProcNo(),
153  recvBuf_.data_bytes(),
154  recvBuf_.size_bytes(),
155  procPatch_.tag(),
156  procPatch_.comm()
157  );
158  }
159 
160  if (doTransform())
161  {
162  const processorCyclicPolyPatch& ppp =
163  procPatch_.procCyclicPolyPatch();
164  const tensor& forwardT = ppp.forwardT()[0];
165 
166  transform(recvBuf_, forwardT, recvBuf_);
167  }
168 
169  // All points are separated
170  this->addToInternalField(pField, recvBuf_);
171  }
172 }
173 
174 
175 // ************************************************************************* //
Foam::processorCyclicPointPatchField.
dictionary dict
processorCyclicPointPatchField(const pointPatch &, const DimensionedField< Type, pointMesh > &)
Construct from patch and internal field.
commsTypes
Communications types.
Definition: UPstream.H:77
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
Processor patch boundary needs to be such that the ordering of points in the patch is the same on bot...
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition: typeInfo.H:172
Foam::pointPatchFieldMapper.
Tensor< scalar > tensor
Definition: symmTensor.H:57
virtual void initSwapAddSeparated(const Pstream::commsTypes commsType, Field< Type > &) const
Initialise swap of non-collocated patch point values.
Spatial transformation functions for primitive fields.
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
Generic templated field type.
Definition: Field.H:62
void read(Istream &, label &val, const dictionary &)
In-place read with dictionary lookup.
Basic pointPatch represents a set of points from the mesh.
Definition: pointPatch.H:61
A Coupled boundary condition for pointField.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const std::string patch
OpenFOAM patch number as a std::string.
volScalarField & p
virtual void swapAddSeparated(const Pstream::commsTypes commsType, Field< Type > &) const
Complete swap of patch point values and add to local values.
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:521