Nektar++
FilterMovingBody.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File FilterMovingBody.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: Output moving body forces during time-stepping.
32 //
33 ///////////////////////////////////////////////////////////////////////////////
34 
37 #include <iomanip>
42 
43 using namespace std;
44 
45 namespace Nektar
46 {
47 
48 std::string FilterMovingBody::className = SolverUtils::GetFilterFactory().
49  RegisterCreatorFunction("MovingBody",
50  FilterMovingBody::create,
51  "Moving Body Filter");
52 /**
53  *
54  */
55 FilterMovingBody::FilterMovingBody(
57  const std::weak_ptr<SolverUtils::EquationSystem> &pEquation,
58  const ParamMap &pParams)
59  : Filter(pSession, pEquation)
60 {
61  // OutputFile
62  auto it = pParams.find("OutputFile");
63  if (it == pParams.end())
64  {
65  m_outputFile_fce = pSession->GetSessionName();
66  m_outputFile_mot = pSession->GetSessionName();
67  }
68  else
69  {
70  ASSERTL0(it->second.length() > 0, "Missing parameter 'OutputFile'.");
71 
72  m_outputFile_fce = it->second;
73  m_outputFile_mot = it->second;
74  }
75  if (!(m_outputFile_fce.length() >= 4 &&
76  m_outputFile_fce.substr(m_outputFile_fce.length() - 4) == ".fce"))
77  {
78  m_outputFile_fce += ".fce";
79  }
80  if (!(m_outputFile_mot.length() >= 4 &&
81  m_outputFile_mot.substr(m_outputFile_mot.length() - 4) == ".mot"))
82  {
83  m_outputFile_mot += ".mot";
84  }
85 
86  // OutputFrequency
87  it = pParams.find("OutputFrequency");
88  if (it == pParams.end())
89  {
91  }
92  else
93  {
94  LibUtilities::Equation equ(m_session->GetInterpreter(), it->second);
95  m_outputFrequency = round(equ.Evaluate());
96  }
97 
98  pSession->MatchSolverInfo("Homogeneous", "1D", m_isHomogeneous1D, false);
99  ASSERTL0(m_isHomogeneous1D, "Moving Body implemented just for 3D "
100  "Homogeneous 1D discetisations.");
101 
102  // Boundary (to calculate the forces)
103  it = pParams.find("Boundary");
104  ASSERTL0(it != pParams.end(), "Missing parameter 'Boundary'.");
105  ASSERTL0(it->second.length() > 0, "Missing parameter 'Boundary'.");
106  m_BoundaryString = it->second;
107 }
108 
109 
110 /**
111  *
112  */
114 {
115 
116 }
117 
118 
119 /**
120  *
121  */
124  const NekDouble &time)
125 {
126  m_index_f = 0;
127  m_index_m = 0;
129  // Parse the boundary regions into a list.
130  std::string::size_type FirstInd = m_BoundaryString.find_first_of('[') + 1;
131  std::string::size_type LastInd = m_BoundaryString.find_last_of(']') - 1;
132 
133  ASSERTL0(FirstInd <= LastInd,
134  (std::string("Error reading boundary region definition:") +
135  m_BoundaryString).c_str());
136 
137  std::string IndString = m_BoundaryString.substr(FirstInd,
138  LastInd - FirstInd + 1);
139 
140  bool parseGood = ParseUtils::GenerateSeqVector(IndString,
142 
143  ASSERTL0(parseGood && !m_boundaryRegionsIdList.empty(),
144  (std::string("Unable to read boundary regions index "
145  "range for FilterAeroForces: ") + IndString).c_str());
146 
147  // determine what boundary regions need to be considered
148  unsigned int numBoundaryRegions
149  = pFields[0]->GetBndConditions().num_elements();
150 
152  numBoundaryRegions, 0);
153 
154  SpatialDomains::BoundaryConditions bcs(m_session,pFields[0]->GetGraph());
155 
157  = bcs.GetBoundaryRegions();
158 
159  int cnt = 0;
160  for (auto &it : bregions)
161  {
162  if ( std::find(m_boundaryRegionsIdList.begin(),
163  m_boundaryRegionsIdList.end(), it.first) !=
165  {
166  m_boundaryRegionIsInList[cnt] = 1;
167  }
168  ++cnt;
169  }
170 
171  LibUtilities::CommSharedPtr vComm = pFields[0]->GetComm();
172 
173  if (vComm->GetRank() == 0)
174  {
175  // Open output stream for cable forces
176  m_outputStream[0].open(m_outputFile_fce.c_str());
177  m_outputStream[0] << "#";
178  m_outputStream[0].width(7);
179  m_outputStream[0] << "Time";
180  m_outputStream[0].width(15);
181  m_outputStream[0] << "z";
182  m_outputStream[0].width(15);
183  m_outputStream[0] << "Fx (press)";
184  m_outputStream[0].width(15);
185  m_outputStream[0] << "Fx (visc)";
186  m_outputStream[0].width(15);
187  m_outputStream[0] << "Fx (tot)";
188  m_outputStream[0].width(15);
189  m_outputStream[0] << "Fy (press)";
190  m_outputStream[0].width(15);
191  m_outputStream[0] << "Fy (visc)";
192  m_outputStream[0].width(15);
193  m_outputStream[0] << "Fy (tot)";
194  m_outputStream[0] << endl;
195 
196  // Open output stream for cable motions
197  m_outputStream[1].open(m_outputFile_mot.c_str());
198  m_outputStream[1] << "#";
199  m_outputStream[1].width(7);
200  m_outputStream[1] << "Time";
201  m_outputStream[1].width(15);
202  m_outputStream[1] << "z";
203  m_outputStream[1].width(15);
204  m_outputStream[1] << "Disp_x";
205  m_outputStream[1].width(15);
206  m_outputStream[1] << "Vel_x";
207  m_outputStream[1].width(15);
208  m_outputStream[1] << "Acel_x";
209  m_outputStream[1].width(15);
210  m_outputStream[1] << "Disp_y";
211  m_outputStream[1].width(15);
212  m_outputStream[1] << "Vel_y";
213  m_outputStream[1].width(15);
214  m_outputStream[1] << "Acel_y";
215  m_outputStream[1] << endl;
216  }
217 }
218 
219 
220 /**
221  *
222  */
224  const LibUtilities::SessionReaderSharedPtr &pSession,
226  Array<OneD, NekDouble> &Aeroforces,
227  const NekDouble &time)
228 {
229  int n, cnt, elmtid, nq, offset, boundary;
230  int nt = pFields[0]->GetNpoints();
231  int dim = pFields.num_elements()-1;
232 
234  Array<OneD, int> BoundarytoElmtID;
235  Array<OneD, int> BoundarytoTraceID;
237 
242 
246 
250 
251  LibUtilities::CommSharedPtr vComm = pFields[0]->GetComm();
252  LibUtilities::CommSharedPtr vRowComm = vComm->GetRowComm();
253  LibUtilities::CommSharedPtr vColComm = vComm->GetColumnComm();
254  LibUtilities::CommSharedPtr vCColComm = vColComm->GetColumnComm();
255 
256  // set up storage space for forces on all the planes (z-positions)
257  // on each processors some of them will remain empty as we may
258  // have just few planes per processor
259  int Num_z_pos = pFields[0]->GetHomogeneousBasis()->GetNumModes();
260  Array<OneD, NekDouble> Fx(Num_z_pos,0.0);
261  Array<OneD, NekDouble> Fxp(Num_z_pos,0.0);
262  Array<OneD, NekDouble> Fxv(Num_z_pos,0.0);
263  Array<OneD, NekDouble> Fy(Num_z_pos,0.0);
264  Array<OneD, NekDouble> Fyp(Num_z_pos,0.0);
265  Array<OneD, NekDouble> Fyv(Num_z_pos,0.0);
266 
267 
268  NekDouble rho = (pSession->DefinesParameter("rho"))
269  ? (pSession->GetParameter("rho"))
270  : 1;
271  NekDouble mu = rho*pSession->GetParameter("Kinvis");
272 
273  for(int i = 0; i < pFields.num_elements(); ++i)
274  {
275  pFields[i]->SetWaveSpace(false);
276  pFields[i]->BwdTrans(pFields[i]->GetCoeffs(),
277  pFields[i]->UpdatePhys());
278  pFields[i]->SetPhysState(true);
279  }
280 
281  // Get the number of local planes on the process and their IDs
282  // to properly locate the forces in the Fx, Fy etc. vectors.
284  ZIDs = pFields[0]->GetZIDs();
285  int local_planes = ZIDs.num_elements();
286 
287  // Homogeneous 1D case Compute forces on all WALL boundaries
288  // This only has to be done on the zero (mean) Fourier mode.
289  for(int plane = 0 ; plane < local_planes; plane++)
290  {
291  pFields[0]->GetPlane(plane)->GetBoundaryToElmtMap(BoundarytoElmtID,
292  BoundarytoTraceID);
293  BndExp = pFields[0]->GetPlane(plane)->GetBndCondExpansions();
295 
296  // loop over the types of boundary conditions
297  for(cnt = n = 0; n < BndExp.num_elements(); ++n)
298  {
299  if(m_boundaryRegionIsInList[n] == 1)
300  {
301  for(int i = 0; i < BndExp[n]->GetExpSize(); ++i, cnt++)
302  {
303  // find element of this expansion.
304  elmtid = BoundarytoElmtID[cnt];
305  elmt = pFields[0]->GetPlane(plane)->GetExp(elmtid);
306  nq = elmt->GetTotPoints();
307  offset = pFields[0]->GetPlane(plane)
308  ->GetPhys_Offset(elmtid);
309 
310  // Initialise local arrays for the velocity
311  // gradients size of total number of quadrature
312  // points for each element (hence local).
313  for(int j = 0; j < dim; ++j)
314  {
315  gradU[j] = Array<OneD, NekDouble>(nq,0.0);
316  gradV[j] = Array<OneD, NekDouble>(nq,0.0);
317  gradW[j] = Array<OneD, NekDouble>(nq,0.0);
318  }
319 
320  // identify boundary of element
321  boundary = BoundarytoTraceID[cnt];
322 
323  // Extract fields
324  U = pFields[0]->GetPlane(plane)->GetPhys() + offset;
325  V = pFields[1]->GetPlane(plane)->GetPhys() + offset;
326  P = pFields[3]->GetPlane(plane)->GetPhys() + offset;
327 
328  // compute the gradients
329  elmt->PhysDeriv(U,gradU[0],gradU[1]);
330  elmt->PhysDeriv(V,gradV[0],gradV[1]);
331 
332  // Get face 1D expansion from element expansion
333  bc = BndExp[n]->GetExp(i)->as<LocalRegions::Expansion1D>();
334 
335  // number of points on the boundary
336  int nbc = bc->GetTotPoints();
337 
338  // several vectors for computing the forces
339  Array<OneD, NekDouble> Pb(nbc,0.0);
340 
341  for(int j = 0; j < dim; ++j)
342  {
343  fgradU[j] = Array<OneD, NekDouble>(nbc,0.0);
344  fgradV[j] = Array<OneD, NekDouble>(nbc,0.0);
345  }
346 
347  Array<OneD, NekDouble> drag_t(nbc,0.0);
348  Array<OneD, NekDouble> lift_t(nbc,0.0);
349  Array<OneD, NekDouble> drag_p(nbc,0.0);
350  Array<OneD, NekDouble> lift_p(nbc,0.0);
351  Array<OneD, NekDouble> temp(nbc,0.0);
352  Array<OneD, NekDouble> temp2(nbc,0.0);
353 
354  // identify boundary of element .
355  boundary = BoundarytoTraceID[cnt];
356 
357  // extraction of the pressure and wss on the
358  // boundary of the element
359  elmt->GetEdgePhysVals(boundary,bc,P,Pb);
360 
361  for(int j = 0; j < dim; ++j)
362  {
363  elmt->GetEdgePhysVals(boundary,bc,gradU[j],fgradU[j]);
364  elmt->GetEdgePhysVals(boundary,bc,gradV[j],fgradV[j]);
365  }
366 
367  //normals of the element
368  const Array<OneD, Array<OneD, NekDouble> > &normals
369  = elmt->GetEdgeNormal(boundary);
370 
371  //
372  // Compute viscous tractive forces on wall from
373  //
374  // t_i = - T_ij * n_j (minus sign for force
375  // exerted BY fluid ON wall),
376  //
377  // where
378  //
379  // T_ij = viscous stress tensor (here in Cartesian
380  // coords)
381  // dU_i dU_j
382  // = RHO * KINVIS * ( ---- + ---- ) .
383  // dx_j dx_i
384 
385  //a) DRAG TERMS
386  //-rho*kinvis*(2*du/dx*nx+(du/dy+dv/dx)*ny)
387 
388  Vmath::Vadd(nbc, fgradU[1], 1, fgradV[0], 1, drag_t, 1);
389  Vmath::Vmul(nbc, drag_t, 1, normals[1], 1, drag_t, 1);
390 
391  Vmath::Smul(nbc, 2.0, fgradU[0], 1, fgradU[0], 1);
392  Vmath::Vmul(nbc, fgradU[0], 1, normals[0], 1, temp2, 1);
393  Vmath::Smul(nbc, 0.5, fgradU[0], 1, fgradU[0], 1);
394 
395  Vmath::Vadd(nbc, temp2, 1, drag_t, 1, drag_t, 1);
396  Vmath::Smul(nbc, -mu, drag_t, 1, drag_t, 1);
397 
398  //zero temporary storage vector
399  Vmath::Zero(nbc, temp, 0.0);
400  Vmath::Zero(nbc, temp2, 0.0);
401 
402 
403  //b) LIFT TERMS
404  //-rho*kinvis*(2*dv/dy*ny+(du/dy+dv/dx)*nx)
405 
406  Vmath::Vadd(nbc, fgradU[1], 1, fgradV[0], 1, lift_t, 1);
407  Vmath::Vmul(nbc, lift_t, 1, normals[0], 1, lift_t, 1);
408 
409  Vmath::Smul(nbc, 2.0, fgradV[1], 1, fgradV[1], 1);
410  Vmath::Vmul(nbc, fgradV[1], 1, normals[1], 1, temp2, 1);
411  Vmath::Smul(nbc, -0.5, fgradV[1], 1, fgradV[1], 1);
412 
413 
414  Vmath::Vadd(nbc, temp2, 1, lift_t, 1, lift_t, 1);
415  Vmath::Smul(nbc, -mu, lift_t, 1, lift_t, 1);
416 
417  // Compute normal tractive forces on all WALL
418  // boundaries
419 
420  Vmath::Vvtvp(nbc, Pb, 1, normals[0], 1, drag_p, 1,
421  drag_p, 1);
422  Vmath::Vvtvp(nbc, Pb, 1, normals[1], 1, lift_p, 1,
423  lift_p, 1);
424 
425  //integration over the boundary
426  Fxv[ZIDs[plane]] += bc->Integral(drag_t);
427  Fyv[ZIDs[plane]] += bc->Integral(lift_t);
428 
429  Fxp[ZIDs[plane]] += bc->Integral(drag_p);
430  Fyp[ZIDs[plane]] += bc->Integral(lift_p);
431  }
432  }
433  else
434  {
435  cnt += BndExp[n]->GetExpSize();
436  }
437  }
438  }
439 
440  for(int i = 0; i < pFields.num_elements(); ++i)
441  {
442  pFields[i]->SetWaveSpace(true);
443  pFields[i]->BwdTrans(pFields[i]->GetCoeffs(),
444  pFields[i]->UpdatePhys());
445  pFields[i]->SetPhysState(false);
446  }
447 
448  // In case we are using an hybrid parallelisation we need
449  // to reduce the forces on the same plane coming from
450  // different mesh partitions.
451  // It is quite an expensive communication, therefore
452  // we check to make sure it is actually required.
453 
454  if(vComm->GetRowComm()->GetSize() > 0)
455  {
456  // NOTE 1: We can eventually sum the viscous and pressure
457  // component before doing the communication, thus
458  // reducing by a factor 2 the communication.
459  // NOTE 2: We may want to set up in the Comm class an AllReduce
460  // routine wich can handle arrays more efficiently
461  for(int plane = 0 ; plane < local_planes; plane++)
462  {
463  vRowComm->AllReduce(Fxp[ZIDs[plane]], LibUtilities::ReduceSum);
464  vRowComm->AllReduce(Fxv[ZIDs[plane]], LibUtilities::ReduceSum);
465  vRowComm->AllReduce(Fyp[ZIDs[plane]], LibUtilities::ReduceSum);
466  vRowComm->AllReduce(Fyv[ZIDs[plane]], LibUtilities::ReduceSum);
467  }
468  }
469 
470  // At this point on rank (0,n) of the Mesh partion communicator we have
471  // the total areo forces on the planes which are on the same column
472  // communicator. Since the planes are scattered on different processors
473  // some of the entries of the vector Fxp, Fxp etc. are still zero.
474  // Now we need to reduce the values on a single vector on rank (0,0) of the
475  // global communicator.
476  if(!pSession->DefinesSolverInfo("HomoStrip"))
477  {
478  if(vComm->GetRowComm()->GetRank() == 0)
479  {
480  for(int z = 0 ; z < Num_z_pos; z++)
481  {
482  vColComm->AllReduce(Fxp[z], LibUtilities::ReduceSum);
483  vColComm->AllReduce(Fxv[z], LibUtilities::ReduceSum);
484  vColComm->AllReduce(Fyp[z], LibUtilities::ReduceSum);
485  vColComm->AllReduce(Fyv[z], LibUtilities::ReduceSum);
486  }
487  }
488  }
489  else
490  {
491  if(vComm->GetRowComm()->GetRank() == 0)
492  {
493  for(int z = 0 ; z < Num_z_pos; z++)
494  {
495  vCColComm->AllReduce(Fxp[z], LibUtilities::ReduceSum);
496  vCColComm->AllReduce(Fxv[z], LibUtilities::ReduceSum);
497  vCColComm->AllReduce(Fyp[z], LibUtilities::ReduceSum);
498  vCColComm->AllReduce(Fyv[z], LibUtilities::ReduceSum);
499  }
500  }
501  }
502 
503  if(!pSession->DefinesSolverInfo("HomoStrip"))
504  {
505  //set the forces imparted on the cable's wall
506  for(int plane = 0 ; plane < local_planes; plane++)
507  {
508  Aeroforces[plane] = Fxp[ZIDs[plane]]
509  + Fxv[ZIDs[plane]];
510  Aeroforces[plane + local_planes] = Fyp[ZIDs[plane]]
511  + Fyv[ZIDs[plane]];
512  }
513 
514  // Only output every m_outputFrequency.
515  if ((m_index_f++) % m_outputFrequency)
516  {
517  return;
518  }
519 
520  // At thi point in rank (0,0) we have the full vectors
521  // containing Fxp,Fxv,Fyp and Fyv where different positions
522  // in the vectors correspond to different planes.
523  // Here we write it to file. We do it just on one porcess
524 
525  Array<OneD, NekDouble> z_coords(Num_z_pos,0.0);
527  = pFields[0]->GetHomogeneousBasis()->GetZ();
528 
529  NekDouble LZ;
530  pSession->LoadParameter("LZ", LZ);
531  Vmath::Smul(Num_z_pos,LZ/2.0,pts,1,z_coords,1);
532  Vmath::Sadd(Num_z_pos,LZ/2.0,z_coords,1,z_coords,1);
533  if (vComm->GetRank() == 0)
534  {
535 
536  Vmath::Vadd(Num_z_pos,Fxp,1,Fxv,1,Fx,1);
537  Vmath::Vadd(Num_z_pos,Fyp,1,Fyv,1,Fy,1);
538 
539  for(int i = 0 ; i < Num_z_pos; i++)
540  {
541  m_outputStream[0].width(8);
542  m_outputStream[0] << setprecision(6) << time;
543 
544  m_outputStream[0].width(15);
545  m_outputStream[0] << setprecision(6) << z_coords[i];
546 
547  m_outputStream[0].width(15);
548  m_outputStream[0] << setprecision(8) << Fxp[i];
549  m_outputStream[0].width(15);
550  m_outputStream[0] << setprecision(8) << Fxv[i];
551  m_outputStream[0].width(15);
552  m_outputStream[0] << setprecision(8) << Fx[i];
553 
554  m_outputStream[0].width(15);
555  m_outputStream[0] << setprecision(8) << Fyp[i];
556  m_outputStream[0].width(15);
557  m_outputStream[0] << setprecision(8) << Fyv[i];
558  m_outputStream[0].width(15);
559  m_outputStream[0] << setprecision(8) << Fy[i];
560  m_outputStream[0] << endl;
561  }
562  }
563  }
564  else
565  {
566  // average the forces over local planes for each processor
567  Array<OneD, NekDouble> fces(6,0.0);
568  for(int plane = 0 ; plane < local_planes; plane++)
569  {
570  fces[0] += Fxp[ZIDs[plane]] + Fxv[ZIDs[plane]];
571  fces[1] += Fyp[ZIDs[plane]] + Fyv[ZIDs[plane]];
572  fces[2] += Fxp[ZIDs[plane]] ;
573  fces[3] += Fyp[ZIDs[plane]] ;
574  fces[4] += Fxv[ZIDs[plane]] ;
575  fces[5] += Fyv[ZIDs[plane]] ;
576  }
577 
578  fces[0] = fces[0]/local_planes;
579  fces[1] = fces[1]/local_planes;
580  fces[2] = fces[2]/local_planes;
581  fces[3] = fces[3]/local_planes;
582  fces[4] = fces[4]/local_planes;
583  fces[5] = fces[5]/local_planes;
584 
585  // average the forces over communicators within each strip
586  vCColComm->AllReduce(fces[0], LibUtilities::ReduceSum);
587  vCColComm->AllReduce(fces[1], LibUtilities::ReduceSum);
588  vCColComm->AllReduce(fces[2], LibUtilities::ReduceSum);
589  vCColComm->AllReduce(fces[3], LibUtilities::ReduceSum);
590  vCColComm->AllReduce(fces[4], LibUtilities::ReduceSum);
591  vCColComm->AllReduce(fces[5], LibUtilities::ReduceSum);
592 
593  int npts = vComm->GetColumnComm()->GetColumnComm()->GetSize();
594 
595  fces[0] = fces[0]/npts;
596  fces[1] = fces[1]/npts;
597  fces[2] = fces[2]/npts;
598  fces[3] = fces[3]/npts;
599  fces[4] = fces[4]/npts;
600  fces[5] = fces[5]/npts;
601 
602  for(int plane = 0 ; plane < local_planes; plane++)
603  {
604  Aeroforces[plane] = fces[0];
605  Aeroforces[plane+local_planes] = fces[1];
606  }
607 
608  // Only output every m_outputFrequency.
610  {
611  return;
612  }
613 
614  int colrank = vColComm->GetRank();
615  int nstrips;
616 
617  NekDouble DistStrip;
618 
619  pSession->LoadParameter("Strip_Z", nstrips);
620  pSession->LoadParameter("DistStrip", DistStrip);
621 
622  Array<OneD, NekDouble> z_coords(nstrips);
623  for(int i = 0; i < nstrips; i++)
624  {
625  z_coords[i] = i * DistStrip;
626  }
627 
628  if(colrank == 0)
629  {
630  m_outputStream[0].width(8);
631  m_outputStream[0] << setprecision(6) << time;
632 
633  m_outputStream[0].width(15);
634  m_outputStream[0] << setprecision(6) << z_coords[0];
635 
636  m_outputStream[0].width(15);
637  m_outputStream[0] << setprecision(8) << fces[2];
638  m_outputStream[0].width(15);
639  m_outputStream[0] << setprecision(8) << fces[4];
640  m_outputStream[0].width(15);
641  m_outputStream[0] << setprecision(8) << fces[0];
642 
643  m_outputStream[0].width(15);
644  m_outputStream[0] << setprecision(8) << fces[3];
645  m_outputStream[0].width(15);
646  m_outputStream[0] << setprecision(8) << fces[5];
647  m_outputStream[0].width(15);
648  m_outputStream[0] << setprecision(8) << fces[1];
649  m_outputStream[0] << endl;
650 
651  for(int i = 1; i < nstrips; i++)
652  {
653  vColComm->Recv(i, fces);
654 
655  m_outputStream[0].width(8);
656  m_outputStream[0] << setprecision(6) << time;
657 
658  m_outputStream[0].width(15);
659  m_outputStream[0] << setprecision(6) << z_coords[i];
660 
661  m_outputStream[0].width(15);
662  m_outputStream[0] << setprecision(8) << fces[2];
663  m_outputStream[0].width(15);
664  m_outputStream[0] << setprecision(8) << fces[4];
665  m_outputStream[0].width(15);
666  m_outputStream[0] << setprecision(8) << fces[0];
667 
668  m_outputStream[0].width(15);
669  m_outputStream[0] << setprecision(8) << fces[3];
670  m_outputStream[0].width(15);
671  m_outputStream[0] << setprecision(8) << fces[5];
672  m_outputStream[0].width(15);
673  m_outputStream[0] << setprecision(8) << fces[1];
674  m_outputStream[0] << endl;
675  }
676  }
677  else
678  {
679  for(int i = 1; i < nstrips; i++)
680  {
681  if(colrank == i)
682  {
683  vColComm->Send(0, fces);
684  }
685  }
686  }
687  }
688 }
689 
690 
691 /**
692  *
693  */
695  const LibUtilities::SessionReaderSharedPtr &pSession,
697  Array<OneD, NekDouble> &MotionVars,
698  const NekDouble &time)
699 {
700  // Only output every m_outputFrequency.
701  if ((m_index_m++) % m_outputFrequency)
702  {
703  return;
704  }
705 
706  int npts;
707  // Length of the cable
708  NekDouble Length;
709 
710  if(!pSession->DefinesSolverInfo("HomoStrip"))
711  {
712  pSession->LoadParameter("LZ", Length);
713  npts = m_session->GetParameter("HomModesZ");
714  }
715  else
716  {
717  pSession->LoadParameter("LC", Length);
718  npts = m_session->GetParameter("HomStructModesZ");
719  }
720 
721  NekDouble z_coords;
722  for(int n = 0; n < npts; n++)
723  {
724  z_coords = Length/npts*n;
725  m_outputStream[1].width(8);
726  m_outputStream[1] << setprecision(6) << time;
727  m_outputStream[1].width(15);
728  m_outputStream[1] << setprecision(6) << z_coords;
729  m_outputStream[1].width(15);
730  m_outputStream[1] << setprecision(8) << MotionVars[n];
731  m_outputStream[1].width(15);
732  m_outputStream[1] << setprecision(8) << MotionVars[npts+n];
733  m_outputStream[1].width(15);
734  m_outputStream[1] << setprecision(8) << MotionVars[2*npts+n];
735  m_outputStream[1].width(15);
736  m_outputStream[1] << setprecision(8) << MotionVars[3*npts+n];
737  m_outputStream[1].width(15);
738  m_outputStream[1] << setprecision(8) << MotionVars[4*npts+n];
739  m_outputStream[1].width(15);
740  m_outputStream[1] << setprecision(8) << MotionVars[5*npts+n];
741  m_outputStream[1] << endl;
742  }
743 }
744 
745 
746 /**
747  *
748  */
751  const NekDouble &time)
752 {
753  if (pFields[0]->GetComm()->GetRank() == 0)
754  {
755  m_outputStream[0].close();
756  m_outputStream[1].close();
757  }
758 }
759 
760 
761 /**
762  *
763  */
765 {
766  return true;
767 }
768 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
Array< OneD, std::ofstream > m_outputStream
std::vector< bool > m_boundaryRegionIsInList
Determines if a given Boundary Region is in m_boundaryRegionsIdList.
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:53
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.cpp:445
STL namespace.
virtual void v_Finalise(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time)
void UpdateMotion(const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, Array< OneD, NekDouble > &MotionVars, const NekDouble &time)
std::shared_ptr< StdExpansion > StdExpansionSharedPtr
virtual void v_Initialise(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time)
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y.
Definition: Vmath.cpp:216
std::map< int, BoundaryRegionShPtr > BoundaryRegionCollection
Definition: Conditions.h:217
double NekDouble
void Sadd(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Add vector y = alpha + x.
Definition: Vmath.cpp:318
std::map< std::string, std::string > ParamMap
Definition: Filter.h:68
LibUtilities::SessionReaderSharedPtr m_session
Definition: Filter.h:86
void UpdateForce(const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, Array< OneD, NekDouble > &Aeroforces, const NekDouble &time)
InputIterator find(InputIterator first, InputIterator last, InputIterator startingpoint, const EqualityComparable &value)
Definition: StdRegions.hpp:358
FilterFactory & GetFilterFactory()
Definition: Filter.cpp:41
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:140
std::vector< unsigned int > m_boundaryRegionsIdList
ID&#39;s of boundary regions where we want the forces.
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:376
std::shared_ptr< SessionReader > SessionReaderSharedPtr
const BoundaryRegionCollection & GetBoundaryRegions(void) const
Definition: Conditions.h:238
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.cpp:302
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:186
static bool GenerateSeqVector(const std::string &str, std::vector< unsigned int > &out)
Takes a comma-separated compressed string and converts it to entries in a vector. ...
Definition: ParseUtils.cpp:108