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 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  }
45  \endverbatim
46 
47 See also
48  Foam::AMIInterpolation
49 
50 SourceFiles
51  cyclicACMIFvPatchField.C
52 
53 \*---------------------------------------------------------------------------*/
54 
55 #ifndef Foam_cyclicACMIFvPatchField_H
56 #define Foam_cyclicACMIFvPatchField_H
57 
58 #include "coupledFvPatchField.H"
60 #include "cyclicACMIFvPatch.H"
61 
62 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
63 
64 namespace Foam
65 {
66 
67 /*---------------------------------------------------------------------------*\
68  Class cyclicACMIFvPatchField Declaration
69 \*---------------------------------------------------------------------------*/
70 
71 template<class Type>
73 :
74  virtual public cyclicACMILduInterfaceField,
75  public coupledFvPatchField<Type>
76 {
77  // Private data
78 
79  //- Local reference cast into the cyclic patch
80  const cyclicACMIFvPatch& cyclicACMIPatch_;
81 
82 
83  // Private Member Functions
84 
85  //- Return neighbour side field given internal fields
86  template<class Type2>
87  tmp<Field<Type2>> neighbourSideField
88  (
89  const Field<Type2>&
90  ) const;
91 
92  //- Return new matrix coeffs
93  tmp<Field<scalar>> coeffs
94  (
95  fvMatrix<Type>& matrix,
96  const Field<scalar>&,
97  const label
98  ) const;
99 
100 
101 
102 public:
103 
104  //- Runtime type information
105  TypeName(cyclicACMIFvPatch::typeName_());
106 
107 
108  // Constructors
109 
110  //- Construct from patch and internal field
112  (
113  const fvPatch&,
115  );
116 
117  //- Construct from patch, internal field and dictionary
119  (
120  const fvPatch&,
122  const dictionary&
123  );
124 
125  //- Construct by mapping given cyclicACMIFvPatchField onto a new patch
127  (
129  const fvPatch&,
131  const fvPatchFieldMapper&
132  );
133 
134  //- Construct as copy
136 
137  //- Construct and return a clone
138  virtual tmp<fvPatchField<Type>> clone() const
139  {
140  return tmp<fvPatchField<Type>>
141  (
143  );
144  }
145 
146  //- Construct as copy setting internal field reference
148  (
151  );
152 
153  //- Construct and return a clone setting internal field reference
155  (
157  ) const
158  {
159  return tmp<fvPatchField<Type>>
160  (
161  new cyclicACMIFvPatchField<Type>(*this, iF)
162  );
163  }
164 
165 
166  // Member functions
167 
168  // Access
169 
170  //- Return local reference cast into the cyclic AMI patch
171  const cyclicACMIFvPatch& cyclicACMIPatch() const
172  {
173  return cyclicACMIPatch_;
174  }
175 
176 
177  // Evaluation functions
178 
179  //- Return true if coupled. Note that the underlying patch
180  // is not coupled() - the points don't align
181  virtual bool coupled() const;
182 
183  //- Return true if this patch field fixes a value
184  // Needed to check if a level has to be specified while solving
185  // Poisson equations
186  virtual bool fixesValue() const
187  {
188  const scalarField& mask =
189  cyclicACMIPatch_.cyclicACMIPatch().mask();
191  if (gMax(mask) > 1e-5)
192  {
193  // regions connected
194  return false;
195  }
196  else
197  {
198  // fully separated
199  return nonOverlapPatchField().fixesValue();
200  }
201  }
202 
203 
204  //- Return neighbour coupled internal cell data
205  virtual tmp<Field<Type>> patchNeighbourField() const;
206 
207  //- Return reference to neighbour patchField
208  const cyclicACMIFvPatchField<Type>& neighbourPatchField() const;
209 
210  //- Return reference to non-overlapping patchField
212 
213  //- Update result field based on interface functionality
214  virtual void updateInterfaceMatrix
215  (
216  solveScalarField& result,
217  const bool add,
218  const lduAddressing& lduAddr,
219  const label patchId,
220  const solveScalarField& psiInternal,
221  const scalarField& coeffs,
222  const direction cmpt,
223  const Pstream::commsTypes commsType
224  ) const;
225 
226  //- Update result field based on interface functionality
227  virtual void updateInterfaceMatrix
228  (
229  Field<Type>&,
230  const bool add,
231  const lduAddressing& lduAddr,
232  const label patchId,
233  const Field<Type>&,
234  const scalarField&,
235  const Pstream::commsTypes commsType
236  ) const;
237 
238  //- Manipulate matrix
239  virtual void manipulateMatrix(fvMatrix<Type>& matrix);
240 
241  //- Manipulate matrix
242  virtual void manipulateMatrix
243  (
244  fvMatrix<Type>& m,
245  const label iMatrix,
246  const direction cmpt
247  );
248 
249  //- Update the coefficients associated with the patch field
250  virtual void updateCoeffs();
251 
252 
253  // Cyclic AMI coupled interface functions
254 
255  //- Does the patch field perform the transformation
256  virtual bool doTransform() const
257  {
258  return (pTraits<Type>::rank && !cyclicACMIPatch_.parallel());
259  }
260 
261  //- Return face transformation tensor
262  virtual const tensorField& forwardT() const
263  {
264  return cyclicACMIPatch_.forwardT();
265  }
266 
267  //- Return neighbour-cell transformation tensor
268  virtual const tensorField& reverseT() const
269  {
270  return cyclicACMIPatch_.reverseT();
271  }
272 
273  //- Return rank of component for transform
274  virtual int rank() const
275  {
276  return pTraits<Type>::rank;
277  }
278 
279 
280  // I-O
281 
282  //- Write
283  virtual void write(Ostream& os) const;
284 };
285 
286 
287 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
288 
289 } // End namespace Foam
290 
291 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
292 
293 #ifdef NoRepository
294  #include "cyclicACMIFvPatchField.C"
295 #endif
296 
297 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
298 
299 #endif
300 
301 // ************************************************************************* //
label patchId(-1)
virtual int rank() const
Return rank of component for transform.
uint8_t direction
Definition: direction.H:48
TypeName(cyclicACMIFvPatch::typeName_())
Runtime type information.
const fvPatchField< Type > & nonOverlapPatchField() const
Return reference to non-overlapping patchField.
commsTypes
Communications types.
Definition: UPstream.H:74
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:120
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 tmp< Field< Type > > patchNeighbourField() const
Return neighbour coupled internal cell data.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
A traits class, which is primarily used for primitives.
Definition: pTraits.H:51
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 void write(Ostream &os) const
Write.
Generic templated field type.
Definition: Field.H:62
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)
Field< tensor > tensorField
Specialisation of Field<T> for tensor.
const scalarField & mask() const
Mask field where 1 = overlap(coupled), 0 = no-overlap.
const cyclicACMIPolyPatch & cyclicACMIPatch() const
Return local reference cast into the cyclic patch.
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
virtual void manipulateMatrix(fvMatrix< Type > &matrix)
Manipulate matrix.
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.