Loading...
Searching...
No Matches
cyclicAMIFvPatchField.H
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-2016 OpenFOAM Foundation
9 Copyright (C) 2019-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
27Class
28 Foam::cyclicAMIFvPatchField
29
30Group
31 grpCoupledBoundaryConditions
32
33Description
34 This boundary condition enforces a cyclic condition between a pair of
35 boundaries, whereby communication between the patches is performed using
36 an arbitrary mesh interface (AMI) interpolation.
37
38Usage
39 Example of the boundary condition specification:
40 \verbatim
41 <patchName>
42 {
43 type cyclicAMI;
44 value <initial value>;
45 neighbourValue <initial value of neighbour patch cells>;
46 }
47 \endverbatim
48
49Note
50 The outer boundary of the patch pairs must be similar, i.e. if the owner
51 patch is transformed to the neighbour patch, the outer perimeter of each
52 patch should be identical (or very similar).
53
54 The \c neighbourValue is only needed when running distributed,
55 i.e. the neighbour cells are on a different processor from the owner cells.
56 It guarantees consistent restart e.g. when doing a snGrad and avoids
57 additional communication.
58
59See also
60 Foam::AMIInterpolation
61 Foam::processorFvPatchField
62
63SourceFiles
64 cyclicAMIFvPatchField.C
65
66\*---------------------------------------------------------------------------*/
67
68#ifndef Foam_cyclicAMIFvPatchField_H
69#define Foam_cyclicAMIFvPatchField_H
70
71#include "cyclicAMIFvPatch.H"
72#include "coupledFvPatchField.H"
74
75// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
76
77namespace Foam
78{
79
80/*---------------------------------------------------------------------------*\
81 Class cyclicAMIFvPatchField Declaration
82\*---------------------------------------------------------------------------*/
83
84template<class Type>
86:
87 virtual public cyclicAMILduInterfaceField,
88 public coupledFvPatchField<Type>
89{
90 // Private Data
91
92 //- Local reference cast into the cyclic patch
93 const cyclicAMIFvPatch& cyclicAMIPatch_;
94
95
96 // Sending and receiving (distributed AMI)
97
98 //- Current range of send requests (non-blocking)
99 mutable labelRange sendRequests_;
100
101 //- Current range of recv requests (non-blocking)
102 mutable labelRange recvRequests_;
103
104 //- Send buffers
105 mutable PtrList<List<Type>> sendBufs_;
106
107 //- Receive buffers
108 mutable PtrList<List<Type>> recvBufs_;
109
110 //- Scalar send buffers
111 mutable PtrList<List<solveScalar>> scalarSendBufs_;
112
113 //- Scalar receive buffers
114 mutable PtrList<List<solveScalar>> scalarRecvBufs_;
115
116
117 // Only used for AMI caching
118
119 //- Current range of send requests (non-blocking)
120 mutable labelRange sendRequests1_;
121
122 //- Current range of recv requests (non-blocking)
123 mutable labelRange recvRequests1_;
124
125 //- Send buffers
126 mutable PtrList<List<Type>> sendBufs1_;
127
128 //- Receive buffers_
129 mutable PtrList<List<Type>> recvBufs1_;
130
131 //- Scalar send buffers
132 mutable PtrList<List<solveScalar>> scalarSendBufs1_;
133
134 //- Scalar receive buffers
135 mutable PtrList<List<solveScalar>> scalarRecvBufs1_;
136
137
138 //- Neighbour coupled internal cell data
139 mutable autoPtr<Field<Type>> patchNeighbourFieldPtr_;
140
141
142 // Private Member Functions
143
144 //- Return the AMI corresponding to the owner side
145 const AMIPatchToPatchInterpolation& ownerAMI() const
146 {
147 return
148 (
149 cyclicAMIPatch_.owner()
150 ? cyclicAMIPatch_.AMI()
151 : cyclicAMIPatch_.neighbPatch().AMI()
152 );
153 }
154
155 //- All receive/send requests have completed
156 virtual bool all_ready() const;
157
158 //- Use neighbour field caching
159 bool cacheNeighbourField() const;
160
161 //- Return neighbour coupled internal cell data
162 tmp<Field<Type>> getNeighbourField(const UList<Type>&) const;
163
164 //- Return neighbour coupled internal cell data (cached or extracted),
165 //- with optional check that AMI.comm() is valid.
166 tmp<Field<Type>> getPatchNeighbourField(bool checkCommunicator) const;
167
168 //- Return new matrix coeffs
169 tmp<Field<scalar>> coeffs
170 (
171 fvMatrix<Type>& matrix,
172 const Field<scalar>&,
173 const label
174 ) const;
175
176 template<class Type2>
177 void collectStencilData
178 (
179 const refPtr<mapDistribute>& mapPtr,
180 const labelListList& stencil,
181 const Type2& data,
182 List<Type2>& expandedData
183 );
184
185
186public:
187
188 //- Runtime type information
189 TypeName(cyclicAMIFvPatch::typeName_());
190
191
192 // Constructors
193
194 //- Construct from patch and internal field
196 (
197 const fvPatch&,
199 );
200
201 //- Construct from patch, internal field and dictionary
203 (
204 const fvPatch&,
206 const dictionary&
207 );
208
209 //- Construct by mapping given cyclicAMIFvPatchField onto a new patch
211 (
213 const fvPatch&,
215 const fvPatchFieldMapper&
216 );
217
218 //- Construct as copy
220
221 //- Construct as copy setting internal field reference
223 (
226 );
227
228 //- Return a clone
229 virtual tmp<fvPatchField<Type>> clone() const
230 {
231 return fvPatchField<Type>::Clone(*this);
232 }
233
234 //- Clone with an internal field reference
236 (
238 ) const
239 {
240 return fvPatchField<Type>::Clone(*this, iF);
241 }
242
243
244 // Member Functions
245
246 //- Return local reference cast into the cyclic AMI patch
247 const cyclicAMIFvPatch& cyclicAMIPatch() const noexcept
248 {
249 return cyclicAMIPatch_;
250 }
251
252
253 // Coupling
254
255 //- Return true if coupled. Note that the underlying patch
256 // is not coupled() - the points don't align.
257 virtual bool coupled() const { return cyclicAMIPatch_.coupled(); }
258
259 //- Are all (receive) data available?
260 virtual bool ready() const;
261
262 //- Return neighbour coupled internal cell data
263 virtual tmp<Field<Type>> patchNeighbourField() const;
264
265 //- Retrieve neighbour coupled internal cell data
266 virtual void patchNeighbourField(UList<Type>& pnf) const;
267
268 //- Return reference to neighbour patchField
270
271
272 // Mapping Functions
273
274 //- Map (and resize as needed) from self given a mapping object
275 virtual void autoMap
276 (
277 const fvPatchFieldMapper&
278 );
279
280 //- Reverse map the given fvPatchField onto this fvPatchField
281 virtual void rmap
282 (
283 const fvPatchField<Type>&,
284 const labelList&
285 );
287
288 // Evaluation
289
290 //- Initialise the evaluation of the patch field
291 virtual void initEvaluate(const Pstream::commsTypes commsType);
292
293 //- Evaluate the patch field
294 virtual void evaluate(const Pstream::commsTypes commsType);
295
296
297 // Coupled interface functionality
298
299 //- Initialise neighbour matrix update
301 (
302 solveScalarField& result,
303 const bool add,
304 const lduAddressing& lduAddr,
305 const label patchId,
306 const solveScalarField& psiInternal,
307 const scalarField& coeffs,
308 const direction cmpt,
309 const Pstream::commsTypes commsType
310 ) const;
311
312 //- Update result field based on interface functionality
314 (
315 solveScalarField& result,
316 const bool add,
317 const lduAddressing& lduAddr,
318 const label patchId,
319 const solveScalarField& psiInternal,
320 const scalarField& coeffs,
321 const direction cmpt,
322 const Pstream::commsTypes commsType
323 ) const;
324
325 //- Initialise neighbour matrix update
326 virtual void initInterfaceMatrixUpdate
327 (
328 Field<Type>& result,
329 const bool add,
330 const lduAddressing& lduAddr,
331 const label patchId,
332 const Field<Type>& psiInternal,
333 const scalarField& coeffs,
334 const Pstream::commsTypes commsType
335 ) const;
336
337 //- Update result field based on interface functionality
338 virtual void updateInterfaceMatrix
339 (
341 const bool add,
342 const lduAddressing& lduAddr,
343 const label patchId,
344 const Field<Type>&,
345 const scalarField&,
346 const Pstream::commsTypes commsType
347 ) const;
348
349
350 //- Manipulate matrix
351 virtual void manipulateMatrix
352 (
354 const label iMatrix,
355 const direction cmpt
356 );
357
358
359 // Coupled interface functions
360
361 //- Does the patch field perform the transformation
362 virtual bool doTransform() const
363 {
364 return
365 (
367 && !cyclicAMIPatch_.parallel()
368 );
369 }
370
371 //- Return face transformation tensor
372 virtual const tensorField& forwardT() const
373 {
374 return cyclicAMIPatch_.forwardT();
375 }
376
377 //- Return neighbour-cell transformation tensor
378 virtual const tensorField& reverseT() const
379 {
380 return cyclicAMIPatch_.reverseT();
381 }
382
383 //- Return rank of component for transform
384 virtual int rank() const
385 {
386 return pTraits<Type>::rank;
387 }
388
389
390 // I-O
391
392 //- Write
393 virtual void write(Ostream& os) const;
394
395
396 // Member Operators
397
398 virtual void operator=(const fvPatchField<Type>&);
399 virtual void operator==(const fvPatchField<Type>&);
400};
401
402
403// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
404
405} // End namespace Foam
406
407// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
408
409#ifdef NoRepository
410 #include "cyclicAMIFvPatchField.C"
411#endif
412
413// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
414
415#endif
416
417// ************************************************************************* //
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
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition UList.H:89
UList(const UList< Type > &) noexcept=default
Type * data() noexcept
Definition UListI.H:274
commsTypes
Communications types.
Definition UPstream.H:81
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
coupledFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
This boundary condition enforces a cyclic condition between a pair of boundaries, whereby communicati...
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
virtual bool doTransform() const
Does the patch field perform the transformation.
virtual void initInterfaceMatrixUpdate(solveScalarField &result, const bool add, const lduAddressing &lduAddr, const label patchId, const solveScalarField &psiInternal, const scalarField &coeffs, const direction cmpt, const Pstream::commsTypes commsType) const
Initialise neighbour matrix update.
cyclicAMIFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
const cyclicAMIFvPatchField< Type > & neighbourPatchField() const
Return reference to neighbour patchField.
virtual int rank() const
Return rank of component for transform.
virtual void manipulateMatrix(fvMatrix< Type > &m, const label iMatrix, const direction cmpt)
Manipulate matrix.
virtual void operator==(const fvPatchField< Type > &)
virtual bool coupled() const
Return true if coupled. Note that the underlying patch.
virtual void updateInterfaceMatrix(solveScalarField &result, const bool add, const lduAddressing &lduAddr, const label patchId, const solveScalarField &psiInternal, const scalarField &coeffs, const direction cmpt, const Pstream::commsTypes commsType) const
Update result field based on interface functionality.
virtual const tensorField & reverseT() const
Return neighbour-cell transformation tensor.
virtual void initEvaluate(const Pstream::commsTypes commsType)
Initialise the evaluation of the patch field.
virtual void operator=(const fvPatchField< Type > &)
virtual void evaluate(const Pstream::commsTypes commsType)
Evaluate the patch field.
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
virtual tmp< fvPatchField< Type > > clone() const
Return a clone.
const cyclicAMIFvPatch & cyclicAMIPatch() const noexcept
Return local reference cast into the cyclic AMI patch.
virtual bool ready() const
Are all (receive) data available?
virtual tmp< Field< Type > > patchNeighbourField() const
Return neighbour coupled internal cell data.
TypeName(cyclicAMIFvPatch::typeName_())
Runtime type information.
virtual tmp< fvPatchField< Type > > clone(const DimensionedField< Type, volMesh > &iF) const
Clone with an internal field reference.
virtual const tensorField & forwardT() const
Return face transformation tensor.
Cyclic patch for Arbitrary Mesh Interface (AMI).
virtual const tensorField & reverseT() const
Return neighbour-cell transformation tensor.
virtual const tensorField & forwardT() const
Return face transformation tensor.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition fvMatrix.H:118
A FieldMapper for finite-volume patch fields.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
static tmp< fvPatchField< Type > > Clone(const DerivedPatchField &pf, Args &&... args)
Clone a patch field, optionally with internal field reference etc.
fvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
A range or interval of labels defined by a start and a size.
Definition labelRange.H:66
The class contains the addressing required by the lduMatrix: upper, lower and losort.
A traits class, which is primarily used for primitives and vector-space.
Definition pTraits.H:64
A class for managing references or pointers (no reference counting).
Definition refPtr.H:54
A class for managing temporary objects.
Definition tmp.H:75
OBJstream os(runTime.globalPath()/outputName)
label patchId(-1)
Namespace for OpenFOAM.
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
void add(DimensionedField< scalar, GeoMesh > &result, const dimensioned< scalar > &dt1, const DimensionedField< scalar, GeoMesh > &f2)
AMIInterpolation AMIPatchToPatchInterpolation
Patch-to-patch interpolation == Foam::AMIInterpolation.
Field< solveScalar > solveScalarField
uint8_t direction
Definition direction.H:49
const direction noexcept
Definition scalarImpl.H:265
Field< tensor > tensorField
Specialisation of Field<T> for tensor.
constexpr bool is_rotational_vectorspace_v
The is_rotational_vectorspace value of Type.
Definition pTraits.H:251
runTime write()
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition typeInfo.H:68