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-2020 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  patchType_(word::null)
50 {}
51 
52 
53 template<class Type>
55 (
56  const fvPatch& p,
58  const Type& value
59 )
60 :
61  Field<Type>(p.size(), value),
62  patch_(p),
63  internalField_(iF),
64  updated_(false),
65  manipulatedMatrix_(false),
66  patchType_(word::null)
67 {}
68 
69 
70 template<class Type>
72 (
73  const fvPatch& p,
75  const word& patchType
76 )
77 :
78  Field<Type>(p.size()),
79  patch_(p),
80  internalField_(iF),
81  updated_(false),
82  manipulatedMatrix_(false),
83  patchType_(patchType)
84 {}
85 
86 
87 template<class Type>
89 (
90  const fvPatch& p,
92  const Field<Type>& f
93 )
94 :
95  Field<Type>(f),
96  patch_(p),
97  internalField_(iF),
98  updated_(false),
99  manipulatedMatrix_(false),
100  patchType_(word::null)
101 {}
102 
103 
104 template<class Type>
106 (
107  const fvPatch& p,
109  const dictionary& dict,
110  const bool valueRequired
111 )
112 :
113  Field<Type>(p.size()),
114  patch_(p),
115  internalField_(iF),
116  updated_(false),
117  manipulatedMatrix_(false),
118  patchType_(dict.getOrDefault<word>("patchType", word::null))
119 {
120  if (valueRequired)
121  {
122  if (dict.found("value"))
123  {
125  (
126  Field<Type>("value", dict, p.size())
127  );
128  }
129  else
130  {
132  << "Essential entry 'value' missing on patch " << p.name() << nl
133  << exit(FatalIOError);
134  }
135  }
136 }
137 
138 
139 template<class Type>
141 (
142  const fvPatchField<Type>& ptf,
143  const fvPatch& p,
145  const fvPatchFieldMapper& mapper
146 )
147 :
148  Field<Type>(p.size()),
149  patch_(p),
150  internalField_(iF),
151  updated_(false),
152  manipulatedMatrix_(false),
153  patchType_(ptf.patchType_)
154 {
155  // For unmapped faces set to internal field value (zero-gradient)
156  if (notNull(iF) && mapper.hasUnmapped())
157  {
158  fvPatchField<Type>::operator=(this->patchInternalField());
159  }
160  this->map(ptf, mapper);
161 }
162 
163 
164 template<class Type>
166 (
167  const fvPatchField<Type>& ptf
168 )
169 :
170  Field<Type>(ptf),
171  patch_(ptf.patch_),
172  internalField_(ptf.internalField_),
173  updated_(false),
174  manipulatedMatrix_(false),
175  patchType_(ptf.patchType_)
176 {}
177 
178 
179 template<class Type>
181 (
182  const fvPatchField<Type>& ptf,
184 )
185 :
186  Field<Type>(ptf),
187  patch_(ptf.patch_),
188  internalField_(iF),
189  updated_(false),
190  manipulatedMatrix_(false),
191  patchType_(ptf.patchType_)
192 {}
193 
194 
195 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
196 
197 template<class Type>
199 {
200  return patch_.boundaryMesh().mesh();
201 }
202 
203 
204 template<class Type>
206 {
207  if (&patch_ != &(ptf.patch_))
208  {
210  << "different patches for fvPatchField<Type>s"
211  << abort(FatalError);
212  }
213 }
214 
215 
216 template<class Type>
218 {
219  return patch_.deltaCoeffs()*(*this - patchInternalField());
220 }
221 
222 
223 template<class Type>
226 {
227  return patch_.patchInternalField(internalField_);
228 }
229 
230 
231 template<class Type>
233 {
234  patch_.patchInternalField(internalField_, pif);
235 }
236 
237 
238 template<class Type>
240 (
241  const fvPatchFieldMapper& mapper
242 )
243 {
244  Field<Type>& f = *this;
245 
246  if (!this->size() && !mapper.distributed())
247  {
248  f.setSize(mapper.size());
249  if (f.size())
250  {
251  f = this->patchInternalField();
252  }
253  }
254  else
255  {
256  // Map all faces provided with mapping data
257  Field<Type>::autoMap(mapper);
258 
259 
260  // For unmapped faces set to internal field value (zero-gradient)
261  if (mapper.hasUnmapped())
262  {
263  Field<Type> pif(this->patchInternalField());
264 
265  if
266  (
267  mapper.direct()
268  && notNull(mapper.directAddressing())
269  && mapper.directAddressing().size()
270  )
271  {
272  const labelList& mapAddressing = mapper.directAddressing();
273 
274  forAll(mapAddressing, i)
275  {
276  if (mapAddressing[i] < 0)
277  {
278  f[i] = pif[i];
279  }
280  }
281  }
282  else if (!mapper.direct() && mapper.addressing().size())
283  {
284  const labelListList& mapAddressing = mapper.addressing();
285 
286  forAll(mapAddressing, i)
287  {
288  const labelList& localAddrs = mapAddressing[i];
289 
290  if (!localAddrs.size())
291  {
292  f[i] = pif[i];
293  }
294  }
295  }
296  }
297  }
298 }
299 
300 
301 template<class Type>
303 (
304  const fvPatchField<Type>& ptf,
305  const labelList& addr
306 )
307 {
308  Field<Type>::rmap(ptf, addr);
309 }
310 
311 
312 template<class Type>
314 {
315  updated_ = true;
316 }
317 
318 
319 template<class Type>
321 {
322  // Default behaviour ignores the weights
323  if (!updated_)
324  {
325  updateCoeffs();
326 
327  updated_ = true;
328  }
329 }
330 
331 
332 template<class Type>
334 {
335  if (!updated_)
336  {
337  updateCoeffs();
338  }
339 
340  updated_ = false;
341  manipulatedMatrix_ = false;
342 }
343 
344 
345 template<class Type>
347 {
348  manipulatedMatrix_ = true;
349 }
350 
351 
352 template<class Type>
354 (
355  fvMatrix<Type>& matrix,
356  const scalarField& weights
357 )
358 {
359  manipulatedMatrix_ = true;
360 }
361 
362 
363 template<class Type>
365 {
366  os.writeEntry("type", type());
367 
368  if (patchType_.size())
369  {
370  os.writeEntry("patchType", patchType_);
371  }
372 }
373 
374 
375 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
376 
377 template<class Type>
379 (
380  const UList<Type>& ul
381 )
382 {
384 }
385 
386 
387 template<class Type>
389 (
390  const fvPatchField<Type>& ptf
391 )
392 {
393  check(ptf);
395 }
396 
397 
398 template<class Type>
400 (
401  const fvPatchField<Type>& ptf
402 )
403 {
404  check(ptf);
406 }
407 
408 
409 template<class Type>
411 (
412  const fvPatchField<Type>& ptf
413 )
414 {
415  check(ptf);
417 }
418 
419 
420 template<class Type>
422 (
423  const fvPatchField<scalar>& ptf
424 )
425 {
426  if (&patch_ != &ptf.patch())
427  {
429  << "incompatible patches for patch fields"
430  << abort(FatalError);
431  }
432 
434 }
435 
436 
437 template<class Type>
439 (
440  const fvPatchField<scalar>& ptf
441 )
442 {
443  if (&patch_ != &ptf.patch())
444  {
446  << abort(FatalError);
447  }
448 
450 }
451 
452 
453 template<class Type>
455 (
456  const Field<Type>& tf
457 )
458 {
460 }
461 
462 
463 template<class Type>
465 (
466  const Field<Type>& tf
467 )
468 {
470 }
471 
472 
473 template<class Type>
475 (
476  const scalarField& tf
477 )
478 {
480 }
481 
482 
483 template<class Type>
485 (
486  const scalarField& tf
487 )
488 {
490 }
491 
492 
493 template<class Type>
495 (
496  const Type& t
497 )
498 {
500 }
501 
502 
503 template<class Type>
505 (
506  const Type& t
507 )
508 {
510 }
511 
512 
513 template<class Type>
515 (
516  const Type& t
517 )
518 {
520 }
521 
522 
523 template<class Type>
525 (
526  const scalar s
527 )
528 {
530 }
531 
532 
533 template<class Type>
535 (
536  const scalar s
537 )
538 {
540 }
541 
542 
543 template<class Type>
545 (
546  const fvPatchField<Type>& ptf
547 )
548 {
550 }
551 
552 
553 template<class Type>
555 (
556  const Field<Type>& tf
557 )
558 {
560 }
561 
562 
563 template<class Type>
565 (
566  const Type& t
567 )
568 {
570 }
571 
572 
573 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
574 
575 template<class Type>
577 {
578  ptf.write(os);
579 
580  os.check(FUNCTION_NAME);
581 
582  return os;
583 }
584 
585 
586 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
587 
588  #include "fvPatchFieldNew.C"
589 
590 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: volSurfaceMapping.H:50
Foam::fvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:364
Foam::fvPatchField::snGrad
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Definition: fvPatchField.C:217
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
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::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:225
Foam::fvPatchField::evaluate
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
Evaluate the patch field, sets Updated to false.
Definition: fvPatchField.C:333
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:303
IOobject.H
Foam::IOstream::check
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:51
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:121
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:346
Foam::fvPatchField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:313
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:381
Foam::FieldMapper::addressing
virtual const labelListList & addressing() const
Definition: FieldMapper.H:98
Foam::nl
constexpr char nl
Definition: Ostream.H:385
Foam::fvPatchField::check
void check(const fvPatchField< Type > &) const
Check fvPatchField<Type> against given fvPatchField<Type>
Definition: fvPatchField.C:205
f
labelList f(nPoints)
Foam::fvPatchField::db
const objectRegistry & db() const
Return local objectRegistry.
Definition: fvPatchField.C:198
Foam::List< label >
Foam::UList< Type >
dictionary.H
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:76
FUNCTION_NAME
#define FUNCTION_NAME
Definition: messageStream.H:270
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:232
Foam::roots::type
type
Types of root.
Definition: Roots.H:54
Foam::UList::size
void size(const label n) noexcept
Override size to be inconsistent with allocated storage.
Definition: UListI.H:360
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:401
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::fvPatchField::patch
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:343
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:320
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:240
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54