Loading...
Searching...
No Matches
coordSetWriterTemplates.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) 2022-2025 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "transformField.H"
29
30// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
31
32template<class Type>
34(
35 const word& fieldName,
36 const tmp<Field<Type>>& tfield
37) const
38{
39 if (verbose_)
40 {
41 Info<< "Writing field " << fieldName;
42 }
43
44 tmp<Field<Type>> tadjusted;
45
46 // Output scaling for the variable, but not for integer types
47 // which are typically ids etc.
48 if constexpr (std::is_integral_v<Type>)
49 {
50 return tfield;
51 }
52 else
53 {
54 scalar value;
55
56 // Remove *uniform* reference level
57 if
58 (
59 fieldLevel_.readIfPresent(fieldName, value, keyType::REGEX)
60 && !equal(value, 0)
61 )
62 {
63 // Could also detect brackets (...) and read accordingly
64 // or automatically scale by 1/sqrt(nComponents) instead ...
65
66 Type refLevel;
67 if constexpr (is_vectorspace_v<Type>)
68 {
69 refLevel.fill(value);
70 }
71 else
72 {
73 refLevel = value;
74 }
75
76 if (verbose_)
77 {
78 Info<< " [level " << refLevel << ']';
79 }
80
81 if (!tadjusted)
82 {
83 // Steal or clone
84 tadjusted.reset(tfield.ptr());
85 }
86
87 // Remove offset level
88 tadjusted.ref() -= refLevel;
89 }
90
91 // Apply scaling
92 if
93 (
94 fieldScale_.readIfPresent(fieldName, value, keyType::REGEX)
95 && !equal(value, 1)
96 )
97 {
98 if (verbose_)
99 {
100 Info<< " [scaling " << value << ']';
101 }
102
103 if (!tadjusted)
104 {
105 // Steal or clone
106 tadjusted.reset(tfield.ptr());
107 }
108
109 // Apply scaling
110 tadjusted.ref() *= value;
111 }
112 }
113
114 // Rotate fields (vector and non-spherical tensors)
116 {
117 if (geometryTransform_.good() && !geometryTransform_.R().is_identity())
118 {
119 if (!tadjusted)
120 {
121 // Steal or clone
122 tadjusted.reset(tfield.ptr());
123 }
124
126 (
127 tadjusted.ref(),
129 tadjusted()
130 );
131 }
132 }
134 return (tadjusted ? tadjusted : tfield);
135}
136
137
138template<class Type>
141{
142 UPtrList<const Field<Type>> fieldPtrs(1);
143 fieldPtrs.set(0, &field);
145 return fieldPtrs;
146}
147
148
149template<class Type>
152{
153 UPtrList<const Field<Type>> fieldPtrs(fieldValues.size());
154 forAll(fieldValues, i)
155 {
156 fieldPtrs.set(i, &(fieldValues[i]));
158
159 return fieldPtrs;
160}
161
162
163template<class Type>
165(
166 Ostream& os,
167 const coordSet& coords,
168 const UList<Type>& values,
169 const char* sep
170)
171{
172 forAll(coords, pointi)
173 {
174 // Output coordinate (point or scalar) with separator
175 if (coords.hasVectorAxis())
176 {
177 const vector& p = coords.vectorCoord(pointi);
178 os << p.x() << sep << p.y() << sep << p.z();
179 }
180 else
181 {
182 os << coords.scalarCoord(pointi);
183 }
184
185 // Output component values with separator
186 const auto& val = values[pointi];
187 for (direction d=0; d < pTraits<Type>::nComponents; ++d)
188 {
189 os << sep << component(val, d);
190 }
191 os << nl;
192 }
193}
194
195
196// ************************************************************************* //
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
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
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
Definition UPtrList.H:101
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie,...
Definition UPtrList.H:366
coordSystem::cartesian geometryTransform_
Local coordinate system transformation.
static UPtrList< const Field< Type > > repackageFields(const Field< Type > &field)
Repackage field into a UPtrList.
dictionary fieldLevel_
Field level to remove (on output).
static void writeTable(Ostream &os, const coordSet &coords, const UList< Type > &values, const char *sep)
Write coordinates and values.
bool verbose_
Additional output verbosity.
tmp< Field< Type > > adjustFieldTemplate(const word &fieldName, const tmp< Field< Type > > &tfield) const
dictionary fieldScale_
Field scaling (on output).
Holds list of sampling positions.
Definition coordSet.H:52
const vector & vectorCoord(const label index) const
Get point according to axis="xyz" specification.
Definition coordSet.C:181
bool hasVectorAxis() const noexcept
True if axis specification is a vector.
Definition coordSet.C:128
scalar scalarCoord(const label index) const
Get coordinate of point according to axis specification.
Definition coordSet.C:134
@ REGEX
Regular expression.
Definition keyType.H:83
A class for managing temporary objects.
Definition tmp.H:75
void reset(tmp< T > &&other) noexcept
Clear existing and transfer ownership.
Definition tmpI.H:338
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
rDeltaTY field()
OBJstream os(runTime.globalPath()/outputName)
List< T > values(const HashTable< T, Key, Hash > &tbl, const bool doSort=false)
List of values from HashTable, optionally sorted.
Definition HashOps.H:164
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
bool equal(const T &a, const T &b)
Compare two values for equality.
Definition label.H:180
refinementData transform(const tensor &, const refinementData val)
No-op rotational transform for base types.
messageStream Info
Information stream (stdout output on master, null elsewhere).
constexpr bool is_vectorspace_v
The is_vectorspace value of Type.
Definition pTraits.H:244
uint8_t direction
Definition direction.H:49
constexpr bool is_rotational_vectorspace_v
The is_rotational_vectorspace value of Type.
Definition pTraits.H:251
Vector< scalar > vector
Definition vector.H:57
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Spatial transformation functions for primitive fields.