rayleigh_instability_soluble_surfactant.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 // Driver for axisymmetric single-layer fluid problem. Plateau-Rayleigh
27 // instability (unstable if H>2*pi*R -> forming drops)
28 // in the presence of a soluble surfactant. The problem is described in
29 // Numerical analaysis of the Rayleigh instability in capillary tubes:
30 // The influence of surfactant solubility by
31 // D. M. Campana & F. A. Saita, Phys. Fluids.
32 // vol 18, 022104 (2006)
33 // The parameters are taken from their Table 3
34 
35 // Generic oomph-lib header
36 #include "generic.h"
37 
38 // Axisymmetric Navier-Stokes headers
39 #include "axisym_navier_stokes.h"
40 
41 // Interface headers
42 #include "fluid_interface.h"
43 //Header for the pesudo-solid mesh updates
44 #include "solid.h"
45 #include "constitutive.h"
46 
47 // The standard rectangular mesh
48 #include "meshes/rectangular_quadmesh.h"
49 
50 // Equations for bulk surfactant transport
52 
53 //Use the oomph and std namespaces
54 using namespace oomph;
55 using namespace std;
56 
57 //==start_of_namespace===================================================
58 /// Namespace for physical parameters
59 /// The parameter values are chosen to be those used in Figures 7, 12
60 /// in Campana & Saita
61 //=======================================================================
63 {
64 
65  //Film thickness parameter
66  double Film_Thickness = 0.18;
67 
68  /// Reynolds number
69  double Re = 1.0;
70 
71  /// Womersley number
72  double ReSt = Re; // (St = 1)
73 
74  /// Product of Reynolds number and inverse of Froude number
75  double ReInvFr = 0.0; // (Fr = 0)
76 
77  /// Capillary number
78  double Ca = pow(Film_Thickness,3.0);
79 
80  /// External pressure
81  double P_ext = 0.0;
82 
83  /// Direction of gravity
84  Vector<double> G(3);
85 
86  /// Wavelength of the domain
87  double Alpha = 0.8537;
88 
89  /// Free surface cosine deformation parameter
90  double Epsilon = 1.0e-3;
91 
92  /// Surface Elasticity number
93  double Beta = 0.1;
94 
95  /// Alpha for absorption kinetics
96  double Alpha_absorption = 1.0;
97 
98  /// K parameter that
99  /// describes solubility number
100  double K = 0.01;
101 
102  // Bulk Peclet number
103  double Pe = 10000.0;
104  // Bulk Peclet Strouhal number
105  double PeSt = Pe;
106 
107  //We shall adopt the same scaling
108  //as Campana & Saita, so divide
109  //the bulk equation through by
110  //Peclet number to give 1 as the
111  //nominal peclet number
112  double Pe_reference_scale=1.0;
113 
114  //Chose our diffusion in the bulk equations
115  //to be 1/Pe, after dividing through by Pe
116  double Diff = 1.0/Pe;
117 
118  /// Surface Peclet number
119  double Peclet_S = 10000.0;
120 
121  /// Sufrace Peclet number multiplied by Strouhal number
122  double Peclet_St_S = 1.0;
123 
124  /// Pvd file -- a wrapper for all the different
125  /// vtu output files plus information about continuous time
126  /// to facilitate animations in paraview
127  ofstream Pvd_file;
128 
129  /// Pseudo-solid Poisson ratio
130  double Nu = 0.1;
131 
132 } // End of namespace
133 
134 
135 
136 //================================================================
137 /// Interface class to handle the mass transport between bulk
138 /// and surface as well as the surfactant transport along the
139 //interface
140 //================================================================
141 template<class ELEMENT>
144 {
145 private:
146  /// Pointer to adsorption number
147  double *Alpha_pt;
148  /// Pointer to solubility number
149  double *K_pt;
150 
151  /// Storage for the index at
152  /// which the bulk concentration is stored
153  unsigned C_bulk_index;
154 
155 protected:
156 
157  /// Get the bulk surfactant concentration
158  double interpolated_C_bulk(const Vector<double> &s)
159  {
160  //Find number of nodes
161  unsigned n_node = this->nnode();
162 
163  //Get the nodal index at which the unknown is stored
164  const unsigned c_index = C_bulk_index;
165 
166  //Local shape function
167  Shape psi(n_node);
168 
169  //Find values of shape function
170  this->shape(s,psi);
171 
172  //Initialise value of C
173  double C = 0.0;
174 
175  //Loop over the local nodes and sum
176  for(unsigned l=0;l<n_node;l++)
177  {
178  C += this->nodal_value(l,c_index)*psi(l);
179  }
180 
181  return(C);
182  }
183 
184  /// The time derivative of the bulk surface concentration
185  double dc_bulk_dt_surface(const unsigned &l) const
186  {
187  // Get the data's timestepper
188  TimeStepper* time_stepper_pt= this->node_pt(l)->time_stepper_pt();
189 
190  //Initialise dudt
191  double dcdt=0.0;
192  //Loop over the timesteps, if there is a non Steady timestepper
193  if (time_stepper_pt->type()!="Steady")
194  {
195  //Find the index at which the variable is stored
196  const unsigned c_index = C_bulk_index;
197 
198  // Number of timsteps (past & present)
199  const unsigned n_time = time_stepper_pt->ntstorage();
200 
201  for(unsigned t=0;t<n_time;t++)
202  {
203  dcdt += time_stepper_pt->weight(1,t)*this->nodal_value(t,l,c_index);
204  }
205  }
206  return dcdt;
207  }
208 
209  /// Overload the Helper function to calculate the residuals and
210  /// jacobian entries. This particular function ensures that the
211  /// additional entries are calculated inside the integration loop
213  Vector<double> &residuals, DenseMatrix<double> &jacobian,
214  const unsigned &flag,const Shape &psif, const DShape &dpsifds,
215  const DShape &dpsifdS, const DShape &dpsifdS_div,
216  const Vector<double> &s,
217  const Vector<double> &interpolated_x, const Vector<double> &interpolated_n,
218  const double &W,const double &J)
219  {
220  //Call the standard transport equations
223  residuals, jacobian, flag,psif,dpsifds,dpsifdS,dpsifdS_div,
224  s, interpolated_x, interpolated_n, W, J);
225 
226  //Add the additional mass transfer terms
227  const double k = this->k();
228  const double alpha = this->alpha();
229  //Find out how many nodes there are
230  unsigned n_node = this->nnode();
231 
232  //Storage for the local equation numbers and unknowns
233  int local_eqn = 0, local_unknown = 0;
234 
235  //Surface advection-diffusion equation
236 
237  //Find the index at which the bulk concentration is stored
238  unsigned c_bulk_index = this->C_bulk_index;
239 
240  //Now calculate the bulk and surface concentrations at this point
241  //Assuming the same shape functions are used (which they are)
242  double interpolated_C = 0.0;
243  double interpolated_C_bulk = 0.0;
244 
245  //Loop over the shape functions
246  for(unsigned l=0;l<n_node;l++)
247  {
248  const double psi = psif(l);
249  const double C_ = this->nodal_value(l,this->C_index[l]);
250 
251  interpolated_C += C_*psi;
252  interpolated_C_bulk += this->nodal_value(l,c_bulk_index)*psi;
253  }
254 
255  //Pre compute the flux between the surface and bulk
256  double flux = alpha*(interpolated_C_bulk - interpolated_C);
257 
258  //Now we add the flux term to the appropriate residuals
259  for(unsigned l=0;l<n_node;l++)
260  {
261  //BULK FLUX
262 
263  //Read out the appropriate local equation
264  local_eqn = this->nodal_local_eqn(l,c_bulk_index);
265 
266  //If not a boundary condition
267  if(local_eqn >= 0)
268  {
269  //Add flux to the bulk
270  residuals[local_eqn] -= k*flux*psif(l)*W*J;
271 
272  if(flag)
273  {
274  for(unsigned l2=0;l2<n_node;l2++)
275  {
276  local_unknown = this->nodal_local_eqn(l2,this->C_index[l2]);
277  if(local_unknown >= 0)
278  {
279  jacobian(local_eqn,local_unknown) += k*alpha*psif(l2)*psif(l)*W*J;
280  }
281 
282  local_unknown = this->nodal_local_eqn(l2,c_bulk_index);
283  if(local_unknown >= 0)
284  {
285  jacobian(local_eqn,local_unknown) -= k*alpha*psif(l2)*psif(l)*W*J;
286  }
287  }
288  }
289  } //End of contribution to bulk
290 
291 
292  //SURFACE FLUX
293  //Read out the appropriate local equation
294  local_eqn = this->nodal_local_eqn(l,this->C_index[l]);
295 
296  //If not a boundary condition
297  if(local_eqn >= 0)
298  {
299  //Add flux to the surface
300  residuals[local_eqn] -= flux*psif(l)*W*J;
301 
302  if(flag)
303  {
304  for(unsigned l2=0;l2<n_node;l2++)
305  {
306  local_unknown = this->nodal_local_eqn(l2,this->C_index[l2]);
307  if(local_unknown >= 0)
308  {
309  jacobian(local_eqn,local_unknown) += alpha*psif(l2)*psif(l)*W*J;
310  }
311 
312  local_unknown = this->nodal_local_eqn(l2,c_bulk_index);
313  if(local_unknown >= 0)
314  {
315  jacobian(local_eqn,local_unknown) -= alpha*psif(l2)*psif(l)*W*J;
316  }
317  }
318  } //End of contribution to surface
319 
320  }
321  }
322 
323  }
324 
325 public:
326  /// Constructor that passes the bulk element and face index down
327  /// to the underlying
329  FiniteElement* const &element_pt, const int &face_index) :
331  (element_pt,face_index)
332  {
333  //Initialise the values
334  Alpha_pt = &this->Default_Physical_Constant_Value;
335  K_pt = &this->Default_Physical_Constant_Value;
336 
337  //Hack because I know the bulk index in this case
338  C_bulk_index = 3;
339  }
340 
341 
342  /// Return the adsorption number
343  double alpha() {return *Alpha_pt;}
344 
345  /// Return the solubility nubmer
346  double k() {return *K_pt;}
347 
348  /// Access function for pointer to adsorption number
349  double* &alpha_pt() {return Alpha_pt;}
350 
351  /// Access function for pointer to solubility number
352  double* &k_pt() {return K_pt;}
353 
354  /// Overload the output function
355  void output(std::ostream &outfile, const unsigned &n_plot)
356  {
357  outfile.precision(16);
358 
359  //Set output Vector
360  Vector<double> s(1);
361 
362  //Tecplot header info
363  outfile << "ZONE I=" << n_plot << std::endl;
364 
365  const unsigned n_node = this->nnode();
366  const unsigned dim = this->dim();
367 
368  Shape psi(n_node);
369  DShape dpsi(n_node,dim);
370 
371  //Loop over plot points
372  for(unsigned l=0;l<n_plot;l++)
373  {
374  s[0] = -1.0 + l*2.0/(n_plot-1);
375 
376  this->dshape_local(s,psi,dpsi);
377  Vector<double> interpolated_tangent(2,0.0);
378  for(unsigned l=0;l<n_node;l++)
379  {
380  const double dpsi_ = dpsi(l,0);
381  for(unsigned i=0;i<2;i++)
382  {
383  interpolated_tangent[i] += this->nodal_position(l,i)*dpsi_;
384  }
385  }
386 
387  //Tangent
388  double t_vel = (this->interpolated_u(s,0)*interpolated_tangent[0] + this->interpolated_u(s,1)*interpolated_tangent[1])/
389  sqrt(interpolated_tangent[0]*interpolated_tangent[0] + interpolated_tangent[1]*interpolated_tangent[1]);
390 
391 
392  //Output the x,y,u,v
393  for(unsigned i=0;i<2;i++) outfile << this->interpolated_x(s,i) << " ";
394  for(unsigned i=0;i<2;i++) outfile << this->interpolated_u(s,i) << " ";
395  //Output a dummy pressure
396  outfile << 0.0 << " ";
397  //Output the concentrations
398  outfile << this->interpolated_C(s) << " ";
399  outfile << interpolated_C_bulk(s) << " ";
400  //Output the interfacial tension
401  outfile << this->sigma(s) << " " << t_vel << std::endl;
402  }
403  outfile << std::endl;
404  }
405 
406 };
407 
408 //==start_of_problem_class===============================================
409 /// Single axisymmetric fluid interface problem including the
410 /// transport of an soluble surfactant.
411 //=======================================================================
412 template<class ELEMENT, class TIMESTEPPER>
413 class InterfaceProblem : public Problem
414 {
415 
416 public:
417 
418  /// Constructor: Pass the number of elements in radial and axial directions
419  /// and the length of the domain in the z direction)
420  InterfaceProblem(const unsigned &n_r, const unsigned &n_z,
421  const double &l_z);
422 
423  /// Destructor (empty)
425 
426  /// Set initial conditions: Set all nodal velocities to zero and
427  /// initialise the previous velocities to correspond to an impulsive
428  /// start
430  {
431  // Determine number of nodes in mesh
432  const unsigned n_node = Bulk_mesh_pt->nnode();
433 
434  // Loop over all nodes in mesh
435  for(unsigned n=0;n<n_node;n++)
436  {
437  // Loop over the three velocity components
438  for(unsigned i=0;i<3;i++)
439  {
440  // Set velocity component i of node n to zero
441  Bulk_mesh_pt->node_pt(n)->set_value(i,0.0);
442  }
443  //Set the bulk concentration to be 1
444  Bulk_mesh_pt->node_pt(n)->set_value(3,1.0);
445  }
446 
447  // Initialise the previous velocity values for timestepping
448  // corresponding to an impulsive start
449  assign_initial_values_impulsive();
450 
451  } // End of set_initial_condition
452 
453 
454  /// The global temporal error norm, based on the movement of the nodes
455  /// in the radial direction only (because that's the only direction
456  /// in which they move!)
458  {
459  //Temp
460  double global_error = 0.0;
461 
462  //Find out how many nodes there are in the problem
463  const unsigned n_node = Bulk_mesh_pt->nnode();
464 
465  //Loop over the nodes and calculate the errors in the positions
466  for(unsigned n=0;n<n_node;n++)
467  {
468  //Set the dimensions to be restricted to the radial direction only
469  const unsigned n_dim = 1;
470  //Set the position error to zero
471  double node_position_error = 0.0;
472  //Loop over the dimensions
473  for(unsigned i=0;i<n_dim;i++)
474  {
475  //Get position error
476  double error =
477  Bulk_mesh_pt->node_pt(n)->position_time_stepper_pt()->
478  temporal_error_in_position(Bulk_mesh_pt->node_pt(n),i);
479 
480  //Add the square of the individual error to the position error
481  node_position_error += error*error;
482  }
483 
484  //Divide the position error by the number of dimensions
485  node_position_error /= n_dim;
486  //Now add to the global error
487  global_error += node_position_error;
488  }
489 
490  //Now the global error must be divided by the number of nodes
491  global_error /= n_node;
492 
493  //Return the square root of the errr
494  return sqrt(global_error);
495  }
496 
497 
498  /// Access function for the specific mesh
499  ElasticRectangularQuadMesh<ELEMENT>* Bulk_mesh_pt;
500 
501  /// Mesh for the free surface (interface) elements
502  Mesh* Interface_mesh_pt;
503 
504  /// Pointer to the constitutive law
505  ConstitutiveLaw* Constitutive_law_pt;
506 
507  /// Node for documentatin
508  Node* Document_node_pt;
509 
510  /// Doc the solution
511  void doc_solution(DocInfo &doc_info);
512 
513  /// Do unsteady run up to maximum time t_max with given timestep dt
514  void unsteady_run(const double &t_max, const double &dt);
515 
516  /// Compute the total mass of the insoluble surfactant
518  {
519  //Initialise to zero
520  double mass = 0.0;
521 
522  // Determine number of 1D interface elements in mesh
523  const unsigned n_interface_element = Interface_mesh_pt->nelement();
524 
525  // Loop over the interface elements
526  for(unsigned e=0;e<n_interface_element;e++)
527  {
528  // Upcast from GeneralisedElement to the present element
531  (Interface_mesh_pt->element_pt(e));
532  //Add contribution from each element
533  mass += el_pt->integrate_c();
534  }
535  return mass;
536  } // End of compute_total_mass
537 
538 
539 private:
540 
541  /// Deform the mesh/free surface to a prescribed function
542  void deform_free_surface(const double &epsilon)
543  {
544  //Loop over all nodes in the mesh
545  const unsigned n_node = Bulk_mesh_pt->nnode();
546  for(unsigned n=0;n<n_node;n++)
547  {
548  //Find out the r value
549  double r = Bulk_mesh_pt->node_pt(n)->x(0);
550  // Determine z coordinate of spine
551  double z_value = Bulk_mesh_pt->node_pt(n)->x(1);
552 
553  //Set the thickess of the filme
554  double thickness =
556  + epsilon*cos(Global_Physical_Variables::Alpha*z_value));
557 
558  //Now scale the position accordingly
559  Bulk_mesh_pt->node_pt(n)->x(0) = 1.0 - (1.0-r)*thickness;
560  }
561 
562  //Reset the lagrangian coordinates
563  Bulk_mesh_pt->set_lagrangian_nodal_coordinates();
564  } // End of deform_free_surface
565 
566  /// Trace file
567 ofstream Trace_file;
568 
569 }; // End of problem class
570 
571 
572 
573 //==start_of_constructor==================================================
574 /// Constructor for single fluid interface problem
575 //========================================================================
576 template<class ELEMENT, class TIMESTEPPER>
578 InterfaceProblem(const unsigned &n_r,
579  const unsigned &n_z,
580  const double &l_z)
581 
582 {
583  // Allocate the timestepper (this constructs the time object as well)
584  add_time_stepper_pt(new TIMESTEPPER(true));
585 
586  // Build and assign mesh (the "false" boolean flag tells the mesh
587  // constructor that the domain is not periodic in r)
588  Bulk_mesh_pt = new ElasticRectangularQuadMesh<ELEMENT>(n_r,n_z,1.0,l_z,time_stepper_pt());
589 
590  //Create "surface mesh" that will only contain the interface elements
591  Interface_mesh_pt = new Mesh;
592  {
593  // How many bulk elements are adjacent to boundary b?
594  // Boundary 1 is the outer boundary
595  // The boundaries are labelled
596  // 2
597  // 3 1
598  // 0
599 
600  unsigned n_element = Bulk_mesh_pt->nboundary_element(3);
601 
602  // Loop over the bulk elements adjacent to boundary b?
603  for(unsigned e=0;e<n_element;e++)
604  {
605  // Get pointer to the bulk element that is adjacent to boundary b
606  ELEMENT* bulk_elem_pt = dynamic_cast<ELEMENT*>(
607  Bulk_mesh_pt->boundary_element_pt(3,e));
608 
609  // Find the index of the face of element e along boundary b
610  int face_index = Bulk_mesh_pt->face_index_at_boundary(3,e);
611 
612  // Build the corresponding free surface element
615 
616  //Add the prescribed-flux element to the surface mesh
617  Interface_mesh_pt->add_element_pt(interface_element_pt);
618 
619  } //end of loop over bulk elements adjacent to boundary b
620  }
621 
622 
623  //Use the first node as our document
624  Document_node_pt = Bulk_mesh_pt->node_pt(0);
625 
626  // Add the two sub meshes to the problem
627  add_sub_mesh(Bulk_mesh_pt);
628  add_sub_mesh(Interface_mesh_pt);
629 
630  // Combine all submeshes into a single Mesh
631  build_global_mesh();
632 
633 
634  // --------------------------------------------
635  // Set the boundary conditions for this problem
636  // --------------------------------------------
637 
638  //Pin all azimuthal velocities
639  //and all vertical positions so the
640  //nodes can only move horizontally
641  {
642  unsigned n_node = this->Bulk_mesh_pt->nnode();
643  for(unsigned n=0;n<n_node;++n)
644  {
645  this->Bulk_mesh_pt->node_pt(n)->pin(2);
646  this->Bulk_mesh_pt->node_pt(n)->pin_position(1);
647  }
648  }
649 
650  // All nodes are free by default -- just pin the ones that have
651  // Dirichlet conditions here
652  unsigned ibound = 3;
653  unsigned n_node = Bulk_mesh_pt->nboundary_node(ibound);
654  for(unsigned n=0;n<n_node;n++)
655  {
656  Bulk_mesh_pt->boundary_node_pt(ibound,n)->set_value(4,1.0);
657  }
658 
659  // Determine number of mesh boundaries
660  const unsigned n_boundary = Bulk_mesh_pt->nboundary();
661 
662  // Loop over mesh boundaries
663  for(unsigned b=0;b<n_boundary;b++)
664  {
665  // Determine number of nodes on boundary b
666  const unsigned n_node = Bulk_mesh_pt->nboundary_node(b);
667 
668  // Loop over nodes on boundary b
669  for(unsigned n=0;n<n_node;n++)
670  {
671  // Pin azimuthal velocity on bounds
672  Bulk_mesh_pt->boundary_node_pt(b,n)->pin(2);
673 
674  // Pin velocity on the outer wall
675  if(b==1)
676  {
677  Bulk_mesh_pt->boundary_node_pt(b,n)->pin(0);
678  Bulk_mesh_pt->boundary_node_pt(b,n)->pin(1);
679  Bulk_mesh_pt->boundary_node_pt(b,n)->pin_position(0);
680  }
681  // Pin axial velocity on bottom and top walls (no penetration)
682  if(b==0 || b==2)
683  {
684  Bulk_mesh_pt->boundary_node_pt(b,n)->pin(1);
685  }
686  } // End of loop over nodes on boundary b
687  } // End of loop over mesh boundaries
688 
689  // ----------------------------------------------------------------
690  // Complete the problem setup to make the elements fully functional
691  // ----------------------------------------------------------------
692  Constitutive_law_pt = new GeneralisedHookean(&Global_Physical_Variables::Nu);
693 
694  // Determine number of bulk elements in mesh
695  const unsigned n_bulk = Bulk_mesh_pt->nelement();
696 
697  // Loop over the bulk elements
698  for(unsigned e=0;e<n_bulk;e++)
699  {
700  // Upcast from GeneralisedElement to the present element
701  ELEMENT *el_pt = dynamic_cast<ELEMENT*>(Bulk_mesh_pt->element_pt(e));
702 
703  // Set the Reynolds number
704  el_pt->re_pt() = &Global_Physical_Variables::Re;
705 
706  // Set the Womersley number
707  el_pt->re_st_pt() = &Global_Physical_Variables::ReSt;
708 
709  // Set the product of the Reynolds number and the inverse of the
710  // Froude number
711  el_pt->re_invfr_pt() = &Global_Physical_Variables::ReInvFr;
712 
713  // Set the direction of gravity
714  el_pt->g_pt() = &Global_Physical_Variables::G;
715 
716  // Set the constitutive law
717  el_pt->constitutive_law_pt() = Constitutive_law_pt;
718 
719  //Set the peclet numbers
721 
723  //Set the diffusion scale
724  el_pt->d_pt() = &Global_Physical_Variables::Diff;
725  } // End of loop over bulk elements
726 
727  // Create a Data object whose single value stores the external pressure
728  Data* external_pressure_data_pt = new Data(1);
729 
730  // Pin and set the external pressure to some arbitrary value
731  double p_ext = Global_Physical_Variables::P_ext;
732 
733  external_pressure_data_pt->pin(0);
734  external_pressure_data_pt->set_value(0,p_ext);
735 
736  // Determine number of 1D interface elements in mesh
737  const unsigned n_interface_element = Interface_mesh_pt->nelement();
738 
739  // Loop over the interface elements
740  for(unsigned e=0;e<n_interface_element;e++)
741  {
742  // Upcast from GeneralisedElement to the present element
745  (Interface_mesh_pt->element_pt(e));
746 
747  // Set the Capillary number
749 
750  // Pass the Data item that contains the single external pressure value
751  el_pt->set_external_pressure_data(external_pressure_data_pt);
752 
753  // Set the surface elasticity number
755 
756  // Set the sorption
758 
759  // Set the K parameter
761 
762  // Set the surface peclect number
764 
765  // Set the surface peclect number multiplied by strouhal number
767 
768  } // End of loop over interface elements
769 
770  // Setup equation numbering scheme
771  cout << "Number of equations: " << assign_eqn_numbers() << std::endl;
772 
773 } // End of constructor
774 
775 
776 
777 //==start_of_doc_solution=================================================
778 /// Doc the solution
779 //========================================================================
780 template<class ELEMENT, class TIMESTEPPER>
782 doc_solution(DocInfo &doc_info)
783 {
784 
785  // Output the time
786  double t= time_pt()->time();
787  cout << "Time is now " << t << std::endl;
788 
789  // Document in trace file
790  Trace_file << time_pt()->time() << " "
791  << Document_node_pt->x(0)
792  << " " << this->compute_total_mass() << std::endl;
793 
794  ofstream some_file;
795  char filename[100];
796 
797  // Set number of plot points (in each coordinate direction)
798  const unsigned npts = 5;
799 
800  // Open solution output file
801 // sprintf(filename,"%s/soln%i.dat",doc_info.directory().c_str(),
802 // doc_info.number());
803 // some_file.open(filename);
804 
805  // Output solution to file
806 // Bulk_mesh_pt->output(some_file,npts);
807 // some_file.close();
808  //Put interface in separate file
809  sprintf(filename,"%s/int%i.dat",doc_info.directory().c_str(),
810  doc_info.number());
811  some_file.open(filename);
812  Interface_mesh_pt->output(some_file,npts);
813  some_file.close();
814 
815  // Write file as a tecplot text object...
816 // some_file << "TEXT X=2.5,Y=93.6,F=HELV,HU=POINT,C=BLUE,H=26,T=\"time = "
817  // << time_pt()->time() << "\"";
818  // ...and draw a horizontal line whose length is proportional
819  // to the elapsed time
820  //some_file << "GEOMETRY X=2.5,Y=98,T=LINE,C=BLUE,LT=0.4" << std::endl;
821  //some_file << "1" << std::endl;
822  //some_file << "2" << std::endl;
823  //some_file << " 0 0" << std::endl;
824  //some_file << time_pt()->time()*20.0 << " 0" << std::endl;
825 
826  // Close solution output file
827 // some_file.close();
828 
829  // Output solution to file in paraview format
830  //sprintf(filename,"%s/soln%i.vtu",doc_info.directory().c_str(),
831  // doc_info.number());
832  //some_file.open(filename);
833  //Bulk_mesh_pt->output_paraview(some_file,npts);
834  //some_file.close();
835 
836  // Write pvd information
837  //string file_name="soln"+StringConversion::to_string(doc_info.number())
838  // +".vtu";
839  //ParaviewHelper::write_pvd_information(Global_Physical_Variables::Pvd_file,
840  // file_name,t);
841 
842 } // End of doc_solution
843 
844 
845 
846 //==start_of_unsteady_run=================================================
847 /// Perform run up to specified time t_max with given timestep dt
848 //========================================================================
849 template<class ELEMENT, class TIMESTEPPER>
851 unsteady_run(const double &t_max, const double &dt)
852 {
853 
854  // Set value of epsilon
855  double epsilon = Global_Physical_Variables::Epsilon;
856 
857  // Deform the mesh/free surface
858  deform_free_surface(epsilon);
859 
860  // Initialise DocInfo object
861  DocInfo doc_info;
862 
863  // Set output directory
864  doc_info.set_directory("RESLT");
865 
866  // Initialise counter for solutions
867  doc_info.number()=0;
868 
869  // Open trace file
870  char filename[100];
871  sprintf(filename,"%s/trace.dat",doc_info.directory().c_str());
872  Trace_file.open(filename);
873 
874  // Initialise trace file
875  Trace_file << "time" << ", "
876  << "edge spine height" << ", "
877  << "mass " << ", " << std::endl;
878 
879  // Initialise timestep
880  initialise_dt(dt);
881 
882  // Set initial conditions
883  set_initial_condition();
884 
885  // Determine number of timesteps
886  const unsigned n_timestep = unsigned(t_max/dt);
887 
888  // Open pvd file -- a wrapper for all the different
889  // vtu output files plus information about continuous time
890  // to facilitate animations in paraview
891  //sprintf(filename,"%s/soln.pvd",doc_info.directory().c_str());
892  //Global_Physical_Variables::Pvd_file.open(filename);
893  //ParaviewHelper::write_pvd_header(Global_Physical_Variables::Pvd_file);
894 
895  // Doc initial solution
896  doc_solution(doc_info);
897 
898  // Increment counter for solutions
899  doc_info.number()++;
900 
901  //double dt_desired = dt;
902 
903  // Timestepping loop
904  for(unsigned t=1;t<=n_timestep;t++)
905  {
906  // Output current timestep to screen
907  cout << "\nTimestep " << t << " of " << n_timestep << std::endl;
908 
909  // Take one fixed timestep
910  unsteady_newton_solve(dt);
911 
912  //double dt_actual =
913  // adaptive_unsteady_newton_solve(dt_desired,1.0e-6);
914  //dt_desired = dt_actual;
915 
916 
917  // Doc solution
918  doc_solution(doc_info);
919 
920  // Increment counter for solutions
921  doc_info.number()++;
922 
923  //Reset the lagrangian coordinates
924  Bulk_mesh_pt->set_lagrangian_nodal_coordinates();
925 
926  } // End of timestepping loop
927 
928  // write footer and close pvd file
929  //ParaviewHelper::write_pvd_footer(Global_Physical_Variables::Pvd_file);
930  //Global_Physical_Variables::Pvd_file.close();
931 
932 } // End of unsteady_run
933 
934 
935 /// ////////////////////////////////////////////////////////////////////
936 /// ////////////////////////////////////////////////////////////////////
937 /// ////////////////////////////////////////////////////////////////////
938 
939 
940 //==start_of_main=========================================================
941 /// Driver code for single fluid axisymmetric horizontal interface problem
942 //========================================================================
943 int main(int argc, char* argv[])
944 {
945 
946  // Store command line arguments
947  CommandLineArgs::setup(argc,argv);
948 
949  /// Maximum time
950  double t_max = 1000.0;
951 
952  /// Duration of timestep
953  const double dt = 0.1;
954 
955  // If we are doing validation run, use smaller number of timesteps
956  if(CommandLineArgs::Argc>1)
957  {
958  t_max = 0.5;
959  }
960 
961  // Number of elements in radial (r) direction
962  const unsigned n_r = 10;
963 
964  // Number of elements in axial (z) direction
965  const unsigned n_z = 80;
966 
967  // Height of domain
968  const double l_z = MathematicalConstants::Pi/Global_Physical_Variables::Alpha;
969 
970  // Set direction of gravity (vertically downwards)
974 
975  // Set up the spine test problem with AxisymmetricQCrouzeixRaviartElements,
976  // using the BDF<2> timestepper
978  problem(n_r,n_z,l_z);
979 
980  // Run the unsteady simulation
981  problem.unsteady_run(t_max,dt);
982 } // End of main
983 
Interface class to handle the mass transport between bulk and surface as well as the surfactant trans...
void add_additional_residual_contributions_interface(Vector< double > &residuals, DenseMatrix< double > &jacobian, const unsigned &flag, const Shape &psif, const DShape &dpsifds, const DShape &dpsifdS, const DShape &dpsifdS_div, const Vector< double > &s, const Vector< double > &interpolated_x, const Vector< double > &interpolated_n, const double &W, const double &J)
Overload the Helper function to calculate the residuals and jacobian entries. This particular functio...
unsigned C_bulk_index
Storage for the index at which the bulk concentration is stored.
ElasticAxisymmetricSolubleSurfactantTransportInterfaceElement(FiniteElement *const &element_pt, const int &face_index)
Constructor that passes the bulk element and face index down to the underlying.
void output(std::ostream &outfile, const unsigned &n_plot)
Overload the output function.
double *& k_pt()
Access function for pointer to solubility number.
double *& alpha_pt()
Access function for pointer to adsorption number.
double dc_bulk_dt_surface(const unsigned &l) const
The time derivative of the bulk surface concentration.
double interpolated_C_bulk(const Vector< double > &s)
Get the bulk surfactant concentration.
Single fluid interface problem including transport of an insoluble surfactant.
void set_initial_condition()
Set initial conditions: Set all nodal velocities to zero and initialise the previous velocities to co...
void unsteady_run(const unsigned &nstep)
Run an unsteady simulation with specified number of steps.
InterfaceProblem(const unsigned &n_r, const unsigned &n_y, const unsigned &n_theta, const double &r_min, const double &r_max, const double &l_y, const double &theta_max)
Constructor: Pass number of elements in x and y directions. Also lengths of the domain in x- and y-di...
void deform_free_surface(const double &epsilon)
Deform the mesh/free surface to a prescribed function.
void doc_solution(DocInfo &doc_info)
Doc the solution.
ConstitutiveLaw * Constitutive_law_pt
Pointer to the constitutive law.
double global_temporal_error_norm()
The global temporal error norm, based on the movement of the nodes in the radial direction only (beca...
double compute_total_mass()
Compute the total mass of the insoluble surfactant.
void unsteady_run(const double &t_max, const double &dt)
Do unsteady run up to maximum time t_max with given timestep dt.
InterfaceProblem(const unsigned &n_r, const unsigned &n_z, const double &l_z)
Constructor: Pass the number of elements in radial and axial directions and the length of the domain ...
ElasticRectangularQuadMesh< ELEMENT > * Bulk_mesh_pt
Access function for the specific mesh.
void add_additional_residual_contributions_interface(Vector< double > &residuals, DenseMatrix< double > &jacobian, const unsigned &flag, const Shape &psif, const DShape &dpsifds, const DShape &dpsifdS, const DShape &dpsifdS_div, const Vector< double > &s, const Vector< double > &interpolated_x, const Vector< double > &interpolated_n, const double &W, const double &J)
Helper function to calculate the additional contributions to be added at each integration point....
double *& ca_pt()
Pointer to the Capillary number.
void set_external_pressure_data(Data *external_pressure_data_pt)
Set the Data that contains the single pressure value that specifies the "external pressure" for the i...
double *& peclet_s_pt()
Access function for pointer to the surface Peclet number.
double *& beta_pt()
Access function for pointer to the Elasticity number.
double integrate_c()
Compute the concentration intergated over the surface area.
double *& peclet_strouhal_s_pt()
Access function for pointer to the surface Peclet x Strouhal number.
Namepspace for global parameters, chosen from Campana et al. as in the axisymmetric problem.
double Alpha_absorption
Alpha for absorption kinetics.
double ReSt
Womersley = Reynolds times Strouhal.
double Epsilon
Free surface cosine deformation parameter.
double Beta
Surface Elasticity number (weak case)
ofstream Pvd_file
Pvd file – a wrapper for all the different vtu output files plus information about continuous time to...
double K
K parameter that describes solubility number.
double Alpha
Wavelength of the domain.
double ReInvFr
Product of Reynolds and Froude number.
double Peclet_St_S
Sufrace Peclet number multiplied by Strouhal number.
Vector< double > G(3)
Direction of gravity.
int main(int argc, char *argv[])
//////////////////////////////////////////////////////////////////// ////////////////////////////////...