Nektar++
InputSemtex.cpp
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////////////
2 //
3 // File: InputSemtex.cpp
4 //
5 // For more information, please see: http://www.nektar.info/
6 //
7 // The MIT License
8 //
9 // Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10 // Department of Aeronautics, Imperial College London (UK), and Scientific
11 // Computing and Imaging Institute, University of Utah (USA).
12 //
13 // Permission is hereby granted, free of charge, to any person obtaining a
14 // copy of this software and associated documentation files (the "Software"),
15 // to deal in the Software without restriction, including without limitation
16 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
17 // and/or sell copies of the Software, and to permit persons to whom the
18 // Software is furnished to do so, subject to the following conditions:
19 //
20 // The above copyright notice and this permission notice shall be included
21 // in all copies or substantial portions of the Software.
22 //
23 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29 // DEALINGS IN THE SOFTWARE.
30 //
31 // Description: Reads a Semtex checkpoint file.
32 //
33 ////////////////////////////////////////////////////////////////////////////////
34 
35 #include <iostream>
36 #include <fstream>
37 #include <string>
38 #include <vector>
39 using namespace std;
40 
41 #include <boost/core/ignore_unused.hpp>
42 #include <boost/algorithm/string.hpp>
43 
45 
46 #include "InputSemtex.h"
47 
48 namespace Nektar
49 {
50 namespace FieldUtils
51 {
52 
53 ModuleKey InputSemtex::m_className[1] = {
55  ModuleKey(eInputModule, "fldsem"), InputSemtex::create,
56  "Reads Semtex field file.")
57 };
58 
59 /**
60  * @brief Set up InputSemtex object.
61  *
62  */
63 InputSemtex::InputSemtex(FieldSharedPtr f) : InputModule(f)
64 {
65  m_allowedFiles.insert("fldsem");
66 }
67 
68 /**
69  *
70  */
72 {
73 }
74 
75 /**
76  * @brief Process Semtex input file.
77  *
78  * This routine reads a binary-format Semtex field file, loads the data into
79  * memory and populates the field definitions to match the data format. Semtex
80  * is a classic nodal-Lagrangian spectral element code at a single polynomial
81  * order, meaning that the field data are set up according to this structure.
82  */
83 void InputSemtex::Process(po::variables_map &vm)
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
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 }
280 
281 }
282 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
Abstract base class for input modules.
Definition: Module.h:251
virtual void Process(po::variables_map &vm)
Process Semtex input file.
Definition: InputSemtex.cpp:83
std::set< std::string > m_allowedFiles
List of allowed file formats.
Definition: Module.h:238
FieldSharedPtr m_f
Field object.
Definition: Module.h:230
std::map< std::string, ConfigOption > m_config
List of configuration values.
Definition: Module.h:233
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:200
General purpose memory allocation routines with the ability to allocate from thread specific memory p...
def copy(self)
Definition: pycml.py:2663
std::shared_ptr< Field > FieldSharedPtr
Definition: Field.hpp:989
std::pair< ModuleType, std::string > ModuleKey
Definition: Module.h:290
void swap_endian(T &u)
Swap endian ordering of the input variable.
Definition: Module.h:102
ModuleFactory & GetModuleFactory()
Definition: Module.cpp:49
EndianType Endianness(void)
std::shared_ptr< FieldDefinitions > FieldDefinitionsSharedPtr
Definition: FieldIO.h:179
@ eGLL_Lagrange
Lagrange for SEM basis .
Definition: BasisType.h:54
@ eFourier
Fourier Expansion .
Definition: BasisType.h:53
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:1
double NekDouble
Metadata that describes the storage properties of field output.
Definition: FieldIO.h:103