fvPatchField.C
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) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2015-2021 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 \*---------------------------------------------------------------------------*/
28 
29 #include "IOobject.H"
30 #include "dictionary.H"
31 #include "fvMesh.H"
32 #include "fvPatchFieldMapper.H"
33 #include "volMesh.H"
34 
35 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
36 
37 template<class Type>
39 (
40  const fvPatch& p,
42 )
43 :
44  Field<Type>(p.size()),
45  patch_(p),
46  internalField_(iF),
47  updated_(false),
48  manipulatedMatrix_(false),
49  useImplicit_(false),
50  patchType_()
51 {}
52 
53 
54 template<class Type>
56 (
57  const fvPatch& p,
59  const Type& value
60 )
61 :
62  Field<Type>(p.size(), value),
63  patch_(p),
64  internalField_(iF),
65  updated_(false),
66  manipulatedMatrix_(false),
67  useImplicit_(false),
68  patchType_()
69 {}
70 
71 
72 template<class Type>
74 (
75  const fvPatch& p,
77  const word& patchType
78 )
79 :
80  Field<Type>(p.size()),
81  patch_(p),
82  internalField_(iF),
83  updated_(false),
84  manipulatedMatrix_(false),
85  useImplicit_(false),
86  patchType_(patchType)
87 {}
88 
89 
90 template<class Type>
92 (
93  const fvPatch& p,
95  const Field<Type>& f
96 )
97 :
98  Field<Type>(f),
99  patch_(p),
100  internalField_(iF),
101  updated_(false),
102  manipulatedMatrix_(false),
103  useImplicit_(false),
104  patchType_()
105 {}
106 
107 
108 template<class Type>
110 (
111  const fvPatch& p,
113  const dictionary& dict,
114  const bool valueRequired
115 )
116 :
117  Field<Type>(p.size()),
118  patch_(p),
119  internalField_(iF),
120  updated_(false),
121  manipulatedMatrix_(false),
122  useImplicit_(dict.getOrDefault<bool>("useImplicit", false)),
123  patchType_(dict.getOrDefault<word>("patchType", word::null))
124 {
125  if (valueRequired)
126  {
127  if (dict.found("value"))
128  {
130  (
131  Field<Type>("value", dict, p.size())
132  );
133  }
134  else
135  {
137  << "Essential entry 'value' missing on patch " << p.name() << nl
138  << exit(FatalIOError);
139  }
140  }
141 }
142 
143 
144 template<class Type>
146 (
147  const fvPatchField<Type>& ptf,
148  const fvPatch& p,
150  const fvPatchFieldMapper& mapper
151 )
152 :
153  Field<Type>(p.size()),
154  patch_(p),
155  internalField_(iF),
156  updated_(false),
157  manipulatedMatrix_(false),
158  useImplicit_(ptf.useImplicit_),
159  patchType_(ptf.patchType_)
160 {
161  // For unmapped faces set to internal field value (zero-gradient)
162  if (notNull(iF) && mapper.hasUnmapped())
163  {
164  fvPatchField<Type>::operator=(this->patchInternalField());
165  }
166  this->map(ptf, mapper);
167 }
168 
169 
170 template<class Type>
172 (
173  const fvPatchField<Type>& ptf
174 )
175 :
176  Field<Type>(ptf),
177  patch_(ptf.patch_),
178  internalField_(ptf.internalField_),
179  updated_(false),
180  manipulatedMatrix_(false),
181  useImplicit_(ptf.useImplicit_),
182  patchType_(ptf.patchType_)
183 {}
184 
185 
186 template<class Type>
188 (
189  const fvPatchField<Type>& ptf,
191 )
192 :
193  Field<Type>(ptf),
194  patch_(ptf.patch_),
195  internalField_(iF),
196  updated_(false),
197  manipulatedMatrix_(false),
198  useImplicit_(ptf.useImplicit_),
199  patchType_(ptf.patchType_)
200 {}
201 
202 
203 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
204 
205 template<class Type>
207 {
208  return patch_.boundaryMesh().mesh();
209 }
210 
211 
212 template<class Type>
214 {
215  if (&patch_ != &(ptf.patch_))
216  {
218  << "different patches for fvPatchField<Type>s"
219  << abort(FatalError);
220  }
221 }
222 
223 
224 template<class Type>
226 {
227  return patch_.deltaCoeffs()*(*this - patchInternalField());
228 }
229 
230 
231 template<class Type>
234 {
235  return patch_.patchInternalField(internalField_);
236 }
237 
238 
239 template<class Type>
241 {
242  patch_.patchInternalField(internalField_, pif);
243 }
244 
245 
246 template<class Type>
248 (
249  const fvPatchFieldMapper& mapper
250 )
251 {
252  Field<Type>& f = *this;
253 
254  if (!this->size() && !mapper.distributed())
255  {
256  f.setSize(mapper.size());
257  if (f.size())
258  {
259  f = this->patchInternalField();
260  }
261  }
262  else
263  {
264  // Map all faces provided with mapping data
265  Field<Type>::autoMap(mapper);
266 
267 
268  // For unmapped faces set to internal field value (zero-gradient)
269  if (mapper.hasUnmapped())
270  {
271  Field<Type> pif(this->patchInternalField());
272 
273  if
274  (
275  mapper.direct()
276  && notNull(mapper.directAddressing())
277  && mapper.directAddressing().size()
278  )
279  {
280  const labelList& mapAddressing = mapper.directAddressing();
281 
282  forAll(mapAddressing, i)
283  {
284  if (mapAddressing[i] < 0)
285  {
286  f[i] = pif[i];
287  }
288  }
289  }
290  else if (!mapper.direct() && mapper.addressing().size())
291  {
292  const labelListList& mapAddressing = mapper.addressing();
293 
294  forAll(mapAddressing, i)
295  {
296  const labelList& localAddrs = mapAddressing[i];
297 
298  if (!localAddrs.size())
299  {
300  f[i] = pif[i];
301  }
302  }
303  }
304  }
305  }
306 }
307 
308 
309 template<class Type>
311 (
312  const fvPatchField<Type>& ptf,
313  const labelList& addr
314 )
315 {
316  Field<Type>::rmap(ptf, addr);
317 }
318 
319 
320 template<class Type>
322 {
323  updated_ = true;
324 }
325 
326 
327 template<class Type>
329 {
330  // Default behaviour ignores the weights
331  if (!updated_)
332  {
333  updateCoeffs();
334 
335  updated_ = true;
336  }
337 }
338 
339 
340 template<class Type>
342 {
343  if (!updated_)
344  {
345  updateCoeffs();
346  }
347 
348  updated_ = false;
349  manipulatedMatrix_ = false;
350 }
351 
352 
353 template<class Type>
355 {
356  manipulatedMatrix_ = true;
357 }
358 
359 
360 template<class Type>
362 (
363  fvMatrix<Type>& matrix,
364  const scalarField& weights
365 )
366 {
367  manipulatedMatrix_ = true;
368 }
369 
370 
371 template<class Type>
373 (
374  fvMatrix<Type>& matrix,
375  const label iMatrix,
376  const direction cmp
377 )
378 {
379  manipulatedMatrix_ = true;
380 }
381 
382 
383 template<class Type>
385 {
386  os.writeEntry("type", type());
387 
388  if (useImplicit_)
389  {
390  os.writeEntry("useImplicit", "true");
391  }
392 
393  if (patchType_.size())
394  {
395  os.writeEntry("patchType", patchType_);
396  }
397 }
398 
399 
400 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
401 
402 template<class Type>
404 (
405  const UList<Type>& ul
406 )
407 {
409 }
410 
411 
412 template<class Type>
414 (
415  const fvPatchField<Type>& ptf
416 )
417 {
418  check(ptf);
420 }
421 
422 
423 template<class Type>
425 (
426  const fvPatchField<Type>& ptf
427 )
428 {
429  check(ptf);
431 }
432 
433 
434 template<class Type>
436 (
437  const fvPatchField<Type>& ptf
438 )
439 {
440  check(ptf);
442 }
443 
444 
445 template<class Type>
447 (
448  const fvPatchField<scalar>& ptf
449 )
450 {
451  if (&patch_ != &ptf.patch())
452  {
454  << "incompatible patches for patch fields"
455  << abort(FatalError);
456  }
457 
459 }
460 
461 
462 template<class Type>
464 (
465  const fvPatchField<scalar>& ptf
466 )
467 {
468  if (&patch_ != &ptf.patch())
469  {
471  << abort(FatalError);
472  }
473 
475 }
476 
477 
478 template<class Type>
480 (
481  const Field<Type>& tf
482 )
483 {
485 }
486 
487 
488 template<class Type>
490 (
491  const Field<Type>& tf
492 )
493 {
495 }
496 
497 
498 template<class Type>
500 (
501  const scalarField& tf
502 )
503 {
505 }
506 
507 
508 template<class Type>
510 (
511  const scalarField& tf
512 )
513 {
515 }
516 
517 
518 template<class Type>
520 (
521  const Type& t
522 )
523 {
525 }
526 
527 
528 template<class Type>
530 (
531  const Type& t
532 )
533 {
535 }
536 
537 
538 template<class Type>
540 (
541  const Type& t
542 )
543 {
545 }
546 
547 
548 template<class Type>
550 (
551  const scalar s
552 )
553 {
555 }
556 
557 
558 template<class Type>
560 (
561  const scalar s
562 )
563 {
565 }
566 
567 
568 template<class Type>
570 (
571  const fvPatchField<Type>& ptf
572 )
573 {
575 }
576 
577 
578 template<class Type>
580 (
581  const Field<Type>& tf
582 )
583 {
585 }
586 
587 
588 template<class Type>
590 (
591  const Type& t
592 )
593 {
595 }
596 
597 
598 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
599 
600 template<class Type>
602 {
603  ptf.write(os);
604 
605  os.check(FUNCTION_NAME);
606 
607  return os;
608 }
609 
610 
611 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
612 
613 #include "fvPatchFieldNew.C"
614 
615 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: volSurfaceMapping.H:51
Foam::fvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:384
Foam::fvPatchField::snGrad
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Definition: fvPatchField.C:225
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
s
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Definition: gmvOutputSpray.H:25
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::FieldMapper::size
virtual label size() const =0
volMesh.H
Foam::FatalIOError
IOerror FatalIOError
fvPatchFieldMapper.H
Foam::FieldMapper::direct
virtual bool direct() const =0
Foam::FieldMapper::hasUnmapped
virtual bool hasUnmapped() const =0
Are there unmapped values? I.e. do all size() elements get.
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::objectRegistry
Registry of regIOobjects.
Definition: objectRegistry.H:60
Foam::operator<<
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
Definition: boundaryPatch.C:83
Foam::check
static void check(const int retVal, const char *what)
Definition: ptscotchDecomp.C:80
Foam::Field
Generic templated field type.
Definition: Field.H:63
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
Foam::fvPatchField::patchInternalField
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:233
Foam::fvPatchField::evaluate
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
Evaluate the patch field, sets Updated to false.
Definition: fvPatchField.C:341
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::notNull
bool notNull(const T *ptr)
True if ptr is not a pointer (of type T) to the nullObject.
Definition: nullObject.H:207
Foam::fvPatchField::rmap
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Definition: fvPatchField.C:311
IOobject.H
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
os
OBJstream os(runTime.globalPath()/outputName)
Foam::FieldMapper::directAddressing
virtual const labelUList & directAddressing() const
Definition: FieldMapper.H:89
fvMesh.H
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
fvPatchFieldNew.C
Foam::fvPatchField::manipulateMatrix
virtual void manipulateMatrix(fvMatrix< Type > &matrix)
Manipulate matrix.
Definition: fvPatchField.C:354
Foam::fvPatchField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:321
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::UPstream::commsTypes
commsTypes
Types of communications.
Definition: UPstream.H:69
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::FieldMapper::addressing
virtual const labelListList & addressing() const
Definition: FieldMapper.H:98
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::fvPatchField::check
void check(const fvPatchField< Type > &) const
Check fvPatchField<Type> against given fvPatchField<Type>
Definition: fvPatchField.C:213
f
labelList f(nPoints)
Foam::fvPatchField::db
const objectRegistry & db() const
Return local objectRegistry.
Definition: fvPatchField.C:206
Foam::List< label >
Foam::UList< Type >
dictionary.H
Foam::direction
uint8_t direction
Definition: direction.H:52
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:68
FUNCTION_NAME
#define FUNCTION_NAME
Definition: messageStream.H:295
Foam::roots::type
type
Types of root.
Definition: Roots.H:54
Foam::fvPatchField::fvPatchField
fvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
Definition: fvPatchField.C:39
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:473
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::fvPatchField::patch
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:357
Foam::UList::size
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::fvPatchField::updateWeightedCoeffs
virtual void updateWeightedCoeffs(const scalarField &weights)
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:328
Foam::FieldMapper::distributed
virtual bool distributed() const
Definition: FieldMapper.H:71
Foam::fvPatchField::autoMap
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Definition: fvPatchField.C:248
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54