cyclicACMIFvPatchField.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) 2013-2016 OpenFOAM Foundation
9  Copyright (C) 2019,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 Class
28  Foam::cyclicACMIFvPatchField
29 
30 Group
31  grpCoupledBoundaryConditions
32 
33 Description
34  This boundary condition enforces a cyclic condition between a pair of
35  boundaries, whereby communication between the patches is performed using
36  an arbitrarily coupled mesh interface (ACMI) interpolation.
37 
38 Usage
39  Example of the boundary condition specification:
40  \verbatim
41  <patchName>
42  {
43  type cyclicACMI;
44  value <initial value>;
45  neighbourValue <initial value of neighbour patch cells>;
46  }
47  \endverbatim
48 
49 See also
50  Foam::AMIInterpolation
51 
52 SourceFiles
53  cyclicACMIFvPatchField.C
54 
55 \*---------------------------------------------------------------------------*/
56 
57 #ifndef Foam_cyclicACMIFvPatchField_H
58 #define Foam_cyclicACMIFvPatchField_H
59 
60 #include "coupledFvPatchField.H"
62 #include "cyclicACMIFvPatch.H"
63 
64 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
65 
66 namespace Foam
67 {
68 
69 /*---------------------------------------------------------------------------*\
70  Class cyclicACMIFvPatchField Declaration
71 \*---------------------------------------------------------------------------*/
72 
73 template<class Type>
75 :
76  virtual public cyclicACMILduInterfaceField,
77  public coupledFvPatchField<Type>
78 {
79  // Private data
80 
81  //- Local reference cast into the cyclic patch
82  const cyclicACMIFvPatch& cyclicACMIPatch_;
83 
84 
85  // Sending and receiving (distributed AMI)
86 
87  //- Current range of send requests (non-blocking)
88  mutable labelRange sendRequests_;
89 
90  //- Current range of recv requests (non-blocking)
91  mutable labelRange recvRequests_;
92 
93  //- Send buffers
94  mutable PtrList<List<Type>> sendBufs_;
95 
96  //- Receive buffers_
97  mutable PtrList<List<Type>> recvBufs_;
98 
99  //- Scalar send buffers
100  mutable PtrList<List<solveScalar>> scalarSendBufs_;
101 
102  //- Scalar receive buffers
103  mutable PtrList<List<solveScalar>> scalarRecvBufs_;
104 
105  //- Neighbour coupled internal cell data
106  mutable autoPtr<Field<Type>> patchNeighbourFieldPtr_;
107 
108 
109  // Private Member Functions
110 
111  //- Return the AMI corresponding to the owner side
112  const AMIPatchToPatchInterpolation& ownerAMI() const
113  {
114  return
115  (
116  cyclicACMIPatch_.owner()
117  ? cyclicACMIPatch_.AMI()
118  : cyclicACMIPatch_.neighbPatch().AMI()
119  );
120  }
121 
122  //- All receive/send requests have completed
123  virtual bool all_ready() const;
124 
125  //- Use neighbour field caching
126  static bool cacheNeighbourField();
127 
128  //- Return neighbour coupled internal cell data
130 
131  //- Return neighbour side field given internal fields
132  template<class Type2>
133  tmp<Field<Type2>> neighbourSideField
134  (
135  const Field<Type2>&
136  ) const;
137 
138  //- Return new matrix coeffs
139  tmp<Field<scalar>> coeffs
140  (
141  fvMatrix<Type>& matrix,
142  const Field<scalar>&,
143  const label
144  ) const;
145 
146 
147 
148 public:
149 
150  //- Runtime type information
151  TypeName(cyclicACMIFvPatch::typeName_());
152 
153 
154  // Constructors
155 
156  //- Construct from patch and internal field
158  (
159  const fvPatch&,
161  );
162 
163  //- Construct from patch, internal field and dictionary
165  (
166  const fvPatch&,
168  const dictionary&
169  );
170 
171  //- Construct by mapping given cyclicACMIFvPatchField onto a new patch
173  (
175  const fvPatch&,
177  const fvPatchFieldMapper&
178  );
179 
180  //- Construct as copy
182 
183  //- Construct as copy setting internal field reference
185  (
188  );
189 
190  //- Return a clone
191  virtual tmp<fvPatchField<Type>> clone() const
192  {
193  return fvPatchField<Type>::Clone(*this);
194  }
195 
196  //- Clone with an internal field reference
198  (
200  ) const
201  {
202  return fvPatchField<Type>::Clone(*this, iF);
203  }
204 
205 
206  // Member functions
207 
208  // Access
209 
210  //- Return local reference cast into the cyclic AMI patch
211  const cyclicACMIFvPatch& cyclicACMIPatch() const
212  {
213  return cyclicACMIPatch_;
214  }
215 
216 
217  // Coupling
218 
219  //- Return true if coupled. Note that the underlying patch
220  // is not coupled() - the points don't align
221  virtual bool coupled() const;
222 
223  //- Are all (receive) data available?
224  virtual bool ready() const;
225 
226  //- Return true if this patch field fixes a value
227  // Needed to check if a level has to be specified while solving
228  // Poisson equations
229  virtual bool fixesValue() const
230  {
231  const scalarField& mask =
232  cyclicACMIPatch_.cyclicACMIPatch().mask();
233 
234  if (gMax(mask) > 1e-5)
235  {
236  // regions connected
237  return false;
238  }
239  else
240  {
241  // fully separated
242  return nonOverlapPatchField().fixesValue();
243  }
244  }
245 
246  //- Return neighbour coupled internal cell data
247  virtual tmp<Field<Type>> patchNeighbourField() const;
248 
249  //- Return reference to neighbour patchField
250  const cyclicACMIFvPatchField<Type>& neighbourPatchField() const;
251 
252  //- Return reference to non-overlapping patchField
254 
255 
256  // Evaluation
257 
258  //- Initialise the evaluation of the patch field
259  virtual void initEvaluate(const Pstream::commsTypes commsType);
260 
261  //- Evaluate the patch field
262  virtual void evaluate(const Pstream::commsTypes commsType);
263 
264 
265  // Coupled interface functionality
266 
267  //- Initialise neighbour matrix update
268  virtual void initInterfaceMatrixUpdate
269  (
270  solveScalarField& result,
271  const bool add,
272  const lduAddressing& lduAddr,
273  const label patchId,
274  const solveScalarField& psiInternal,
275  const scalarField& coeffs,
276  const direction cmpt,
277  const Pstream::commsTypes commsType
278  ) const;
279 
280  //- Update result field based on interface functionality
281  virtual void updateInterfaceMatrix
282  (
283  solveScalarField& result,
284  const bool add,
285  const lduAddressing& lduAddr,
286  const label patchId,
287  const solveScalarField& psiInternal,
288  const scalarField& coeffs,
289  const direction cmpt,
290  const Pstream::commsTypes commsType
291  ) const;
292 
293  //- Initialise neighbour matrix update
294  virtual void initInterfaceMatrixUpdate
295  (
296  Field<Type>& result,
297  const bool add,
298  const lduAddressing& lduAddr,
299  const label patchId,
300  const Field<Type>& psiInternal,
301  const scalarField& coeffs,
302  const Pstream::commsTypes commsType
303  ) const;
304 
305  //- Update result field based on interface functionality
306  virtual void updateInterfaceMatrix
307  (
308  Field<Type>&,
309  const bool add,
310  const lduAddressing& lduAddr,
311  const label patchId,
312  const Field<Type>&,
313  const scalarField&,
314  const Pstream::commsTypes commsType
315  ) const;
316 
317  //- Manipulate matrix
318  virtual void manipulateMatrix(fvMatrix<Type>& matrix);
319 
320  //- Manipulate matrix
321  virtual void manipulateMatrix
322  (
323  fvMatrix<Type>& m,
324  const label iMatrix,
325  const direction cmpt
326  );
327 
328  //- Update the coefficients associated with the patch field
329  virtual void updateCoeffs();
330 
331 
332  // Cyclic AMI coupled interface functions
333 
334  //- Does the patch field perform the transformation
335  virtual bool doTransform() const
336  {
337  return (pTraits<Type>::rank && !cyclicACMIPatch_.parallel());
338  }
339 
340  //- Return face transformation tensor
341  virtual const tensorField& forwardT() const
342  {
343  return cyclicACMIPatch_.forwardT();
344  }
345 
346  //- Return neighbour-cell transformation tensor
347  virtual const tensorField& reverseT() const
348  {
349  return cyclicACMIPatch_.reverseT();
350  }
351 
352  //- Return rank of component for transform
353  virtual int rank() const
354  {
355  return pTraits<Type>::rank;
356  }
357 
358 
359  // I-O
360 
361  //- Write
362  virtual void write(Ostream& os) const;
363 };
364 
365 
366 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
367 
368 } // End namespace Foam
369 
370 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
371 
372 #ifdef NoRepository
373  #include "cyclicACMIFvPatchField.C"
374 #endif
375 
376 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
377 
378 #endif
379 
380 // ************************************************************************* //
label patchId(-1)
virtual void evaluate(const Pstream::commsTypes commsType)
Evaluate the patch field.
virtual int rank() const
Return rank of component for transform.
uint8_t direction
Definition: direction.H:46
TypeName(cyclicACMIFvPatch::typeName_())
Runtime type information.
const fvPatchField< Type > & nonOverlapPatchField() const
Return reference to non-overlapping patchField.
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
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 bool owner() const
virtual tmp< Field< Type > > patchNeighbourField() const
Return neighbour coupled internal cell data.
A range or interval of labels defined by a start and a size.
Definition: labelRange.H:52
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
A traits class, which is primarily used for primitives and vector-space.
Definition: pTraits.H:75
const cyclicACMIFvPatch & cyclicACMIPatch() const
Return local reference cast into the cyclic AMI patch.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
virtual tmp< fvPatchField< Type > > clone() const
Return a clone.
This boundary condition enforces a cyclic condition between a pair of boundaries, whereby communicati...
virtual bool coupled() const
Return true if coupled. Note that the underlying patch.
virtual bool fixesValue() const
Return true if this patch field fixes a value.
virtual const tensorField & reverseT() const
Return neighbour-cell transformation tensor.
static tmp< fvPatchField< Type > > Clone(const DerivedPatchField &pf, Args &&... args)
Clone a patch field, optionally with internal field reference etc.
Definition: fvPatchField.H:597
virtual bool doTransform() const
Does the patch field perform the transformation.
Cyclic patch for Arbitrarily Coupled Mesh Interface (ACMI)
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
virtual bool ready() const
Are all (receive) data available?
virtual void write(Ostream &os) const
Write.
Generic templated field type.
Definition: Field.H:62
virtual const cyclicACMIFvPatch & neighbPatch() const
Return neighbour fvPatch.
A FieldMapper for finite-volume patch fields.
virtual const tensorField & forwardT() const
Return face transformation tensor.
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:64
virtual bool parallel() const
Are the cyclic planes parallel.
const cyclicACMIFvPatchField< Type > & neighbourPatchField() const
Return reference to neighbour patchField.
cyclicACMIFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
Abstract base class for coupled patches.
Abstract base class for cyclic ACMI coupled interfaces.
void add(FieldField< Field1, typename typeOfSum< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Type gMax(const FieldField< Field, Type > &f)
OBJstream os(runTime.globalPath()/outputName)
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.
Field< tensor > tensorField
Specialisation of Field<T> for tensor.
const scalarField & mask() const
Mask field where 1 = overlap(coupled), 0 = no-overlap.
virtual void initEvaluate(const Pstream::commsTypes commsType)
Initialise the evaluation of the patch field.
const cyclicACMIPolyPatch & cyclicACMIPatch() const
Return local reference cast into the cyclic patch.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
virtual const tensorField & forwardT() const
Return face transformation tensor.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Interpolation class dealing with transfer of data between two primitive patches with an arbitrary mes...
virtual void manipulateMatrix(fvMatrix< Type > &matrix)
Manipulate matrix.
::Foam::direction rank(const expressions::valueTypeCode) noexcept
The vector-space rank associated with given valueTypeCode.
Definition: exprTraits.C:70
virtual const AMIPatchToPatchInterpolation & AMI() const
Return a reference to the AMI interpolator.
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
The class contains the addressing required by the lduMatrix: upper, lower and losort.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
virtual const tensorField & reverseT() const
Return neighbour-cell transformation tensor.
Namespace for OpenFOAM.