Loading...
Searching...
No Matches
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-2025 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
31#include "processorPolyPatch.H"
32
33
34// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35
36template<class Type>
38(
39 const pointPatch& p,
41)
43 coupledPointPatchField<Type>(p, iF),
44 procPatch_(refCast<const processorCyclicPointPatch>(p))
45{}
46
47
48template<class Type>
50(
51 const pointPatch& p,
53 const dictionary& dict
54)
67 const pointPatchFieldMapper& mapper
68)
70 coupledPointPatchField<Type>(ptf, p, iF, mapper),
71 procPatch_(refCast<const processorCyclicPointPatch>(ptf.patch()))
72{}
73
74
75template<class Type>
77(
80)
81:
82 coupledPointPatchField<Type>(ptf, iF),
83 procPatch_(refCast<const processorCyclicPointPatch>(ptf.patch()))
84{}
85
86
87// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
88
89template<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
102 sendBuf_.resize_nocopy(this->size()); // <- presize buffer
103 recvBuf_.resize_nocopy(this->size());
104
105 this->patchInternalField
106 (
107 pField,
108 procPatch_.reverseMeshPoints(),
109 sendBuf_
110 );
111
112 if (commsType == Pstream::commsTypes::nonBlocking)
113 {
115 (
116 commsType,
117 procPatch_.neighbProcNo(),
118 recvBuf_,
119 procPatch_.tag(),
120 procPatch_.comm()
121 );
122 }
123
124 UOPstream::write
125 (
126 commsType,
127 procPatch_.neighbProcNo(),
128 sendBuf_,
129 procPatch_.tag(),
130 procPatch_.comm()
131 );
132 }
133}
134
135
136template<class Type>
138(
139 const Pstream::commsTypes commsType,
140 Field<Type>& pField
141) const
142{
143 if (UPstream::parRun())
144 {
145 // If nonblocking, data is already in receive buffer
146
147 if (commsType != Pstream::commsTypes::nonBlocking)
148 {
149 recvBuf_.resize_nocopy(this->size()); // In general a no-op
150
151 UIPstream::read
152 (
153 commsType,
154 procPatch_.neighbProcNo(),
155 recvBuf_,
156 procPatch_.tag(),
157 procPatch_.comm()
158 );
159 }
160
161 if (doTransform())
162 {
163 const auto& ppp = 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// ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
static std::streamsize read(const UPstream::commsTypes commsType, const int fromProcNo, Type *buffer, std::streamsize count, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm, UPstream::Request *req=nullptr)
Receive buffer contents (contiguous types) from given processor.
static bool write(const UPstream::commsTypes commsType, const int toProcNo, const Type *buffer, std::streamsize count, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm, UPstream::Request *req=nullptr, const UPstream::sendModes sendMode=UPstream::sendModes::normal)
Write buffer contents (contiguous types only) to given processor.
commsTypes
Communications types.
Definition UPstream.H:81
@ nonBlocking
"nonBlocking" (immediate) : (MPI_Isend, MPI_Irecv)
Definition UPstream.H:84
static bool parRun(const bool on) noexcept
Set as parallel run on/off.
Definition UPstream.H:1669
A Coupled boundary condition for pointField.
coupledPointPatchField(const pointPatch &, const DimensionedField< Type, pointMesh > &)
Construct from patch and internal field.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
const pointPatch & patch() const noexcept
Return the patch.
Foam::pointPatchFieldMapper.
void patchInternalField(const UList< Type1 > &internalData, const labelUList &addressing, UList< Type1 > &pfld) const
Extract field using specified addressing.
void addToInternalField(Field< Type1 > &iF, const Field< Type1 > &pF) const
Given the internal field and a patch field, add the patch field to the internal field.
label size() const noexcept
Return the patch size.
Basic pointPatch represents a set of points from the mesh.
Definition pointPatch.H:67
Foam::processorCyclicPointPatchField.
processorCyclicPointPatchField(const pointPatch &, const DimensionedField< Type, pointMesh > &)
Construct from patch and internal field.
virtual bool doTransform() const
Does the patch field perform the transformation.
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...
Tensor of scalars, i.e. Tensor<scalar>.
volScalarField & p
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
refinementData transform(const tensor &, const refinementData val)
No-op rotational transform for base types.
Tensor< scalar > tensor
Definition symmTensor.H:57
dictionary dict
Spatial transformation functions for primitive fields.