33 #include "navier_stokes.h"
36 #include "fluid_interface.h"
39 #include "constitutive.h"
45 #include "meshes/rectangular_quadmesh.h"
49 using namespace oomph;
101 template <
class ELEMENT>
103 public virtual ElasticRefineableRectangularQuadMesh<ELEMENT>
119 const bool& periodic_in_x,
120 TimeStepper* time_stepper_pt=
121 &Mesh::Default_TimeStepper)
122 : RectangularQuadMesh<ELEMENT>(nx,ny1+ny2,lx,h1+h2,
123 periodic_in_x,time_stepper_pt),
124 ElasticRectangularQuadMesh<ELEMENT>(nx,ny1+ny2,lx,h1+h2,
125 periodic_in_x,time_stepper_pt),
126 ElasticRefineableRectangularQuadMesh<ELEMENT>(nx,ny1+ny2,lx,h1+h2,
135 this->set_nboundary(5);
138 for(
unsigned e=0;e<nx;e++)
141 FiniteElement* el_pt = this->finite_element_pt(nx*(ny1-1)+e);
144 const unsigned n_node = el_pt->nnode();
148 for(
unsigned n=0;n<3;n++)
150 Node* nod_pt = el_pt->node_pt(n_node-3+n);
151 this->convert_to_boundary_node(nod_pt);
152 this->add_boundary_node(4,nod_pt);
157 this->setup_boundary_element_info();
172 template <
class ELEMENT,
class TIMESTEPPER>
185 void set_initial_condition();
188 void set_boundary_conditions();
191 void doc_solution(DocInfo &doc_info);
194 void unsteady_run(
const double &t_max,
const double &dt);
208 Bulk_mesh_pt->set_lagrangian_nodal_coordinates();
212 void actions_before_adapt();
215 void actions_after_adapt();
218 void create_interface_elements();
221 void delete_interface_elements();
224 void deform_free_surface(
const double &epsilon,
const unsigned &n_periods);
228 const unsigned &pdof,
229 const double &pvalue)
232 dynamic_cast<ELEMENT*
>(mesh_pt()->element_pt(e))->
233 fix_pressure(pdof,pvalue);
258 template <
class ELEMENT,
class TIMESTEPPER>
263 add_time_stepper_pt(
new TIMESTEPPER);
266 const unsigned n_x = 3;
269 const unsigned n_y1 = 3;
272 const unsigned n_y2 = 3;
275 const double l_x = 1.0;
279 const double h1 = 1.0;
282 const double h2 = 1.0;
287 (n_x,n_y1,n_y2,l_x,h1,h2,
true,time_stepper_pt());
290 Bulk_mesh_pt->spatial_error_estimator_pt() =
new Z2ErrorEstimator;
293 Bulk_mesh_pt->max_refinement_level() = 4;
298 Surface_mesh_pt =
new Mesh;
302 create_interface_elements();
305 add_sub_mesh(Bulk_mesh_pt);
306 add_sub_mesh(Surface_mesh_pt);
319 const unsigned n_boundary = Bulk_mesh_pt->nboundary();
322 for(
unsigned b=0;b<n_boundary;b++)
325 const unsigned n_node = Bulk_mesh_pt->nboundary_node(b);
328 for(
unsigned n=0;n<n_node;n++)
335 if(b!=4) { Bulk_mesh_pt->boundary_node_pt(b,n)->pin(0); }
338 if(b==0 || b==2) { Bulk_mesh_pt->boundary_node_pt(b,n)->pin(1); }
344 if(b==0 || b==2) { Bulk_mesh_pt->boundary_node_pt(b,n)->pin_position(1); }
350 const unsigned n_node = Bulk_mesh_pt->nnode();
351 for(
unsigned n=0;n<n_node;n++) { Bulk_mesh_pt->node_pt(n)->pin_position(0); }
361 const unsigned n_lower = n_x*n_y1;
362 const unsigned n_upper = n_x*n_y2;
365 for(
unsigned e=0;e<n_lower;e++)
368 ELEMENT *el_pt =
dynamic_cast<ELEMENT*
>(Bulk_mesh_pt->element_pt(e));
389 for(
unsigned e=n_lower;e<(n_lower+n_upper);e++)
392 ELEMENT *el_pt =
dynamic_cast<ELEMENT*
>(Bulk_mesh_pt->element_pt(e));
419 fix_pressure(0,0,0.0);
422 set_boundary_conditions();
425 cout <<
"Number of equations: " << assign_eqn_numbers() << std::endl;
436 template <
class ELEMENT,
class TIMESTEPPER>
440 const unsigned n_node = mesh_pt()->nnode();
443 for(
unsigned n=0;n<n_node;n++)
446 for(
unsigned i=0;i<2;i++)
449 mesh_pt()->node_pt(n)->set_value(i,0.0);
455 assign_initial_values_impulsive();
466 template <
class ELEMENT,
class TIMESTEPPER>
470 const unsigned n_boundary = Bulk_mesh_pt->nboundary();
473 for(
unsigned b=0;b<n_boundary;b++)
476 const unsigned n_node = Bulk_mesh_pt->nboundary_node(b);
479 for(
unsigned n=0;n<n_node;n++)
483 if(b!=4) { Bulk_mesh_pt->boundary_node_pt(b,n)->set_value(0,0.0); }
489 Bulk_mesh_pt->boundary_node_pt(b,n)->set_value(1,0.0);
501 template<
class ELEMENT,
class TIMESTEPPER>
505 delete_interface_elements();
508 rebuild_global_mesh();
517 template<
class ELEMENT,
class TIMESTEPPER>
521 this->create_interface_elements();
524 rebuild_global_mesh();
527 const unsigned n_node = Bulk_mesh_pt->nnode();
528 for(
unsigned n=0;n<n_node;n++) { Bulk_mesh_pt->node_pt(n)->pin_position(0); }
531 RefineableNavierStokesEquations<2>::
532 unpin_all_pressure_dofs(Bulk_mesh_pt->element_pt());
535 RefineableNavierStokesEquations<2>::
536 pin_redundant_nodal_pressures(Bulk_mesh_pt->element_pt());
539 fix_pressure(0,0,0.0);
542 PVDEquationsBase<2>::pin_redundant_nodal_solid_pressures(
543 Bulk_mesh_pt->element_pt());
546 set_boundary_conditions();
557 template <
class ELEMENT,
class TIMESTEPPER>
561 const unsigned n_element = this->Bulk_mesh_pt->nboundary_element(4);
564 for(
unsigned e=0;e<n_element;e++)
567 ELEMENT* bulk_elem_pt =
dynamic_cast<ELEMENT*
>(
568 this->Bulk_mesh_pt->boundary_element_pt(4,e));
573 if(bulk_elem_pt->viscosity_ratio_pt()
577 const int face_index = this->Bulk_mesh_pt->face_index_at_boundary(4,e);
580 FiniteElement* interface_element_element_pt =
581 new ElasticLineFluidInterfaceElement<ELEMENT>(bulk_elem_pt,face_index);
584 this->Surface_mesh_pt->add_element_pt(interface_element_element_pt);
593 const unsigned n_interface_element = this->Surface_mesh_pt->nelement();
596 for(
unsigned e=0;e<n_interface_element;e++)
599 ElasticLineFluidInterfaceElement<ELEMENT>* el_pt =
600 dynamic_cast<ElasticLineFluidInterfaceElement<ELEMENT>*
>
601 (Surface_mesh_pt->element_pt(e));
618 template <
class ELEMENT,
class TIMESTEPPER>
622 const unsigned n_interface_element = Surface_mesh_pt->nelement();
625 for(
unsigned e=0;e<n_interface_element;e++)
627 delete Surface_mesh_pt->element_pt(e);
631 Surface_mesh_pt->flush_element_and_node_storage();
640 template <
class ELEMENT,
class TIMESTEPPER>
645 const unsigned n_node = Bulk_mesh_pt->nnode();
648 for(
unsigned n=0;n<n_node;n++)
651 const double current_x_pos = Bulk_mesh_pt->node_pt(n)->x(0);
652 const double current_y_pos = Bulk_mesh_pt->node_pt(n)->x(1);
655 const double new_y_pos = current_y_pos
656 + (1.0-fabs(1.0-current_y_pos))*epsilon
657 *(cos(2.0*n_periods*MathematicalConstants::Pi*current_x_pos/Lx));
660 Bulk_mesh_pt->node_pt(n)->x(1) = new_y_pos;
669 template <
class ELEMENT,
class TIMESTEPPER>
674 cout <<
"Time is now " << time_pt()->time() << std::endl;
677 ElasticLineFluidInterfaceElement<ELEMENT>* el_pt =
678 dynamic_cast<ElasticLineFluidInterfaceElement<ELEMENT>*
>
679 (Surface_mesh_pt->element_pt(0));
684 << el_pt->node_pt(0)->x(1) << std::endl;
690 const unsigned npts = 5;
693 sprintf(filename,
"%s/soln%i.dat",
694 doc_info.directory().c_str(),doc_info.number());
695 some_file.open(filename);
698 Bulk_mesh_pt->output(some_file,npts);
704 sprintf(filename,
"%s/interface_soln%i.dat",
705 doc_info.directory().c_str(),doc_info.number());
706 some_file.open(filename);
709 Surface_mesh_pt->output(some_file,npts);
721 template <
class ELEMENT,
class TIMESTEPPER>
727 const double epsilon = 0.1;
730 const unsigned n_periods = 1;
733 deform_free_surface(epsilon,n_periods);
739 doc_info.set_directory(
"RESLT");
746 sprintf(filename,
"%s/trace.dat",doc_info.directory().c_str());
750 Trace_file <<
"time, interface height" << std::endl;
756 set_initial_condition();
759 unsigned max_adapt = 2;
762 for(
unsigned i=0;i<2;i++) { refine_uniformly(); }
765 doc_solution(doc_info);
771 const unsigned n_timestep = unsigned(t_max/dt);
774 bool first_timestep =
true;
777 for(
unsigned t=1;t<=n_timestep;t++)
780 cout <<
"\nTimestep " << t <<
" of " << n_timestep << std::endl;
783 unsteady_newton_solve(dt,max_adapt,first_timestep);
786 first_timestep =
false;
792 doc_solution(doc_info);
810 int main(
int argc,
char* argv[])
813 CommandLineArgs::setup(argc,argv);
821 CommandLineArgs::specify_command_line_flag(
"--validation");
824 CommandLineArgs::parse_and_assign();
827 CommandLineArgs::doc_specified_flags();
834 std::ostringstream error_stream;
835 error_stream <<
"Definition of Global_Physical_Variables::ReSt is "
836 <<
"inconsistant with those\n"
837 <<
"of Global_Physical_Variables::Re and "
838 <<
"Global_Physical_Variables::St." << std::endl;
839 throw OomphLibError(error_stream.str(),
840 OOMPH_CURRENT_FUNCTION,OOMPH_EXCEPTION_LOCATION);
847 const double dt = 0.0025;
850 if(CommandLineArgs::command_line_flag_has_been_set(
"--validation"))
862 RefineableQCrouzeixRaviartElement<2>,RefineableQPVDElement<2,3> >, BDF<2> >
/////////////////////////////////////////////////////////////////////// /////////////////////////////...
ElasticRefineableTwoLayerMesh(const unsigned &nx, const unsigned &ny1, const unsigned &ny2, const double &lx, const double &h1, const double &h2, const bool &periodic_in_x, TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
Constructor: Pass number of elements in x-direction, number of elements in y-direction in bottom and ...
/////////////////////////////////////////////////////////////////////// /////////////////////////////...
Mesh * Surface_mesh_pt
Pointer to the "surface" mesh.
void set_initial_condition()
Set initial conditions.
void actions_after_adapt()
Rebuild the mesh of interface elements after adapting the bulk mesh.
void deform_free_surface(const double &epsilon, const unsigned &n_periods)
Deform the mesh/free surface to a prescribed function.
ofstream Trace_file
Trace file.
void doc_solution(DocInfo &doc_info)
Document the solution.
InterfaceProblem()
Constructor.
ConstitutiveLaw * Constitutive_law_pt
void set_boundary_conditions()
Set boundary conditions.
~InterfaceProblem()
Destructor (empty)
void actions_before_adapt()
Strip off the interface elements before adapting the bulk mesh.
void create_interface_elements()
Create the 1d interface elements.
double Lx
Width of domain.
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 delete_interface_elements()
Delete the 1d interface elements.
ElasticRefineableTwoLayerMesh< ELEMENT > * Bulk_mesh_pt
Pointer to the (specific) "bulk" mesh.
void actions_before_newton_solve()
No actions required before solve step.
void unsteady_run(const double &t_max, const double &dt)
Do unsteady run up to maximum time t_max with given timestep dt.
void actions_before_implicit_timestep()
Actions before the timestep: For maximum stability, reset the current nodal positions to be the "stre...
void actions_after_newton_solve()
No actions required after solve step.
int main(int argc, char *argv[])
/////////////////////////////////////////////////////////////////////// /////////////////////////////...
Namespace for physical parameters.
double ReSt
Womersley number (Reynolds x Strouhal)
Vector< double > G(2)
Direction of gravity.
double Nu
Pseudo-solid Poisson ratio.
double St
Strouhal number.
double Density_Ratio
Ratio of density in upper fluid to density in lower fluid. Reynolds number etc. is based on density i...
double Ca
Capillary number.
double ReInvFr
Product of Reynolds number and inverse of Froude number.
double Re
Reynolds number.
double Viscosity_Ratio
Ratio of viscosity in upper fluid to viscosity in lower fluid. Reynolds number etc....
ofstream Trace_file
Trace file.
ConstitutiveLaw * Constitutive_law_pt
Constitutive law used to determine the mesh deformation.