fvPatchField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2022 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 {}
47 
48 
49 template<class Type>
51 (
52  const fvPatch& p,
54  const Field<Type>& f
55 )
56 :
57  Field<Type>(f),
58  patch_(p),
59  internalField_(iF),
60  updated_(false),
61  manipulatedMatrix_(false)
62 {}
63 
64 
65 template<class Type>
67 (
68  const fvPatch& p,
70  const dictionary& dict,
71  const bool valueRequired
72 )
73 :
74  Field<Type>(p.size()),
75  patch_(p),
76  internalField_(iF),
77  updated_(false),
78  manipulatedMatrix_(false)
79 {
80  if (valueRequired)
81  {
82  if (dict.found("value"))
83  {
85  (
86  Field<Type>("value", dict, p.size())
87  );
88  }
89  else
90  {
92  (
93  dict
94  ) << "Essential entry 'value' missing"
95  << exit(FatalIOError);
96  }
97  }
98 }
99 
100 
101 template<class Type>
103 (
104  const fvPatchField<Type>& ptf,
105  const fvPatch& p,
107  const fvPatchFieldMapper& mapper,
108  const bool mappingRequired
109 )
110 :
111  Field<Type>(p.size()),
112  patch_(p),
113  internalField_(iF),
114  updated_(false),
115  manipulatedMatrix_(false)
116 {
117  if (mappingRequired)
118  {
119  // For unmapped faces set to internal field value (zero-gradient)
120  if (notNull(iF) && mapper.hasUnmapped())
121  {
122  fvPatchField<Type>::operator=(this->patchInternalField());
123  }
124  mapper(*this, ptf);
125  }
126 }
127 
128 
129 template<class Type>
131 (
132  const fvPatchField<Type>& ptf,
134 )
135 :
136  Field<Type>(ptf),
137  patch_(ptf.patch_),
138  internalField_(iF),
139  updated_(false),
140  manipulatedMatrix_(false)
141 {}
142 
143 
144 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
145 
146 template<class Type>
148 {
149  return patch_.boundaryMesh().mesh();
150 }
151 
152 
153 template<class Type>
155 {
156  if (&patch_ != &(ptf.patch_))
157  {
159  << "different patches for fvPatchField<Type>s"
160  << abort(FatalError);
161  }
162 }
163 
164 
165 template<class Type>
167 {
168  return patch_.deltaCoeffs()*(*this - patchInternalField());
169 }
170 
171 
172 template<class Type>
175 {
176  return patch_.patchInternalField(internalField_);
177 }
178 
179 
180 template<class Type>
182 {
183  patch_.patchInternalField(internalField_, pif);
184 }
185 
186 
187 template<class Type>
189 (
190  const fvPatchFieldMapper& mapper
191 )
192 {
193  mapper(*this, *this);
194 }
195 
196 
197 template<class Type>
199 (
200  const fvPatchField<Type>& ptf,
201  const labelList& addr
202 )
203 {
204  Field<Type>::rmap(ptf, addr);
205 }
206 
207 
208 template<class Type>
210 {
211  Field<Type>::reset(ptf);
212 }
213 
214 
215 template<class Type>
217 {
218  updated_ = true;
219 }
220 
221 
222 template<class Type>
224 {
225  if (!updated_)
226  {
227  updateCoeffs();
228  }
229 
230  updated_ = false;
231  manipulatedMatrix_ = false;
232 }
233 
234 
235 template<class Type>
237 {
238  manipulatedMatrix_ = true;
239 }
240 
241 
242 template<class Type>
244 {
245  writeEntry(os, "type", type());
246 
247  if (overridesConstraint())
248  {
249  writeEntry(os, "patchType", patch().type());
250  }
251 }
252 
253 
254 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
255 
256 template<class Type>
257 void Foam::fvPatchField<Type>::operator=
258 (
259  const UList<Type>& ul
260 )
261 {
263 }
264 
265 
266 template<class Type>
267 void Foam::fvPatchField<Type>::operator=
268 (
269  const fvPatchField<Type>& ptf
270 )
271 {
272  check(ptf);
274 }
275 
276 
277 template<class Type>
278 void Foam::fvPatchField<Type>::operator+=
279 (
280  const fvPatchField<Type>& ptf
281 )
282 {
283  check(ptf);
285 }
286 
287 
288 template<class Type>
289 void Foam::fvPatchField<Type>::operator-=
290 (
291  const fvPatchField<Type>& ptf
292 )
293 {
294  check(ptf);
296 }
297 
298 
299 template<class Type>
300 void Foam::fvPatchField<Type>::operator*=
301 (
302  const fvPatchField<scalar>& ptf
303 )
304 {
305  if (&patch_ != &ptf.patch())
306  {
308  << "incompatible patches for patch fields"
309  << abort(FatalError);
310  }
311 
313 }
314 
315 
316 template<class Type>
317 void Foam::fvPatchField<Type>::operator/=
318 (
319  const fvPatchField<scalar>& ptf
320 )
321 {
322  if (&patch_ != &ptf.patch())
323  {
325  << abort(FatalError);
326  }
327 
329 }
330 
331 
332 template<class Type>
333 void Foam::fvPatchField<Type>::operator+=
334 (
335  const Field<Type>& tf
336 )
337 {
339 }
340 
341 
342 template<class Type>
343 void Foam::fvPatchField<Type>::operator-=
344 (
345  const Field<Type>& tf
346 )
347 {
349 }
350 
351 
352 template<class Type>
353 void Foam::fvPatchField<Type>::operator*=
354 (
355  const scalarField& tf
356 )
357 {
359 }
360 
361 
362 template<class Type>
363 void Foam::fvPatchField<Type>::operator/=
364 (
365  const scalarField& tf
366 )
367 {
369 }
370 
371 
372 template<class Type>
373 void Foam::fvPatchField<Type>::operator=
374 (
375  const Type& t
376 )
377 {
379 }
380 
381 
382 template<class Type>
383 void Foam::fvPatchField<Type>::operator+=
384 (
385  const Type& t
386 )
387 {
389 }
390 
391 
392 template<class Type>
393 void Foam::fvPatchField<Type>::operator-=
394 (
395  const Type& t
396 )
397 {
399 }
400 
401 
402 template<class Type>
403 void Foam::fvPatchField<Type>::operator*=
404 (
405  const scalar s
406 )
407 {
409 }
410 
411 
412 template<class Type>
413 void Foam::fvPatchField<Type>::operator/=
414 (
415  const scalar s
416 )
417 {
419 }
420 
421 
422 template<class Type>
423 void Foam::fvPatchField<Type>::operator==
424 (
425  const fvPatchField<Type>& ptf
426 )
427 {
429 }
430 
431 
432 template<class Type>
433 void Foam::fvPatchField<Type>::operator==
434 (
435  const Field<Type>& tf
436 )
437 {
439 }
440 
441 
442 template<class Type>
443 void Foam::fvPatchField<Type>::operator==
444 (
445  const Type& t
446 )
447 {
449 }
450 
451 
452 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
453 
454 template<class Type>
455 Foam::Ostream& Foam::operator<<(Ostream& os, const fvPatchField<Type>& ptf)
456 {
457  ptf.write(os);
458 
459  os.check("Ostream& operator<<(Ostream&, const fvPatchField<Type>&");
460 
461  return os;
462 }
463 
464 
465 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
466 
467 #include "fvPatchFieldNew.C"
468 
469 // ************************************************************************* //
dictionary dict
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:663
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Definition: fvPatchField.C:189
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Definition: fvPatchField.C:166
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:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
commsTypes
Types of communications.
Definition: UPstream.H:64
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:243
const tensorField & tf
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))
virtual void reset(const fvPatchField< Type > &)
Reset the fvPatchField to the given fvPatchField.
Definition: fvPatchField.C:209
Pre-declare SubField and related Field type.
Definition: Field.H:56
Foam::fvPatchFieldMapper.
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Definition: fvPatchField.C:199
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:72
virtual bool hasUnmapped() const =0
Are there unmapped values? I.e. do all size() elements get.
virtual void manipulateMatrix(fvMatrix< Type > &matrix)
Manipulate matrix.
Definition: fvPatchField.C:236
void check(const fvPatchField< Type > &) const
Check fvPatchField<Type> against given fvPatchField<Type>
Definition: fvPatchField.C:154
errorManip< error > abort(error &err)
Definition: errorManip.H:131
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:351
virtual label size() const
Return size.
Definition: fvPatch.H:157
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:174
labelList f(nPoints)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
bool notNull(const T &t)
Return true if t is not a reference to the nullObject of type T.
Definition: nullObjectI.H:52
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
Evaluate the patch field, sets Updated to false.
Definition: fvPatchField.C:223
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:318
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:216
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const objectRegistry & db() const
Return local objectRegistry.
Definition: fvPatchField.C:147
A class for managing temporary objects.
Definition: PtrList.H:53
Registry of regIOobjects.
fvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
Definition: fvPatchField.C:36
IOerror FatalIOError