tensioned_string.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-2024 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 function for a simple beam proble
27 
28 //OOMPH-LIB includes
29 #include "generic.h"
30 #include "beam.h"
31 #include "meshes/one_d_lagrangian_mesh.h"
32 
33 using namespace std;
34 
35 using namespace oomph;
36 
37 //========start_of_namespace========================
38 /// Namespace for physical parameters
39 //==================================================
41 {
42  /// Non-dimensional thickness
43  double H;
44 
45  /// 2nd Piola Kirchhoff pre-stress
46  double Sigma0;
47 
48  /// Pressure load
49  double P_ext;
50 
51  /// Load function: Apply a constant external pressure to the beam
52  void load(const Vector<double>& xi, const Vector<double> &x,
53  const Vector<double>& N, Vector<double>& load)
54  {
55  for(unsigned i=0;i<2;i++) {load[i] = -P_ext*N[i];}
56  }
57 
58 } // end of namespace
59 
60 //======start_of_problem_class==========================================
61 /// Beam problem object
62 //======================================================================
63 class ElasticBeamProblem : public Problem
64 {
65 public:
66 
67  /// Constructor: The arguments are the number of elements,
68  /// the length of domain
69  ElasticBeamProblem(const unsigned &n_elem, const double &length);
70 
71  /// Conduct a parameter study
72  void parameter_study();
73 
74  /// Return pointer to the mesh
75  OneDLagrangianMesh<HermiteBeamElement>* mesh_pt()
76  {return dynamic_cast<OneDLagrangianMesh<HermiteBeamElement>*>
77  (Problem::mesh_pt());}
78 
79  /// No actions need to be performed after a solve
81 
82  /// No actions need to be performed before a solve
84 
85 private:
86 
87  /// Pointer to the node whose displacement is documented
88  Node* Doc_node_pt;
89 
90  /// Length of domain (in terms of the Lagrangian coordinates)
91  double Length;
92 
93  /// Pointer to geometric object that represents the beam's undeformed shape
94  GeomObject* Undef_beam_pt;
95 
96 }; // end of problem class
97 
98 
99 //=============start_of_constructor=====================================
100 /// Constructor for elastic beam problem
101 //======================================================================
103  const double &length) : Length(length)
104 {
105  // Set the undeformed beam to be a straight line at y=0
106  Undef_beam_pt=new StraightLine(0.0);
107 
108  // Create the (Lagrangian!) mesh, using the geometric object
109  // Undef_beam_pt to specify the initial (Eulerian) position of the
110  // nodes.
111  Problem::mesh_pt() =
112  new OneDLagrangianMesh<HermiteBeamElement>(n_elem,length,Undef_beam_pt);
113 
114  // Set the boundary conditions: Each end of the beam is fixed in space
115  // Loop over the boundaries (ends of the beam)
116  for(unsigned b=0;b<2;b++)
117  {
118  // Pin displacements in both x and y directions
119  // [Note: The mesh_pt() function has been overloaded
120  // to return a pointer to the actual mesh, rather than
121  // a pointer to the Mesh base class. The current mesh is derived
122  // from the SolidMesh class. In such meshes, all access functions
123  // to the nodes, such as boundary_node_pt(...), are overloaded
124  // to return pointers to SolidNodes (whose position can be
125  // pinned) rather than "normal" Nodes.]
126  mesh_pt()->boundary_node_pt(b,0)->pin_position(0);
127  mesh_pt()->boundary_node_pt(b,0)->pin_position(1);
128  }
129 
130  //Find number of elements in the mesh
131  unsigned n_element = mesh_pt()->nelement();
132 
133  //Loop over the elements to set physical parameters etc.
134  for(unsigned e=0;e<n_element;e++)
135  {
136  // Upcast to the specific element type
137  HermiteBeamElement *elem_pt =
138  dynamic_cast<HermiteBeamElement*>(mesh_pt()->element_pt(e));
139 
140  // Set physical parameters for each element:
141  elem_pt->sigma0_pt() = &Global_Physical_Variables::Sigma0;
142  elem_pt->h_pt() = &Global_Physical_Variables::H;
143 
144  // Set the load Vector for each element
145  elem_pt->load_vector_fct_pt() = &Global_Physical_Variables::load;
146 
147  // Set the undeformed shape for each element
148  elem_pt->undeformed_beam_pt() = Undef_beam_pt;
149  } // end of loop over elements
150 
151  // Choose node at which displacement is documented (halfway along -- provided
152  // we have an odd number of nodes; complain if this is not the
153  // case because the comparison with the exact solution will be wrong
154  // otherwise!)
155  unsigned n_nod=mesh_pt()->nnode();
156  if (n_nod%2!=1)
157  {
158  cout << "Warning: Even number of nodes " << n_nod << std::endl;
159  cout << "Comparison with exact solution will be misleading..." << std::endl;
160  }
161  Doc_node_pt=mesh_pt()->node_pt((n_nod+1)/2-1);
162 
163  // Assign the global and local equation numbers
164  cout << "# of dofs " << assign_eqn_numbers() << std::endl;
165 
166 } // end of constructor
167 
168 
169 //=======start_of_parameter_study==========================================
170 /// Solver loop to perform parameter study
171 //=========================================================================
173 {
174  // Over-ride the default maximum value for the residuals
175  Problem::Max_residuals = 1.0e10;
176 
177  // Set the increments in control parameters
178  double pext_increment = 0.001;
179 
180  // Set initial values for control parameters
181  Global_Physical_Variables::P_ext = 0.0 - pext_increment;
182 
183  // Create label for output
184  DocInfo doc_info;
185 
186  // Set output directory -- this function checks if the output
187  // directory exists and issues a warning if it doesn't.
188  doc_info.set_directory("RESLT");
189 
190  // Open a trace file
191  ofstream trace("RESLT/trace_beam.dat");
192 
193  // Write a header for the trace file
194  trace <<
195  "VARIABLES=\"p_e_x_t\",\"d\"" <<
196  ", \"p_e_x_t_(_e_x_a_c_t_)\"" << std::endl;
197 
198  // Output file stream used for writing results
199  ofstream file;
200  // String used for the filename
201  char filename[100];
202 
203  // Loop over parameter increments
204  unsigned nstep=10;
205  for(unsigned i=1;i<=nstep;i++)
206  {
207  // Increment pressure
208  Global_Physical_Variables::P_ext += pext_increment;
209 
210  // Solve the system
211  newton_solve();
212 
213  // Calculate exact solution for `string under tension' (applicable for
214  // small wall thickness and pinned ends)
215 
216  // The tangent of the angle beta
217  double tanbeta =-2.0*Doc_node_pt->x(1)/Length;
218 
219  double exact_pressure = 0.0;
220  //If the beam has deformed, calculate the pressure required
221  if(tanbeta!=0)
222  {
223 
224  //Calculate the opening angle alpha
225  double alpha = 2.0*atan(2.0*tanbeta/(1.0-tanbeta*tanbeta));
226 
227  // Jump back onto the main branch if alpha>180 degrees
228  if (alpha<0) alpha+=2.0*MathematicalConstants::Pi;
229 
230  // Green strain:
231  double gamma=0.5*(0.25*alpha*alpha/(sin(0.5*alpha)*sin(0.5*alpha))-1.0);
232 
233  //Calculate the exact pressure
234  exact_pressure=Global_Physical_Variables::H*
236  }
237 
238  // Document the solution
239  sprintf(filename,"RESLT/beam%i.dat",i);
240  file.open(filename);
241  mesh_pt()->output(file,5);
242  file.close();
243 
244  // Write trace file: Pressure, displacement and exact solution
245  // (for string under tension)
246  trace << Global_Physical_Variables::P_ext << " "
247  << abs(Doc_node_pt->x(1))
248  << " " << exact_pressure
249  << std::endl;
250  }
251 
252 } // end of parameter study
253 
254 //========start_of_main================================================
255 /// Driver for beam (string under tension) test problem
256 //=====================================================================
257 int main()
258 {
259 
260  // Set the non-dimensional thickness
262 
263  // Set the 2nd Piola Kirchhoff prestress
265 
266  // Set the length of domain
267  double L = 10.0;
268 
269  // Number of elements (choose an even number if you want the control point
270  // to be located at the centre of the beam)
271  unsigned n_element = 10;
272 
273  // Construst the problem
274  ElasticBeamProblem problem(n_element,L);
275 
276  // Check that we're ready to go:
277  cout << "\n\n\nProblem self-test ";
278  if (problem.self_test()==0)
279  {
280  cout << "passed: Problem can be solved." << std::endl;
281  }
282  else
283  {
284  throw OomphLibError("Self test failed",
285  OOMPH_CURRENT_FUNCTION,
286  OOMPH_EXCEPTION_LOCATION);
287  }
288 
289  // Conduct parameter study
290  problem.parameter_study();
291 
292 } // end of main
293 
Beam problem object.
GeomObject * Undef_beam_pt
Pointer to geometric object that represents the beam's undeformed shape.
ElasticBeamProblem(const unsigned &n_elem, const double &length)
Constructor: The arguments are the number of elements, the length of domain.
void parameter_study()
Conduct a parameter study.
void actions_after_newton_solve()
No actions need to be performed after a solve.
void actions_before_newton_solve()
No actions need to be performed before a solve.
Node * Doc_node_pt
Pointer to the node whose displacement is documented.
double Length
Length of domain (in terms of the Lagrangian coordinates)
OneDLagrangianMesh< HermiteBeamElement > * mesh_pt()
Return pointer to the mesh.
Namespace for physical parameters.
double P_ext
Pressure load.
void load(const Vector< double > &xi, const Vector< double > &x, const Vector< double > &N, Vector< double > &load)
Load function: Apply a constant external pressure to the beam.
double Sigma0
2nd Piola Kirchhoff pre-stress
double H
Non-dimensional thickness.
int main()
Driver for beam (string under tension) test problem.