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-2021 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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"
31#include "processorPolyPatch.H"
32
33
34// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35
36template<class Type>
38(
39 const pointPatch& p,
41)
42:
43 coupledPointPatchField<Type>(p, iF),
44 procPatch_(refCast<const processorCyclicPointPatch>(p)),
45 receiveBuf_(0)
46{}
47
48
49template<class Type>
51(
52 const pointPatch& p,
54 const dictionary& dict
55)
56:
57 coupledPointPatchField<Type>(p, iF, dict),
58 procPatch_(refCast<const processorCyclicPointPatch>(p, dict)),
59 receiveBuf_(0)
60{}
61
62
63template<class Type>
65(
67 const pointPatch& p,
69 const pointPatchFieldMapper& mapper
70)
71:
72 coupledPointPatchField<Type>(ptf, p, iF, mapper),
73 procPatch_(refCast<const processorCyclicPointPatch>(ptf.patch())),
74 receiveBuf_(0)
75{}
76
77
78template<class Type>
80(
83)
84:
85 coupledPointPatchField<Type>(ptf, iF),
86 procPatch_(refCast<const processorCyclicPointPatch>(ptf.patch())),
87 receiveBuf_(0)
88{}
89
90
91// * * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * //
92
93template<class Type>
95{}
96
97
98// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
99
100template<class Type>
102(
103 const Pstream::commsTypes commsType,
104 Field<Type>& pField
105) const
106{
107 if (Pstream::parRun())
108 {
109 // Get internal field into correct order for opposite side
110 Field<Type> pf
111 (
112 this->patchInternalField
113 (
114 pField,
115 procPatch_.reverseMeshPoints()
116 )
117 );
118
119 if (commsType == Pstream::commsTypes::nonBlocking)
120 {
121 receiveBuf_.setSize(pf.size());
123 (
124 commsType,
125 procPatch_.neighbProcNo(),
126 receiveBuf_.data_bytes(),
127 receiveBuf_.size_bytes(),
128 procPatch_.tag(),
129 procPatch_.comm()
130 );
131 }
133 (
134 commsType,
135 procPatch_.neighbProcNo(),
136 pf.cdata_bytes(),
137 pf.size_bytes(),
138 procPatch_.tag(),
139 procPatch_.comm()
140 );
141 }
142}
143
144
145template<class Type>
147(
148 const Pstream::commsTypes commsType,
149 Field<Type>& pField
150) const
151{
152 if (Pstream::parRun())
153 {
154 // If nonblocking data has already been received into receiveBuf_
155 if (commsType != Pstream::commsTypes::nonBlocking)
156 {
157 receiveBuf_.setSize(this->size());
159 (
160 commsType,
161 procPatch_.neighbProcNo(),
162 receiveBuf_.data_bytes(),
163 receiveBuf_.size_bytes(),
164 procPatch_.tag(),
165 procPatch_.comm()
166 );
167 }
168
169 if (doTransform())
170 {
171 const processorCyclicPolyPatch& ppp =
172 procPatch_.procCyclicPolyPatch();
173 const tensor& forwardT = ppp.forwardT()[0];
174
175 transform(receiveBuf_, forwardT, receiveBuf_);
176 }
177
178 // All points are separated
179 this->addToInternalField(pField, receiveBuf_);
180 }
181}
182
183
184// ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic templated field type.
Definition: Field.H:82
virtual bool read()
Re-read model coefficients if they have changed.
const char * cdata_bytes() const noexcept
Return pointer to the underlying array serving as data storage,.
Definition: UListI.H:244
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
std::streamsize size_bytes() const noexcept
Number of contiguous bytes for the List data.
Definition: UListI.H:258
commsTypes
Types of communications.
Definition: UPstream.H:67
@ nonBlocking
"nonBlocking"
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:433
A Coupled boundary condition for pointField.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
virtual bool write()
Write the output fields.
Foam::pointPatchFieldMapper.
Basic pointPatch represents a set of points from the mesh.
Definition: pointPatch.H:64
Foam::processorCyclicPointPatchField.
virtual void swapAddSeparated(const Pstream::commsTypes commsType, Field< Type > &) const
Complete swap of patch point values and add to local values.
virtual void initSwapAddSeparated(const Pstream::commsTypes commsType, Field< Type > &) const
Initialise swap of non-collocated patch point values.
Processor patch boundary needs to be such that the ordering of points in the patch is the same on bot...
virtual const tensorField & forwardT() const
Return face transformation tensor.
volScalarField & p
To & refCast(From &r)
Reference type cast template function.
Definition: typeInfo.H:131
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:536
dictionary dict
Spatial transformation functions for primitive fields.