fsi_channel_with_leaflet.cc
Go to the documentation of this file.
1 //LIC// ====================================================================
2 //LIC// This file forms part of oomph-lib, the object-oriented,
3 //LIC// multi-physics finite-element library, available
4 //LIC// at http://www.oomph-lib.org.
5 //LIC//
6 //LIC// Copyright (C) 2006-2024 Matthias Heil and Andrew Hazel
7 //LIC//
8 //LIC// This library is free software; you can redistribute it and/or
9 //LIC// modify it under the terms of the GNU Lesser General Public
10 //LIC// License as published by the Free Software Foundation; either
11 //LIC// version 2.1 of the License, or (at your option) any later version.
12 //LIC//
13 //LIC// This library is distributed in the hope that it will be useful,
14 //LIC// but WITHOUT ANY WARRANTY; without even the implied warranty of
15 //LIC// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
16 //LIC// Lesser General Public License for more details.
17 //LIC//
18 //LIC// You should have received a copy of the GNU Lesser General Public
19 //LIC// License along with this library; if not, write to the Free Software
20 //LIC// Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
21 //LIC// 02110-1301 USA.
22 //LIC//
23 //LIC// The authors may be contacted at oomph-lib@maths.man.ac.uk.
24 //LIC//
25 //LIC//====================================================================
26 // Generic oomph-lib includes
27 #include "generic.h"
28 #include "beam.h"
29 #include "navier_stokes.h"
30 #include "multi_physics.h"
31 
32 //Include the mesh
33 #include "meshes/channel_with_leaflet_mesh.h"
34 
35 // The wall mesh
36 #include "meshes/one_d_lagrangian_mesh.h"
37 
38 
39 using namespace std;
40 using namespace oomph;
41 
42 
43 //==== start_of_global_parameters================================
44 /// Global parameters
45 //===============================================================
47 {
48  /// Reynolds number
49  double Re=50.0;
50 
51  /// Womersley number: Product of Reynolds and Strouhal numbers
52  double ReSt=50.0;
53 
54  /// Non-dimensional wall thickness.
55  double H=0.05;
56 
57  /// Fluid structure interaction parameter: Ratio of stresses used for
58  /// non-dimensionalisation of fluid to solid stresses.
59  double Q=1.0e-6;
60 
61  /// Period for fluctuations in flux
62  double Period=2.0;
63 
64  /// Min. flux
65  double Min_flux=1.0;
66 
67  /// Max. flux
68  double Max_flux=2.0;
69 
70  /// Flux: Pulsatile flow fluctuating between Min_flux and Max_flux
71  /// with period Period
72  double flux(const double& t)
73  {
74  return Min_flux+
75  (Max_flux-Min_flux)*0.5*(1.0-cos(2.0*MathematicalConstants::Pi*t/Period));
76  }
77 
78 } // end_of_namespace
79 
80 
81 
82 /// ////////////////////////////////////////////////////////////////////
83 /// ///////////////////////////////////////////////////////////////////
84 /// ////////////////////////////////////////////////////////////////////
85 
86 
87 //=====start_of_undeformed_leaflet=================================
88 /// GeomObject: Undeformed straight, vertical leaflet
89 //=================================================================
90 class UndeformedLeaflet : public GeomObject
91 {
92 
93 public:
94 
95  /// Constructor: argument is the x-coordinate of the leaflet
96  UndeformedLeaflet(const double& x0): GeomObject(1,2)
97  {
98  X0=x0;
99  }
100 
101  /// Position vector at Lagrangian coordinate zeta
102  void position(const Vector<double>& zeta, Vector<double>& r) const
103  {
104  // Position Vector
105  r[0] = X0;
106  r[1] = zeta[0];
107  }
108 
109 
110  /// Parametrised position on object: r(zeta). Evaluated at
111  /// previous timestep. t=0: current time; t>0: previous
112  /// timestep. Calls steady version.
113  void position(const unsigned& t, const Vector<double>& zeta,
114  Vector<double>& r) const
115  {
116  // Use the steady version
117  position(zeta,r);
118  } // end of position
119 
120 
121  /// Posn vector and its 1st & 2nd derivatives
122  /// w.r.t. to coordinates:
123  /// \f$ \frac{dR_i}{d \zeta_\alpha}\f$ = drdzeta(alpha,i).
124  /// \f$ \frac{d^2R_i}{d \zeta_\alpha d \zeta_\beta}\f$ =
125  /// ddrdzeta(alpha,beta,i). Evaluated at current time.
126  void d2position(const Vector<double>& zeta,
127  Vector<double>& r,
128  DenseMatrix<double> &drdzeta,
129  RankThreeTensor<double> &ddrdzeta) const
130  {
131  // Position vector
132  r[0] = X0;
133  r[1] = zeta[0];
134 
135  // Tangent vector
136  drdzeta(0,0)=0.0;
137  drdzeta(0,1)=1.0;
138 
139  // Derivative of tangent vector
140  ddrdzeta(0,0,0)=0.0;
141  ddrdzeta(0,0,1)=0.0;
142  } // end of d2position
143 
144  /// Number of geometric Data in GeomObject: None.
145  unsigned ngeom_data() const {return 0;}
146 
147  private :
148 
149  /// x position of the undeformed leaflet's origin.
150  double X0;
151 
152 }; //end_of_undeformed_wall
153 
154 
155 /// ////////////////////////////////////////////////////////////////////
156 /// ///////////////////////////////////////////////////////////////////
157 /// ////////////////////////////////////////////////////////////////////
158 
159 
160 //=====start_of_problem_class========================================
161 /// FSI leaflet in channel
162 //===================================================================
163 template<class ELEMENT>
164 class FSIChannelWithLeafletProblem : public Problem
165 {
166 
167 public:
168 
169  /// Constructor: Pass the lenght of the domain at the left
170  /// of the leaflet lleft,the lenght of the domain at the right of the
171  /// leaflet lright,the height of the leaflet hleaflet, the total height
172  /// of the domain htot, the number of macro-elements at the left of the
173  /// leaflet nleft, the number of macro-elements at the right of the
174  /// leaflet nright, the number of macro-elements under hleaflet ny1,
175  /// the number of macro-elements above hleaflet ny2, the abscissa
176  /// of the origin of the leaflet x_0.
177  FSIChannelWithLeafletProblem(const double& lleft,
178  const double& lright, const double& hleaflet,
179  const double& htot,
180  const unsigned& nleft, const unsigned& nright,
181  const unsigned& ny1, const unsigned& ny2,
182  const double& x_0);
183 
184  /// Destructor empty
186 
187  /// Actions after solve (empty)
189 
190  /// Actions before solve (empty)
192 
193  /// Actions after adaptation
194  void actions_after_adapt();
195 
196  /// Access function to the wall mesh
197  OneDLagrangianMesh<FSIHermiteBeamElement>* wall_mesh_pt()
198  {
199  return Wall_mesh_pt;
200  }
201 
202  /// Access function to fluid mesh
203  RefineableAlgebraicChannelWithLeafletMesh<ELEMENT>* fluid_mesh_pt()
204  {
205  return Fluid_mesh_pt;
206  }
207 
208  /// Doc the solution
209  void doc_solution(DocInfo& doc_info, ofstream& trace);
210 
211 
212 /// Update the inflow velocity
214  {
215  // Actual time
216  double t=time_pt()->time();
217 
218  // Amplitude of flow
219  double ampl=Global_Physical_Variables::flux(t);
220 
221  // Update parabolic flow along boundary 3
222  unsigned ibound=3;
223  unsigned num_nod= Fluid_mesh_pt->nboundary_node(ibound);
224  for (unsigned inod=0;inod<num_nod;inod++)
225  {
226  double ycoord = Fluid_mesh_pt->boundary_node_pt(ibound,inod)->x(1);
227  double uy = ampl*6.0*ycoord/Htot*(1.0-ycoord/Htot);
228  Fluid_mesh_pt->boundary_node_pt(ibound,inod)->set_value(0,uy);
229  Fluid_mesh_pt->boundary_node_pt(ibound,inod)->set_value(1,0.0);
230  }
231  } // end of actions_before_implicit_timestep
232 
233  /// Update before checking Newton convergence: Update the
234  /// nodal positions in the fluid mesh in response to possible
235  /// changes in the wall shape
237  {
238  Fluid_mesh_pt->node_update();
239  }
240 
241 private:
242 
243  /// Pointer to the fluid mesh
244  RefineableAlgebraicChannelWithLeafletMesh<ELEMENT>* Fluid_mesh_pt;
245 
246  /// Pointer to the "wall" mesh
247  OneDLagrangianMesh<FSIHermiteBeamElement>* Wall_mesh_pt;
248 
249  /// Pointer to the GeomObject that represents the wall
250  GeomObject* Leaflet_pt;
251 
252  /// Total height of the domain
253  double Htot;
254 
255 };
256 
257 
258 
259 
260 
261 //=====start_of_constructor==============================================
262 /// Constructor
263 //=======================================================================
264 template <class ELEMENT>
266  const double& lleft,
267  const double& lright,
268  const double& hleaflet,
269  const double& htot,
270  const unsigned& nleft,
271  const unsigned& nright,
272  const unsigned& ny1,
273  const unsigned& ny2,
274  const double& x_0) : Htot(htot)
275 {
276  // Timesteppers:
277  //--------------
278 
279  // Allocate the timestepper
280  BDF<2>* fluid_time_stepper_pt=new BDF<2>;
281  add_time_stepper_pt(fluid_time_stepper_pt);
282 
283  // Allocate the wall timestepper
284  Steady<2>* wall_time_stepper_pt=new Steady<2>;
285  add_time_stepper_pt(wall_time_stepper_pt);
286 
287 
288  // Discretise leaflet
289  //-------------------
290 
291  // Geometric object that represents the undeformed leaflet
292  UndeformedLeaflet* undeformed_wall_pt=new UndeformedLeaflet(x_0);
293 
294  //Create the "wall" mesh with FSI Hermite beam elements
295  unsigned n_wall_el=5;
296  Wall_mesh_pt = new OneDLagrangianMesh<FSIHermiteBeamElement>
297  (n_wall_el,hleaflet,undeformed_wall_pt,wall_time_stepper_pt);
298 
299 
300  // Provide GeomObject representation of leaflet mesh and build fluid mesh
301  //-----------------------------------------------------------------------
302 
303  // Build a geometric object (one Lagrangian, two Eulerian coordinates)
304  // from the wall mesh
305  MeshAsGeomObject* wall_geom_object_pt=
306  new MeshAsGeomObject(Wall_mesh_pt);
307 
308 //Build the mesh
309  Fluid_mesh_pt =new RefineableAlgebraicChannelWithLeafletMesh<ELEMENT>(
310  wall_geom_object_pt,
311  lleft, lright,
312  hleaflet,
313  htot,nleft,
314  nright,ny1,ny2,
315  fluid_time_stepper_pt);
316 
317  // Set error estimator
318  Z2ErrorEstimator* error_estimator_pt=new Z2ErrorEstimator;
319  Fluid_mesh_pt->spatial_error_estimator_pt()=error_estimator_pt;
320 
321 
322 
323  // Build global mesh
324  //------------------
325 
326  // Add the sub meshes to the problem
327  add_sub_mesh(Fluid_mesh_pt);
328  add_sub_mesh(Wall_mesh_pt);
329 
330  // Combine all submeshes into a single Mesh
331  build_global_mesh();
332 
333 
334 
335  // Fluid boundary conditions
336  //--------------------------
337 
338  //Pin the boundary nodes of the fluid mesh
339  unsigned num_bound = Fluid_mesh_pt->nboundary();
340  for(unsigned ibound=0;ibound<num_bound;ibound++)
341  {
342  unsigned num_nod= Fluid_mesh_pt->nboundary_node(ibound);
343  for (unsigned inod=0;inod<num_nod;inod++)
344  {
345  Fluid_mesh_pt->boundary_node_pt(ibound,inod)->pin(1);
346 
347  // Do not pin the x velocity of the outflow
348  if( ibound != 1)
349  {
350  Fluid_mesh_pt->boundary_node_pt(ibound,inod)->pin(0);
351  }
352  }
353  }
354  // end loop over boundaries
355 
356 
357  // Setup parabolic flow along boundary 3 (everything else that's
358  // pinned has homogenous boundary conditions so no action is required
359  // as that's the default assignment). Inflow profile is parabolic
360  // and this is interpolated correctly during mesh refinement so
361  // no re-assignment necessary after adaptation.
362  unsigned ibound=3;
363  unsigned num_nod= Fluid_mesh_pt->nboundary_node(ibound);
364  for (unsigned inod=0;inod<num_nod;inod++)
365  {
366  double ycoord = Fluid_mesh_pt->boundary_node_pt(ibound,inod)->x(1);
367  double uy = 6.0*ycoord/htot*(1.0-ycoord/htot);
368  Fluid_mesh_pt->boundary_node_pt(ibound,inod)->set_value(0,uy);
369  Fluid_mesh_pt->boundary_node_pt(ibound,inod)->set_value(1,0.0);
370  }// end of setup boundary condition
371 
372 
373 
374  // Boundary conditions for wall mesh
375  //----------------------------------
376 
377  // Set the boundary conditions: the lower end of the beam is fixed in space
378  unsigned b=0;
379 
380  // Pin displacements in both x and y directions
381  wall_mesh_pt()->boundary_node_pt(b,0)->pin_position(0);
382  wall_mesh_pt()->boundary_node_pt(b,0)->pin_position(1);
383 
384  // Infinite slope: Pin type 1 (slope) dof for displacement direction 0
385  wall_mesh_pt()->boundary_node_pt(b,0)->pin_position(1,0);
386 
387 
388 
389  // Complete build of fluid elements
390  //---------------------------------
391  unsigned n_element = Fluid_mesh_pt->nelement();
392 
393  // Loop over the elements to set up element-specific
394  // things that cannot be handled by constructor
395  for(unsigned e=0;e<n_element;e++)
396  {
397  // Upcast from GeneralisedElement to the present element
398  ELEMENT* el_pt = dynamic_cast<ELEMENT*>(Fluid_mesh_pt->element_pt(e));
399 
400  //Set the Reynolds number
401  el_pt->re_pt() = &Global_Physical_Variables::Re;
402 
403  //Set the Womersley number
404  el_pt->re_st_pt() = &Global_Physical_Variables::ReSt;
405 
406  }// end loop over elements
407 
408 
409  // Pin redudant pressure dofs
410  RefineableNavierStokesEquations<2>::
411  pin_redundant_nodal_pressures(Fluid_mesh_pt->element_pt());
412 
413 
414  // Complete build of wall elements
415  //--------------------------------
416  n_element = wall_mesh_pt()->nelement();
417  for(unsigned e=0;e<n_element;e++)
418  {
419  // Upcast to the specific element type
420  FSIHermiteBeamElement *elem_pt =
421  dynamic_cast<FSIHermiteBeamElement*>(wall_mesh_pt()->element_pt(e));
422 
423  // Set physical parameters for each element:
424  elem_pt->h_pt() = &Global_Physical_Variables::H;
425 
426  // Function that specifies the load ratios
427  elem_pt->q_pt() = &Global_Physical_Variables::Q;
428 
429  // Set the undeformed shape for each element
430  elem_pt->undeformed_beam_pt() = undeformed_wall_pt;
431 
432  // Leaflet is immersed and loaded by fluid on both sides
433  elem_pt->enable_fluid_loading_on_both_sides();
434 
435  // The normal to the leaflet, as computed by the
436  // FSIHermiteElements points away from the fluid rather than
437  // into the fluid (as assumed by default) when viewed from
438  // the "front" (face 0).
439  elem_pt->set_normal_pointing_out_of_fluid();
440 
441  } // end of loop over elements
442 
443 
444  // Setup FSI
445  //----------
446 
447  // The velocity of the fluid nodes on the wall (fluid mesh boundary 4,5)
448  // is set by the wall motion -- hence the no-slip condition must be
449  // re-applied whenever a node update is performed for these nodes.
450  // Such tasks may be performed automatically by the auxiliary node update
451  // function specified by a function pointer:
452  for(unsigned ibound=4;ibound<6;ibound++ )
453  {
454  unsigned num_nod= Fluid_mesh_pt->nboundary_node(ibound);
455  for (unsigned inod=0;inod<num_nod;inod++)
456  {
457  Fluid_mesh_pt->boundary_node_pt(ibound, inod)->
458  set_auxiliary_node_update_fct_pt(
459  FSI_functions::apply_no_slip_on_moving_wall);
460  }
461  }// aux node update fct has been set
462 
463  // Work out which fluid dofs affect the residuals of the wall elements:
464  // We pass the boundary between the fluid and solid meshes and
465  // pointers to the meshes. The interaction boundary is boundary 4 and 5
466  // of the 2D fluid mesh.
467 
468  // Front of leaflet: Set face=0 (which is also the default so this argument
469  // could be omitted)
470  unsigned face=0;
471  FSI_functions::setup_fluid_load_info_for_solid_elements<ELEMENT,2>
472  (this,4,Fluid_mesh_pt,Wall_mesh_pt,face);
473 
474  // Back of leaflet: face 1, needs to be specified explicitly
475  face=1;
476  FSI_functions::setup_fluid_load_info_for_solid_elements<ELEMENT,2>
477  (this,5,Fluid_mesh_pt,Wall_mesh_pt,face);
478 
479  // Setup equation numbering scheme
480  cout <<"Number of equations: " << assign_eqn_numbers() << std::endl;
481 
482 
483 
484  // Choose iterative solvers with FSI preconditioner (only if not test)
485  //====================================================================
486  if (CommandLineArgs::Argc==1)
487  {
488 
489  // Build iterative linear solver
490  //------------------------------
491  GMRES<CRDoubleMatrix>* iterative_linear_solver_pt =
492  new GMRES<CRDoubleMatrix>;
493 
494  // Set maximum number of iterations
495  iterative_linear_solver_pt->max_iter() = 200;
496 
497  // Pass solver to problem:
498  linear_solver_pt()=iterative_linear_solver_pt;
499 
500 
501  // Build preconditioner
502  //---------------------
503  FSIPreconditioner* prec_pt=new FSIPreconditioner(this);
504 
505  // Set Navier Stokes mesh:
506  prec_pt->set_navier_stokes_mesh(Fluid_mesh_pt);
507 
508  // Set solid mesh:
509  prec_pt->set_wall_mesh(Wall_mesh_pt);
510 
511  // By default, the Schur complement preconditioner uses SuperLU as
512  // an exact preconditioner (i.e. a solver) for the
513  // momentum and Schur complement blocks.
514  // Can overwrite this by passing pointers to
515  // other preconditioners that perform the (approximate)
516  // solves of these blocks
517 
518 
519 #ifdef OOMPH_HAS_HYPRE
520 
521  // Create internal preconditioner used on Schur block
522  //---------------------------------------------------
523  HyprePreconditioner* p_preconditioner_pt = new HyprePreconditioner;
524 
525  // Shut up!
526  p_preconditioner_pt->disable_doc_time();
527 
528  // Set defaults parameters for use as preconditioner on Poisson-type problem
529  Hypre_default_settings::
530  set_defaults_for_2D_poisson_problem(p_preconditioner_pt);
531 
532  // Pass to preconditioner
533  //prec_pt->set_p_preconditioner(p_preconditioner_pt);
534 
535 
536  // Create internal preconditioner used on momentum block
537  //------------------------------------------------------
538  HyprePreconditioner* f_preconditioner_pt = new HyprePreconditioner;
539 
540  // Shut up!
541  f_preconditioner_pt->disable_doc_time();
542 
543  // Set default parameters for use as preconditioner in for momentum
544  // block in Navier-Stokes problem
545  Hypre_default_settings::
546  set_defaults_for_navier_stokes_momentum_block(f_preconditioner_pt);
547 
548  // Use Hypre for momentum block
549  //prec_pt->set_f_preconditioner(f_preconditioner_pt);
550 
551 #endif
552 
553  // Retain fluid onto solid terms in FSI preconditioner
554  prec_pt->use_block_triangular_version_with_fluid_on_solid();
555 
556  // Pass preconditioner to iterative linear solver
557  iterative_linear_solver_pt->preconditioner_pt()= prec_pt;
558 
559  }
560 
561 
562 }//end of constructor
563 
564 
565 
566 
567 
568 //==== start_of_actions_after_adapt=================================
569 /// Actions_after_adapt()
570 //==================================================================
571 template<class ELEMENT>
573 {
574  // Unpin all pressure dofs
575  RefineableNavierStokesEquations<2>::
576  unpin_all_pressure_dofs(Fluid_mesh_pt->element_pt());
577 
578  // Pin redundant pressure dofs
579  RefineableNavierStokesEquations<2>::
580  pin_redundant_nodal_pressures(Fluid_mesh_pt->element_pt());
581 
582 
583  // (Re-)apply the no slip condition on the moving wall
584  //-----------------------------------------------------
585 
586  // The velocity of the fluid nodes on the wall (fluid mesh boundary 4,5)
587  // is set by the wall motion -- hence the no-slip condition must be
588  // re-applied whenever a node update is performed for these nodes.
589  // Such tasks may be performed automatically by the auxiliary node update
590  // function specified by a function pointer:
591  for(unsigned ibound=4;ibound<6;ibound++ )
592  {
593  unsigned num_nod= Fluid_mesh_pt->nboundary_node(ibound);
594  for (unsigned inod=0;inod<num_nod;inod++)
595  {
596  Fluid_mesh_pt->boundary_node_pt(ibound, inod)->
597  set_auxiliary_node_update_fct_pt(
598  FSI_functions::apply_no_slip_on_moving_wall);
599  }
600  } // aux node update fct has been (re-)set
601 
602 
603 
604  // Re-setup FSI
605  //-------------
606 
607  // Work out which fluid dofs affect the residuals of the wall elements:
608  // We pass the boundary between the fluid and solid meshes and
609  // pointers to the meshes. The interaction boundary is boundary 4 and 5
610  // of the 2D fluid mesh.
611 
612  // Front of leaflet: Set face=0 (which is also the default so this argument
613  // could be omitted)
614  unsigned face=0;
615  FSI_functions::setup_fluid_load_info_for_solid_elements<ELEMENT,2>
616  (this,4,Fluid_mesh_pt,Wall_mesh_pt,face);
617 
618  // Back of leaflet: face 1, needs to be specified explicitly
619  face=1;
620  FSI_functions::setup_fluid_load_info_for_solid_elements<ELEMENT,2>
621  (this,5,Fluid_mesh_pt,Wall_mesh_pt,face);
622 
623 } // end_of_actions_after_adapt
624 
625 
626 
627 
628 
629 //==start_of_doc_solution=================================================
630 /// Doc the solution
631 //========================================================================
632 template<class ELEMENT>
634  ofstream& trace)
635 {
636  ofstream some_file;
637  char filename[100];
638 
639  // Number of plot points
640  unsigned npts;
641  npts=5;
642 
643  // Output fluid solution
644  sprintf(filename,"%s/soln%i.dat",doc_info.directory().c_str(),
645  doc_info.number());
646  some_file.open(filename);
647  Fluid_mesh_pt->output(some_file,npts);
648  some_file.close();
649 
650  // Output wall solution
651  sprintf(filename,"%s/wall_soln%i.dat",doc_info.directory().c_str(),
652  doc_info.number());
653  some_file.open(filename);
654  Wall_mesh_pt->output(some_file,npts);
655  some_file.close();
656 
657  // Get node at tip of leaflet
658  unsigned n_el_wall=Wall_mesh_pt->nelement();
659  Node* tip_node_pt=Wall_mesh_pt->finite_element_pt(n_el_wall-1)->node_pt(1);
660 
661  // Get time:
662  double time=time_pt()->time();
663 
664  // Write trace file
665  trace << time << " "
666  << Global_Physical_Variables::flux(time) << " "
667  << tip_node_pt->x(0) << " "
668  << tip_node_pt->x(1) << " "
669  << tip_node_pt->dposition_dt(0) << " "
670  << tip_node_pt->dposition_dt(1) << " "
671  << doc_info.number() << " "
672  << std::endl;
673 
674 
675  // Help me figure out what the "front" and "back" faces of the leaflet are
676  //------------------------------------------------------------------------
677 
678  // Output fluid elements on fluid mesh boundary 4 (associated with
679  // the "front")
680  unsigned bound=4;
681  sprintf(filename,"%s/fluid_boundary_elements_front_%i.dat",
682  doc_info.directory().c_str(),
683  doc_info.number());
684  some_file.open(filename);
685  unsigned nel= Fluid_mesh_pt->nboundary_element(bound);
686  for (unsigned e=0;e<nel;e++)
687  {
688  dynamic_cast<ELEMENT*>(Fluid_mesh_pt->boundary_element_pt(bound,e))
689  ->output(some_file,npts);
690  }
691  some_file.close();
692 
693 
694  // Output fluid elements on fluid mesh boundary 5 (associated with
695  // the "back")
696  bound=5;
697  sprintf(filename,"%s/fluid_boundary_elements_back_%i.dat",
698  doc_info.directory().c_str(),
699  doc_info.number());
700  some_file.open(filename);
701  nel= Fluid_mesh_pt->nboundary_element(bound);
702  for (unsigned e=0;e<nel;e++)
703  {
704  dynamic_cast<ELEMENT*>(Fluid_mesh_pt->boundary_element_pt(bound,e))
705  ->output(some_file,npts);
706  }
707  some_file.close();
708 
709 
710  // Output normal vector on wall elements
711  sprintf(filename,"%s/wall_normal_%i.dat",
712  doc_info.directory().c_str(),
713  doc_info.number());
714  some_file.open(filename);
715  nel=Wall_mesh_pt->nelement();
716  Vector<double> s(1);
717  Vector<double> x(2);
718  Vector<double> xi(1);
719  Vector<double> N(2);
720  for (unsigned e=0;e<nel;e++)
721  {
722  // Get pointer to element
723  FSIHermiteBeamElement* el_pt=
724  dynamic_cast<FSIHermiteBeamElement*>(Wall_mesh_pt->element_pt(e));
725 
726  // Loop over plot points
727  for (unsigned i=0;i<npts;i++)
728  {
729  s[0]=-1.0+2.0*double(i)/double(npts-1);
730 
731  // Get Eulerian position
732  el_pt->interpolated_x(s,x);
733 
734  // Get unit normal
735  el_pt->get_normal(s,N);
736 
737  // Get Lagrangian coordinate
738  el_pt->interpolated_xi(s,xi);
739 
740  some_file << x[0] << " " << x[1] << " "
741  << N[0] << " " << N[1] << " "
742  << xi[0] << std::endl;
743  }
744  }
745  some_file.close();
746 
747 } // end_of_doc_solution
748 
749 
750 
751 
752 //======= start_of_main================================================
753 /// Driver code -- pass a command line argument if you want to run
754 /// the code in validation mode where it only performs a few steps
755 //=====================================================================
756 int main(int argc, char* argv[])
757 {
758  // Store command line arguments
759  CommandLineArgs::setup(argc,argv);
760 
761  //Parameters for the leaflet: x-position of root and height
762  double x_0 = 1.0;
763  double hleaflet=0.5;
764 
765  // Number of elements in various regions of mesh
766  unsigned nleft=6;
767  unsigned nright=18;
768  unsigned ny1=3;
769  unsigned ny2=3;
770 
771  // Dimensions of fluid mesh: length to the left and right of leaflet
772  // and total height
773  double lleft =1.0;
774  double lright=3.0;
775  double htot=1.0;
776 
777  //Build the problem
779  AlgebraicElement<RefineableQTaylorHoodElement<2> > >
780  problem(lleft,lright,hleaflet,
781  htot,nleft,nright,ny1,ny2,x_0);
782 
783  // Set up doc info
784  DocInfo doc_info;
785  doc_info.set_directory("RESLT");
786 
787  // Trace file
788  ofstream trace;
789  char filename[100];
790  sprintf(filename,"%s/trace.dat",doc_info.directory().c_str());
791  trace.open(filename);
792 
793 
794  // Number of timesteps (reduced for validation)
795  unsigned nstep=200;
796  if (CommandLineArgs::Argc>1)
797  {
798  nstep=2;
799  }
800 
801  //Timestep:
802  double dt=0.05;
803 
804  // Initialise timestep
805  problem.initialise_dt(dt);
806 
807  // Doc initial guess for steady solve
808  problem.doc_solution(doc_info,trace);
809  doc_info.number()++;
810 
811 
812  // Initial loop to increment the Reynolds number in sequence of steady solves
813  //---------------------------------------------------------------------------
814  unsigned n_increment=4;
815  // Just to one step for validation run
816  if (CommandLineArgs::Argc>1)
817  {
818  n_increment=1;
819  }
820 
821  // Set max. number of adaptations
822  unsigned max_adapt=3;
823 
825  for (unsigned i=0;i<n_increment;i++)
826  {
827  // Increase Re and ReSt (for St=1)
830 
831  // Solve the steady problem
832  std::cout << "Computing a steady solution for Re="
833  << Global_Physical_Variables::Re << std::endl;
834  problem.steady_newton_solve(max_adapt);
835  problem.doc_solution(doc_info,trace);
836  doc_info.number()++;
837  } // reached final Reynolds number
838 
839 
840 
841  // Proper time-dependent run
842  //--------------------------
843 
844  // Limit the number of adaptations during unsteady run to one per timestep
845  max_adapt=1;
846 
847  // Don't re-set the initial conditions when adapting the mesh
848  bool first = false;
849 
850  // Timestepping loop
851  for (unsigned istep=0;istep<nstep;istep++)
852  {
853  // Solve the problem
854  problem.unsteady_newton_solve(dt,max_adapt,first);
855 
856  // Output the solution
857  problem.doc_solution(doc_info,trace);
858 
859  // Step number
860  doc_info.number()++;
861 
862  }
863 
864 }//end of main
865 
866 
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
void doc_solution(DocInfo &doc_info, ofstream &trace)
Doc the solution.
RefineableAlgebraicChannelWithLeafletMesh< ELEMENT > * Fluid_mesh_pt
Pointer to the fluid mesh.
~FSIChannelWithLeafletProblem()
Destructor empty.
double Htot
Total height of the domain.
void actions_before_newton_solve()
Actions before solve (empty)
OneDLagrangianMesh< FSIHermiteBeamElement > * Wall_mesh_pt
Pointer to the "wall" mesh.
FSIChannelWithLeafletProblem(const double &lleft, const double &lright, const double &hleaflet, const double &htot, const unsigned &nleft, const unsigned &nright, const unsigned &ny1, const unsigned &ny2, const double &x_0)
Constructor: Pass the lenght of the domain at the left of the leaflet lleft,the lenght of the domain ...
OneDLagrangianMesh< FSIHermiteBeamElement > * wall_mesh_pt()
Access function to the wall mesh.
void actions_after_newton_solve()
Actions after solve (empty)
void actions_before_implicit_timestep()
Update the inflow velocity.
void actions_after_adapt()
Actions after adaptation.
GeomObject * Leaflet_pt
Pointer to the GeomObject that represents the wall.
RefineableAlgebraicChannelWithLeafletMesh< ELEMENT > * fluid_mesh_pt()
Access function to fluid mesh.
void actions_before_newton_convergence_check()
Update before checking Newton convergence: Update the nodal positions in the fluid mesh in response t...
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
void d2position(const Vector< double > &zeta, Vector< double > &r, DenseMatrix< double > &drdzeta, RankThreeTensor< double > &ddrdzeta) const
Posn vector and its 1st & 2nd derivatives w.r.t. to coordinates: = drdzeta(alpha,...
unsigned ngeom_data() const
Number of geometric Data in GeomObject: None.
void position(const unsigned &t, const Vector< double > &zeta, Vector< double > &r) const
Parametrised position on object: r(zeta). Evaluated at previous timestep. t=0: current time; t>0: pre...
void position(const Vector< double > &zeta, Vector< double > &r) const
Position vector at Lagrangian coordinate zeta.
double X0
x position of the undeformed leaflet's origin.
UndeformedLeaflet(const double &x0)
Constructor: argument is the x-coordinate of the leaflet.
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 on...
double ReSt
Womersley number: Product of Reynolds and Strouhal numbers.
double Period
Period for fluctuations in flux.
double Q
Fluid structure interaction parameter: Ratio of stresses used for non-dimensionalisation of fluid to ...
double flux(const double &t)
Flux: Pulsatile flow fluctuating between Min_flux and Max_flux with period Period.
double H
Non-dimensional wall thickness.