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 <fstream>
36#include <iostream>
37#include <string>
38#include <vector>
39using namespace std;
40
41#include <boost/algorithm/string.hpp>
42
44
45#include "InputSemtex.h"
46
47namespace Nektar::FieldUtils
48{
49
53 "Reads Semtex field file.")};
54
55/**
56 * @brief Set up InputSemtex object.
57 *
58 */
60{
61 m_allowedFiles.insert("fldsem");
62}
63
64/**
65 *
66 */
68{
69}
70
71/**
72 * @brief Process Semtex input file.
73 *
74 * This routine reads a binary-format Semtex field file, loads the data into
75 * memory and populates the field definitions to match the data format. Semtex
76 * is a classic nodal-Lagrangian spectral element code at a single polynomial
77 * order, meaning that the field data are set up according to this structure.
78 */
79void InputSemtex::v_Process([[maybe_unused]] po::variables_map &vm)
80{
81 // Variables to be read from session file
82 string sessionName, date, fields, endian;
83 int nr, ns, nz, nelmt, step;
84 NekDouble time, dt, kinvis, beta;
85
86 ifstream file(m_config["infile"].as<string>().c_str(), ios::binary);
87
88 // -- Read header information.
89 char buf[25];
90 string line;
91
92 // Session name
93 file.read(buf, 25);
94 sessionName = string(buf, 25);
95 boost::trim(sessionName);
96 getline(file, line);
97 m_f->m_fieldMetaDataMap["SessionName0"] = sessionName;
98
99 // Date
100 file.read(buf, 25);
101 date = string(buf, 25);
102 boost::trim(date);
103 getline(file, line);
104
105 // nP, nZ, nElmt
106 file >> nr >> ns >> nz >> nelmt;
107 getline(file, line);
108
109 // Step
110 file >> step;
111 getline(file, line);
112
113 // Time
114 file >> time;
115 getline(file, line);
116 m_f->m_fieldMetaDataMap["Time"] = boost::lexical_cast<string>(time);
117
118 // Timestep
119 file >> dt;
120 getline(file, line);
121 m_f->m_fieldMetaDataMap["TimeStep"] = boost::lexical_cast<string>(dt);
122
123 // Viscosity
124 file >> kinvis;
125 getline(file, line);
126 m_f->m_fieldMetaDataMap["Kinvis"] = boost::lexical_cast<string>(kinvis);
127
128 // Beta
129 file >> beta;
130 getline(file, line);
131
132 // Fields
133 file.read(buf, 25);
134 fields = string(buf, 25);
135 boost::trim(fields);
136 getline(file, line);
137
138 // Endian-ness
140 std::string endianSearch;
141 if (systemEndian == LibUtilities::eEndianBig)
142 {
143 endianSearch = "big";
144 }
145 else if (systemEndian == LibUtilities::eEndianLittle)
146 {
147 endianSearch = "little";
148 }
149 else
150 {
151 ASSERTL0(false, "Only little- or big-endian systems are supported");
152 }
153
154 file.read(buf, 25);
155 endian = string(buf, 25);
156 bool byteSwap = endian.find(endianSearch) == string::npos;
157 getline(file, line);
158
159 // Print some basic information for input if in verbose mode.
160 if (m_f->m_verbose)
161 {
162 cout << "Found header information:" << endl;
163 cout << " -- From session : " << sessionName << endl;
164 cout << " -- File generated : " << date << endl;
165 cout << " -- Polynomial order : " << nr - 1 << endl;
166 cout << " -- Number of planes : " << nz << endl;
167 cout << " -- Number of elements : " << nelmt << endl;
168 cout << " -- Simulation time : " << time << endl;
169 cout << " -- Timestep : " << dt << endl;
170 cout << " -- Viscosity : " << kinvis << endl;
171 cout << " -- Fields : " << fields << " (" << fields.size()
172 << " total)" << endl;
173
174 if (nz > 1)
175 {
176 cout << " -- Homogeneous length : " << 2 * M_PI / beta << endl;
177 }
178
179 cout << " -- " << (byteSwap ? "" : "do not ") << "need to swap endian"
180 << endl;
181 }
182
183 ASSERTL0(nr == ns, "Semtex reader assumes values of nr and ns are equal");
184
185 // Set up a field definition
188 fielddef->m_shapeType = LibUtilities::eQuadrilateral;
189 fielddef->m_homoStrips = false;
190 fielddef->m_pointsDef = false;
191 fielddef->m_uniOrder = true;
192 fielddef->m_numPointsDef = false;
193
194 // Set up basis
195 fielddef->m_basis.push_back(LibUtilities::eGLL_Lagrange);
196 fielddef->m_basis.push_back(LibUtilities::eGLL_Lagrange);
197 fielddef->m_numModes.push_back(nr);
198 fielddef->m_numModes.push_back(nr);
199
200 // Set up elements
201 fielddef->m_elementIDs.resize(nelmt);
202 for (int i = 0; i < nelmt; ++i)
203 {
204 fielddef->m_elementIDs[i] = i;
205 }
206
207 // Deal with homogeneous direction.
208 if (nz > 1)
209 {
210 fielddef->m_numHomogeneousDir = 1;
211 fielddef->m_homogeneousLengths.push_back(2 * M_PI / beta);
212 fielddef->m_numModes.push_back(nz);
213 fielddef->m_basis.push_back(LibUtilities::eFourier);
214
215 for (int i = 0; i < nz; ++i)
216 {
217 fielddef->m_homogeneousZIDs.push_back(i);
218 }
219 }
220 else
221 {
222 fielddef->m_numHomogeneousDir = 0;
223 }
224
225 for (string::size_type i = 0; i < fields.size(); ++i)
226 {
227 fielddef->m_fields.push_back(string(&fields[i], 1));
228 }
229
230 // Size of data to read.
231 size_t elmtSize = nr * ns;
232 size_t planeSize = elmtSize * nelmt;
233 size_t fieldSize = planeSize * nz;
234 size_t dataSize = fieldSize * fields.size();
235
236 // Allocate our storage.
237 m_f->m_data.resize(1);
238 m_f->m_data[0].resize(dataSize);
239
240 // Temporary storage for one plane of data.
241 vector<NekDouble> tmp(planeSize);
242 size_t offset = nz * nr * ns;
243
244 // Now reorder data; Semtex ordering traverses memory fastest over planes,
245 // whereas Nektar++ expects it over elements
246 for (int i = 0; i < fields.size(); ++i)
247 {
248 NekDouble *data = &m_f->m_data[0][i * fieldSize];
249 for (int j = 0; j < nz; ++j)
250 {
251 size_t elSizeJ = j * elmtSize;
252 file.read((char *)&tmp[0], planeSize * sizeof(NekDouble));
253
254 if (byteSwap)
255 {
256 swap_endian(tmp);
257 }
258
259 double *x = &tmp[0];
260 for (int k = 0; k < nelmt; ++k)
261 {
262 std::copy(x, x + elmtSize, data + k * offset + elSizeJ);
263 x += elmtSize;
264 }
265 }
266 }
267
268 m_f->m_fielddef.push_back(fielddef);
269
270 // save field names
271 m_f->m_variables = m_f->m_fielddef[0]->m_fields;
272}
273
274} // namespace Nektar::FieldUtils
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
Abstract base class for input modules.
Definition: Module.h:287
InputSemtex(FieldSharedPtr f)
Set up InputSemtex object.
Definition: InputSemtex.cpp:59
static ModuleSharedPtr create(FieldSharedPtr f)
Creates an instance of this class.
Definition: InputSemtex.h:53
void v_Process(po::variables_map &vm) override
Process Semtex input file.
Definition: InputSemtex.cpp:79
static ModuleKey m_className[]
ModuleKey for class.
Definition: InputSemtex.h:58
std::set< std::string > m_allowedFiles
List of allowed file formats.
Definition: Module.h:274
FieldSharedPtr m_f
Field object.
Definition: Module.h:239
std::map< std::string, ConfigOption > m_config
List of configuration values.
Definition: Module.h:272
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:197
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
def copy(self)
Definition: pycml.py:2663
std::shared_ptr< Field > FieldSharedPtr
Definition: Field.hpp:990
std::pair< ModuleType, std::string > ModuleKey
Definition: Module.h:180
void swap_endian(T &u)
Swap endian ordering of the input variable.
Definition: Module.h:98
ModuleFactory & GetModuleFactory()
Definition: Module.cpp:47
EndianType Endianness(void)
std::shared_ptr< FieldDefinitions > FieldDefinitionsSharedPtr
Definition: FieldIO.h:184
@ beta
Gauss Radau pinned at x=-1,.
Definition: PointsType.h:59
@ eGLL_Lagrange
Lagrange for SEM basis .
Definition: BasisType.h:56
@ eFourier
Fourier Expansion .
Definition: BasisType.h:55
double NekDouble