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 and return a clone
184  virtual tmp<fvPatchField<Type>> clone() const
185  {
186  return tmp<fvPatchField<Type>>
187  (
189  );
190  }
191 
192  //- Construct as copy setting internal field reference
194  (
197  );
198 
199  //- Construct and return a clone setting internal field reference
201  (
203  ) const
204  {
205  return tmp<fvPatchField<Type>>
206  (
207  new cyclicACMIFvPatchField<Type>(*this, iF)
208  );
209  }
210 
211 
212  // Member functions
213 
214  // Access
215 
216  //- Return local reference cast into the cyclic AMI patch
217  const cyclicACMIFvPatch& cyclicACMIPatch() const
218  {
219  return cyclicACMIPatch_;
220  }
221 
222 
223  // Coupling
224 
225  //- Return true if coupled. Note that the underlying patch
226  // is not coupled() - the points don't align
227  virtual bool coupled() const;
228 
229  //- Are all (receive) data available?
230  virtual bool ready() const;
231 
232  //- Return true if this patch field fixes a value
233  // Needed to check if a level has to be specified while solving
234  // Poisson equations
235  virtual bool fixesValue() const
236  {
237  const scalarField& mask =
238  cyclicACMIPatch_.cyclicACMIPatch().mask();
239 
240  if (gMax(mask) > 1e-5)
241  {
242  // regions connected
243  return false;
244  }
245  else
246  {
247  // fully separated
248  return nonOverlapPatchField().fixesValue();
249  }
250  }
251 
252  //- Return neighbour coupled internal cell data
253  virtual tmp<Field<Type>> patchNeighbourField() const;
254 
255  //- Return reference to neighbour patchField
256  const cyclicACMIFvPatchField<Type>& neighbourPatchField() const;
257 
258  //- Return reference to non-overlapping patchField
260 
261 
262  // Evaluation
263 
264  //- Initialise the evaluation of the patch field
265  virtual void initEvaluate(const Pstream::commsTypes commsType);
266 
267  //- Evaluate the patch field
268  virtual void evaluate(const Pstream::commsTypes commsType);
269 
270 
271  // Coupled interface functionality
272 
273  //- Initialise neighbour matrix update
274  virtual void initInterfaceMatrixUpdate
275  (
276  solveScalarField& result,
277  const bool add,
278  const lduAddressing& lduAddr,
279  const label patchId,
280  const solveScalarField& psiInternal,
281  const scalarField& coeffs,
282  const direction cmpt,
283  const Pstream::commsTypes commsType
284  ) const;
285 
286  //- Update result field based on interface functionality
287  virtual void updateInterfaceMatrix
288  (
289  solveScalarField& result,
290  const bool add,
291  const lduAddressing& lduAddr,
292  const label patchId,
293  const solveScalarField& psiInternal,
294  const scalarField& coeffs,
295  const direction cmpt,
296  const Pstream::commsTypes commsType
297  ) const;
298 
299  //- Initialise neighbour matrix update
300  virtual void initInterfaceMatrixUpdate
301  (
302  Field<Type>& result,
303  const bool add,
304  const lduAddressing& lduAddr,
305  const label patchId,
306  const Field<Type>& psiInternal,
307  const scalarField& coeffs,
308  const Pstream::commsTypes commsType
309  ) const;
310 
311  //- Update result field based on interface functionality
312  virtual void updateInterfaceMatrix
313  (
314  Field<Type>&,
315  const bool add,
316  const lduAddressing& lduAddr,
317  const label patchId,
318  const Field<Type>&,
319  const scalarField&,
320  const Pstream::commsTypes commsType
321  ) const;
322 
323  //- Manipulate matrix
324  virtual void manipulateMatrix(fvMatrix<Type>& matrix);
325 
326  //- Manipulate matrix
327  virtual void manipulateMatrix
328  (
329  fvMatrix<Type>& m,
330  const label iMatrix,
331  const direction cmpt
332  );
333 
334  //- Update the coefficients associated with the patch field
335  virtual void updateCoeffs();
336 
337 
338  // Cyclic AMI coupled interface functions
339 
340  //- Does the patch field perform the transformation
341  virtual bool doTransform() const
342  {
343  return (pTraits<Type>::rank && !cyclicACMIPatch_.parallel());
344  }
345 
346  //- Return face transformation tensor
347  virtual const tensorField& forwardT() const
348  {
349  return cyclicACMIPatch_.forwardT();
350  }
351 
352  //- Return neighbour-cell transformation tensor
353  virtual const tensorField& reverseT() const
354  {
355  return cyclicACMIPatch_.reverseT();
356  }
357 
358  //- Return rank of component for transform
359  virtual int rank() const
360  {
361  return pTraits<Type>::rank;
362  }
363 
364 
365  // I-O
366 
367  //- Write
368  virtual void write(Ostream& os) const;
369 };
370 
371 
372 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
373 
374 } // End namespace Foam
375 
376 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
377 
378 #ifdef NoRepository
379  #include "cyclicACMIFvPatchField.C"
380 #endif
381 
382 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
383 
384 #endif
385 
386 // ************************************************************************* //
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:72
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
Construct and 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.
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...
Definition: areaFieldsFwd.H:42
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.