31#include "oomph_crbond_bessel.h"
34#include "meshes/quarter_circle_sector_mesh.h"
93 TimeStepper* time_stepper_pt);
100 void position(
const Vector<double>& xi, Vector<double>& r)
const;
103 void veloc(
const Vector<double>& xi, Vector<double>&
veloc);
107 void accel(
const Vector<double>& xi, Vector<double>&
accel);
112 Vector<double>& drdt)
132 std::ostringstream error_message;
133 error_message << j <<
"th derivative not implemented\n";
135 throw OomphLibError(error_message.str(),
136 OOMPH_CURRENT_FUNCTION,
137 OOMPH_EXCEPTION_LOCATION);
146 const Vector<double>& omega,
147 Vector<double>& residual);
175 TimeStepper* time_stepper_pt) :
176 GeomObject(2,2,time_stepper_pt), Ampl(ampl), Nu(nu)
179 Vector<double> param(1);
183 Vector<double> omega(1);
194 T=2.0*MathematicalConstants::Pi/
Omega;
196 std::cout <<
"Period of oscillation: " <<
T << std::endl;
207 Vector<double>& r)
const
210 double time=Geom_object_time_stepper_pt->time_pt()->time();
213 double lagr_radius=sqrt( pow(xi[0],2) + pow(xi[1],2) );
215 if (lagr_radius<1.0e-12)
224 double j0,j1,y0,y1,j0p,j1p,y0p,y1p;
225 CRBond_Bessel::bessjy01a(
Omega*lagr_radius,j0,j1,y0,y1,j0p,j1p,y0p,y1p);
228 double u=
Ampl*j1*sin(2.0*MathematicalConstants::Pi*time/
T);
231 r[0]=(xi[0]+xi[0]/lagr_radius*u);
232 r[1]=(xi[1]+xi[1]/lagr_radius*u);
245 Vector<double>& veloc)
248 double time=Geom_object_time_stepper_pt->time_pt()->time();
251 double lagr_radius=sqrt(pow(xi[0],2)+pow(xi[1],2));
253 if (lagr_radius<1.0e-12)
262 double j0,j1,y0,y1,j0p,j1p,y0p,y1p;
263 CRBond_Bessel::bessjy01a(
Omega*lagr_radius,j0,j1,y0,y1,j0p,j1p,y0p,y1p);
266 double udot=2.0*MathematicalConstants::Pi/
T*
Ampl*j1*
267 cos(2.0*MathematicalConstants::Pi*time/
T);
270 veloc[0]=(xi[0]/lagr_radius*udot);
271 veloc[1]=(xi[1]/lagr_radius*udot);
284 Vector<double>& accel)
287 double time=Geom_object_time_stepper_pt->time_pt()->time();
290 double lagr_radius=sqrt(pow(xi[0],2)+pow(xi[1],2));
293 if (lagr_radius<1.0e-12)
302 double j0,j1,y0,y1,j0p,j1p,y0p,y1p;
303 CRBond_Bessel::bessjy01a(
Omega*lagr_radius,j0,j1,y0,y1,j0p,j1p,y0p,y1p);
306 double udotdot=-pow(2.0*MathematicalConstants::Pi/
T,2)*
Ampl*j1*
307 sin(2.0*MathematicalConstants::Pi*time/
T);
310 accel[0]=(xi[0]/lagr_radius*udotdot);
311 accel[1]=(xi[1]/lagr_radius*udotdot);
324 const Vector<double>& param,
const Vector<double>& omega,
325 Vector<double>& residual)
334 double j0,j1,y0,y1,j0p,j1p,y0p,y1p;
335 CRBond_Bessel::bessjy01a(arg,j0,j1,y0,y1,j0p,j1p,y0p,y1p);
338 residual[0]=nu/(1.0-2.0*nu)*(j1+(j0-j1/omega[0])*omega[0])+
339 (j0-j1/omega[0])*omega[0];
358template <
class ELEMENT>
360 public virtual RefineableQuarterCircleSectorMesh<ELEMENT>,
361 public virtual SolidMesh
385 "Element needs to be derived from SolidFiniteElement\n",
413template<
class ELEMENT>
461template<
class ELEMENT>
491 unsigned nnod0=mesh_pt()->nboundary_node(0);
492 unsigned nnod1=mesh_pt()->nboundary_node(1);
497 Trace_node_pt[
j]=mesh_pt()->boundary_node_pt(0,
j);
503 Trace_node_pt[
j+
nnod0]=mesh_pt()->boundary_node_pt(1,
j);
509 unsigned n_hor = mesh_pt()->nboundary_node(0);
512 mesh_pt()->boundary_node_pt(0,
i)->pin_position(1);
516 unsigned n_vert = mesh_pt()->nboundary_node(2);
520 mesh_pt()->boundary_node_pt(2,
i)->pin_position(0);
526 unsigned n_element =mesh_pt()->nelement();
533 el_pt->constitutive_law_pt() =
541 mesh_pt()->refine_uniformly();
553template<
class ELEMENT>
588 Trace_file <<
time_pt()->time() <<
" "
595 Trace_file <<
sqrt(
pow(Trace_node_pt[
j]->
x(0),2)+
596 pow(Trace_node_pt[
j]->
x(1),2)) <<
" ";
598 Trace_file << std::endl;
607 unsigned nelem=mesh_pt()->nboundary_element(0);
640 (mesh_pt()->boundary_element_pt(0,
e));
659 <<
dxdt[0] << std::endl;
676 unsigned nelem=mesh_pt()->nelement();
709 (mesh_pt()->element_pt(
e));
728 <<
x[0]-
xi[0] <<
" " <<
x[1]-
xi[1] <<
" "
729 <<
dxdt[0] <<
" " <<
dxdt[1] << std::endl;
752template<
class ELEMENT>
789 SolidMesh::Solid_IC_problem.set_newmark_initial_condition_consistently(
Axisymmetrially oscillating disk with displacement field according to linear elasticity.
void dposition_dt(const Vector< double > &xi, const unsigned &j, Vector< double > &drdt)
Parametrised j-th time-derivative on object at current time: .
AxisymOscillatingDisk(const double &l, const double &nu, TimeStepper *time_stepper_pt)
Constructor: 2 Lagrangian coordinate, 2 Eulerian coords. Pass amplitude of oscillation,...
double T
Period of oscillation.
void veloc(const Vector< double > &xi, Vector< double > &veloc)
Parametrised velocity on object at current time: veloc = d r(xi)/dt.
void accel(const Vector< double > &xi, Vector< double > &accel)
Parametrised acceleration on object at current time: accel = d^2 r(xi)/dt^2.
double Ampl
Amplitude of oscillation.
double Omega
Eigenfrequency.
void position(const Vector< double > &xi, Vector< double > &r) const
Position vector at Lagrangian coordinate xi at present time.
~AxisymOscillatingDisk()
Destructor (empty)
double Nu
Poisson ratio nu.
static void residual_for_dispersion(const Vector< double > ¶m, const Vector< double > &omega, Vector< double > &residual)
Residual of dispersion relation for use in black-box Newton method which requires global (or static) ...
Problem class to simulate small-amplitude oscillations of a circular disk.
Vector< Node * > Trace_node_pt
Vector of pointers to nodes whose position we're tracing.
AxisymOscillatingDisk * IC_geom_object_pt
Geometric object that specifies the initial conditions.
ofstream Trace_file
Trace file.
ElasticRefineableQuarterCircleSectorMesh< ELEMENT > * mesh_pt()
Access function for the solid mesh.
DiskOscillationProblem()
Constructor.
void actions_after_newton_solve()
Update function (empty)
void run(const unsigned &nstep)
Run the problem: Pass number of timesteps to be performed.
void doc_solution(DocInfo &doc_info)
Doc the solution.
void actions_before_newton_solve()
Update function (empty)
Elastic quarter circle sector mesh: We "upgrade" the RefineableQuarterCircleSectorMesh to become an S...
ElasticRefineableQuarterCircleSectorMesh(GeomObject *wall_pt, const double &xi_lo, const double &fract_mid, const double &xi_hi, TimeStepper *time_stepper_pt=&Mesh::Default_TimeStepper)
Constructor: Build mesh and copy Eulerian coords to Lagrangian ones so that the initial configuration...
int main(int argc, char *argv[])
Driver for disk oscillation problem.
double multiplier(const Vector< double > &xi)
Multiplier for inertia terms (needed for consistent assignment of initial conditions in Newmark schem...
ConstitutiveLaw * Constitutive_law_pt
Pointer to constitutive law.
double Nu
Poisson's ratio.
double Lambda_sq
Timescale ratio.