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