Nektar++
Public Member Functions | Static Public Member Functions | Static Public Attributes | Protected Member Functions | 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 ()
 
- 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 ()
 
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

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

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 48 of file InputSemtex.h.

Constructor & Destructor Documentation

◆ InputSemtex()

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

Set up InputSemtex object.

Definition at line 62 of file InputSemtex.cpp.

62  : InputModule(f)
63 {
64  m_allowedFiles.insert("fldsem");
65 }
InputModule(FieldSharedPtr p_m)
Definition: Module.cpp:63
std::set< std::string > m_allowedFiles
List of allowed file formats.
Definition: Module.h:265

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

◆ ~InputSemtex()

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

Definition at line 70 of file InputSemtex.cpp.

71 {
72 }

Member Function Documentation

◆ create()

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

Creates an instance of this class.

Definition at line 55 of file InputSemtex.h.

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

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

◆ v_GetModuleDescription()

virtual std::string Nektar::FieldUtils::InputSemtex::v_GetModuleDescription ( )
inlineoverrideprotectedvirtual

Reimplemented from Nektar::FieldUtils::Module.

Definition at line 70 of file InputSemtex.h.

71  {
72  return "Processing Semtex field file";
73  }

◆ v_GetModuleName()

virtual std::string Nektar::FieldUtils::InputSemtex::v_GetModuleName ( )
inlineoverrideprotectedvirtual

Reimplemented from Nektar::FieldUtils::Module.

Definition at line 65 of file InputSemtex.h.

66  {
67  return "InputSemtex";
68  }

◆ v_GetModulePriority()

virtual ModulePriority Nektar::FieldUtils::InputSemtex::v_GetModulePriority ( )
inlineoverrideprotectedvirtual

Reimplemented from Nektar::FieldUtils::Module.

Definition at line 75 of file InputSemtex.h.

76  {
77  return eCreateFieldData;
78  }

References Nektar::FieldUtils::eCreateFieldData.

◆ v_Process()

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

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.

Reimplemented from Nektar::FieldUtils::Module.

Definition at line 82 of file InputSemtex.cpp.

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

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL0, Nektar::LibUtilities::beta, 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().

Member Data Documentation

◆ m_className

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

ModuleKey for class.

Definition at line 60 of file InputSemtex.h.