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