Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 // License for the specific language governing rights and limitations under
14 // Permission is hereby granted, free of charge, to any person obtaining a
15 // copy of this software and associated documentation files (the "Software"),
16 // to deal in the Software without restriction, including without limitation
17 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
18 // and/or sell copies of the Software, and to permit persons to whom the
19 // Software is furnished to do so, subject to the following conditions:
20 //
21 // The above copyright notice and this permission notice shall be included
22 // in all copies or substantial portions of the Software.
23 //
24 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
25 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
26 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
27 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
28 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
29 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
30 // DEALINGS IN THE SOFTWARE.
31 //
32 // Description: Load a file of interpolated vales at physical
33 // quadrature points and project to a fld file
34 //
35 ////////////////////////////////////////////////////////////////////////////////
36 #include <iostream>
37 #include <string>
38 
39 using namespace std;
40 
41 #include "ProcessPointDataToFld.h"
42 
45 #include <boost/math/special_functions/fpclassify.hpp>
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 .pts .xml and .fld files.");
58 
59 ProcessPointDataToFld::ProcessPointDataToFld(FieldSharedPtr f)
60  : ProcessModule(f)
61 {
62  m_requireEquiSpaced = true;
63 
64  m_config["setnantovalue"] = ConfigOption(
65  false, "NotSet", "reset any nan value to prescribed value");
66 
67  if ((f->m_inputfiles.count("pts") == 0))
68  {
69  cout << endl
70  << "A pts input file must be specified for the boundary "
71  "extraction module"
72  << endl;
73 
74  cout
75  << "Usage: Fieldconvert -m pointdatatofld file.pts file.xml out.fld"
76  << endl;
77  exit(3);
78  }
79 
80  if ((f->m_inputfiles.count("xml") == 0) &&
81  (f->m_inputfiles.count("xml.gz") == 0))
82  {
83  cout << endl
84  << "An xml or xml.gz input file must be specified for the "
85  "boundary extraction module"
86  << endl;
87  cout
88  << "Usage: Fieldconvert -m pointdatatofld file.pts file.xml out.fld"
89  << endl;
90  exit(3);
91  }
92 }
93 
95 {
96 }
97 
98 void ProcessPointDataToFld::Process(po::variables_map &vm)
99 {
100  if (m_f->m_verbose)
101  {
102  if (m_f->m_comm->TreatAsRankZero())
103  {
104  cout << "ProcessPointDataToFld: projecting data to expansion..."
105  << endl;
106  }
107  }
108 
109  int i, j;
110  bool setnantovalue = false;
111  NekDouble defvalue=0.0;
112 
113  if (!boost::iequals(m_config["setnantovalue"].as<string>(), "NotSet"))
114  {
115  setnantovalue = true;
116  defvalue = m_config["setnantovalue"].as<NekDouble>();
117  }
118 
119  // Check for command line point specification if no .pts file specified
120  ASSERTL0(m_f->m_fieldPts != LibUtilities::NullPtsField,
121  "No input points found");
122 
123  int nFields = m_f->m_fieldPts->GetNFields();
124  ASSERTL0(nFields > 0, "No field values provided in input");
125 
126  int dim = m_f->m_fieldPts->GetDim();
127 
128  // assume one field is already defined from input file.
129  m_f->m_exp.resize(nFields);
130  for (i = 1; i < nFields; ++i)
131  {
132  m_f->m_exp[i] = m_f->AppendExpList(0);
133  }
134 
136  m_f->m_fieldPts->GetPts(pts);
137 
138  // set any nan values to default value if requested
139  if (setnantovalue)
140  {
141  for (int i = 0; i < pts[0].num_elements(); ++i)
142  {
143  for (int j = 0; j < nFields; ++j)
144  {
145  if ((boost::math::isnan)(pts[j + dim][i]))
146  {
147  pts[j + dim][i] = defvalue;
148  }
149  }
150  }
151  }
152 
153  if (m_f->m_fieldPts->m_ptsInfo.count(LibUtilities::eIsEquiSpacedData) != 0)
154  {
155  int totcoeffs = m_f->m_exp[0]->GetNcoeffs();
156 
157  ASSERTL0(pts[0].num_elements() != totcoeffs,
158  "length of points in .pts file is different "
159  "to the number of coefficients in expansion ");
160 
161  for (int i = 0; i < nFields; ++i)
162  {
163  Array<OneD, NekDouble> coeffs = m_f->m_exp[i]->UpdateCoeffs(), tmp;
164  int cnt = 0;
165  for (int e = 0; e < m_f->m_exp[0]->GetNumElmts(); ++e)
166  {
167  int ncoeffs = m_f->m_exp[i]->GetExp(e)->GetNcoeffs();
168  int offset = m_f->m_exp[i]->GetCoeff_Offset(e);
169  Vmath::Vcopy(ncoeffs, &pts[i + dim][cnt], 1, &coeffs[offset],
170  1);
171 
172  m_f->m_exp[i]->GetExp(e)->EquiSpacedToCoeffs(
173  coeffs + offset, tmp = coeffs + offset);
174  cnt += ncoeffs;
175  }
176  }
177  }
178  else
179  {
180  int totpoints = m_f->m_exp[0]->GetTotPoints();
182 
183  coords[0] = Array<OneD, NekDouble>(totpoints);
184  coords[1] = Array<OneD, NekDouble>(totpoints);
185  coords[2] = Array<OneD, NekDouble>(totpoints);
186 
187  m_f->m_exp[0]->GetCoords(coords[0], coords[1], coords[2]);
188 
189  if (pts[0].num_elements() != totpoints)
190  {
191  WARNINGL0(false, "length of points in .pts file is different to "
192  "the number of quadrature points in xml file");
193  totpoints = min(totpoints, (int)pts[0].num_elements());
194  }
195 
196  int mismatch = 0;
197  for (i = 0; i < totpoints; ++i)
198  {
199  for (j = 0; j < dim; ++j)
200  {
201  if (fabs(coords[j][i] - pts[j][i]) > 1e-4)
202  {
203  string outstring =
204  "Coordinates do not match within 1e-4: " +
205  boost::lexical_cast<string>(coords[j][i]) + " versus " +
206  boost::lexical_cast<string>(pts[j][i]) + " diff " +
207  boost::lexical_cast<string>(
208  fabs(coords[j][i] - pts[j][i]));
209  ;
210  WARNINGL0(false, outstring);
211  mismatch += 1;
212  }
213  }
214 
215  for (j = 0; j < nFields; ++j)
216  {
217  m_f->m_exp[j]->SetPhys(i, pts[j + dim][i]);
218  }
219  }
220 
221  if (m_f->m_session->GetComm()->GetRank() == 0)
222  {
223  cout << endl;
224  }
225 
226  // forward transform fields
227  for (i = 0; i < nFields; ++i)
228  {
229  m_f->m_exp[i]->FwdTrans_IterPerExp(m_f->m_exp[i]->GetPhys(),
230  m_f->m_exp[i]->UpdateCoeffs());
231  }
232  }
233 
234  // set up output fld file.
235  std::vector<LibUtilities::FieldDefinitionsSharedPtr> FieldDef =
236  m_f->m_exp[0]->GetFieldDefinitions();
237  std::vector<std::vector<NekDouble> > FieldData(FieldDef.size());
238 
239  for (j = 0; j < nFields; ++j)
240  {
241  for (i = 0; i < FieldDef.size(); ++i)
242  {
243  FieldDef[i]->m_fields.push_back(m_f->m_fieldPts->GetFieldName(j));
244 
245  m_f->m_exp[j]->AppendFieldData(FieldDef[i], FieldData[i]);
246  }
247  }
248 
249  m_f->m_fielddef = FieldDef;
250  m_f->m_data = FieldData;
251 }
252 }
253 }
map< string, ConfigOption > m_config
List of configuration values.
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
Represents a command-line configuration option.
virtual void Process(po::variables_map &vm)
Write mesh to output file.
STL namespace.
pair< ModuleType, string > ModuleKey
boost::shared_ptr< Field > FieldSharedPtr
Definition: Field.hpp:767
#define WARNINGL0(condition, msg)
Definition: ErrorUtil.hpp:204
double NekDouble
static PtsFieldSharedPtr NullPtsField
Definition: PtsField.h:179
Abstract base class for processing modules.
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1061
ModuleFactory & GetModuleFactory()
FieldSharedPtr m_f
Field object.
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, tDescription pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:215