Nektar++
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 // Permission is hereby granted, free of charge, to any person obtaining a
14 // copy of this software and associated documentation files (the "Software"),
15 // to deal in the Software without restriction, including without limitation
16 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
17 // and/or sell copies of the Software, and to permit persons to whom the
18 // Software is furnished to do so, subject to the following conditions:
19 //
20 // The above copyright notice and this permission notice shall be included
21 // in all copies or substantial portions of the Software.
22 //
23 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29 // DEALINGS IN THE SOFTWARE.
30 //
31 // Description: Interpolate one field to another.
32 //
33 ////////////////////////////////////////////////////////////////////////////////
34 #include <iostream>
35 #include <string>
36 using namespace std;
37 
38 #include <boost/core/ignore_unused.hpp>
39 #include <boost/geometry.hpp>
40 #include <boost/math/special_functions/fpclassify.hpp>
41 
46 
47 #include "ProcessInterpField.h"
48 
49 namespace bg = boost::geometry;
50 namespace bgi = boost::geometry::index;
51 
52 namespace Nektar
53 {
54 namespace FieldUtils
55 {
56 
57 ModuleKey ProcessInterpField::className =
59  ModuleKey(eProcessModule, "interpfield"),
60  ProcessInterpField::create,
61  "Interpolates one field to another, requires fromxml, "
62  "fromfld to be defined");
63 
64 ProcessInterpField::ProcessInterpField(FieldSharedPtr f) : ProcessModule(f)
65 {
66 
67  m_config["fromxml"] = ConfigOption(
68  false, "NotSet", "Xml file from which to interpolate field");
69  m_config["fromfld"] = ConfigOption(
70  false, "NotSet", "Fld file from which to interpolate field");
71 
72  m_config["clamptolowervalue"] =
73  ConfigOption(false, "-10000000", "Lower bound for interpolation value");
74  m_config["clamptouppervalue"] =
75  ConfigOption(false, "10000000", "Upper bound for interpolation value");
76  m_config["defaultvalue"] =
77  ConfigOption(false, "0", "Default value if point is outside domain");
78 }
79 
81 {
82 }
83 
84 void ProcessInterpField::Process(po::variables_map &vm)
85 {
86  boost::ignore_unused(vm);
87 
88  FieldSharedPtr fromField = std::shared_ptr<Field>(new Field());
89 
90  std::vector<std::string> files;
91 
92  // set up session file for from field
93  char *argv[] = { const_cast<char *>("FieldConvert"), nullptr };
94  ParseUtils::GenerateVector(m_config["fromxml"].as<string>(), files);
95  fromField->m_session =
97  1, argv, files,
98  LibUtilities::GetCommFactory().CreateInstance("Serial", 0, 0));
99 
100  // Set up range based on min and max of local parallel partition
103 
104  int coordim = m_f->m_exp[0]->GetCoordim(0);
105  int npts = m_f->m_exp[0]->GetTotPoints();
107 
108  for (int i = 0; i < coordim; ++i)
109  {
110  coords[i] = Array<OneD, NekDouble>(npts);
111  }
112 
113  for (int i = coordim; i < 3; ++i)
114  {
115  coords[i] = NullNekDouble1DArray;
116  }
117 
118  m_f->m_exp[0]->GetCoords(coords[0], coords[1], coords[2]);
119 
120  rng->m_checkShape = false;
121  switch (coordim)
122  {
123  case 3:
124  rng->m_doZrange = true;
125  rng->m_zmin = Vmath::Vmin(npts, coords[2], 1);
126  rng->m_zmax = Vmath::Vmax(npts, coords[2], 1);
127  /* Falls through. */
128  case 2:
129  rng->m_doYrange = true;
130  rng->m_ymin = Vmath::Vmin(npts, coords[1], 1);
131  rng->m_ymax = Vmath::Vmax(npts, coords[1], 1);
132  /* Falls through. */
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  NEKERROR(ErrorUtil::efatal, "coordim should be <= 3");
140  }
141 
142  // setup rng parameters.
143  fromField->m_graph =
144  SpatialDomains::MeshGraph::Read(fromField->m_session, rng);
145 
146  // Read in local from field partitions
147  const SpatialDomains::ExpansionMap &expansions =
148  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 
159  int i = 0;
160  for (auto &expIt : expansions)
161  {
162  ElementGIDs[i++] = expIt.second->m_geomShPtr->GetGlobalID();
163  }
164 
165  string fromfld = m_config["fromfld"].as<string>();
166  m_f->FieldIOForFile(fromfld)->Import(
167  fromfld, fromField->m_fielddef, fromField->m_data,
169 
170  int NumHomogeneousDir = fromField->m_fielddef[0]->m_numHomogeneousDir;
171 
172  //----------------------------------------------
173  // Set up Expansion information to use mode order from field
174  fromField->m_graph->SetExpansions(fromField->m_fielddef);
175 
176  int nfields = fromField->m_fielddef[0]->m_fields.size();
177 
178  fromField->m_exp.resize(nfields);
179  fromField->m_exp[0] =
180  fromField->SetUpFirstExpList(NumHomogeneousDir, true);
181 
182  m_f->m_exp.resize(nfields);
183 
184  // declare auxiliary fields.
185  for (i = 1; i < nfields; ++i)
186  {
187  m_f->m_exp[i] = m_f->AppendExpList(NumHomogeneousDir);
188  fromField->m_exp[i] = fromField->AppendExpList(NumHomogeneousDir);
189  }
190 
191  // load field into expansion in fromfield.
192  for (int j = 0; j < nfields; ++j)
193  {
194  for (i = 0; i < fromField->m_fielddef.size(); i++)
195  {
196  fromField->m_exp[j]->ExtractDataToCoeffs(
197  fromField->m_fielddef[i], fromField->m_data[i],
198  fromField->m_fielddef[0]->m_fields[j],
199  fromField->m_exp[j]->UpdateCoeffs());
200  }
201  fromField->m_exp[j]->BwdTrans(fromField->m_exp[j]->GetCoeffs(),
202  fromField->m_exp[j]->UpdatePhys());
203  }
204 
205  int nq1 = m_f->m_exp[0]->GetTotPoints();
206 
207  NekDouble clamp_low = m_config["clamptolowervalue"].as<NekDouble>();
208  NekDouble clamp_up = m_config["clamptouppervalue"].as<NekDouble>();
209  NekDouble def_value = m_config["defaultvalue"].as<NekDouble>();
210 
211  for (int i = 0; i < nfields; i++)
212  {
213  for (int j = 0; j < nq1; ++j)
214  {
215  m_f->m_exp[i]->UpdatePhys()[j] = def_value;
216  }
217  }
218 
219  Interpolator interp;
220  if (m_f->m_verbose && m_f->m_comm->TreatAsRankZero())
221  {
223  }
224  interp.Interpolate(fromField->m_exp, m_f->m_exp);
225  if (m_f->m_verbose && m_f->m_comm->TreatAsRankZero())
226  {
227  cout << endl;
228  }
229 
230  for (int i = 0; i < nfields; ++i)
231  {
232  for (int j = 0; j < nq1; ++j)
233  {
234  if (m_f->m_exp[i]->GetPhys()[j] > clamp_up)
235  {
236  m_f->m_exp[i]->UpdatePhys()[j] = clamp_up;
237  }
238  else if (m_f->m_exp[i]->GetPhys()[j] < clamp_low)
239  {
240  m_f->m_exp[i]->UpdatePhys()[j] = clamp_low;
241  }
242  }
243  m_f->m_exp[i]->FwdTrans_IterPerExp(
244  m_f->m_exp[i]->GetPhys(), m_f->m_exp[i]->UpdateCoeffs());
245  }
246  // save field names
247  m_f->m_variables = fromField->m_fielddef[0]->m_fields;
248 }
249 
250 void ProcessInterpField::PrintProgressbar(const int position,
251  const int goal) const
252 {
253  LibUtilities::PrintProgressbar(position, goal, "Interpolating");
254 }
255 }
256 }
A class that contains algorithms for interpolation between pts fields, expansions and different meshe...
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mod...
Definition: ErrorUtil.hpp:209
int PrintProgressbar(const int position, const int goal, const std::string message, int lastprogress=-1)
Prints a progressbar.
Definition: Progressbar.hpp:67
std::map< std::string, ConfigOption > m_config
List of configuration values.
static Array< OneD, NekDouble > NullNekDouble1DArray
void SetProgressCallback(FuncPointerT func, ObjectPointerT obj)
sets a callback funtion which gets called every time the interpolation progresses ...
Represents a command-line configuration option.
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:782
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:874
void PrintProgressbar(const int position, const int goal) const
STL namespace.
std::shared_ptr< Field > FieldSharedPtr
Definition: Field.hpp:762
std::shared_ptr< DomainRange > DomainRangeShPtr
Definition: MeshGraph.h:128
CommFactory & GetCommFactory()
static SessionReaderSharedPtr CreateInstance(int argc, char *argv[])
Creates an instance of the SessionReader class.
std::pair< ModuleType, std::string > ModuleKey
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
FIELD_UTILS_EXPORT void Interpolate(const std::vector< MultiRegions::ExpListSharedPtr > expInField, std::vector< MultiRegions::ExpListSharedPtr > &expOutField, NekDouble def_value=0.0)
Interpolate from an expansion to an expansion.
double NekDouble
static bool GenerateVector(const std::string &str, std::vector< T > &out)
Takes a comma-separated string and converts it to entries in a vector.
Definition: ParseUtils.cpp:135
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:199
virtual void Process(po::variables_map &vm)
Write mesh to output file.
Abstract base class for processing modules.
static FieldMetaDataMap NullFieldMetaDataMap
Definition: FieldIO.h:53
static MeshGraphSharedPtr Read(const LibUtilities::SessionReaderSharedPtr pSession, DomainRangeShPtr rng=NullDomainRangeShPtr, bool fillGraph=true)
Definition: MeshGraph.cpp:113
std::map< int, ExpansionShPtr > ExpansionMap
Definition: MeshGraph.h:152
ModuleFactory & GetModuleFactory()
FieldSharedPtr m_f
Field object.