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 "ProcessInterpField.h"
40 
45 #include <boost/math/special_functions/fpclassify.hpp>
46 namespace Nektar
47 {
48 namespace FieldUtils
49 {
50 
51 ModuleKey ProcessInterpField::className =
53  ModuleKey(eProcessModule, "interpfield"),
54  ProcessInterpField::create,
55  "Interpolates one field to another, requires fromxml, "
56  "fromfld to be defined");
57 
58 ProcessInterpField::ProcessInterpField(FieldSharedPtr f) : ProcessModule(f)
59 {
60 
61  m_config["fromxml"] = ConfigOption(
62  false, "NotSet", "Xml file form which to interpolate field");
63  m_config["fromfld"] = ConfigOption(
64  false, "NotSet", "Fld file form which to interpolate field");
65 
66  m_config["clamptolowervalue"] =
67  ConfigOption(false, "-10000000", "Lower bound for interpolation value");
68  m_config["clamptouppervalue"] =
69  ConfigOption(false, "10000000", "Upper bound for interpolation value");
70  m_config["defaultvalue"] =
71  ConfigOption(false, "0", "Default value if point is outside domain");
72 }
73 
75 {
76 }
77 
78 void ProcessInterpField::Process(po::variables_map &vm)
79 {
80  if (m_f->m_verbose)
81  {
82  if (m_f->m_comm->TreatAsRankZero())
83  {
84  cout << "ProcessInterpField: Interpolating field..." << endl;
85  }
86  }
87 
88  m_fromField = boost::shared_ptr<Field>(new Field());
89 
90  std::vector<std::string> files;
91 
92  // set up session file for from field
93  ParseUtils::GenerateOrderedStringVector(m_config["fromxml"].as<string>().c_str(), files);
94  m_fromField->m_session =
96 
97  // Set up range based on min and max of local parallel partition
100 
101  int coordim = m_f->m_exp[0]->GetCoordim(0);
102  int npts = m_f->m_exp[0]->GetTotPoints();
104 
105  for (int i = 0; i < coordim; ++i)
106  {
107  coords[i] = Array<OneD, NekDouble>(npts);
108  }
109 
110  for (int i = coordim; i < 3; ++i)
111  {
112  coords[i] = NullNekDouble1DArray;
113  }
114 
115  m_f->m_exp[0]->GetCoords(coords[0], coords[1], coords[2]);
116 
117  rng->m_checkShape = false;
118  switch (coordim)
119  {
120  case 3:
121  rng->m_doZrange = true;
122  rng->m_zmin = Vmath::Vmin(npts, coords[2], 1);
123  rng->m_zmax = Vmath::Vmax(npts, coords[2], 1);
124  case 2:
125  rng->m_doYrange = true;
126  rng->m_ymin = Vmath::Vmin(npts, coords[1], 1);
127  rng->m_ymax = Vmath::Vmax(npts, coords[1], 1);
128  case 1:
129  rng->m_doXrange = true;
130  rng->m_xmin = Vmath::Vmin(npts, coords[0], 1);
131  rng->m_xmax = Vmath::Vmax(npts, coords[0], 1);
132  break;
133  default:
134  ASSERTL0(false, "too many values specfied in range");
135  }
136 
137  // setup rng parameters.
138  m_fromField->m_graph =
140 
141  // Read in local from field partitions
142  const SpatialDomains::ExpansionMap &expansions =
143  m_fromField->m_graph->GetExpansions();
144 
145  // check for case where no elements are specified on this
146  // parallel partition
147  if (!expansions.size())
148  {
149  return;
150  }
151 
152  Array<OneD, int> ElementGIDs(expansions.size());
153  SpatialDomains::ExpansionMap::const_iterator expIt;
154 
155  int i = 0;
156  for (expIt = expansions.begin(); expIt != expansions.end(); ++expIt)
157  {
158  ElementGIDs[i++] = expIt->second->m_geomShPtr->GetGlobalID();
159  }
160 
161  string fromfld = m_config["fromfld"].as<string>();
162  m_f->FieldIOForFile(fromfld)->Import(
163  fromfld, m_fromField->m_fielddef, m_fromField->m_data,
165 
166  int NumHomogeneousDir = m_fromField->m_fielddef[0]->m_numHomogeneousDir;
167 
168  //----------------------------------------------
169  // Set up Expansion information to use mode order from field
170  m_fromField->m_graph->SetExpansions(m_fromField->m_fielddef);
171 
172  int nfields = m_fromField->m_fielddef[0]->m_fields.size();
173 
174  m_fromField->m_exp.resize(nfields);
175  m_fromField->m_exp[0] =
176  m_fromField->SetUpFirstExpList(NumHomogeneousDir, true);
177 
178  m_f->m_exp.resize(nfields);
179 
180  // declare auxiliary fields.
181  for (i = 1; i < nfields; ++i)
182  {
183  m_f->m_exp[i] = m_f->AppendExpList(NumHomogeneousDir);
184  m_fromField->m_exp[i] = m_fromField->AppendExpList(NumHomogeneousDir);
185  }
186 
187  // load field into expansion in fromfield.
188  for (int j = 0; j < nfields; ++j)
189  {
190  for (i = 0; i < m_fromField->m_fielddef.size(); i++)
191  {
192  m_fromField->m_exp[j]->ExtractDataToCoeffs(
193  m_fromField->m_fielddef[i], m_fromField->m_data[i],
194  m_fromField->m_fielddef[0]->m_fields[j],
195  m_fromField->m_exp[j]->UpdateCoeffs());
196  }
197  m_fromField->m_exp[j]->BwdTrans(m_fromField->m_exp[j]->GetCoeffs(),
198  m_fromField->m_exp[j]->UpdatePhys());
199  }
200 
201  int nq1 = m_f->m_exp[0]->GetTotPoints();
202 
203  Array<OneD, NekDouble> x1(nq1);
204  Array<OneD, NekDouble> y1(nq1);
205  Array<OneD, NekDouble> z1(nq1);
206 
207  if (coordim == 2)
208  {
209  m_f->m_exp[0]->GetCoords(x1, y1);
210  }
211  else if (coordim == 3)
212  {
213  m_f->m_exp[0]->GetCoords(x1, y1, z1);
214  }
215 
216  NekDouble clamp_low = m_config["clamptolowervalue"].as<NekDouble>();
217  NekDouble clamp_up = m_config["clamptouppervalue"].as<NekDouble>();
218  NekDouble def_value = m_config["defaultvalue"].as<NekDouble>();
219 
220  for (int i = 0; i < nfields; i++)
221  {
222  for (int j = 0; j < nq1; ++j)
223  {
224  m_f->m_exp[i]->UpdatePhys()[j] = def_value;
225  }
226  }
227 
228  Interpolator interp;
229  if (m_f->m_comm->GetRank() == 0)
230  {
232  }
233  interp.Interpolate(m_fromField->m_exp, m_f->m_exp);
234  if (m_f->m_comm->GetRank() == 0)
235  {
236  cout << endl;
237  }
238 
239  for (int i = 0; i < nfields; ++i)
240  {
241  for (int j = 0; j < nq1; ++j)
242  {
243  if (m_f->m_exp[i]->GetPhys()[j] > clamp_up)
244  {
245  m_f->m_exp[i]->UpdatePhys()[j] = clamp_up;
246  }
247  else if (m_f->m_exp[i]->GetPhys()[j] < clamp_low)
248  {
249  m_f->m_exp[i]->UpdatePhys()[j] = clamp_low;
250  }
251  }
252  }
253 
254  // put field into field data for output
255  std::vector<LibUtilities::FieldDefinitionsSharedPtr> FieldDef =
256  m_f->m_exp[0]->GetFieldDefinitions();
257  std::vector<std::vector<NekDouble> > FieldData(FieldDef.size());
258 
259  for (int j = 0; j < nfields; ++j)
260  {
261  m_f->m_exp[j]->FwdTrans(m_f->m_exp[j]->GetPhys(),
262  m_f->m_exp[j]->UpdateCoeffs());
263  for (i = 0; i < FieldDef.size(); ++i)
264  {
265  FieldDef[i]->m_fields.push_back(
266  m_fromField->m_fielddef[0]->m_fields[j]);
267  m_f->m_exp[j]->AppendFieldData(FieldDef[i], FieldData[i]);
268  }
269  }
270 
271  m_f->m_fielddef = FieldDef;
272  m_f->m_data = FieldData;
273 }
274 
275 void ProcessInterpField::PrintProgressbar(const int position,
276  const int goal) const
277 {
278  LibUtilities::PrintProgressbar(position, goal, "Interpolating");
279 }
280 }
281 }
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:78
void SetProgressCallback(FuncPointerT func, ObjectPointerT obj)
sets a callback funtion which gets called every time the interpolation progresses ...
Definition: Interpolator.h:159
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:740
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