refineable_two_layer_interface.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 
27 //Generic routines
28 #include "generic.h"
29 
30 
31 // The equations
32 #include "navier_stokes.h"
33 #include "solid.h"
34 #include "constitutive.h"
35 #include "fluid_interface.h"
36 
37 // The mesh
38 #include "meshes/triangle_mesh.h"
39 
40 using namespace std;
41 using namespace oomph;
42 
43 
44 namespace oomph
45 {
46 
47 
48 //==============================================================
49 /// Overload CrouzeixRaviart element to modify output
50 //==============================================================
52  public virtual PseudoSolidNodeUpdateElement<TCrouzeixRaviartElement<2>,
53  TPVDBubbleEnrichedElement<2,3> >
54  {
55 
56  private:
57 
58  /// Storage for elemental error estimate -- used for post-processing
59  double Error;
60 
61  public:
62 
63  /// Constructor initialise error
64  MyCrouzeixRaviartElement() : PseudoSolidNodeUpdateElement<TCrouzeixRaviartElement<2>,
65  TPVDBubbleEnrichedElement<2,3> >()
66  {
67  Error=0.0;
68  }
69 
70  /// Set error value for post-processing
71  void set_error(const double& error){Error=error;}
72 
73  /// Return variable identifier
74  std::string variable_identifier()
75  {
76  std::string txt="VARIABLES=";
77  txt+="\"x\",";
78  txt+="\"y\",";
79  txt+="\"u\",";
80  txt+="\"v\",";
81  txt+="\"p\",";
82  txt+="\"du/dt\",";
83  txt+="\"dv/dt\",";
84  txt+="\"u_m\",";
85  txt+="\"v_m\",";
86  txt+="\"x_h1\",";
87  txt+="\"y_h1\",";
88  txt+="\"x_h2\",";
89  txt+="\"y_h2\",";
90  txt+="\"u_h1\",";
91  txt+="\"v_h1\",";
92  txt+="\"u_h2\",";
93  txt+="\"v_h2\",";
94  txt+="\"error\",";
95  txt+="\"size\",";
96  txt+="\n";
97  return txt;
98  }
99 
100 
101  /// Overload output function
102  void output(std::ostream &outfile,
103  const unsigned &nplot)
104  {
105 
106  // Assign dimension
107  unsigned el_dim=2;
108 
109  // Vector of local coordinates
110  Vector<double> s(el_dim);
111 
112  // Acceleration
113  Vector<double> dudt(el_dim);
114 
115  // Mesh elocity
116  Vector<double> mesh_veloc(el_dim,0.0);
117 
118  // Tecplot header info
119  outfile << tecplot_zone_string(nplot);
120 
121  // Find out how many nodes there are
122  unsigned n_node = nnode();
123 
124  //Set up memory for the shape functions
125  Shape psif(n_node);
126  DShape dpsifdx(n_node,el_dim);
127 
128  // Loop over plot points
129  unsigned num_plot_points=nplot_points(nplot);
130  for (unsigned iplot=0;iplot<num_plot_points;iplot++)
131  {
132 
133  // Get local coordinates of plot point
134  get_s_plot(iplot,nplot,s);
135 
136  //Call the derivatives of the shape and test functions
137  dshape_eulerian(s,psif,dpsifdx);
138 
139  //Allocate storage
140  Vector<double> mesh_veloc(el_dim);
141  Vector<double> dudt(el_dim);
142  Vector<double> dudt_ALE(el_dim);
143  DenseMatrix<double> interpolated_dudx(el_dim,el_dim);
144 
145  //Initialise everything to zero
146  for(unsigned i=0;i<el_dim;i++)
147  {
148  mesh_veloc[i]=0.0;
149  dudt[i]=0.0;
150  dudt_ALE[i]=0.0;
151  for(unsigned j=0;j<el_dim;j++)
152  {
153  interpolated_dudx(i,j) = 0.0;
154  }
155  }
156 
157  //Calculate velocities and derivatives
158 
159  //Loop over directions
160  for(unsigned i=0;i<el_dim;i++)
161  {
162  //Get the index at which velocity i is stored
163  unsigned u_nodal_index = u_index_nst(i);
164  // Loop over nodes
165  for(unsigned l=0;l<n_node;l++)
166  {
167  dudt[i]+=du_dt_nst(l,u_nodal_index)*psif[l];
168  mesh_veloc[i]+=dnodal_position_dt(l,i)*psif[l];
169 
170  //Loop over derivative directions for velocity gradients
171  for(unsigned j=0;j<el_dim;j++)
172  {
173  interpolated_dudx(i,j) += nodal_value(l,u_nodal_index)*
174  dpsifdx(l,j);
175  }
176  }
177  }
178 
179 
180  // Get dudt in ALE form (incl mesh veloc)
181  for(unsigned i=0;i<el_dim;i++)
182  {
183  dudt_ALE[i]=dudt[i];
184  for (unsigned k=0;k<el_dim;k++)
185  {
186  dudt_ALE[i]-=mesh_veloc[k]*interpolated_dudx(i,k);
187  }
188  }
189 
190 
191  // Coordinates
192  for(unsigned i=0;i<el_dim;i++)
193  {
194  outfile << interpolated_x(s,i) << " ";
195  }
196 
197  // Velocities
198  for(unsigned i=0;i<el_dim;i++)
199  {
200  outfile << interpolated_u_nst(s,i) << " ";
201  }
202 
203  // Pressure
204  outfile << interpolated_p_nst(s) << " ";
205 
206  // Accelerations
207  for(unsigned i=0;i<el_dim;i++)
208  {
209  outfile << dudt_ALE[i] << " ";
210  }
211 
212  // Mesh velocity
213  for(unsigned i=0;i<el_dim;i++)
214  {
215  outfile << mesh_veloc[i] << " ";
216  }
217 
218  // History values of coordinates
219  unsigned n_prev=node_pt(0)->position_time_stepper_pt()->ntstorage();
220  for (unsigned t=1;t<n_prev;t++)
221  {
222  for(unsigned i=0;i<el_dim;i++)
223  {
224  outfile << interpolated_x(t,s,i) << " ";
225  }
226  }
227 
228  // History values of velocities
229  n_prev=node_pt(0)->time_stepper_pt()->ntstorage();
230  for (unsigned t=1;t<n_prev;t++)
231  {
232  for(unsigned i=0;i<el_dim;i++)
233  {
234  outfile << interpolated_u_nst(t,s,i) << " ";
235  }
236  }
237 
238  outfile << Error << " "
239  << size() << std::endl;
240  }
241 
242  // Write tecplot footer (e.g. FE connectivity lists)
243  write_tecplot_zone_footer(outfile,nplot);
244  }
245 
246  };
247 
248 
249 //=======================================================================
250 /// Face geometry for element is the same as that for the underlying
251 /// wrapped element
252 //=======================================================================
253  template<>
254  class FaceGeometry<MyCrouzeixRaviartElement>
255  : public virtual SolidTElement<1,3>
256  {
257  public:
258  FaceGeometry() : SolidTElement<1,3>() {}
259  };
260 
261 
262 //=======================================================================
263 /// Face geometry of Face geometry for element is the same
264 /// as that for the underlying
265 /// wrapped element
266 //=======================================================================
267  template<>
268  class FaceGeometry<FaceGeometry<MyCrouzeixRaviartElement> >
269  : public virtual SolidPointElement
270  {
271  public:
272  FaceGeometry() : SolidPointElement() {}
273  };
274 
275 
276 } //End of namespace extension
277 
278 
279 
280 /// ////////////////////////////////////////////////////////
281 /// ////////////////////////////////////////////////////////
282 /// ////////////////////////////////////////////////////////
283 
284 
285 //==start_of_namespace==============================
286 /// Namespace for Problem Parameter
287 //==================================================
289  {
290  /// Doc info
291  DocInfo Doc_info;
292 
293  /// Reynolds number
294  double Re=5.0;
295 
296  /// Strouhal number
297  double St = 1.0;
298 
299  /// Womersley number (Reynolds x Strouhal)
300  double ReSt = 5.0;
301 
302  /// Product of Reynolds number and inverse of Froude number
303  double ReInvFr = 5.0;
304 
305  /// Ratio of viscosity in upper fluid to viscosity in lower
306  /// fluid. Reynolds number etc. is based on viscosity in lower fluid.
307  double Viscosity_Ratio = 0.1;
308 
309  /// Ratio of density in upper fluid to density in lower
310  /// fluid. Reynolds number etc. is based on density in lower fluid.
311  double Density_Ratio = 0.5;
312 
313  /// Capillary number
314  double Ca = 0.01;
315 
316  /// Direction of gravity
317  Vector<double> G(2);
318 
319  /// Pseudo-solid Poisson ratio
320  double Nu=0.1;
321 
322  /// Length of the channel
323  double Length = 1.0;
324 
325  // Height of the channel
326  double Height = 2.0;
327 
328  // Relative Height of the interface
329  double Interface_Height = 0.5;
330 
331  /// Constitutive law used to determine the mesh deformation
332  ConstitutiveLaw *Constitutive_law_pt=0;
333 
334  /// Trace file
335  ofstream Trace_file;
336 
337  /// File to document the norm of the solution (for validation
338  /// purposes -- triangle doesn't give fully reproducible results so
339  /// mesh generation/adaptation may generate slightly different numbers
340  /// of elements on different machines!)
341  ofstream Norm_file;
342 
343  } // end_of_namespace
344 
345 
346 
347 //==start_of_problem_class============================================
348 /// Problem class to simulate viscous drop propagating along 2D channel
349 //====================================================================
350 template<class ELEMENT>
351 class TwoLayerInterfaceProblem : public Problem
352 {
353 
354 public:
355 
356  /// Constructor
358 
359  /// Destructor
361  {
362  // Fluid timestepper
363  delete this->time_stepper_pt(0);
364 
365  // Kill data associated with outer boundary
366  unsigned n=Outer_boundary_polyline_pt->npolyline();
367  for (unsigned j=0;j<n;j++)
368  {
369  delete Outer_boundary_polyline_pt->polyline_pt(j);
370  }
371  delete Outer_boundary_polyline_pt;
372 
373  // Flush element of free surface elements
374  delete_free_surface_elements();
375  delete Free_surface_mesh_pt;
376 
377  // Delete error estimator
378  delete Fluid_mesh_pt->spatial_error_estimator_pt();
379 
380  // Delete fluid mesh
381  delete Fluid_mesh_pt;
382 
383  // Delete the global pressure drop data
384  delete Drop_pressure_data_pt;
385 
386  // Kill const eqn
388 
389  }
390 
391 
392  /// Actions before adapt: Wipe the mesh of free surface elements
394  {
395  // Kill the elements and wipe surface mesh
396  delete_free_surface_elements();
397 
398  // Rebuild the Problem's global mesh from its various sub-meshes
399  this->rebuild_global_mesh();
400 
401  }// end of actions_before_adapt
402 
403 
404  /// Actions after adapt: Rebuild the mesh of free surface elements
406  {
407  // Create the elements that impose the displacement constraint
408  create_free_surface_elements();
409 
410  // Rebuild the Problem's global mesh from its various sub-meshes
411  this->rebuild_global_mesh();
412 
413  // Setup the problem again -- remember that fluid mesh has been
414  // completely rebuilt and its element's don't have any
415  // pointers to Re etc. yet
416  complete_problem_setup();
417 
418  }// end of actions_after_adapt
419 
420 
421  /// Update the after solve (empty)
423 
424  /// Update the problem specs before solve
426  {
427  //Reset the Lagrangian coordinates of the nodes to be the current
428  //Eulerian coordinates -- this makes the current configuration
429  //stress free
430  Fluid_mesh_pt->set_lagrangian_nodal_coordinates();
431  }
432 
433  /// Set boundary conditions and complete the build of all elements
434  void complete_problem_setup();
435 
436  /// Doc the solution
437  void doc_solution(const std::string& comment="");
438 
439  /// Compute the error estimates and assign to elements for plotting
440  void compute_error_estimate(double& max_err,
441  double& min_err);
442 
443 /// Set the initial conditions
445 {
446  // Determine number of nodes in mesh
447  const unsigned n_node = Fluid_mesh_pt->nnode();
448 
449  // Loop over all nodes in mesh
450  for(unsigned n=0;n<n_node;n++)
451  {
452  // Loop over the two velocity components
453  for(unsigned i=0;i<2;i++)
454  {
455  // Set velocity component i of node n to zero
456  Fluid_mesh_pt->node_pt(n)->set_value(i,0.0);
457  }
458  }
459 
460  // Initialise the previous velocity values for timestepping
461  // corresponding to an impulsive start
462  assign_initial_values_impulsive();
463 
464 } // End of set_initial_condition
465 
466 
467 
468 
469  /// Function to deform the interface
470 void deform_interface(const double &epsilon,const unsigned &n_periods)
471  {
472  // Determine number of nodes in the "bulk" mesh
473  const unsigned n_node = Fluid_mesh_pt->nnode();
474 
475  // Loop over all nodes in mesh
476  for(unsigned n=0;n<n_node;n++)
477  {
478  // Determine eulerian position of node
479  const double current_x_pos = Fluid_mesh_pt->node_pt(n)->x(0);
480  const double current_y_pos = Fluid_mesh_pt->node_pt(n)->x(1);
481 
482  // Determine new vertical position of node
483  const double new_y_pos = current_y_pos
484  + (1.0-fabs(1.0-current_y_pos))*epsilon
485  *(cos(2.0*n_periods*MathematicalConstants::Pi*current_x_pos/Problem_Parameter::Length));
486 
487  // Set new position
488  Fluid_mesh_pt->node_pt(n)->x(1) = new_y_pos;
489  }
490 } // End of deform_interface
491 
492 
493 
494  /// Fix pressure in element e at pressure dof pdof and set to pvalue
495  void fix_pressure(const unsigned &e,
496  const unsigned &pdof,
497  const double &pvalue)
498  {
499  // Fix the pressure at that element
500  dynamic_cast<ELEMENT*>(Fluid_mesh_pt->element_pt(e))->
501  fix_pressure(pdof,pvalue);
502  }
503 
504 
505 private:
506 
507 
508  /// Create free surface elements
509  void create_free_surface_elements();
510 
511  /// Delete free surface elements
513  {
514  // How many surface elements are in the surface mesh
515  unsigned n_element = Free_surface_mesh_pt->nelement();
516 
517  // Loop over the surface elements
518  for(unsigned e=0;e<n_element;e++)
519  {
520  // Kill surface element
521  delete Free_surface_mesh_pt->element_pt(e);
522  }
523 
524 
525  // Wipe the mesh
526  Free_surface_mesh_pt->flush_element_and_node_storage();
527 
528 
529  } // end of delete_free_surface_elements
530 
531 
532  /// Pointers to mesh of free surface elements
534 
535 
536  /// Pointer to Fluid_mesh
537  RefineableSolidTriangleMesh<ELEMENT>* Fluid_mesh_pt;
538 
539  /// Vector storing pointer to the drop polygons
540  Vector<TriangleMeshPolygon*> Drop_polygon_pt;
541 
542  /// Triangle mesh polygon for outer boundary
543  TriangleMeshPolygon* Outer_boundary_polyline_pt;
544 
545  /// Pointer to a global drop pressure datum
547 
548  /// Backed up drop pressure between adaptations
550 
551  /// Pointer to hijacked element
553 
554  /// Bool to indicate if volume constraint is applied (only for steady run)
556 
557  /// Enumeration of channel boundaries
558  enum
559  {
560  Inflow_boundary_id=0,
561  Upper_wall_boundary_id=1,
562  Outflow_boundary_id=2,
563  Bottom_wall_boundary_id=3,
564  Interface_boundary_id=4,
565  };
566 
567 
568 }; // end_of_problem_class
569 
570 
571 //==start_constructor=====================================================
572 /// Constructor
573 //========================================================================
574 template<class ELEMENT>
576 {
577  // Output directory
578  Problem_Parameter::Doc_info.set_directory("RESLT");
579 
580  // Allocate the timestepper -- this constructs the Problem's
581  // time object with a sufficient amount of storage to store the
582  // previous timsteps.
583  this->add_time_stepper_pt(new BDF<2>);
584 
585  // Build the boundary segments for outer boundary, consisting of
586  //--------------------------------------------------------------
587  // four separate polylines
588  //------------------------
589  Vector<TriangleMeshCurveSection*> boundary_polyline_pt(4);
590 
591  // Each polyline only has three vertices -- provide storage for their
592  // coordinates
593  Vector<Vector<double> > vertex_coord(3);
594  for(unsigned i=0;i<3;i++)
595  {
596  vertex_coord[i].resize(2);
597  }
598 
599  const double H = Problem_Parameter::Height;
600  const double h = Problem_Parameter::Interface_Height*H;
601 
602  // First polyline: Inflow
603  vertex_coord[0][0]=0.0;
604  vertex_coord[0][1]=0.0;
605  vertex_coord[1][0]=0.0;
606  vertex_coord[1][1]= h;
607  vertex_coord[2][0]=0.0;
608  vertex_coord[2][1]= H;
609 
610  // Build the 1st boundary polyline
611  boundary_polyline_pt[0] = new TriangleMeshPolyLine(vertex_coord,
612  Inflow_boundary_id);
613 
614  // Second boundary polyline: Upper wall
615  vertex_coord[0][0]=0.0;
616  vertex_coord[0][1]=H;
617  vertex_coord[1][0]=0.5*Problem_Parameter::Length;
618  vertex_coord[1][1]=H;
619  vertex_coord[2][0]=Problem_Parameter::Length;
620  vertex_coord[2][1]=H;
621 
622  // Build the 2nd boundary polyline
623  boundary_polyline_pt[1] = new TriangleMeshPolyLine(vertex_coord,
624  Upper_wall_boundary_id);
625 
626  // Third boundary polyline: Outflow
627  vertex_coord[0][0]=Problem_Parameter::Length;
628  vertex_coord[0][1]=H;
629  vertex_coord[1][0]=Problem_Parameter::Length;
630  vertex_coord[1][1]=h;
631  vertex_coord[2][0]=Problem_Parameter::Length;
632  vertex_coord[2][1]=0.0;
633 
634  // Build the 3rd boundary polyline
635  boundary_polyline_pt[2] = new TriangleMeshPolyLine(vertex_coord,
636  Outflow_boundary_id);
637 
638  // Fourth boundary polyline: Bottom wall
639  vertex_coord[0][0]=Problem_Parameter::Length;
640  vertex_coord[0][1]=0.0;
641  vertex_coord[1][0]=0.5*Problem_Parameter::Length;
642  vertex_coord[1][1]=0.0;
643  vertex_coord[2][0]=0.0;
644  vertex_coord[2][1]=0.0;
645 
646  // Build the 4th boundary polyline
647  boundary_polyline_pt[3] = new TriangleMeshPolyLine(vertex_coord,
648  Bottom_wall_boundary_id);
649 
650  // Create the triangle mesh polygon for outer boundary
651  Outer_boundary_polyline_pt = new TriangleMeshPolygon(boundary_polyline_pt);
652 
653  //Here we need to put the dividing internal line in
654  Vector<TriangleMeshOpenCurve *> interface_pt(1);
655  //Set the vertex coordinates
656  vertex_coord[0][0]=0.0;
657  vertex_coord[0][1]=h;
658  vertex_coord[1][0]=0.5*Problem_Parameter::Length;
659  vertex_coord[1][1]=h;
660  vertex_coord[2][0]=Problem_Parameter::Length;
661  vertex_coord[2][1]=h;
662 
663 //Create the internal line
664  TriangleMeshPolyLine* interface_polyline_pt =
665  new TriangleMeshPolyLine(vertex_coord,
666  Interface_boundary_id);
667 
668  // Do the connection with the destination boundary, in this case
669  // the connection is done with the inflow boundary
670  interface_polyline_pt->connect_initial_vertex_to_polyline(
671  dynamic_cast<TriangleMeshPolyLine*>(boundary_polyline_pt[0]),1);
672 
673  // Do the connection with the destination boundary, in this case
674  // the connection is done with the outflow boundary
675  interface_polyline_pt->connect_final_vertex_to_polyline(
676  dynamic_cast<TriangleMeshPolyLine*>(boundary_polyline_pt[2]),1);
677 
678  Vector<TriangleMeshCurveSection*> interface_curve_pt(1);
679  interface_curve_pt[0] = interface_polyline_pt;
680 
681  interface_pt[0] = new TriangleMeshOpenCurve(interface_curve_pt);
682 
683  // Now build the mesh, based on the boundaries specified by
684  //---------------------------------------------------------
685  // polygons just created
686  //----------------------
687 
688  // Convert to "closed curve" objects
689  TriangleMeshClosedCurve* outer_closed_curve_pt=Outer_boundary_polyline_pt;
690 
691  unsigned n_internal_closed_boundaries = 0;
692  Vector<TriangleMeshClosedCurve *>
693  inner_boundaries_pt(n_internal_closed_boundaries);
694 
695  // Target area for initial mesh
696  double uniform_element_area=0.01;
697 
698  // Use the TriangleMeshParameter object for gathering all
699  // the necessary arguments for the TriangleMesh object
700  TriangleMeshParameters triangle_mesh_parameters(
701  outer_closed_curve_pt);
702 
703  //Define the inner boundaries
704  triangle_mesh_parameters.internal_closed_curve_pt() = inner_boundaries_pt;
705 
706  // Define the holes on the boundary
707  triangle_mesh_parameters.internal_open_curves_pt() = interface_pt;
708 
709  // Define the maximum element area
710  triangle_mesh_parameters.element_area() =
711  uniform_element_area;
712 
713  Vector<double> lower_region(2);
714  lower_region[0] = 0.5*Problem_Parameter::Length;
715  lower_region[1] = 0.5*Problem_Parameter::Interface_Height;
716 
717  // Define the region
718  triangle_mesh_parameters.add_region_coordinates(1, lower_region);
719 
720  // Create the mesh
721  Fluid_mesh_pt =
722  new RefineableSolidTriangleMesh<ELEMENT>(
723  triangle_mesh_parameters, this->time_stepper_pt());
724 
725  // Set error estimator for bulk mesh
726  Z2ErrorEstimator* error_estimator_pt=new Z2ErrorEstimator;
727  Fluid_mesh_pt->spatial_error_estimator_pt()=error_estimator_pt;
728 
729  // Set targets for spatial adaptivity
730  Fluid_mesh_pt->max_permitted_error()=0.005;
731  Fluid_mesh_pt->min_permitted_error()=0.001;
732  Fluid_mesh_pt->max_element_size()=0.2;
733  Fluid_mesh_pt->min_element_size()=0.001;
734 
735  // Use coarser mesh during validation
736  if (CommandLineArgs::command_line_flag_has_been_set("--validation"))
737  {
738  Fluid_mesh_pt->min_element_size()=0.01;
739  }
740 
741  // Output boundary and mesh initial mesh for information
742  this->Fluid_mesh_pt->output_boundaries("boundaries.dat");
743  this->Fluid_mesh_pt->output("mesh.dat");
744 
745  // Set boundary condition and complete the build of all elements
746  complete_problem_setup();
747 
748  // Construct the mesh of free surface elements
749  Free_surface_mesh_pt=new Mesh;
750  create_free_surface_elements();
751 
752  // Combine meshes
753  //---------------
754 
755  // Add Fluid_mesh_pt sub meshes
756  this->add_sub_mesh(Fluid_mesh_pt);
757 
758  // Add Free_surface sub meshes
759  this->add_sub_mesh(this->Free_surface_mesh_pt);
760 
761  // Build global mesh
762  this->build_global_mesh();
763 
764  // Setup equation numbering scheme
765  cout <<"Number of equations: " << this->assign_eqn_numbers() << std::endl;
766 
767 } // end_of_constructor
768 
769 
770 //============start_of_create_free_surface_elements===============
771 /// Create elements that impose the kinematic and dynamic bcs
772 /// for the pseudo-solid fluid mesh
773 //=======================================================================
774 template<class ELEMENT>
776 {
777 
778  //Loop over the free surface boundaries
779  unsigned nb=Fluid_mesh_pt->nboundary();
780  for(unsigned b=4;b<nb;b++)
781  {
782  // Note: region is important
783  // How many bulk fluid elements are adjacent to boundary b in region 0?
784  unsigned n_element = Fluid_mesh_pt->nboundary_element_in_region(b,0);
785 
786  // Loop over the bulk fluid elements adjacent to boundary b?
787  for(unsigned e=0;e<n_element;e++)
788  {
789  // Get pointer to the bulk fluid element that is
790  // adjacent to boundary b in region 0
791  ELEMENT* bulk_elem_pt = dynamic_cast<ELEMENT*>(
792  Fluid_mesh_pt->boundary_element_in_region_pt(b,0,e));
793 
794  //Find the index of the face of element e along boundary b in region 0
795  int face_index = Fluid_mesh_pt->face_index_at_boundary_in_region(b,0,e);
796 
797  // Create new element
798  ElasticLineFluidInterfaceElement<ELEMENT>* el_pt =
799  new ElasticLineFluidInterfaceElement<ELEMENT>(
800  bulk_elem_pt,face_index);
801 
802  // Add it to the mesh
803  Free_surface_mesh_pt->add_element_pt(el_pt);
804 
805  //Add the appropriate boundary number
806  el_pt->set_boundary_number_in_bulk_mesh(b);
807 
808  //Specify the Strouhal number
809  el_pt->st_pt() = &Problem_Parameter::St;
810 
811  //Specify the capillary number
812  el_pt->ca_pt() = &Problem_Parameter::Ca;
813  }
814  }
815 }
816 // end of create_free_surface_elements
817 
818 
819 
820 
821 //==start_of_complete_problem_setup=======================================
822 /// Set boundary conditions and complete the build of all elements
823 //========================================================================
824 template<class ELEMENT>
826  {
827  // Re-set the boundary conditions for fluid problem: All nodes are
828  // free by default -- just pin the ones that have Dirichlet conditions
829  // here.
830  unsigned nbound=Fluid_mesh_pt->nboundary();
831  for(unsigned ibound=0;ibound<nbound;ibound++)
832  {
833  unsigned num_nod=Fluid_mesh_pt->nboundary_node(ibound);
834  for (unsigned inod=0;inod<num_nod;inod++)
835  {
836  // Get node
837  Node* nod_pt=Fluid_mesh_pt->boundary_node_pt(ibound,inod);
838 
839  //Pin both velocities on side boundaries (1 and 3)
840  if((ibound==1) || (ibound==3))
841  {
842  nod_pt->pin(0);
843  nod_pt->pin(1);
844  }
845 
846  //If it's the inflow or outflow pin only the horizontal velocity
847  if((ibound==0) || (ibound==2)) {nod_pt->pin(0);}
848 
849  // Pin horizontal pseudo-solid positions
850  SolidNode* solid_node_pt = dynamic_cast<SolidNode*>(nod_pt);
851  solid_node_pt->pin_position(0);
852 
853  //Pin the vertical position on the upper and lower walls
854  if((ibound==1) || (ibound==3))
855  {
856  solid_node_pt->pin_position(1);
857  }
858  else
859  {
860  solid_node_pt->unpin_position(1);
861  }
862  } //End of loop over nodes
863  } // end loop over boundaries
864 
865  // Complete the build of all elements so they are fully functional
866  // Remember that adaptation for triangle meshes involves a complete
867  // regneration of the mesh (rather than splitting as in tree-based
868  // meshes where such parameters can be passed down from the father
869  // element!)
870  unsigned n_element = Fluid_mesh_pt->nelement();
871  for(unsigned e=0;e<n_element;e++)
872  {
873  // Upcast from GeneralisedElement to the present element
874  ELEMENT* el_pt = dynamic_cast<ELEMENT*>(Fluid_mesh_pt->element_pt(e));
875 
876  // Set the Reynolds number
877  el_pt->re_pt() = &Problem_Parameter::Re;
878 
879  // Set the Womersley number (same as Re since St=1)
880  el_pt->re_st_pt() = &Problem_Parameter::ReSt;
881 
882  // Set the product of the Reynolds number and the inverse of the
883  // Froude number
884  el_pt->re_invfr_pt() = &Problem_Parameter::ReInvFr;
885 
886  // Set the direction of gravity
887  el_pt->g_pt() = &Problem_Parameter::G;
888 
889  // Set the constitutive law for pseudo-elastic mesh deformation
890  el_pt->constitutive_law_pt()=Problem_Parameter::Constitutive_law_pt;
891  }
892 
893  //For the elements in the upper region (region 0),
894  //set the viscosity ratio
895  n_element = Fluid_mesh_pt->nregion_element(1);
896  for(unsigned e=0;e<n_element;e++)
897  {
898  // Upcast from GeneralisedElement to the present element
899  ELEMENT* el_pt =
900  dynamic_cast<ELEMENT*>(Fluid_mesh_pt->region_element_pt(1,e));
901 
902  el_pt->viscosity_ratio_pt() = &Problem_Parameter::Viscosity_Ratio;
903 
904  el_pt->density_ratio_pt() = &Problem_Parameter::Density_Ratio;
905  }
906 
907 
908  // Re-apply boundary values on Dirichlet boundary conditions
909  // (Boundary conditions are ignored when the solution is transferred
910  // from the old to the new mesh by projection; this leads to a slight
911  // change in the boundary values (which are, of course, never changed,
912  // unlike the actual unknowns for which the projected values only
913  // serve as an initial guess)
914 
915  // Set velocity and history values of velocity on walls
916  nbound=this->Fluid_mesh_pt->nboundary();
917  for(unsigned ibound=0;ibound<nbound;++ibound)
918  {
919  // Loop over nodes on this boundary
920  unsigned num_nod=this->Fluid_mesh_pt->nboundary_node(ibound);
921  for (unsigned inod=0;inod<num_nod;inod++)
922  {
923  // Get node
924  Node* nod_pt=this->Fluid_mesh_pt->boundary_node_pt(ibound,inod);
925 
926  // Get number of previous (history) values
927  unsigned n_prev=nod_pt->time_stepper_pt()->nprev_values();
928 
929  //If we are on the upper or lower walls
930  // Velocity is and was zero at all previous times
931  if((ibound==1) || (ibound==3))
932  {
933  //Loop over time history values
934  for (unsigned t=0;t<=n_prev;t++)
935  {
936  nod_pt->set_value(t,0,0.0);
937  nod_pt->set_value(t,1,0.0);
938 
939  // Nodes have always been there...
940  nod_pt->x(t,0)=nod_pt->x(0,0);
941  nod_pt->x(t,1)=nod_pt->x(0,1);
942  }
943  }
944  //If we are on the side wals there is no horizontal velocity or
945  //change in horizontal position
946  if((ibound==0) || (ibound==2))
947  {
948  for (unsigned t=0;t<=n_prev;t++)
949  {
950  nod_pt->set_value(t,0,0.0);
951  // Nodes have always been there...
952  nod_pt->x(t,0)=nod_pt->x(0,0);
953  }
954  }
955  //But there is a change in vertical position!
956  }
957  } //End of loop over boundaries
958 
959  fix_pressure(0,0,0.0);
960 
961  }
962 
963 
964 
965 //==start_of_doc_solution=================================================
966 /// Doc the solution
967 //========================================================================
968 template<class ELEMENT>
969 void TwoLayerInterfaceProblem<ELEMENT>::doc_solution(const std::string& comment)
970 {
971 
972  // Output the time
973  cout << "Time is now " << time_pt()->time() << std::endl;
974 
975  double min_x_coordinate = 2.0;
976  unsigned min_boundary_node = 0;
977  //Find the left-most node on the boundary
978  unsigned n_bound = Fluid_mesh_pt->nboundary_node(4);
979  for(unsigned n=0;n<n_bound;++n)
980  {
981  Node* nod_pt = Fluid_mesh_pt->boundary_node_pt(4,n);
982  if(nod_pt->x(0) < min_x_coordinate)
983  {
984  min_x_coordinate = nod_pt->x(0);
985  min_boundary_node = n;
986  }
987  }
988 
989  // Document time and vertical position of left hand side of interface
990  // in trace file
992  << time_pt()->time() << " "
993  << Fluid_mesh_pt->boundary_node_pt(4,min_boundary_node)->x(1) << std::endl;
994 
995  ofstream some_file;
996  char filename[100];
997 
998  // Set number of plot points (in each coordinate direction)
999  const unsigned npts = 5;
1000 
1001  // Open solution output file
1002  sprintf(filename,"%s/soln%i.dat",
1003  Problem_Parameter::Doc_info.directory().c_str(),
1004  Problem_Parameter::Doc_info.number());
1005  some_file.open(filename);
1006 
1007  // Output solution to file
1008  Fluid_mesh_pt->output(some_file,npts);
1009 
1010  // Close solution output file
1011  some_file.close();
1012 
1013  // Open interface solution output file
1014  sprintf(filename,"%s/interface_soln%i.dat",
1015  Problem_Parameter::Doc_info.directory().c_str(),
1016  Problem_Parameter::Doc_info.number());
1017  some_file.open(filename);
1018 
1019  // Output solution to file
1020  Free_surface_mesh_pt->output(some_file,npts);
1021 
1022  // Close solution output file
1023  some_file.close();
1024 
1025  // Increment the doc_info number
1026  Problem_Parameter::Doc_info.number()++;
1027 }
1028 
1029 //========================================================================
1030 /// Compute error estimates and assign to elements for plotting
1031 //========================================================================
1032 template<class ELEMENT>
1034  double& min_err)
1035 {
1036  // Get error estimator
1037  ErrorEstimator* err_est_pt=Fluid_mesh_pt->spatial_error_estimator_pt();
1038 
1039  // Get/output error estimates
1040  unsigned nel=Fluid_mesh_pt->nelement();
1041  Vector<double> elemental_error(nel);
1042 
1043  // We need a dynamic cast, get_element_errors from the Fluid_mesh_pt
1044  // Dynamic cast is used because get_element_errors require a Mesh* ans
1045  // not a SolidMesh*
1046  Mesh* fluid_mesh_pt=dynamic_cast<Mesh*>(Fluid_mesh_pt);
1047  err_est_pt->get_element_errors(fluid_mesh_pt,
1048  elemental_error);
1049 
1050  // Set errors for post-processing and find extrema
1051  max_err=0.0;
1052  min_err=DBL_MAX;
1053  for (unsigned e=0;e<nel;e++)
1054  {
1055  dynamic_cast<MyCrouzeixRaviartElement*>(Fluid_mesh_pt->element_pt(e))->
1056  set_error(elemental_error[e]);
1057 
1058  max_err=std::max(max_err,elemental_error[e]);
1059  min_err=std::min(min_err,elemental_error[e]);
1060  }
1061 
1062 }
1063 
1064 
1065 //============================================================
1066 /// Driver code for moving block problem
1067 //============================================================
1068 int main(int argc, char **argv)
1069 {
1070 
1071  // Store command line arguments
1072  CommandLineArgs::setup(argc,argv);
1073 
1074  // Define possible command line arguments and parse the ones that
1075  // were actually specified
1076 
1077  // Validation?
1078  CommandLineArgs::specify_command_line_flag("--validation");
1079 
1080  // Parse command line
1081  CommandLineArgs::parse_and_assign();
1082 
1083  // Doc what has actually been specified on the command line
1084  CommandLineArgs::doc_specified_flags();
1085 
1086  // Create generalised Hookean constitutive equations
1088  new GeneralisedHookean(&Problem_Parameter::Nu);
1089 
1090  Problem_Parameter::G[0] = 0.0;
1091  Problem_Parameter::G[1] = -1.0;
1092 
1093  // Create problem in initial configuration
1094  TwoLayerInterfaceProblem<Hijacked<ProjectableCrouzeixRaviartElement<
1095  MyCrouzeixRaviartElement> > > problem;
1096 
1097  // Set value of epsilon
1098  const double epsilon = 0.1;
1099 
1100  // Set number of periods for cosine term
1101  const unsigned n_periods = 1;
1102 
1103  const double dt = 0.0025;
1104 
1105  double t_max = 0.6;
1106 
1107  // Deform the mesh/free surface
1108  problem.deform_interface(epsilon,n_periods);
1109 
1110  // Open trace file
1111  char filename[100];
1112  sprintf(filename,"%s/trace.dat",Problem_Parameter::Doc_info.directory().c_str());
1113  Problem_Parameter::Trace_file.open(filename);
1114 
1115  // Initialise trace file
1116  Problem_Parameter::Trace_file << "time, interface height" << std::endl;
1117 
1118  // Initialise timestep
1119  problem.initialise_dt(dt);
1120 
1121  // Set initial condition
1122  problem.set_initial_condition();
1123 
1124  // Maximum number of spatial adaptations per timestep
1125  unsigned max_adapt = 2;
1126 
1127  // Doc initial solution
1128  problem.doc_solution();
1129 
1130  // Determine number of timesteps
1131  const unsigned n_timestep = unsigned(t_max/dt);
1132 
1133  // Are we on the first timestep? At this point, yes!
1134  bool first_timestep = true;
1135 
1136  // Timestepping loop
1137  for(unsigned t=1;t<=n_timestep;t++)
1138  {
1139  // Output current timestep to screen
1140  cout << "\nTimestep " << t << " of " << n_timestep << std::endl;
1141 
1142  // Take one fixed timestep with spatial adaptivity
1143  problem.unsteady_newton_solve(dt,max_adapt,first_timestep);
1144 
1145  // No longer on first timestep, so set first_timestep flag to false
1146  first_timestep = false;
1147 
1148  // Reset maximum number of adaptations for all future timesteps
1149  max_adapt = 1;
1150 
1151  // Doc solution
1152  problem.doc_solution();
1153 
1154  } // End of timestepping loop
1155 
1156 
1157 } //End of main
Problem class to simulate viscous drop propagating along 2D channel.
ELEMENT * Hijacked_element_pt
Pointer to hijacked element.
void doc_solution(const std::string &comment="")
Doc the solution.
Mesh * Free_surface_mesh_pt
Pointers to mesh of free surface elements.
void delete_free_surface_elements()
Delete free surface elements.
void actions_after_adapt()
Actions after adapt: Rebuild the mesh of free surface elements.
void fix_pressure(const unsigned &e, const unsigned &pdof, const double &pvalue)
Fix pressure in element e at pressure dof pdof and set to pvalue.
void compute_error_estimate(double &max_err, double &min_err)
Compute the error estimates and assign to elements for plotting.
void actions_before_adapt()
Actions before adapt: Wipe the mesh of free surface elements.
void create_free_surface_elements()
Create free surface elements.
void complete_problem_setup()
Set boundary conditions and complete the build of all elements.
void actions_before_newton_solve()
Update the problem specs before solve.
void actions_after_newton_solve()
Update the after solve (empty)
void set_initial_condition()
Set the initial conditions.
Vector< TriangleMeshPolygon * > Drop_polygon_pt
Vector storing pointer to the drop polygons.
RefineableSolidTriangleMesh< ELEMENT > * Fluid_mesh_pt
Pointer to Fluid_mesh.
void deform_interface(const double &epsilon, const unsigned &n_periods)
Function to deform the interface.
bool Use_volume_constraint
Bool to indicate if volume constraint is applied (only for steady run)
TriangleMeshPolygon * Outer_boundary_polyline_pt
Triangle mesh polygon for outer boundary.
Data * Drop_pressure_data_pt
Pointer to a global drop pressure datum.
double Initial_value_for_drop_pressure
Backed up drop pressure between adaptations.
Overload CrouzeixRaviart element to modify output.
void set_error(const double &error)
Set error value for post-processing.
MyCrouzeixRaviartElement()
Constructor initialise error.
double Error
Storage for elemental error estimate – used for post-processing.
void output(std::ostream &outfile, const unsigned &nplot)
Overload output function.
std::string variable_identifier()
Return variable identifier.
//////////////////////////////////////////////////////// ////////////////////////////////////////////...
double ReSt
Womersley number (Reynolds x Strouhal)
ofstream Norm_file
File to document the norm of the solution (for validation purposes – triangle doesn't give fully repr...
double Length
Length of the channel.
Vector< double > G(2)
Direction of gravity.
double Density_Ratio
Ratio of density in upper fluid to density in lower fluid. Reynolds number etc. is based on density i...
ConstitutiveLaw * Constitutive_law_pt
Constitutive law used to determine the mesh deformation.
double ReInvFr
Product of Reynolds number and inverse of Froude number.
double Nu
Pseudo-solid Poisson ratio.
double Viscosity_Ratio
Ratio of viscosity in upper fluid to viscosity in lower fluid. Reynolds number etc....
int main(int argc, char **argv)
Driver code for moving block problem.