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-2015 OpenFOAM Foundation
6  \\/ M anipulation | Copyright (C) 2015 OpenCFD Ltd.
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,
38  const DimensionedField<Type, volMesh>& iF
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,
54  const DimensionedField<Type, volMesh>& iF,
55  const word& patchType
56 )
57 :
58  Field<Type>(p.size()),
59  patch_(p),
60  internalField_(iF),
61  updated_(false),
62  manipulatedMatrix_(false),
63  patchType_(patchType)
64 {}
65 
66 
67 template<class Type>
69 (
70  const fvPatch& p,
71  const DimensionedField<Type, volMesh>& iF,
72  const Field<Type>& f
73 )
74 :
75  Field<Type>(f),
76  patch_(p),
77  internalField_(iF),
78  updated_(false),
79  manipulatedMatrix_(false),
80  patchType_(word::null)
81 {}
82 
83 
84 template<class Type>
86 (
87  const fvPatchField<Type>& ptf,
88  const fvPatch& p,
89  const DimensionedField<Type, volMesh>& iF,
90  const fvPatchFieldMapper& mapper
91 )
92 :
93  Field<Type>(p.size()),
94  patch_(p),
95  internalField_(iF),
96  updated_(false),
97  manipulatedMatrix_(false),
98  patchType_(ptf.patchType_)
99 {
100  // For unmapped faces set to internal field value (zero-gradient)
101  if (notNull(iF) && mapper.hasUnmapped())
102  {
103  fvPatchField<Type>::operator=(this->patchInternalField());
104  }
105  this->map(ptf, mapper);
106 }
107 
108 
109 template<class Type>
111 (
112  const fvPatch& p,
113  const DimensionedField<Type, volMesh>& iF,
114  const dictionary& dict,
115  const bool valueRequired
116 )
117 :
118  Field<Type>(p.size()),
119  patch_(p),
120  internalField_(iF),
121  updated_(false),
122  manipulatedMatrix_(false),
123  patchType_(dict.lookupOrDefault<word>("patchType", word::null))
124 {
125  if (dict.found("value"))
126  {
128  (
129  Field<Type>("value", dict, p.size())
130  );
131  }
132  else if (!valueRequired)
133  {
134  fvPatchField<Type>::operator=(pTraits<Type>::zero);
135  }
136  else
137  {
139  (
140  dict
141  ) << "Essential entry 'value' missing"
142  << exit(FatalIOError);
143  }
144 }
145 
146 
147 template<class Type>
149 (
150  const fvPatchField<Type>& ptf
151 )
152 :
153  Field<Type>(ptf),
154  patch_(ptf.patch_),
155  internalField_(ptf.internalField_),
156  updated_(false),
157  manipulatedMatrix_(false),
158  patchType_(ptf.patchType_)
159 {}
160 
161 
162 template<class Type>
164 (
165  const fvPatchField<Type>& ptf,
166  const DimensionedField<Type, volMesh>& iF
167 )
168 :
169  Field<Type>(ptf),
170  patch_(ptf.patch_),
171  internalField_(iF),
172  updated_(false),
173  manipulatedMatrix_(false),
174  patchType_(ptf.patchType_)
175 {}
176 
177 
178 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
179 
180 template<class Type>
182 {
183  return patch_.boundaryMesh().mesh();
184 }
185 
186 
187 template<class Type>
189 {
190  if (&patch_ != &(ptf.patch_))
191  {
193  << "different patches for fvPatchField<Type>s"
194  << abort(FatalError);
195  }
196 }
197 
198 
199 template<class Type>
201 {
202  return patch_.deltaCoeffs()*(*this - patchInternalField());
203 }
204 
205 
206 template<class Type>
209 {
210  return patch_.patchInternalField(internalField_);
211 }
212 
213 
214 template<class Type>
216 {
217  patch_.patchInternalField(internalField_, pif);
218 }
219 
220 
221 template<class Type>
223 (
224  const fvPatchFieldMapper& mapper
225 )
226 {
227  Field<Type>& f = *this;
228 
229  if (!this->size() && !mapper.distributed())
230  {
231  f.setSize(mapper.size());
232  if (f.size())
233  {
234  f = this->patchInternalField();
235  }
236  }
237  else
238  {
239  // Map all faces provided with mapping data
240  Field<Type>::autoMap(mapper);
241 
242 
243  // For unmapped faces set to internal field value (zero-gradient)
244  if (mapper.hasUnmapped())
245  {
246  Field<Type> pif(this->patchInternalField());
247 
248  if
249  (
250  mapper.direct()
251  && notNull(mapper.directAddressing())
252  && mapper.directAddressing().size()
253  )
254  {
255  const labelList& mapAddressing = mapper.directAddressing();
256 
257  forAll(mapAddressing, i)
258  {
259  if (mapAddressing[i] < 0)
260  {
261  f[i] = pif[i];
262  }
263  }
264  }
265  else if (!mapper.direct() && mapper.addressing().size())
266  {
267  const labelListList& mapAddressing = mapper.addressing();
268 
269  forAll(mapAddressing, i)
270  {
271  const labelList& localAddrs = mapAddressing[i];
272 
273  if (!localAddrs.size())
274  {
275  f[i] = pif[i];
276  }
277  }
278  }
279  }
280  }
281 }
282 
283 
284 template<class Type>
286 (
287  const fvPatchField<Type>& ptf,
288  const labelList& addr
289 )
290 {
291  Field<Type>::rmap(ptf, addr);
292 }
293 
294 
295 template<class Type>
297 {
298  updated_ = true;
299 }
300 
301 
302 template<class Type>
304 {
305  if (!updated_)
306  {
307  updateCoeffs();
308 
309  Field<Type>& fld = *this;
310  fld *= weights;
311 
312  updated_ = true;
313  }
314 }
315 
316 
317 template<class Type>
319 {
320  if (!updated_)
321  {
322  updateCoeffs();
323  }
324 
325  updated_ = false;
326  manipulatedMatrix_ = false;
327 }
328 
329 
330 template<class Type>
332 {
333  manipulatedMatrix_ = true;
334 }
335 
336 
337 template<class Type>
339 (
340  fvMatrix<Type>& matrix,
341  const scalarField& weights
342 )
343 {
344  manipulatedMatrix_ = true;
345 }
346 
347 
348 template<class Type>
350 {
351  os.writeKeyword("type") << type() << token::END_STATEMENT << nl;
352 
353  if (patchType_.size())
354  {
355  os.writeKeyword("patchType") << patchType_
356  << token::END_STATEMENT << nl;
357  }
358 }
359 
360 
361 template<class Type>
362 template<class EntryType>
364 (
365  Ostream& os,
366  const word& entryName,
367  const EntryType& value1,
368  const EntryType& value2
369 ) const
370 {
371  if (value1 != value2)
372  {
373  os.writeKeyword(entryName) << value2 << token::END_STATEMENT << nl;
374  }
375 }
376 
377 
378 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
379 
380 template<class Type>
382 (
383  const UList<Type>& ul
384 )
385 {
387 }
388 
389 
390 template<class Type>
392 (
393  const fvPatchField<Type>& ptf
394 )
395 {
396  check(ptf);
398 }
399 
400 
401 template<class Type>
403 (
404  const fvPatchField<Type>& ptf
405 )
406 {
407  check(ptf);
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<scalar>& ptf
427 )
428 {
429  if (&patch_ != &ptf.patch())
430  {
432  << "incompatible patches for patch fields"
433  << abort(FatalError);
434  }
435 
437 }
438 
439 
440 template<class Type>
442 (
443  const fvPatchField<scalar>& ptf
444 )
445 {
446  if (&patch_ != &ptf.patch())
447  {
449  << abort(FatalError);
450  }
451 
453 }
454 
455 
456 template<class Type>
458 (
459  const Field<Type>& tf
460 )
461 {
463 }
464 
465 
466 template<class Type>
468 (
469  const Field<Type>& tf
470 )
471 {
473 }
474 
475 
476 template<class Type>
478 (
479  const scalarField& tf
480 )
481 {
483 }
484 
485 
486 template<class Type>
488 (
489  const scalarField& tf
490 )
491 {
493 }
494 
495 
496 template<class Type>
498 (
499  const Type& t
500 )
501 {
503 }
504 
505 
506 template<class Type>
508 (
509  const Type& t
510 )
511 {
513 }
514 
515 
516 template<class Type>
518 (
519  const Type& t
520 )
521 {
523 }
524 
525 
526 template<class Type>
528 (
529  const scalar s
530 )
531 {
533 }
534 
535 
536 template<class Type>
538 (
539  const scalar s
540 )
541 {
543 }
544 
545 
546 // Force an assignment, overriding fixedValue status
547 template<class Type>
549 (
550  const fvPatchField<Type>& ptf
551 )
552 {
554 }
555 
556 
557 template<class Type>
559 (
560  const Field<Type>& tf
561 )
562 {
564 }
565 
566 
567 template<class Type>
569 (
570  const Type& t
571 )
572 {
574 }
575 
576 
577 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
578 
579 template<class Type>
580 Foam::Ostream& Foam::operator<<(Ostream& os, const fvPatchField<Type>& ptf)
581 {
582  ptf.write(os);
583 
584  os.check("Ostream& operator<<(Ostream&, const fvPatchField<Type>&");
585 
586  return os;
587 }
588 
589 
590 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
591 
592 # include "fvPatchFieldNew.C"
593 
594 // ************************************************************************* //
Foam::fvPatchField< Type >
Foam::fvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:349
Foam::fvPatchField::operator=
virtual void operator=(const UList< Type > &)
Foam::fvPatchField::snGrad
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Definition: fvPatchField.C:200
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:48
Foam::token::END_STATEMENT
@ END_STATEMENT
Definition: token.H:99
Foam::fvPatchField::evaluate
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::blocking)
Evaluate the patch field, sets Updated to false.
Definition: fvPatchField.C:318
Foam::Field::operator-=
void operator-=(const UList< Type > &)
p
p
Definition: pEqn.H:62
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::Field::autoMap
void autoMap(const FieldMapper &map, const bool applyFlip=true)
Map from self.
Definition: Field.C:508
Foam::notNull
bool notNull(const T &t)
Return true if t is not a reference to the nullObject of type T.
Definition: nullObjectI.H:46
Foam::FieldMapper::size
virtual label size() const =0
volMesh.H
Foam::fvPatchField< Type >::operator
friend Ostream & operator(Ostream &, const fvPatchField< Type > &)
Foam::Field::operator+=
void operator+=(const UList< Type > &)
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.
Foam::fvPatchField::writeEntryIfDifferent
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:364
Foam::objectRegistry
Registry of regIOobjects.
Definition: objectRegistry.H:50
tf
const tensorField & tf
Definition: getPatchFieldTensor.H:36
Foam::Field< Type >
Foam::Ostream::write
virtual Ostream & write(const token &)=0
Write next token to stream.
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::fvPatchField::patchInternalField
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:208
Foam::operator<<
Ostream & operator<<(Ostream &, const edgeMesh &)
Definition: edgeMeshIO.C:130
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::Field::rmap
void rmap(const UList< Type > &mapF, const labelUList &mapAddressing)
1 to 1 reverse-map from the given field
Definition: Field.C:571
Foam::fvPatchField::rmap
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Definition: fvPatchField.C:286
IOobject.H
Foam::IOstream::check
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:92
Foam::FatalError
error FatalError
Foam::fvPatchField::patch_
const fvPatch & patch_
Reference to patch.
Definition: fvPatchField.H:86
fld
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){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Foam::FieldMapper::directAddressing
virtual const labelUList & directAddressing() const
Definition: FieldMapper.H:85
fvMesh.H
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
fvPatchFieldNew.C
Foam::fvPatchField::manipulateMatrix
virtual void manipulateMatrix(fvMatrix< Type > &matrix)
Manipulate matrix.
Definition: fvPatchField.C:331
s
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))
Foam::fvPatchField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:296
Foam::Field::operator=
void operator=(const Field< Type > &)
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::List::setSize
void setSize(const label)
Reset size of List.
Foam::Field::operator/=
void operator/=(const UList< scalar > &)
Foam::UPstream::commsTypes
commsTypes
Types of communications.
Definition: UPstream.H:64
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::FieldMapper::addressing
virtual const labelListList & addressing() const
Definition: FieldMapper.H:94
Foam::fvPatchField::check
void check(const fvPatchField< Type > &) const
Check fvPatchField<Type> against given fvPatchField<Type>
Definition: fvPatchField.C:188
f
labelList f(nPoints)
Foam::fvPatchField::db
const objectRegistry & db() const
Return local objectRegistry.
Definition: fvPatchField.C:181
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
scalarField
volScalarField scalarField(fieldObject, mesh)
Foam::Ostream::writeKeyword
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:59
Foam::UList< Type >
dictionary.H
Foam::Field::operator*=
void operator*=(const UList< scalar > &)
Foam::fvMatrix< Type >
Foam::fvPatchField::fvPatchField
fvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:330
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:45
Foam::fvPatchField::patch
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:300
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::UList::size
label size() const
Return the number of elements in the UList.
Definition: UListI.H:299
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
Foam::FieldMapper::distributed
virtual bool distributed() const
Definition: FieldMapper.H:68
Foam::fvPatchField::autoMap
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Definition: fvPatchField.C:223