Nektar++
Public Member Functions | Static Public Member Functions | Static Public Attributes | Protected Member Functions | List of all members
Nektar::FieldUtils::InputNek5000 Class Reference

#include <InputNek5000.h>

Inheritance diagram for Nektar::FieldUtils::InputNek5000:
[legend]

Public Member Functions

 InputNek5000 (FieldSharedPtr f)
 Set up InputNek5000 object. More...
 
 ~InputNek5000 () override
 
- Public Member Functions inherited from Nektar::FieldUtils::InputModule
 InputModule (FieldSharedPtr p_m)
 
void PrintSummary ()
 Print a brief summary of information. More...
 
- Public Member Functions inherited from Nektar::FieldUtils::Module
FIELD_UTILS_EXPORT Module (FieldSharedPtr p_f)
 
virtual ~Module ()=default
 
void Process (po::variables_map &vm)
 
std::string GetModuleName ()
 
std::string GetModuleDescription ()
 
const ConfigOptionGetConfigOption (const std::string &key) const
 
ModulePriority GetModulePriority ()
 
std::vector< ModuleKeyGetModulePrerequisites ()
 
FIELD_UTILS_EXPORT void RegisterConfig (std::string key, std::string value="")
 Register a configuration option with a module. More...
 
FIELD_UTILS_EXPORT void PrintConfig ()
 Print out all configuration options for a module. More...
 
FIELD_UTILS_EXPORT void SetDefaults ()
 Sets default configuration options for those which have not been set. More...
 
FIELD_UTILS_EXPORT void AddFile (std::string fileType, std::string fileName)
 
FIELD_UTILS_EXPORT void EvaluateTriFieldAtEquiSpacedPts (LocalRegions::ExpansionSharedPtr &exp, const Array< OneD, const NekDouble > &infield, Array< OneD, NekDouble > &outfield)
 

Static Public Member Functions

static ModuleSharedPtr create (FieldSharedPtr f)
 Creates an instance of this class. More...
 
- Static Public Member Functions inherited from Nektar::FieldUtils::InputModule
static FIELD_UTILS_EXPORT std::string GuessFormat (std::string fileName)
 Tries to guess the format of the input file. More...
 

Static Public Attributes

static ModuleKey m_className []
 ModuleKey for class. More...
 

Protected Member Functions

void v_Process (po::variables_map &vm) override
 Process Nek5000 input file. More...
 
std::string v_GetModuleName () override
 
std::string v_GetModuleDescription () override
 
ModulePriority v_GetModulePriority () override
 
- Protected Member Functions inherited from Nektar::FieldUtils::Module
 Module ()
 
virtual void v_Process (po::variables_map &vm)
 
virtual std::string v_GetModuleName ()
 
virtual std::string v_GetModuleDescription ()
 
virtual ModulePriority v_GetModulePriority ()
 
virtual std::vector< ModuleKeyv_GetModulePrerequisites ()
 

Additional Inherited Members

- Public Attributes inherited from Nektar::FieldUtils::Module
FieldSharedPtr m_f
 Field object. More...
 
- Protected Attributes inherited from Nektar::FieldUtils::Module
std::map< std::string, ConfigOptionm_config
 List of configuration values. More...
 
std::set< std::string > m_allowedFiles
 List of allowed file formats. More...
 

Detailed Description

Converter for Fld files.

Definition at line 46 of file InputNek5000.h.

Constructor & Destructor Documentation

◆ InputNek5000()

Nektar::FieldUtils::InputNek5000::InputNek5000 ( FieldSharedPtr  f)

Set up InputNek5000 object.

Definition at line 57 of file InputNek5000.cpp.

57 : InputModule(f)
58{
59 m_allowedFiles.insert("fld5000");
60}
InputModule(FieldSharedPtr p_m)
Definition: Module.cpp:61
std::set< std::string > m_allowedFiles
List of allowed file formats.
Definition: Module.h:274

References Nektar::FieldUtils::Module::m_allowedFiles.

◆ ~InputNek5000()

Nektar::FieldUtils::InputNek5000::~InputNek5000 ( )
override

Definition at line 65 of file InputNek5000.cpp.

66{
67}

Member Function Documentation

◆ create()

static ModuleSharedPtr Nektar::FieldUtils::InputNek5000::create ( FieldSharedPtr  f)
inlinestatic

Creates an instance of this class.

Definition at line 53 of file InputNek5000.h.

54 {
56 }
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr().

◆ v_GetModuleDescription()

std::string Nektar::FieldUtils::InputNek5000::v_GetModuleDescription ( )
inlineoverrideprotectedvirtual

Reimplemented from Nektar::FieldUtils::Module.

Definition at line 68 of file InputNek5000.h.

69 {
70 return "Processing Nek5000 field file";
71 }

◆ v_GetModuleName()

std::string Nektar::FieldUtils::InputNek5000::v_GetModuleName ( )
inlineoverrideprotectedvirtual

Reimplemented from Nektar::FieldUtils::Module.

Definition at line 63 of file InputNek5000.h.

64 {
65 return "InputNek5000";
66 }

◆ v_GetModulePriority()

ModulePriority Nektar::FieldUtils::InputNek5000::v_GetModulePriority ( )
inlineoverrideprotectedvirtual

Reimplemented from Nektar::FieldUtils::Module.

Definition at line 73 of file InputNek5000.h.

74 {
75 return eCreateFieldData;
76 }

References Nektar::FieldUtils::eCreateFieldData.

◆ v_Process()

void Nektar::FieldUtils::InputNek5000::v_Process ( po::variables_map &  vm)
overrideprotectedvirtual

Process Nek5000 input file.

This routine reads a binary-format Nek5000 field file, loads the data into memory and populates the field definitions to match the data format. Nek5000 is a classic nodal-Lagrangian spectral element code at a single polynomial order, meaning that the field data are set up according to this structure.

This module is adapted from the VisIt visualisation software, which supports a number of Nek5000 inputs.

Reimplemented from Nektar::FieldUtils::Module.

Definition at line 80 of file InputNek5000.cpp.

81{
82 ifstream file(m_config["infile"].as<string>().c_str(), ios::binary);
83
84 // Header: 132 bytes for binary.
85 vector<char> data(132);
86 file.read(&data[0], 132);
87
88 // Check header: should be the four characters #std
89 string check(&data[0], 4);
90 string header(&data[4], 128);
91
92 ASSERTL0(check == "#std", "Unable to read file");
93
94 // Determine whether we need to byte-swap data: 4-byte float at byte 80
95 // should be 6.54321.
96 bool byteSwap = false;
97 float test;
98
99 file.read((char *)(&test), 4);
100 if (test > 6.5 && test < 6.6)
101 {
102 byteSwap = false;
103 }
104 else
105 {
107 ASSERTL0(test > 6.5 && test < 6.6,
108 "Unable to determine endian-ness of input file");
109 byteSwap = true;
110 }
111
112 stringstream ss;
113 ss.str(header);
114
115 int nBytes, nBlocksXYZ[3], nBlocks, nTotBlocks, dir, nDirs, nCycle, nDim;
116 NekDouble time;
117
118 // Read header information (this is written as ASCII)
119 string remain;
120 ss >> nBytes >> nBlocksXYZ[0] >> nBlocksXYZ[1] >> nBlocksXYZ[2] >>
121 nBlocks >> nTotBlocks >> time >> nCycle >> dir >> nDirs >> remain;
122 boost::trim(remain);
123
124 nDim = nBlocksXYZ[2] == 1 ? 2 : 3;
125
126 // Print some basic information for input if in verbose mode.
127 if (m_f->m_verbose)
128 {
129 cout << "Found header information:" << endl;
130 cout << " -- " << (byteSwap ? "" : "do not ") << "need to swap endian"
131 << endl;
132 cout << " -- Data byte size : " << nBytes << endl;
133 cout << " -- Number of xyz blocks : " << nBlocksXYZ[0] << "x"
134 << nBlocksXYZ[1] << "x" << nBlocksXYZ[2] << endl;
135 cout << " -- Blocks in file/total : " << nBlocks << "/" << nTotBlocks
136 << endl;
137 cout << " -- Simulation time : " << time << endl;
138 cout << " -- Number of cycles : " << nCycle << endl;
139 cout << " -- Number of dirs : " << dir << "/" << nDirs << endl;
140 cout << " -- Remaining header : " << remain << endl;
141 }
142
143 // Major limitation: we don't read out of multiple directories
144 ASSERTL0(nDirs == 1, "Number of directories must be one");
145
146 // We also don't read partial files.
147 ASSERTL0(nBlocks == nTotBlocks, "Partial field output not supported");
148
149 // We don't support non-double data
150 ASSERTL0(nBytes == 8, "Data file must contain double-precision data");
151
152 // Set up a field definition
155 fielddef->m_shapeType = LibUtilities::eHexahedron;
156 fielddef->m_numHomogeneousDir = 0;
157 fielddef->m_homoStrips = false;
158 fielddef->m_pointsDef = false;
159 fielddef->m_uniOrder = true;
160 fielddef->m_numPointsDef = false;
161
162 for (int i = 0; i < nDim; ++i)
163 {
164 fielddef->m_basis.push_back(LibUtilities::eGLL_Lagrange);
165 fielddef->m_numModes.push_back(nBlocksXYZ[i]);
166 }
167
168 // Read element IDs
169 NekUInt32 maxID = 0, minID = numeric_limits<NekUInt32>::max();
170 for (NekUInt32 i = 0; i < nBlocks; ++i)
171 {
172 NekUInt32 blockNum;
173 file.read((char *)&blockNum, 4);
174 if (byteSwap)
175 {
176 swap_endian(blockNum);
177 }
178 fielddef->m_elementIDs.push_back(blockNum - 1);
179
180 maxID = maxID > blockNum ? maxID : blockNum;
181 minID = minID < blockNum ? minID : blockNum;
182 }
183
184 // Determine how many fields we have to extract
185 size_t blockSize = nBlocksXYZ[0] * nBlocksXYZ[1] * nBlocksXYZ[2];
186 size_t dataSize = blockSize * nBlocks;
187
188 for (string::size_type i = 0; i < remain.size(); ++i)
189 {
190 switch (remain[i])
191 {
192 case 'U':
193 fielddef->m_fields.push_back("u");
194 fielddef->m_fields.push_back("v");
195 if (nDim == 3)
196 {
197 fielddef->m_fields.push_back("w");
198 }
199 break;
200 case 'P':
201 fielddef->m_fields.push_back("p");
202 break;
203 case 'T':
204 fielddef->m_fields.push_back("T");
205 break;
206 case '1':
207 case '2':
208 case '3':
209 fielddef->m_fields.push_back(string("scalar") + remain[i]);
210 break;
211 case ' ':
212 continue;
213 default:
214 cerr << "Field contains unknown variable: " << remain[i]
215 << endl;
216 abort();
217 }
218 }
219
220 m_f->m_data.resize(1);
221 m_f->m_data[0].resize(fielddef->m_fields.size() * dataSize);
222
223 // Now reprocess since different fields need different logic
224 for (size_t i = 0, cnt = 0; i < remain.size(); ++i)
225 {
226 switch (remain[i])
227 {
228 case 'U':
229 {
230 size_t cntVel[3] = {cnt, cnt + dataSize, cnt + 2 * dataSize};
231
232 for (size_t j = 0; j < nBlocks; ++j)
233 {
234 for (size_t k = 0; k < nDim; ++k)
235 {
236 file.read((char *)&m_f->m_data[0][cntVel[k]],
237 blockSize * sizeof(NekDouble));
238 cntVel[k] += blockSize;
239 }
240 }
241
242 cnt += nDim * dataSize;
243 break;
244 }
245 case 'P':
246 {
247 file.read((char *)&m_f->m_data[0][cnt],
248 dataSize * sizeof(NekDouble));
249 cnt += dataSize;
250 break;
251 }
252 case '1':
253 case '2':
254 case '3':
255 {
256 file.read((char *)&m_f->m_data[0][cnt],
257 dataSize * sizeof(NekDouble));
258 cnt += dataSize;
259 break;
260 }
261 case ' ':
262 continue;
263 }
264 }
265
266 m_f->m_fielddef.push_back(fielddef);
267
268 // save field names
269 m_f->m_variables = m_f->m_fielddef[0]->m_fields;
270}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
FieldSharedPtr m_f
Field object.
Definition: Module.h:239
std::map< std::string, ConfigOption > m_config
List of configuration values.
Definition: Module.h:272
void swap_endian(T &u)
Swap endian ordering of the input variable.
Definition: Module.h:98
std::shared_ptr< FieldDefinitions > FieldDefinitionsSharedPtr
Definition: FieldIO.h:184
@ eGLL_Lagrange
Lagrange for SEM basis .
Definition: BasisType.h:56
std::uint32_t NekUInt32
double NekDouble

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL0, Nektar::LibUtilities::eGLL_Lagrange, Nektar::LibUtilities::eHexahedron, Nektar::FieldUtils::Module::m_config, Nektar::FieldUtils::Module::m_f, Nektar::FieldUtils::swap_endian(), and Nektar::UnitTests::test().

Member Data Documentation

◆ m_className

ModuleKey Nektar::FieldUtils::InputNek5000::m_className
static
Initial value:
= {
"Reads Nek5000 field file.")}
static ModuleSharedPtr create(FieldSharedPtr f)
Creates an instance of this class.
Definition: InputNek5000.h:53
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:197
std::pair< ModuleType, std::string > ModuleKey
Definition: Module.h:180
ModuleFactory & GetModuleFactory()
Definition: Module.cpp:47

ModuleKey for class.

Definition at line 58 of file InputNek5000.h.