Nektar++
Functions | Variables
CompressibleBL.cpp File Reference
#include <cstdio>
#include <cstdlib>
#include <string>
#include <iostream>
#include <iomanip>
#include <cmath>
#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/ExpList3D.h>
#include <MultiRegions/AssemblyMap/AssemblyMapDG.h>
#include <MultiRegions/ContField2D.h>
#include <MultiRegions/ContField3D.h>
#include <MultiRegions/DisContField2D.h>
#include <MultiRegions/DisContField3D.h>
#include <LocalRegions/MatrixKey.h>
#include <LocalRegions/Expansion2D.h>
#include <LocalRegions/Expansion3D.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 <SpatialDomains/MeshGraph2D.h>
#include <SpatialDomains/MeshGraph3D.h>
#include <SolverUtils/SolverUtilsDeclspec.h>
Include dependency graph for CompressibleBL.cpp:

Go to the source code of this file.

Functions

void COMPBL (Array< OneD, NekDouble > v, Array< OneD, NekDouble > dv)
 
void RK4 (Array< OneD, NekDouble > y, Array< OneD, NekDouble > dydx, int n, NekDouble x, NekDouble h, Array< OneD, NekDouble > yout)
 
void RKDUMB (Array< OneD, NekDouble > vstart, int nvar, NekDouble x1, NekDouble x2, int m_xpoints, Array< OneD, NekDouble > xx, Array< OneD, Array< OneD, NekDouble > > y)
 
void OUTPUT (int m_xpoints, Array< OneD, NekDouble > xx, Array< OneD, Array< OneD, NekDouble > > ff, int nQuadraturePts, Array< OneD, NekDouble > x_QuadraturePts, Array< OneD, NekDouble > y_QuadraturePts, Array< OneD, NekDouble > u_QuadraturePts, Array< OneD, NekDouble > v_QuadraturePts, Array< OneD, NekDouble > rho_QuadraturePts, Array< OneD, NekDouble > T_QuadraturePts)
 
int main (int argc, char *argv[])
 

Variables

NekDouble m_Re
 
NekDouble m_Mach
 
NekDouble L
 
NekDouble m_Tinf
 
NekDouble m_Suth
 
NekDouble m_Tw
 
NekDouble m_Twall
 
NekDouble m_Gamma
 
NekDouble m_Pr
 
NekDouble m_long
 
NekDouble m_uInf
 
NekDouble m_rhoInf
 
NekDouble m_R
 
NekDouble m_vInf
 
NekDouble m_mu
 
NekDouble m_To = 273.11
 
const int m_xpoints = 1000001
 
const NekDouble Nvisc = 1
 
const NekDouble Omega = 1
 
const NekDouble etamax = 10.0
 
const NekDouble errtol = 1e-5
 

Function Documentation

void COMPBL ( Array< OneD, NekDouble v,
Array< OneD, NekDouble dv 
)

Calculate the compressible boundary layer using the similarity solution

Definition at line 105 of file CompressibleBL.cpp.

References m_Gamma, m_Mach, m_Pr, m_Suth, m_Twall, Nvisc, and Omega.

Referenced by OUTPUT(), RK4(), and RKDUMB().

107 {
108  NekDouble c, dcdg, cp;
109 
110  if (Nvisc == 1)
111  {
112  c = sqrt(v[3]) * (1.0 + m_Suth) / (v[3] + m_Suth);
113  dcdg = 1.0 / (2. * sqrt(v[3])) - sqrt(v[3]) / (v[3]+m_Suth);
114  dcdg = dcdg * (1.0 + m_Suth) / (v[3] + m_Suth);
115  cp = dcdg * v[4];
116  }
117  if (Nvisc == 2)
118  {
119  c = pow(v[3], (Omega-1.0));
120  dcdg = (Omega - 1.0) * pow(v[3], (Omega - 2.0));
121  cp = dcdg * v[4];
122  }
123  if (Nvisc == 3)
124  {
125  c = sqrt(m_Twall) * (1.0 + m_Suth) / (m_Suth + m_Twall);
126  cp = 0.0;
127  }
128 
129  dv[0] = v[1];
130  dv[1] = v[2];
131  dv[2] = - v[2] * (cp + v[0]) / c;
132  dv[3] = v[4];
133  dv[4] = - v[4] * (cp + m_Pr * v[0]) / c -
134  m_Pr * (m_Gamma - 1.0) * pow(m_Mach, 2.0) *
135  pow(v[2], 2);
136 }
NekDouble m_Mach
const NekDouble Omega
NekDouble m_Twall
NekDouble m_Pr
double NekDouble
NekDouble m_Suth
const NekDouble Nvisc
NekDouble m_Gamma
int main ( int  argc,
char *  argv[] 
)

Calculate the compressible boundary layer solution for a given 3D mesh and dump the solution into a .rst file.

Definition at line 360 of file CompressibleBL.cpp.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL0, Nektar::LibUtilities::SessionReader::CreateInstance(), errtol, etamax, m_Gamma, m_long, m_Mach, m_mu, m_Pr, m_R, m_Re, m_rhoInf, m_Suth, m_Tinf, m_Tw, m_Twall, m_uInf, m_vInf, m_xpoints, OUTPUT(), Nektar::SpatialDomains::MeshGraph::Read(), RKDUMB(), Vmath::Vcopy(), and Nektar::LibUtilities::Write().

361 {
362  // Variable initialisation
363  int nmax = 5;
364  int maxit = 10;
365 
366  string opt;
367 
368  int i, j, numModes;
369 
372  Array<OneD, NekDouble> parameter(9, 0.0);
373 
374  for (int i = 0; i < nmax; i++)
375  {
377  }
378 
379  Array<OneD, NekDouble > vstart(nmax, 0.0);
385 
386  NekDouble al11, al21, al12, al22, det;
387 
388  // Reading the session file
390  = LibUtilities::SessionReader::CreateInstance(argc, argv);
391 
392  // Read in mesh from input file and create an object of class MeshGraph
394  = SpatialDomains::MeshGraph::Read(vSession);
395 
396  int expdim = graphShPt->GetMeshDimension();
397 
398  int nElements, nQuadraturePts = 0;
399  Array<OneD, NekDouble> x_QuadraturePts;
400  Array<OneD, NekDouble> y_QuadraturePts;
401  Array<OneD, NekDouble> z_QuadraturePts;
402 
403  if (expdim == 2)
404  {
406  ::AllocateSharedPtr(vSession);
407 
410  ::AllocateSharedPtr(vSession, graphShPt,
411  vSession->GetVariable(0));
412 
413  // Get the total number of elements
414  nElements = Domain->GetExpSize();
415  std::cout << "Number of elements = "
416  << nElements << std::endl;
417 
418  // Get the total number of quadrature points (depends on n. modes)
419  nQuadraturePts = Domain->GetTotPoints();
420  std::cout << "Number of quadrature points = "
421  << nQuadraturePts << std::endl;
422 
423  // Coordinates of the quadrature points
424  x_QuadraturePts = Array<OneD, NekDouble>(nQuadraturePts);
425  y_QuadraturePts = Array<OneD, NekDouble>(nQuadraturePts);
426  z_QuadraturePts = Array<OneD, NekDouble>(nQuadraturePts);
427  Domain->GetCoords(x_QuadraturePts, y_QuadraturePts, z_QuadraturePts);
428  }
429  else if (expdim == 3)
430  {
432  ::AllocateSharedPtr(vSession);
433 
436  ::AllocateSharedPtr(vSession, graphShPt, vSession->GetVariable(0));
437 
438  // Get the total number of elements
439  nElements = Domain->GetExpSize();
440  std::cout << "Number of elements = "
441  << nElements << std::endl;
442 
443  // Get the total number of quadrature points (depends on n. modes)
444  nQuadraturePts = Domain->GetTotPoints();
445  std::cout << "Number of quadrature points = "
446  << nQuadraturePts << std::endl;
447 
448  // Coordinates of the quadrature points
449  x_QuadraturePts = Array<OneD, NekDouble>(nQuadraturePts);
450  y_QuadraturePts = Array<OneD, NekDouble>(nQuadraturePts);
451  z_QuadraturePts = Array<OneD, NekDouble>(nQuadraturePts);
452  Domain->GetCoords(x_QuadraturePts, y_QuadraturePts, z_QuadraturePts);
453  }
454  else
455  {
456  ASSERTL0(false, "Routine available for 2D and 3D problem only.")
457  }
458 
459  // Loading parameters from session file
460  vSession->LoadParameter("Re", m_Re, 1.0);
461  vSession->LoadParameter("Mach", m_Mach, 1.0);
462  vSession->LoadParameter("TInf", m_Tinf, 1.0);
463  vSession->LoadParameter("Twall", m_Twall, 1.0);
464  vSession->LoadParameter("Gamma", m_Gamma, 1.0);
465  vSession->LoadParameter("Pr", m_Pr, 1.0);
466  vSession->LoadParameter("L", m_long, 1.0);
467  vSession->LoadParameter("rhoInf", m_rhoInf, 1.0);
468  vSession->LoadParameter("uInf", m_uInf, 1.0);
469  vSession->LoadParameter("GasConstant", m_R, 1.0);
470  vSession->LoadParameter("vInf", m_vInf, 1.0);
471  vSession->LoadParameter("mu", m_mu, 1.0);
472 
473  // Rescaling factors
474  m_Suth = 110.4 / m_Tinf;
475  m_Tw = m_Twall / m_Tinf;
476  m_Re = m_Re / m_long;
477 
478  cout << "Number of points" << " " << m_xpoints << endl;
479 
480  // Defining the solution arrays
481  Array<OneD, NekDouble> u_QuadraturePts (nQuadraturePts, 0.0);
482  Array<OneD, NekDouble> v_QuadraturePts (nQuadraturePts, 0.0);
483  Array<OneD, NekDouble> rho_QuadraturePts(nQuadraturePts, 0.0);
484  Array<OneD, NekDouble> T_QuadraturePts (nQuadraturePts, 0.0);
485 
486  // Calculation of the similarity variables
487  if (m_Tw > 0)
488  {
489  vstart[3] = m_Tw;
490  }
491  if (m_Tw < 0.0)
492  {
493  v[1] = 1.0 + 0.5 * 0.84 * (m_Gamma - 1) * (m_Mach * m_Mach);
494  v[0] = 0.47 * pow(v[1], 0.21);
495  }
496  else
497  {
498  v[1] = 0.062 * pow(m_Mach, 2) - 0.1 * (m_Tw - 1.0) *
499  (10 + m_Mach) / (0.2 + m_Mach);
500  v[0] = 0.45 - 0.01 * m_Mach + (m_Tw - 1.0) * 0.06;
501  m_Twall = m_Tw;
502  }
503 
504  dv[0] = v[0] * 0.01;
505 
506  if (m_Tw < 0.0)
507  {
508  dv[1] = v[1] * 0.01;
509  }
510  else
511  {
512  dv[1] = 0.1;
513  }
514 
515  vstart[2] = v[0];
516 
517  if (m_Tw < 0)
518  {
519  vstart[3] = v[1];
520  m_Twall = vstart[3];
521  }
522  else
523  {
524  vstart[4] = v[1];
525  }
526 
527  RKDUMB(vstart, 5, 0.0, etamax, m_xpoints, xx, ff);
528 
529  for (int k = 0; k < maxit; k++)
530  {
531  vstart[2] = v[0];
532 
533  if (m_Tw < 0)
534  {
535  vstart[3] = v[1];
536  m_Twall = vstart[3];
537  }
538  else
539  {
540  vstart[4] = v[1];
541  }
542 
543  RKDUMB(vstart, 5, 0.0, etamax, m_xpoints, xx, ff);
544 
545  NekDouble err = fabs(ff[1][m_xpoints] - 1) +
546  fabs(ff[3][m_xpoints] - 1);
547 
548  cout << "err" << scientific << setprecision(9) << " " << err << endl;
549 
550  if (expdim == 2)
551  {
552  if (err < errtol)
553  {
554  cout << "Calculating" << endl;
555  OUTPUT(m_xpoints, xx, ff, nQuadraturePts, x_QuadraturePts,
556  y_QuadraturePts, u_QuadraturePts, v_QuadraturePts,
557  rho_QuadraturePts, T_QuadraturePts);
558  break;
559  }
560  else
561  {
562  f[0] = ff[1][m_xpoints] - 1;
563  f[1] = ff[3][m_xpoints] - 1;
564  vstart[2] = v[0] + dv[0];
565 
566  if (m_Tw < 0)
567  {
568  vstart[3] = v[1];
569  m_Twall = vstart[3];
570  }
571  else
572  {
573  vstart[4] = v[1];
574  }
575 
576  RKDUMB(vstart, 5, 0.0, etamax, m_xpoints, xx, ff);
577 
578  f1[0] = ff[1][m_xpoints] - 1;
579  f1[1] = ff[3][m_xpoints] - 1;
580 
581  vstart[2] = v[0];
582 
583  if (m_Tw < 0)
584  {
585  vstart[3] = v[1] + dv[1];
586  m_Twall = vstart[3];
587  }
588  else
589  {
590  vstart[4] = v[1] + dv[1];
591  }
592 
593  RKDUMB(vstart, 5, 0.0, etamax, m_xpoints, xx, ff);
594 
595  f2[0] = ff[1][m_xpoints] - 1;
596  f2[1] = ff[3][m_xpoints] - 1;
597 
598  al11 = (f1[0] - f[0]) / dv[0];
599  al21 = (f1[1] - f[1]) / dv[0];
600  al12 = (f2[0] - f[0]) / dv[1];
601  al22 = (f2[1] - f[1]) / dv[1];
602  det = al11 * al22 - al21 * al12;
603 
604  dv[0] = ( - al22 * f[0] + al12 * f[1]) / det;
605  dv[1] = (al21 * f[0] - al11 * f[1]) / det;
606  v[0] = v[0] + dv[0];
607  v[1] = v[1] + dv[1];
608  }
609  }
610  else if (expdim == 3)
611  {
612  {
613  if (err < errtol)
614  {
615  cout << "Calculating" << endl;
616  OUTPUT(m_xpoints, xx, ff, nQuadraturePts, x_QuadraturePts,
617  z_QuadraturePts, u_QuadraturePts, v_QuadraturePts,
618  rho_QuadraturePts, T_QuadraturePts);
619  break;
620  }
621  else
622  {
623  f[0] = ff[1][m_xpoints] - 1;
624  f[1] = ff[3][m_xpoints] - 1;
625  vstart[2] = v[0] + dv[0];
626 
627  if (m_Tw < 0)
628  {
629  vstart[3] = v[1];
630  m_Twall = vstart[3];
631  }
632  else
633  {
634  vstart[4] = v[1];
635  }
636 
637  RKDUMB(vstart, 5, 0.0, etamax, m_xpoints, xx, ff);
638 
639  f1[0] = ff[1][m_xpoints] - 1;
640  f1[1] = ff[3][m_xpoints] - 1;
641 
642  vstart[2] = v[0];
643 
644  if (m_Tw < 0)
645  {
646  vstart[3] = v[1] + dv[1];
647  m_Twall = vstart[3];
648  }
649  else
650  {
651  vstart[4] = v[1] + dv[1];
652  }
653 
654  RKDUMB(vstart, 5, 0.0, etamax, m_xpoints, xx, ff);
655 
656  f2[0] = ff[1][m_xpoints] - 1;
657  f2[1] = ff[3][m_xpoints] - 1;
658 
659  al11 = (f1[0] - f[0]) / dv[0];
660  al21 = (f1[1] - f[1]) / dv[0];
661  al12 = (f2[0] - f[0]) / dv[1];
662  al22 = (f2[1] - f[1]) / dv[1];
663  det = al11 * al22 - al21 * al12;
664 
665  dv[0] = ( - al22 * f[0] + al12 * f[1]) / det;
666  dv[1] = (al21 * f[0] - al11 * f[1]) / det;
667  v[0] = v[0] + dv[0];
668  v[1] = v[1] + dv[1];
669  }
670  }
671  }
672  }
673 
674  // Verification of the compressible similarity solution
675  ofstream verif;
676  verif.open("similarity_solution.dat");
677  for (int i=0; i< nQuadraturePts; i++)
678  {
679  verif << scientific << setprecision(9) << x_QuadraturePts[i]
680  << " \t " << y_QuadraturePts[i] << " \t " ;
681  verif << scientific << setprecision(9) << u_QuadraturePts[i]
682  << " \t " << v_QuadraturePts[i] << " \t " ;
683  verif << scientific << setprecision(9) << rho_QuadraturePts[i]
684  << " \t " << T_QuadraturePts[i] << endl;
685  }
686  verif.close();
687 
688  // Calculation of the physical variables
689  for (int i = 0; i < nQuadraturePts; i++)
690  {
691  rho_QuadraturePts[i] = rho_QuadraturePts[i] * m_rhoInf;
692  u_QuadraturePts[i] = u_QuadraturePts[i] * m_uInf;
693  v_QuadraturePts[i] = v_QuadraturePts[i] * m_uInf;
694  T_QuadraturePts[i] = T_QuadraturePts[i] * m_Tinf;
695 
696  T_QuadraturePts[i] = T_QuadraturePts[i] * rho_QuadraturePts[i] * m_R;
697  T_QuadraturePts[i] = T_QuadraturePts[i] / (m_Gamma-1);
698  T_QuadraturePts[i] = T_QuadraturePts[i] + 0.5 * rho_QuadraturePts[i] * (
699  pow(u_QuadraturePts[i], 2.0) + pow(v_QuadraturePts[i], 2.0));
700 
701  u_QuadraturePts[i] = u_QuadraturePts[i] * rho_QuadraturePts[i];
702  v_QuadraturePts[i] = v_QuadraturePts[i] * rho_QuadraturePts[i];
703  }
704  string file_name;
705  if (expdim == 2)
706  {
708  ::AllocateSharedPtr(vSession);
709 
712  ::AllocateSharedPtr(vSession, graphShPt, vSession->GetVariable(0));
713 
717  ::AllocateSharedPtr(vSession, graphShPt);
718 
721  ::AllocateSharedPtr(vSession, graphShPt);
722 
725  ::AllocateSharedPtr(vSession, graphShPt);
726 
729  ::AllocateSharedPtr(vSession, graphShPt);
730 
731  // Filling the 2D expansion using a recursive algorithm based on the
732  // mesh ordering
734  Basis = Domain->GetExp(0)->GetBasis(0);
735  numModes = Basis->GetNumModes();
736 
737  std::cout << "Number of modes = " << numModes << std::endl;
738 
739  // Copying the ukGlobal vector in m_phys (with the same pattern of
740  // m_phys)
741  Vmath::Vcopy(nQuadraturePts, u_QuadraturePts, 1,
742  Exp2D_uk->UpdatePhys(), 1);
743  Vmath::Vcopy(nQuadraturePts, v_QuadraturePts, 1,
744  Exp2D_vk->UpdatePhys(), 1);
745  Vmath::Vcopy(nQuadraturePts, rho_QuadraturePts, 1,
746  Exp2D_rhok->UpdatePhys(), 1);
747  Vmath::Vcopy(nQuadraturePts, T_QuadraturePts, 1,
748  Exp2D_Tk->UpdatePhys(), 1);
749 
750  // Initialisation of the ExpList Exp
751  Exp[0] = Exp2D_rhok;
752  Exp[1] = Exp2D_uk;
753  Exp[2] = Exp2D_vk;
754  Exp[3] = Exp2D_Tk;
755 
756  // Expansion coefficient extraction (necessary to write the .fld file)
757  Exp[0]->FwdTrans(Exp2D_rhok->GetPhys(), Exp[0]->UpdateCoeffs());
758  Exp[1]->FwdTrans(Exp2D_uk->GetPhys(), Exp[1]->UpdateCoeffs());
759  Exp[2]->FwdTrans(Exp2D_vk->GetPhys(), Exp[2]->UpdateCoeffs());
760  Exp[3]->FwdTrans(Exp2D_Tk->GetPhys(), Exp[3]->UpdateCoeffs());
761 
762  // Definition of the name of the .fld file
763  cout << argv[1] << endl;
764  string tmp = argv[1];
765  int len = tmp.size();
766  for (int i = 0; i < len-4; ++i)
767  {
768  file_name += argv[1][i];
769  }
770  file_name = file_name+".rst";
771 
772  // Definition of the Field
773  std::vector<LibUtilities::FieldDefinitionsSharedPtr>
774  FieldDef = Exp[0]->GetFieldDefinitions();
775  std::vector<std::vector<NekDouble> > FieldData(FieldDef.size());
776 
777  for (j = 0; j < 4; j++)
778  {
779  for (i = 0; i < FieldDef.size(); i++)
780  {
781  if (j == 0)
782  {
783  FieldDef[i]->m_fields.push_back("rho");
784  }
785  else if (j == 1)
786  {
787  FieldDef[i]->m_fields.push_back("rhou");
788  }
789  else if (j == 2 )
790  {
791  FieldDef[i]->m_fields.push_back("rhov");
792  }
793  else if (j == 3 )
794  {
795  FieldDef[i]->m_fields.push_back("E");
796  }
797  Exp[j]->AppendFieldData(FieldDef[i], FieldData[i]);
798  }
799  }
800 
801  LibUtilities::Write(file_name, FieldDef, FieldData);
802  }
803  else if (expdim == 3)
804  {
806  ::AllocateSharedPtr(vSession);
807 
810  ::AllocateSharedPtr(vSession, graphShPt, vSession->GetVariable(0));
811 
812  Array<OneD,NekDouble> w_QuadraturePts;
813  w_QuadraturePts = Array<OneD,NekDouble>(nQuadraturePts, 0.0);
815 
818  ::AllocateSharedPtr(vSession, graphShPt);
819 
822  ::AllocateSharedPtr(vSession, graphShPt);
823 
826  ::AllocateSharedPtr(vSession, graphShPt);
827 
830  ::AllocateSharedPtr(vSession, graphShPt);
831 
834  ::AllocateSharedPtr(vSession, graphShPt);
835 
836  // Filling the 3D expansion using a recursive algorithm based
837  // on the mesh ordering
839  Basis = Domain->GetExp(0)->GetBasis(0);
840  numModes = Basis->GetNumModes();
841 
842  std::cout<< "Number of modes = " << numModes << std::endl;
843 
844  // Copying the ukGlobal vector in m_phys (with the same pattern
845  // of m_phys)
846  Vmath::Vcopy(nQuadraturePts, rho_QuadraturePts, 1,
847  Exp3D_rhok->UpdatePhys(), 1);
848  Vmath::Vcopy(nQuadraturePts, u_QuadraturePts, 1,
849  Exp3D_uk->UpdatePhys(), 1);
850  Vmath::Vcopy(nQuadraturePts, w_QuadraturePts, 1,
851  Exp3D_vk->UpdatePhys(), 1);
852  Vmath::Vcopy(nQuadraturePts, v_QuadraturePts, 1,
853  Exp3D_wk->UpdatePhys(), 1);
854  Vmath::Vcopy(nQuadraturePts, T_QuadraturePts, 1,
855  Exp3D_Tk->UpdatePhys(), 1);
856 
857  // Initialisation of the ExpList Exp
858  Exp[0] = Exp3D_rhok;
859  Exp[1] = Exp3D_uk;
860  Exp[2] = Exp3D_vk;
861  Exp[3] = Exp3D_wk;
862  Exp[4] = Exp3D_Tk;
863 
864  // Expansion coefficient extraction (necessary to write the .fld file)
865  Exp[0]->FwdTrans(Exp3D_rhok->GetPhys(), Exp[0]->UpdateCoeffs());
866  Exp[1]->FwdTrans(Exp3D_uk->GetPhys(), Exp[1]->UpdateCoeffs());
867  Exp[2]->FwdTrans(Exp3D_vk->GetPhys(), Exp[2]->UpdateCoeffs());
868  Exp[3]->FwdTrans(Exp3D_wk->GetPhys(), Exp[3]->UpdateCoeffs());
869  Exp[4]->FwdTrans(Exp3D_Tk->GetPhys(), Exp[4]->UpdateCoeffs());
870 
871  // Definition of the name of the .fld file
872  cout << argv[1] << endl;
873  string tmp = argv[1];
874  int len = tmp.size();
875  for (int i = 0; i < len-4; ++i)
876  {
877  file_name += argv[1][i];
878  }
879  file_name = file_name+".rst";
880 
881  // Definition of the Field
882  std::vector<LibUtilities::FieldDefinitionsSharedPtr>
883  FieldDef = Exp[0]->GetFieldDefinitions();
884  std::vector<std::vector<NekDouble> > FieldData(FieldDef.size());
885 
886  for (j = 0; j < 5; j++)
887  {
888  for (i = 0; i < FieldDef.size(); i++)
889  {
890  if (j == 0)
891  {
892  FieldDef[i]->m_fields.push_back("rho");
893  }
894  else if (j == 1)
895  {
896  FieldDef[i]->m_fields.push_back("rhou");
897  }
898  else if (j == 2 )
899  {
900  FieldDef[i]->m_fields.push_back("rhov");
901  }
902  else if (j == 3 )
903  {
904  FieldDef[i]->m_fields.push_back("rhow");
905  }
906  else if (j == 4 )
907  {
908  FieldDef[i]->m_fields.push_back("E");
909  }
910  Exp[j]->AppendFieldData(FieldDef[i], FieldData[i]);
911  }
912  }
913 
914  LibUtilities::Write(file_name, FieldDef, FieldData);
915  }
916 
917  std::cout <<"----------------------------------------------------\n";
918  std::cout <<"\n=================================================\n";
919  std::cout <<"Similarity solution \n";
920  std::cout <<"===================================================\n";
921  std::cout <<"***************************************************\n";
922  std::cout <<"DATA FROM THE SESSION FILE:\n";
923  std::cout << "Reynolds number = " << m_Re
924  << "\t[-]" << std::endl;
925  std::cout << "Mach number = " << m_Mach
926  << "\t[-]" << std::endl;
927  std::cout << "Characteristic length = " << m_long
928  << "\t[m]" << std::endl;
929  std::cout << "U_infinity = " << m_uInf
930  << "\t[m/s]" << std::endl;
931  std::cout <<"***************************************************\n";
932  std::cout <<"---------------------------------------------------\n";
933  std::cout <<"MESH and EXPANSION DATA:\n";
934  std::cout << "Done." << std::endl;
935 
936  return 0;
937 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
NekDouble m_Mach
void OUTPUT(int m_xpoints, Array< OneD, NekDouble > xx, Array< OneD, Array< OneD, NekDouble > > ff, int nQuadraturePts, Array< OneD, NekDouble > x_QuadraturePts, Array< OneD, NekDouble > y_QuadraturePts, Array< OneD, NekDouble > u_QuadraturePts, Array< OneD, NekDouble > v_QuadraturePts, Array< OneD, NekDouble > rho_QuadraturePts, Array< OneD, NekDouble > T_QuadraturePts)
General purpose memory allocation routines with the ability to allocate from thread specific memory p...
NekDouble m_Tw
NekDouble m_Tinf
NekDouble m_vInf
NekDouble m_Twall
NekDouble m_uInf
boost::shared_ptr< ContField2D > ContField2DSharedPtr
Definition: ContField2D.h:291
NekDouble m_mu
boost::shared_ptr< SessionReader > SessionReaderSharedPtr
Definition: MeshPartition.h:50
NekDouble m_rhoInf
const NekDouble etamax
NekDouble m_Pr
NekDouble m_Re
NekDouble m_long
void RKDUMB(Array< OneD, NekDouble > vstart, int nvar, NekDouble x1, NekDouble x2, int m_xpoints, Array< OneD, NekDouble > xx, Array< OneD, Array< OneD, NekDouble > > y)
double NekDouble
NekDouble m_Suth
boost::shared_ptr< ExpList2D > ExpList2DSharedPtr
Shared pointer to an ExpList2D object.
Definition: ExpList2D.h:49
const int m_xpoints
boost::shared_ptr< ExpList3D > ExpList3DSharedPtr
Shared pointer to an ExpList3D object.
Definition: ExpList3D.h:110
NekDouble m_R
void Write(const std::string &outFile, std::vector< FieldDefinitionsSharedPtr > &fielddefs, std::vector< std::vector< NekDouble > > &fielddata, const FieldMetaDataMap &fieldinfomap)
Write a field file in serial only.
Definition: FieldIO.cpp:72
boost::shared_ptr< Basis > BasisSharedPtr
const NekDouble errtol
boost::shared_ptr< ContField3D > ContField3DSharedPtr
Definition: ContField3D.h:190
boost::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:432
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1038
NekDouble m_Gamma
void OUTPUT ( int  m_xpoints,
Array< OneD, NekDouble xx,
Array< OneD, Array< OneD, NekDouble > >  ff,
int  nQuadraturePts,
Array< OneD, NekDouble x_QuadraturePts,
Array< OneD, NekDouble y_QuadraturePts,
Array< OneD, NekDouble u_QuadraturePts,
Array< OneD, NekDouble v_QuadraturePts,
Array< OneD, NekDouble rho_QuadraturePts,
Array< OneD, NekDouble T_QuadraturePts 
)

Create the output file

Definition at line 234 of file CompressibleBL.cpp.

References COMPBL(), etamax, m_mu, m_Pr, m_Re, m_rhoInf, m_Suth, m_uInf, and m_xpoints.

Referenced by main().

244 {
253  Array <OneD, NekDouble > velocity(m_xpoints, 0.0);
255 
256 
257  NekDouble dd, dm, scale, flg, dlta;
258  NekDouble xcher, ycher;
259  int index;
260 
261  z[0] = 0.0;
262  NekDouble sumd = 0.0;
263 
264  for (int i = 1; i < m_xpoints ; i++)
265  {
266  z[i] = z[i-1] + 0.5 * (xx[i] - xx[i-1]) * (ff[3][i] + ff[3][i-1]);
267  dm = ff[3][i-1] - ff[1][i-1];
268  dd = ff[3][i] - ff[1][i];
269  sumd = sumd + 0.5 * (xx[i] - xx[i-1]) * (dd + dm);
270 
271  if ((ff[1][i] > 0.999) && (flg < 1.0))
272  {
273  dlta = z[i];
274  flg = 2.0;
275  }
276  }
277 
278  scale = sumd;
279 
280  ofstream file3;
281  file3.open("physical_data.dat");
282 
283  NekDouble xin, rex, delsx, delta;
284 
285  for (int i = 0; i < m_xpoints; i++)
286  {
287  for (int k = 0; k < 5; k++)
288  {
289  v[k] = ff[k][i];
290  }
291  COMPBL(v, dv);
292  u[i] = ff[1][i];
293  t[i] = ff[3][i];
294  rho[i] = (1.0 / ff[3][i]);
295  vv[i] = -ff[0][i]/sqrt(m_uInf);
296  mu[i] = pow(t[i], 1.5) * (1 + m_Suth) / (t[i] + m_Suth) / (m_Re);
297  velocity[i] = ff[0][i] ;
298  }
299 
300  NekDouble scale2, coeff;
301 
302  for (int i = 0; i < nQuadraturePts; i++)
303  {
304  if (i%100000 == 0)
305  {
306  cout << "i" << " " << i << "/" << nQuadraturePts << endl;
307  }
308 
309  xcher = x_QuadraturePts[i];
310  ycher = y_QuadraturePts[i];
311 
312  scale = sumd;
313  xin = xcher;
314  rex = 0.5 * pow(((m_Re) / scale), 2) + (m_Re) * xin;
315  delsx = sqrt(2.0 / rex) * scale * (xin)* m_Pr;
316  scale = scale / delsx;
317  delta = 4.91 * sqrt((xin * m_mu) / (m_rhoInf * m_uInf));
318  scale2 = ycher * (scale * delta) / sqrt(etamax) ;
319  coeff = 0.5 * sqrt( 2 / (xcher*m_Re)) ;
320 
321  if (scale2 > z[m_xpoints-3])
322  {
323  u_QuadraturePts[i] = 1;
324  rho_QuadraturePts[i] = 1;
325  T_QuadraturePts[i] = 1.0 / rho_QuadraturePts[i];
326  v_QuadraturePts[i] = coeff * (z[m_xpoints-3] -
327  velocity[m_xpoints-3]);
328 
329  file3 << xcher << " "
330  << ycher << " "
331  << velocity[m_xpoints-3] << " "
332  << z[m_xpoints-3] << " "
333  << u[m_xpoints-3]
334  << endl;
335  }
336  else
337  {
338  for (int j = 0 ; j< m_xpoints-1; j++)
339  {
340  if ((z[j] <= scale2) && (z[j+1] > scale2))
341  {
342  index = j;
343  break;
344  }
345  }
346 
347  u_QuadraturePts[i] = u[index];
348  rho_QuadraturePts[i] = rho[index];
349  T_QuadraturePts[i] = 1.0/rho_QuadraturePts[i];
350  v_QuadraturePts[i] = coeff * (u[index]*scale2 - velocity[index]);
351  }
352  }
353 }
NekDouble m_uInf
NekDouble m_mu
NekDouble m_rhoInf
const NekDouble etamax
NekDouble m_Pr
NekDouble m_Re
double NekDouble
NekDouble m_Suth
const int m_xpoints
void COMPBL(Array< OneD, NekDouble > v, Array< OneD, NekDouble > dv)
void RK4 ( Array< OneD, NekDouble y,
Array< OneD, NekDouble dydx,
int  n,
NekDouble  x,
NekDouble  h,
Array< OneD, NekDouble yout 
)

Perform the RK4 integration

Definition at line 141 of file CompressibleBL.cpp.

References COMPBL().

Referenced by RKDUMB().

147 {
148  int nmax = 5;
149 
150  Array<OneD, NekDouble> yt (nmax, 0.0);
151  Array<OneD, NekDouble> dyt(nmax, 0.0);
152  Array<OneD, NekDouble> dym(nmax, 0.0);
153  NekDouble hh = h * 0.5;
154  NekDouble h6 = h / 6;
155 
156  for (int i = 0; i < n ; i++)
157  {
158  yt[i] = y[i] + hh * dydx[i];
159  }
160 
161  COMPBL(yt, dyt);
162 
163  for (int i = 0; i < n; i++)
164  {
165  yt[i] = y[i] + hh * dyt[i];
166  }
167 
168  COMPBL(yt, dym);
169 
170  for (int i = 0; i < n; i++)
171  {
172  yt[i] = y[i] + h * dym[i];
173  dym[i] = dyt[i] + dym[i];
174  }
175 
176  COMPBL(yt, dyt);
177 
178  for (int i = 0; i < n; i++)
179  {
180  yout[i] = y[i] + h6 * (dydx[i] + dyt[i] + 2 * dym[i]);
181  }
182 }
double NekDouble
void COMPBL(Array< OneD, NekDouble > v, Array< OneD, NekDouble > dv)
void RKDUMB ( Array< OneD, NekDouble vstart,
int  nvar,
NekDouble  x1,
NekDouble  x2,
int  m_xpoints,
Array< OneD, NekDouble xx,
Array< OneD, Array< OneD, NekDouble > >  y 
)

Calculate initial guess for RK4

Definition at line 188 of file CompressibleBL.cpp.

References COMPBL(), m_xpoints, and RK4().

Referenced by main().

195 {
196  int nmax = 5;
197  NekDouble x, h;
198  Array<OneD, NekDouble> v (nmax, 0.0);
199  Array<OneD, NekDouble> dv(nmax, 0.0);
200 
201  for (int i = 0; i < nvar; i++)
202  {
203  v[i] = vstart[i];
204  y[i][0] = v[i];
205  }
206 
207  xx[0] = x1;
208  x = x1;
209  h = (x2-x1) / m_xpoints;
210 
211  for (int k = 0; k < m_xpoints; k++)
212  {
213  COMPBL(v, dv);
214  RK4 (v, dv, nvar, x, h, v);
215 
216  if (x + h == x)
217  {
218  cout << "bug" << endl;
219  }
220 
221  x = x + h;
222  xx[k+1] = x;
223 
224  for (int i = 0; i < nvar; i++)
225  {
226  y[i][k+1] = v[i];
227  }
228  }
229 }
void RK4(Array< OneD, NekDouble > y, Array< OneD, NekDouble > dydx, int n, NekDouble x, NekDouble h, Array< OneD, NekDouble > yout)
double NekDouble
const int m_xpoints
void COMPBL(Array< OneD, NekDouble > v, Array< OneD, NekDouble > dv)

Variable Documentation

const NekDouble errtol = 1e-5

Definition at line 99 of file CompressibleBL.cpp.

Referenced by main().

const NekDouble etamax = 10.0

Definition at line 98 of file CompressibleBL.cpp.

Referenced by main(), and OUTPUT().

Definition at line 79 of file CompressibleBL.cpp.

Referenced by main().

NekDouble m_Gamma

Definition at line 84 of file CompressibleBL.cpp.

Referenced by COMPBL(), and main().

NekDouble m_long

Definition at line 86 of file CompressibleBL.cpp.

Referenced by main().

NekDouble m_Mach

Definition at line 78 of file CompressibleBL.cpp.

Referenced by COMPBL(), and main().

NekDouble m_mu

Definition at line 91 of file CompressibleBL.cpp.

Referenced by main(), and OUTPUT().

NekDouble m_Pr

Definition at line 85 of file CompressibleBL.cpp.

Referenced by COMPBL(), main(), and OUTPUT().

NekDouble m_R

Definition at line 89 of file CompressibleBL.cpp.

Referenced by main().

NekDouble m_Re

Definition at line 77 of file CompressibleBL.cpp.

Referenced by main(), and OUTPUT().

NekDouble m_rhoInf

Definition at line 88 of file CompressibleBL.cpp.

Referenced by main(), and OUTPUT().

NekDouble m_Suth

Definition at line 81 of file CompressibleBL.cpp.

Referenced by COMPBL(), main(), and OUTPUT().

NekDouble m_Tinf

Definition at line 80 of file CompressibleBL.cpp.

Referenced by main().

NekDouble m_To = 273.11

Definition at line 92 of file CompressibleBL.cpp.

NekDouble m_Tw

Definition at line 82 of file CompressibleBL.cpp.

Referenced by main().

NekDouble m_Twall

Definition at line 83 of file CompressibleBL.cpp.

Referenced by COMPBL(), and main().

NekDouble m_uInf

Definition at line 87 of file CompressibleBL.cpp.

Referenced by main(), and OUTPUT().

NekDouble m_vInf

Definition at line 90 of file CompressibleBL.cpp.

Referenced by main().

const int m_xpoints = 1000001

Definition at line 94 of file CompressibleBL.cpp.

Referenced by main(), OUTPUT(), and RKDUMB().

const NekDouble Nvisc = 1

Definition at line 96 of file CompressibleBL.cpp.

Referenced by COMPBL().

const NekDouble Omega = 1

Definition at line 97 of file CompressibleBL.cpp.

Referenced by COMPBL().