Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator 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 <string>
37 #include <iostream>
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 Utilities
49 {
50 
51 ModuleKey ProcessPointDataToFld::className =
53  ModuleKey(eProcessModule, "pointdatatofld"),
54  ProcessPointDataToFld::create,
55  "Given discrete data at quadrature points project them onto an expansion"
56  "basis and output fld file. Requires .pts .xml and .fld files.");
57 
58 
59 ProcessPointDataToFld::ProcessPointDataToFld(FieldSharedPtr f)
60  : ProcessModule(f)
61 {
62  m_requireEquiSpaced = true;
63 
64  m_config["setnantovalue"] = ConfigOption(false,"NotSet",
65  "reset any nan value to prescribed value");
66 
67  if((f->m_inputfiles.count("pts") == 0))
68  {
69  cout << endl << "A pts input file must be specified for the boundary extraction module" << endl;
70 
71  cout << "Usage: Fieldconvert -m pointdatatofld file.pts file.xml out.fld" << endl;
72  exit(3);
73  }
74 
75  if((f->m_inputfiles.count("xml") == 0)&&(f->m_inputfiles.count("xml.gz") == 0))
76  {
77  cout << endl << "An xml or xml.gz input file must be specified for the boundary extraction module" << endl;
78  cout << "Usage: Fieldconvert -m pointdatatofld file.pts file.xml out.fld" << endl;
79  exit(3);
80  }
81 
82 }
83 
85 {
86 }
87 
88 void ProcessPointDataToFld::Process(po::variables_map &vm)
89 {
90  if (m_f->m_verbose)
91  {
92  if(m_f->m_comm->GetRank() == 0)
93  {
94  cout << "ProcessPointDataToFld: projecting data to expansion..."
95  << endl;
96  }
97  }
98 
99  int i,j;
100  bool setnantovalue = false;
101  NekDouble defvalue;
102 
103  if(!boost::iequals(m_config["setnantovalue"].as<string>(),"NotSet"))
104  {
105  setnantovalue = true;
106  defvalue = m_config["setnantovalue"].as<NekDouble>();
107  }
108 
109  // Check for command line point specification if no .pts file specified
110  ASSERTL0(m_f->m_fieldPts != LibUtilities::NullPtsField,
111  "No input points found");
112 
113  int nFields = m_f->m_fieldPts->GetNFields();
114  ASSERTL0(nFields > 0,
115  "No field values provided in input");
116 
117  int dim = m_f->m_fieldPts->GetDim();
118 
119  // assume one field is already defined from input file.
120  m_f->m_exp.resize(nFields);
121  for(i = 1; i < nFields; ++i)
122  {
123  m_f->m_exp[i] = m_f->AppendExpList(0);
124  }
125 
127  m_f->m_fieldPts->GetPts(pts);
128 
129  // set any nan values to default value if requested
130  if(setnantovalue)
131  {
132  for(int i = 0; i < pts[0].num_elements(); ++i)
133  {
134  for(int j = 0; j < nFields; ++j)
135  {
136  if ((boost::math::isnan)(pts[j+dim][i]))
137  {
138  pts[j+dim][i] = defvalue;
139  }
140  }
141  }
142  }
143 
144  if(m_f->m_fieldPts->m_ptsInfo.count(LibUtilities::eIsEquiSpacedData) != 0)
145  {
146  int totcoeffs = m_f->m_exp[0]->GetNcoeffs();
147 
148  ASSERTL0(pts[0].num_elements() != totcoeffs,"length of points in .pts file is different "
149  "to the number of coefficients in expansion ");
150 
151  for(int i = 0; i < nFields; ++i)
152  {
153  Array<OneD, NekDouble> coeffs = m_f->m_exp[i]->UpdateCoeffs(), tmp;
154  int cnt = 0;
155  for(int e = 0; e < m_f->m_exp[0]->GetNumElmts(); ++e)
156  {
157  int ncoeffs = m_f->m_exp[i]->GetExp(e)->GetNcoeffs();
158  int offset = m_f->m_exp[i]->GetCoeff_Offset(e);
159  Vmath::Vcopy(ncoeffs,&pts[i+dim][cnt],1,&coeffs[offset],1);
160 
161  m_f->m_exp[i]->GetExp(e)->EquiSpacedToCoeffs(coeffs+offset,tmp = coeffs+offset);
162  cnt += ncoeffs;
163  }
164  }
165  }
166  else
167  {
168  int totpoints = m_f->m_exp[0]->GetTotPoints();
170 
171  coords[0] = Array<OneD,NekDouble>(totpoints);
172  coords[1] = Array<OneD,NekDouble>(totpoints);
173  coords[2] = Array<OneD,NekDouble>(totpoints);
174 
175  m_f->m_exp[0]->GetCoords(coords[0],coords[1],coords[2]);
176 
177 
178  if(pts[0].num_elements() != totpoints)
179  {
180  WARNINGL0(false,"length of points in .pts file is different to the number of quadrature points in xml file");
181  totpoints = min(totpoints,(int)pts[0].num_elements());
182  }
183 
184  int mismatch = 0;
185  for(i = 0; i < totpoints; ++i)
186  {
187  for(j = 0; j < dim; ++j)
188  {
189  if(fabs(coords[j][i]-pts[j][i]) > 1e-4)
190  {
191  string outstring = "Coordinates do not match within 1e-4: "
192  + boost::lexical_cast<string>(coords[j][i])
193  + " versus "
194  + boost::lexical_cast<string>(pts[j][i])
195  + " diff " + boost::lexical_cast<string>(fabs(coords[j][i]-pts[j][i]));;
196  WARNINGL0(false,outstring);
197  mismatch +=1;
198  }
199  }
200 
201  for(j = 0;j < nFields; ++j)
202  {
203  m_f->m_exp[j]->SetPhys(i, pts[j+dim][i]);
204  }
205  }
206 
207  if(m_f->m_session->GetComm()->GetRank() == 0)
208  {
209  cout << endl;
210  }
211 
212  // forward transform fields
213  for(i = 0; i < nFields; ++i)
214  {
215  m_f->m_exp[i]->FwdTrans_IterPerExp(m_f->m_exp[i]->GetPhys(),
216  m_f->m_exp[i]->UpdateCoeffs());
217  }
218  }
219 
220  // set up output fld file.
221  std::vector<LibUtilities::FieldDefinitionsSharedPtr> FieldDef
222  = m_f->m_exp[0]->GetFieldDefinitions();
223  std::vector<std::vector<NekDouble> > FieldData(FieldDef.size());
224 
225  for (j = 0; j < nFields; ++j)
226  {
227  for (i = 0; i < FieldDef.size(); ++i)
228  {
229  FieldDef[i]->m_fields.push_back(m_f->m_fieldPts->GetFieldName(j));
230 
231  m_f->m_exp[j]->AppendFieldData(FieldDef[i], FieldData[i]);
232  }
233  }
234 
235  m_f->m_fielddef = FieldDef;
236  m_f->m_data = FieldData;
237 }
238 
239 }
240 }
241 
242 
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:188
pair< ModuleType, string > ModuleKey
virtual void Process()=0
map< string, ConfigOption > m_config
List of configuration values.
STL namespace.
FieldSharedPtr m_f
Field object.
#define WARNINGL0(condition, msg)
Definition: ErrorUtil.hpp:194
double NekDouble
boost::shared_ptr< Field > FieldSharedPtr
Definition: Field.hpp:695
Represents a command-line configuration option.
static PtsFieldSharedPtr NullPtsField
Definition: PtsField.h:263
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1047
ModuleFactory & GetModuleFactory()
Abstract base class for processing modules.
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, tDescription pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:215