immersedBoundaryFvPatchField.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | foam-extend: Open Source CFD
4  \\ / O peration | Version: 3.2
5  \\ / A nd | Web: http://www.foam-extend.org
6  \\/ M anipulation | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
8 License
9  This file is part of foam-extend.
10 
11  foam-extend is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by the
13  Free Software Foundation, either version 3 of the License, or (at your
14  option) any later version.
15 
16  foam-extend is distributed in the hope that it will be useful, but
17  WITHOUT ANY WARRANTY; without even the implied warranty of
18  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19  General Public License for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with foam-extend. If not, see <http://www.gnu.org/licenses/>.
23 
24 Class
25  Foam::immersedBoundaryFvPatchField
26 
27 Description
28  Foam::immersedBoundaryFvPatchField
29 
30 Author
31  Hrvoje Jasak
32 
33 SourceFiles
34  immersedBoundaryFvPatchField.C
35 
36 \*---------------------------------------------------------------------------*/
37 
38 #ifndef immersedBoundaryFvPatchField_H
39 #define immersedBoundaryFvPatchField_H
40 
41 #include "fvPatchField.H"
43 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
46 namespace Foam
47 {
48 
49 /*---------------------------------------------------------------------------*\
50  Class immersedBoundaryFvPatchField Declaration
51 \*---------------------------------------------------------------------------*/
52 
53 template<class Type>
55 :
56  public fvPatchField<Type>
57 {
58  // Private data
59 
60  //- Local reference cast into the processor patch
62 
63  //- Local reference to fvMesh
64  const fvMesh& mesh_;
65 
66  // Defining fields
67  // Note: defining fields carry values on faces of the IB patch
68  // represented as a triangulated surface
69 
70  //- Defining value field
72 
73  //- Defining normal gradient field
75 
76 
77  //- Does the boundary condition fix the value
79 
80 
81  // Dead cell controls
82 
83  //- Set dead cell value
85 
86  //- Dead cell value
87  Type deadCellValue_;
88 
89 
90  // Fields on the IB intersection points
91  // Field size equals the number of cells intesected with the IB
92 
93  //- Time index for last update of mesh or moving boundary
94  mutable label ibUpdateTimeIndex_;
95 
96  //- Value field
97  mutable Field<Type> ibValue_;
98 
99  //- Normal gradient field
100  mutable Field<Type> ibGrad_;
101 
102 
103  // Private Member Functions
104 
105  //- Update IB value and gradient
106  void updateIbValues() const;
107 
108  //- Impose Dirichlet BC at IB cells and return corrected cells values
109  // Calculate value and gradient on IB intersection points
111 
112  //- Impose Neumann BC at IB cells and return corrected cells values
113  // Calculate value and gradient on IB intersection points
115 
116  //- Impose condition in dead cells
117  void imposeDeadCondition();
118 
119 
120  // Manipulate matrix
121 
122  //- Check and correct zero diagonal in dead cells
123  void correctDiag
124  (
125  fvMatrix<Type>& eqn
126  ) const;
127 
128  //- Impose fixed gradient condition by manipulating matrix
129  // Note: reconsider pre-factor for diffusivity
130  // HJ, 16/Apr/2012
131  void correctOffDiag
132  (
133  fvMatrix<Type>& eqn
134  ) const;
135 
136 
137 protected:
138 
139  // Protected Member Functions
140 
141  //- Does immersed boundary patch field need motion update?
142  bool motionUpdateRequired() const;
143 
144  //- Execute immersed boundary patch field motion update
145  virtual void motionUpdate() const;
146 
147  //- Set IB cell values
148  virtual void setIbCellValues(const Field<Type>&) const;
149 
150 
151 public:
152 
153  //- Runtime type information
154  TypeName(immersedBoundaryFvPatch::typeName_());
155 
156 
157  // Static algorithm control data
158 
159  //- Number of boundary condition update iterations
160  static const debug::optimisationSwitch nBcIter_;
161 
162  //- Boundary condition update tolerance
163  static const debug::tolerancesSwitch bcTolerance_;
164 
165 
166  // Constructors
167 
168  //- Construct from patch and internal field
170  (
171  const fvPatch&,
173  );
174 
175  //- Construct from patch, internal field and dictionary
177  (
178  const fvPatch&,
180  const dictionary&
181  );
182 
183  //- Construct by mapping given immersedBoundaryFvPatchField
184  // onto a new patch
186  (
188  const fvPatch&,
190  const fvPatchFieldMapper&
191  );
192 
193  //- Construct as copy
195  (
197  );
198 
199  //- Construct and return a clone
200  virtual tmp<fvPatchField<Type> > clone() const
201  {
202  return tmp<fvPatchField<Type> >
203  (
205  );
206  }
207 
208  //- Construct as copy setting internal field reference
210  (
213  );
214 
215  //- Construct and return a clone setting internal field reference
217  (
219  ) const
220  {
221  return tmp<fvPatchField<Type> >
222  (
224  );
225  }
226 
227 
228  //- Destructor
230  {}
231 
232 
233  // Member functions
234 
235  //- Return reference to immersed boundary patch
236  const immersedBoundaryFvPatch& ibPatch() const
237  {
238  return ibPatch_;
239  }
240 
241 
242  // Return defining fields
243  // Note: defining fields carry values on faces of the IB patch
244  // represented as a triangulated surface
245 
246  //- Return access to reference value
247  virtual Field<Type>& refValue()
248  {
249  return refValue_;
250  }
251 
252  //- Return reference value
253  virtual const Field<Type>& refValue() const
254  {
255  return refValue_;
256  }
257 
258  //- Return access to reference gradient
259  virtual Field<Type>& refGrad()
260  {
261  return refGrad_;
262  }
263 
264  //- Return reference gradient
265  virtual const Field<Type>& refGrad() const
266  {
267  return refGrad_;
268  }
269 
270  //- Return true if this patch field fixes a value
271  virtual bool fixesValue() const
272  {
273  return fixesValue_;
274  }
275 
276  //- Return access to fixes value switch
277  Switch& fixesValue()
278  {
279  return fixesValue_;
280  }
281 
282 
283  //- Return fields on intersection points interacting with the IB
284 
285  //- Return IB field
286  const Field<Type>& ibValue() const;
287 
288  //- Return IB surface-normal gradient
289  const Field<Type>& ibGrad() const;
290 
291  //- Return IB cell field
292  tmp<Field<Type> > ibCellValue() const;
293 
294  //- Return IB sample point field
296 
297 
298  //- Return fields on triangular faces
299 
300  //- Return IB field
301  tmp<Field<Type> > triValue() const;
302 
303  //- Return IB surface-normal gradient
304  tmp<Field<Type> > triGrad() const;
305 
306 
307  // Mapping functions
308 
309  //- Map (and resize as needed) from self given a mapping object
310  virtual void autoMap
311  (
312  const fvPatchFieldMapper&
313  )
314  {}
315 
316  //- Reverse map the given fvPatchField onto this fvPatchField
317  virtual void rmap
318  (
319  const fvPatchField<Type>&,
320  const labelList&
321  )
322  {}
323 
324 
325  // Evaluation functions
326 
327  //- Update the coefficients associated with the patch field
328  void updateCoeffs();
329 
330  //- Initialise the evaluation of the patch field
331  virtual void initEvaluate
332  (
333  const Pstream::commsTypes commsType = Pstream::blocking
334  );
335 
336  //- Evaluate the patch field
337  virtual void evaluate
338  (
339  const Pstream::commsTypes commsType = Pstream::blocking
340  );
341 
342  //- Return the matrix diagonal coefficients corresponding to the
343  // evaluation of the value of this patchField with given weights
345  (
346  const tmp<scalarField>&
347  ) const
348  {
349  return tmp<Field<Type> >(new Field<Type>(0));
350  }
351 
352  //- Return the matrix source coefficients corresponding to the
353  // evaluation of the value of this patchField with given weights
355  (
356  const tmp<scalarField>&
357  ) const
358  {
359  return tmp<Field<Type> >(new Field<Type>(0));
360  }
361 
362  //- Return the matrix diagonal coefficients corresponding to the
363  // evaluation of the gradient of this patchField
365  {
366  return tmp<Field<Type> >(new Field<Type>(0));
367  }
368 
369  //- Return the matrix source coefficients corresponding to the
370  // evaluation of the gradient of this patchField
372  {
373  return tmp<Field<Type> >(new Field<Type>(0));
374  }
375 
376  //- Manipulate matrix
377  virtual void manipulateMatrix(fvMatrix<Type>& matrix);
378 
379 
380  // I-O
381 
382  //- Write
383  virtual void write(Ostream&) const;
384 };
385 
386 
387 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
388 
389 } // End namespace Foam
390 
391 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
392 
393 #ifdef NoRepository
395 #endif
396 
397 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
398 
399 #endif
400 
401 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
Foam::immersedBoundaryFvPatchField::manipulateMatrix
virtual void manipulateMatrix(fvMatrix< Type > &matrix)
Manipulate matrix.
Definition: immersedBoundaryFvPatchField.C:1034
Foam::Switch
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:60
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::immersedBoundaryFvPatchField::fixesValue_
Switch fixesValue_
Does the boundary condition fix the value.
Definition: immersedBoundaryFvPatchField.H:77
Foam::immersedBoundaryFvPatchField::gradientInternalCoeffs
tmp< Field< Type > > gradientInternalCoeffs() const
Return the matrix diagonal coefficients corresponding to the.
Definition: immersedBoundaryFvPatchField.H:363
Foam::immersedBoundaryFvPatchField::evaluate
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::blocking)
Evaluate the patch field.
Definition: immersedBoundaryFvPatchField.C:1017
Foam::immersedBoundaryFvPatchField::refValue_
Field< Type > refValue_
Defining value field.
Definition: immersedBoundaryFvPatchField.H:70
Foam::immersedBoundaryFvPatchField::initEvaluate
virtual void initEvaluate(const Pstream::commsTypes commsType=Pstream::blocking)
Initialise the evaluation of the patch field.
Definition: immersedBoundaryFvPatchField.C:1007
Foam::immersedBoundaryFvPatchField::refGrad
virtual Field< Type > & refGrad()
Return access to reference gradient.
Definition: immersedBoundaryFvPatchField.H:258
Foam::immersedBoundaryFvPatchField::mesh_
const fvMesh & mesh_
Local reference to fvMesh.
Definition: immersedBoundaryFvPatchField.H:63
Foam::immersedBoundaryFvPatchField::~immersedBoundaryFvPatchField
virtual ~immersedBoundaryFvPatchField()
Destructor.
Definition: immersedBoundaryFvPatchField.H:228
Foam::immersedBoundaryFvPatchField::fixesValue
Switch & fixesValue()
Return access to fixes value switch.
Definition: immersedBoundaryFvPatchField.H:276
Foam::immersedBoundaryFvPatchField::ibValue
const Field< Type > & ibValue() const
Return fields on intersection points interacting with the IB.
Definition: immersedBoundaryFvPatchField.C:905
Foam::immersedBoundaryFvPatchField::valueInternalCoeffs
virtual tmp< Field< Type > > valueInternalCoeffs(const tmp< scalarField > &) const
Return the matrix diagonal coefficients corresponding to the.
Definition: immersedBoundaryFvPatchField.H:344
Foam::immersedBoundaryFvPatchField::motionUpdateRequired
bool motionUpdateRequired() const
Does immersed boundary patch field need motion update?
Definition: immersedBoundaryFvPatchField.C:688
Foam::immersedBoundaryFvPatchField::immersedBoundaryFvPatchField
immersedBoundaryFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
Definition: immersedBoundaryFvPatchField.C:760
Foam::immersedBoundaryFvPatchField::refValue
virtual Field< Type > & refValue()
Return access to reference value.
Definition: immersedBoundaryFvPatchField.H:246
Foam::UPstream::blocking
@ blocking
Definition: UPstream.H:66
Foam::immersedBoundaryFvPatchField::imposeDeadCondition
void imposeDeadCondition()
Impose condition in dead cells.
Definition: immersedBoundaryFvPatchField.C:474
Foam::immersedBoundaryFvPatchField::updateCoeffs
void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: immersedBoundaryFvPatchField.C:984
Foam::immersedBoundaryFvPatchField::updateIbValues
void updateIbValues() const
Update IB value and gradient.
Definition: immersedBoundaryFvPatchField.C:59
Foam::debug::optimisationSwitch
int optimisationSwitch(const char *name, const int defaultValue=0)
Lookup optimisation switch or add default value.
Definition: debug.C:182
Foam::immersedBoundaryFvPatchField::refValue
virtual const Field< Type > & refValue() const
Return reference value.
Definition: immersedBoundaryFvPatchField.H:252
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::Field< Type >
Foam::immersedBoundaryFvPatch
Immersed boundary FV patch.
Definition: immersedBoundaryFvPatch.H:66
Foam::immersedBoundaryFvPatchField::ibPatch
const immersedBoundaryFvPatch & ibPatch() const
Return reference to immersed boundary patch.
Definition: immersedBoundaryFvPatchField.H:235
Foam::immersedBoundaryFvPatchField::ibGrad
const Field< Type > & ibGrad() const
Return IB surface-normal gradient.
Definition: immersedBoundaryFvPatchField.C:925
Foam::immersedBoundaryFvPatchField::fixesValue
virtual bool fixesValue() const
Return true if this patch field fixes a value.
Definition: immersedBoundaryFvPatchField.H:270
Foam::immersedBoundaryFvPatchField::correctOffDiag
void correctOffDiag(fvMatrix< Type > &eqn) const
Impose fixed gradient condition by manipulating matrix.
Definition: immersedBoundaryFvPatchField.C:526
Foam::immersedBoundaryFvPatchField::imposeDirichletCondition
tmp< Field< Type > > imposeDirichletCondition() const
Impose Dirichlet BC at IB cells and return corrected cells values.
Definition: immersedBoundaryFvPatchField.C:90
Foam::immersedBoundaryFvPatchField::ibCellValue
tmp< Field< Type > > ibCellValue() const
Return IB cell field.
Definition: immersedBoundaryFvPatchField.C:945
Foam::immersedBoundaryFvPatchField::deadCellValue_
Type deadCellValue_
Dead cell value.
Definition: immersedBoundaryFvPatchField.H:86
Foam::immersedBoundaryFvPatchField::ibSamplingPointValue
tmp< Field< Type > > ibSamplingPointValue() const
Return IB sample point field.
Definition: immersedBoundaryFvPatchField.C:963
Foam::immersedBoundaryFvPatchField::ibGrad_
Field< Type > ibGrad_
Normal gradient field.
Definition: immersedBoundaryFvPatchField.H:99
Foam::immersedBoundaryFvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: immersedBoundaryFvPatchField.C:1066
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
immersedBoundaryFvPatchField.C
Foam::immersedBoundaryFvPatchField::triGrad
tmp< Field< Type > > triGrad() const
Return IB surface-normal gradient.
Definition: immersedBoundaryFvPatchField.C:977
Foam::immersedBoundaryFvPatchField::rmap
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Definition: immersedBoundaryFvPatchField.H:317
Foam::immersedBoundaryFvPatchField::refGrad_
Field< Type > refGrad_
Defining normal gradient field.
Definition: immersedBoundaryFvPatchField.H:73
Foam::immersedBoundaryFvPatchField::clone
virtual tmp< fvPatchField< Type > > clone() const
Construct and return a clone.
Definition: immersedBoundaryFvPatchField.H:199
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::immersedBoundaryFvPatchField::autoMap
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Definition: immersedBoundaryFvPatchField.H:310
Foam::immersedBoundaryFvPatchField::TypeName
TypeName(immersedBoundaryFvPatch::typeName_())
Runtime type information.
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::immersedBoundaryFvPatchField::bcTolerance_
static const debug::tolerancesSwitch bcTolerance_
Boundary condition update tolerance.
Definition: immersedBoundaryFvPatchField.H:162
Foam::immersedBoundaryFvPatchField::setDeadCellValue_
Switch setDeadCellValue_
Set dead cell value.
Definition: immersedBoundaryFvPatchField.H:83
Foam::immersedBoundaryFvPatchField::nBcIter_
static const debug::optimisationSwitch nBcIter_
Number of boundary condition update iterations.
Definition: immersedBoundaryFvPatchField.H:159
Foam::immersedBoundaryFvPatchField::correctDiag
void correctDiag(fvMatrix< Type > &eqn) const
Check and correct zero diagonal in dead cells.
Definition: immersedBoundaryFvPatchField.C:495
Foam::UPstream::commsTypes
commsTypes
Types of communications.
Definition: UPstream.H:64
Foam::immersedBoundaryFvPatchField::gradientBoundaryCoeffs
tmp< Field< Type > > gradientBoundaryCoeffs() const
Return the matrix source coefficients corresponding to the.
Definition: immersedBoundaryFvPatchField.H:370
Foam::immersedBoundaryFvPatchField::setIbCellValues
virtual void setIbCellValues(const Field< Type > &) const
Set IB cell values.
Definition: immersedBoundaryFvPatchField.C:725
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::immersedBoundaryFvPatchField::ibValue_
Field< Type > ibValue_
Value field.
Definition: immersedBoundaryFvPatchField.H:96
Foam::immersedBoundaryFvPatchField::valueBoundaryCoeffs
virtual tmp< Field< Type > > valueBoundaryCoeffs(const tmp< scalarField > &) const
Return the matrix source coefficients corresponding to the.
Definition: immersedBoundaryFvPatchField.H:354
Foam::immersedBoundaryFvPatchField
Foam::immersedBoundaryFvPatchField.
Definition: immersedBoundaryFvPatchField.H:53
Foam::fvMatrix< Type >
Foam::immersedBoundaryFvPatchField::imposeNeumannCondition
tmp< Field< Type > > imposeNeumannCondition() const
Impose Neumann BC at IB cells and return corrected cells values.
Definition: immersedBoundaryFvPatchField.C:284
Foam::immersedBoundaryFvPatchField::ibPatch_
const immersedBoundaryFvPatch & ibPatch_
Local reference cast into the processor patch.
Definition: immersedBoundaryFvPatchField.H:60
immersedBoundaryFvPatch.H
Foam::immersedBoundaryFvPatchField::ibUpdateTimeIndex_
label ibUpdateTimeIndex_
Time index for last update of mesh or moving boundary.
Definition: immersedBoundaryFvPatchField.H:93
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:45
Foam::immersedBoundaryFvPatchField::triValue
tmp< Field< Type > > triValue() const
Return fields on triangular faces.
Definition: immersedBoundaryFvPatchField.C:970
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::immersedBoundaryFvPatchField::refGrad
virtual const Field< Type > & refGrad() const
Return reference gradient.
Definition: immersedBoundaryFvPatchField.H:264
Foam::immersedBoundaryFvPatchField::motionUpdate
virtual void motionUpdate() const
Execute immersed boundary patch field motion update.
Definition: immersedBoundaryFvPatchField.C:712
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:51
fvPatchField.H