Nektar++
ProcessPointDataToFld.cpp
Go to the documentation of this file.
1////////////////////////////////////////////////////////////////////////////////
2//
3// File: ProcessPointDataToFld.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: Load a file of interpolated vales at physical
32// quadrature points and project to a fld file
33//
34////////////////////////////////////////////////////////////////////////////////
35#include <iostream>
36#include <string>
37using namespace std;
38
42
43namespace Nektar::FieldUtils
44{
45
48 ModuleKey(eProcessModule, "pointdatatofld"),
50 "Given discrete data at quadrature points project them onto an "
51 "expansion"
52 "basis and output fld file. Requires frompts and .xml and .fld files.");
53
55 : ProcessModule(f)
56{
57 m_config["setnantovalue"] = ConfigOption(
58 false, "NotSet", "reset any nan value to prescribed value");
59
60 m_config["frompts"] = ConfigOption(
61 false, "NotSet", "Pts file from which to interpolate field");
62}
63
65{
66}
67
68void ProcessPointDataToFld::v_Process(po::variables_map &vm)
69{
70 m_f->SetUpExp(vm);
71
72 int i, j;
73 bool setnantovalue = false;
74 NekDouble defvalue = 0.0;
75
76 if (!boost::iequals(m_config["setnantovalue"].as<string>(), "NotSet"))
77 {
78 setnantovalue = true;
79 defvalue = m_config["setnantovalue"].as<NekDouble>();
80 }
81
82 // Check for command line point specification if no .pts file specified
83 // Load pts file
85 ASSERTL0(m_config["frompts"].as<string>().compare("NotSet") != 0,
86 "ProcessInterpPointDataToFld requires frompts parameter");
87 string inFile = m_config["frompts"].as<string>().c_str();
90
91 // Determine file format from file extension
92 if (fs::path(inFile).extension() == ".pts")
93 {
94 LibUtilities::PtsIO(c).Import(inFile, fieldPts);
95 }
96 else if (fs::path(inFile).extension() == ".csv")
97 {
98 LibUtilities::CsvIO(c).Import(inFile, fieldPts);
99 }
100 else
101 {
103 "Unsupported file format for the \"frompts\" parameter. "
104 "Supported formats: \".pts\" and \".csv\"");
105 }
106
107 int nFields = fieldPts->GetNFields();
108 ASSERTL0(nFields > 0, "No field values provided in input");
109
110 int dim = fieldPts->GetDim();
111
112 // assume one field is already defined from input file.
113 ASSERTL0(
114 m_f->m_numHomogeneousDir == 0,
115 "ProcessInterpPointDataToFld does not support homogeneous expansion");
116
117 m_f->m_exp.resize(nFields);
118 for (i = 1; i < nFields; ++i)
119 {
120 m_f->m_exp[i] = m_f->AppendExpList(m_f->m_numHomogeneousDir);
121 }
122
124 fieldPts->GetPts(pts);
125
126 // set any nan values to default value if requested
127 if (setnantovalue)
128 {
129 for (int i = 0; i < pts[0].size(); ++i)
130 {
131 for (int j = 0; j < nFields; ++j)
132 {
133 if ((boost::math::isnan)(pts[j + dim][i]))
134 {
135 pts[j + dim][i] = defvalue;
136 }
137 }
138 }
139 }
140
141 if (fieldPts->m_ptsInfo.count(LibUtilities::eIsEquiSpacedData) != 0)
142 {
143 int totcoeffs = m_f->m_exp[0]->GetNcoeffs();
144
145 ASSERTL0(pts[0].size() != totcoeffs,
146 "length of points in .pts file is different "
147 "to the number of coefficients in expansion ");
148
149 for (int i = 0; i < nFields; ++i)
150 {
151 Array<OneD, NekDouble> coeffs = m_f->m_exp[i]->UpdateCoeffs(), tmp;
152 int cnt = 0;
153 for (int e = 0; e < m_f->m_exp[0]->GetNumElmts(); ++e)
154 {
155 int ncoeffs = m_f->m_exp[i]->GetExp(e)->GetNcoeffs();
156 int offset = m_f->m_exp[i]->GetCoeff_Offset(e);
157 Vmath::Vcopy(ncoeffs, &pts[i + dim][cnt], 1, &coeffs[offset],
158 1);
159
160 m_f->m_exp[i]->GetExp(e)->EquiSpacedToCoeffs(
161 coeffs + offset, tmp = coeffs + offset);
162 cnt += ncoeffs;
163 }
164 m_f->m_exp[i]->BwdTrans(m_f->m_exp[i]->GetCoeffs(),
165 m_f->m_exp[i]->UpdatePhys());
166 }
167 }
168 else
169 {
170 int totpoints = m_f->m_exp[0]->GetTotPoints();
172
173 coords[0] = Array<OneD, NekDouble>(totpoints);
174 coords[1] = Array<OneD, NekDouble>(totpoints);
175 coords[2] = Array<OneD, NekDouble>(totpoints);
176
177 m_f->m_exp[0]->GetCoords(coords[0], coords[1], coords[2]);
178
179 if (pts[0].size() != totpoints)
180 {
181 WARNINGL0(false, "length of points in .pts file is different to "
182 "the number of quadrature points in xml file");
183 totpoints = min(totpoints, (int)pts[0].size());
184 }
185
186 for (j = 0; j < dim; ++j)
187 {
188 for (i = 0; i < totpoints; ++i)
189 {
190 if (fabs(coords[j][i] - pts[j][i]) > 1e-4)
191 {
192 string outstring =
193 "Coordinates do not match within 1e-4: " +
194 boost::lexical_cast<string>(coords[j][i]) + " versus " +
195 boost::lexical_cast<string>(pts[j][i]) + " diff " +
196 boost::lexical_cast<string>(
197 fabs(coords[j][i] - pts[j][i]));
198 ;
199 WARNINGL0(false, outstring);
200 }
201 }
202 }
203
204 for (j = 0; j < nFields; ++j)
205 {
206 Array<OneD, NekDouble> phys = m_f->m_exp[j]->UpdatePhys();
207
208 for (i = 0; i < totpoints; ++i)
209 {
210 phys[i] = pts[j + dim][i];
211 }
212
213 // forward transform fields
214 m_f->m_exp[j]->FwdTransLocalElmt(phys,
215 m_f->m_exp[j]->UpdateCoeffs());
216 }
217 }
218
219 // save field names
220 for (int j = 0; j < fieldPts->GetNFields(); ++j)
221 {
222 m_f->m_variables.push_back(fieldPts->GetFieldName(j));
223 }
224}
225} // namespace Nektar::FieldUtils
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:202
#define WARNINGL0(condition, msg)
Definition: ErrorUtil.hpp:215
FieldSharedPtr m_f
Field object.
Definition: Module.h:239
std::map< std::string, ConfigOption > m_config
List of configuration values.
Definition: Module.h:272
Abstract base class for processing modules.
Definition: Module.h:301
void v_Process(po::variables_map &vm) override
Write mesh to output file.
static std::shared_ptr< Module > create(FieldSharedPtr f)
Creates an instance of this class.
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:197
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:143
void Import(const std::string &inFile, PtsFieldSharedPtr &ptsField, FieldMetaDataMap &fieldmetadatamap=NullFieldMetaDataMap, DomainRangeShPtr &Range=NullDomainRangeShPtr)
Import a pts field from file.
Definition: PtsIO.cpp:66
std::shared_ptr< Field > FieldSharedPtr
Definition: Field.hpp:990
std::pair< ModuleType, std::string > ModuleKey
Definition: Module.h:180
ModuleFactory & GetModuleFactory()
Definition: Module.cpp:47
std::shared_ptr< PtsField > PtsFieldSharedPtr
Definition: PtsField.h:184
CommFactory & GetCommFactory()
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:55
double NekDouble
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825
Represents a command-line configuration option.
Definition: Module.h:129