Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties 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 <iostream>
36 #include <string>
37 using namespace std;
38 
39 #include <boost/geometry.hpp>
40 #include "ProcessInterpField.h"
41 
46 #include <boost/math/special_functions/fpclassify.hpp>
47 
48 namespace bg = boost::geometry;
49 namespace bgi = boost::geometry::index;
50 
51 namespace Nektar
52 {
53 namespace FieldUtils
54 {
55 
56 ModuleKey ProcessInterpField::className =
58  ModuleKey(eProcessModule, "interpfield"),
59  ProcessInterpField::create,
60  "Interpolates one field to another, requires fromxml, "
61  "fromfld to be defined");
62 
63 ProcessInterpField::ProcessInterpField(FieldSharedPtr f) : ProcessModule(f)
64 {
65 
66  m_config["fromxml"] = ConfigOption(
67  false, "NotSet", "Xml file form which to interpolate field");
68  m_config["fromfld"] = ConfigOption(
69  false, "NotSet", "Fld file form which to interpolate field");
70 
71  m_config["clamptolowervalue"] =
72  ConfigOption(false, "-10000000", "Lower bound for interpolation value");
73  m_config["clamptouppervalue"] =
74  ConfigOption(false, "10000000", "Upper bound for interpolation value");
75  m_config["defaultvalue"] =
76  ConfigOption(false, "0", "Default value if point is outside domain");
77 }
78 
80 {
81 }
82 
83 void ProcessInterpField::Process(po::variables_map &vm)
84 {
85  if (m_f->m_verbose)
86  {
87  if (m_f->m_comm->TreatAsRankZero())
88  {
89  cout << "ProcessInterpField: Interpolating field..." << endl;
90  }
91  }
92 
93  m_fromField = boost::shared_ptr<Field>(new Field());
94 
95  std::vector<std::string> files;
96 
97  // set up session file for from field
98  ParseUtils::GenerateOrderedStringVector(m_config["fromxml"].as<string>().c_str(), files);
99  m_fromField->m_session =
101 
102  // Set up range based on min and max of local parallel partition
105 
106  int coordim = m_f->m_exp[0]->GetCoordim(0);
107  int npts = m_f->m_exp[0]->GetTotPoints();
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  rng->m_checkShape = false;
123  switch (coordim)
124  {
125  case 3:
126  rng->m_doZrange = true;
127  rng->m_zmin = Vmath::Vmin(npts, coords[2], 1);
128  rng->m_zmax = Vmath::Vmax(npts, coords[2], 1);
129  case 2:
130  rng->m_doYrange = true;
131  rng->m_ymin = Vmath::Vmin(npts, coords[1], 1);
132  rng->m_ymax = Vmath::Vmax(npts, coords[1], 1);
133  case 1:
134  rng->m_doXrange = true;
135  rng->m_xmin = Vmath::Vmin(npts, coords[0], 1);
136  rng->m_xmax = Vmath::Vmax(npts, coords[0], 1);
137  break;
138  default:
139  ASSERTL0(false, "too many values specfied in range");
140  }
141 
142  // setup rng parameters.
143  m_fromField->m_graph =
145 
146  // Read in local from field partitions
147  const SpatialDomains::ExpansionMap &expansions =
148  m_fromField->m_graph->GetExpansions();
149 
150  // check for case where no elements are specified on this
151  // parallel partition
152  if (!expansions.size())
153  {
154  return;
155  }
156 
157  Array<OneD, int> ElementGIDs(expansions.size());
158  SpatialDomains::ExpansionMap::const_iterator expIt;
159 
160  int i = 0;
161  for (expIt = expansions.begin(); expIt != expansions.end(); ++expIt)
162  {
163  ElementGIDs[i++] = expIt->second->m_geomShPtr->GetGlobalID();
164  }
165 
166  string fromfld = m_config["fromfld"].as<string>();
167  m_f->FieldIOForFile(fromfld)->Import(
168  fromfld, m_fromField->m_fielddef, m_fromField->m_data,
170 
171  int NumHomogeneousDir = m_fromField->m_fielddef[0]->m_numHomogeneousDir;
172 
173  //----------------------------------------------
174  // Set up Expansion information to use mode order from field
175  m_fromField->m_graph->SetExpansions(m_fromField->m_fielddef);
176 
177  int nfields = m_fromField->m_fielddef[0]->m_fields.size();
178 
179  m_fromField->m_exp.resize(nfields);
180  m_fromField->m_exp[0] =
181  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  // load field into expansion in fromfield.
193  for (int j = 0; j < nfields; ++j)
194  {
195  for (i = 0; i < m_fromField->m_fielddef.size(); i++)
196  {
197  m_fromField->m_exp[j]->ExtractDataToCoeffs(
198  m_fromField->m_fielddef[i], m_fromField->m_data[i],
199  m_fromField->m_fielddef[0]->m_fields[j],
200  m_fromField->m_exp[j]->UpdateCoeffs());
201  }
202  m_fromField->m_exp[j]->BwdTrans(m_fromField->m_exp[j]->GetCoeffs(),
203  m_fromField->m_exp[j]->UpdatePhys());
204  }
205 
206  int nq1 = m_f->m_exp[0]->GetTotPoints();
207 
208  Array<OneD, NekDouble> x1(nq1);
209  Array<OneD, NekDouble> y1(nq1);
210  Array<OneD, NekDouble> z1(nq1);
211 
212  if (coordim == 2)
213  {
214  m_f->m_exp[0]->GetCoords(x1, y1);
215  }
216  else if (coordim == 3)
217  {
218  m_f->m_exp[0]->GetCoords(x1, y1, z1);
219  }
220 
221  NekDouble clamp_low = m_config["clamptolowervalue"].as<NekDouble>();
222  NekDouble clamp_up = m_config["clamptouppervalue"].as<NekDouble>();
223  NekDouble def_value = m_config["defaultvalue"].as<NekDouble>();
224 
225  for (int i = 0; i < nfields; i++)
226  {
227  for (int j = 0; j < nq1; ++j)
228  {
229  m_f->m_exp[i]->UpdatePhys()[j] = def_value;
230  }
231  }
232 
233  Interpolator interp;
234  if (m_f->m_comm->GetRank() == 0)
235  {
237  }
238  interp.Interpolate(m_fromField->m_exp, m_f->m_exp);
239  if (m_f->m_comm->GetRank() == 0)
240  {
241  cout << endl;
242  }
243 
244  for (int i = 0; i < nfields; ++i)
245  {
246  for (int j = 0; j < nq1; ++j)
247  {
248  if (m_f->m_exp[i]->GetPhys()[j] > clamp_up)
249  {
250  m_f->m_exp[i]->UpdatePhys()[j] = clamp_up;
251  }
252  else if (m_f->m_exp[i]->GetPhys()[j] < clamp_low)
253  {
254  m_f->m_exp[i]->UpdatePhys()[j] = clamp_low;
255  }
256  }
257  }
258 
259  // put field into field data for output
260  std::vector<LibUtilities::FieldDefinitionsSharedPtr> FieldDef =
261  m_f->m_exp[0]->GetFieldDefinitions();
262  std::vector<std::vector<NekDouble> > FieldData(FieldDef.size());
263 
264  for (int j = 0; j < nfields; ++j)
265  {
266  m_f->m_exp[j]->FwdTrans(m_f->m_exp[j]->GetPhys(),
267  m_f->m_exp[j]->UpdateCoeffs());
268  for (i = 0; i < FieldDef.size(); ++i)
269  {
270  FieldDef[i]->m_fields.push_back(
271  m_fromField->m_fielddef[0]->m_fields[j]);
272  m_f->m_exp[j]->AppendFieldData(FieldDef[i], FieldData[i]);
273  }
274  }
275 
276  m_f->m_fielddef = FieldDef;
277  m_f->m_data = FieldData;
278 }
279 
280 void ProcessInterpField::PrintProgressbar(const int position,
281  const int goal) const
282 {
283  LibUtilities::PrintProgressbar(position, goal, "Interpolating");
284 }
285 }
286 }
map< string, ConfigOption > m_config
List of configuration values.
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
static bool GenerateOrderedStringVector(const char *const str, std::vector< std::string > &vec)
Definition: ParseUtils.hpp:143
A class that contains algorithms for interpolation between pts fields, expansions and different meshe...
Definition: Interpolator.h:75
void SetProgressCallback(FuncPointerT func, ObjectPointerT obj)
sets a callback funtion which gets called every time the interpolation progresses ...
Definition: Interpolator.h:156
static Array< OneD, NekDouble > NullNekDouble1DArray
int PrintProgressbar(const int position, const int goal, const string message, int lastprogress=-1)
Prints a progressbar.
Definition: Progressbar.hpp:69
static boost::shared_ptr< MeshGraph > Read(const LibUtilities::SessionReaderSharedPtr &pSession, DomainRangeShPtr &rng=NullDomainRangeShPtr)
Definition: MeshGraph.cpp:124
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
Represents a command-line configuration option.
FIELD_UTILS_EXPORT void Interpolate(const LibUtilities::PtsFieldSharedPtr ptsInField, LibUtilities::PtsFieldSharedPtr &ptsOutField)
Interpolate from a pts field to a pts field.
T Vmax(int n, const T *x, const int incx)
Return the maximum element in x – called vmax to avoid conflict with max.
Definition: Vmath.cpp:779
T Vmin(int n, const T *x, const int incx)
Return the minimum element in x - called vmin to avoid conflict with min.
Definition: Vmath.cpp:871
STL namespace.
pair< ModuleType, string > ModuleKey
void PrintProgressbar(const int position, const int goal) const
static SessionReaderSharedPtr CreateInstance(int argc, char *argv[])
Creates an instance of the SessionReader class.
boost::shared_ptr< Field > FieldSharedPtr
Definition: Field.hpp:767
static std::string npts
Definition: InputFld.cpp:43
boost::shared_ptr< DomainRange > DomainRangeShPtr
Definition: MeshGraph.h:157
double NekDouble
virtual void Process(po::variables_map &vm)
Write mesh to output file.
Abstract base class for processing modules.
static FieldMetaDataMap NullFieldMetaDataMap
Definition: FieldIO.h:55
std::map< int, ExpansionShPtr > ExpansionMap
Definition: MeshGraph.h:174
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