mpi/solvers/fsi_channel_with_leaflet.cc File Reference

Classes

class  UndeformedLeaflet
 GeomObject: Undeformed straight, vertical leaflet. More...
 
class  FSIChannelWithLeafletProblem< ELEMENT >
 FSI leaflet in channel. More...
 
class  SimpleFSIPreconditioner< MATRIX >
 

Namespaces

 Global_Physical_Variables
 Global variables.
 

Functions

double Global_Physical_Variables::flux (const double &t)
 
int main (int argc, char *argv[])
 

Function Documentation

◆ main()

int main ( int argc  ,
char argv[] 
)

//////////////////////////////////////////////////////////////////////// //////////////////////////////////////////////////////////////////////// //////////////////////////////////////////////////////////////////////// Driver code – pass a command line argument if you want to run the code in validation mode where it only performs a few steps

896 {
897 #ifdef OOMPH_HAS_MPI
898  MPI_Helpers::init(argc,argv);
899 #endif
900 
901  //Parameters for the leaflet: x-position of root and height
902  double x_0 = 1.0;
903  double hleaflet=0.5;
904 
905  // Number of elements in various regions of mesh
906  unsigned nleft=6;
907  unsigned nright=18;
908  unsigned ny1=3;
909  unsigned ny2=3;
910 
911  // Dimensions of fluid mesh: length to the left and right of leaflet
912  // and total height
913  double lleft =1.0;
914  double lright=3.0;
915  double htot=1.0;
916 
917  //Build the problem
920  problem(lleft,lright,hleaflet,
921  htot,nleft,nright,ny1,ny2,x_0);
922 
923  // Create the solver.
924 #ifdef OOMPH_HAS_TRILINOS
926  solver_pt->solver_type() = TrilinosAztecOOSolver::GMRES;
927 #else
929 #endif
930 
931  // Pass the solver to the problem.
932  problem.linear_solver_pt() = solver_pt;
933 
934  // Create the preconditioner
935  SimpleFSIPreconditioner<CRDoubleMatrix>* preconditioner_pt
937 
938  // Pass the meshes to the preconditioner.
939  preconditioner_pt->set_navier_stokes_mesh(problem.fluid_mesh_pt());
940  preconditioner_pt->set_solid_mesh(problem.solid_mesh_pt());
941 
942  // Pass the preconditioner to the solver
943  solver_pt->preconditioner_pt() = preconditioner_pt;
944 
945  // Set up doc info
946  DocInfo doc_info;
947  doc_info.set_directory("RESLT");
948 
949  // Set max. number of adaptations
950  unsigned max_adapt=3;
951 
952  // solve and document
953  problem.steady_newton_solve(max_adapt);
954  problem.doc_solution(doc_info);
955 
956 #ifdef OOMPH_HAS_MPI
957  MPI_Helpers::finalize();
958 #endif
959 }//end of main
FSI leaflet in channel.
Definition: interaction/fsi_channel_with_leaflet/fsi_channel_with_leaflet.cc:165
Definition: algebraic_elements.h:540
Definition: oomph_utilities.h:499
void set_directory(const std::string &directory)
Definition: oomph_utilities.cc:298
The GMRES method.
Definition: iterative_linear_solver.h:1227
Preconditioner *& preconditioner_pt()
Access function to preconditioner.
Definition: iterative_linear_solver.h:95
Definition: fsi_preconditioners.h:485
void set_navier_stokes_mesh(Mesh *mesh_pt, const bool &allow_multiple_element_type_in_navier_stokes_mesh=false)
Broken assignment operator.
Definition: fsi_preconditioners.h:531
Definition: trilinos_solver.h:267
unsigned & solver_type()
Access function to Solver_type.
Definition: trilinos_solver.h:442
Constructor for SteadyAxisymAdvectionDiffusion problem
Definition: steady_axisym_advection_diffusion.cc:213

References oomph::IterativeLinearSolver::preconditioner_pt(), problem, oomph::DocInfo::set_directory(), oomph::SimpleFSIPreconditioner< MATRIX >::set_navier_stokes_mesh(), and oomph::TrilinosAztecOOSolver::solver_type().