inclined_plane.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-2025 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//A demo driver that solves the classic fluid flow problem of flow
27//of a fluid film along an inclined plane. Stability analysis performed by
28//Yih (1963), Benjamin (1957), ... and Blyth & Pozrikidis (2004).
29
30//This is an example of the subtleties involved in even a seemingly simple
31//free surface problem.
32
33//Standard C++ library includes
34#include <iostream>
35#include <fstream>
36#include <string>
37#include <cmath>
38
39//Finite-Element library routines
40#include "generic.h"
41#include "navier_stokes.h"
42#include "solid.h"
43#include "fluid_interface.h"
44#include "meshes/simple_rectangular_quadmesh.h"
45
46using namespace std;
47
48using namespace oomph;
49
50//The global physical variables
52{
53 /// Reynolds number, based on the average velocity within the fluid film
54 double Re=0.0;
55
56 /// The product of Reynolds number and inverse Froude number
57 /// is set to two in this problem, which gives the free surface velocity
58 /// to be sin(alpha). [Set to three in order to get the same scale as
59 /// used by Yih, Benjamin, etc]
60 double ReInvFr=2.0;
61
62 /// Angle of incline of the slope (45 degrees)
63 double Alpha = 1.0*atan(1.0);
64
65 /// The Vector direction of gravity, set in main()
66 Vector<double> G(2,0.0);
67
68 /// The Capillary number
69 double Ca= 1.0;
70
71 /// Set the wavenumber
72 double K = 0.1;
73
74 /// Set the number of waves desired in the domain
75 double N_wave = 3;
76
77 /// The length of the domain to fit the desired number of waves
78 double Length = 2*N_wave*4.0*atan(1.0)/K;
79
80 /// Direction of the wall normal vector (at the inlet)
81 Vector<double> Wall_normal;
82
83 /// Function that specifies the wall unit normal at the inlet
84 void wall_unit_normal_inlet_fct(const Vector<double> &x,
85 Vector<double> &normal)
86 {
87 normal=Wall_normal;
88 }
89
90 /// Function that specified the wall unit normal at the outlet
91 void wall_unit_normal_outlet_fct(const Vector<double> &x,
92 Vector<double> &normal)
93 {
94 //Set the normal
95 normal = Wall_normal;
96 //and flip the sign
97 unsigned n_dim = normal.size();
98 for(unsigned i=0;i<n_dim;++i) {normal[i] *= -1.0;}
99 }
100
101 /// The contact angle that is imposed at the inlet (pi)
102 double Inlet_Angle = 2.0*atan(1.0);
103
104
105 /// Function that prescribes the hydrostatic pressure field at the outlet
106 void hydrostatic_pressure_outlet(const double& time, const Vector<double> &x,
107 const Vector<double> &n,
108 Vector<double> &traction)
109 {
110 traction[0] = ReInvFr*G[1]*(1.0 - x[1]);
111 traction[1] = 0.0;
112 }
113
114 /// Function that prescribes hydrostatic pressure field at the inlet
115 void hydrostatic_pressure_inlet(const double& time, const Vector<double> &x,
116 const Vector<double> &n,
117 Vector<double> &traction)
118 {
119 traction[0] = -ReInvFr*G[1]*(1.0 - x[1]);
120 traction[1] = 0.0;
121 }
122 //end of traction functions
123
124 /// Constitutive law used to determine the mesh deformation
125 ConstitutiveLaw *Constitutive_law_pt;
126
127 /// Pseudo-solid Poisson ratio
128 double Nu=0.1;
129
130}
131
132//=====================================================================
133/// Generic problem class that will form the base class for both
134/// spine and elastic mesh-updates of the problem.
135/// Templated by the bulk element and interface element types
136//====================================================================
137template<class ELEMENT, class INTERFACE_ELEMENT>
138class InclinedPlaneProblem : public Problem
139{
140
141protected:
142
143 /// Bulk fluid mesh
145
146 /// Mesh for the traction elements that are added at inlet and outlet
148
149 /// Mesh for the free surface elements
151
152 /// Mesh for the point elements at each end of the free surface
154
155 /// Prefix for output files
156 std::string Output_prefix;
157
158public:
159
160 /// Generic Constructor (empty)
161 InclinedPlaneProblem(const unsigned &nx, const unsigned &ny,
162 const double &length) :
163 Output_prefix("Unset") { }
164
165 /// Solve the steady problem
167
168 /// Take n_tsteps timesteps of size dt
169 void timestep(const double &dt, const unsigned &n_tsteps);
170
171 /// Actions before the timestep
172 /// (update the time-dependent boundary conditions)
174 {
175 //Read out the current time
176 double time = this->time_pt()->time();
177 //Now add a temporary sinusoidal suction and blowing to the base
178 //Amplitude of the perturbation
179 double epsilon = 0.01;
180 //Loop over the nodes on the base
181 unsigned n_node = this->Bulk_mesh_pt->nboundary_node(0);
182 for(unsigned n=0;n<n_node;n++)
183 {
184 Node* nod_pt = this->Bulk_mesh_pt->boundary_node_pt(0,n);
185 double arg = Global_Physical_Variables::K*nod_pt->x(0);
186 double value = sin(arg)*epsilon*time*exp(-time);
187 nod_pt->set_value(1,value);
188 }
189 } //end_of_actions_before_implicit_timestep
190
191 /// Function to add the traction boundary elements to boundaries
192 /// 3(inlet) and 1(outlet) of the mesh
194 {
195 //Create a new (empty mesh)
196 Traction_mesh_pt = new Mesh;
197 //Inlet boundary conditions (boundary 3)
198 {
199 unsigned b = 3;
200 //Find the number of elements adjacent to mesh boundary
201 unsigned n_boundary_element = Bulk_mesh_pt->nboundary_element(b);
202 //Loop over these elements and create the traction elements
203 for(unsigned e=0;e<n_boundary_element;e++)
204 {
205 NavierStokesTractionElement<ELEMENT> *surface_element_pt =
206 new NavierStokesTractionElement<ELEMENT>
207 (Bulk_mesh_pt->boundary_element_pt(b,e),
208 Bulk_mesh_pt->face_index_at_boundary(b,e));
209 //Add the elements to the mesh
210 Traction_mesh_pt->add_element_pt(surface_element_pt);
211 //Set the traction function
212 surface_element_pt->traction_fct_pt() =
214 }
215 }
216
217 //Outlet boundary conditions (boundary 1)
218 {
219 unsigned b=1;
220 //Find the number of elements adjacent to mesh boundary
221 unsigned n_boundary_element = Bulk_mesh_pt->nboundary_element(b);
222 //Loop over these elements and create the traction elements
223 for(unsigned e=0;e<n_boundary_element;e++)
224 {
225 NavierStokesTractionElement<ELEMENT> *surface_element_pt =
226 new NavierStokesTractionElement<ELEMENT>
227 (Bulk_mesh_pt->boundary_element_pt(b,e),
228 Bulk_mesh_pt->face_index_at_boundary(b,e));
229 //Add the elements to the mesh
230 Traction_mesh_pt->add_element_pt(surface_element_pt);
231 //Set the traction function
232 surface_element_pt->traction_fct_pt() =
234 }
235 }
236 } //end of make_traction_elements
237
238 //Make the free surface elements on the top surface
240 {
241 //Create the (empty) meshes
242 Surface_mesh_pt = new Mesh;
243 Point_mesh_pt = new Mesh;
244
245 //The free surface is on the boundary 2
246 unsigned b = 2;
247 unsigned n_boundary_element = Bulk_mesh_pt->nboundary_element(b);
248 //Loop over the elements and create the appropriate interface elements
249 for(unsigned e=0;e<n_boundary_element;e++)
250 {
251 INTERFACE_ELEMENT *surface_element_pt =
252 new INTERFACE_ELEMENT
253 (Bulk_mesh_pt->boundary_element_pt(b,e),
254 Bulk_mesh_pt->face_index_at_boundary(b,e));
255 //Add elements to the mesh
256 Surface_mesh_pt->add_element_pt(surface_element_pt);
257 //Assign the capillary number to the free surface
258 surface_element_pt->ca_pt() =
260
261 //Make a point element from left-hand side of the
262 //first surface element (note that this relies on knowledge of
263 //the element order within the mesh)
264 if(e==0)
265 {
266 FluidInterfaceBoundingElement* point_element_pt =
267 surface_element_pt->make_bounding_element(-1);
268 //Add element to the point mesh
269 Point_mesh_pt->add_element_pt(point_element_pt);
270 //Set the capillary number
271 point_element_pt->ca_pt() = &Global_Physical_Variables::Ca;
272 //Set the wall normal
273 point_element_pt->wall_unit_normal_fct_pt() =
275 //Set the contact angle (using the strong version of the constraint)
276 point_element_pt->set_contact_angle(
278 }
279
280 //Make another point element from the right-hand side of the
281 //last surface element (note that this relies on knowledge of
282 //the element order within the mesh)
283 if(e==n_boundary_element-1)
284 {
285 FluidInterfaceBoundingElement* point_element_pt =
286 surface_element_pt->make_bounding_element(1);
287 //Add element to the mesh
288 Point_mesh_pt->add_element_pt(point_element_pt);
289 //Set the capillary number
290 point_element_pt->ca_pt() = &Global_Physical_Variables::Ca;
291 // Set the function that specifies the wall normal
292 point_element_pt->wall_unit_normal_fct_pt() =
294 }
295 }
296 } //end of make_free_surface_elements
297
298 /// Complete the build of the problem setting all standard
299 /// parameters and boundary conditions
301 {
302 using namespace Global_Physical_Variables;
303
304 //Complete the build of the fluid elements by passing physical parameters
305 //Find the number of bulk elements
306 unsigned n_element = Bulk_mesh_pt->nelement();
307 //Loop over all the fluid elements
308 for(unsigned e=0;e<n_element;e++)
309 {
310 //Cast to a fluid element
311 ELEMENT *temp_pt = dynamic_cast<ELEMENT*>(Bulk_mesh_pt->element_pt(e));
312
313 //Set the Reynolds number
314 temp_pt->re_pt() = &Re;
315 //The Strouhal number is 1, so ReSt = Re
316 temp_pt->re_st_pt() = &Re;
317 //Set the Reynolds number / Froude number
318 temp_pt->re_invfr_pt() = &ReInvFr;
319 //Set the direction of gravity
320 temp_pt->g_pt() = &G;
321 }
322
323 //------------Set the boundary conditions for this problem----------
324
325 {
326 //Determine whether we are solving an elastic problem or not
327 bool elastic = false;
328 if(dynamic_cast<SolidNode*>(Bulk_mesh_pt->node_pt(0))) {elastic=true;}
329
330 //Loop over the bottom of the mesh (the wall of the channel)
331 unsigned n_node = Bulk_mesh_pt->nboundary_node(0);
332 for(unsigned j=0;j<n_node;j++)
333 {
334 //Pin the u- and v- velocities
335 Bulk_mesh_pt->boundary_node_pt(0,j)->pin(0);
336 Bulk_mesh_pt->boundary_node_pt(0,j)->pin(1);
337
338 //If we are formulating the elastic problem pin both positions
339 //of nodes
340 if(elastic)
341 {
342 static_cast<SolidNode*>(Bulk_mesh_pt->boundary_node_pt(0,j))
343 ->pin_position(0);
344 static_cast<SolidNode*>(Bulk_mesh_pt->boundary_node_pt(0,j))
345 ->pin_position(1);
346 }
347 }
348
349 //Loop over the inlet and set the Dirichlet condition
350 //of no vertical velocity
351 n_node = Bulk_mesh_pt->nboundary_node(3);
352 for(unsigned j=0;j<n_node;j++)
353 {
354 Bulk_mesh_pt->boundary_node_pt(3,j)->pin(1);
355
356 //If elastic pin horizontal position of nodes
357 if(elastic)
358 {
359 static_cast<SolidNode*>(Bulk_mesh_pt->boundary_node_pt(3,j))
360 ->pin_position(0);
361 }
362 }
363
364 //Loop over the outlet and set the Dirichlet condition
365 //of no vertical velocity
366 n_node = Bulk_mesh_pt->nboundary_node(1);
367 for(unsigned j=0;j<n_node;j++)
368 {
369 Bulk_mesh_pt->boundary_node_pt(1,j)->pin(1);
370
371 //If elastic pin horizontal position
372 if(elastic)
373 {
374 static_cast<SolidNode*>(Bulk_mesh_pt->boundary_node_pt(1,j))
375 ->pin_position(0);
376 }
377 }
378 }
379
380 //Attach the boundary conditions to the mesh
381 std::cout << assign_eqn_numbers() << " in the main problem" << std::endl;
382 } //end of complete_build
383
384 /// Generic desructor to clean up the memory allocated in the problem
386 {
387 //Clear node storage before the mesh can be deleted,
388 //to avoid double deletion
389 this->Point_mesh_pt->flush_node_storage();
390 //Then delete the mesh
391 delete this->Point_mesh_pt;
392 //Clear node storage and then delete mesh
393 this->Surface_mesh_pt->flush_node_storage();
394 delete this->Surface_mesh_pt;
395 //Clear node storage and then delete mesh
396 this->Traction_mesh_pt->flush_node_storage();
397 delete this->Traction_mesh_pt;
398 //Delete the bulk mesh (no need to clear node storage)
399 delete this->Bulk_mesh_pt;
400 //Delete the time stepper
401 delete this->time_stepper_pt();
402 }
403
404};
405
406
407//-------------------------------------------------------------------------
408/// Solve the steady problem
409//-------------------------------------------------------------------------
410template<class ELEMENT,class INTERFACE_ELEMENT>
412{
413 //Load the namespace
414 using namespace Global_Physical_Variables;
415
416 //Initially set all nodes to the Nusselt flat-film solution
417 {
418 unsigned n_node = Bulk_mesh_pt->nnode();
419 for(unsigned n=0;n<n_node;n++)
420 {
421 double y = Bulk_mesh_pt->node_pt(n)->x(1);
422 //Top row
423 Bulk_mesh_pt->node_pt(n)->set_value(0,0.5*ReInvFr*sin(Alpha)*(2.0*y - y*y));
424 }
425 }
426
427 //Do one steady solve
428 steady_newton_solve();
429
430 //Output the full flow field
431 std::string filename = Output_prefix;;
432 filename.append("_output.dat");
433 ofstream file(filename.c_str());
434 Bulk_mesh_pt->output(file,5);
435 file.close();
436} //end of solve_steady
437
438
439//----------------------------------------------------------------------
440/// Perform n_tsteps timesteps of size dt
441//----------------------------------------------------------------------
442template<class ELEMENT, class INTERFACE_ELEMENT>
444timestep(const double &dt, const unsigned &n_tsteps)
445{
446 //Need to use the Global variables here
447 using namespace Global_Physical_Variables;
448
449 //Open an output file
450 std::string filename = Output_prefix;
451 filename.append("_time_trace.dat");
452 ofstream trace(filename.c_str());
453 //Counter that will be used to output the full flowfield
454 //at certain timesteps
455 int counter=0;
456
457 //Initial output of the time and the value of the vertical position at the
458 //left and right-hand end of the free surface
459 trace << time_pt()->time() << " "
460 << Bulk_mesh_pt->boundary_node_pt(2,0)->value(1)
461 << " "
462 << Bulk_mesh_pt->
463 boundary_node_pt(2, Bulk_mesh_pt->nboundary_node(2)-1)->x(1)
464 << " "
465 << std::endl;
466
467 //Loop over the desired number of timesteps
468 for(unsigned t=1;t<=n_tsteps;t++)
469 {
470 //Increase the counter
471 counter++;
472 cout << std::endl;
473 cout << "--------------TIMESTEP " << t<< " ------------------" << std::endl;
474
475 //Take a timestep of size dt
476 unsteady_newton_solve(dt);
477
478 //Uncomment to get full solution output
479 if(counter==2) //Change this number to get output every n steps
480 {
481 std::ofstream file;
482 std::ostringstream filename;
483 filename << Output_prefix << "_step" << Re << "_" << t << ".dat";
484 file.open(filename.str().c_str());
485 Bulk_mesh_pt->output(file,5);
486 file.close();
487
488 counter=0;
489 }
490
491 //Always output the interface
492 {
493 std::ofstream file;
494 std::ostringstream filename;
495 filename << Output_prefix << "_interface_" << Re << "_" << t << ".dat";
496 file.open(filename.str().c_str());
497 Surface_mesh_pt->output(file,5);
498 file.close();
499 }
500
501 //Output the time and value of the vertical position of the free surface
502 //at the left- and right-hand ends
503 trace << time_pt()->time() << " "
504 << Bulk_mesh_pt->boundary_node_pt(2,0)->x(1) << " "
505 <<
506 Bulk_mesh_pt->
507 boundary_node_pt(2,Bulk_mesh_pt->nboundary_node(2)-1)->x(1) << " "
508 << std::endl;
509 }
510} //end of timestep
511
512
513//====================================================================
514// Spine formulation of the problem
515//===================================================================
516
517
518//======================================================================
519/// Create a spine mesh for the problem
520//======================================================================
521template <class ELEMENT>
523 public SimpleRectangularQuadMesh<ELEMENT>,
524 public SpineMesh
525{
526public:
527 SpineInclinedPlaneMesh(const unsigned &nx, const unsigned &ny,
528 const double &lx, const double &ly,
529 TimeStepper* time_stepper_pt) :
530 SimpleRectangularQuadMesh<ELEMENT>
531 (nx,ny,lx,ly,time_stepper_pt), SpineMesh()
532 {
533 //Find the number of linear points in the element
534 unsigned n_p = dynamic_cast<ELEMENT*>(finite_element_pt(0))->nnode_1d();
535 //Reserve storage for the number of spines
536 Spine_pt.reserve((n_p-1)*nx + 1);
537
538 //Create single pointer to a spine
539 Spine* new_spine_pt=0;
540
541 //Now loop over the elements horizontally
542 for(unsigned long j=0;j<nx;j++)
543 {
544 //In most elements, we don't assign a spine to the last column,
545 //beacuse that will be done by the next element
546 unsigned n_pmax = n_p-1;
547 //In the last element, however, we must assign the final spine
548 if(j==nx-1) {n_pmax = n_p;}
549
550 //Loop over all nodes horizontally
551 for(unsigned l2=0;l2<n_pmax;l2++)
552 {
553 //Create a new spine with unit height and add to the mesh
554 new_spine_pt=new Spine(1.0);
555 Spine_pt.push_back(new_spine_pt);
556
557 // Get the node
558 SpineNode* nod_pt=element_node_pt(j,l2);
559 //Set the pointer to spine
560 nod_pt->spine_pt() = new_spine_pt;
561 //Set the fraction
562 nod_pt->fraction() = 0.0;
563 // Pointer to the mesh that implements the update fct
564 nod_pt->spine_mesh_pt() = this;
565
566 //Loop vertically along the spine
567 //Loop over the elements
568 for(unsigned long i=0;i<ny;i++)
569 {
570 //Loop over the vertical nodes, apart from the first
571 for(unsigned l1=1;l1<n_p;l1++)
572 {
573 // Get the node
574 SpineNode* nod_pt=element_node_pt(i*nx+j,l1*n_p+l2);
575 //Set the pointer to the spine
576 nod_pt->spine_pt() = new_spine_pt;
577 //Set the fraction
578 nod_pt->fraction()=(double(i)+double(l1)/double(n_p-1))/double(ny);
579 // Pointer to the mesh that implements the update fct
580 nod_pt->spine_mesh_pt() = this;
581 }
582 }
583 }
584 } //End of horizontal loop over elements
585 } //end of constructor
586
587 /// General node update function implements pure virtual function
588 /// defined in SpineMesh base class and performs specific node update
589 /// actions: along vertical spines
590 virtual void spine_node_update(SpineNode* spine_node_pt)
591 {
592 //Get fraction along the spine
593 double W = spine_node_pt->fraction();
594 //Get spine height
595 double H = spine_node_pt->h();
596 //Set the value of y
597 spine_node_pt->x(1) = W*H;
598 }
599};
600
601
602
603//============================================================================
604//Specific class for inclined plane problem using spines
605//============================================================================
606template<class ELEMENT, class TIMESTEPPER>
608 public InclinedPlaneProblem<ELEMENT,SpineLineFluidInterfaceElement<ELEMENT> >
609{
610public:
611
612 //Constructor
613 SpineInclinedPlaneProblem(const unsigned &nx, const unsigned &ny,
614 const double &length):
615 InclinedPlaneProblem<ELEMENT,SpineLineFluidInterfaceElement<ELEMENT> >
616 (nx,ny,length)
617 {
618 //Set the name
619 this->Output_prefix = "spine";
620
621 //Create our one and only timestepper, with adaptive timestepping
622 this->add_time_stepper_pt(new TIMESTEPPER);
623
624 //Create the bulk mesh
626 nx,ny,length,1.0,this->time_stepper_pt());
627
628 //Create the traction elements
630 //Create the free surface elements
632
633 //Add all sub meshes to the problem
634 this->add_sub_mesh(this->Bulk_mesh_pt);
635 this->add_sub_mesh(this->Traction_mesh_pt);
636 this->add_sub_mesh(this->Surface_mesh_pt);
637 this->add_sub_mesh(this->Point_mesh_pt);
638 //Create the global mesh
639 this->build_global_mesh();
640
641 //Complete the build of the problem
642 this->complete_build();
643 }
644
645 /// Spine heights/lengths are unknowns in the problem so their
646 /// values get corrected during each Newton step. However,
647 /// changing their value does not automatically change the
648 /// nodal positions, so we need to update all of them
650 {this->Bulk_mesh_pt->node_update();}
651
652};
653
654
655//====================================================================
656// Elastic formulation of the problem
657//===================================================================
658
659
660//======================================================================
661/// Create an Elastic mesh for the problem
662//======================================================================
663template <class ELEMENT>
665 public SimpleRectangularQuadMesh<ELEMENT>,
666 public SolidMesh
667{
668 //Public functions
669 public:
670 ElasticInclinedPlaneMesh(const unsigned &nx, const unsigned &ny,
671 const double &lx, const double &ly,
673 SimpleRectangularQuadMesh<ELEMENT>(nx,ny,lx,ly,time_stepper_pt), SolidMesh()
674 {
675 //Make the current configuration the undeformed one
677 }
678};
679
680
681
682
683//============================================================================
684//Specific class for inclined plane problem using pseudo-elastic formulation
685//============================================================================
686template<class ELEMENT, class TIMESTEPPER>
688 public InclinedPlaneProblem<ELEMENT,ElasticLineFluidInterfaceElement<ELEMENT> >
689{
690public:
691 //Constructor
692 ElasticInclinedPlaneProblem(const unsigned &nx, const unsigned &ny,
693 const double &length) :
695 (nx,ny,length)
696 {
697 //Set the name
698 this->Output_prefix = "elastic";
699
700 //Create our one and only timestepper, with adaptive timestepping
702
703 //Create the bulk mesh
705 nx,ny,length,1.0,this->time_stepper_pt());
706
707 //Set the consititutive law for the elements
708 unsigned n_element = this->Bulk_mesh_pt->nelement();
709 //Loop over all the fluid elements
710 for(unsigned e=0;e<n_element;e++)
711 {
712 //Cast to a fluid element
713 ELEMENT *temp_pt = dynamic_cast<ELEMENT*>(
714 this->Bulk_mesh_pt->element_pt(e));
715 //Set the constitutive law
716 temp_pt->constitutive_law_pt() =
718 }
719
720 //Create the traction elements
722 //Create the free surface element
724
725 //Add all sub meshes to the problem
726 this->add_sub_mesh(this->Bulk_mesh_pt);
727 this->add_sub_mesh(this->Traction_mesh_pt);
728 this->add_sub_mesh(this->Surface_mesh_pt);
729 this->add_sub_mesh(this->Point_mesh_pt);
730 //Create the global mesh
731 this->build_global_mesh();
732
733 //Complete the rest of the build
734 this->complete_build();
735 } //end of constructor
736
737 /// Update Lagrangian positions after each timestep
738 /// (updated-lagrangian approach)
740 {
741 //Now loop over all the nodes and reset their Lagrangian coordinates
742 unsigned n_node = this->Bulk_mesh_pt->nnode();
743 for(unsigned n=0;n<n_node;n++)
744 {
745 //Cast node to an elastic node
747 static_cast<SolidNode*>(this->Bulk_mesh_pt->node_pt(n));
748 for(unsigned j=0;j<2;j++) {temp_pt->xi(j) = temp_pt->x(j);}
749 }
750 } //end of actions_after_implicit_timestep
751
752};
753
754
755//start of main
756int main(int argc, char **argv)
757{
758 using namespace Global_Physical_Variables;
759
760 //Set the constitutive law for the mesh deformation
761 Constitutive_law_pt = new GeneralisedHookean(&Global_Physical_Variables::Nu);
762
763#ifdef CR_ELEMENT
764#define FLUID_ELEMENT QCrouzeixRaviartElement<2>
765#else
766#define FLUID_ELEMENT QTaylorHoodElement<2>
767#endif
768
769 //Initialise physical parameters
770 //Scale Reynolds number to be independent of alpha.
771 Re = 4.0/sin(Alpha);
772
773 //Set the direction of gravity
774 G[0] = sin(Alpha);
775 G[1] = -cos(Alpha);
776
777 //The wall normal to the inlet is in the negative x direction
778 Wall_normal.resize(2);
779 Wall_normal[0] = -1.0;
780 Wall_normal[1] = 0.0;
781
782 //Spine problem
783 {
784 //Create the problem
787
788 //Solve the steady problem
789 problem.solve_steady();
790
791 //Prepare the problem for timestepping
792 //(assume that it's been at the flat-film solution for all previous time)
793 double dt = 0.1;
794 problem.assign_initial_values_impulsive(dt);
795
796 //Timestep it
797 problem.timestep(dt,2);
798 } //End of spine problem
799
800
801 //Elastic problem
802 {
803 //Create the problem
807
808 //Solve the steady problem
809 problem.solve_steady();
810
811 //Prepare the problem for timestepping
812 //(assume that it's been at the flat-film solution for all previous time)
813 double dt = 0.1;
814 problem.assign_initial_values_impulsive(dt);
815
816 //Timestep it
817 problem.timestep(dt,2);
818 } //End of elastic problem
819}
Create an Elastic mesh for the problem.
ElasticInclinedPlaneMesh(const unsigned &nx, const unsigned &ny, const double &lx, const double &ly, TimeStepper *time_stepper_pt)
ElasticInclinedPlaneProblem(const unsigned &nx, const unsigned &ny, const double &length)
void actions_after_implicit_timestep()
Update Lagrangian positions after each timestep (updated-lagrangian approach)
Generic problem class that will form the base class for both spine and elastic mesh-updates of the pr...
void solve_steady()
Solve the steady problem.
Mesh * Bulk_mesh_pt
Bulk fluid mesh.
void make_traction_elements()
Function to add the traction boundary elements to boundaries 3(inlet) and 1(outlet) of the mesh.
InclinedPlaneProblem(const unsigned &nx, const unsigned &ny, const double &length)
Generic Constructor (empty)
std::string Output_prefix
Prefix for output files.
void timestep(const double &dt, const unsigned &n_tsteps)
Take n_tsteps timesteps of size dt.
Mesh * Traction_mesh_pt
Mesh for the traction elements that are added at inlet and outlet.
Mesh * Surface_mesh_pt
Mesh for the free surface elements.
~InclinedPlaneProblem()
Generic desructor to clean up the memory allocated in the problem.
Mesh * Point_mesh_pt
Mesh for the point elements at each end of the free surface.
void complete_build()
Complete the build of the problem setting all standard parameters and boundary conditions.
void actions_before_implicit_timestep()
Actions before the timestep (update the time-dependent boundary conditions)
Create a spine mesh for the problem.
virtual void spine_node_update(SpineNode *spine_node_pt)
General node update function implements pure virtual function defined in SpineMesh base class and per...
SpineInclinedPlaneMesh(const unsigned &nx, const unsigned &ny, const double &lx, const double &ly, TimeStepper *time_stepper_pt)
SpineInclinedPlaneProblem(const unsigned &nx, const unsigned &ny, const double &length)
void actions_before_newton_convergence_check()
Spine heights/lengths are unknowns in the problem so their values get corrected during each Newton st...
int main(int argc, char **argv)
double Inlet_Angle
The contact angle that is imposed at the inlet (pi)
ConstitutiveLaw * Constitutive_law_pt
Constitutive law used to determine the mesh deformation.
double Nu
Pseudo-solid Poisson ratio.
Vector< double > Wall_normal
Direction of the wall normal vector (at the inlet)
void wall_unit_normal_outlet_fct(const Vector< double > &x, Vector< double > &normal)
Function that specified the wall unit normal at the outlet.
double Ca
The Capillary number.
double Length
The length of the domain to fit the desired number of waves.
double K
Set the wavenumber.
void wall_unit_normal_inlet_fct(const Vector< double > &x, Vector< double > &normal)
Function that specifies the wall unit normal at the inlet.
double Alpha
Angle of incline of the slope (45 degrees)
double ReInvFr
The product of Reynolds number and inverse Froude number is set to two in this problem,...
void hydrostatic_pressure_outlet(const double &time, const Vector< double > &x, const Vector< double > &n, Vector< double > &traction)
Function that prescribes the hydrostatic pressure field at the outlet.
double Re
Reynolds number, based on the average velocity within the fluid film.
Vector< double > G(2, 0.0)
The Vector direction of gravity, set in main()
double N_wave
Set the number of waves desired in the domain.
void hydrostatic_pressure_inlet(const double &time, const Vector< double > &x, const Vector< double > &n, Vector< double > &traction)
Function that prescribes hydrostatic pressure field at the inlet.