Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
FilterHistoryPoints.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File FilterHistoryPoints.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: Outputs values at specific points during time-stepping.
33 //
34 ///////////////////////////////////////////////////////////////////////////////
35 
37 #include <iomanip>
39 
40 #include <boost/format.hpp>
41 
42 using namespace std;
43 
44 namespace Nektar
45 {
46 namespace SolverUtils
47 {
48 std::string FilterHistoryPoints::className = GetFilterFactory().RegisterCreatorFunction("HistoryPoints", FilterHistoryPoints::create);
49 
50 /**
51  *
52  */
53 FilterHistoryPoints::FilterHistoryPoints(
55  const ParamMap &pParams) :
56  Filter(pSession)
57 {
58  ParamMap::const_iterator it;
59 
60  // OutputFile
61  it = pParams.find("OutputFile");
62  if (it == pParams.end())
63  {
64  m_outputFile = m_session->GetSessionName();
65  }
66  else
67  {
68  ASSERTL0(it->second.length() > 0, "Missing parameter 'OutputFile'.");
69  m_outputFile = it->second;
70  }
71  if (!(m_outputFile.length() >= 4
72  && m_outputFile.substr(m_outputFile.length() - 4) == ".his"))
73  {
74  m_outputFile += ".his";
75  }
76 
77  // OutputFrequency
78  it = pParams.find("OutputFrequency");
79  if (it == pParams.end())
80  {
82  }
83  else
84  {
85  LibUtilities::Equation equ(m_session, it->second);
86  m_outputFrequency = floor(equ.Evaluate());
87  }
88 
89  // OutputPlane
90  m_session->MatchSolverInfo("Homogeneous", "1D", m_isHomogeneous1D, false);
92  {
93  it = pParams.find("OutputPlane");
94  if (it == pParams.end())
95  {
96  m_outputPlane = 0;
97  }
98  else
99  {
100  LibUtilities::Equation equ(m_session, it->second);
101  m_outputPlane = floor(equ.Evaluate());
102  }
103  }
104 
105  // Points
106  it = pParams.find("Points");
107  ASSERTL0(it != pParams.end(), "Missing parameter 'Points'.");
108  m_historyPointStream.str(it->second);
109  m_index = 0;
110 }
111 
112 
113 /**
114  *
115  */
117 {
118 
119 }
120 
121 
122 /**
123  *
124  */
127  const NekDouble &time)
128 {
130  "No history points in stream.");
131 
132  m_index = 0;
133  m_historyList.clear();
134 
135  // Read history points
136  Array<OneD, NekDouble> gloCoord(3,0.0);
137  int dim = pFields[0]->GetGraph()->GetSpaceDimension();
138  int i = 0;
139  while (!m_historyPointStream.fail())
140  {
141  m_historyPointStream >> gloCoord[0]
142  >> gloCoord[1]
143  >> gloCoord[2];
144  if(m_isHomogeneous1D) // overwrite with plane z
145  {
146  NekDouble Z = (pFields[0]->GetHomogeneousBasis()
147  ->GetZ())[m_outputPlane];
148  if(fabs(gloCoord[2]-Z) > NekConstants::kVertexTheSameDouble)
149  {
150  cout << "Reseting History point from " << gloCoord[2]
151  << " to " << Z << endl;
152  }
153  gloCoord[2] = Z;
154  }
155 
156  if (!m_historyPointStream.fail())
157  {
160  ::AllocateSharedPtr(dim, i, gloCoord[0],
161  gloCoord[1], gloCoord[2]);
162 
163  m_historyPoints.push_back(vert);
164  ++i;
165  }
166  }
167 
168 
169  // Determine the unique process responsible for each history point
170  // For points on a partition boundary, must select a single process
171  LibUtilities::CommSharedPtr vComm = pFields[0]->GetComm();
172  int vRank = vComm->GetRank();
173  int vHP = m_historyPoints.size();
174  Array<OneD, int> procList(vHP, -1 );
175  Array<OneD, int> idList (vHP, -1 );
176  Array<OneD, NekDouble> dist (vHP, 1e16);
177  Array<OneD, NekDouble> dist_loc(vHP, 1e16);
178  std::vector<Array<OneD, NekDouble> > LocCoords;
179 
180  // Find the nearest element on this process to which the history
181  // point could belong and note down the distance from the element
182  // and the process ID.
183  for (i = 0; i < vHP; ++i)
184  {
185  Array<OneD, NekDouble> locCoords(3);
186  m_historyPoints[i]->GetCoords( gloCoord[0],
187  gloCoord[1],
188  gloCoord[2]);
189 
190  // Determine the expansion and local coordinates
191  idList[i] = pFields[0]->GetExpIndex(gloCoord,locCoords,
193 
194  // Save Local coordinates for later
195  LocCoords.push_back(locCoords);
196 
197  // For those points for which a potential nearby element exists
198  // compute the perp. distance from the point to the element and
199  // store in the distances array.
200  if (idList[i] != -1)
201  {
203  pFields[0]->GetExp(idList[i])->GetGeom();
204  StdRegions::StdExpansionSharedPtr e = g->GetXmap();
205  Array<OneD, NekDouble> coordVals(e->GetTotPoints());
206  dist_loc[i] = 0.0;
207  for (int j = 0; j < g->GetCoordim(); ++j)
208  {
209  e->BwdTrans(g->GetCoeffs(j), coordVals);
210  NekDouble x = e->PhysEvaluate(locCoords, coordVals)
211  - gloCoord[j];
212  dist_loc[i] += x*x;
213  }
214  }
215  }
216 
217  // Reduce distances of points from elements, keeping the smallest
218  // distance.
219  Vmath::Vcopy(vHP, dist_loc, 1, dist, 1);
220  vComm->AllReduce(dist, LibUtilities::ReduceMin);
221 
222  // If multiple processes find they are the nearest (e.g. point lies
223  // on a partition boundary, we will choose the process of highest
224  // rank.
225  for (i = 0; i < vHP; ++i)
226  {
227  if (dist_loc[i] == dist[i])
228  {
229  // Set element id to Vid of m_history point for later use
230  m_historyPoints[i]->SetVid(idList[i]);
231  }
232  else
233  {
234  // This history point is not handled by this process
235  idList[i] = -1;
236  }
237 
238  // If a matching element is found on this process, note the
239  // process ID
240  if (idList[i] != -1)
241  {
243  {
244  int j;
246  = pFields[0]->GetZIDs();
247  for(j = 0; j < IDs.num_elements(); ++j)
248  {
249  if(IDs[j] == m_outputPlane)
250  {
251  break;
252  }
253  }
254 
255  if(j != IDs.num_elements())
256  {
257  m_outputPlane = j;
258  procList[i] = vRank;
259  }
260  }
261  else
262  {
263  procList[i] = vRank;
264  }
265  }
266  }
267 
268  // Reduce process IDs for all history points. The process with
269  // largest rank will handle the history point in the case where the
270  // distance was the same.
271  vComm->AllReduce(procList, LibUtilities::ReduceMax);
272 
273  // Determine the element in which each history point resides.
274  // If point is not in mesh (on this process), id is -1.
275  for (i = 0; i < vHP; ++i)
276  {
277  // If point lies on partition boundary, only the proc with max
278  // rank retains possession.
279  if (procList[i] != vRank)
280  {
281  idList[i] = -1;
282  }
283 
284  // If the current process owns this history point, add it to its
285  // local list of history points.
286  if (idList[i] != -1)
287  {
289  m_historyList.push_back(
292  (m_historyPoints[i], LocCoords[i]));
293  }
294  }
295 
296  // Collate the element ID list across processes and check each
297  // history point is allocated to a process
298  vComm->AllReduce(idList, LibUtilities::ReduceMax);
299  if (vComm->GetRank() == 0)
300  {
301  for (i = 0; i < vHP; ++i)
302  {
303  m_historyPoints[i]->GetCoords( gloCoord[0],
304  gloCoord[1],
305  gloCoord[2]);
306 
307  // Write an error if no process owns history point
308  ASSERTL0(idList[i] != -1,
309  "History point "
310  + boost::lexical_cast<std::string>(gloCoord[0])
311  + ", "
312  + boost::lexical_cast<std::string>(gloCoord[1])
313  + ", "
314  + boost::lexical_cast<std::string>(gloCoord[2])
315  + " cannot be found in the mesh.");
316 
317  // Print a warning if a process owns it but it is not close
318  // enough to the element.
319  if (dist[i] > NekConstants::kGeomFactorsTol)
320  {
321  cout << "Warning: History point " << i << " at ("
322  << gloCoord[0] << "," << gloCoord[1] << ","
323  << gloCoord[2] << ") lies a distance of "
324  << sqrt(dist[i]) << " from the manifold." << endl;
325  }
326  }
327 
328  // Open output stream
329  m_outputStream.open(m_outputFile.c_str());
330  m_outputStream << "# History data for variables (:";
331 
332  for (i = 0; i < pFields.num_elements(); ++i)
333  {
334  m_outputStream << m_session->GetVariable(i) <<",";
335  }
336 
338  {
339  m_outputStream << ") at points:";
340  }
341  else
342  {
343  m_outputStream << ") at points:" << endl;
344  }
345 
346  for (i = 0; i < vHP; ++i)
347  {
348  m_historyPoints[i]->GetCoords( gloCoord[0],
349  gloCoord[1],
350  gloCoord[2]);
351 
352  m_outputStream << "# " << boost::format("%6.0f") % i;
353  m_outputStream << " " << boost::format("%25e") % gloCoord[0];
354  m_outputStream << " " << boost::format("%25e") % gloCoord[1];
355  m_outputStream << " " << boost::format("%25e") % gloCoord[2];
356  m_outputStream << endl;
357  }
358 
360  {
361  m_outputStream << "(in Wavespace)" << endl;
362  }
363  }
364  v_Update(pFields, time);
365 }
366 
367 
368 /**
369  *
370  */
372 {
373  // Only output every m_outputFrequency.
374  if ((m_index++) % m_outputFrequency)
375  {
376  return;
377  }
378 
379  int j = 0;
380  int k = 0;
381  int numPoints = m_historyPoints.size();
382  int numFields = pFields.num_elements();
383  LibUtilities::CommSharedPtr vComm = pFields[0]->GetComm();
384  Array<OneD, NekDouble> data(numPoints*numFields, 0.0);
385  Array<OneD, NekDouble> gloCoord(3, 0.0);
386  std::list<std::pair<SpatialDomains::PointGeomSharedPtr, Array<OneD, NekDouble> > >::iterator x;
387 
388  Array<OneD, NekDouble> physvals;
389  Array<OneD, NekDouble> locCoord;
390  int expId;
391 
392  // Pull out data values field by field
393  for (j = 0; j < numFields; ++j)
394  {
396  {
397  for (k = 0, x = m_historyList.begin(); x != m_historyList.end();
398  ++x, ++k)
399  {
400  locCoord = (*x).second;
401  expId = (*x).first->GetVid();
402 
403  physvals = pFields[j]->GetPlane(m_outputPlane)->UpdatePhys() + pFields[j]->GetPhys_Offset(expId);
404 
405  // transform elemental data if required.
406  if(pFields[j]->GetPhysState() == false)
407  {
408  pFields[j]->GetPlane(m_outputPlane)->GetExp(expId)->BwdTrans(pFields[j]->GetPlane(m_outputPlane)->GetCoeffs() + pFields[j]->GetCoeff_Offset(expId),physvals);
409  }
410 
411  // interpolate point can do with zero plane methods
412  data[m_historyLocalPointMap[k]*numFields+j] = pFields[j]->GetExp(expId)->StdPhysEvaluate(locCoord,physvals);
413 
414  }
415  }
416  else
417  {
418  for (k = 0, x = m_historyList.begin(); x != m_historyList.end(); ++x, ++k)
419  {
420  locCoord = (*x).second;
421  expId = (*x).first->GetVid();
422 
423  physvals = pFields[j]->UpdatePhys() + pFields[j]->GetPhys_Offset(expId);
424 
425  // transform elemental data if required.
426  if(pFields[j]->GetPhysState() == false)
427  {
428  pFields[j]->GetExp(expId)->BwdTrans(pFields[j]->GetCoeffs() + pFields[j]->GetCoeff_Offset(expId),physvals);
429  }
430 
431  // interpolate point
432  data[m_historyLocalPointMap[k]*numFields+j] = pFields[j]->GetExp(expId)->StdPhysEvaluate(locCoord,physvals);
433  }
434  }
435  }
436 
437  // Exchange history data
438  // This could be improved to reduce communication but works for now
439  vComm->AllReduce(data, LibUtilities::ReduceSum);
440 
441  // Only the root process writes out history data
442  if (vComm->GetRank() == 0)
443  {
444 
445  // Write data values point by point
446  for (k = 0; k < m_historyPoints.size(); ++k)
447  {
448  m_outputStream << boost::format("%25e") % time;
449  for (int j = 0; j < numFields; ++j)
450  {
451  m_outputStream << " " << boost::format("%25e") % data[k*numFields+j];
452  }
453  m_outputStream << endl;
454  }
455  }
456 }
457 
458 
459 /**
460  *
461  */
463 {
464  if (pFields[0]->GetComm()->GetRank() == 0)
465  {
466  m_outputStream.close();
467  }
468 }
469 
470 
471 /**
472  *
473  */
475 {
476  return true;
477 }
478 }
479 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
STL namespace.
static const NekDouble kVertexTheSameDouble
boost::shared_ptr< SessionReader > SessionReaderSharedPtr
Definition: MeshPartition.h:51
boost::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:53
unsigned int m_outputPlane
plane to take history point from if using a homogeneous1D expansion
NekDouble Evaluate() const
Definition: Equation.h:102
SpatialDomains::PointGeomVector m_historyPoints
double NekDouble
std::map< std::string, std::string > ParamMap
Definition: Filter.h:67
LibUtilities::SessionReaderSharedPtr m_session
Definition: Filter.h:84
virtual SOLVER_UTILS_EXPORT void v_Finalise(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time)
virtual SOLVER_UTILS_EXPORT void v_Initialise(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time)
StandardMatrixTag boost::call_traits< LhsDataType >::const_reference rhs typedef NekMatrix< LhsDataType, StandardMatrixTag >::iterator iterator
virtual SOLVER_UTILS_EXPORT bool v_IsTimeDependent()
static const NekDouble kGeomFactorsTol
FilterFactory & GetFilterFactory()
Definition: Filter.cpp:42
boost::shared_ptr< StdExpansion > StdExpansionSharedPtr
virtual SOLVER_UTILS_EXPORT void v_Update(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time)
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1047
boost::shared_ptr< Geometry > GeometrySharedPtr
Definition: Geometry.h:53
std::list< std::pair< SpatialDomains::PointGeomSharedPtr, Array< OneD, NekDouble > > > m_historyList
boost::shared_ptr< PointGeom > PointGeomSharedPtr
Definition: Geometry.h:60
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, tDescription pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:215