33 #include "constitutive.h"
36 #include "meshes/rectangular_quadmesh.h"
40 using namespace oomph;
53 bool operator()(FiniteElement*
const& el1_pt, FiniteElement*
const& el2_pt)
56 return el1_pt->node_pt(0)->x(0) < el2_pt->node_pt(0)->x(0);
80 BrokenCopy::broken_copy(
"WarpedLine");
86 BrokenCopy::broken_assign(
"WarpedLine");
94 void position(
const Vector<double>& zeta, Vector<double>& r)
const
97 r[0] = zeta[0]+5.0*Ampl*zeta[0]*(zeta[0]-1.0)*(zeta[0]-0.7);
98 r[1] = 1.0+Ampl*0.5*(1.0-cos(2.0*MathematicalConstants::Pi*zeta[0]));
104 void position(
const unsigned& t,
const Vector<double>& zeta,
105 Vector<double>& r)
const
157 template<
class ELEMENT>
174 {
return Solid_mesh_pt;}
177 void actions_before_adapt();
180 void actions_after_adapt();
189 void create_lagrange_multiplier_elements();
193 void delete_lagrange_multiplier_elements();
210 template<
class ELEMENT>
229 solid_mesh_pt() =
new ElasticRefineableRectangularQuadMesh<ELEMENT>(
233 solid_mesh_pt()->spatial_error_estimator_pt()=
new Z2ErrorEstimator;
237 unsigned n_element =solid_mesh_pt()->nelement();
238 for(
unsigned i=0;i<n_element;i++)
241 ELEMENT *el_pt =
dynamic_cast<ELEMENT*
>(solid_mesh_pt()->element_pt(i));
244 el_pt->constitutive_law_pt()=&Global_Physical_Variables::Constitutive_law;
248 solid_mesh_pt()->refine_uniformly();
252 Lagrange_multiplier_mesh_pt=
new SolidMesh;
253 create_lagrange_multiplier_elements();
256 add_sub_mesh(solid_mesh_pt());
259 add_sub_mesh(Lagrange_multiplier_mesh_pt);
265 for (
unsigned b=0;b<4;b++)
269 unsigned n_side = solid_mesh_pt()->nboundary_node(b);
272 for(
unsigned i=0;i<n_side;i++)
274 solid_mesh_pt()->boundary_node_pt(b,i)->pin_position(0);
275 solid_mesh_pt()->boundary_node_pt(b,i)->pin_position(1);
281 PVDEquationsBase<2>::pin_redundant_nodal_solid_pressures(
282 solid_mesh_pt()->element_pt());
285 cout <<
"Number of dofs: " << assign_eqn_numbers() << std::endl;
288 Doc_info.set_directory(
"RESLT");
297 template<
class ELEMENT>
301 delete_lagrange_multiplier_elements();
304 rebuild_global_mesh();
314 template<
class ELEMENT>
320 create_lagrange_multiplier_elements();
323 rebuild_global_mesh();
326 PVDEquationsBase<2>::pin_redundant_nodal_solid_pressures(
327 solid_mesh_pt()->element_pt());
336 template<
class ELEMENT>
344 unsigned n_element = solid_mesh_pt()->nboundary_element(b);
347 for(
unsigned e=0;e<n_element;e++)
350 ELEMENT* bulk_elem_pt =
dynamic_cast<ELEMENT*
>(
351 solid_mesh_pt()->boundary_element_pt(b,e));
354 int face_index = solid_mesh_pt()->face_index_at_boundary(b,e);
357 Lagrange_multiplier_mesh_pt->add_element_pt(
358 new ImposeDisplacementByLagrangeMultiplierElement<ELEMENT>(
359 bulk_elem_pt,face_index));
365 n_element=Lagrange_multiplier_mesh_pt->nelement();
366 for(
unsigned i=0;i<n_element;i++)
369 ImposeDisplacementByLagrangeMultiplierElement<ELEMENT> *el_pt =
370 dynamic_cast<ImposeDisplacementByLagrangeMultiplierElement<ELEMENT>*
>
371 (Lagrange_multiplier_mesh_pt->element_pt(i));
376 el_pt->set_boundary_shape_geom_object_pt(
380 unsigned nnod=el_pt->nnode();
381 for (
unsigned j=0;j<nnod;j++)
383 Node* nod_pt = el_pt->node_pt(j);
386 if ((nod_pt->is_on_boundary(1))||(nod_pt->is_on_boundary(3)))
390 unsigned n_bulk_value=el_pt->nbulk_value(j);
393 unsigned nval=nod_pt->nvalue();
394 for (
unsigned j=n_bulk_value;j<nval;j++)
411 template<
class ELEMENT>
415 unsigned n_element = Lagrange_multiplier_mesh_pt->nelement();
418 for(
unsigned e=0;e<n_element;e++)
421 delete Lagrange_multiplier_mesh_pt->element_pt(e);
425 Lagrange_multiplier_mesh_pt->flush_element_and_node_storage();
434 template<
class ELEMENT>
447 sprintf(filename,
"%s/soln%i.dat",Doc_info.directory().c_str(),
449 some_file.open(filename);
450 solid_mesh_pt()->output(some_file,n_plot);
455 sprintf(filename,
"%s/lagr%i.dat",Doc_info.directory().c_str(),
457 some_file.open(filename);
461 std::vector<FiniteElement*> el_pt;
462 unsigned nelem=Lagrange_multiplier_mesh_pt->nelement();
463 for (
unsigned e=0;e<nelem;e++)
465 el_pt.push_back(Lagrange_multiplier_mesh_pt->finite_element_pt(e));
468 for (
unsigned e=0;e<nelem;e++)
470 el_pt[e]->output(some_file);
494 unsigned max_adapt=1;
498 for(
unsigned i=0;i<nstep;i++)
505 problem.newton_solve(max_adapt);
Function-type-object to compare finite elements based on their x coordinate.
bool operator()(FiniteElement *const &el1_pt, FiniteElement *const &el2_pt) const
Comparison. Is x coordinate of el1_pt less than that of el2_pt?
Problem class for deformation of elastic block by prescribed boundary motion.
void delete_lagrange_multiplier_elements()
Delete elements that enforce prescribed boundary motion by Lagrange multiplilers.
void actions_before_adapt()
Actions before adapt: Wipe the mesh of Lagrange multiplier elements.
void actions_after_newton_solve()
Update function (empty)
ElasticRefineableRectangularQuadMesh< ELEMENT > * Solid_mesh_pt
Pointer to solid mesh.
void actions_after_adapt()
Actions after adapt: Rebuild the mesh of Lagrange multiplier elements.
SolidMesh * Lagrange_multiplier_mesh_pt
Pointers to meshes of Lagrange multiplier elements.
void doc_solution()
Doc the solution.
PrescribedBoundaryDisplacementProblem()
Constructor:
ElasticRefineableRectangularQuadMesh< ELEMENT > *& solid_mesh_pt()
Access function for the solid mesh.
void actions_before_newton_solve()
Update function (empty)
void create_lagrange_multiplier_elements()
Create elements that enforce prescribed boundary motion by Lagrange multiplilers.
DocInfo Doc_info
DocInfo object for output.
void position(const unsigned &t, const Vector< double > &zeta, Vector< double > &r) const
Parametrised position on object: r(zeta). Evaluated at previous timestep. t=0: current time; t>0: pre...
~WarpedLine()
Empty Destructor.
WarpedLine(const WarpedLine &dummy)
Broken copy constructor.
WarpedLine(const double &l)
Constructor: Specify amplitude of deflection from straight horizontal line.
unsigned ngeom_data() const
How many items of Data does the shape of the object depend on? None.
void position(const Vector< double > &zeta, Vector< double > &r) const
Position vector at Lagrangian coordinate zeta.
double Ampl
Amplitude of perturbation.
void operator=(const WarpedLine &)
Broken assignment operator.
double & ampl()
Access to amplitude.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
double Nu
Poisson's ratio.
WarpedLine Boundary_geom_object(0.0)
GeomObject specifying the shape of the boundary: Initially it's flat.