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