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 | Copyright (C) 2011-2016 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 \*---------------------------------------------------------------------------*/
25 
26 #include "IOobject.H"
27 #include "dictionary.H"
28 #include "fvMesh.H"
29 #include "fvPatchFieldMapper.H"
30 #include "volMesh.H"
31 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
34 template<class Type>
36 (
37  const fvPatch& p,
39 )
40 :
41  Field<Type>(p.size()),
42  patch_(p),
43  internalField_(iF),
44  updated_(false),
45  manipulatedMatrix_(false),
46  patchType_(word::null)
47 {}
48 
49 
50 template<class Type>
52 (
53  const fvPatch& p,
55  const Type& value
56 )
57 :
58  Field<Type>(p.size(), value),
59  patch_(p),
60  internalField_(iF),
61  updated_(false),
62  manipulatedMatrix_(false),
63  patchType_(word::null)
64 {}
65 
66 
67 template<class Type>
69 (
70  const fvPatch& p,
72  const word& patchType
73 )
74 :
75  Field<Type>(p.size()),
76  patch_(p),
77  internalField_(iF),
78  updated_(false),
79  manipulatedMatrix_(false),
80  patchType_(patchType)
81 {}
82 
83 
84 template<class Type>
86 (
87  const fvPatch& p,
89  const Field<Type>& f
90 )
91 :
92  Field<Type>(f),
93  patch_(p),
94  internalField_(iF),
95  updated_(false),
96  manipulatedMatrix_(false),
97  patchType_(word::null)
98 {}
99 
100 
101 template<class Type>
103 (
104  const fvPatch& p,
106  const dictionary& dict,
107  const bool valueRequired
108 )
109 :
110  Field<Type>(p.size()),
111  patch_(p),
112  internalField_(iF),
113  updated_(false),
114  manipulatedMatrix_(false),
115  patchType_(dict.lookupOrDefault<word>("patchType", word::null))
116 {
117  if (dict.found("value"))
118  {
120  (
121  Field<Type>("value", dict, p.size())
122  );
123  }
124  else if (!valueRequired)
125  {
127  }
128  else
129  {
131  (
132  dict
133  ) << "Essential entry 'value' missing"
134  << exit(FatalIOError);
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  // For unmapped faces set to internal field value (zero-gradient)
260  if (mapper.hasUnmapped())
261  {
262  Field<Type> pif(this->patchInternalField());
263 
264  if
265  (
266  mapper.direct()
267  && notNull(mapper.directAddressing())
268  && mapper.directAddressing().size()
269  )
270  {
271  const labelList& mapAddressing = mapper.directAddressing();
272 
273  forAll(mapAddressing, i)
274  {
275  if (mapAddressing[i] < 0)
276  {
277  f[i] = pif[i];
278  }
279  }
280  }
281  else if (!mapper.direct() && mapper.addressing().size())
282  {
283  const labelListList& mapAddressing = mapper.addressing();
284 
285  forAll(mapAddressing, i)
286  {
287  const labelList& localAddrs = mapAddressing[i];
288 
289  if (!localAddrs.size())
290  {
291  f[i] = pif[i];
292  }
293  }
294  }
295  }
296  }
297 }
298 
299 
300 template<class Type>
302 (
303  const fvPatchField<Type>& ptf,
304  const labelList& addr
305 )
306 {
307  Field<Type>::rmap(ptf, addr);
308 }
309 
310 
311 template<class Type>
313 {
314  updated_ = true;
315 }
316 
317 
318 template<class Type>
320 {
321  // Default behaviour ignores the weights
322  if (!updated_)
323  {
324  updateCoeffs();
325 
326  updated_ = true;
327  }
328 }
329 
330 
331 template<class Type>
333 {
334  if (!updated_)
335  {
336  updateCoeffs();
337  }
338 
339  updated_ = false;
340  manipulatedMatrix_ = false;
341 }
342 
343 
344 template<class Type>
346 {
347  manipulatedMatrix_ = true;
348 }
349 
350 
351 template<class Type>
353 (
354  fvMatrix<Type>& matrix,
355  const scalarField& weights
356 )
357 {
358  manipulatedMatrix_ = true;
359 }
360 
361 
362 template<class Type>
364 {
365  os.writeKeyword("type") << type() << token::END_STATEMENT << nl;
366 
367  if (patchType_.size())
368  {
369  os.writeKeyword("patchType") << patchType_
370  << token::END_STATEMENT << nl;
371  }
372 }
373 
374 
375 template<class Type>
376 template<class EntryType>
378 (
379  Ostream& os,
380  const word& entryName,
381  const EntryType& value1,
382  const EntryType& value2
383 ) const
384 {
385  if (value1 != value2)
386  {
387  os.writeKeyword(entryName) << value2 << token::END_STATEMENT << nl;
388  }
389 }
390 
391 
392 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
393 
394 template<class Type>
395 void Foam::fvPatchField<Type>::operator=
396 (
397  const UList<Type>& ul
398 )
399 {
401 }
402 
403 
404 template<class Type>
405 void Foam::fvPatchField<Type>::operator=
406 (
407  const fvPatchField<Type>& ptf
408 )
409 {
410  check(ptf);
412 }
413 
414 
415 template<class Type>
416 void Foam::fvPatchField<Type>::operator+=
417 (
418  const fvPatchField<Type>& ptf
419 )
420 {
421  check(ptf);
423 }
424 
425 
426 template<class Type>
427 void Foam::fvPatchField<Type>::operator-=
428 (
429  const fvPatchField<Type>& ptf
430 )
431 {
432  check(ptf);
434 }
435 
436 
437 template<class Type>
438 void Foam::fvPatchField<Type>::operator*=
439 (
440  const fvPatchField<scalar>& ptf
441 )
442 {
443  if (&patch_ != &ptf.patch())
444  {
446  << "incompatible patches for patch fields"
447  << abort(FatalError);
448  }
449 
451 }
452 
453 
454 template<class Type>
455 void Foam::fvPatchField<Type>::operator/=
456 (
457  const fvPatchField<scalar>& ptf
458 )
459 {
460  if (&patch_ != &ptf.patch())
461  {
463  << abort(FatalError);
464  }
465 
467 }
468 
469 
470 template<class Type>
471 void Foam::fvPatchField<Type>::operator+=
472 (
473  const Field<Type>& tf
474 )
475 {
477 }
478 
479 
480 template<class Type>
481 void Foam::fvPatchField<Type>::operator-=
482 (
483  const Field<Type>& tf
484 )
485 {
487 }
488 
489 
490 template<class Type>
491 void Foam::fvPatchField<Type>::operator*=
492 (
493  const scalarField& tf
494 )
495 {
497 }
498 
499 
500 template<class Type>
501 void Foam::fvPatchField<Type>::operator/=
502 (
503  const scalarField& tf
504 )
505 {
507 }
508 
509 
510 template<class Type>
511 void Foam::fvPatchField<Type>::operator=
512 (
513  const Type& t
514 )
515 {
517 }
518 
519 
520 template<class Type>
521 void Foam::fvPatchField<Type>::operator+=
522 (
523  const Type& t
524 )
525 {
527 }
528 
529 
530 template<class Type>
531 void Foam::fvPatchField<Type>::operator-=
532 (
533  const Type& t
534 )
535 {
537 }
538 
539 
540 template<class Type>
541 void Foam::fvPatchField<Type>::operator*=
542 (
543  const scalar s
544 )
545 {
547 }
548 
549 
550 template<class Type>
551 void Foam::fvPatchField<Type>::operator/=
552 (
553  const scalar s
554 )
555 {
557 }
558 
559 
560 template<class Type>
561 void Foam::fvPatchField<Type>::operator==
562 (
563  const fvPatchField<Type>& ptf
564 )
565 {
567 }
568 
569 
570 template<class Type>
571 void Foam::fvPatchField<Type>::operator==
572 (
573  const Field<Type>& tf
574 )
575 {
577 }
578 
579 
580 template<class Type>
581 void Foam::fvPatchField<Type>::operator==
582 (
583  const Type& t
584 )
585 {
587 }
588 
589 
590 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
591 
592 template<class Type>
593 Foam::Ostream& Foam::operator<<(Ostream& os, const fvPatchField<Type>& ptf)
594 {
595  ptf.write(os);
596 
597  os.check("Ostream& operator<<(Ostream&, const fvPatchField<Type>&");
598 
599  return os;
600 }
601 
602 
603 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
604 
605  #include "fvPatchFieldNew.C"
606 
607 // ************************************************************************* //
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:339
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::blocking)
Evaluate the patch field, sets Updated to false.
Definition: fvPatchField.C:332
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Definition: fvPatchField.C:240
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:225
commsTypes
Types of communications.
Definition: UPstream.H:64
void size(const label)
Override size to be inconsistent with allocated storage.
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:92
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
virtual bool direct() const =0
virtual bool distributed() const
Definition: FieldMapper.H:68
const tensorField & tf
virtual const labelUList & directAddressing() const
Definition: FieldMapper.H:85
virtual const labelListList & addressing() const
Definition: FieldMapper.H:94
virtual label size() const =0
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){gmvFile<< iter().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))
Pre-declare SubField and related Field type.
Definition: Field.H:57
A class for handling words, derived from string.
Definition: word.H:59
Foam::fvPatchFieldMapper.
virtual label size() const
Return size.
Definition: fvPatch.H:161
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Definition: fvPatchField.C:302
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:71
virtual void manipulateMatrix(fvMatrix< Type > &matrix)
Manipulate matrix.
Definition: fvPatchField.C:345
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:306
static const zero Zero
Definition: zero.H:91
errorManip< error > abort(error &err)
Definition: errorManip.H:131
void check(const fvPatchField< Type > &) const
Check fvPatchField<Type> against given fvPatchField<Type>
Definition: fvPatchField.C:205
virtual bool hasUnmapped() const =0
Are there unmapped values? I.e. do all size() elements get.
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Definition: fvPatchField.C:217
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
static const char nl
Definition: Ostream.H:262
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:54
labelList f(nPoints)
const objectRegistry & db() const
Return local objectRegistry.
Definition: fvPatchField.C:198
bool notNull(const T &t)
Return true if t is not a reference to the nullObject of type T.
Definition: nullObjectI.H:46
label size() const
Return the number of elements in the UList.
Definition: UListI.H:299
void setSize(const label)
Reset size of List.
void writeEntryIfDifferent(Ostream &os, const word &entryName, const EntryType &value1, const EntryType &value2) const
Helper function to write the keyword and entry only if the.
Definition: fvPatchField.C:378
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:461
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:331
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:312
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
virtual void updateWeightedCoeffs(const scalarField &weights)
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:319
A class for managing temporary objects.
Definition: PtrList.H:54
Registry of regIOobjects.
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
fvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
Definition: fvPatchField.C:36
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:363
IOerror FatalIOError