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

#include <InputSemtex.h>

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

Public Member Functions

 InputSemtex (FieldSharedPtr f)
 Set up InputSemtex object. More...
 
virtual ~InputSemtex ()
 
virtual void Process (po::variables_map &vm)
 Process Semtex 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 InputSemtex.h.

Constructor & Destructor Documentation

◆ InputSemtex()

Nektar::FieldUtils::InputSemtex::InputSemtex ( FieldSharedPtr  f)

Set up InputSemtex object.

Definition at line 63 of file InputSemtex.cpp.

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

63  : InputModule(f)
64 {
65  m_allowedFiles.insert("fldsem");
66 }
std::set< std::string > m_allowedFiles

◆ ~InputSemtex()

Nektar::FieldUtils::InputSemtex::~InputSemtex ( )
virtual

Definition at line 71 of file InputSemtex.cpp.

72 {
73 }

Member Function Documentation

◆ create()

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

Creates an instance of this class.

Definition at line 56 of file InputSemtex.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::InputSemtex::GetModuleDescription ( )
inlinevirtual

Reimplemented from Nektar::FieldUtils::Module.

Definition at line 68 of file InputSemtex.h.

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

◆ GetModuleName()

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

Implements Nektar::FieldUtils::Module.

Definition at line 63 of file InputSemtex.h.

64  {
65  return "InputSemtex";
66  }

◆ GetModulePriority()

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

◆ Process()

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

Process Semtex input file.

This routine reads a binary-format Semtex field file, loads the data into memory and populates the field definitions to match the data format. Semtex 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.

Implements Nektar::FieldUtils::Module.

Definition at line 83 of file InputSemtex.cpp.

References ASSERTL0, CellMLToNektar.pycml::copy(), Nektar::LibUtilities::eEndianBig, Nektar::LibUtilities::eEndianLittle, Nektar::LibUtilities::eFourier, Nektar::LibUtilities::eGLL_Lagrange, Nektar::LibUtilities::Endianness(), Nektar::LibUtilities::eQuadrilateral, Nektar::FieldUtils::Module::m_config, Nektar::FieldUtils::Module::m_f, and Nektar::FieldUtils::swap_endian().

84 {
85  boost::ignore_unused(vm);
86 
87  // Variables to be read from session file
88  string sessionName, date, fields, endian;
89  int nr, ns, nz, nelmt, step;
90  NekDouble time, dt, kinvis, beta;
91 
92  ifstream file(m_config["infile"].as<string>().c_str(), ios::binary);
93 
94  // -- Read header information.
95  char buf[25];
96  string line;
97 
98  // Session name
99  file.read(buf, 25);
100  sessionName = string(buf, 25);
101  boost::trim(sessionName);
102  getline(file, line);
103  m_f->m_fieldMetaDataMap["SessionName0"] = sessionName;
104 
105  // Date
106  file.read(buf, 25);
107  date = string(buf, 25);
108  boost::trim(date);
109  getline(file, line);
110 
111  // nP, nZ, nElmt
112  file >> nr >> ns >> nz >> nelmt;
113  getline(file, line);
114 
115  // Step
116  file >> step;
117  getline(file, line);
118 
119  // Time
120  file >> time;
121  getline(file, line);
122  m_f->m_fieldMetaDataMap["Time"] = boost::lexical_cast<string>(time);
123 
124  // Timestep
125  file >> dt;
126  getline(file, line);
127  m_f->m_fieldMetaDataMap["TimeStep"] = boost::lexical_cast<string>(dt);
128 
129  // Viscosity
130  file >> kinvis;
131  getline(file, line);
132  m_f->m_fieldMetaDataMap["Kinvis"] = boost::lexical_cast<string>(kinvis);
133 
134  // Beta
135  file >> beta;
136  getline(file, line);
137 
138  // Fields
139  file.read(buf, 25);
140  fields = string(buf, 25);
141  boost::trim(fields);
142  getline(file, line);
143 
144  // Endian-ness
146  std::string endianSearch;
147  if (systemEndian == LibUtilities::eEndianBig)
148  {
149  endianSearch = "big";
150  }
151  else if (systemEndian == LibUtilities::eEndianLittle)
152  {
153  endianSearch = "little";
154  }
155  else
156  {
157  ASSERTL0(false, "Only little- or big-endian systems are supported");
158  }
159 
160  file.read(buf, 25);
161  endian = string(buf, 25);
162  bool byteSwap = endian.find(endianSearch) == string::npos;
163  getline(file, line);
164 
165  // Print some basic information for input if in verbose mode.
166  if (m_f->m_verbose)
167  {
168  cout << "Found header information:" << endl;
169  cout << " -- From session : " << sessionName << endl;
170  cout << " -- File generated : " << date << endl;
171  cout << " -- Polynomial order : " << nr-1 << endl;
172  cout << " -- Number of planes : " << nz << endl;
173  cout << " -- Number of elements : " << nelmt << endl;
174  cout << " -- Simulation time : " << time << endl;
175  cout << " -- Timestep : " << dt << endl;
176  cout << " -- Viscosity : " << kinvis << endl;
177  cout << " -- Fields : " << fields
178  << " (" << fields.size() << " total)" << endl;
179 
180  if (nz > 1)
181  {
182  cout << " -- Homogeneous length : " << 2*M_PI/beta << endl;
183  }
184 
185  cout << " -- " << (byteSwap ? "" : "do not ") << "need to swap endian"
186  << endl;
187  }
188 
189  ASSERTL0(nr == ns, "Semtex reader assumes values of nr and ns are equal");
190 
191  // Set up a field definition
192  LibUtilities::FieldDefinitionsSharedPtr fielddef = MemoryManager<
193  LibUtilities::FieldDefinitions>::AllocateSharedPtr();
194  fielddef->m_shapeType = LibUtilities::eQuadrilateral;
195  fielddef->m_homoStrips = false;
196  fielddef->m_pointsDef = false;
197  fielddef->m_uniOrder = true;
198  fielddef->m_numPointsDef = false;
199 
200  // Set up basis
201  fielddef->m_basis.push_back(LibUtilities::eGLL_Lagrange);
202  fielddef->m_basis.push_back(LibUtilities::eGLL_Lagrange);
203  fielddef->m_numModes.push_back(nr);
204  fielddef->m_numModes.push_back(nr);
205 
206  // Set up elements
207  fielddef->m_elementIDs.resize(nelmt);
208  for (int i = 0; i < nelmt; ++i)
209  {
210  fielddef->m_elementIDs[i] = i;
211  }
212 
213  // Deal with homogeneous direction.
214  if (nz > 1)
215  {
216  fielddef->m_numHomogeneousDir = 1;
217  fielddef->m_homogeneousLengths.push_back(2 * M_PI / beta);
218  fielddef->m_numModes.push_back(nz);
219  fielddef->m_basis.push_back(LibUtilities::eFourier);
220 
221  for (int i = 0; i < nz; ++i)
222  {
223  fielddef->m_homogeneousZIDs.push_back(i);
224  }
225  }
226  else
227  {
228  fielddef->m_numHomogeneousDir = 0;
229  }
230 
231  for (string::size_type i = 0; i < fields.size(); ++i)
232  {
233  fielddef->m_fields.push_back(string(&fields[i], 1));
234  }
235 
236  // Size of data to read.
237  size_t elmtSize = nr * ns;
238  size_t planeSize = elmtSize * nelmt;
239  size_t fieldSize = planeSize * nz;
240  size_t dataSize = fieldSize * fields.size();
241 
242  // Allocate our storage.
243  m_f->m_data.resize(1);
244  m_f->m_data[0].resize(dataSize);
245 
246  // Temporary storage for one plane of data.
247  vector<NekDouble> tmp(planeSize);
248  size_t offset = nz * nr * ns;
249 
250  // Now reorder data; Semtex ordering traverses memory fastest over planes,
251  // whereas Nektar++ expects it over elements
252  for (int i = 0; i < fields.size(); ++i)
253  {
254  NekDouble *data = &m_f->m_data[0][i * fieldSize];
255  for (int j = 0; j < nz; ++j)
256  {
257  size_t elSizeJ = j * elmtSize;
258  file.read((char *)&tmp[0], planeSize * sizeof(NekDouble));
259 
260  if (byteSwap)
261  {
262  swap_endian(tmp);
263  }
264 
265  double* x = &tmp[0];
266  for (int k = 0; k < nelmt; ++k)
267  {
268  std::copy(x, x + elmtSize, data + k * offset + elSizeJ);
269  x += elmtSize;
270 
271  }
272  }
273  }
274 
275  m_f->m_fielddef.push_back(fielddef);
276 
277  // save field names
278  m_f->m_variables = m_f->m_fielddef[0]->m_fields;
279 }
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.
Fourier Expansion .
Definition: BasisType.h:53
def copy(self)
Definition: pycml.py:2663
double NekDouble
EndianType Endianness(void)
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::InputSemtex::m_className
static
Initial value:
= {
"Reads Semtex field file.")
}

ModuleKey for class.

Definition at line 61 of file InputSemtex.h.