Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Functions
ExtractSurface2DCFS.cpp File Reference
#include <cstdio>
#include <cstdlib>
#include <string>
#include <iostream>
#include <iomanip>
#include <MultiRegions/ExpList.h>
#include <MultiRegions/ExpList1D.h>
#include <MultiRegions/ExpList2D.h>
#include <MultiRegions/ExpList2DHomogeneous1D.h>
#include <MultiRegions/ExpList3DHomogeneous1D.h>
#include <MultiRegions/ExpList3DHomogeneous2D.h>
#include <MultiRegions/AssemblyMap/AssemblyMapDG.h>
#include <MultiRegions/DisContField2D.h>
#include <LocalRegions/MatrixKey.h>
#include <LocalRegions/Expansion2D.h>
#include <LocalRegions/Expansion.h>
#include <LibUtilities/BasicUtils/FieldIO.h>
#include <LibUtilities/BasicUtils/NekFactory.hpp>
#include <LibUtilities/BasicUtils/SessionReader.h>
#include <LibUtilities/BasicUtils/SharedArray.hpp>
#include <LibUtilities/Communication/Comm.h>
#include <LibUtilities/Memory/NekMemoryManager.hpp>
#include <MultiRegions/ContField2D.h>
#include <SpatialDomains/MeshGraph2D.h>
#include <SolverUtils/SolverUtilsDeclspec.h>
Include dependency graph for ExtractSurface2DCFS.cpp:

Go to the source code of this file.

Functions

int main (int argc, char *argv[])
 

Function Documentation

int main ( int  argc,
char *  argv[] 
)

Evaluation of the velocity gradient in the cartesian directions Du_x: traceFieldsAdded[10] Du_y: traceFieldsAdded[11] Dv_x: traceFieldsAdded[12] Dv_y: traceFieldsAdded[13]

Definition at line 71 of file ExtractSurface2DCFS.cpp.

References ASSERTL0, Nektar::LibUtilities::ePolyEvenlySpaced, Vmath::Fill(), Nektar::StdRegions::StdExpansion::GetTotPoints(), Nektar::LibUtilities::Import(), m_mu, m_rhoInf, m_Twall, m_uInf, m_vInf, Vmath::Neg(), Vmath::Smul(), Vmath::Vadd(), Vmath::Vcopy(), Vmath::Vdiv(), Vmath::Vmul(), Vmath::Vsqrt(), Vmath::Vsub(), and Vmath::Vvtvp().

72 {
73  string fname = std::string(argv[2]);
74  int fdot = fname.find_last_of('.');
75  if (fdot != std::string::npos)
76  {
77  string ending = fname.substr(fdot);
78 
79  // If .chk or .fld we exchange the extension in the output file.
80  // For all other files (e.g. .bse) we append the extension to avoid
81  // conflicts.
82  if (ending == ".chk" || ending == ".fld")
83  {
84  fname = fname.substr(0,fdot);
85  }
86  }
87 
88  fname = fname + ".txt";
89 
90  int cnt;
91  int id1 = 0;
92  int id2 = 0;
93  int i, j, n, e, b;
94  Array<OneD, NekDouble> auxArray;
95 
96  int nBndEdgePts, nBndEdges, nBndRegions;
97 
98  if (argc < 3)
99  {
100  fprintf(stderr,
101  "Usage: ExtractSurface2DCFS meshfile fieldFile\n");
102  fprintf(stderr,
103  "Extracts a surface from a 2D fld file"
104  "(only for CompressibleFlowSolver and purely 2D .fld files)\n");
105  exit(1);
106  }
107 
109  = LibUtilities::SessionReader::CreateInstance(3, argv);
110 
111  std::string m_ViscosityType;
112 
113  NekDouble m_gamma;
114  NekDouble m_pInf;
118  NekDouble m_wInf;
119  NekDouble m_gasConstant;
121  NekDouble m_mu;
122  NekDouble m_thermalConductivity;
123 
124  int m_spacedim = 2;
125  int nDimensions = m_spacedim;
126  int phys_offset;
127 
128  // Get gamma parameter from session file.
129  ASSERTL0(vSession->DefinesParameter("Gamma"),
130  "Compressible flow sessions must define a Gamma parameter.");
131  vSession->LoadParameter("Gamma", m_gamma, 1.4);
132 
133  // Get E0 parameter from session file.
134  ASSERTL0(vSession->DefinesParameter("pInf"),
135  "Compressible flow sessions must define a pInf parameter.");
136  vSession->LoadParameter("pInf", m_pInf, 101325);
137 
138  // Get rhoInf parameter from session file.
139  ASSERTL0(vSession->DefinesParameter("rhoInf"),
140  "Compressible flow sessions must define a rhoInf parameter.");
141  vSession->LoadParameter("rhoInf", m_rhoInf, 1.225);
142 
143  // Get uInf parameter from session file.
144  ASSERTL0(vSession->DefinesParameter("uInf"),
145  "Compressible flow sessions must define a uInf parameter.");
146  vSession->LoadParameter("uInf", m_uInf, 0.1);
147 
148  // Get vInf parameter from session file.
149  if (m_spacedim == 2 || m_spacedim == 3)
150  {
151  ASSERTL0(vSession->DefinesParameter("vInf"),
152  "Compressible flow sessions must define a vInf parameter"
153  "for 2D/3D problems.");
154  vSession->LoadParameter("vInf", m_vInf, 0.0);
155  }
156 
157  // Get wInf parameter from session file.
158  if (m_spacedim == 3)
159  {
160  ASSERTL0(vSession->DefinesParameter("wInf"),
161  "Compressible flow sessions must define a wInf parameter"
162  "for 3D problems.");
163  vSession->LoadParameter("wInf", m_wInf, 0.0);
164  }
165 
166  vSession->LoadParameter ("GasConstant", m_gasConstant, 287.058);
167  vSession->LoadParameter ("Twall", m_Twall, 300.15);
168  vSession->LoadSolverInfo("ViscosityType", m_ViscosityType, "Constant");
169  vSession->LoadParameter ("mu", m_mu, 1.78e-05);
170  vSession->LoadParameter ("thermalConductivity",
171  m_thermalConductivity, 0.0257);
172 
173  //--------------------------------------------------------------------------
174  // Read in mesh from input file
175  string meshfile(argv[1]);
177  SpatialDomains::MeshGraph::Read(vSession);
178  //--------------------------------------------------------------------------
179 
180  //--------------------------------------------------------------------------
181  // Import field file
182  string fieldFile(argv[2]);
183  vector<LibUtilities::FieldDefinitionsSharedPtr> fieldDef;
184  vector<vector<NekDouble> > fieldData;
185 
186  LibUtilities::Import(fieldFile, fieldDef, fieldData);
187  //--------------------------------------------------------------------------
188 
189  //--------------------------------------------------------------------------
190  // Set up Expansion information
191  vector< vector<LibUtilities::PointsType> > pointsType;
192  for (i = 0; i < fieldDef.size(); ++i)
193  {
194  vector<LibUtilities::PointsType> ptype;
195  for (j = 0; j < 2; ++j)
196  {
197  ptype.push_back(LibUtilities::ePolyEvenlySpaced);
198  }
199  pointsType.push_back(ptype);
200  }
201  graphShPt->SetExpansions(fieldDef, pointsType);
202 
203  //--------------------------------------------------------------------------
204 
205 
206  //--------------------------------------------------------------------------
207  // Define Expansion
208  int nfields = fieldDef[0]->m_fields.size();
211 
212  for(i = 0; i < pFields.num_elements(); i++)
213  {
214  pFields[i] = MemoryManager<
215  MultiRegions::DisContField2D>::AllocateSharedPtr(
216  vSession, graphShPt, vSession->GetVariable(i));
217  }
218 
221  ::AllocateSharedPtr(vSession, graphShPt);
222 
223  Exp[0] = Exp2D;
224 
225  for (i = 1; i < nfields; ++i)
226  {
228  ::AllocateSharedPtr(*Exp2D);
229  }
230 
231  int nSolutionPts = pFields[0]->GetNpoints();
232  int nTracePts = pFields[0]->GetTrace()->GetTotPoints();
233  int nElements = pFields[0]->GetExpSize();
234 
235  Array<OneD, NekDouble> tmp(nSolutionPts, 0.0);
236 
237  Array<OneD, NekDouble> x(nSolutionPts);
238  Array<OneD, NekDouble> y(nSolutionPts);
239  Array<OneD, NekDouble> z(nSolutionPts);
240 
241  Array<OneD, NekDouble> traceX(nTracePts);
242  Array<OneD, NekDouble> traceY(nTracePts);
243  Array<OneD, NekDouble> traceZ(nTracePts);
244 
245  Array<OneD, NekDouble> surfaceX(nTracePts);
246  Array<OneD, NekDouble> surfaceY(nTracePts);
247  Array<OneD, NekDouble> surfaceZ(nTracePts);
248 
249  pFields[0]->GetCoords(x, y, z);
250 
251  pFields[0]->ExtractTracePhys(x, traceX);
252  pFields[0]->ExtractTracePhys(y, traceY);
253  pFields[0]->ExtractTracePhys(z, traceZ);
254  //--------------------------------------------------------------------------
255 
256  //--------------------------------------------------------------------------
257  // Copy data from field file
258  Array<OneD, Array<OneD, NekDouble> > uFields(nfields);
259  Array<OneD, Array<OneD, NekDouble> > traceFields(nfields);
260  Array<OneD, Array<OneD, NekDouble> > surfaceFields(nfields);
261 
262  // Extract the physical values of the solution at the boundaries
263  for (j = 0; j < nfields; ++j)
264  {
265  uFields[j] = Array<OneD, NekDouble>(nSolutionPts, 0.0);
266  traceFields[j] = Array<OneD, NekDouble>(nTracePts, 0.0);
267  surfaceFields[j] = Array<OneD, NekDouble>(nTracePts, 0.0);
268 
269 
270  for (i = 0; i < fieldData.size(); ++i)
271  {
272  Exp[j]->ExtractDataToCoeffs(fieldDef[i], fieldData[i],
273  fieldDef[i]->m_fields[j],
274  Exp[j]->UpdateCoeffs());
275  }
276  Exp[j]->BwdTrans(Exp[j]->GetCoeffs(), Exp[j]->UpdatePhys());
277  Vmath::Vcopy(nSolutionPts, Exp[j]->GetPhys(), 1, uFields[j], 1);
278  pFields[0]->ExtractTracePhys(uFields[j], traceFields[j]);
279  }
280 
281  //Fields to add in the output file
282 
283  int nfieldsAdded = 20;
284  Array<OneD, Array<OneD, NekDouble> > traceFieldsAdded(nfieldsAdded);
285  Array<OneD, Array<OneD, NekDouble> > surfaceFieldsAdded(nfieldsAdded);
286 
287  for (j = 0; j < nfieldsAdded; ++j)
288  {
289  traceFieldsAdded[j] = Array<OneD, NekDouble>(nTracePts, 0.0);
290  surfaceFieldsAdded[j] = Array<OneD, NekDouble>(nTracePts, 0.0);
291  }
292 
293  /******** Evaluation of normals and tangents on the trace *****************
294  * nx -> traceFieldsAdded[0];
295  * ny -> traceFieldsAdded[1];
296  * tx -> traceFieldsAdded[2];
297  * ty -> traceFieldsAdded[3];
298  ***************************************************************************/
299 
300  Array<OneD, Array<OneD, NekDouble> > m_traceNormals (nDimensions);
301  for(i = 0; i < nDimensions; ++i)
302  {
303  m_traceNormals[i] = Array<OneD, NekDouble> (nTracePts, 0.0);
304  }
305  pFields[0]->GetTrace()->GetNormals(m_traceNormals);
306 
307  Array<OneD, Array<OneD, NekDouble> > m_traceTangents (nDimensions);
308  for(i = 0; i < nDimensions; ++i)
309  {
310  m_traceTangents[i] = Array<OneD, NekDouble> (nTracePts, 0.0);
311  }
312 
313 
314  // nx
315  Vmath::Vcopy(nTracePts,
316  &m_traceNormals[0][0], 1,
317  &traceFieldsAdded[0][0], 1);
318 
319  // ny
320  Vmath::Vcopy(nTracePts,
321  &m_traceNormals[1][0], 1,
322  &traceFieldsAdded[1][0], 1);
323 
324  // t_x = - n_y
325  Vmath::Vcopy(nTracePts,
326  &m_traceNormals[1][0], 1,
327  &m_traceTangents[0][0], 1);
328  Vmath::Neg(nTracePts, &m_traceTangents[0][0], 1);
329 
330  Vmath::Vcopy(nTracePts,
331  &m_traceTangents[0][0], 1,
332  &traceFieldsAdded[2][0], 1);
333 
334  // t_y = n_x
335  Vmath::Vcopy(nTracePts,
336  &m_traceNormals[0][0], 1,
337  &m_traceTangents[1][0], 1);
338 
339  Vmath::Vcopy(nTracePts,
340  &m_traceTangents[1][0], 1,
341  &traceFieldsAdded[3][0], 1);
342 
343  /******** Evaluation of the pressure ***************************************
344  * P = (E-1/2.*rho.*((rhou./rho).^2+(rhov./rho).^2))*(gamma - 1);
345  * P -> traceFieldsAdded[4];
346  ***************************************************************************/
347 
348  Array<OneD, NekDouble> pressure(nSolutionPts, 0.0);
349  NekDouble gammaMinusOne = m_gamma - 1.0;
350 
351  for (i = 0; i < m_spacedim; i++)
352  {
353  Vmath::Vmul(nSolutionPts,
354  &uFields[i + 1][0], 1,
355  &uFields[i + 1][0], 1,
356  &tmp[0],1);
357 
358 
359  Vmath::Smul(nSolutionPts, 0.5,
360  &tmp[0], 1,
361  &tmp[0], 1);
362 
363  Vmath::Vadd(nSolutionPts,
364  &pressure[0], 1,
365  &tmp[0], 1,
366  &pressure[0], 1);
367  }
368 
369  Vmath::Vdiv(nSolutionPts,
370  &pressure[0], 1,
371  &uFields[0][0], 1,
372  &pressure[0],1);
373 
374  Vmath::Vsub(nSolutionPts,
375  &uFields[nfields - 1][0], 1,
376  &pressure[0], 1,
377  &pressure[0],1);
378 
379  Vmath::Smul(nSolutionPts, gammaMinusOne,
380  &pressure[0], 1,
381  &pressure[0], 1);
382 
383  // Extract trace
384  pFields[0]->ExtractTracePhys(pressure, traceFieldsAdded[4]);
385 
386  /******** Evaluation of the temperature ************************************
387  * T = P/(R*rho);
388  * T -> traceFieldsAdded[5];
389  ***************************************************************************/
390 
391  Array<OneD, NekDouble> temperature(nSolutionPts, 0.0);
392 
393  Vmath::Vdiv(nSolutionPts,
394  &pressure[0], 1,
395  &uFields[0][0], 1,
396  &temperature[0],1);
397 
398  NekDouble GasConstantInv = 1.0/m_gasConstant;
399  Vmath::Smul(nSolutionPts, GasConstantInv,
400  &temperature[0], 1,
401  &temperature[0], 1);
402 
403  // Extract trace
404  pFields[0]->ExtractTracePhys(temperature, traceFieldsAdded[5]);
405 
406  /*** Evaluation of the temperature gradient in the normal direction ********
407  * DT_n -> traceFieldsAdded[6]
408  ***************************************************************************/
409 
410  Array<OneD, Array<OneD, NekDouble> > Dtemperature(nDimensions);
411  Array<OneD, Array<OneD, NekDouble> > traceDtemperature(nDimensions);
412 
413  for (i = 0; i < nDimensions; ++ i)
414  {
415  Dtemperature[i] = Array<OneD, NekDouble>(nSolutionPts, 0.0);
416  traceDtemperature[i] = Array<OneD, NekDouble>(nTracePts, 0.0);
417  }
418 
419  for (i = 0; i < nDimensions; ++ i)
420  {
421  for (n = 0; n < nElements; n++)
422  {
423  phys_offset = pFields[0]->GetPhys_Offset(n);
424 
425  pFields[i]->GetExp(n)->PhysDeriv(
426  i, temperature + phys_offset,
427  auxArray = Dtemperature[i] + phys_offset);
428  }
429  // Extract trace
430  pFields[0]->ExtractTracePhys(Dtemperature[i], traceDtemperature[i]);
431  }
432 
433  for(i = 0; i < nDimensions; ++i)
434  {
435  Vmath::Vmul(nTracePts,
436  &m_traceNormals[i][0], 1,
437  &traceDtemperature[i][0], 1,
438  &tmp[0],1);
439 
440  Vmath::Vadd(nTracePts,
441  &traceFieldsAdded[6][0], 1,
442  &tmp[0], 1,
443  &traceFieldsAdded[6][0], 1);
444  }
445 
446  /*** Evaluation of the pressure gradient ***********************************
447  * DP_t -> traceFieldsAdded[7] tangent direction
448  * DP_x -> traceFieldsAdded[8]
449  * DP_y -> traceFieldsAdded[9]
450  ***************************************************************************/
451 
452  Array<OneD, Array<OneD, NekDouble> > Dpressure(nDimensions);
453  Array<OneD, Array<OneD, NekDouble> > traceDpressure(nDimensions);
454 
455  for (i = 0; i < nDimensions; ++ i)
456  {
457  Dpressure[i] = Array<OneD, NekDouble>(nSolutionPts, 0.0);
458  traceDpressure[i] = Array<OneD, NekDouble>(nTracePts, 0.0);
459  }
460 
461  for (i = 0; i < nDimensions; ++ i)
462  {
463  for (n = 0; n < nElements; n++)
464  {
465  phys_offset = pFields[0]->GetPhys_Offset(n);
466 
467  pFields[i]->GetExp(n)->PhysDeriv(
468  i, pressure + phys_offset,
469  auxArray = Dpressure[i] + phys_offset);
470  }
471  // Extract trace
472  pFields[0]->ExtractTracePhys(Dpressure[i], traceDpressure[i]);
473  }
474 
475  // Dp_t
476  for(i = 0; i < nDimensions; ++i)
477  {
478  Vmath::Vmul(nTracePts,
479  &m_traceTangents[i][0], 1,
480  &traceDpressure[i][0], 1,
481  &tmp[0],1);
482 
483  Vmath::Vadd(nTracePts,
484  &traceFieldsAdded[7][0], 1,
485  &tmp[0], 1,
486  &traceFieldsAdded[7][0], 1);
487  }
488 
489  // Dp_x
490  Vmath::Vcopy(nTracePts,
491  &traceDpressure[0][0], 1,
492  &traceFieldsAdded[8][0], 1);
493 
494  // Dp_y
495  Vmath::Vcopy(nTracePts,
496  &traceDpressure[1][0], 1,
497  &traceFieldsAdded[9][0], 1);
498 
499  /** Evaluation of the velocity gradient in the cartesian directions
500  * Du_x: traceFieldsAdded[10]
501  * Du_y: traceFieldsAdded[11]
502  * Dv_x: traceFieldsAdded[12]
503  * Dv_y: traceFieldsAdded[13]
504  **/
505  Array<OneD, Array<OneD, Array<OneD, NekDouble> > > Dvelocity(nDimensions);
506  Array<OneD, Array<OneD, Array<OneD, NekDouble> > > traceDvelocity(nDimensions);
507  Array<OneD, Array<OneD, NekDouble> > velocity(nDimensions);
508 
509  for (i = 0; i < nDimensions; ++ i)
510  {
511  Dvelocity[i] = Array<OneD, Array<OneD, NekDouble> >(nDimensions);
512  traceDvelocity[i] = Array<OneD, Array<OneD, NekDouble> >(nDimensions);
513  velocity[i] = Array<OneD, NekDouble>(nSolutionPts, 0.0);
514 
515  Vmath::Vdiv(nSolutionPts, uFields[i+1], 1, uFields[0], 1,
516  velocity[i], 1);
517 
518  for (j = 0; j < nDimensions; ++j)
519  {
520  Dvelocity[i][j] = Array<OneD, NekDouble>(nSolutionPts, 0.0);
521  traceDvelocity[i][j] = Array<OneD, NekDouble>(nTracePts, 0.0);
522  }
523  }
524 
525  for (i = 0; i < nDimensions; ++i)
526  {
527  for (j = 0; j < nDimensions; ++j)
528  {
529  for (n = 0; n < nElements; n++)
530  {
531  phys_offset = pFields[0]->GetPhys_Offset(n);
532 
533  pFields[i]->GetExp(n)->PhysDeriv(
534  j, velocity[i] + phys_offset,
535  auxArray = Dvelocity[i][j] + phys_offset);
536  }
537 
538  // Extract trace
539  pFields[0]->ExtractTracePhys(Dvelocity[i][j], traceDvelocity[i][j]);
540  }
541  }
542 
543  Vmath::Vcopy(nTracePts,
544  &traceDvelocity[0][0][0], 1,
545  &traceFieldsAdded[10][0], 1);
546  Vmath::Vcopy(nTracePts,
547  &traceDvelocity[0][1][0], 1,
548  &traceFieldsAdded[11][0], 1);
549  Vmath::Vcopy(nTracePts,
550  &traceDvelocity[1][0][0], 1,
551  &traceFieldsAdded[12][0], 1);
552  Vmath::Vcopy(nTracePts,
553  &traceDvelocity[1][1][0], 1,
554  &traceFieldsAdded[13][0], 1);
555 
556 
557  /*** Evaluation of shear stresses ******************************************
558  * tau_xx -> traceFieldsAdded[14]
559  * tau_yy -> traceFieldsAdded[15]
560  * tau_xy -> traceFieldsAdded[16]
561  ***************************************************************************/
562 
563  // Stokes hypotesis
564  const NekDouble lambda = -2.0/3.0;
565 
566  // Auxiliary variables
567  Array<OneD, NekDouble > mu (nSolutionPts, 0.0);
568  Array<OneD, NekDouble > mu2 (nSolutionPts, 0.0);
569  Array<OneD, NekDouble > divVel(nSolutionPts, 0.0);
570 
571  // Variable viscosity through the Sutherland's law
572  if (m_ViscosityType == "Variable")
573  {
574  NekDouble mu_star = m_mu;
575  NekDouble T_star = m_pInf / (m_rhoInf * m_gasConstant);
576  NekDouble ratio;
577 
578  for (int i = 0; i < nSolutionPts; ++i)
579  {
580  ratio = temperature[i] / T_star;
581  mu[i] = mu_star * ratio * sqrt(ratio) *
582  (T_star + 110.0) / (temperature[i] + 110.0);
583  }
584  }
585  else
586  {
587  Vmath::Fill(nSolutionPts, m_mu, &mu[0], 1);
588  }
589 
590  // Computing diagonal terms of viscous stress tensor
591  Array<OneD, Array<OneD, NekDouble> > temp(m_spacedim);
592  Array<OneD, Array<OneD, NekDouble> > Sgg(m_spacedim);
593 
594  // mu2 = 2 * mu
595  Vmath::Smul(nSolutionPts, 2.0, &mu[0], 1, &mu2[0], 1);
596 
597  // Velocity divergence
598  Vmath::Vadd(nSolutionPts, &divVel[0], 1,
599  &Dvelocity[0][0][0], 1, &divVel[0], 1);
600  Vmath::Vadd(nSolutionPts, &divVel[0], 1,
601  &Dvelocity[1][1][0], 1, &divVel[0], 1);
602 
603  // Velocity divergence scaled by lambda * mu
604  Vmath::Smul(nSolutionPts, lambda, &divVel[0], 1, &divVel[0], 1);
605  Vmath::Vmul(nSolutionPts, &mu[0], 1, &divVel[0], 1, &divVel[0], 1);
606 
607  // Diagonal terms of viscous stress tensor (Sxx, Syy)
608  // Sjj = 2 * mu * du_j/dx_j - (2 / 3) * mu * sum_j(du_j/dx_j)
609  for (j = 0; j < m_spacedim; ++j)
610  {
611  temp[j] = Array<OneD, NekDouble>(nSolutionPts, 0.0);
612  Sgg[j] = Array<OneD, NekDouble>(nSolutionPts, 0.0);
613 
614  Vmath::Vmul(nSolutionPts, &mu2[0], 1, &Dvelocity[j][j][0], 1,
615  &temp[j][0], 1);
616 
617  Vmath::Vadd(nSolutionPts, &temp[j][0], 1, &divVel[0], 1, &Sgg[j][0], 1);
618  }
619 
620  // Extra diagonal terms of viscous stress tensor (Sxy = Syx)
621  Array<OneD, NekDouble > Sxy(nSolutionPts, 0.0);
622 
623  // Sxy = (du/dy + dv/dx)
624  Vmath::Vadd(nSolutionPts, &Dvelocity[0][1][0], 1,
625  &Dvelocity[1][0][0], 1, &Sxy[0], 1);
626 
627  // Sxy = mu * (du/dy + dv/dx)
628  Vmath::Vmul(nSolutionPts, &mu[0], 1, &Sxy[0], 1, &Sxy[0], 1);
629 
630  pFields[0]->ExtractTracePhys(Sgg[0], traceFieldsAdded[14]);
631  pFields[0]->ExtractTracePhys(Sgg[1], traceFieldsAdded[15]);
632  pFields[0]->ExtractTracePhys(Sxy, traceFieldsAdded[16]);
633 
634  /*** Evaluation of the shear stress in tangent direction *******************
635  * tau_t -> traceFieldsAdded[17]
636  ***************************************************************************/
637  Array<OneD, NekDouble > sigma_diff (nTracePts, 0.0);
638  Array<OneD, NekDouble > cosTeta (nTracePts, 0.0);
639  Array<OneD, NekDouble > sinTeta (nTracePts, 0.0);
640  Array<OneD, NekDouble > cos2Teta (nTracePts, 0.0);
641  Array<OneD, NekDouble > sin2Teta (nTracePts, 0.0);
642  Array<OneD, NekDouble > tau_t (nTracePts, 0.0);
643 
644  Array<OneD, NekDouble > tmpTeta (nTracePts, 0.0);
645 
646  // cos(teta) = nx
647  Vmath::Vcopy(nTracePts, &m_traceNormals[0][0], 1, &cosTeta[0], 1);
648 
649  // sin(teta) = ny
650  Vmath::Vcopy(nTracePts, &m_traceNormals[1][0], 1, &sinTeta[0], 1);
651 
652  // sigma_diff = sigma_x - sigma_y
653  Vmath::Vsub(nTracePts, &traceFieldsAdded[14][0], 1,
654  &traceFieldsAdded[15][0], 1, &sigma_diff[0], 1);
655 
656  // sin(2*teta)
657  Vmath::Vmul(nTracePts, &cosTeta[0], 1, &sinTeta[0], 1, &tmpTeta[0], 1);
658  Vmath::Smul(nTracePts, 2.0, &tmpTeta[0], 1, &sin2Teta[0], 1);
659 
660  // cos(2*teta)
661  Vmath::Vmul(nTracePts, &cosTeta[0], 1, &cosTeta[0], 1, &cos2Teta[0], 1);
662  Vmath::Vmul(nTracePts, &sinTeta[0], 1, &sinTeta[0], 1, &tmpTeta[0], 1);
663  Vmath::Vsub(nTracePts, &cos2Teta[0], 1, &tmpTeta[0], 1, &cos2Teta[0], 1);
664 
665  // tau_t = -0.5*sigma_diff * sin(2*teta) + tau_xy * cos(2*teta)
666  Vmath::Smul(nTracePts, -0.5, &sigma_diff[0], 1, &sigma_diff[0], 1);
667  Vmath::Vmul(nTracePts, &sigma_diff[0], 1, &sin2Teta[0], 1, &tau_t[0], 1);
668  Vmath::Vmul(nTracePts, &traceFieldsAdded[16][0], 1, &cos2Teta[0], 1,
669  &tmpTeta[0], 1);
670  Vmath::Vadd(nTracePts, &tau_t[0], 1, &tmpTeta[0], 1, &tau_t[0], 1);
671 
672  Vmath::Vcopy(nTracePts, &tau_t[0], 1, &traceFieldsAdded[17][0], 1);
673 
674  /*** Evaluation of dinamic viscosity ***************************************
675  * mu -> traceFieldsAdded[18]
676  ***************************************************************************/
677 
678  pFields[0]->ExtractTracePhys(mu, traceFieldsAdded[18]);
679 
680  /*** Evaluation of Mach number *********************************************
681  * M -> traceFieldsAdded[18]
682  ***************************************************************************/
683  NekDouble gamma = m_gamma;
684 
685  // Speed of sound
686  Array<OneD, NekDouble> soundspeed(nSolutionPts, 0.0);
687 
688  Vmath::Vdiv (nSolutionPts, pressure, 1, uFields[0], 1, soundspeed, 1);
689  Vmath::Smul (nSolutionPts, gamma, soundspeed, 1, soundspeed, 1);
690  Vmath::Vsqrt(nSolutionPts, soundspeed, 1, soundspeed, 1);
691 
692  // Mach
693  Array<OneD, NekDouble> mach(nSolutionPts, 0.0);
694 
695  for (int i = 0; i < m_spacedim; ++i)
696  {
697  Vmath::Vvtvp(nSolutionPts, uFields[i + 1], 1, uFields[i + 1], 1,
698  mach, 1, mach, 1);
699  }
700 
701  Vmath::Vdiv(nSolutionPts, mach, 1, uFields[0], 1, mach, 1);
702  Vmath::Vdiv(nSolutionPts, mach, 1, uFields[0], 1, mach, 1);
703  Vmath::Vsqrt(nSolutionPts, mach, 1, mach, 1);
704  Vmath::Vdiv(nSolutionPts, mach, 1, soundspeed, 1, mach, 1);
705 
706  pFields[0]->ExtractTracePhys(mach, traceFieldsAdded[19]);
707 
708  /**************************************************************************/
709  // Extract coordinates
710 
711  if (pFields[0]->GetBndCondExpansions().num_elements())
712  {
713  id1 = 0;
714  cnt = 0;
715  nBndRegions = pFields[0]->GetBndCondExpansions().num_elements();
716  for (b = 0; b < nBndRegions; ++b)
717  {
718  nBndEdges = pFields[0]->GetBndCondExpansions()[b]->GetExpSize();
719  for (e = 0; e < nBndEdges; ++e)
720  {
721  nBndEdgePts = pFields[0]->
722  GetBndCondExpansions()[b]->GetExp(e)->GetNumPoints(0);
723 
724  id2 = pFields[0]->GetTrace()->
725  GetPhys_Offset(pFields[0]->GetTraceMap()->
726  GetBndCondTraceToGlobalTraceMap(cnt++));
727 
728  if (pFields[0]->GetBndConditions()[b]->
729  GetUserDefined() == "WallViscous" ||
730  pFields[0]->GetBndConditions()[b]->
731  GetUserDefined() == "WallAdiabatic" ||
732  pFields[0]->GetBndConditions()[b]->
733  GetUserDefined() == "Wall")
734  {
735  Vmath::Vcopy(nBndEdgePts, &traceX[id2], 1,
736  &surfaceX[id1], 1);
737 
738  Vmath::Vcopy(nBndEdgePts, &traceY[id2], 1,
739  &surfaceY[id1], 1);
740 
741  Vmath::Vcopy(nBndEdgePts, &traceZ[id2], 1,
742  &surfaceZ[id1], 1);
743 
744  id1 += nBndEdgePts;
745  }
746  }
747  }
748  }
749 
750  // Extract fields
751  if (pFields[0]->GetBndCondExpansions().num_elements())
752  {
753  for (j = 0; j < nfields; ++j)
754  {
755  id1 = 0;
756  cnt = 0;
757  nBndRegions = pFields[j]->GetBndCondExpansions().num_elements();
758  for (b = 0; b < nBndRegions; ++b)
759  {
760  nBndEdges = pFields[j]->GetBndCondExpansions()[b]->GetExpSize();
761  for (e = 0; e < nBndEdges; ++e)
762  {
763  nBndEdgePts = pFields[j]->
764  GetBndCondExpansions()[b]->GetExp(e)->GetNumPoints(0);
765 
766  id2 = pFields[j]->GetTrace()->
767  GetPhys_Offset(pFields[j]->GetTraceMap()->
768  GetBndCondTraceToGlobalTraceMap(cnt++));
769 
770  if (pFields[j]->GetBndConditions()[b]->
771  GetUserDefined() == "WallViscous" ||
772  pFields[j]->GetBndConditions()[b]->
773  GetUserDefined() == "WallAdiabatic" ||
774  pFields[j]->GetBndConditions()[b]->
775  GetUserDefined() == "Wall")
776  {
777  Vmath::Vcopy(nBndEdgePts, &traceFields[j][id2], 1,
778  &surfaceFields[j][id1], 1);
779 
780  id1 += nBndEdgePts;
781  }
782  }
783  }
784  }
785  }
786 
787  // Extract fields added
788  if (pFields[0]->GetBndCondExpansions().num_elements())
789  {
790  for (j = 0; j < nfieldsAdded; ++j)
791  {
792  id1 = 0;
793  cnt = 0;
794  nBndRegions = pFields[0]->GetBndCondExpansions().num_elements();
795  for (b = 0; b < nBndRegions; ++b)
796  {
797  nBndEdges = pFields[0]->GetBndCondExpansions()[b]->GetExpSize();
798  for (e = 0; e < nBndEdges; ++e)
799  {
800  nBndEdgePts = pFields[0]->
801  GetBndCondExpansions()[b]->GetExp(e)->GetNumPoints(0);
802 
803  id2 = pFields[0]->GetTrace()->
804  GetPhys_Offset(pFields[0]->GetTraceMap()->
805  GetBndCondTraceToGlobalTraceMap(cnt++));
806 
807  if (pFields[0]->GetBndConditions()[b]->
808  GetUserDefined() == "WallViscous" ||
809  pFields[0]->GetBndConditions()[b]->
810  GetUserDefined() == "WallAdiabatic" ||
811  pFields[0]->GetBndConditions()[b]->
812  GetUserDefined() == "Wall")
813  {
814  Vmath::Vcopy(nBndEdgePts, &traceFieldsAdded[j][id2], 1,
815  &surfaceFieldsAdded[j][id1], 1);
816 
817  id1 += nBndEdgePts;
818  }
819  }
820  }
821  }
822  }
823  //===================================================================================================
824  //===================================================================================================
825  //===================================================================================================
826  std::string vEquation = vSession->GetSolverInfo("EQType");
827 
829  BndExp = pFields[0]->GetBndCondExpansions();
831 
832  NekDouble Fxp(0.0);
833  NekDouble Fyp(0.0);
834  NekDouble Fxv(0.0);
835  NekDouble Fyv(0.0);
836  NekDouble Sref(0.0);
837 
838  int GlobalIndex(0);
839 
840  for(int i = 0; i < BndExp[0]->GetExpSize(); ++i)
841  {
842  bc = BndExp[0]->GetExp(i)->as<LocalRegions::Expansion1D> ();
843 
844  int nbc = bc->GetTotPoints();
845 
846  Array<OneD, NekDouble> nxOnBnd(nbc,0.0);
847  Array<OneD, NekDouble> nyOnBnd(nbc,0.0);
848  Array<OneD, NekDouble> txOnBnd(nbc,0.0);
849  Array<OneD, NekDouble> tyOnBnd(nbc,0.0);
850  // Array<OneD, NekDouble> CoeffAero(nbc,0.0);
851  // Array<OneD, NekDouble> tmp(nbc,0.0);
852 
853  Array<OneD, NekDouble> drag_p(nbc,0.0);
854  Array<OneD, NekDouble> lift_p(nbc,0.0);
855  Array<OneD, NekDouble> PressurOnBnd(nbc,0.0);
856 
857  Array<OneD, NekDouble> drag_v(nbc,0.0);
858  Array<OneD, NekDouble> lift_v(nbc,0.0);
859  Array<OneD, NekDouble> ShearStressOnBnd(nbc,0.0);
860 
861  Array<OneD, NekDouble> Unity(nbc,1.0);
862 
863  for(int j = 0; j < nbc; ++j)
864  {
865 
866  nxOnBnd[j] = surfaceFieldsAdded[0][GlobalIndex];
867  nyOnBnd[j] = surfaceFieldsAdded[1][GlobalIndex];
868  txOnBnd[j] = surfaceFieldsAdded[2][GlobalIndex];
869  tyOnBnd[j] = surfaceFieldsAdded[3][GlobalIndex];
870 
871  PressurOnBnd[j] = surfaceFieldsAdded[4][GlobalIndex];
872 
873  if (vEquation == "NavierStokesCFE")
874  {
875  ShearStressOnBnd[j] = surfaceFieldsAdded[17][GlobalIndex];
876  }
877 
878  // CoeffAero[j] = surfaceFields[0][GlobalIndex];
879  // tmp[j] = surfaceFields[1][GlobalIndex]*surfaceFields[1][GlobalIndex];
880  // tmp[j] = tmp[j] + surfaceFields[2][GlobalIndex]*surfaceFields[2][GlobalIndex];
881  // tmp[j] = sqrt(tmp[j]);
882  // CoeffAero[j] = CoeffAero[j]*tmp[j];
883  // CoeffAero[j] = 1.0/CoeffAero[j];
884  //
885  // PressurOnBnd[j] = CoeffAero[j]*surfaceFieldsAdded[4][GlobalIndex];
886  //
887  // cout << "CoeffAero = " << CoeffAero[j] << endl;
888 
889  GlobalIndex++;
890  }
891 
892  Vmath::Vmul(nbc,PressurOnBnd,1,nxOnBnd,1, drag_p,1);
893  Vmath::Vmul(nbc,PressurOnBnd,1,nyOnBnd,1, lift_p,1);
894 
895  // Vmath::Vmul(nbc,drag_p,1,CoeffAero,1, drag_p,1);
896  // Vmath::Vmul(nbc,lift_p,1,CoeffAero,1, lift_p,1);
897 
898  Fxp += bc->Integral(drag_p);
899  Fyp += bc->Integral(lift_p);
900 
901  if (vEquation == "NavierStokesCFE")
902  {
903  Vmath::Vmul(nbc,ShearStressOnBnd,1,txOnBnd,1, drag_v,1);
904  Vmath::Vmul(nbc,ShearStressOnBnd,1,tyOnBnd,1, lift_v,1);
905 
906  // Vmath::Vdiv(nbc,drag_v,1,CoeffAero,1, drag_v,1);
907  // Vmath::Vdiv(nbc,lift_v,1,CoeffAero,1, lift_v,1);
908 
909  Fxv += bc->Integral(drag_v);
910  Fyv += bc->Integral(lift_v);
911  }
912 
913  Sref += bc->Integral(Unity);
914 
915  }
916 
917  cout << "\n Sref = " << Sref << endl;
918  Fxp = Fxp/Sref;
919  Fyp = Fyp/Sref;
920  Fxv = Fxv/Sref;
921  Fyv = Fyv/Sref;
922  cout << " Pressure drag (Fxp) = " << Fxp << endl;
923  cout << " Pressure lift (Fyp) = " << Fyp << endl;
924  cout << " Viscous drag (Fxv) = " << Fxv << endl;
925  cout << " Viscous lift (Fyv) = " << Fyv << endl;
926  cout << "\n ==> Total drag = " << Fxp+Fxv << endl;
927  cout << " ==> Total lift = " << Fyp+Fyv << "\n" << endl;
928 
929  //===================================================================================================
930  //===================================================================================================
931  //===================================================================================================
932 
933 
934  // Print the surface coordinates and the surface solution in a .txt file
935  ofstream outfile;
936  outfile.open(fname.c_str());
937  outfile << "% x[m] " << " \t"
938  << "y[m] " << " \t"
939  << "z[m] " << " \t"
940  << "nx[] " << " \t"
941  << "ny[] " << " \t"
942  << "tx[] " << " \t"
943  << "ty[] " << " \t"
944  << "rho[kg/m^3] " << " \t"
945  << "rhou[kg/(m^2 s)] " << " \t"
946  << "rhov[kg/(m^2 s)] " << " \t"
947  << "E[Pa] " << " \t"
948  << "p[Pa] " << " \t"
949  << "T[k] " << " \t"
950  << "dT/dn[k/m] " << " \t"
951  << "dp/dT[Pa/m] " << " \t"
952  << "dp/dx[Pa/m] " << " \t"
953  << "dp/dy[Pa/m] " << " \t"
954  << "du/dx[s^-1] " << " \t"
955  << "du/dy[s^-1] " << " \t"
956  << "dv/dx[s^-1] " << " \t"
957  << "dv/dy[s^-1] " << " \t"
958  << "tau_xx[Pa] " << " \t"
959  << "tau_yy[Pa] " << " \t"
960  << "tau_xy[Pa] " << " \t"
961  << "tau_t[Pa] " << " \t"
962  << "mu[Pa s] " << " \t"
963  << "M[] " << " \t"
964  // << "Fxp " << " \t"
965  << endl;
966  for (i = 0; i < id1; ++i)
967  {
968  outfile << scientific
969  << setw (17)
970  << setprecision(16)
971  << surfaceX[i] << " \t "
972  << surfaceY[i] << " \t "
973  << surfaceZ[i] << " \t "
974  << surfaceFieldsAdded[0][i] << " \t "
975  << surfaceFieldsAdded[1][i] << " \t "
976  << surfaceFieldsAdded[2][i] << " \t "
977  << surfaceFieldsAdded[3][i] << " \t "
978  << surfaceFields[0][i] << " \t "
979  << surfaceFields[1][i] << " \t "
980  << surfaceFields[2][i] << " \t "
981  << surfaceFields[3][i] << " \t "
982  << surfaceFieldsAdded[4][i] << " \t "
983  << surfaceFieldsAdded[5][i] << " \t "
984  << surfaceFieldsAdded[6][i] << " \t "
985  << surfaceFieldsAdded[7][i] << " \t "
986  << surfaceFieldsAdded[8][i] << " \t "
987  << surfaceFieldsAdded[9][i] << " \t "
988  << surfaceFieldsAdded[10][i] << " \t "
989  << surfaceFieldsAdded[11][i] << " \t "
990  << surfaceFieldsAdded[12][i] << " \t "
991  << surfaceFieldsAdded[13][i] << " \t "
992  << surfaceFieldsAdded[14][i] << " \t "
993  << surfaceFieldsAdded[15][i] << " \t "
994  << surfaceFieldsAdded[16][i] << " \t "
995  << surfaceFieldsAdded[17][i] << " \t "
996  << surfaceFieldsAdded[18][i] << " \t "
997  << surfaceFieldsAdded[19][i] << " \t "
998  // << Fxp << " \t "
999  << endl;
1000  }
1001  outfile << endl << endl;
1002  outfile.close();
1003 
1004  return 0;
1005 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
void Vsqrt(int n, const T *x, const int incx, T *y, const int incy)
sqrt y = sqrt(x)
Definition: Vmath.cpp:394
General purpose memory allocation routines with the ability to allocate from thread specific memory p...
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.cpp:46
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:428
NekDouble m_vInf
NekDouble m_Twall
NekDouble m_uInf
void Vdiv(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:227
NekDouble m_mu
boost::shared_ptr< SessionReader > SessionReaderSharedPtr
Definition: MeshPartition.h:51
void Import(const std::string &infilename, std::vector< FieldDefinitionsSharedPtr > &fielddefs, std::vector< std::vector< NekDouble > > &fielddata, FieldMetaDataMap &fieldinfomap, const Array< OneD, int > ElementiDs)
Imports an FLD file.
Definition: FieldIO.cpp:115
NekDouble m_rhoInf
1D Evenly-spaced points using Lagrange polynomial
Definition: PointsType.h:63
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:141
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:199
boost::shared_ptr< StdExpansion1D > StdExpansion1DSharedPtr
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:382
double NekDouble
boost::shared_ptr< ExpList2D > ExpList2DSharedPtr
Shared pointer to an ExpList2D object.
Definition: ExpList2D.h:49
void Vsub(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Subtract vector z = x-y.
Definition: Vmath.cpp:329
boost::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:442
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1047
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:285
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:169