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>
37 using namespace std;
38 
39 #include <boost/core/ignore_unused.hpp>
40 #include <boost/math/special_functions/fpclassify.hpp>
41 
43 
44 #include "ProcessPointDataToFld.h"
45 
46 namespace Nektar
47 {
48 namespace FieldUtils
49 {
50 
51 ModuleKey ProcessPointDataToFld::className =
53  ModuleKey(eProcessModule, "pointdatatofld"),
54  ProcessPointDataToFld::create,
55  "Given discrete data at quadrature points project them onto an "
56  "expansion"
57  "basis and output fld file. Requires frompts and .xml and .fld files.");
58 
59 ProcessPointDataToFld::ProcessPointDataToFld(FieldSharedPtr f)
60  : ProcessModule(f)
61 {
62  m_config["setnantovalue"] = ConfigOption(
63  false, "NotSet", "reset any nan value to prescribed value");
64 
65  m_config["frompts"] = ConfigOption(
66  false, "NotSet", "Pts file from which to interpolate field");
67 }
68 
70 {
71 }
72 
73 void ProcessPointDataToFld::Process(po::variables_map &vm)
74 {
75  m_f->SetUpExp(vm);
76 
77  int i, j;
78  bool setnantovalue = false;
79  NekDouble defvalue=0.0;
80 
81  if (!boost::iequals(m_config["setnantovalue"].as<string>(), "NotSet"))
82  {
83  setnantovalue = true;
84  defvalue = m_config["setnantovalue"].as<NekDouble>();
85  }
86 
87  // Check for command line point specification if no .pts file specified
88  // Load pts file
90  ASSERTL0( m_config["frompts"].as<string>().compare("NotSet") != 0,
91  "ProcessInterpPointDataToFld requires frompts parameter");
92  string inFile = m_config["frompts"].as<string>().c_str();
97  ptsIO->Import(inFile, fieldPts);
98 
99  int nFields = fieldPts->GetNFields();
100  ASSERTL0(nFields > 0, "No field values provided in input");
101 
102  int dim = fieldPts->GetDim();
103 
104  // assume one field is already defined from input file.
105  ASSERTL0(m_f->m_numHomogeneousDir == 0,
106  "ProcessInterpPointDataToFld does not support homogeneous expansion");
107 
108  m_f->m_exp.resize(nFields);
109  for (i = 1; i < nFields; ++i)
110  {
111  m_f->m_exp[i] = m_f->AppendExpList(m_f->m_numHomogeneousDir);
112  }
114  fieldPts->GetPts(pts);
115 
116  // set any nan values to default value if requested
117  if (setnantovalue)
118  {
119  for (int i = 0; i < pts[0].size(); ++i)
120  {
121  for (int j = 0; j < nFields; ++j)
122  {
123  if ((boost::math::isnan)(pts[j + dim][i]))
124  {
125  pts[j + dim][i] = defvalue;
126  }
127  }
128  }
129  }
130 
131  if (fieldPts->m_ptsInfo.count(LibUtilities::eIsEquiSpacedData) != 0)
132  {
133  int totcoeffs = m_f->m_exp[0]->GetNcoeffs();
134 
135  ASSERTL0(pts[0].size() != totcoeffs,
136  "length of points in .pts file is different "
137  "to the number of coefficients in expansion ");
138 
139  for (int i = 0; i < nFields; ++i)
140  {
141  Array<OneD, NekDouble> coeffs = m_f->m_exp[i]->UpdateCoeffs(), tmp;
142  int cnt = 0;
143  for (int e = 0; e < m_f->m_exp[0]->GetNumElmts(); ++e)
144  {
145  int ncoeffs = m_f->m_exp[i]->GetExp(e)->GetNcoeffs();
146  int offset = m_f->m_exp[i]->GetCoeff_Offset(e);
147  Vmath::Vcopy(ncoeffs, &pts[i + dim][cnt], 1, &coeffs[offset],
148  1);
149 
150  m_f->m_exp[i]->GetExp(e)->EquiSpacedToCoeffs(
151  coeffs + offset, tmp = coeffs + offset);
152  cnt += ncoeffs;
153  }
154  m_f->m_exp[i]->BwdTrans(m_f->m_exp[i]->GetCoeffs(),
155  m_f->m_exp[i]->UpdatePhys());
156  }
157  }
158  else
159  {
160  int totpoints = m_f->m_exp[0]->GetTotPoints();
162 
163  coords[0] = Array<OneD, NekDouble>(totpoints);
164  coords[1] = Array<OneD, NekDouble>(totpoints);
165  coords[2] = Array<OneD, NekDouble>(totpoints);
166 
167  m_f->m_exp[0]->GetCoords(coords[0], coords[1], coords[2]);
168 
169  if (pts[0].size() != totpoints)
170  {
171  WARNINGL0(false, "length of points in .pts file is different to "
172  "the number of quadrature points in xml file");
173  totpoints = min(totpoints, (int)pts[0].size());
174  }
175 
176  int mismatch = 0;
177  for (i = 0; i < totpoints; ++i)
178  {
179  for (j = 0; j < dim; ++j)
180  {
181  if (fabs(coords[j][i] - pts[j][i]) > 1e-4)
182  {
183  string outstring =
184  "Coordinates do not match within 1e-4: " +
185  boost::lexical_cast<string>(coords[j][i]) + " versus " +
186  boost::lexical_cast<string>(pts[j][i]) + " diff " +
187  boost::lexical_cast<string>(
188  fabs(coords[j][i] - pts[j][i]));
189  ;
190  WARNINGL0(false, outstring);
191  mismatch += 1;
192  }
193  }
194 
195  for (j = 0; j < nFields; ++j)
196  {
197  m_f->m_exp[j]->SetPhys(i, pts[j + dim][i]);
198  }
199  }
200 
201  // forward transform fields
202  for (i = 0; i < nFields; ++i)
203  {
204  m_f->m_exp[i]->FwdTrans_IterPerExp(m_f->m_exp[i]->GetPhys(),
205  m_f->m_exp[i]->UpdateCoeffs());
206  }
207  }
208 
209  // save field names
210  for (int j = 0; j < fieldPts->GetNFields(); ++j)
211  {
212  m_f->m_variables.push_back(fieldPts->GetFieldName(j));
213  }
214 }
215 }
216 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
#define WARNINGL0(condition, msg)
Definition: ErrorUtil.hpp:223
FieldSharedPtr m_f
Field object.
Definition: Module.h:230
std::map< std::string, ConfigOption > m_config
List of configuration values.
Definition: Module.h:233
Abstract base class for processing modules.
Definition: Module.h:265
virtual void Process(po::variables_map &vm)
Write mesh to output file.
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:200
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:145
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
std::shared_ptr< Field > FieldSharedPtr
Definition: Field.hpp:989
std::pair< ModuleType, std::string > ModuleKey
Definition: Module.h:290
ModuleFactory & GetModuleFactory()
Definition: Module.cpp:49
std::shared_ptr< PtsField > PtsFieldSharedPtr
Definition: PtsField.h:182
std::shared_ptr< PtsIO > PtsIOSharedPtr
Definition: PtsIO.h:101
CommFactory & GetCommFactory()
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:54
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:1
double NekDouble
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1199
Represents a command-line configuration option.
Definition: Module.h:134