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