Loading...
Searching...
No Matches
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-2025 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
36
37template<class Type>
39(
40 const dictionary& dict,
42)
43{
44 if (!IOobjectOption::isAnyRead(readOpt)) return false;
45 const auto& p = fvPatchFieldBase::patch();
46
47
48 const auto* eptr = dict.findEntry("value", keyType::LITERAL);
49
50 if (eptr)
51 {
52 Field<Type>::assign(*eptr, p.size());
53 return true;
54 }
55
57 {
59 << "Required entry 'value' : missing for patch " << p.name()
60 << " in dictionary " << dict.relativeName() << nl
62 }
63
64 return false;
65}
66
67
68template<class Type>
70{
71 const auto& p = fvPatchFieldBase::patch();
72 this->resize_nocopy(p.size()); // In general this is a no-op
73 p.patchInternalField(internalField_, *this);
74}
75
76
77// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
78
79template<class Type>
81(
82 const fvPatch& p,
84)
85:
87 Field<Type>(p.size()),
88 internalField_(iF)
89{}
90
91
92template<class Type>
94(
95 const fvPatch& p,
97 const word& patchType
98)
99:
101 Field<Type>(p.size()),
102 internalField_(iF)
103{}
104
105
106template<class Type>
108(
109 const fvPatch& p,
111 const Type& value
112)
113:
115 Field<Type>(p.size(), value),
116 internalField_(iF)
117{}
118
119
120template<class Type>
122(
123 const fvPatch& p,
125 const Field<Type>& pfld
126)
127:
129 Field<Type>(pfld),
130 internalField_(iF)
131{}
132
133
134template<class Type>
136(
137 const fvPatch& p,
139 Field<Type>&& pfld
140)
141:
143 Field<Type>(std::move(pfld)),
144 internalField_(iF)
145{}
146
147
148template<class Type>
150(
151 const fvPatch& p,
153 const dictionary& dict,
154 IOobjectOption::readOption requireValue
155)
156:
158 Field<Type>(p.size()),
159 internalField_(iF)
160{
161 readValueEntry(dict, requireValue);
162}
163
164
165template<class Type>
167(
168 const fvPatchField<Type>& ptf,
169 const fvPatch& p,
171 const fvPatchFieldMapper& mapper
172)
173:
174 fvPatchFieldBase(ptf, p),
175 Field<Type>(p.size()),
176 internalField_(iF)
177{
178 // For unmapped faces set to internal field value (zero-gradient)
179 if (notNull(iF) && mapper.hasUnmapped())
180 {
181 this->extrapolateInternal();
182 }
183 this->map(ptf, mapper);
184}
185
186
187template<class Type>
189(
190 const fvPatchField<Type>& pfld,
191 const fvPatch& p,
193 const Type& value
194)
195:
197 Field<Type>(p.size(), value),
198 internalField_(iF)
199{}
200
201
202template<class Type>
204(
205 const fvPatchField<Type>& ptf,
207)
208:
209 fvPatchFieldBase(ptf),
210 Field<Type>(ptf),
211 internalField_(iF)
212{}
213
214
215// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
216
217template<class Type>
219{
221}
222
223
224template<class Type>
226{
227 // Get patch internal field, store temporarily in result
228 this->patchInternalField(result);
229 const auto& pif = result;
230
231 const Field<Type>& pfld = *this;
232 const auto& dc = patch().deltaCoeffs();
233
234 const label len = result.size();
235
236 for (label i = 0; i < len; ++i)
238 result[i] = dc[i]*(pfld[i] - pif[i]);
239 }
240}
241
242
243template<class Type>
245{
246 auto tfld = tmp<Field<Type>>::New(this->size());
247 this->snGrad(static_cast<UList<Type>&>(tfld.ref()));
248 return tfld;
249}
250
251
252template<class Type>
255{
256 return patch().patchInternalField(internalField_);
257}
258
259
260template<class Type>
262{
263 patch().patchInternalField(internalField_, pfld);
264}
265
266
267template<class Type>
269(
270 const fvPatchFieldMapper& mapper
271)
272{
273 Field<Type>& f = *this;
274
275 if (!this->size() && !mapper.distributed())
276 {
277 f.resize_nocopy(mapper.size());
278 if (f.size())
279 {
280 f = this->patchInternalField();
281 }
282 }
283 else
284 {
285 // Map all faces provided with mapping data
286 Field<Type>::autoMap(mapper);
287
288
289 // For unmapped faces set to internal field value (zero-gradient)
290 if (mapper.hasUnmapped())
291 {
292 Field<Type> pif(this->patchInternalField());
293
294 if
295 (
296 mapper.direct()
297 && notNull(mapper.directAddressing())
298 && mapper.directAddressing().size()
299 )
300 {
301 const labelUList& mapAddressing = mapper.directAddressing();
302
303 forAll(mapAddressing, i)
304 {
305 if (mapAddressing[i] < 0)
306 {
307 f[i] = pif[i];
308 }
309 }
310 }
311 else if (!mapper.direct() && mapper.addressing().size())
312 {
313 const labelListList& mapAddressing = mapper.addressing();
314
315 forAll(mapAddressing, i)
316 {
317 const labelList& localAddrs = mapAddressing[i];
318
319 if (!localAddrs.size())
320 {
321 f[i] = pif[i];
322 }
323 }
325 }
326 }
327}
328
329
330template<class Type>
332(
333 const fvPatchField<Type>& ptf,
334 const labelList& addr
336{
337 Field<Type>::rmap(ptf, addr);
338}
339
340
341template<class Type>
343{
345}
346
347
348template<class Type>
350{
351 // Default behaviour ignores the weights
352 if (!updated())
353 {
354 updateCoeffs();
357 }
358}
359
360
361template<class Type>
363{
364 if (!updated())
365 {
366 updateCoeffs();
367 }
372
373
374template<class Type>
383(
384 fvMatrix<Type>& matrix,
385 const scalarField& weights
387{
389}
390
391
392template<class Type>
394(
395 fvMatrix<Type>& matrix,
396 const label iMatrix,
397 const direction cmp
399{
401}
402
403
404template<class Type>
406{
407 os.writeEntry("type", type());
408
409 if (!patchType().empty())
410 {
411 os.writeEntry("patchType", patchType());
413 if (useImplicit())
414 {
415 os.writeEntry("useImplicit", "true");
417}
418
419
420// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
421
422template<class Type>
424(
425 const UList<Type>& ul
427{
430
431
432template<class Type>
434(
435 const fvPatchField<Type>& ptf
436)
440}
441
442
443template<class Type>
445(
446 const fvPatchField<Type>& ptf
447)
451}
452
453
454template<class Type>
456(
457 const fvPatchField<Type>& ptf
458)
462}
463
464
465template<class Type>
467(
468 const fvPatchField<scalar>& ptf
469)
473}
474
475
476template<class Type>
478(
479 const fvPatchField<scalar>& ptf
480)
484}
485
486
487template<class Type>
489(
490 const Field<Type>& tf
492{
494}
495
496
497template<class Type>
499(
500 const Field<Type>& tf
501)
502{
504}
505
506
507template<class Type>
509(
510 const scalarField& tf
511)
512{
514}
515
516
517template<class Type>
518void Foam::fvPatchField<Type>::operator/=
519(
520 const scalarField& tf
522{
524}
525
526
527template<class Type>
529(
530 const Type& t
532{
534}
535
536
537template<class Type>
539(
540 const Type& t
542{
544}
545
546
547template<class Type>
549(
550 const Type& t
552{
554}
555
556
557template<class Type>
559(
560 const scalar s
562{
564}
565
566
567template<class Type>
569(
570 const scalar s
572{
574}
575
576
577template<class Type>
579(
580 const fvPatchField<Type>& ptf
582{
584}
585
586
587template<class Type>
589(
590 const Field<Type>& tf
592{
594}
595
596
597template<class Type>
599(
600 const Type& t
601)
602{
604}
605
606
607// * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
608
609template<class Type>
610Foam::Ostream& Foam::operator<<(Ostream& os, const fvPatchField<Type>& ptf)
611{
612 ptf.write(os);
613
614 os.check(FUNCTION_NAME);
615
616 return os;
617}
618
619
620// ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
virtual bool direct() const =0
Is it a direct (non-interpolating) mapper?
virtual const labelUList & directAddressing() const
Return the direct addressing values.
Definition FieldMapper.H:91
virtual label size() const =0
The size of the mapper.
virtual bool distributed() const
Does the mapper have remote contributions?
Definition FieldMapper.H:76
virtual const labelListList & addressing() const
Return the interpolation addressing.
virtual bool hasUnmapped() const =0
Any unmapped values?
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
void operator=(const Field< Type > &)
Copy assignment.
Definition Field.C:781
void autoMap(const FieldMapper &map, const bool applyFlip=true)
Map from self.
Definition Field.C:465
void operator+=(const UList< Type > &)
Definition Field.C:833
void operator-=(const UList< Type > &)
Definition Field.C:834
constexpr Field() noexcept
Default construct.
Definition FieldI.H:24
void operator*=(const UList< scalar > &)
Definition Field.C:835
void assign(const entry &e, const label len)
Assign from a primitive dictionary entry with the following behaviour:
Definition Field.C:206
void operator/=(const UList< scalar > &)
Definition Field.C:836
void map(const UList< Type > &mapF, const labelUList &mapAddressing)
1 to 1 map from the given field
Definition Field.C:302
void rmap(const UList< Type > &mapF, const labelUList &mapAddressing)
1 to 1 reverse-map from the given field
Definition Field.C:528
static bool isReadRequired(readOption opt) noexcept
True if (MUST_READ | READ_MODIFIED) bits are set.
readOption
Enumeration defining read preferences.
static bool isAnyRead(readOption opt) noexcept
True if any reading may be required (ie, != NO_READ).
void resize_nocopy(const label len)
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition UList.H:89
bool empty() const noexcept
Definition UList.H:701
void size(const label n)
Definition UList.H:118
commsTypes
Communications types.
Definition UPstream.H:81
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition fvMatrix.H:118
Template invariant parts for fvPatchField.
fvPatchFieldBase(const fvPatch &p)
Construct from patch.
void checkPatch(const fvPatchFieldBase &rhs) const
Check that patches are identical.
const fvPatch & patch() const noexcept
Return the patch.
void setManipulated(bool state) noexcept
Set matrix manipulated state.
bool useImplicit() const noexcept
Use implicit formulation for coupled patches only.
void setUpdated(bool state) noexcept
Set updated state.
const word & patchType() const noexcept
The optional patch type.
bool updated() const noexcept
True if the boundary condition has already been updated.
A FieldMapper for finite-volume patch fields.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
virtual void write(Ostream &) const
Write.
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch.
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void manipulateMatrix(fvMatrix< Type > &matrix)
Manipulate matrix.
virtual void updateWeightedCoeffs(const scalarField &weights)
Update the coefficients associated with the patch field.
virtual void rmap(const fvPatchField< Type > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
bool readValueEntry(const dictionary &dict, IOobjectOption::readOption readOpt=IOobjectOption::LAZY_READ)
Read the "value" entry into *this.
void extrapolateInternal()
Assign the patch field from the internal field.
void check(const fvPatchField< Type > &) const
Check against given patch field.
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::buffered)
Evaluate the patch field, sets updated() to false.
fvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
@ LITERAL
String literal.
Definition keyType.H:82
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition error.H:629
OBJstream os(runTime.globalPath()/outputName)
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))
#define FUNCTION_NAME
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces).
uint8_t direction
Definition direction.H:49
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
void rhs(fvMatrix< typename Expr::value_type > &m, const Expr &expression)
bool notNull(const T *ptr) noexcept
True if ptr is not a pointer (of type T) to the nullObject.
Definition nullObject.H:267
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
UList< label > labelUList
A UList of labels.
Definition UList.H:75
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
labelList f(nPoints)
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299