ensightSurfaceReaderTemplates.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) 2015-2022 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
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 "SpanStream.H"
29 #include "ensightPTraits.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 template<class Type>
35 (
36  const label nSkip,
37  Istream& is,
38  Type& value
39 ) const
40 {
41  skip(nSkip, is);
42 
43  is >> value;
44 }
45 
46 
47 template<class Type>
49 (
50  const label nSkip,
51  const string& buffer,
52  Type& value
53 ) const
54 {
55  ISpanStream is(buffer.data(), buffer.length());
56 
57  readFromLine(nSkip, is, value);
58 }
59 
60 
61 template<class Type>
63 (
64  const fileName& dataFile,
65  const word& fieldName
66 ) const
67 {
68  auto tfield = tmp<Field<Type>>::New(surfPtr_->nFaces(), Zero);
69  auto& field = tfield.ref();
70 
71  if (!masterOnly_ || UPstream::master(UPstream::worldComm))
72  {
73  // Use previously detected ascii/binary format
74  ensightReadFile is(dataFile, readFormat_);
75 
76  if (!is.good())
77  {
79  << "Cannot read file " << is.name()
80  << " for field " << fieldName
81  << exit(FatalError);
82  }
83 
84  // Check that data type is as expected
85  // (assuming OpenFOAM generated the data set)
86  string primitiveType;
87  is.read(primitiveType);
88 
89  DebugInfo << "primitiveType: " << primitiveType << endl;
90 
91  if
92  (
93  debug
94  && primitiveType != ensightPTraits<Type>::typeName
95  && primitiveType != pTraits<Type>::typeName
96  )
97  {
99  << "Expected <" << ensightPTraits<Type>::typeName
100  << "> values for <" << pTraits<Type>::typeName
101  << "> but found " << primitiveType << nl
102  << " This may be okay, but could indicate an error"
103  << nl << nl;
104  }
105 
106  string strValue;
107  label iValue;
108 
109  // Read header info: part index, e.g. part 1
110  is.read(strValue);
111  is.read(iValue);
112 
113  label begFace = 0;
114 
115  // Loop through different element types when reading the field values
116  for (const faceInfoTuple& facesInfo : faceTypeInfo_)
117  {
118  // [faceType, faceCount]
119  const label endFace = begFace + facesInfo.second();
120 
121  DebugInfo
122  << "Reading <" << pTraits<Type>::typeName << "> face type "
123  << ensightFaces::elemNames[facesInfo.first()]
124  << " data:" << facesInfo.second() << endl;
125 
126  if (begFace < endFace)
127  {
128  // The element type, optionally with 'undef'
129  is.read(strValue);
130 
131  if (strValue.contains("undef"))
132  {
133  // Skip undef entry
134  scalar value;
135  is.read(value);
136  }
137 
138  // Ensight fields are written component-wise
139  // (can be in different order than OpenFOAM uses)
140 
141  for (direction d = 0; d < pTraits<Type>::nComponents; ++d)
142  {
143  const direction cmpt =
145 
146  for (label facei = begFace; facei < endFace; ++facei)
147  {
148  scalar value;
149  is.read(value);
150  setComponent(field[facei], cmpt) = value;
151  }
152  }
153 
154  begFace = endFace;
155  }
156  }
157  }
158 
159  if (masterOnly_ && UPstream::parRun())
160  {
162  }
164  return tfield;
165 }
166 
167 
168 template<class Type>
170 (
171  const label timeIndex,
172  const label fieldIndex
173 ) const
174 {
175  if (fieldIndex < 0 || fieldIndex >= fieldNames_.size())
176  {
178  << "Invalid timeIndex:" << timeIndex
179  << " should be in range [0.." << fieldNames_.size() << ')' << nl
180  << "Possibly used incorrect field lookup name. Known field names: "
181  << flatOutput(fieldNames_) << nl
182  << exit(FatalError);
183  }
184 
185  const word& fieldName = fieldNames_[fieldIndex];
186  const label fileIndex = timeStartIndex_ + timeIndex*timeIncrement_;
187 
188  const fileName dataFile
189  (
190  baseDir_/replaceMask(fieldFileNames_[fieldIndex], fileIndex)
191  );
192 
193  if (debug)
194  {
195  Pout<< "Read <" << pTraits<Type>::typeName << "> field, file="
196  << dataFile << endl;
197  }
198 
199  return readField<Type>(dataFile, fieldName);
200 }
201 
202 
203 // ************************************************************************* //
static std::string::size_type length(const char *s)
Length of the character sequence (with nullptr protection)
Definition: string.H:258
rDeltaTY field()
uint8_t direction
Definition: direction.H:46
Input/output streams with (internal or external) character storage.
A class for handling file names.
Definition: fileName.H:72
virtual const fileName & name() const
The name of the stream.
Definition: IOstream.C:33
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
static const char * elemNames[nTypes]
The ensight &#39;Face&#39; element type names.
Definition: ensightFaces.H:94
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:1049
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.
::Foam::direction nComponents(const expressions::valueTypeCode) noexcept
The number of components associated with given valueTypeCode.
Definition: exprTraits.C:40
static label worldComm
Communicator for all ranks. May differ from commGlobal() if local worlds are in use.
Definition: UPstream.H:409
tmp< Field< Type > > readField(const fileName &dataFile, const word &fieldName) const
Helper function to return a field.
static void broadcast(Type &value, const label comm=UPstream::worldComm)
Broadcast content (contiguous or non-contiguous) to all communicator ranks. Does nothing in non-paral...
virtual Istream & read(token &)=0
Return next token from stream.
A class for handling words, derived from Foam::string.
Definition: word.H:63
A variant of IFstream with specialised read() for strings, integers and floats. Correctly handles bin...
void skip(const label n, Istream &is) const
Helper function to skip forward n steps in stream.
#define DebugInfo
Report an information message using Foam::Info.
const char *const typeName
int debug
Static debugging option.
bool good() const noexcept
True if next operation might succeed.
Definition: IOstream.H:281
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition: UPstream.H:1082
#define IOWarningInFunction(ios)
Report an IO warning using Foam::Warning.
const Foam::direction componentOrder[]
void readFromLine(const label nSkip, Istream &is, Type &value) const
Helper function to return Type after skipping n tokens.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
label & setComponent(label &val, const direction) noexcept
Non-const access to integer-type (has no components)
Definition: label.H:144
Similar to IStringStream but using an externally managed buffer for its input. This allows the input ...
Definition: ISpanStream.H:286
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
label timeIndex
Definition: getTimeIndex.H:24
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition: FlatOutput.H:225
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127