Field.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-2021 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 "FieldMapper.H"
30#include "FieldM.H"
31#include "dictionary.H"
32#include "contiguous.H"
33#include "mapDistributeBase.H"
34
35// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
36
37template<class Type>
39(
40 const UList<Type>& mapF,
41 const labelUList& mapAddressing
42)
43:
44 List<Type>(mapAddressing.size())
45{
46 map(mapF, mapAddressing);
47}
48
49
50template<class Type>
52(
53 const tmp<Field<Type>>& tmapF,
54 const labelUList& mapAddressing
55)
56:
57 List<Type>(mapAddressing.size())
58{
59 map(tmapF, mapAddressing);
60}
61
62
63template<class Type>
65(
66 const UList<Type>& mapF,
67 const labelListList& mapAddressing,
68 const scalarListList& mapWeights
69)
70:
71 List<Type>(mapAddressing.size())
72{
73 map(mapF, mapAddressing, mapWeights);
74}
75
76
77template<class Type>
79(
80 const tmp<Field<Type>>& tmapF,
81 const labelListList& mapAddressing,
82 const scalarListList& mapWeights
83)
84:
85 List<Type>(mapAddressing.size())
86{
87 map(tmapF, mapAddressing, mapWeights);
88}
89
90
91template<class Type>
93(
94 const UList<Type>& mapF,
95 const FieldMapper& mapper,
96 const bool applyFlip
97)
98:
99 List<Type>(mapper.size())
100{
101 map(mapF, mapper, applyFlip);
102}
103
104
105template<class Type>
107(
108 const UList<Type>& mapF,
109 const FieldMapper& mapper,
110 const Type& defaultValue,
111 const bool applyFlip
112)
113:
114 List<Type>(mapper.size(), defaultValue)
115{
116 map(mapF, mapper, applyFlip);
117}
118
119
120template<class Type>
122(
123 const UList<Type>& mapF,
124 const FieldMapper& mapper,
125 const UList<Type>& defaultValues,
126 const bool applyFlip
127)
128:
129 List<Type>(defaultValues)
130{
131 map(mapF, mapper, applyFlip);
132}
133
134
135template<class Type>
137(
138 const tmp<Field<Type>>& tmapF,
139 const FieldMapper& mapper,
140 const bool applyFlip
141)
142:
143 List<Type>(mapper.size())
145 map(tmapF, mapper, applyFlip);
146}
147
148
149template<class Type>
152 const tmp<Field<Type>>& tmapF,
153 const FieldMapper& mapper,
154 const Type& defaultValue,
155 const bool applyFlip
156)
157:
158 List<Type>(mapper.size(), defaultValue)
159{
160 map(tmapF, mapper, applyFlip);
161}
162
163
164template<class Type>
167 const tmp<Field<Type>>& tmapF,
168 const FieldMapper& mapper,
169 const UList<Type>& defaultValues,
170 const bool applyFlip
171)
172:
173 List<Type>(defaultValues)
175 map(tmapF, mapper, applyFlip);
176}
177
178
179template<class Type>
181(
182 const word& keyword,
183 const dictionary& dict,
184 const label len
185)
186{
187 if (len)
188 {
189 ITstream& is = dict.lookup(keyword);
190
191 // Read first token
192 token firstToken(is);
193
194 if (firstToken.isWord("uniform"))
195 {
196 this->resize(len);
197 operator=(pTraits<Type>(is));
198 }
199 else if (firstToken.isWord("nonuniform"))
201 is >> static_cast<List<Type>&>(*this);
202 const label lenRead = this->size();
203 if (len != lenRead)
204 {
205 if (len < lenRead && allowConstructFromLargerSize)
206 {
207 #ifdef FULLDEBUG
209 << "Sizes do not match. Truncating " << lenRead
210 << " entries to " << len << endl;
211 #endif
212
213 // Truncate the data
214 this->resize(len);
215 }
216 else
217 {
219 << "size " << lenRead
220 << " is not equal to the expected length " << len
221 << exit(FatalIOError);
222 }
223 }
224 }
225 else
226 {
228 << "Expected keyword 'uniform' or 'nonuniform', found "
229 << firstToken.info() << nl
230 << exit(FatalIOError);
231 }
232 }
233}
234
235
236// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
238template<class Type>
240(
241 const UList<Type>& mapF,
242 const labelUList& mapAddressing
243)
244{
245 Field<Type>& f = *this;
246
247 if (f.size() != mapAddressing.size())
248 {
249 f.setSize(mapAddressing.size());
250 }
251
252 if (mapF.size() > 0)
253 {
254 forAll(f, i)
255 {
256 const label mapI = mapAddressing[i];
257
258 if (mapI >= 0)
259 {
260 f[i] = mapF[mapI];
261 }
262 }
263 }
264}
265
266
267template<class Type>
269(
270 const tmp<Field<Type>>& tmapF,
271 const labelUList& mapAddressing
272)
273{
274 map(tmapF(), mapAddressing);
275 tmapF.clear();
276}
277
278
279template<class Type>
281(
282 const UList<Type>& mapF,
283 const labelListList& mapAddressing,
284 const scalarListList& mapWeights
285)
286{
287 Field<Type>& f = *this;
288
289 if (f.size() != mapAddressing.size())
291 f.setSize(mapAddressing.size());
292 }
293
294 if (mapWeights.size() != mapAddressing.size())
295 {
297 << mapWeights.size() << " map size: " << mapAddressing.size()
299 }
300
301 forAll(f, i)
302 {
303 const labelList& localAddrs = mapAddressing[i];
304 const scalarList& localWeights = mapWeights[i];
305
306 f[i] = Zero;
307
308 forAll(localAddrs, j)
309 {
310 f[i] += localWeights[j]*mapF[localAddrs[j]];
311 }
312 }
314
315
316template<class Type>
318(
319 const tmp<Field<Type>>& tmapF,
320 const labelListList& mapAddressing,
321 const scalarListList& mapWeights
322)
323{
324 map(tmapF(), mapAddressing, mapWeights);
325 tmapF.clear();
326}
328
329template<class Type>
331(
332 const UList<Type>& mapF,
333 const FieldMapper& mapper,
334 const bool applyFlip
336{
337 if (mapper.distributed())
338 {
339 // Fetch remote parts of mapF
340 const mapDistributeBase& distMap = mapper.distributeMap();
341 Field<Type> newMapF(mapF);
342
343 if (applyFlip)
345 distMap.distribute(newMapF);
346 }
347 else
348 {
349 distMap.distribute(newMapF, identityOp());
351
352 if (mapper.direct() && notNull(mapper.directAddressing()))
354 map(newMapF, mapper.directAddressing());
355 }
356 else if (!mapper.direct())
357 {
358 map(newMapF, mapper.addressing(), mapper.weights());
360 else if (mapper.direct() && isNull(mapper.directAddressing()))
361 {
362 // Special case, no local mapper. Assume ordering already correct
363 // from distribution. Note: this behaviour is different compared
364 // to local mapper.
365 this->transfer(newMapF);
366 this->setSize(mapper.size());
367 }
369 else
370 {
371 if
372 (
373 mapper.direct()
374 && notNull(mapper.directAddressing())
375 && mapper.directAddressing().size()
376 )
377 {
378 map(mapF, mapper.directAddressing());
379 }
380 else if (!mapper.direct() && mapper.addressing().size())
381 {
382 map(mapF, mapper.addressing(), mapper.weights());
383 }
384 }
385}
386
387
388template<class Type>
390(
391 const tmp<Field<Type>>& tmapF,
392 const FieldMapper& mapper,
393 const bool applyFlip
394)
395{
396 map(tmapF(), mapper, applyFlip);
397 tmapF.clear();
398}
399
400
401template<class Type>
403(
404 const FieldMapper& mapper,
405 const bool applyFlip
406)
407{
408 if (mapper.distributed())
409 {
410 // Fetch remote parts of *this
411 const mapDistributeBase& distMap = mapper.distributeMap();
412 Field<Type> fCpy(*this);
413
414 if (applyFlip)
416 distMap.distribute(fCpy);
418 else
419 {
420 distMap.distribute(fCpy, identityOp());
422
423 if
425 (mapper.direct()
427 || !mapper.direct()
428 )
430 this->map(fCpy, mapper);
431 }
432 else if (mapper.direct() && isNull(mapper.directAddressing()))
434 // Special case, no local mapper. Assume ordering already correct
435 // from distribution. Note: this behaviour is different compared
436 // to local mapper.
437 this->transfer(fCpy);
438 this->setSize(mapper.size());
439 }
440 }
441 else
442 {
443 if
444 (
445 (
446 mapper.direct()
447 && notNull(mapper.directAddressing())
448 && mapper.directAddressing().size()
449 )
450 || (!mapper.direct() && mapper.addressing().size())
451 )
452 {
453 Field<Type> fCpy(*this);
454 map(fCpy, mapper);
455 }
456 else
457 {
458 this->setSize(mapper.size());
459 }
460 }
461}
462
463
464template<class Type>
466(
467 const UList<Type>& mapF,
468 const labelUList& mapAddressing
469)
470{
471 Field<Type>& f = *this;
472
473 forAll(mapF, i)
474 {
475 label mapI = mapAddressing[i];
476
477 if (mapI >= 0)
478 {
479 f[mapI] = mapF[i];
480 }
481 }
482}
483
484
485template<class Type>
487(
488 const tmp<Field<Type>>& tmapF,
489 const labelUList& mapAddressing
490)
491{
492 rmap(tmapF(), mapAddressing);
493 tmapF.clear();
494}
495
496
497template<class Type>
499(
500 const UList<Type>& mapF,
501 const labelUList& mapAddressing,
502 const UList<scalar>& mapWeights
503)
504{
505 Field<Type>& f = *this;
506
507 f = Zero;
508
509 forAll(mapF, i)
510 {
511 f[mapAddressing[i]] += mapF[i]*mapWeights[i];
512 }
513}
514
515
516template<class Type>
518(
519 const tmp<Field<Type>>& tmapF,
520 const labelUList& mapAddressing,
521 const UList<scalar>& mapWeights
522)
523{
524 rmap(tmapF(), mapAddressing, mapWeights);
525 tmapF.clear();
526}
527
528
529template<class Type>
531{
532 TFOR_ALL_F_OP_OP_F(Type, *this, =, -, Type, *this)
533}
534
535
536// A no-op except for vector specialization
537template<class Type>
539{}
540
541
542template<class Type>
545(
546 const direction d
547) const
548{
549 auto tres = tmp<Field<cmptType>>::New(this->size());
550 ::Foam::component(tres.ref(), *this, d);
551 return tres;
552}
553
554
555template<class Type>
557(
558 const direction d,
559 const UList<cmptType>& sf
560)
561{
562 TFOR_ALL_F_OP_FUNC_S_F(Type, *this, ., replace, const direction, d,
563 cmptType, sf)
564}
565
566
567template<class Type>
569(
570 const direction d,
571 const tmp<Field<cmptType>>& tsf
572)
573{
574 replace(d, tsf());
575 tsf.clear();
576}
577
578
579template<class Type>
581(
582 const direction d,
583 const cmptType& c
584)
585{
586 TFOR_ALL_F_OP_FUNC_S_S(Type, *this, ., replace, const direction, d,
587 cmptType, c)
588}
589
590
591template<class Type>
592template<class VSForm>
593VSForm Foam::Field<Type>::block(const label start) const
594{
595 VSForm vs;
596 for (direction i=0; i<VSForm::nComponents; i++)
597 {
598 vs[i] = this->operator[](start + i);
599 }
600 return vs;
601}
602
603
604template<class Type>
606{
607 auto tres = tmp<Field<Type>>::New(this->size());
608 ::Foam::T(tres.ref(), *this);
609 return tres;
610}
611
612
613template<class Type>
614void Foam::Field<Type>::writeEntry(const word& keyword, Ostream& os) const
615{
616 if (keyword.size())
617 {
618 os.writeKeyword(keyword);
619 }
620
621 // The contents are 'uniform' if the list is non-empty
622 // and all entries have identical values.
623
625 {
626 os << word("uniform") << token::SPACE << this->first();
627 }
628 else
629 {
630 os << word("nonuniform") << token::SPACE;
632 }
633
635}
636
637
638// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
639
640template<class Type>
642{
643 if (this == &rhs)
644 {
645 return; // Self-assignment is a no-op
646 }
647
649}
650
651
652template<class Type>
654{
655 if (this == &(rhs()))
656 {
657 return; // Self-assignment is a no-op
658 }
659
661}
662
663
664template<class Type>
665template<class Form, class Cmpt, Foam::direction nCmpt>
667{
668 TFOR_ALL_F_OP_S(Type, *this, =, VSType, vs)
669}
670
671
672#define COMPUTED_ASSIGNMENT(TYPE, op) \
673 \
674template<class Type> \
675void Foam::Field<Type>::operator op(const UList<TYPE>& f) \
676{ \
677 TFOR_ALL_F_OP_F(Type, *this, op, TYPE, f) \
678} \
679 \
680template<class Type> \
681void Foam::Field<Type>::operator op(const tmp<Field<TYPE>>& tf) \
682{ \
683 operator op(tf()); \
684 tf.clear(); \
685} \
686 \
687template<class Type> \
688void Foam::Field<Type>::operator op(const TYPE& t) \
689{ \
690 TFOR_ALL_F_OP_S(Type, *this, op, TYPE, t) \
691}
692
697
698#undef COMPUTED_ASSIGNMENT
699
700
701// * * * * * * * * * * * * * * * Ostream Operator * * * * * * * * * * * * * //
702
703template<class Type>
705{
706 os << static_cast<const List<Type>&>(f);
707 return os;
708}
709
710
711template<class Type>
713{
714 os << tf();
715 tf.clear();
716 return os;
717}
718
719
720// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
721
722#include "FieldFunctions.C"
723
724// ************************************************************************* //
#define COMPUTED_ASSIGNMENT(TYPE, op)
High performance macro functions for Field<Type> algebra. These expand using either array element acc...
#define TFOR_ALL_F_OP_OP_F(typeF1, f1, OP1, OP2, typeF2, f2)
Definition: FieldM.H:341
#define TFOR_ALL_F_OP_S(typeF, f, OP, typeS, s)
Definition: FieldM.H:359
#define TFOR_ALL_F_OP_FUNC_S_S(typeF1, f1, OP, FUNC, typeS1, s1, typeS2, s2)
Definition: FieldM.H:237
#define TFOR_ALL_F_OP_FUNC_S_F(typeF1, f1, OP, FUNC, typeS, s, typeF2, f2)
Definition: FieldM.H:219
Abstract base class to hold the Field mapping addressing and weights.
Definition: FieldMapper.H:50
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:82
virtual const scalarListList & weights() const
Return the interpolation weights.
Definition: FieldMapper.H:112
virtual label size() const =0
The size of the mapper.
virtual bool distributed() const
Does the mapper have remote contributions?
Definition: FieldMapper.H:72
virtual const labelListList & addressing() const
Return the interpolation addressing.
Definition: FieldMapper.H:102
virtual const mapDistributeBase & distributeMap() const
Return the distribution map.
Definition: FieldMapper.H:92
Generic templated field type.
Definition: Field.H:82
tmp< Field< Type > > T() const
Return the field transpose (only defined for second rank tensors)
Definition: Field.C:605
void operator=(const Field< Type > &)
Copy assignment.
Definition: Field.C:641
void autoMap(const FieldMapper &map, const bool applyFlip=true)
Map from self.
Definition: Field.C:403
void writeEntry(const word &keyword, Ostream &os) const
Write the field as a dictionary entry.
Definition: Field.C:614
void replace(const direction, const UList< cmptType > &)
Replace a component field of the field.
Definition: Field.C:557
constexpr Field() noexcept
Default construct.
Definition: FieldI.H:40
void map(const UList< Type > &mapF, const labelListList &mapAddressing, const scalarListList &weights)
Interpolative map from the given field.
Definition: Field.C:281
void negate()
Inplace negate this field (negative).
Definition: Field.C:530
void map(const UList< Type > &mapF, const labelUList &mapAddressing)
1 to 1 map from the given field
Definition: Field.C:240
void normalise()
Definition: Field.C:538
void rmap(const UList< Type > &mapF, const labelUList &mapAddressing)
1 to 1 reverse-map from the given field
Definition: Field.C:466
An input stream of tokens.
Definition: ITstream.H:56
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:77
void setSize(const label n)
Alias for resize()
Definition: List.H:218
void operator=(const UList< T > &a)
Assignment to UList operator. Takes linear time.
Definition: List.C:480
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
virtual Ostream & writeKeyword(const keyType &kw)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:57
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:94
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
Templated vector space.
Definition: VectorSpace.H:79
Creates a single block of cells from point coordinates, numbers of cells in each direction and an exp...
Definition: block.H:61
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
ITstream & lookup(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionary.C:386
Particle-size distribution model wherein random samples are drawn from the doubly-truncated uniform p...
Definition: uniform.H:164
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:66
Class containing processor-to-processor mapping information.
static void distribute(const Pstream::commsTypes commsType, const List< labelPair > &schedule, const label constructSize, const labelListList &subMap, const bool subHasFlip, const labelListList &constructMap, const bool constructHasFlip, List< T > &field, const NegateOp &negOp, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Distribute data with specified negate operator (for flips).
A traits class, which is primarily used for primitives.
Definition: pTraits.H:59
A class for managing temporary objects.
Definition: tmp.H:65
A token holds an item read from Istream.
Definition: token.H:69
@ END_STATEMENT
End entry [isseparator].
Definition: token.H:154
@ SPACE
Space [isspace].
Definition: token.H:125
A class for handling words, derived from Foam::string.
Definition: word.H:68
patchWriters resize(patchIds.size())
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:473
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
#define IOWarningInFunction(ios)
Report an IO warning using Foam::Warning.
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
bool isNull(const T *ptr)
True if ptr is a pointer (of type T) to the nullObject.
Definition: nullObject.H:192
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
Definition: boundaryPatch.C:83
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
errorManip< error > abort(error &err)
Definition: errorManip.H:144
uint8_t direction
Definition: direction.H:56
bool notNull(const T *ptr)
True if ptr is not a pointer (of type T) to the nullObject.
Definition: nullObject.H:207
IOerror FatalIOError
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
error FatalError
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
points setSize(newPointi)
labelList f(nPoints)
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
A template class to specify that a data type can be considered as being contiguous in memory.
Definition: contiguous.H:78