Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
ProcessInterpField.cpp
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////////////
2 //
3 // File: ProcessInterpField.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: Interpolate one field to another.
33 //
34 ////////////////////////////////////////////////////////////////////////////////
35 #include <string>
36 #include <iostream>
37 using namespace std;
38 
39 #include "ProcessInterpField.h"
40 
43 #include <boost/math/special_functions/fpclassify.hpp>
44 namespace Nektar
45 {
46 namespace Utilities
47 {
48 
49 ModuleKey ProcessInterpField::className =
51  ModuleKey(eProcessModule, "interpfield"),
52  ProcessInterpField::create,
53  "Interpolates one field to another, requires fromxml, "
54  "fromfld to be defined");
55 
56 
57 ProcessInterpField::ProcessInterpField(FieldSharedPtr f) : ProcessModule(f)
58 {
59 
60  m_config["fromxml"] = ConfigOption(false, "NotSet",
61  "Xml file form which to interpolate field");
62 
63  ASSERTL0(m_config["fromxml"].as<string>().compare("NotSet") != 0,
64  "Need to specify fromxml=file.xml");
65 
66  m_config["fromfld"] = ConfigOption(false, "NotSet",
67  "Fld file form which to interpolate field");
68 
69  ASSERTL0(m_config["fromfld"].as<string>().compare("NotSet") != 0,
70  "Need to specify fromfld=file.fld ");
71 
72  m_config["clamptolowervalue"] = ConfigOption(false,"-10000000",
73  "Lower bound for interpolation value");
74  m_config["clamptouppervalue"] = ConfigOption(false,"10000000",
75  "Upper bound for interpolation value");
76  m_config["defaultvalue"] = ConfigOption(false,"0",
77  "Default value if point is outside domain");
78 }
79 
81 {
82 }
83 
84 void ProcessInterpField::Process(po::variables_map &vm)
85 {
86 
87  if(m_f->m_verbose)
88  {
89  cout << "Processing interpolation" << endl;
90  }
91 
92  m_fromField = boost::shared_ptr<Field>(new Field());
93 
94  std::vector<std::string> files;
95 
96  // set up session file for from field
97  files.push_back(m_config["fromxml"].as<string>());
99  CreateInstance(0, 0, files);
100 
101  // Set up range based on min and max of local parallel partition
104 
105  int coordim = m_f->m_exp[0]->GetCoordim(0);
106  int npts = m_f->m_exp[0]->GetTotPoints();
107  Array<OneD, Array<OneD, NekDouble> > coords(3);
108 
109 
110  for(int i = 0; i < coordim; ++i)
111  {
112  coords[i] = Array<OneD, NekDouble>(npts);
113  }
114 
115  for(int i = coordim; i < 3; ++i)
116  {
117  coords[i] = NullNekDouble1DArray;
118  }
119 
120  m_f->m_exp[0]->GetCoords(coords[0],coords[1],coords[2]);
121 
122  switch(coordim)
123  {
124  case 3:
125  rng->doZrange = true;
126  rng->zmin = Vmath::Vmin(npts,coords[2],1);
127  rng->zmax = Vmath::Vmax(npts,coords[2],1);
128  case 2:
129  rng->doYrange = true;
130  rng->ymin = Vmath::Vmin(npts,coords[1],1);
131  rng->ymax = Vmath::Vmax(npts,coords[1],1);
132  case 1:
133  rng->doXrange = true;
134  rng->xmin = Vmath::Vmin(npts,coords[0],1);
135  rng->xmax = Vmath::Vmax(npts,coords[0],1);
136  break;
137  default:
138  ASSERTL0(false,"too many values specfied in range");
139  }
140 
141  // setup rng parameters.
142  m_fromField->m_graph =
144 
145  // Read in local from field partitions
146  const SpatialDomains::ExpansionMap &expansions =
147  m_fromField->m_graph->GetExpansions();
148 
149  // check for case where no elements are specified on this
150  // parallel partition
151  if(!expansions.size())
152  {
153  return;
154  }
155 
156  Array<OneD,int> ElementGIDs(expansions.size());
157  SpatialDomains::ExpansionMap::const_iterator expIt;
158 
159  int i = 0;
160  for (expIt = expansions.begin(); expIt != expansions.end();
161  ++expIt)
162  {
163  ElementGIDs[i++] = expIt->second->m_geomShPtr->GetGlobalID();
164  }
165 
166  string fromfld = m_config["fromfld"].as<string>();
167  m_f->m_fld->Import(fromfld,m_fromField->m_fielddef,
168  m_fromField->m_data,
170  ElementGIDs);
171 
172  int NumHomogeneousDir = m_fromField->m_fielddef[0]->m_numHomogeneousDir;
173 
174  //----------------------------------------------
175  // Set up Expansion information to use mode order from field
176  m_fromField->m_graph->SetExpansions(m_fromField->m_fielddef);
177 
178  int nfields = m_fromField->m_fielddef[0]->m_fields.size();
179 
180  m_fromField->m_exp.resize(nfields);
181  m_fromField->m_exp[0] = m_fromField->SetUpFirstExpList(NumHomogeneousDir,true);
182 
183  m_f->m_exp.resize(nfields);
184 
185  // declare auxiliary fields.
186  for(i = 1; i < nfields; ++i)
187  {
188  m_f->m_exp[i] = m_f->AppendExpList(NumHomogeneousDir);
189  m_fromField->m_exp[i] = m_fromField->AppendExpList(NumHomogeneousDir);
190  }
191 
192 
193  // load field into expansion in fromfield.
194  for(int j = 0; j < nfields; ++j)
195  {
196  for (i = 0; i < m_fromField->m_fielddef.size(); i++)
197  {
198  m_fromField->m_exp[j]->ExtractDataToCoeffs(
199  m_fromField->m_fielddef[i],
200  m_fromField->m_data[i],
201  m_fromField->m_fielddef[0]->m_fields[j],
202  m_fromField->m_exp[j]->UpdateCoeffs());
203  }
204  m_fromField->m_exp[j]->BwdTrans(m_fromField->m_exp[j]->GetCoeffs(),
205  m_fromField->m_exp[j]->UpdatePhys());
206 
207  }
208 
209  int nq1 = m_f->m_exp[0]->GetTotPoints();
210 
211  Array<OneD, NekDouble> x1(nq1);
212  Array<OneD, NekDouble> y1(nq1);
213  Array<OneD, NekDouble> z1(nq1);
214 
215  if (coordim == 2)
216  {
217  m_f->m_exp[0]->GetCoords(x1, y1);
218  }
219  else if (coordim == 3)
220  {
221  m_f->m_exp[0]->GetCoords(x1, y1, z1);
222  }
223 
224  if(m_f->m_session->GetComm()->TreatAsRankZero())
225  {
226  cout << "Interpolating [" << flush;
227  }
228 
229  NekDouble clamp_low = m_config["clamptolowervalue"].as<NekDouble>();
230  NekDouble clamp_up = m_config["clamptouppervalue"].as<NekDouble>();
231  NekDouble def_value = m_config["defaultvalue"].as<NekDouble>();
232 
233  InterpolateField(m_fromField->m_exp, m_f->m_exp,
234  x1, y1, z1, clamp_low, clamp_up,def_value);
235 
236  if(m_f->m_session->GetComm()->TreatAsRankZero())
237  {
238  cout << "]" << endl;
239  }
240 
241 
242  // put field into field data for output
243  std::vector<LibUtilities::FieldDefinitionsSharedPtr> FieldDef
244  = m_f->m_exp[0]->GetFieldDefinitions();
245  std::vector<std::vector<NekDouble> > FieldData(FieldDef.size());
246 
247  for (int j = 0; j < nfields; ++j)
248  {
249  m_f->m_exp[j]->FwdTrans(m_f->m_exp[j]->GetPhys(),
250  m_f->m_exp[j]->UpdateCoeffs());
251  for (i = 0; i < FieldDef.size(); ++i)
252  {
253  FieldDef[i]->m_fields.push_back(m_fromField->m_fielddef[0]->m_fields[j]);
254  m_f->m_exp[j]->AppendFieldData(FieldDef[i], FieldData[i]);
255  }
256  }
257 
258  m_f->m_fielddef = FieldDef;
259  m_f->m_data = FieldData;
260 }
261 
263  vector<MultiRegions::ExpListSharedPtr> &field0,
264  vector<MultiRegions::ExpListSharedPtr> &field1,
265  Array<OneD, NekDouble> x,
266  Array<OneD, NekDouble> y,
267  Array<OneD, NekDouble> z,
268  NekDouble clamp_low,
269  NekDouble clamp_up,
270  NekDouble def_value)
271 {
272  int expdim = field0[0]->GetCoordim(0);
273 
274  Array<OneD, NekDouble> coords(expdim), Lcoords(expdim);
275  int nq1 = field1[0]->GetTotPoints();
276  int elmtid, offset;
277  int r, f;
278  static int intpts = 0;
279 
280  ASSERTL0(field0.size() == field1.size(),
281  "Input field dimension must be same as output dimension");
282 
283  for (r = 0; r < nq1; r++)
284  {
285  coords[0] = x[r];
286  coords[1] = y[r];
287  if (expdim == 3)
288  {
289  coords[2] = z[r];
290  }
291 
292  // Obtain nearest Element and LocalCoordinate to interpolate
293  elmtid = field0[0]->GetExpIndex(coords, Lcoords, 1e-3,true);
294 
295  if(elmtid >= 0)
296  {
297  offset = field0[0]->GetPhys_Offset(field0[0]->
298  GetOffset_Elmt_Id(elmtid));
299 
300  for (f = 0; f < field1.size(); ++f)
301  {
302  NekDouble value;
303  value = field0[f]->GetExp(elmtid)->
304  StdPhysEvaluate(Lcoords, field0[f]->GetPhys() +offset);
305 
306  if ((boost::math::isnan)(value))
307  {
308  ASSERTL0(false, "new value is not a number");
309  }
310  else
311  {
312  value = (value > clamp_up)? clamp_up :
313  ((value < clamp_low)? clamp_low :
314  value);
315 
316  field1[f]->UpdatePhys()[r] = value;
317  }
318  }
319  }
320  else
321  {
322  for (f = 0; f < field1.size(); ++f)
323  {
324  field1[f]->UpdatePhys()[r] = def_value;
325  }
326  }
327 
328  if (intpts%1000 == 0)
329  {
330  cout <<"." << flush;
331  }
332  intpts ++;
333  }
334 }
335 
336 }
337 }
338 
339