#include "generic.h"
#include "navier_stokes.h"
#include "meshes/quarter_circle_sector_mesh.h"
using namespace std;
using namespace MathematicalConstants;
{
}
template<class ELEMENT>
{
public:
FiniteElement::UnsteadyExactSolutionFctPt IC_fct_pt);
GeomObject* wall_pt()
{
return Wall_pt;
}
void actions_after_newton_solve(){}
void actions_before_newton_solve(){}
void actions_before_newton_convergence_check()
{
fluid_mesh_pt()->node_update();
}
void actions_after_adapt()
{
unsigned ibound=1;
{
unsigned num_nod= fluid_mesh_pt()->nboundary_node(ibound);
for (unsigned inod=0;inod<num_nod;inod++)
{
fluid_mesh_pt()->boundary_node_pt(ibound,inod)->
set_auxiliary_node_update_fct_pt(
FSI_functions::apply_no_slip_on_moving_wall);
}
}
dynamic_cast<PseudoBucklingRingElement*>(Wall_pt)
->set_reference_pressure_pt(fluid_mesh_pt()->element_pt(0)
->internal_data_pt(0));
}
void unsteady_run(const unsigned &ntsteps, const bool& restarted,
DocInfo& doc_info);
void set_initial_condition();
void doc_solution(DocInfo& doc_info);
AlgebraicRefineableQuarterCircleSectorMesh<ELEMENT>* fluid_mesh_pt()
{
return Fluid_mesh_pt;
}
void dump_it(ofstream& dump_file, DocInfo doc_info);
void restart(ifstream& restart_file);
private:
void write_trace_file_header();
FiniteElement::UnsteadyExactSolutionFctPt IC_Fct_pt;
GeomObject* Wall_pt;
AlgebraicRefineableQuarterCircleSectorMesh<ELEMENT>* Fluid_mesh_pt;
Mesh* Wall_mesh_pt;
ofstream Trace_file;
Node* Veloc_trace_node_pt;
Node* Sarah_veloc_trace_node_pt;
};
template<class ELEMENT>
FiniteElement::UnsteadyExactSolutionFctPt IC_fct_pt) : IC_Fct_pt(IC_fct_pt)
{
double T=1.0;
add_time_stepper_pt(new BDF<4>);
initialise_dt(dt);
double eps_buckl=0.1;
double ampl_ratio=-0.5;
unsigned n_buckl=2;
double r_0=1.0;
Wall_pt=new PseudoBucklingRingElement(eps_buckl,ampl_ratio,n_buckl,r_0,T,
time_stepper_pt());
double xi_lo=0.0;
double xi_hi=2.0*atan(1.0);
double fract_mid=0.5;
Fluid_mesh_pt=new AlgebraicRefineableQuarterCircleSectorMesh<ELEMENT>(
Wall_pt,xi_lo,fract_mid,xi_hi,time_stepper_pt());
Z2ErrorEstimator* error_estimator_pt=new Z2ErrorEstimator;
Fluid_mesh_pt->spatial_error_estimator_pt()=error_estimator_pt;
add_sub_mesh(Fluid_mesh_pt);
Wall_mesh_pt=new Mesh;
Wall_mesh_pt->add_element_pt(dynamic_cast<GeneralisedElement*>(Wall_pt));
add_sub_mesh(Wall_mesh_pt);
build_global_mesh();
unsigned nnode=fluid_mesh_pt()->finite_element_pt(0)->nnode();
Veloc_trace_node_pt=fluid_mesh_pt()->finite_element_pt(0)->node_pt(nnode-1);
unsigned nnode_1d=dynamic_cast<ELEMENT*>(
fluid_mesh_pt()->finite_element_pt(0))->nnode_1d();
Sarah_veloc_trace_node_pt=fluid_mesh_pt()->
finite_element_pt(0)->node_pt(nnode_1d-1);
dynamic_cast<PseudoBucklingRingElement*>(Wall_pt)
->set_reference_pressure_pt(fluid_mesh_pt()
->element_pt(0)->internal_data_pt(0));
unsigned ibound=0;
{
unsigned num_nod= fluid_mesh_pt()->nboundary_node(ibound);
for (unsigned inod=0;inod<num_nod;inod++)
{
{
fluid_mesh_pt()->boundary_node_pt(ibound,inod)->pin(1);
}
}
}
ibound=1;
{
unsigned num_nod=fluid_mesh_pt()->nboundary_node(ibound);
for (unsigned inod=0;inod<num_nod;inod++)
{
Node* node_pt=fluid_mesh_pt()->boundary_node_pt(ibound,inod);
node_pt->set_auxiliary_node_update_fct_pt(
FSI_functions::apply_no_slip_on_moving_wall);
for (unsigned i=0;i<2;i++)
{
node_pt->pin(i);
}
}
}
ibound=2;
{
unsigned num_nod=fluid_mesh_pt()->nboundary_node(ibound);
for (unsigned inod=0;inod<num_nod;inod++)
{
{
fluid_mesh_pt()->boundary_node_pt(ibound,inod)->pin(0);
}
}
}
unsigned n_element = fluid_mesh_pt()->nelement();
for(unsigned i=0;i<n_element;i++)
{
ELEMENT *el_pt = dynamic_cast<ELEMENT*>(fluid_mesh_pt()->element_pt(i));
}
cout <<"Number of equations: " << assign_eqn_numbers() << std::endl;
eps_buckl();
SarahBL::A=
static_cast<PseudoBucklingRingElement*
>(Wall_pt)->ampl_ratio();
SarahBL::N=
static_cast<PseudoBucklingRingElement*
>(Wall_pt)->n_buckl_float();
}
template<class ELEMENT>
{
dynamic_cast<PseudoBucklingRingElement*>(Wall_pt)->set_R_0(1.0);
double backed_up_time=time_pt()->time();
Vector<double> soln(3);
Vector<double> x(2);
unsigned num_nod = fluid_mesh_pt()->nnode();
int nprev_steps=time_stepper_pt()->nprev_values();
Vector<double> prev_time(nprev_steps+1);
for (int itime=nprev_steps;itime>=0;itime--)
{
prev_time[itime]=time_pt()->time(unsigned(itime));
}
for (int itime=nprev_steps;itime>=0;itime--)
{
double time=prev_time[itime];
time_pt()->time()=time;
cout << "setting IC at time =" << time << std::endl;
fluid_mesh_pt()->node_update();
for (unsigned jnod=0;jnod<num_nod;jnod++)
{
x[0]=fluid_mesh_pt()->node_pt(jnod)->x(0);
x[1]=fluid_mesh_pt()->node_pt(jnod)->x(1);
(*IC_Fct_pt)(time,x,soln);
for (unsigned i=0;i<2;i++)
{
fluid_mesh_pt()->node_pt(jnod)->set_value(itime,i,soln[i]);
}
for (unsigned i=0;i<2;i++)
{
fluid_mesh_pt()->node_pt(jnod)->x(itime,i)=x[i];
}
}
}
time_pt()->time()=backed_up_time;
}
template<class ELEMENT>
{
cout << "Doc-ing step " << doc_info.number()
<< " for time " << time_stepper_pt()->time_pt()->time() << std::endl;
ofstream some_file;
char filename[100];
unsigned npts;
npts=5;
sprintf(filename,"%s/soln%i.dat",doc_info.directory().c_str(),
doc_info.number());
some_file.open(filename);
unsigned nelem=fluid_mesh_pt()->nelement();
for (unsigned ielem=0;ielem<nelem;ielem++)
{
dynamic_cast<ELEMENT*>(fluid_mesh_pt()->element_pt(ielem))->
full_output(some_file,npts);
}
some_file.close();
sprintf(filename,"%s/Wall%i.dat",doc_info.directory().c_str(),
doc_info.number());
some_file.open(filename);
unsigned nplot=100;
Vector<double > xi_wall(1);
Vector<double > r_wall(2);
for (unsigned iplot=0;iplot<nplot;iplot++)
{
xi_wall[0]=0.5*Pi*double(iplot)/double(nplot-1);
Wall_pt->position(xi_wall,r_wall);
some_file << r_wall[0] << " " << r_wall[1] << std::endl;
}
some_file.close();
sprintf(filename,"%s/exact_soln%i.dat",doc_info.directory().c_str(),
doc_info.number());
some_file.open(filename);
fluid_mesh_pt()->output_fct(some_file,npts,
time_stepper_pt()->time_pt()->time(),
some_file.close();
Vector<double> r(2);
Vector<double> xi(1);
xi[0]=MathematicalConstants::Pi/2.0;
wall_pt()->position(xi,r);
double area=0.0;
double press_int=0.0;
double diss=0.0;
double kin_en=0.0;
nelem=fluid_mesh_pt()->nelement();
for (unsigned ielem=0;ielem<nelem;ielem++)
{
area+=fluid_mesh_pt()->finite_element_pt(ielem)->size();
press_int+=dynamic_cast<ELEMENT*>(fluid_mesh_pt()->element_pt(ielem))
->pressure_integral();
diss+=dynamic_cast<ELEMENT*>(fluid_mesh_pt()->element_pt(ielem))->
dissipation();
kin_en+=dynamic_cast<ELEMENT*>(fluid_mesh_pt()->element_pt(ielem))->
kin_energy();
}
double global_kin=4.0*kin_en;
unsigned min_level;
unsigned max_level;
fluid_mesh_pt()->get_refinement_levels(min_level,max_level);
double time=time_pt()->time();
Vector<double> x_sarah(2);
Vector<double> soln_sarah(3);
x_sarah[0]=Sarah_veloc_trace_node_pt->x(0);
x_sarah[1]=Sarah_veloc_trace_node_pt->x(1);
Trace_file << time_pt()->time()
<< " " << r[1]
<< " " << global_kin
<< " " << static_cast<PseudoBucklingRingElement*>(Wall_pt)->r_0()
<< " " << area
<< " " << press_int/area
<< " " << diss
<< " " << diss_asympt
<< " " << Veloc_trace_node_pt->x(0)
<< " " << Veloc_trace_node_pt->x(1)
<< " " << Veloc_trace_node_pt->value(0)
<< " " << Veloc_trace_node_pt->value(1)
<< " " << fluid_mesh_pt()->nelement()
<< " " << ndof()
<< " " << min_level
<< " " << max_level
<< " " << fluid_mesh_pt()->nrefinement_overruled()
<< " " << fluid_mesh_pt()->max_error()
<< " " << fluid_mesh_pt()->min_error()
<< " " << fluid_mesh_pt()->max_permitted_error()
<< " " << fluid_mesh_pt()->min_permitted_error()
<< " " << fluid_mesh_pt()->max_keep_unrefined()
<< " " << doc_info.number()
<< " " << Sarah_veloc_trace_node_pt->x(0)
<< " " << Sarah_veloc_trace_node_pt->x(1)
<< " " << Sarah_veloc_trace_node_pt->value(0)
<< " " << Sarah_veloc_trace_node_pt->value(1)
<< " " << x_sarah[0]
<< " " << x_sarah[1]
<< " " << soln_sarah[0]
<< " " << soln_sarah[1]
<< " "
<< static_cast<PseudoBucklingRingElement*>(Wall_pt)->r_0()-1.0
<< std::endl;
Vector<Tree*> all_element_pt;
fluid_mesh_pt()->forest_pt()->
stick_all_tree_nodes_into_vector(all_element_pt);
Mesh* coarse_mesh_pt = new Mesh();
nelem=all_element_pt.size();
for (unsigned ielem=0;ielem<nelem;ielem++)
{
Tree* el_pt=all_element_pt[ielem];
if (el_pt->level()==min_level)
{
coarse_mesh_pt->add_element_pt(el_pt->object_pt());
}
}
sprintf(filename,"%s/coarse_soln%i.dat",doc_info.directory().c_str(),
doc_info.number());
some_file.open(filename);
nelem=coarse_mesh_pt->nelement();
for (unsigned ielem=0;ielem<nelem;ielem++)
{
dynamic_cast<ELEMENT*>(coarse_mesh_pt->element_pt(ielem))->
full_output(some_file,npts);
}
some_file.close();
sprintf(filename,"%s/restart%i.dat",doc_info.directory().c_str(),
doc_info.number());
some_file.open(filename);
dump_it(some_file,doc_info);
some_file.close();
}
template<class ELEMENT>
{
Problem::dump(dump_file);
}
template<class ELEMENT>
{
Problem::read(restart_file);
}
template<class ELEMENT>
const bool& restarted,
DocInfo& doc_info)
{
char filename[100];
sprintf(filename,"%s/trace.dat",doc_info.directory().c_str());
Trace_file.open(filename);
unsigned max_adapt;
if (restarted)
{
max_adapt=0;
}
else
{
max_adapt=1;
}
fluid_mesh_pt()->max_permitted_error()= 0.5e-2;
fluid_mesh_pt()->min_permitted_error()= 0.5e-3;
fluid_mesh_pt()->min_refinement_level()=2;
fluid_mesh_pt()->max_refinement_level()=6;
fluid_mesh_pt()->max_keep_unrefined()=20;
unsigned min_refinement_level;
unsigned max_refinement_level;
fluid_mesh_pt()->get_refinement_levels(min_refinement_level,
max_refinement_level);
cout << "\n Initial mesh: min/max refinement levels: "
<< min_refinement_level << " " << max_refinement_level << std::endl << std::endl;
fluid_mesh_pt()->doc_adaptivity_targets(cout);
write_trace_file_header();
doc_solution(doc_info);
doc_info.number()++;
doc_info.disable_doc();
bool first;
bool shift;
if (restarted)
{
first=false;
shift=false;
time_pt()->time()-=time_pt()->dt();
}
else
{
first=true;
shift=true;
}
double dt=time_pt()->dt();
unsteady_newton_solve(dt,max_adapt,first,shift);
doc_info.enable_doc();
doc_solution(doc_info);
doc_info.number()++;
max_adapt=1;
for(unsigned t=2;t<=ntsteps;t++)
{
doc_info.disable_doc();
first=false;
unsteady_newton_solve(dt,max_adapt,first);
doc_info.enable_doc();
{
doc_solution(doc_info);
doc_info.number()++;
}
}
Trace_file.close();
}
template<class ELEMENT>
{
Trace_file << "# err_max " << fluid_mesh_pt()->max_permitted_error() << std::endl;
Trace_file << "# err_min " << fluid_mesh_pt()->min_permitted_error() << std::endl;
Trace_file << "# dt " << time_stepper_pt()->time_pt()->dt() << std::endl;
Trace_file << "# initial # elements " << mesh_pt()->nelement() << std::endl;
Trace_file << "# min_refinement_level "
<< fluid_mesh_pt()->min_refinement_level() << std::endl;
Trace_file << "# max_refinement_level "
<< fluid_mesh_pt()->max_refinement_level() << std::endl;
Trace_file << "VARIABLES=\"time\",\"V_c_t_r_l\",\"e_k_i_n\"";
Trace_file << ",\"e_k_i_n_(_a_s_y_m_p_t_)\",\"R_0\",\"Area\"" ;
Trace_file << ",\"Average pressure\",\"Total dissipation (quarter domain)\"";
Trace_file << ",\"Asymptotic dissipation (quarter domain)\"";
Trace_file << ",\"x<sub>1</sub><sup>(track)</sup>\"";
Trace_file << ",\"x<sub>2</sub><sup>(track)</sup>\"";
Trace_file << ",\"u<sub>1</sub><sup>(track)</sup>\"";
Trace_file << ",\"u<sub>2</sub><sup>(track)</sup>\"";
Trace_file << ",\"N<sub>element</sub>\"";
Trace_file << ",\"N<sub>dof</sub>\"";
Trace_file << ",\"max. refinement level\"";
Trace_file << ",\"min. refinement level\"";
Trace_file << ",\"# of elements whose refinement was over-ruled\"";
Trace_file << ",\"max. error\"";
Trace_file << ",\"min. error\"";
Trace_file << ",\"max. permitted error\"";
Trace_file << ",\"min. permitted error\"";
Trace_file << ",\"max. permitted # of unrefined elements\"";
Trace_file << ",\"doc number\"";
Trace_file << ",\"x<sub>1</sub><sup>(track2 FE)</sup>\"";
Trace_file << ",\"x<sub>2</sub><sup>(track2 FE)</sup>\"";
Trace_file << ",\"u<sub>1</sub><sup>(track2 FE)</sup>\"";
Trace_file << ",\"u<sub>2</sub><sup>(track2 FE)</sup>\"";
Trace_file << ",\"x<sub>1</sub><sup>(track2 Sarah)</sup>\"";
Trace_file << ",\"x<sub>2</sub><sup>(track2 Sarah)</sup>\"";
Trace_file << ",\"u<sub>1</sub><sup>(track2 Sarah)</sup>\"";
Trace_file << ",\"u<sub>2</sub><sup>(track2 Sarah)</sup>\"";
Trace_file << ",\"R<sub>0</sub>-1\"";
Trace_file << std::endl;
}
int main(
int argc,
char* argv[])
{
CommandLineArgs::setup(argc,argv);
unsigned nstep_per_period=40;
unsigned nperiod=3;
unsigned nstep=nstep_per_period*nperiod;
double dt=1.0/double(nstep_per_period);
bool restarted=false;
ifstream* restart_file_pt=0;
if (CommandLineArgs::Argc!=2)
{
cout << "No restart" << std::endl;
restarted=false;
problem.refine_uniformly();
problem.refine_uniformly();
problem.refine_uniformly();
problem.set_initial_condition();
}
else if (CommandLineArgs::Argc==2)
{
restarted=true;
restart_file_pt=new ifstream(CommandLineArgs::Argv[1],ios_base::in);
if (restart_file_pt!=0)
{
cout << "Have opened " << CommandLineArgs::Argv[1] <<
" for restart. " << std::endl;
}
else
{
std::ostringstream error_stream;
error_stream << "ERROR while trying to open "
<< CommandLineArgs::Argv[2]
<< " for restart." << std::endl;
throw OomphLibError(error_stream.str(),
OOMPH_CURRENT_FUNCTION,
OOMPH_EXCEPTION_LOCATION);
}
problem.restart(*restart_file_pt);
}
if (CommandLineArgs::Argc==3)
{
nstep=3;
cout << "Only doing nstep steps for validation: " << nstep << std::endl;
}
DocInfo doc_info;
doc_info.set_directory("RESLT");
problem.unsteady_run(nstep,restarted,doc_info);
if (CommandLineArgs::Argc==3)
{
DocInfo restarted_doc_info;
restarted_doc_info.set_directory("RESLT_restarted");
restarted_doc_info.number()=0;
restart_file_pt=new ifstream("RESLT/restart2.dat");
restarted_problem.restart(*restart_file_pt);
unsigned nstep=2;
bool restarted=true;
restarted_problem.unsteady_run(nstep,restarted,restarted_doc_info);
}
}
///////////////////////////////////////////////////////////////////// ///////////////////////////////...
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.
void doc_solution(DocInfo &doc_info)
Doc the solution.
void write_trace_file_header()
Write header for trace file.
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 ...
void dump_it(ofstream &dump_file, DocInfo doc_info)
Dump problem data.
int main(int argc, char *argv[])
Driver for fsi ring test problem.
Namespace for physical parameters.
double ReSt
Reynolds x Strouhal number.
double Re
Reynolds number.
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)