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 | Copyright (C) 2013-2014 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Class
25  Foam::cyclicACMIFvPatchField
26 
27 Group
28  grpCoupledBoundaryConditions
29 
30 Description
31  This boundary condition enforces a cyclic condition between a pair of
32  boundaries, whereby communication between the patches is performed using
33  an arbitrarily coupled mesh interface (ACMI) interpolation.
34 
35  \heading Patch usage
36 
37  Example of the boundary condition specification:
38  \verbatim
39  myPatch
40  {
41  type cyclicACMI;
42  }
43  \endverbatim
44 
45 SeeAlso
46  Foam::AMIInterpolation
47 
48 SourceFiles
49  cyclicACMIFvPatchField.C
50 
51 \*---------------------------------------------------------------------------*/
52 
53 #ifndef cyclicACMIFvPatchField_H
54 #define cyclicACMIFvPatchField_H
55 
56 #include "coupledFvPatchField.H"
58 #include "cyclicACMIFvPatch.H"
59 
60 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
61 
62 namespace Foam
63 {
64 
65 /*---------------------------------------------------------------------------*\
66  Class cyclicACMIFvPatchField Declaration
67 \*---------------------------------------------------------------------------*/
68 
69 template<class Type>
71 :
72  virtual public cyclicACMILduInterfaceField,
73  public coupledFvPatchField<Type>
74 {
75  // Private data
76 
77  //- Local reference cast into the cyclic patch
79 
80 
81  // Private Member Functions
82 
83  //- Return neighbour side field given internal fields
84  template<class Type2>
86  (
87  const Field<Type2>&
88  ) const;
89 
90 
91 public:
92 
93  //- Runtime type information
94  TypeName(cyclicACMIFvPatch::typeName_());
95 
96 
97  // Constructors
98 
99  //- Construct from patch and internal field
101  (
102  const fvPatch&,
104  );
105 
106  //- Construct from patch, internal field and dictionary
108  (
109  const fvPatch&,
111  const dictionary&
112  );
113 
114  //- Construct by mapping given cyclicACMIFvPatchField onto a new patch
116  (
118  const fvPatch&,
120  const fvPatchFieldMapper&
121  );
122 
123  //- Construct as copy
125 
126  //- Construct and return a clone
127  virtual tmp<fvPatchField<Type> > clone() const
128  {
129  return tmp<fvPatchField<Type> >
130  (
132  );
133  }
134 
135  //- Construct as copy setting internal field reference
137  (
140  );
141 
142  //- Construct and return a clone setting internal field reference
144  (
146  ) const
147  {
148  return tmp<fvPatchField<Type> >
149  (
150  new cyclicACMIFvPatchField<Type>(*this, iF)
151  );
152  }
153 
154 
155  // Member functions
156 
157  // Access
158 
159  //- Return local reference cast into the cyclic AMI patch
160  const cyclicACMIFvPatch& cyclicACMIPatch() const
161  {
162  return cyclicACMIPatch_;
163  }
164 
165 
166  // Evaluation functions
167 
168  //- Return true if coupled. Note that the underlying patch
169  // is not coupled() - the points don't align
170  virtual bool coupled() const;
171 
172  //- Return true if this patch field fixes a value
173  // Needed to check if a level has to be specified while solving
174  // Poissons equations
175  virtual bool fixesValue() const
176  {
177  const scalarField& mask =
179 
180  if (gMax(mask) > 1e-5)
181  {
182  // regions connected
183  return false;
184  }
185  else
186  {
187  // fully separated
188  return nonOverlapPatchField().fixesValue();
189  }
190  }
191 
192 
193  //- Return neighbour coupled internal cell data
194  virtual tmp<Field<Type> > patchNeighbourField() const;
195 
196  //- Return reference to neighbour patchField
198 
199  //- Return reference to non-overlapping patchField
201 
202  //- Update result field based on interface functionality
203  virtual void updateInterfaceMatrix
204  (
205  scalarField& result,
206  const scalarField& psiInternal,
207  const scalarField& coeffs,
208  const direction cmpt,
209  const Pstream::commsTypes commsType
210  ) const;
211 
212  //- Update result field based on interface functionality
213  virtual void updateInterfaceMatrix
214  (
215  Field<Type>&,
216  const Field<Type>&,
217  const scalarField&,
218  const Pstream::commsTypes commsType
219  ) const;
220 
221  //- Return patch-normal gradient
222  virtual tmp<Field<Type> > snGrad
223  (
224  const scalarField& deltaCoeffs
225  ) const;
226 
227  //- Update the coefficients associated with the patch field
228  void updateCoeffs();
229 
230  //- Initialise the evaluation of the patch field
231  virtual void initEvaluate
232  (
233  const Pstream::commsTypes commsType
234  );
235 
236  //- Evaluate the patch field
237  virtual void evaluate
238  (
239  const Pstream::commsTypes commsType
240  );
241 
242  //- Return the matrix diagonal coefficients corresponding to the
243  // evaluation of the value of this patchField with given weights
245  (
246  const tmp<scalarField>&
247  ) const;
248 
249  //- Return the matrix source coefficients corresponding to the
250  // evaluation of the value of this patchField with given weights
252  (
253  const tmp<scalarField>&
254  ) const;
255 
256  //- Return the matrix diagonal coefficients corresponding to the
257  // evaluation of the gradient of this patchField
259  (
260  const scalarField& deltaCoeffs
261  ) const;
262 
263  //- Return the matrix diagonal coefficients corresponding to the
264  // evaluation of the gradient of this patchField
265  virtual tmp<Field<Type> > gradientInternalCoeffs() const;
266 
267  //- Return the matrix source coefficients corresponding to the
268  // evaluation of the gradient of this patchField
270  (
271  const scalarField& deltaCoeffs
272  ) const;
273 
274  //- Return the matrix source coefficients corresponding to the
275  // evaluation of the gradient of this patchField
276  virtual tmp<Field<Type> > gradientBoundaryCoeffs() const;
277 
278  //- Manipulate matrix
279  virtual void manipulateMatrix(fvMatrix<Type>& matrix);
280 
281 
282  // Cyclic AMI coupled interface functions
283 
284  //- Does the patch field perform the transformation
285  virtual bool doTransform() const
286  {
287  return
289  }
290 
291  //- Return face transformation tensor
292  virtual const tensorField& forwardT() const
293  {
294  return cyclicACMIPatch_.forwardT();
295  }
296 
297  //- Return neighbour-cell transformation tensor
298  virtual const tensorField& reverseT() const
299  {
300  return cyclicACMIPatch_.reverseT();
301  }
302 
303  //- Return rank of component for transform
304  virtual int rank() const
305  {
306  return pTraits<Type>::rank;
307  }
308 
309 
310  // I-O
311 
312  //- Write
313  virtual void write(Ostream& os) const;
314 };
315 
316 
317 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
318 
319 } // End namespace Foam
320 
321 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
322 
323 #ifdef NoRepository
324 # include "cyclicACMIFvPatchField.C"
325 #endif
326 
327 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
328 
329 #endif
330 
331 // ************************************************************************* //
Foam::fvPatchField< Type >
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:48
Foam::cyclicACMIFvPatchField::valueInternalCoeffs
virtual tmp< Field< Type > > valueInternalCoeffs(const tmp< scalarField > &) const
Return the matrix diagonal coefficients corresponding to the.
Definition: cyclicACMIFvPatchField.C:335
Foam::cyclicACMIFvPatchField::initEvaluate
virtual void initEvaluate(const Pstream::commsTypes commsType)
Initialise the evaluation of the patch field.
Definition: cyclicACMIFvPatchField.C:288
Foam::cyclicACMIFvPatchField::reverseT
virtual const tensorField & reverseT() const
Return neighbour-cell transformation tensor.
Definition: cyclicACMIFvPatchField.H:297
Foam::cyclicACMIFvPatchField::nonOverlapPatchField
const fvPatchField< Type > & nonOverlapPatchField() const
Return reference to non-overlapping patchField.
Definition: cyclicACMIFvPatchField.C:188
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::cyclicACMIFvPatchField::doTransform
virtual bool doTransform() const
Does the patch field perform the transformation.
Definition: cyclicACMIFvPatchField.H:284
Foam::cyclicACMIFvPatchField::neighbourSideField
tmp< Field< Type2 > > neighbourSideField(const Field< Type2 > &) const
Return neighbour side field given internal fields.
Foam::cyclicACMIFvPatchField::evaluate
virtual void evaluate(const Pstream::commsTypes commsType)
Evaluate the patch field.
Definition: cyclicACMIFvPatchField.C:312
Foam::cyclicACMIFvPatchField
This boundary condition enforces a cyclic condition between a pair of boundaries, whereby communicati...
Definition: cyclicACMIFvPatchField.H:69
Foam::cyclicACMIFvPatchField::gradientInternalCoeffs
virtual tmp< Field< Type > > gradientInternalCoeffs() const
Return the matrix diagonal coefficients corresponding to the.
Definition: cyclicACMIFvPatchField.C:373
Foam::cyclicACMIFvPatch::parallel
virtual bool parallel() const
Are the cyclic planes parallel.
Definition: cyclicACMIFvPatch.H:143
cyclicACMIFvPatchField.C
Foam::fvPatchField::fixesValue
virtual bool fixesValue() const
Return true if this patch field fixes a value.
Definition: fvPatchField.H:336
Foam::cyclicACMILduInterfaceField
Abstract base class for cyclic ACMI coupled interfaces.
Definition: cyclicACMILduInterfaceField.H:48
Foam::cyclicACMIFvPatchField::updateCoeffs
void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: cyclicACMIFvPatchField.C:273
Foam::cyclicACMIFvPatchField::neighbourPatchField
const cyclicACMIFvPatchField< Type > & neighbourPatchField() const
Return reference to neighbour patchField.
Definition: cyclicACMIFvPatchField.C:171
Foam::cyclicACMIFvPatchField::TypeName
TypeName(cyclicACMIFvPatch::typeName_())
Runtime type information.
Foam::cyclicACMIFvPatchField::valueBoundaryCoeffs
virtual tmp< Field< Type > > valueBoundaryCoeffs(const tmp< scalarField > &) const
Return the matrix source coefficients corresponding to the.
Definition: cyclicACMIFvPatchField.C:348
Foam::cyclicACMIFvPatchField::updateInterfaceMatrix
virtual void updateInterfaceMatrix(scalarField &result, const scalarField &psiInternal, const scalarField &coeffs, const direction cmpt, const Pstream::commsTypes commsType) const
Update result field based on interface functionality.
Definition: cyclicACMIFvPatchField.C:202
Foam::cyclicACMIFvPatch::forwardT
virtual const tensorField & forwardT() const
Return face transformation tensor.
Definition: cyclicACMIFvPatch.H:149
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::cyclicACMIFvPatch::reverseT
virtual const tensorField & reverseT() const
Return neighbour-cell transformation tensor.
Definition: cyclicACMIFvPatch.H:155
Foam::cyclicACMIFvPatchField::coupled
virtual bool coupled() const
Return true if coupled. Note that the underlying patch.
Definition: cyclicACMIFvPatchField.C:129
Foam::cyclicACMIFvPatchField::manipulateMatrix
virtual void manipulateMatrix(fvMatrix< Type > &matrix)
Manipulate matrix.
Definition: cyclicACMIFvPatchField.C:406
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Foam::cyclicACMIFvPatchField::cyclicACMIPatch_
const cyclicACMIFvPatch & cyclicACMIPatch_
Local reference cast into the cyclic patch.
Definition: cyclicACMIFvPatchField.H:77
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
cyclicACMILduInterfaceField.H
Foam::cyclicACMIFvPatchField::clone
virtual tmp< fvPatchField< Type > > clone() const
Construct and return a clone.
Definition: cyclicACMIFvPatchField.H:126
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::e
const double e
Elementary charge.
Definition: doubleFloat.H:94
Foam::UPstream::commsTypes
commsTypes
Types of communications.
Definition: UPstream.H:64
Foam::cyclicACMIFvPatchField::rank
virtual int rank() const
Return rank of component for transform.
Definition: cyclicACMIFvPatchField.H:303
Foam::coupledFvPatchField
Abstract base class for coupled patches.
Definition: coupledFvPatchField.H:54
Foam::pTraits
Traits class for primitives.
Definition: pTraits.H:50
Foam::cyclicACMIFvPatchField::write
virtual void write(Ostream &os) const
Write.
Definition: cyclicACMIFvPatchField.C:421
Foam::cyclicACMIFvPatch::cyclicACMIPatch
const cyclicACMIPolyPatch & cyclicACMIPatch() const
Return local reference cast into the cyclic patch.
Definition: cyclicACMIFvPatch.H:94
Foam::cyclicACMIFvPatchField::cyclicACMIFvPatchField
cyclicACMIFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
Definition: cyclicACMIFvPatchField.C:33
Foam::cyclicACMIFvPatchField::patchNeighbourField
virtual tmp< Field< Type > > patchNeighbourField() const
Return neighbour coupled internal cell data.
Definition: cyclicACMIFvPatchField.C:137
Foam::fvMatrix< Type >
Foam::direction
unsigned char direction
Definition: direction.H:43
cyclicACMIFvPatch.H
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:45
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::cyclicACMIFvPatchField::cyclicACMIPatch
const cyclicACMIFvPatch & cyclicACMIPatch() const
Return local reference cast into the cyclic AMI patch.
Definition: cyclicACMIFvPatchField.H:159
coupledFvPatchField.H
Foam::coupledFvPatchField::snGrad
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Definition: coupledFvPatchField.H:147
Foam::cyclicACMIFvPatchField::gradientBoundaryCoeffs
virtual tmp< Field< Type > > gradientBoundaryCoeffs() const
Return the matrix source coefficients corresponding to the.
Definition: cyclicACMIFvPatchField.C:396
Foam::cyclicACMIFvPatchField::forwardT
virtual const tensorField & forwardT() const
Return face transformation tensor.
Definition: cyclicACMIFvPatchField.H:291
Foam::cyclicACMIPolyPatch::mask
const scalarField & mask() const
Mask field where 1 = overlap, 0 = no-overlap.
Definition: cyclicACMIPolyPatchI.H:70
Foam::gMax
Type gMax(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:562
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:51
Foam::cyclicACMIFvPatchField::fixesValue
virtual bool fixesValue() const
Return true if this patch field fixes a value.
Definition: cyclicACMIFvPatchField.H:174
Foam::cyclicACMIFvPatch
Cyclic patch for Arbitrarily Coupled Mesh Interface (ACMI)
Definition: cyclicACMIFvPatch.H:51