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