osc_ring_alg.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 2D Navier Stokes flow, driven by oscillating ring
27 // with pseudo-elasticity: The mean radius of ring is adjusted to
28 // allow conservation of volume (area).
29 
30 // Oomph-lib includes
31 #include "generic.h"
32 #include "navier_stokes.h"
33 
34 
35 //Need to instantiate templated mesh
36 #include "meshes/quarter_circle_sector_mesh.h"
37 
38 //Include namespace containing Sarah's asymptotics
40 
41 
42 using namespace std;
43 
44 using namespace oomph;
45 
46 using namespace MathematicalConstants;
47 
48 
49 
50 /// /////////////////////////////////////////////////////////////////////
51 /// /////////////////////////////////////////////////////////////////////
52 /// /////////////////////////////////////////////////////////////////////
53 
54 
55 //==================================================
56 /// Namespace for physical parameters
57 //==================================================
59 {
60 
61  /// Reynolds number
62  double Re=100.0; // ADJUST_PRIORITY
63 
64  /// Reynolds x Strouhal number
65  double ReSt=100.0; // ADJUST_PRIORITY
66 
67 }
68 
69 
70 
71 
72 /// /////////////////////////////////////////////////////////////////////
73 /// /////////////////////////////////////////////////////////////////////
74 /// /////////////////////////////////////////////////////////////////////
75 
76 
77 //====================================================================
78 /// Driver for oscillating ring problem: Wall performs oscillations
79 /// that resemble eigenmodes of freely oscillating ring and drives
80 /// viscous fluid flow. Mean radius of wall is adjustable and
81 /// responds to a pressure value in the fluid to allow for
82 /// mass conservation.
83 //====================================================================
84 template<class ELEMENT>
85 class OscRingNStProblem : public Problem
86 {
87 
88 public:
89 
90  /// Constructor: Pass timestep and function pointer to the
91  /// solution that provides the initial conditions for the fluid
92  OscRingNStProblem(const double& dt,
93  FiniteElement::UnsteadyExactSolutionFctPt IC_fct_pt);
94 
95  /// Destructor (empty)
97 
98  /// Get pointer to wall as geometric object
99  GeomObject* wall_pt()
100  {
101  return Wall_pt;
102  }
103 
104  /// Update after solve (empty)
106 
107  /// Update the problem specs before solve (empty)
109 
110  /// Update the problem specs before checking Newton
111  /// convergence: Update the fluid mesh and re-set velocity
112  /// boundary conditions -- no slip velocity on the wall means
113  /// that the velocity on the wall is dependent.
115  {
116  // Update the fluid mesh -- auxiliary update function for algebraic
117  // nodes automatically updates no slip condition.
118  fluid_mesh_pt()->node_update();
119  }
120 
121 
122  /// Update the problem specs after adaptation:
123  /// Set auxiliary update function that applies no slip on all
124  /// boundary nodes and choose fluid pressure dof that drives
125  /// the wall deformation
127  {
128  // Ring boundary: No slip; this also implies that the velocity needs
129  // to be updated in response to wall motion. This needs to be reset
130  // every time the mesh is changed -- there's no mechanism by which
131  // auxiliary update functions are copied to newly created nodes.
132  // (that's because unlike boundary conditions, they don't
133  // occur exclusively at boundaries)
134  unsigned ibound=1;
135  {
136  unsigned num_nod= fluid_mesh_pt()->nboundary_node(ibound);
137  for (unsigned inod=0;inod<num_nod;inod++)
138  {
139  fluid_mesh_pt()->boundary_node_pt(ibound,inod)->
140  set_auxiliary_node_update_fct_pt(
141  FSI_functions::apply_no_slip_on_moving_wall);
142  }
143  }
144 
145  // Set the reference pressure as the constant pressure in element 0
146  dynamic_cast<PseudoBucklingRingElement*>(Wall_pt)
147  ->set_reference_pressure_pt(fluid_mesh_pt()->element_pt(0)
148  ->internal_data_pt(0));
149  }
150 
151  /// Run the time integration for ntsteps steps
152  void unsteady_run(const unsigned &ntsteps, const bool& restarted,
153  DocInfo& doc_info);
154 
155  /// Set initial condition (incl previous timesteps) according
156  /// to specified function.
157  void set_initial_condition();
158 
159  /// Doc the solution
160  void doc_solution(DocInfo& doc_info);
161 
162  /// Access function for the fluid mesh
163  AlgebraicRefineableQuarterCircleSectorMesh<ELEMENT>* fluid_mesh_pt()
164  {
165  return Fluid_mesh_pt;
166  }
167 
168  /// Dump problem data.
169  void dump_it(ofstream& dump_file, DocInfo doc_info);
170 
171  /// Read problem data.
172  void restart(ifstream& restart_file);
173 
174 private:
175 
176  /// Write header for trace file
177  void write_trace_file_header();
178 
179  /// Function pointer to set the intial condition
180  FiniteElement::UnsteadyExactSolutionFctPt IC_Fct_pt;
181 
182  /// Pointer to wall
183  GeomObject* Wall_pt;
184 
185  /// Pointer to fluid mesh
186  AlgebraicRefineableQuarterCircleSectorMesh<ELEMENT>* Fluid_mesh_pt;
187 
188  /// Pointer to wall mesh (contains only a single GeneralisedElement)
190 
191  /// Trace file
192  ofstream Trace_file;
193 
194  /// Pointer to node on coarsest mesh on which velocity is traced
196 
197  /// Pointer to node in symmetry plane on coarsest mesh at
198  /// which velocity is traced
200 
201 };
202 
203 
204 //========================================================================
205 /// Constructor: Pass (constant) timestep and function pointer to the solution
206 /// that provides the initial conditions for the fluid.
207 //========================================================================
208 template<class ELEMENT>
210 FiniteElement::UnsteadyExactSolutionFctPt IC_fct_pt) : IC_Fct_pt(IC_fct_pt)
211 {
212 
213  // Period of oscillation
214  double T=1.0;
215 
216  //Allocate the timestepper
217  add_time_stepper_pt(new BDF<4>);
218 
219  // Initialise timestep -- also sets the weights for all timesteppers
220  // in the problem.
221  initialise_dt(dt);
222 
223  // Parameters for pseudo-buckling ring
224  double eps_buckl=0.1; // ADJUST_PRIORITY
225  double ampl_ratio=-0.5; // ADJUST_PRIORITY
226  unsigned n_buckl=2; // ADJUST_PRIORITY
227  double r_0=1.0;
228 
229  // Build wall geometric element
230  Wall_pt=new PseudoBucklingRingElement(eps_buckl,ampl_ratio,n_buckl,r_0,T,
231  time_stepper_pt());
232 
233  // Fluid mesh is suspended from wall between these two Lagrangian
234  // coordinates:
235  double xi_lo=0.0;
236  double xi_hi=2.0*atan(1.0);
237 
238  // Fractional position of dividing line for two outer blocks in fluid mesh
239  double fract_mid=0.5;
240 
241  // Build fluid mesh
242  Fluid_mesh_pt=new AlgebraicRefineableQuarterCircleSectorMesh<ELEMENT>(
243  Wall_pt,xi_lo,fract_mid,xi_hi,time_stepper_pt());
244 
245  // Set error estimator
246  Z2ErrorEstimator* error_estimator_pt=new Z2ErrorEstimator;
247  Fluid_mesh_pt->spatial_error_estimator_pt()=error_estimator_pt;
248 
249  // Fluid mesh is first sub-mesh
250  add_sub_mesh(Fluid_mesh_pt);
251 
252  // Build wall mesh
253  Wall_mesh_pt=new Mesh;
254 
255  // Wall mesh is completely empty: Add Wall element in its GeneralisedElement
256  // incarnation
257  Wall_mesh_pt->add_element_pt(dynamic_cast<GeneralisedElement*>(Wall_pt));
258 
259  // Wall mesh is second sub-mesh
260  add_sub_mesh(Wall_mesh_pt);
261 
262  // Combine all submeshes into a single Mesh
263  build_global_mesh();
264 
265  // Extract pointer to node at center of mesh (this node exists
266  // at all refinement levels and can be used to doc continuous timetrace
267  // of velocities)
268  unsigned nnode=fluid_mesh_pt()->finite_element_pt(0)->nnode();
269  Veloc_trace_node_pt=fluid_mesh_pt()->finite_element_pt(0)->node_pt(nnode-1);
270 
271  // Extract pointer to node in symmetry plane (this node exists
272  // at all refinement levels and can be used to doc continuous timetrace
273  // of velocities)
274  unsigned nnode_1d=dynamic_cast<ELEMENT*>(
275  fluid_mesh_pt()->finite_element_pt(0))->nnode_1d();
277  finite_element_pt(0)->node_pt(nnode_1d-1);
278 
279  // The "pseudo-elastic" wall element is "loaded" by a pressure.
280  // Use the "constant" pressure component in Crouzeix Raviart
281  // fluid element as that pressure.
282  dynamic_cast<PseudoBucklingRingElement*>(Wall_pt)
283  ->set_reference_pressure_pt(fluid_mesh_pt()
284  ->element_pt(0)->internal_data_pt(0));
285 
286 
287  // Set the boundary conditions for this problem:
288  //----------------------------------------------
289 
290  // All nodes are free by default -- just pin the ones that have
291  // Dirichlet conditions here.
292 
293  // Bottom boundary:
294  unsigned ibound=0;
295  {
296  unsigned num_nod= fluid_mesh_pt()->nboundary_node(ibound);
297  for (unsigned inod=0;inod<num_nod;inod++)
298  {
299  // Pin vertical velocity
300  {
301  fluid_mesh_pt()->boundary_node_pt(ibound,inod)->pin(1);
302  }
303  }
304  }
305 
306  // Ring boundary: No slip; this also implies that the velocity needs
307  // to be updated in response to wall motion
308  ibound=1;
309  {
310  unsigned num_nod=fluid_mesh_pt()->nboundary_node(ibound);
311  for (unsigned inod=0;inod<num_nod;inod++)
312  {
313  // Which node are we dealing with?
314  Node* node_pt=fluid_mesh_pt()->boundary_node_pt(ibound,inod);
315 
316  // Set auxiliary update function pointer to apply no-slip condition
317  // to velocities whenever nodal position is updated
318  node_pt->set_auxiliary_node_update_fct_pt(
319  FSI_functions::apply_no_slip_on_moving_wall);
320 
321  // Pin both velocities
322  for (unsigned i=0;i<2;i++)
323  {
324  node_pt->pin(i);
325  }
326  }
327  }
328 
329  // Left boundary:
330  ibound=2;
331  {
332  unsigned num_nod=fluid_mesh_pt()->nboundary_node(ibound);
333  for (unsigned inod=0;inod<num_nod;inod++)
334  {
335  // Pin horizontal velocity
336  {
337  fluid_mesh_pt()->boundary_node_pt(ibound,inod)->pin(0);
338  }
339  }
340  }
341 
342 
343  // Complete the build of all elements so they are fully functional
344  //----------------------------------------------------------------
345 
346  //Find number of elements in mesh
347  unsigned n_element = fluid_mesh_pt()->nelement();
348 
349  // Loop over the elements to set up element-specific
350  // things that cannot be handled by constructor
351  for(unsigned i=0;i<n_element;i++)
352  {
353  // Upcast from FiniteElement to the present element
354  ELEMENT *el_pt = dynamic_cast<ELEMENT*>(fluid_mesh_pt()->element_pt(i));
355 
356  //Set the Reynolds number, etc
357  el_pt->re_pt() = &Global_Physical_Variables::Re;
358  el_pt->re_st_pt() = &Global_Physical_Variables::ReSt;
359  }
360 
361 
362  //Attach the boundary conditions to the mesh
363  cout <<"Number of equations: " << assign_eqn_numbers() << std::endl;
364 
365 
366  // Set parameters for Sarah's asymptotic solution
367  //-----------------------------------------------
368 
369  // Amplitude of the oscillation
370  SarahBL::epsilon=static_cast<PseudoBucklingRingElement*>(Wall_pt)->
371  eps_buckl();
372 
373  // Womersley number is the same as square root of Reynolds number
375 
376  // Amplitude ratio
377  SarahBL::A=static_cast<PseudoBucklingRingElement*>(Wall_pt)->ampl_ratio();
378 
379  // Buckling wavenumber
380  SarahBL::N=static_cast<PseudoBucklingRingElement*>(Wall_pt)->n_buckl_float();
381 
382  // Frequency of oscillation (period is one)
383  SarahBL::Omega=2.0*MathematicalConstants::Pi;
384 
385 }
386 
387 
388 
389 
390 
391 //========================================================================
392 /// Set initial condition: Assign previous and current values
393 /// of the velocity from the velocity field specified via
394 /// the function pointer.
395 ///
396 /// Values are assigned so that the velocities and accelerations
397 /// are correct for current time.
398 //========================================================================
399 template<class ELEMENT>
401 {
402 
403  // Elastic wall: We're starting from a given initial state in which
404  // the wall is undeformed. If set_initial_condition() is called again
405  // after mesh refinement for first timestep, this needs to be reset.
406  dynamic_cast<PseudoBucklingRingElement*>(Wall_pt)->set_R_0(1.0);
407 
408  // Backup time in global timestepper
409  double backed_up_time=time_pt()->time();
410 
411  // Past history for velocities needs to be established for t=time0-deltat, ...
412  // Then provide current values (at t=time0) which will also form
413  // the initial guess for first solve at t=time0+deltat
414 
415  // Vector of exact solution values (includes pressure)
416  Vector<double> soln(3);
417  Vector<double> x(2);
418 
419  //Find number of nodes in mesh
420  unsigned num_nod = fluid_mesh_pt()->nnode();
421 
422  // Get continuous times at previous timesteps
423  int nprev_steps=time_stepper_pt()->nprev_values();
424  Vector<double> prev_time(nprev_steps+1);
425  for (int itime=nprev_steps;itime>=0;itime--)
426  {
427  prev_time[itime]=time_pt()->time(unsigned(itime));
428  }
429 
430  // Loop over current & previous timesteps (in outer loop because
431  // the mesh might also move)
432  for (int itime=nprev_steps;itime>=0;itime--)
433  {
434  double time=prev_time[itime];
435 
436  // Set global time (because this is how the geometric object refers
437  // to continous time
438  time_pt()->time()=time;
439 
440  cout << "setting IC at time =" << time << std::endl;
441 
442  // Update the fluid mesh for this value of the continuous time
443  // (The wall object reads the continous time from the same
444  // global time object)
445  fluid_mesh_pt()->node_update();
446 
447  // Loop over the nodes to set initial guess everywhere
448  for (unsigned jnod=0;jnod<num_nod;jnod++)
449  {
450 
451  // Get nodal coordinates
452  x[0]=fluid_mesh_pt()->node_pt(jnod)->x(0);
453  x[1]=fluid_mesh_pt()->node_pt(jnod)->x(1);
454 
455  // Get intial solution
456  (*IC_Fct_pt)(time,x,soln);
457 
458  // Loop over velocity directions (velocities are in soln[0] and soln[1]).
459  for (unsigned i=0;i<2;i++)
460  {
461  fluid_mesh_pt()->node_pt(jnod)->set_value(itime,i,soln[i]);
462  }
463 
464  // Loop over coordinate directions
465  for (unsigned i=0;i<2;i++)
466  {
467  fluid_mesh_pt()->node_pt(jnod)->x(itime,i)=x[i];
468  }
469  }
470  }
471 
472  // Reset backed up time for global timestepper
473  time_pt()->time()=backed_up_time;
474 
475 }
476 
477 
478 
479 
480 
481 //========================================================================
482 /// Doc the solution
483 ///
484 //========================================================================
485 template<class ELEMENT>
487 {
488 
489  cout << "Doc-ing step " << doc_info.number()
490  << " for time " << time_stepper_pt()->time_pt()->time() << std::endl;
491 
492  ofstream some_file;
493  char filename[100];
494 
495  // Number of plot points
496  unsigned npts;
497  npts=5;
498 
499 
500  // Output solution on fluid mesh
501  //-------------------------------
502  sprintf(filename,"%s/soln%i.dat",doc_info.directory().c_str(),
503  doc_info.number());
504  //some_file.precision(20);
505  some_file.open(filename);
506  unsigned nelem=fluid_mesh_pt()->nelement();
507  for (unsigned ielem=0;ielem<nelem;ielem++)
508  {
509  dynamic_cast<ELEMENT*>(fluid_mesh_pt()->element_pt(ielem))->
510  full_output(some_file,npts);
511  }
512  some_file.close();
513 
514 
515  // Plot wall posn
516  //---------------
517  sprintf(filename,"%s/Wall%i.dat",doc_info.directory().c_str(),
518  doc_info.number());
519  some_file.open(filename);
520 
521  unsigned nplot=100;
522  Vector<double > xi_wall(1);
523  Vector<double > r_wall(2);
524  for (unsigned iplot=0;iplot<nplot;iplot++)
525  {
526  xi_wall[0]=0.5*Pi*double(iplot)/double(nplot-1);
527  Wall_pt->position(xi_wall,r_wall);
528  some_file << r_wall[0] << " " << r_wall[1] << std::endl;
529  }
530  some_file.close();
531 
532 
533  // Doc Sarah's asymptotic solution
534  //--------------------------------
535  sprintf(filename,"%s/exact_soln%i.dat",doc_info.directory().c_str(),
536  doc_info.number());
537  some_file.open(filename);
538  fluid_mesh_pt()->output_fct(some_file,npts,
539  time_stepper_pt()->time_pt()->time(),
541  some_file.close();
542 
543 
544 
545 
546  // Get position of control point
547  //------------------------------
548  Vector<double> r(2);
549  Vector<double> xi(1);
550  xi[0]=MathematicalConstants::Pi/2.0;
551  wall_pt()->position(xi,r);
552 
553 
554 
555  // Get total volume (area) of computational domain, energies and average
556  //----------------------------------------------------------------------
557  // pressure
558  //---------
559  double area=0.0;
560  double press_int=0.0;
561  double diss=0.0;
562  double kin_en=0.0;
563  nelem=fluid_mesh_pt()->nelement();
564 
565  for (unsigned ielem=0;ielem<nelem;ielem++)
566  {
567  area+=fluid_mesh_pt()->finite_element_pt(ielem)->size();
568  press_int+=dynamic_cast<ELEMENT*>(fluid_mesh_pt()->element_pt(ielem))
569  ->pressure_integral();
570  diss+=dynamic_cast<ELEMENT*>(fluid_mesh_pt()->element_pt(ielem))->
571  dissipation();
572  kin_en+=dynamic_cast<ELEMENT*>(fluid_mesh_pt()->element_pt(ielem))->
573  kin_energy();
574  }
575 
576  // Total kinetic energy in entire domain
577  double global_kin=4.0*kin_en;
578 
579  // Max/min refinement level
580  unsigned min_level;
581  unsigned max_level;
582  fluid_mesh_pt()->get_refinement_levels(min_level,max_level);
583 
584 
585  // Total dissipation for quarter domain
586  double time=time_pt()->time();
587  double diss_asympt=SarahBL::Total_Diss_sarah(time)/4.0;
588 
589  // Asymptotic predictions for velocities at control point
590  Vector<double> x_sarah(2);
591  Vector<double> soln_sarah(3);
592  x_sarah[0]=Sarah_veloc_trace_node_pt->x(0);
593  x_sarah[1]=Sarah_veloc_trace_node_pt->x(1);
594  SarahBL::exact_soln(time,x_sarah,soln_sarah);
595 
596  // Doc
597  Trace_file << time_pt()->time()
598  << " " << r[1]
599  << " " << global_kin
600  << " " << SarahBL::Kin_energy_sarah(time_pt()->time())
601  << " " << static_cast<PseudoBucklingRingElement*>(Wall_pt)->r_0()
602  << " " << area
603  << " " << press_int/area
604  << " " << diss
605  << " " << diss_asympt
606  << " " << Veloc_trace_node_pt->x(0)
607  << " " << Veloc_trace_node_pt->x(1)
608  << " " << Veloc_trace_node_pt->value(0)
609  << " " << Veloc_trace_node_pt->value(1)
610  << " " << fluid_mesh_pt()->nelement()
611  << " " << ndof()
612  << " " << min_level
613  << " " << max_level
614  << " " << fluid_mesh_pt()->nrefinement_overruled()
615  << " " << fluid_mesh_pt()->max_error()
616  << " " << fluid_mesh_pt()->min_error()
617  << " " << fluid_mesh_pt()->max_permitted_error()
618  << " " << fluid_mesh_pt()->min_permitted_error()
619  << " " << fluid_mesh_pt()->max_keep_unrefined()
620  << " " << doc_info.number()
621  << " " << Sarah_veloc_trace_node_pt->x(0)
622  << " " << Sarah_veloc_trace_node_pt->x(1)
623  << " " << Sarah_veloc_trace_node_pt->value(0)
624  << " " << Sarah_veloc_trace_node_pt->value(1)
625  << " " << x_sarah[0]
626  << " " << x_sarah[1]
627  << " " << soln_sarah[0]
628  << " " << soln_sarah[1]
629  << " "
630  << static_cast<PseudoBucklingRingElement*>(Wall_pt)->r_0()-1.0
631  << std::endl;
632 
633 
634  // Output fluid solution on coarse-ish mesh
635  //-----------------------------------------
636 
637  // Extract all elements from quadtree representation
638  Vector<Tree*> all_element_pt;
639  fluid_mesh_pt()->forest_pt()->
640  stick_all_tree_nodes_into_vector(all_element_pt);
641 
642  // Build a coarse mesh
643  Mesh* coarse_mesh_pt = new Mesh();
644 
645  //Loop over all elements and check if their refinement level is
646  //equal to the mesh's minimum refinement level
647  nelem=all_element_pt.size();
648  for (unsigned ielem=0;ielem<nelem;ielem++)
649  {
650  Tree* el_pt=all_element_pt[ielem];
651  if (el_pt->level()==min_level)
652  {
653  coarse_mesh_pt->add_element_pt(el_pt->object_pt());
654  }
655  }
656 
657  // Output fluid solution on coarse mesh
658  sprintf(filename,"%s/coarse_soln%i.dat",doc_info.directory().c_str(),
659  doc_info.number());
660  some_file.open(filename);
661  nelem=coarse_mesh_pt->nelement();
662  for (unsigned ielem=0;ielem<nelem;ielem++)
663  {
664  dynamic_cast<ELEMENT*>(coarse_mesh_pt->element_pt(ielem))->
665  full_output(some_file,npts);
666  }
667  some_file.close();
668 
669  // Write restart file
670  sprintf(filename,"%s/restart%i.dat",doc_info.directory().c_str(),
671  doc_info.number());
672  some_file.open(filename);
673  dump_it(some_file,doc_info);
674  some_file.close();
675 
676 }
677 
678 
679 
680 
681 
682 
683 //========================================================================
684 /// Dump the solution
685 //========================================================================
686 template<class ELEMENT>
687 void OscRingNStProblem<ELEMENT>::dump_it(ofstream& dump_file,DocInfo doc_info)
688 {
689 
690  // Dump refinement status of mesh
691  //fluid_mesh_pt()->dump_refinement(dump_file);
692 
693  // Call generic dump()
694  Problem::dump(dump_file);
695 
696 }
697 
698 
699 
700 //========================================================================
701 /// Read solution from disk
702 //========================================================================
703 template<class ELEMENT>
704 void OscRingNStProblem<ELEMENT>::restart(ifstream& restart_file)
705 {
706  // Refine fluid mesh according to the instructions in restart file
707  //fluid_mesh_pt()->refine(restart_file);
708 
709  // Rebuild the global mesh
710  //rebuild_global_mesh();
711 
712  // Read generic problem data
713  Problem::read(restart_file);
714 
715 // // Complete build of all elements so they are fully functional
716 // finish_problem_setup();
717 
718  //Assign equation numbers
719  //cout <<"\nNumber of equations: " << assign_eqn_numbers()
720  // << std::endl<< std::endl;
721 }
722 
723 
724 
725 //========================================================================
726 /// Driver for timestepping the problem: Fixed timestep but
727 /// guaranteed spatial accuracy. Beautiful, innit?
728 ///
729 //========================================================================
730 template<class ELEMENT>
731 void OscRingNStProblem<ELEMENT>::unsteady_run(const unsigned& ntsteps,
732  const bool& restarted,
733  DocInfo& doc_info)
734 {
735 
736  // Open trace file
737  char filename[100];
738  sprintf(filename,"%s/trace.dat",doc_info.directory().c_str());
739  Trace_file.open(filename);
740 
741  // Max. number of adaptations per solve
742  unsigned max_adapt;
743 
744  // Max. number of adaptations per solve
745  if (restarted)
746  {
747  max_adapt=0;
748  }
749  else
750  {
751  max_adapt=1;
752  }
753 
754  // Max. and min. error for adaptive refinement/unrefinement
755  fluid_mesh_pt()->max_permitted_error()= 0.5e-2;
756  fluid_mesh_pt()->min_permitted_error()= 0.5e-3;
757 
758  // Don't allow refinement to drop under given level
759  fluid_mesh_pt()->min_refinement_level()=2;
760 
761  // Don't allow refinement beyond given level
762  fluid_mesh_pt()->max_refinement_level()=6;
763 
764  // Don't bother adapting the mesh if no refinement is required
765  // and if less than ... elements are to be merged.
766  fluid_mesh_pt()->max_keep_unrefined()=20;
767 
768 
769  // Get max/min refinement levels in mesh
770  unsigned min_refinement_level;
771  unsigned max_refinement_level;
772  fluid_mesh_pt()->get_refinement_levels(min_refinement_level,
773  max_refinement_level);
774 
775  cout << "\n Initial mesh: min/max refinement levels: "
776  << min_refinement_level << " " << max_refinement_level << std::endl << std::endl;
777 
778  // Doc refinement targets
779  fluid_mesh_pt()->doc_adaptivity_targets(cout);
780 
781  // Write header for trace file
782  write_trace_file_header();
783 
784  // Doc initial condition
785  doc_solution(doc_info);
786  doc_info.number()++;
787 
788  // Switch off doc during solve
789  doc_info.disable_doc();
790 
791  // If we set first to true, then initial guess will be re-assigned
792  // after mesh adaptation. Don't want this if we've done a restart.
793  bool first;
794  bool shift;
795  if (restarted)
796  {
797  first=false;
798  shift=false;
799  // Move time back by dt to make sure we're re-solving the read-in solution
800  time_pt()->time()-=time_pt()->dt();
801  }
802  else
803  {
804  first=true;
805  shift=true;
806  }
807 
808  //Take the first fixed timestep with specified tolerance for Newton solver
809  double dt=time_pt()->dt();
810  unsteady_newton_solve(dt,max_adapt,first,shift);
811 
812  // Switch doc back on
813  doc_info.enable_doc();
814 
815  // Doc initial solution
816  doc_solution(doc_info);
817  doc_info.number()++;
818 
819  // Now do normal run; allow for one mesh adaptation per timestep
820  max_adapt=1;
821 
822  //Loop over the remaining timesteps
823  for(unsigned t=2;t<=ntsteps;t++)
824  {
825 
826  // Switch off doc during solve
827  doc_info.disable_doc();
828 
829  //Take fixed timestep
830  first=false;
831  unsteady_newton_solve(dt,max_adapt,first);
832 
833  // Switch doc back on
834  doc_info.enable_doc();
835 
836  // Doc solution
837  //if (icount%10==0)
838  {
839  doc_solution(doc_info);
840  doc_info.number()++;
841  }
842  }
843 
844  // Close trace file
845  Trace_file.close();
846 
847 }
848 
849 
850 
851 
852 //========================================================================
853 /// Write trace file header
854 //========================================================================
855 template<class ELEMENT>
857 {
858 
859  // Doc parameters in trace file
860  Trace_file << "# err_max " << fluid_mesh_pt()->max_permitted_error() << std::endl;
861  Trace_file << "# err_min " << fluid_mesh_pt()->min_permitted_error() << std::endl;
862  Trace_file << "# Re " << Global_Physical_Variables::Re << std::endl;
863  Trace_file << "# St " << Global_Physical_Variables::ReSt/
864  Global_Physical_Variables::Re << std::endl;
865  Trace_file << "# dt " << time_stepper_pt()->time_pt()->dt() << std::endl;
866  Trace_file << "# initial # elements " << mesh_pt()->nelement() << std::endl;
867  Trace_file << "# min_refinement_level "
868  << fluid_mesh_pt()->min_refinement_level() << std::endl;
869  Trace_file << "# max_refinement_level "
870  << fluid_mesh_pt()->max_refinement_level() << std::endl;
871 
872 
873 
874  Trace_file << "VARIABLES=\"time\",\"V_c_t_r_l\",\"e_k_i_n\"";
875  Trace_file << ",\"e_k_i_n_(_a_s_y_m_p_t_)\",\"R_0\",\"Area\"" ;
876  Trace_file << ",\"Average pressure\",\"Total dissipation (quarter domain)\"";
877  Trace_file << ",\"Asymptotic dissipation (quarter domain)\"";
878  Trace_file << ",\"x<sub>1</sub><sup>(track)</sup>\"";
879  Trace_file << ",\"x<sub>2</sub><sup>(track)</sup>\"";
880  Trace_file << ",\"u<sub>1</sub><sup>(track)</sup>\"";
881  Trace_file << ",\"u<sub>2</sub><sup>(track)</sup>\"";
882  Trace_file << ",\"N<sub>element</sub>\"";
883  Trace_file << ",\"N<sub>dof</sub>\"";
884  Trace_file << ",\"max. refinement level\"";
885  Trace_file << ",\"min. refinement level\"";
886  Trace_file << ",\"# of elements whose refinement was over-ruled\"";
887  Trace_file << ",\"max. error\"";
888  Trace_file << ",\"min. error\"";
889  Trace_file << ",\"max. permitted error\"";
890  Trace_file << ",\"min. permitted error\"";
891  Trace_file << ",\"max. permitted # of unrefined elements\"";
892  Trace_file << ",\"doc number\"";
893  Trace_file << ",\"x<sub>1</sub><sup>(track2 FE)</sup>\"";
894  Trace_file << ",\"x<sub>2</sub><sup>(track2 FE)</sup>\"";
895  Trace_file << ",\"u<sub>1</sub><sup>(track2 FE)</sup>\"";
896  Trace_file << ",\"u<sub>2</sub><sup>(track2 FE)</sup>\"";
897  Trace_file << ",\"x<sub>1</sub><sup>(track2 Sarah)</sup>\"";
898  Trace_file << ",\"x<sub>2</sub><sup>(track2 Sarah)</sup>\"";
899  Trace_file << ",\"u<sub>1</sub><sup>(track2 Sarah)</sup>\"";
900  Trace_file << ",\"u<sub>2</sub><sup>(track2 Sarah)</sup>\"";
901  Trace_file << ",\"R<sub>0</sub>-1\"";
902  Trace_file << std::endl;
903 
904 }
905 
906 
907 
908 
909 /// /////////////////////////////////////////////////////////////////////
910 /// /////////////////////////////////////////////////////////////////////
911 /// /////////////////////////////////////////////////////////////////////
912 
913 
914 
915 
916 
917 
918 
919 //========================================================================
920 /// Demonstrate how to solve OscRingNSt problem in deformable domain
921 /// with mesh adaptation
922 //========================================================================
923 int main(int argc, char* argv[])
924 {
925 
926  // Store command line arguments
927  CommandLineArgs::setup(argc,argv);
928 
929  //Do a certain number of timesteps per period
930  unsigned nstep_per_period=40; // 80; // ADJUST_PRIORITY
931  unsigned nperiod=3;
932 
933  // Work out total number of steps and timestep
934  unsigned nstep=nstep_per_period*nperiod;
935  double dt=1.0/double(nstep_per_period);
936 
937  // Set up the problem: Pass timestep and Sarah's asymptotic solution for
938  // generation of initial condition
940  problem(dt,&SarahBL::full_exact_soln);
941 
942 
943  // Restart?
944  //---------
945  bool restarted=false;
946 
947  // Pointer to restart file
948  ifstream* restart_file_pt=0;
949 
950  // No restart
951  //-----------
952  if (CommandLineArgs::Argc!=2)
953  {
954  cout << "No restart" << std::endl;
955  restarted=false;
956 
957  // Refine uniformly
958  problem.refine_uniformly();
959  problem.refine_uniformly();
960  problem.refine_uniformly();
961 
962  // Set initial condition on uniformly refined domain (if this isn't done
963  // then the mesh contains the interpolation of the initial guess
964  // on the original coarse mesh -- if the nodal values were zero in
965  // the interior and nonzero on the boundary, then the the waviness of
966  // of the interpolated i.g. between the nodes on the coarse mesh
967  // gets transferred onto the fine mesh where we can do better
968  problem.set_initial_condition();
969 
970  }
971 
972  // Restart
973  //--------
974  else if (CommandLineArgs::Argc==2)
975  {
976  restarted=true;
977 
978  // Open restart file
979  restart_file_pt=new ifstream(CommandLineArgs::Argv[1],ios_base::in);
980  if (restart_file_pt!=0)
981  {
982  cout << "Have opened " << CommandLineArgs::Argv[1] <<
983  " for restart. " << std::endl;
984  }
985  else
986  {
987  std::ostringstream error_stream;
988  error_stream << "ERROR while trying to open "
989  << CommandLineArgs::Argv[2]
990  << " for restart." << std::endl;
991 
992  throw OomphLibError(error_stream.str(),
993  OOMPH_CURRENT_FUNCTION,
994  OOMPH_EXCEPTION_LOCATION);
995  }
996  // Do the actual restart
997  problem.restart(*restart_file_pt);
998  }
999 
1000 
1001  // Two command line arguments: do validation run
1002  if (CommandLineArgs::Argc==3)
1003  {
1004  nstep=3;
1005  cout << "Only doing nstep steps for validation: " << nstep << std::endl;
1006  }
1007 
1008 
1009 
1010 
1011  // Setup labels for output
1012  DocInfo doc_info;
1013 
1014  // Output directory
1015  doc_info.set_directory("RESLT");
1016 
1017 
1018  // Do unsteady run of the problem for nstep steps
1019  //-----------------------------------------------
1020  problem.unsteady_run(nstep,restarted,doc_info);
1021 
1022 
1023  // Validate the restart procedure
1024  //-------------------------------
1025  if (CommandLineArgs::Argc==3)
1026  {
1027 
1028  // Build problem and restart
1029 
1030  // Set up the problem: Pass timestep and Sarah's asymptotic solution for
1031  // generation of initial condition
1033  restarted_problem(dt,&SarahBL::full_exact_soln);
1034 
1035  // Setup labels for output
1036  DocInfo restarted_doc_info;
1037 
1038  // Output directory
1039  restarted_doc_info.set_directory("RESLT_restarted");
1040 
1041  // Step number
1042  restarted_doc_info.number()=0;
1043 
1044  // Copy by performing a restart from old problem
1045  restart_file_pt=new ifstream("RESLT/restart2.dat");
1046 
1047  // Do the actual restart
1048  restarted_problem.restart(*restart_file_pt);
1049 
1050  // Do unsteady run of the problem for one step
1051  unsigned nstep=2;
1052  bool restarted=true;
1053  restarted_problem.unsteady_run(nstep,restarted,restarted_doc_info);
1054 
1055  }
1056 
1057 }
1058 
1059 
1060 
1061 
1062 
1063 
1064 
1065 
1066 
1067 
///////////////////////////////////////////////////////////////////// ///////////////////////////////...
Definition: osc_ring_alg.cc:86
void unsteady_run(const unsigned &ntsteps, const bool &restarted, DocInfo &doc_info)
Run the time integration for ntsteps steps.
void restart(ifstream &restart_file)
Read problem data.
GeomObject * Wall_pt
Pointer to wall.
void actions_after_adapt()
Update the problem specs after adaptation: Set auxiliary update function that applies no slip on all ...
ofstream Trace_file
Trace file.
Node * Sarah_veloc_trace_node_pt
Pointer to node in symmetry plane on coarsest mesh at which velocity is traced.
void actions_after_newton_solve()
Update after solve (empty)
void actions_before_newton_solve()
Update the problem specs before solve (empty)
void doc_solution(DocInfo &doc_info)
Doc the solution.
GeomObject * wall_pt()
Get pointer to wall as geometric object.
Definition: osc_ring_alg.cc:99
void write_trace_file_header()
Write header for trace file.
~OscRingNStProblem()
Destructor (empty)
Definition: osc_ring_alg.cc:96
void actions_before_newton_convergence_check()
Update the problem specs before checking Newton convergence: Update the fluid mesh and re-set velocit...
Mesh * Wall_mesh_pt
Pointer to wall mesh (contains only a single GeneralisedElement)
void set_initial_condition()
Set initial condition (incl previous timesteps) according to specified function.
OscRingNStProblem(const double &dt, FiniteElement::UnsteadyExactSolutionFctPt IC_fct_pt)
Constructor: Pass timestep and function pointer to the solution that provides the initial conditions ...
FiniteElement::UnsteadyExactSolutionFctPt IC_Fct_pt
Function pointer to set the intial condition.
AlgebraicRefineableQuarterCircleSectorMesh< ELEMENT > * fluid_mesh_pt()
Access function for the fluid mesh.
Node * Veloc_trace_node_pt
Pointer to node on coarsest mesh on which velocity is traced.
void dump_it(ofstream &dump_file, DocInfo doc_info)
Dump problem data.
AlgebraicRefineableQuarterCircleSectorMesh< ELEMENT > * Fluid_mesh_pt
Pointer to fluid mesh.
Namespace for physical parameters.
Definition: fsi_osc_ring.cc:46
double ReSt
Reynolds x Strouhal number.
Definition: fsi_osc_ring.cc:77
double Re
Reynolds number.
Definition: fsi_osc_ring.cc:74
double Kin_energy_sarah(double t)
void full_exact_soln(const double &time, const Vector< double > &x, Vector< double > &soln)
Full exact solution: x,y,u,v,p,du/dt,dv/dt,diss.
void exact_soln(const double &time, const Vector< double > &x, Vector< double > &soln)
Exact solution: x,y,u,v,p.
double Total_Diss_sarah(double t)
int main(int argc, char *argv[])
///////////////////////////////////////////////////////////////////// ///////////////////////////////...