flux_control.cc File Reference
#include <iostream>
#include <cstring>
#include "generic.h"
#include "navier_stokes.h"
#include "meshes/collapsible_channel_mesh.h"
#include "womersley.h"
#include "multi_physics.h"

Classes

class  OscillatingWall
 
class  CollapsibleChannelProblem< ELEMENT >
 Problem class. More...
 

Namespaces

 BL_Squash
 
 Global_Physical_Variables
 Global variables.
 

Functions

double BL_Squash::squash_fct (const double &s)
 
void Global_Physical_Variables::prescribed_traction (const double &t, const Vector< double > &x, const Vector< double > &n, Vector< double > &traction)
 Traction applied on the fluid at the left (inflow) boundary. More...
 
int main (int argc, char *argv[])
 

Variables

DataGlobal_Physical_Variables::Pout_data_pt
 Pointer to Data holding downstream pressure load. More...
 

Function Documentation

◆ main()

int main ( int argc  ,
char argv[] 
)

Driver code for an unsteady adaptive collapsible channel problem with prescribed wall motion. Presence of command line arguments indicates validation run with coarse resolution and small number of timesteps.

1528 {
1529  // Store command line arguments
1530  CommandLineArgs::setup(argc,argv);
1531 
1532  //Uncomment this lot if you want to turn this into a parallel code
1533 /*#ifdef OOMPH_HAS_MPI
1534  // Set up MPI_Helpers
1535  MPI_Helpers::init(argc,argv);
1536 
1537  // Swtich off output modifier
1538  oomph_info.output_modifier_pt() = &default_output_modifier;
1539 
1540  // switch off oomph_info output for all processors but rank 0
1541  if (MPI_Helpers::My_rank!=0)
1542  {
1543  oomph_info.stream_pt() = &oomph_nullstream;
1544  OomphLibWarning::set_stream_pt(&oomph_nullstream);
1545  OomphLibError::set_stream_pt(&oomph_nullstream);
1546  }
1547  else
1548  {
1549  oomph_info << "\n\n=====================================================\n";
1550  oomph_info << "Number of processors: "
1551  << MPI_Helpers::Nproc << "\n";
1552  }
1553 #endif
1554 */
1555 
1556  // Set default values
1557  string directory_for_data = "RESLT";
1558  unsigned outflow = 0;
1560  unsigned solver=1;
1561  unsigned p_solver=2;
1562  unsigned f_solver=2;
1563  unsigned f_bamg_smoother=1;
1564  double f_bamg_damping=1.0;
1565  double f_bamg_strength=0.25;
1566  unsigned f_ml_settings=1;
1567  unsigned refinement=0;
1568  unsigned nsteps=200;
1569  double amplitude=1.0e-2;
1570  double period=0.45;
1571  bool validation_run=false;
1572 
1573  // Parse command line
1574  int arg_index = 0;
1575  unsigned print_help = 0;
1576 
1577  while (arg_index < argc)
1578  {
1579  if ( strcmp(argv[arg_index], "-dir") == 0 )
1580  {
1581  arg_index++;
1582  string directory_number = argv[arg_index++];
1583  directory_for_data += directory_number;
1584  }
1585  else if ( strcmp(argv[arg_index], "-outflow") == 0 )
1586  {
1587  arg_index++;
1588  outflow = atoi(argv[arg_index++]);
1589  }
1590  else if ( strcmp(argv[arg_index], "-Re") == 0 )
1591  {
1592  arg_index++;
1593  Global_Physical_Variables::Re = atof(argv[arg_index++]);
1594  }
1595  else if ( strcmp(argv[arg_index], "-solver") == 0 )
1596  {
1597  arg_index++;
1598  solver = atoi(argv[arg_index++]);
1599  }
1600  else if ( strcmp(argv[arg_index], "-p_solver") == 0 )
1601  {
1602  arg_index++;
1603  p_solver = atoi(argv[arg_index++]);
1604  }
1605  else if ( strcmp(argv[arg_index], "-f_solver") == 0 )
1606  {
1607  arg_index++;
1608  f_solver = atoi(argv[arg_index++]);
1609  }
1610  else if ( strcmp(argv[arg_index], "-f_bamg_smoother") == 0 )
1611  {
1612  arg_index++;
1613  f_bamg_smoother = atoi(argv[arg_index++]);
1614  }
1615  else if ( strcmp(argv[arg_index], "-f_bamg_damping") == 0 )
1616  {
1617  arg_index++;
1618  f_bamg_damping = atof(argv[arg_index++]);
1619  }
1620  else if ( strcmp(argv[arg_index], "-f_bamg_strength") == 0 )
1621  {
1622  arg_index++;
1623  f_bamg_strength = atof(argv[arg_index++]);
1624  }
1625  else if ( strcmp(argv[arg_index], "-f_ml_settings") == 0 )
1626  {
1627  arg_index++;
1628  f_ml_settings = atoi(argv[arg_index++]);
1629  }
1630  else if ( strcmp(argv[arg_index], "-refine") == 0 )
1631  {
1632  arg_index++;
1633  refinement = atoi(argv[arg_index++]);
1634  }
1635  else if ( strcmp(argv[arg_index], "-amplitude") == 0 )
1636  {
1637  arg_index++;
1638  amplitude = atof(argv[arg_index++]);
1639  }
1640  else if ( strcmp(argv[arg_index], "-nsteps") == 0 )
1641  {
1642  arg_index++;
1643  nsteps = atoi(argv[arg_index++]);
1644  }
1645  else if ( strcmp(argv[arg_index], "-period") == 0 )
1646  {
1647  arg_index++;
1648  period = atof(argv[arg_index++]);
1649  }
1650  else if ( strcmp(argv[arg_index], "-validation_run") == 0 )
1651  {
1652  arg_index++;
1653  validation_run=true;
1654  }
1655  else if ( strcmp(argv[arg_index], "-help") == 0 )
1656  {
1657  print_help = 1;
1658  break;
1659  }
1660  else
1661  {
1662  arg_index++;
1663  }
1664  }
1665 
1666  if (print_help)
1667  {
1668  //if (MPI_Helpers::My_rank==0)
1669  {
1670  oomph_info << "\n\nOption flags:\n";
1671  oomph_info << "-dir <n> Data saved to /RESLTn\n";
1672  oomph_info << "-outflow <ID> Outflow type:\n"
1673  << " ID=0 flux control\n"
1674  << " ID=1 pressure driven with impedance"
1675  << " with LSC preconditioner\n"
1676  << " ID=2 pressure driven with impedance"
1677  << " with block triangluar and LSC preconditioner\n"
1678  << " ID=3 prescribed outflow\n"
1679 
1680  << " otherwise pressure driven flow\n";
1681  oomph_info << "-Re <value> Sets Re=value\n";
1682  oomph_info << "-nsteps <value> Number of timesteps\n";
1683  oomph_info << "-solver <ID> Linear solver type:\n"
1684  << " ID=0 GMRES\n"
1685  << " otherwise SuperLU is used\n";
1686  oomph_info << "-p_solver <ID> P matrix solve type:\n"
1687  << " ID=0 BoomerAMG\n"
1688  << " ID=1 ML SA-AMG\n"
1689  << " otherwise SuperLU is used\n";
1690  oomph_info << "-f_solver <ID> F matrix solve type\n"
1691  << " ID=0 BoomerAMG\n"
1692  << " ID=1 ML SA-AMG\n"
1693  << " otherwise SuperLU is used\n";
1694  oomph_info << "-f_bamg_smoother <ID> BoomerAMG smoother for F matrix\n"
1695  << " ID=0 damped Jacobi\n"
1696  << " ID=1 Gauss-Seidel\n";
1697  oomph_info << "-f_bamg_damping <value> BoomerAMG Jacobi damping for F matrix\n";
1698  oomph_info << "-f_bamg_strength <value> BoomerAMG strength paramater F matrix\n";
1699  oomph_info << "-f_ml_settings <ID> Settings for ML on F matrix\n"
1700  << " ID=0 NSSA\n"
1701  << " otherwise SA\n";
1702  oomph_info << "-refine <n> Refines mesh n times\n";
1703  oomph_info << "-period <value> Set period to value\n";
1704  oomph_info << "-amplitude <value> Set amplitude to value\n";
1705  oomph_info << "-validation_run Generate validation data for 1D Womersley\n";
1706  }
1707 
1708  return (0);
1709  }
1710 
1711  // Set validation values
1712  if (validation_run)
1713  {
1714  if (outflow==0)
1715  {
1716  directory_for_data = "RESLT_flux_control";
1717  }
1718  else if (outflow==1)
1719  {
1720  directory_for_data = "RESLT_impedance_tube";
1721  }
1722  else if (outflow==2)
1723  {
1724  directory_for_data = "RESLT_impedance_tube_with_flux_control";
1725  }
1728  solver=0;
1729  refinement=0;
1730  amplitude=0.5;
1731  period=10.0;
1732  }
1733 
1734  // Let St=1
1736 
1737  // Output settings
1738  oomph_info << "Storing data in /"
1739  << directory_for_data << "\n";
1740  oomph_info << "Outflow type: "
1741  << outflow << "\n";
1742  oomph_info << "Re="
1744  << "\n";
1745  if (solver==0)
1746  {
1747  oomph_info << "Using GMRES solver\n";
1748  if (p_solver==0)
1749  {
1750 #ifdef OOMPH_HAS_HYPRE
1751  oomph_info << "Using BoomerAMG on P matrix\n";
1752 #else
1753  oomph_info
1754  << "Warning: Hypre not available. Using SuperLU on P matrix\n";
1755 #endif
1756  }
1757  else if (p_solver==1)
1758  {
1759 #ifdef OOMPH_HAS_TRILINOS
1760  oomph_info << "Using ML SA-AMG on P matrix\n";
1761 #else
1762  oomph_info
1763  << "Warning: Trilinos not available. Using SuperLU on P matrix\n";
1764 #endif
1765  }
1766  else
1767  {
1768  oomph_info << "Using SuperLU on P matrix\n";
1769  }
1770  if (f_solver==0)
1771  {
1772 #ifdef OOMPH_HAS_HYPRE
1773  oomph_info << "Using BoomerAMG on F matrix\n"
1774  << "smoother: " << f_bamg_smoother
1775  << "\ndamping: " << f_bamg_damping
1776  << "\nstrength: " << f_bamg_strength
1777  << "\n";
1778 #else
1779  oomph_info
1780  << "Warning: Hypre not available. Using SuperLU on F matrix\n";
1781 #endif
1782  }
1783  else if (f_solver==1)
1784  {
1785 #ifdef OOMPH_HAS_TRILINOS
1786  oomph_info << "Using ML SA-AMG on F matrix ";
1787  if (f_ml_settings==0)
1788  {
1789  oomph_info << "with NSSA settings\n";
1790  }
1791  else
1792  {
1793  oomph_info << "with SA settings\n";
1794  }
1795 #else
1796  oomph_info
1797  << "Warning: Trilinos not available. Using SuperLU on F matrix\n";
1798 #endif
1799  }
1800  else
1801  {
1802  oomph_info << "Using SuperLU on F matrix\n";
1803  }
1804  }
1805  else
1806  {
1807  oomph_info << "Using SuperLU solver\n";
1808  }
1809 
1810  oomph_info << "Refining mesh " << refinement << " times\n";
1811 
1812  // Number of elements in the domain
1813  unsigned nup=20;
1814  unsigned ncollapsible=8;
1815  unsigned ndown=40;
1816  unsigned ny=16;
1817 
1818  // Length of the domain
1819  double lup=5.0;
1820  double lcollapsible=1.0;
1821  double ldown=10.0;
1822  double ly=1.0;
1823  double limpedance=0.0;
1824 
1825  // Modify values if using impedance tube
1826  if ((outflow==1) || (outflow==2))
1827  {
1828  limpedance = 5.0;
1829  ldown-=limpedance;
1830  ndown=20;
1831  }
1832 
1833  // Modify values for validation run
1834  if (validation_run)
1835  {
1836  limpedance = 8.0;
1837  ldown=2.0;
1838  nup=10;
1839  ncollapsible=4;
1840  ny=8;
1841  ndown=4;
1842  }
1843 
1844  oomph_info << "Amplitude=" << amplitude
1845  << "\nPeriod=" << period << "\n";
1846 
1847  // Solver stuff (definitions needed if hypre and trilinos are around)
1848  Preconditioner* p_preconditioner_pt;
1849  p_preconditioner_pt=0;
1850  Preconditioner* f_preconditioner_pt;
1851  f_preconditioner_pt=0;
1852  GMRES<CRDoubleMatrix>* iterative_solver_pt=0;
1853  NavierStokesSchurComplementPreconditioner* ns_preconditioner_pt=0;
1854 
1855  // Don't build the problem with Crouzeix Raviart Elements if using the
1856  // LSC preconditioner!
1858  problem(nup, ncollapsible, ndown, ny,
1859  lup, lcollapsible, ldown, ly, limpedance,
1860  amplitude, period, outflow);
1861 
1862  //problem.switch_off_wall_oscillations();
1863 
1864  // Refine problem
1865  for (unsigned i=0;i<refinement;i++)
1866  {
1867  problem.refine_uniformly();
1868  }
1869 
1870  if (solver==0)
1871  {
1872  // Set up the solver and pass to the problem
1873  iterative_solver_pt = new GMRES<CRDoubleMatrix>;
1874  // iterative_solver_pt->set_preconditioner_RHS();
1875  iterative_solver_pt->max_iter() = 200;
1876  iterative_solver_pt->tolerance() = 1.0e-8;
1877  problem.linear_solver_pt() = iterative_solver_pt;
1878 
1879  // Set up the preconditioner
1880  NavierStokesSchurComplementPreconditioner* ns_preconditioner_pt = 0;
1881  FSIPreconditioner* fsi_preconditioner_pt = 0;
1882 
1883  if (outflow==2)
1884  {
1885  // Construct the FSI preconditioner
1886  fsi_preconditioner_pt = new FSIPreconditioner(&problem);
1887  fsi_preconditioner_pt->enable_doc_time();
1888 
1889  // Get a pointer to the LSC preconditioner
1890  ns_preconditioner_pt =
1892  (fsi_preconditioner_pt->navier_stokes_preconditioner_pt());
1893 
1894  // Suppress warning about elements that don't make contribution
1895  // to diagonal mass matrix
1896  ns_preconditioner_pt->
1897  enable_accept_non_NavierStokesElementWithDiagonalMassMatrices_elements();
1898 
1899  // Set solid mesh with "true" to tolerate multiple element types
1900  // in the mesh.
1901  fsi_preconditioner_pt->set_wall_mesh(
1902  problem.outflow_impedance_master_mesh_pt(),true);
1903 
1904  // Set fluid mesh
1905  fsi_preconditioner_pt->set_navier_stokes_mesh(
1906  problem.create_mesh_for_navier_stokes_preconditioner());
1907 
1908  // Pass the FSI preconditioner to the solver
1909  iterative_solver_pt->preconditioner_pt() = fsi_preconditioner_pt;
1910  }
1911  else
1912  {
1913  // Construct the LSC preconditioner
1914  ns_preconditioner_pt = new
1916 
1917 
1918  // Suppress warning about elements that don't make contribution
1919  // to diagonal mass matrix
1920  ns_preconditioner_pt->
1921  enable_accept_non_NavierStokesElementWithDiagonalMassMatrices_elements();
1922 
1923  // Setup the fluid mesh with "true" to tolerate multiple element types
1924  // in the mesh
1925  ns_preconditioner_pt->set_navier_stokes_mesh(
1926  problem.create_mesh_for_navier_stokes_preconditioner(),true);
1927 
1928  // Pass the LSC preconditioner to the solver
1929  iterative_solver_pt->preconditioner_pt() = ns_preconditioner_pt;
1930  }
1931 
1932  // Settings for the LCS preconditioner, tolerate multiple element types
1933  // in this mesh.
1934  ns_preconditioner_pt->set_navier_stokes_mesh(
1935  problem.create_mesh_for_navier_stokes_preconditioner(),true);
1936  ns_preconditioner_pt->enable_doc_time();
1937 
1938  // Set up the P sub block preconditioner
1939 #ifdef OOMPH_HAS_HYPRE
1940  if (p_solver==0)
1941  {
1942  p_preconditioner_pt = new HyprePreconditioner;
1943  HyprePreconditioner* hypre_preconditioner_pt =
1944  static_cast<HyprePreconditioner*>(p_preconditioner_pt);
1946  set_defaults_for_2D_poisson_problem(hypre_preconditioner_pt);
1947  ns_preconditioner_pt->set_p_preconditioner(p_preconditioner_pt);
1948  }
1949 #endif
1950 
1951 #ifdef OOMPH_HAS_TRILINOS
1952  if (p_solver==1)
1953  {
1954  p_preconditioner_pt = new TrilinosMLPreconditioner;
1955  ns_preconditioner_pt->set_p_preconditioner(p_preconditioner_pt);
1956  }
1957 #endif
1958 
1959  // Set up the F sub block precondioner
1960 #ifdef OOMPH_HAS_HYPRE
1961  if (f_solver==0)
1962  {
1963  f_preconditioner_pt = new HyprePreconditioner;
1964  HyprePreconditioner* hypre_preconditioner_pt =
1965  static_cast<HyprePreconditioner*>(f_preconditioner_pt);
1966  hypre_preconditioner_pt->use_BoomerAMG();
1967  hypre_preconditioner_pt->set_amg_iterations(1);
1968  hypre_preconditioner_pt->amg_using_simple_smoothing();
1969  hypre_preconditioner_pt->amg_simple_smoother()=f_bamg_smoother;
1970  hypre_preconditioner_pt->amg_damping()=f_bamg_damping;
1971  hypre_preconditioner_pt->amg_strength()=f_bamg_strength;
1972  ns_preconditioner_pt->set_f_preconditioner(f_preconditioner_pt);
1973  }
1974 #endif
1975 
1976 #ifdef OOMPH_HAS_TRILINOS
1977  if (f_solver==1)
1978  {
1979  f_preconditioner_pt = new TrilinosMLPreconditioner;
1980  TrilinosMLPreconditioner* trilinos_preconditioner_pt =
1981  static_cast<TrilinosMLPreconditioner*>(f_preconditioner_pt);
1982  if (f_ml_settings==1)
1983  {
1984  trilinos_preconditioner_pt->set_NSSA_default_values();
1985  }
1986  else
1987  {
1988  trilinos_preconditioner_pt->set_SA_default_values();
1989  }
1990  ns_preconditioner_pt->set_f_preconditioner(f_preconditioner_pt);
1991  }
1992 #endif
1993  }
1994 
1995  // run the problem
1996  problem.unsteady_run(directory_for_data, nsteps, validation_run);
1997 
1998  // delete stuff - remembering P and F matrix preconditioners are deleted in the
1999  // LSC preconditioner
2000  delete iterative_solver_pt;
2001  iterative_solver_pt=0;
2002  delete ns_preconditioner_pt;
2003  ns_preconditioner_pt=0;
2004 
2005 
2006  std::cout << "dummy output to suppress compiler warnings: "
2007  << f_bamg_smoother+f_bamg_damping+f_bamg_strength+f_ml_settings
2008  << " " << p_preconditioner_pt
2009  << " " << f_preconditioner_pt
2010  << std::endl;
2011 
2012 
2013 } //end of driver code
BiCGSTAB< SparseMatrix< double > > solver
Definition: BiCGSTAB_simple.cpp:5
int i
Definition: BiCGSTAB_step_by_step.cpp:9
Problem class.
Definition: collapsible_channel.cc:236
Definition: fsi_preconditioners.h:51
void enable_doc_time()
Enable documentation of time.
Definition: fsi_preconditioners.h:205
void set_navier_stokes_mesh(Mesh *mesh_pt, const bool &allow_multiple_element_type_in_navier_stokes_mesh=false)
Definition: fsi_preconditioners.h:164
void set_wall_mesh(Mesh *mesh_pt, const bool &allow_multiple_element_type_in_wall_mesh=false)
Definition: fsi_preconditioners.h:179
NavierStokesSchurComplementPreconditioner * navier_stokes_preconditioner_pt() const
Access function to the Navier Stokes preconditioner (inexact solver)
Definition: fsi_preconditioners.h:198
The GMRES method.
Definition: iterative_linear_solver.h:1227
Definition: hypre_solver.h:826
unsigned & amg_simple_smoother()
Access function to AMG_simple_smoother flag.
Definition: hypre_solver.h:983
double & amg_strength()
Access function to AMG_strength.
Definition: hypre_solver.h:1038
double & amg_damping()
Access function to AMG_damping parameter.
Definition: hypre_solver.h:1032
void set_amg_iterations(const unsigned &amg_iterations)
Function to set the number of times to apply BoomerAMG.
Definition: hypre_solver.h:964
void amg_using_simple_smoothing()
Definition: hypre_solver.h:977
void use_BoomerAMG()
Function to select BoomerAMG as the preconditioner.
Definition: hypre_solver.h:958
Preconditioner *& preconditioner_pt()
Access function to preconditioner.
Definition: iterative_linear_solver.h:95
double & tolerance()
Access to convergence tolerance.
Definition: iterative_linear_solver.h:107
unsigned & max_iter()
Access to max. number of iterations.
Definition: iterative_linear_solver.h:113
Definition: navier_stokes_preconditioners.h:607
void set_p_preconditioner(Preconditioner *new_p_preconditioner_pt)
Function to set a new pressure matrix preconditioner (inexact solver)
Definition: navier_stokes_preconditioners.h:745
void set_f_preconditioner(Preconditioner *new_f_preconditioner_pt)
Function to set a new momentum matrix preconditioner (inexact solver)
Definition: navier_stokes_preconditioners.h:769
void set_navier_stokes_mesh(Mesh *mesh_pt, const bool &allow_multiple_element_type_in_navier_stokes_mesh=false)
Definition: navier_stokes_preconditioners.h:732
void enable_doc_time()
Enable documentation of time.
Definition: navier_stokes_preconditioners.h:805
Definition: preconditioner.h:54
Definition: trilinos_preconditioners.h:152
void set_NSSA_default_values()
Broken assignment operator.
Definition: trilinos_preconditioners.h:177
void set_SA_default_values()
Definition: trilinos_preconditioners.h:200
void setup(Time *time_pt)
Create all GeomObjects needed to define the cylinder and the flag.
Definition: turek_flag_non_fsi.cc:277
double ReSt
Womersley number.
Definition: rayleigh_instability.cc:56
double Re
Reynolds number.
Definition: fibre.cc:55
const double ly
Definition: ConstraintElementsUnitTest.cpp:34
const unsigned ny
Definition: ConstraintElementsUnitTest.cpp:31
void set_defaults_for_2D_poisson_problem(HyprePreconditioner *hypre_preconditioner_pt)
Definition: hypre_solver.cc:45
OomphInfo oomph_info
Definition: oomph_definitions.cc:319
Constructor for SteadyAxisymAdvectionDiffusion problem
Definition: steady_axisym_advection_diffusion.cc:213

References oomph::HyprePreconditioner::amg_damping(), oomph::HyprePreconditioner::amg_simple_smoother(), oomph::HyprePreconditioner::amg_strength(), oomph::HyprePreconditioner::amg_using_simple_smoothing(), oomph::FSIPreconditioner::enable_doc_time(), oomph::NavierStokesSchurComplementPreconditioner::enable_doc_time(), i, Mesh_Parameters::ly, oomph::IterativeLinearSolver::max_iter(), oomph::FSIPreconditioner::navier_stokes_preconditioner_pt(), Mesh_Parameters::ny, oomph::oomph_info, oomph::IterativeLinearSolver::preconditioner_pt(), problem, Global_Physical_Variables::Re, Global_Physical_Variables::ReSt, oomph::HyprePreconditioner::set_amg_iterations(), oomph::Hypre_default_settings::set_defaults_for_2D_poisson_problem(), oomph::NavierStokesSchurComplementPreconditioner::set_f_preconditioner(), oomph::FSIPreconditioner::set_navier_stokes_mesh(), oomph::NavierStokesSchurComplementPreconditioner::set_navier_stokes_mesh(), oomph::TrilinosMLPreconditioner::set_NSSA_default_values(), oomph::NavierStokesSchurComplementPreconditioner::set_p_preconditioner(), oomph::TrilinosMLPreconditioner::set_SA_default_values(), oomph::FSIPreconditioner::set_wall_mesh(), Flag_definition::setup(), solver, oomph::IterativeLinearSolver::tolerance(), and oomph::HyprePreconditioner::use_BoomerAMG().