fsi_collapsible_channel_adapt.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-2023 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 #include <iostream>
27 
28 // Generic oomph-lib includes
29 #include "generic.h"
30 #include "navier_stokes.h"
31 #include "beam.h"
32 
33 // The wall mesh
34 #include "meshes/one_d_lagrangian_mesh.h"
35 
36 //Include the fluid mesh
37 #include "meshes/collapsible_channel_mesh.h"
38 
39 using namespace std;
40 
41 using namespace oomph;
42 
43 
44 
45 
46 //==========start_of_BL_Squash =========================================
47 /// Namespace to define the mapping [0,1] -> [0,1] that re-distributes
48 /// nodal points across the channel width.
49 //======================================================================
50 namespace BL_Squash
51 {
52 
53  /// Boundary layer width
54  double Delta=0.1;
55 
56  /// Fraction of points in boundary layer
57  double Fract_in_BL=0.5;
58 
59  /// Mapping [0,1] -> [0,1] that re-distributes
60  /// nodal points across the channel width
61  double squash_fct(const double& s)
62  {
63  // Default return
64  double y=s;
65  if (s<0.5*Fract_in_BL)
66  {
67  y=Delta*2.0*s/Fract_in_BL;
68  }
69  else if (s>1.0-0.5*Fract_in_BL)
70  {
71  y=2.0*Delta/Fract_in_BL*s+1.0-2.0*Delta/Fract_in_BL;
72  }
73  else
74  {
75  y=(1.0-2.0*Delta)/(1.0-Fract_in_BL)*s+
76  (Delta-0.5*Fract_in_BL)/(1.0-Fract_in_BL);
77  }
78 
79  return y;
80  }
81 }// end of BL_Squash
82 
83 
84 
85 /// ///////////////////////////////////////////////////////////////////////
86 /// ///////////////////////////////////////////////////////////////////////
87 /// ///////////////////////////////////////////////////////////////////////
88 
89 
90 //====start_of_underformed_wall============================================
91 /// Undeformed wall is a steady, straight 1D line in 2D space
92 /// \f[ x = X_0 + \zeta \f]
93 /// \f[ y = H \f]
94 //=========================================================================
95 class UndeformedWall : public GeomObject
96 {
97 
98 public:
99 
100  /// Constructor: arguments are the starting point and the height
101  /// above y=0.
102  UndeformedWall(const double& x0, const double& h): GeomObject(1,2)
103  {
104  X0=x0;
105  H=h;
106  }
107 
108 
109  /// Position vector at Lagrangian coordinate zeta
110  void position(const Vector<double>& zeta, Vector<double>& r) const
111  {
112  // Position Vector
113  r[0] = zeta[0]+X0;
114  r[1] = H;
115  }
116 
117 
118  /// Parametrised position on object: r(zeta). Evaluated at
119  /// previous timestep. t=0: current time; t>0: previous
120  /// timestep. Calls steady version.
121  void position(const unsigned& t, const Vector<double>& zeta,
122  Vector<double>& r) const
123  {
124  // Use the steady version
125  position(zeta,r);
126 
127  } // end of position
128 
129 
130  /// Posn vector and its 1st & 2nd derivatives
131  /// w.r.t. to coordinates:
132  /// \f$ \frac{dR_i}{d \zeta_\alpha}\f$ = drdzeta(alpha,i).
133  /// \f$ \frac{d^2R_i}{d \zeta_\alpha d \zeta_\beta}\f$ =
134  /// ddrdzeta(alpha,beta,i). Evaluated at current time.
135  virtual void d2position(const Vector<double>& zeta,
136  Vector<double>& r,
137  DenseMatrix<double> &drdzeta,
138  RankThreeTensor<double> &ddrdzeta) const
139  {
140  // Position vector
141  r[0] = zeta[0]+X0;
142  r[1] = H;
143 
144  // Tangent vector
145  drdzeta(0,0)=1.0;
146  drdzeta(0,1)=0.0;
147 
148  // Derivative of tangent vector
149  ddrdzeta(0,0,0)=0.0;
150  ddrdzeta(0,0,1)=0.0;
151 
152  } // end of d2position
153 
154  private :
155 
156  /// x position of the undeformed beam's left end.
157  double X0;
158 
159  /// Height of the undeformed wall above y=0.
160  double H;
161 
162 }; //end_of_undeformed_wall
163 
164 
165 
166 
167 
168 
169 
170 //====start_of_physical_parameters=====================
171 /// Namespace for phyical parameters
172 //======================================================
174 {
175  /// Reynolds number
176  double Re=50.0;
177 
178  /// Womersley = Reynolds times Strouhal
179  double ReSt=50.0;
180 
181  /// Default pressure on the left boundary
182  double P_up=0.0;
183 
184  /// Traction applied on the fluid at the left (inflow) boundary
185  void prescribed_traction(const double& t,
186  const Vector<double>& x,
187  const Vector<double>& n,
188  Vector<double>& traction)
189  {
190  traction.resize(2);
191  traction[0]=P_up;
192  traction[1]=0.0;
193 
194  } //end traction
195 
196  /// Non-dimensional wall thickness. As in Jensen & Heil (2003) paper.
197  double H=1.0e-2;
198 
199  /// 2nd Piola Kirchhoff pre-stress. As in Jensen & Heil (2003) paper.
200  double Sigma0=1.0e3;
201 
202  /// External pressure
203  double P_ext=0.0;
204 
205  /// Load function: Apply a constant external pressure to the wall.
206  /// Note: This is the load without the fluid contribution!
207  /// Fluid load gets added on by FSIWallElement.
208  void load(const Vector<double>& xi, const Vector<double>& x,
209  const Vector<double>& N, Vector<double>& load)
210  {
211  for(unsigned i=0;i<2;i++)
212  {
213  load[i] = -P_ext*N[i];
214  }
215  } //end of load
216 
217 
218  /// Fluid structure interaction parameter: Ratio of stresses used for
219  /// non-dimensionalisation of fluid to solid stresses.
220  double Q=1.0e-5;
221 
222 
223 } // end of namespace
224 
225 
226 
227 
228 //====start_of_problem_class==========================================
229 /// Problem class
230 //====================================================================
231 template <class ELEMENT>
232 class FSICollapsibleChannelProblem : public Problem
233 {
234 
235  public :
236 
237 /// Constructor: The arguments are the number of elements and
238 /// the lengths of the domain.
239  FSICollapsibleChannelProblem(const unsigned& nup,
240  const unsigned& ncollapsible,
241  const unsigned& ndown,
242  const unsigned& ny,
243  const double& lup,
244  const double& lcollapsible,
245  const double& ldown,
246  const double& ly);
247 
248  /// Destructor (empty)
250 
251 #ifdef MACRO_ELEMENT_NODE_UPDATE
252 
253  /// Access function for the specific bulk (fluid) mesh
255  {
256  // Upcast from pointer to the Mesh base class to the specific
257  // element type that we're using here.
258  return dynamic_cast<
260  (Bulk_mesh_pt);
261  }
262 
263 #else
264 
265  /// Access function for the specific bulk (fluid) mesh
267  {
268  // Upcast from pointer to the Mesh base class to the specific
269  // element type that we're using here.
270  return dynamic_cast<
272  (Bulk_mesh_pt);
273  }
274 
275 #endif
276 
277 
278  /// Access function for the wall mesh
280  {
281  return Wall_mesh_pt;
282 
283  } // end of access to wall mesh
284 
285 
286  /// Actions before adapt: Wipe the mesh of prescribed traction elements
287  void actions_before_adapt();
288 
289  /// Actions after adapt: Rebuild the mesh of prescribed traction elements
290  /// and reset FSI
291  void actions_after_adapt();
292 
293  /// Update the problem specs before solve (empty)
295 
296  /// Update the problem after solve (empty)
298 
299  /// Update before checking Newton convergence: Update the
300  /// nodal positions in the fluid mesh in response to possible
301  /// changes in the wall shape
303  {
304  Bulk_mesh_pt->node_update();
305  }
306 
307  /// Doc the solution
308  void doc_solution(DocInfo& doc_info,ofstream& trace_file);
309 
310  /// Apply initial conditions
312 
313 private :
314 
315  /// Create the prescribed traction elements on boundary b
316  void create_traction_elements(const unsigned &b,
317  Mesh* const &bulk_mesh_pt,
318  Mesh* const &traction_mesh_pt);
319 
320  /// Delete prescribed traction elements from the surface mesh
321  void delete_traction_elements(Mesh* const &traction_mesh_pt);
322 
323  /// Number of elements in the x direction in the upstream part of the channel
324  unsigned Nup;
325 
326  /// Number of elements in the x direction in the collapsible part of
327  /// the channel
328  unsigned Ncollapsible;
329 
330  /// Number of elements in the x direction in the downstream part of the channel
331  unsigned Ndown;
332 
333  /// Number of elements across the channel
334  unsigned Ny;
335 
336  /// x-length in the upstream part of the channel
337  double Lup;
338 
339  /// x-length in the collapsible part of the channel
340  double Lcollapsible;
341 
342  /// x-length in the downstream part of the channel
343  double Ldown;
344 
345  /// Transverse length
346  double Ly;
347 
348 #ifdef MACRO_ELEMENT_NODE_UPDATE
349 
350  /// Pointer to the "bulk" mesh
352 
353 #else
354 
355  /// Pointer to the "bulk" mesh
357 
358 #endif
359 
360  /// Pointer to the "surface" mesh that applies the traction at the
361  /// inflow
362  Mesh* Applied_fluid_traction_mesh_pt;
363 
364  /// Pointer to the "wall" mesh
366 
367  /// Pointer to the left control node
368  Node* Left_node_pt;
369 
370  /// Pointer to right control node
371  Node* Right_node_pt;
372 
373  /// Pointer to control node on the wall
374  Node* Wall_node_pt;
375 
376 };//end of problem class
377 
378 
379 
380 
381 //=====start_of_constructor======================================
382 /// Constructor for the collapsible channel problem
383 //===============================================================
384 template <class ELEMENT>
386  const unsigned& nup,
387  const unsigned& ncollapsible,
388  const unsigned& ndown,
389  const unsigned& ny,
390  const double& lup,
391  const double& lcollapsible,
392  const double& ldown,
393  const double& ly)
394 {
395  // Store problem parameters
396  Nup=nup;
397  Ncollapsible=ncollapsible;
398  Ndown=ndown;
399  Ny=ny;
400  Lup=lup;
401  Lcollapsible=lcollapsible;
402  Ldown=ldown;
403  Ly=ly;
404 
405 
406  // Overwrite maximum allowed residual to accomodate bad initial guesses
407  Problem::Max_residuals=1000.0;
408 
409  // Allocate the timestepper for the Navier-Stokes equations
410  BDF<2>* fluid_time_stepper_pt=new BDF<2>;
411 
412  // Add the fluid timestepper to the Problem's collection of timesteppers.
413  add_time_stepper_pt(fluid_time_stepper_pt);
414 
415  // Create a dummy Steady timestepper that stores two history values
416  Steady<2>* wall_time_stepper_pt = new Steady<2>;
417 
418  // Add the wall timestepper to the Problem's collection of timesteppers.
419  add_time_stepper_pt(wall_time_stepper_pt);
420 
421  // Geometric object that represents the undeformed wall:
422  // A straight line at height y=ly; starting at x=lup.
423  UndeformedWall* undeformed_wall_pt=new UndeformedWall(lup,ly);
424 
425  //Create the "wall" mesh with FSI Hermite beam elements, passing the
426  //dummy wall timestepper to the constructor
428  (Ncollapsible,Lcollapsible,undeformed_wall_pt,wall_time_stepper_pt);
429 
430  // Build a geometric object (one Lagrangian, two Eulerian coordinates)
431  // from the wall mesh
432  MeshAsGeomObject* wall_geom_object_pt=
433  new MeshAsGeomObject(Wall_mesh_pt);
434 
435 #ifdef MACRO_ELEMENT_NODE_UPDATE
436 
437  //Build bulk (fluid) mesh
438  Bulk_mesh_pt =
440  (nup, ncollapsible, ndown, ny,
441  lup, lcollapsible, ldown, ly,
442  wall_geom_object_pt,
443  fluid_time_stepper_pt);
444 
445  // Set a non-trivial boundary-layer-squash function...
446  Bulk_mesh_pt->bl_squash_fct_pt() = &BL_Squash::squash_fct;
447 
448  // ... and update the nodal positions accordingly
449  Bulk_mesh_pt->node_update();
450 
451 #else
452 
453  //Build bulk (fluid) mesh
454  Bulk_mesh_pt =
456  (nup, ncollapsible, ndown, ny,
457  lup, lcollapsible, ldown, ly,
458  wall_geom_object_pt,
460  fluid_time_stepper_pt);
461 
462 #endif
463 
464  if (CommandLineArgs::Argc>1)
465  {
466  // One round of uniform refinement
467  Bulk_mesh_pt->refine_uniformly();
468  }
469 
470  // Create "surface mesh" that will contain only the prescribed-traction
471  // elements. The constructor just creates the mesh without
472  // giving it any elements, nodes, etc.
473  Applied_fluid_traction_mesh_pt = new Mesh;
474 
475  // Create prescribed-traction elements from all elements that are
476  // adjacent to boundary 5 (left boundary), but add them to a separate mesh.
477  create_traction_elements(5,Bulk_mesh_pt,Applied_fluid_traction_mesh_pt);
478 
479  // Add the sub meshes to the problem
480  add_sub_mesh(Bulk_mesh_pt);
481  add_sub_mesh(Applied_fluid_traction_mesh_pt);
482  add_sub_mesh(Wall_mesh_pt);
483 
484  // Combine all submeshes into a single Mesh
485  build_global_mesh();
486 
487  //Set errror estimator
488  Z2ErrorEstimator* error_estimator_pt=new Z2ErrorEstimator;
489  bulk_mesh_pt()->spatial_error_estimator_pt()=error_estimator_pt;
490 
491 
492  // Complete build of fluid mesh
493  //-----------------------------
494 
495  // Loop over the elements to set up element-specific
496  // things that cannot be handled by constructor
497  unsigned n_element=Bulk_mesh_pt->nelement();
498  for(unsigned e=0;e<n_element;e++)
499  {
500  // Upcast from GeneralisedElement to the present element
501  ELEMENT* el_pt = dynamic_cast<ELEMENT*>(Bulk_mesh_pt->element_pt(e));
502 
503  //Set the Reynolds number
504  el_pt->re_pt() = &Global_Physical_Variables::Re;
505 
506  // Set the Womersley number
507  el_pt->re_st_pt() = &Global_Physical_Variables::ReSt;
508 
509  } // end loop over elements
510 
511 
512  // Pin redudant pressure dofs
513  RefineableNavierStokesEquations<2>::
514  pin_redundant_nodal_pressures(Bulk_mesh_pt->element_pt());
515 
516 
517  // Apply boundary conditions for fluid
518  //------------------------------------
519 
520 
521  //Pin the velocity on the boundaries
522  //x and y-velocities pinned along boundary 0 (bottom boundary) :
523  unsigned ibound=0;
524  unsigned num_nod= bulk_mesh_pt()->nboundary_node(ibound);
525  for (unsigned inod=0;inod<num_nod;inod++)
526  {
527  for(unsigned i=0;i<2;i++)
528  {
529  bulk_mesh_pt()->boundary_node_pt(ibound, inod)->pin(i);
530  }
531  }
532 
533  //x and y-velocities pinned along boundaries 2, 3, 4 (top boundaries) :
534  for(ibound=2;ibound<5;ibound++)
535  {
536  num_nod= bulk_mesh_pt()->nboundary_node(ibound);
537  for (unsigned inod=0;inod<num_nod;inod++)
538  {
539  for(unsigned i=0;i<2;i++)
540  {
541  bulk_mesh_pt()->boundary_node_pt(ibound, inod)->pin(i);
542  }
543  }
544  }
545 
546  //y-velocity pinned along boundary 1 (right boundary):
547  ibound=1;
548  num_nod= bulk_mesh_pt()->nboundary_node(ibound);
549  for (unsigned inod=0;inod<num_nod;inod++)
550  {
551  bulk_mesh_pt()->boundary_node_pt(ibound, inod)->pin(1);
552  }
553 
554 
555  //y-velocity pinned along boundary 5 (left boundary):
556  ibound=5;
557  num_nod= bulk_mesh_pt()->nboundary_node(ibound);
558  for (unsigned inod=0;inod<num_nod;inod++)
559  {
560 
561  bulk_mesh_pt()->boundary_node_pt(ibound, inod)->pin(1);
562 
563  }//end of pin_velocity
564 
565 
566 
567  // Complete build of applied traction elements
568  //--------------------------------------------
569 
570  // Loop over the traction elements to pass pointer to prescribed
571  // traction function
572  unsigned n_el=Applied_fluid_traction_mesh_pt->nelement();
573  for(unsigned e=0;e<n_el;e++)
574  {
575  // Upcast from GeneralisedElement to NavierStokes traction element
576  NavierStokesTractionElement<ELEMENT> *el_pt =
577  dynamic_cast< NavierStokesTractionElement<ELEMENT>*>(
578  Applied_fluid_traction_mesh_pt->element_pt(e));
579 
580  // Set the pointer to the prescribed traction function
581  el_pt->traction_fct_pt() = &Global_Physical_Variables::prescribed_traction;
582 
583  }
584 
585 
586 
587 
588  // Complete build of wall elements
589  //--------------------------------
590 
591  //Loop over the elements to set physical parameters etc.
592  n_element = wall_mesh_pt()->nelement();
593  for(unsigned e=0;e<n_element;e++)
594  {
595  // Upcast to the specific element type
596  FSIHermiteBeamElement *elem_pt =
597  dynamic_cast<FSIHermiteBeamElement*>(wall_mesh_pt()->element_pt(e));
598 
599  // Set physical parameters for each element:
600  elem_pt->sigma0_pt() = &Global_Physical_Variables::Sigma0;
601  elem_pt->h_pt() = &Global_Physical_Variables::H;
602 
603  // Set the load vector for each element
604  elem_pt->load_vector_fct_pt() = &Global_Physical_Variables::load;
605 
606  // Function that specifies the load ratios
607  elem_pt->q_pt() = &Global_Physical_Variables::Q;
608 
609  // Set the undeformed shape for each element
610  elem_pt->undeformed_beam_pt() = undeformed_wall_pt;
611 
612  // The normal on the wall elements as computed by the FSIHermiteElements
613  // points away from the fluid rather than into the fluid (as assumed
614  // by default)
615  elem_pt->set_normal_pointing_out_of_fluid();
616 
617  } // end of loop over elements
618 
619 
620 
621  // Boundary conditions for wall mesh
622  //----------------------------------
623 
624  // Set the boundary conditions: Each end of the beam is fixed in space
625  // Loop over the boundaries (ends of the beam)
626  for(unsigned b=0;b<2;b++)
627  {
628  // Pin displacements in both x and y directions
629  wall_mesh_pt()->boundary_node_pt(b,0)->pin_position(0);
630  wall_mesh_pt()->boundary_node_pt(b,0)->pin_position(1);
631  }
632 
633 
634 
635  //Choose control nodes
636  //---------------------
637 
638  // Left boundary
639  ibound=5;
640  num_nod= bulk_mesh_pt()->nboundary_node(ibound);
641  unsigned control_nod=num_nod/2;
642  Left_node_pt= bulk_mesh_pt()->boundary_node_pt(ibound, control_nod);
643 
644  // Right boundary
645  ibound=1;
646  num_nod= bulk_mesh_pt()->nboundary_node(ibound);
647  control_nod=num_nod/2;
648  Right_node_pt= bulk_mesh_pt()->boundary_node_pt(ibound, control_nod);
649 
650 
651  // Set the pointer to the control node on the wall
652  num_nod= wall_mesh_pt()->nnode();
653  Wall_node_pt=wall_mesh_pt()->node_pt(Ncollapsible/2);
654 
655 
656  // Setup FSI
657  //----------
658 
659  // The velocity of the fluid nodes on the wall (fluid mesh boundary 3)
660  // is set by the wall motion -- hence the no-slip condition needs to be
661  // re-applied whenever a node update is performed for these nodes.
662  // Such tasks may be performed automatically by the auxiliary node update
663  // function specified by a function pointer:
664  ibound=3;
665  num_nod= bulk_mesh_pt()->nboundary_node(ibound);
666  for (unsigned inod=0;inod<num_nod;inod++)
667  {
668  bulk_mesh_pt()->boundary_node_pt(ibound, inod)->
669  set_auxiliary_node_update_fct_pt(
670  FSI_functions::apply_no_slip_on_moving_wall);
671  }
672 
673 
674 
675  // Work out which fluid dofs affect the residuals of the wall elements:
676  // We pass the boundary between the fluid and solid meshes and
677  // pointers to the meshes. The interaction boundary is boundary 3 of the
678  // 2D fluid mesh.
679  FSI_functions::setup_fluid_load_info_for_solid_elements<ELEMENT,2>
680  (this,3,Bulk_mesh_pt,Wall_mesh_pt);
681 
682  // Setup equation numbering scheme
683  cout <<"Number of equations: " << assign_eqn_numbers() << std::endl;
684 
685 
686 }//end of constructor
687 
688 
689 
690 
691 //====start_of_doc_solution===================================================
692 /// Doc the solution
693 //============================================================================
694 template <class ELEMENT>
696  ofstream& trace_file)
697 {
698 
699 #ifdef MACRO_ELEMENT_NODE_UPDATE
700 
701  // Doc fsi
702  if (CommandLineArgs::Argc>1)
703  {
704  FSI_functions::doc_fsi<MacroElementNodeUpdateNode>
705  (Bulk_mesh_pt,Wall_mesh_pt,doc_info);
706  }
707 
708 #else
709 
710  // Doc fsi
711  if (CommandLineArgs::Argc>1)
712  {
713  FSI_functions::doc_fsi<AlgebraicNode>(Bulk_mesh_pt,Wall_mesh_pt,doc_info);
714  }
715 
716 #endif
717 
718  ofstream some_file;
719  char filename[100];
720 
721  // Number of plot points
722  unsigned npts;
723  npts=5;
724 
725  // Output fluid solution
726  sprintf(filename,"%s/soln%i.dat",doc_info.directory().c_str(),
727  doc_info.number());
728  some_file.open(filename);
729  bulk_mesh_pt()->output(some_file,npts);
730  some_file.close();
731 
732  // Document the wall shape
733  sprintf(filename,"%s/beam%i.dat",doc_info.directory().c_str(),
734  doc_info.number());
735  some_file.open(filename);
736  wall_mesh_pt()->output(some_file,npts);
737  some_file.close();
738 
739  // Loop over all elements do dump out previous solutions
740  // (get the number of previous timesteps available from the wall
741  // timestepper)
742  unsigned nsteps=time_stepper_pt(1)->nprev_values();
743  for (unsigned t=0;t<=nsteps;t++)
744  {
745  sprintf(filename,"%s/wall%i-%i.dat",doc_info.directory().c_str(),
746  doc_info.number(),t);
747  some_file.open(filename);
748  unsigned n_elem=wall_mesh_pt()->nelement();
749  for (unsigned ielem=0;ielem<n_elem;ielem++)
750  {
751  dynamic_cast<FSIHermiteBeamElement*>(wall_mesh_pt()->element_pt(ielem))->
752  output(t,some_file,npts);
753  }
754  some_file.close();
755  } // end of output of previous solutions
756 
757 
758 
759  // Write trace file
760  trace_file << time_pt()->time() << " "
761  << Wall_node_pt->x(1) << " "
762  << Left_node_pt->value(0) << " "
763  << Right_node_pt->value(0) << " "
765  << std::endl;
766 
767 } // end_of_doc_solution
768 
769 
770 
771 
772 
773 //=====start_of_create_traction_elements======================================
774 /// Create the traction elements
775 //============================================================================
776 template <class ELEMENT>
778  const unsigned &b, Mesh* const &bulk_mesh_pt, Mesh* const &traction_mesh_pt)
779 {
780 
781  // How many bulk elements are adjacent to boundary b?
782  unsigned n_element = bulk_mesh_pt->nboundary_element(b);
783 
784  // Loop over the bulk elements adjacent to boundary b?
785  for(unsigned e=0;e<n_element;e++)
786  {
787  // Get pointer to the bulk element that is adjacent to boundary b
788  ELEMENT* bulk_elem_pt = dynamic_cast<ELEMENT*>
789  (bulk_mesh_pt->boundary_element_pt(b,e));
790 
791  //What is the index of the face of element e along boundary
792  int face_index = bulk_mesh_pt->face_index_at_boundary(b,e);
793 
794  // Build the corresponding prescribed-traction element
795  NavierStokesTractionElement<ELEMENT>* flux_element_pt =
796  new NavierStokesTractionElement<ELEMENT>(bulk_elem_pt,face_index);
797 
798  //Add the prescribed-traction element to the surface mesh
799  traction_mesh_pt->add_element_pt(flux_element_pt);
800 
801  } //end of loop over bulk elements adjacent to boundary b
802 
803 } // end of create_traction_elements
804 
805 
806 
807 //============start_of_delete_traction_elements==============================
808 /// Delete traction elements and wipe the surface mesh
809 //=======================================================================
810 template<class ELEMENT>
812 delete_traction_elements(Mesh* const &surface_mesh_pt)
813 {
814  // How many surface elements are in the surface mesh
815  unsigned n_element = surface_mesh_pt->nelement();
816 
817  // Loop over the surface elements
818  for(unsigned e=0;e<n_element;e++)
819  {
820  // Kill surface element
821  delete surface_mesh_pt->element_pt(e);
822  }
823 
824  // Wipe the mesh
825  surface_mesh_pt->flush_element_and_node_storage();
826 
827 } // end of delete_traction_elements
828 
829 
830 
831 //====start_of_apply_initial_condition========================================
832 /// Apply initial conditions
833 //============================================================================
834 template <class ELEMENT>
836 {
837  // Check that timestepper is from the BDF family
838  if (time_stepper_pt()->type()!="BDF")
839  {
840  std::ostringstream error_stream;
841  error_stream << "Timestepper has to be from the BDF family!\n"
842  << "You have specified a timestepper from the "
843  << time_stepper_pt()->type() << " family" << std::endl;
844 
845  throw OomphLibError(error_stream.str(),
846  OOMPH_CURRENT_FUNCTION,
847  OOMPH_EXCEPTION_LOCATION);
848  }
849 
850  // Update the mesh
851  bulk_mesh_pt()->node_update();
852 
853  // Loop over the nodes to set initial guess everywhere
854  unsigned num_nod = bulk_mesh_pt()->nnode();
855  for (unsigned n=0;n<num_nod;n++)
856  {
857  // Get nodal coordinates
858  Vector<double> x(2);
859  x[0]=bulk_mesh_pt()->node_pt(n)->x(0);
860  x[1]=bulk_mesh_pt()->node_pt(n)->x(1);
861 
862  // Assign initial condition: Steady Poiseuille flow
863  bulk_mesh_pt()->node_pt(n)->set_value(0,6.0*(x[1]/Ly)*(1.0-(x[1]/Ly)));
864  bulk_mesh_pt()->node_pt(n)->set_value(1,0.0);
865  }
866 
867  // Assign initial values for an impulsive start
868  bulk_mesh_pt()->assign_initial_values_impulsive();
869  wall_mesh_pt()->assign_initial_values_impulsive();
870 
871 } // end of set_initial_condition
872 
873 
874 
875 
876 
877 //=========start_of_actions_before_adapt==================================
878 /// Actions before adapt: Wipe the mesh of prescribed traction elements
879 //========================================================================
880 template<class ELEMENT>
882 {
883  // Kill the traction elements and wipe surface mesh
884  delete_traction_elements(Applied_fluid_traction_mesh_pt);
885 
886  // Rebuild the global mesh.
887  rebuild_global_mesh();
888 
889 } // end of actions_before_adapt
890 
891 
892 
893 //==========start_of_actions_after_adapt==================================
894 /// Actions after adapt: Rebuild the mesh of prescribed traction elements
895 //========================================================================
896 template<class ELEMENT>
898 {
899  // Create prescribed-flux elements from all elements that are
900  // adjacent to boundary 5 and add them to surface mesh
901  create_traction_elements(5,Bulk_mesh_pt,Applied_fluid_traction_mesh_pt);
902 
903  // Rebuild the global mesh
904  rebuild_global_mesh();
905 
906  // Unpin all pressure dofs
907  RefineableNavierStokesEquations<2>::
908  unpin_all_pressure_dofs(Bulk_mesh_pt->element_pt());
909 
910  // Pin redundant pressure dofs
911  RefineableNavierStokesEquations<2>::
912  pin_redundant_nodal_pressures(Bulk_mesh_pt->element_pt());
913 
914  // Loop over the traction elements to pass pointer to prescribed
915  // traction function
916  unsigned n_element=Applied_fluid_traction_mesh_pt->nelement();
917  for(unsigned e=0;e<n_element;e++)
918  {
919  // Upcast from GeneralisedElement to NavierStokesTractionElement element
920  NavierStokesTractionElement<ELEMENT> *el_pt =
921  dynamic_cast<NavierStokesTractionElement<ELEMENT>*>(
922  Applied_fluid_traction_mesh_pt->element_pt(e));
923 
924  // Set the pointer to the prescribed traction function
925  el_pt->traction_fct_pt() = &Global_Physical_Variables::prescribed_traction;
926  }
927 
928  // The functions used to update the no slip boundary conditions
929  // must be set on any new nodes that have been created during the
930  // mesh adaptation process.
931  // There is no mechanism by which auxiliary update functions
932  // are copied to newly created nodes.
933  // (because, unlike boundary conditions, they don't occur exclusively
934  // at boundaries)
935 
936  // The velocity of the fluid nodes on the wall (fluid mesh boundary 3)
937  // is set by the wall motion -- hence the no-slip condition needs to be
938  // re-applied whenever a node update is performed for these nodes.
939  // Such tasks may be performed automatically by the auxiliary node update
940  // function specified by a function pointer:
941  unsigned ibound=3;
942  unsigned num_nod= bulk_mesh_pt()->nboundary_node(ibound);
943  for (unsigned inod=0;inod<num_nod;inod++)
944  {
945  bulk_mesh_pt()->boundary_node_pt(ibound, inod)->
946  set_auxiliary_node_update_fct_pt(
947  FSI_functions::apply_no_slip_on_moving_wall);
948  }
949 
950  // (Re-)setup fsi: Work out which fluid dofs affect wall elements
951  // the correspondance between wall dofs and fluid elements is handled
952  // during the remeshing, but the "reverse" association must be done
953  // separately. We need to set up the interaction every time because the fluid
954  // element adjacent to a given solid element's integration point may have
955  // changed.We pass the boundary between the fluid and solid meshes and
956  // pointers to the meshes. The interaction boundary is boundary 3 of
957  // the Fluid mesh.
958  FSI_functions::setup_fluid_load_info_for_solid_elements<ELEMENT,2>
959  (this,3,Bulk_mesh_pt,Wall_mesh_pt);
960 
961 
962 } // end of actions_after_adapt
963 
964 
965 
966 //============start_of_main====================================================
967 /// Driver code for a collapsible channel problem with FSI.
968 /// Presence of command line arguments indicates validation run with
969 /// coarse resolution and small number of timesteps.
970 //=============================================================================
971 int main(int argc, char* argv[])
972 {
973 
974  // Store command line arguments
975  CommandLineArgs::setup(argc,argv);
976 
977  // Reduction in resolution for validation run?
978  unsigned coarsening_factor=4;
979  if (CommandLineArgs::Argc>1)
980  {
981  coarsening_factor=4;
982  }
983 
984  // Number of elements in the domain
985  unsigned nup=20/coarsening_factor;
986  unsigned ncollapsible=40/coarsening_factor;
987  unsigned ndown=40/coarsening_factor;
988  unsigned ny=16/coarsening_factor;
989 
990  // Length of the domain
991  double lup=5.0;
992  double lcollapsible=10.0;
993  double ldown=10.0;
994  double ly=1.0;
995 
996  // Set external pressure (on the wall stiffness scale).
998 
999  // Pressure on the left boundary: This is consistent with steady
1000  // Poiseuille flow
1001  Global_Physical_Variables::P_up=12.0*(lup+lcollapsible+ldown);
1002 
1003 
1004 #ifdef MACRO_ELEMENT_NODE_UPDATE
1005 
1006 #ifdef TAYLOR_HOOD
1007 
1008  // Build the problem with QTaylorHoodElements
1010  <MacroElementNodeUpdateElement<RefineableQTaylorHoodElement<2> > >
1011  problem(nup, ncollapsible, ndown, ny,
1012  lup, lcollapsible, ldown, ly);
1013 
1014 #else
1015 
1016  // Build the problem with QCrouzeixRaviartElements
1018  <MacroElementNodeUpdateElement<RefineableQCrouzeixRaviartElement<2> > >
1019  problem(nup, ncollapsible, ndown, ny,
1020  lup, lcollapsible, ldown, ly);
1021 
1022 #endif
1023 
1024 #else
1025 
1026 #ifdef TAYLOR_HOOD
1027 
1028  // Build the problem with QTaylorHoodElements
1030  <AlgebraicElement<RefineableQTaylorHoodElement<2> > >
1031  problem(nup, ncollapsible, ndown, ny,
1032  lup, lcollapsible, ldown, ly);
1033 
1034 #else
1035 
1036  // Build the problem with QCrouzeixRaviartElements
1038  <AlgebraicElement<RefineableQCrouzeixRaviartElement<2> > >
1039  problem(nup, ncollapsible, ndown, ny,
1040  lup, lcollapsible, ldown, ly);
1041 
1042 #endif
1043 
1044 #endif
1045 
1046 
1047  // Timestep. Note: Preliminary runs indicate that the period of
1048  // the oscillation is about 1 so this gives us 40 steps per period.
1049  double dt=1.0/40.0;
1050 
1051  // Initial time for the simulation
1052  double t_min=0.0;
1053 
1054  // Maximum time for simulation
1055  double t_max=3.5;
1056 
1057  // Initialise timestep
1058  problem.time_pt()->time()=t_min;
1059  problem.initialise_dt(dt);
1060 
1061  // Apply initial condition
1062  problem.set_initial_condition();
1063 
1064  //Set output directory
1065  DocInfo doc_info;
1066  doc_info.set_directory("RESLT");
1067 
1068  // Open a trace file
1069  ofstream trace_file;
1070  char filename[100];
1071  sprintf(filename,"%s/trace.dat",doc_info.directory().c_str());
1072  trace_file.open(filename);
1073 
1074  // Output the initial condition
1075  problem.doc_solution(doc_info, trace_file);
1076 
1077  // Increment step number
1078  doc_info.number()++;
1079 
1080  // Find number of timesteps (reduced for validation)
1081  unsigned nstep = unsigned((t_max-t_min)/dt);
1082  if (CommandLineArgs::Argc>1)
1083  {
1084  nstep=3;
1085  }
1086 
1087 
1088  // Set targets for spatial adaptivity
1089  problem.bulk_mesh_pt()->max_permitted_error()=1.0e-3;
1090  problem.bulk_mesh_pt()->min_permitted_error()=1.0e-5;
1091 
1092  // Overwrite for validation run
1093  if (CommandLineArgs::Argc>1)
1094  {
1095  // Set targets for spatial adaptivity
1096  problem.bulk_mesh_pt()->max_permitted_error()=0.5e-2;
1097  problem.bulk_mesh_pt()->min_permitted_error()=0.5e-4;
1098  }
1099 
1100  // When performing the first timestep, we can adapt the mesh as many times
1101  // as we want because the initial condition can be re-set
1102  unsigned max_adapt=3;
1103  bool first=true;
1104 
1105  // Timestepping loop
1106  for (unsigned istep=0;istep<nstep;istep++)
1107  {
1108  // Solve the problem
1109  problem.unsteady_newton_solve(dt,max_adapt,first);
1110 
1111  // Outpt the solution
1112  problem.doc_solution(doc_info, trace_file);
1113 
1114  // Step number
1115  doc_info.number()++;
1116 
1117  // We've done the first step
1118  first=false;
1119  max_adapt=1;
1120  }
1121 
1122 
1123  // Close trace file.
1124  trace_file.close();
1125 
1126 }//end of main
1127 
MacroElementNodeUpdateRefineableCollapsibleChannelMesh< ELEMENT > * Bulk_mesh_pt
Pointer to the "bulk" mesh.
void actions_after_newton_solve()
Update the problem after solve (empty)
MacroElementNodeUpdateRefineableCollapsibleChannelMesh< ELEMENT > * bulk_mesh_pt()
Access function for the specific bulk (fluid) mesh.
void actions_before_adapt()
Actions before adapt: Wipe the mesh of prescribed traction elements.
RefineableAlgebraicCollapsibleChannelMesh< ELEMENT > * Bulk_mesh_pt
Pointer to the "bulk" mesh.
void actions_before_newton_convergence_check()
Update before checking Newton convergence: Update the nodal positions in the fluid mesh in response t...
void actions_before_newton_solve()
Update the problem specs before solve (empty)
void actions_after_adapt()
Actions after adapt: Rebuild the mesh of prescribed traction elements and reset FSI.
OneDLagrangianMesh< FSIHermiteBeamElement > * wall_mesh_pt()
Access function for the wall mesh.
void delete_traction_elements(Mesh *const &traction_mesh_pt)
Delete prescribed traction elements from the surface mesh.
void create_traction_elements(const unsigned &b, Mesh *const &bulk_mesh_pt, Mesh *const &traction_mesh_pt)
Create the prescribed traction elements on boundary b.
RefineableAlgebraicCollapsibleChannelMesh< ELEMENT > * bulk_mesh_pt()
Access function for the specific bulk (fluid) mesh.
void set_initial_condition()
Apply initial conditions.
FSICollapsibleChannelProblem(const unsigned &nup, const unsigned &ncollapsible, const unsigned &ndown, const unsigned &ny, const double &lup, const double &lcollapsible, const double &ldown, const double &ly)
Constructor: The arguments are the number of elements and the lengths of the domain.
void doc_solution(DocInfo &doc_info, ofstream &trace_file)
Doc the solution.
/////////////////////////////////////////////////////////////////////// /////////////////////////////...
virtual 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,...
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.
UndeformedWall(const double &x0, const double &h)
Constructor: arguments are the starting point and the height above y=0.
///////////////////////////////////////////////////////////////////////// ///////////////////////////...
//////////////////////////////////////////////////////////////////////// ////////////////////////////...
int main(int argc, char *argv[])
Driver code for a collapsible channel problem with FSI. Presence of command line arguments indicates ...
Namespace to define the mapping [0,1] -> [0,1] that re-distributes nodal points across the channel wi...
double squash_fct(const double &s)
Mapping [0,1] -> [0,1] that re-distributes nodal points across the channel width.
double Delta
Boundary layer width.
double Fract_in_BL
Fraction of points in boundary layer.
Namespace for phyical parameters.
double P_ext
External pressure.
double ReSt
Womersley = Reynolds times Strouhal.
void 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.
void load(const Vector< double > &xi, const Vector< double > &x, const Vector< double > &N, Vector< double > &load)
Load function: Apply a constant external pressure to the wall. Note: This is the load without the flu...
double Sigma0
2nd Piola Kirchhoff pre-stress. As in Jensen & Heil (2003) paper.
double Q
Fluid structure interaction parameter: Ratio of stresses used for non-dimensionalisation of fluid to ...
double P_up
Default pressure on the left boundary.
double H
Non-dimensional wall thickness. As in Jensen & Heil (2003) paper.
////////////////////////////////////////////////////////////////////// //////////////////////////////...