Nektar++
Public Member Functions | Static Public Member Functions | Static Public Attributes | 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...
 
virtual ~InputNek5000 ()
 
virtual void Process (po::variables_map &vm)
 Process Nek5000 input file. More...
 
virtual std::string GetModuleName ()
 
virtual std::string GetModuleDescription ()
 
virtual ModulePriority GetModulePriority ()
 
- Public Member Functions inherited from Nektar::FieldUtils::InputModule
 InputModule (FieldSharedPtr p_m)
 
FIELD_UTILS_EXPORT void AddFile (std::string fileType, std::string fileName)
 
- Public Member Functions inherited from Nektar::FieldUtils::Module
FIELD_UTILS_EXPORT Module (FieldSharedPtr p_f)
 
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 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...
 

Additional Inherited Members

- Protected Member Functions inherited from Nektar::FieldUtils::InputModule
void PrintSummary ()
 Print summary of elements. More...
 
- Protected Member Functions inherited from Nektar::FieldUtils::Module
 Module ()
 
- Protected Attributes inherited from Nektar::FieldUtils::InputModule
std::set< std::string > m_allowedFiles
 
- Protected Attributes inherited from Nektar::FieldUtils::Module
FieldSharedPtr m_f
 Field object. More...
 
std::map< std::string, ConfigOptionm_config
 List of configuration values. More...
 

Detailed Description

Converter for Fld files.

Definition at line 48 of file library/FieldUtils/InputModules/InputNek5000.h.

Constructor & Destructor Documentation

◆ InputNek5000()

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

Set up InputNek5000 object.

Definition at line 61 of file library/FieldUtils/InputModules/InputNek5000.cpp.

References Nektar::FieldUtils::InputModule::m_allowedFiles.

61  : InputModule(f)
62 {
63  m_allowedFiles.insert("fld5000");
64 }
std::set< std::string > m_allowedFiles

◆ ~InputNek5000()

Nektar::FieldUtils::InputNek5000::~InputNek5000 ( )
virtual

Definition at line 69 of file library/FieldUtils/InputModules/InputNek5000.cpp.

70 {
71 }

Member Function Documentation

◆ create()

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

Creates an instance of this class.

Definition at line 56 of file library/FieldUtils/InputModules/InputNek5000.h.

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

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

◆ GetModuleDescription()

virtual std::string Nektar::FieldUtils::InputNek5000::GetModuleDescription ( )
inlinevirtual

Reimplemented from Nektar::FieldUtils::Module.

Definition at line 68 of file library/FieldUtils/InputModules/InputNek5000.h.

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

◆ GetModuleName()

virtual std::string Nektar::FieldUtils::InputNek5000::GetModuleName ( )
inlinevirtual

Implements Nektar::FieldUtils::Module.

Definition at line 63 of file library/FieldUtils/InputModules/InputNek5000.h.

64  {
65  return "InputNek5000";
66  }

◆ GetModulePriority()

virtual ModulePriority Nektar::FieldUtils::InputNek5000::GetModulePriority ( )
inlinevirtual

◆ Process()

void Nektar::FieldUtils::InputNek5000::Process ( po::variables_map &  vm)
virtual

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.

Implements Nektar::FieldUtils::Module.

Definition at line 84 of file library/FieldUtils/InputModules/InputNek5000.cpp.

References ASSERTL0, Nektar::LibUtilities::eGLL_Lagrange, Nektar::LibUtilities::eHexahedron, Nektar::FieldUtils::Module::m_config, Nektar::FieldUtils::Module::m_f, and Nektar::FieldUtils::swap_endian().

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

Member Data Documentation

◆ m_className

ModuleKey Nektar::FieldUtils::InputNek5000::m_className
static
Initial value:
= {
"Reads Nek5000 field file.")
}

ModuleKey for class.

Definition at line 61 of file library/FieldUtils/InputModules/InputNek5000.h.