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#include <boost/core/ignore_unused.hpp>
43
45
46#include "InputSemtex.h"
47
48namespace Nektar
49{
50namespace FieldUtils
51{
52
56 "Reads Semtex field file.")};
57
58/**
59 * @brief Set up InputSemtex object.
60 *
61 */
63{
64 m_allowedFiles.insert("fldsem");
65}
66
67/**
68 *
69 */
71{
72}
73
74/**
75 * @brief Process Semtex input file.
76 *
77 * This routine reads a binary-format Semtex field file, loads the data into
78 * memory and populates the field definitions to match the data format. Semtex
79 * is a classic nodal-Lagrangian spectral element code at a single polynomial
80 * order, meaning that the field data are set up according to this structure.
81 */
82void InputSemtex::v_Process(po::variables_map &vm)
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}
278
279} // namespace FieldUtils
280} // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
Abstract base class for input modules.
Definition: Module.h:278
InputSemtex(FieldSharedPtr f)
Set up InputSemtex object.
Definition: InputSemtex.cpp:62
static ModuleSharedPtr create(FieldSharedPtr f)
Creates an instance of this class.
Definition: InputSemtex.h:55
virtual void v_Process(po::variables_map &vm) override
Process Semtex input file.
Definition: InputSemtex.cpp:82
static ModuleKey m_className[]
ModuleKey for class.
Definition: InputSemtex.h:60
std::set< std::string > m_allowedFiles
List of allowed file formats.
Definition: Module.h:265
FieldSharedPtr m_f
Field object.
Definition: Module.h:234
std::map< std::string, ConfigOption > m_config
List of configuration values.
Definition: Module.h:263
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:198
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:991
std::pair< ModuleType, std::string > ModuleKey
Definition: Module.h:317
void swap_endian(T &u)
Swap endian ordering of the input variable.
Definition: Module.h:100
ModuleFactory & GetModuleFactory()
Definition: Module.cpp:49
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
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
double NekDouble