tensioned_string2.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-2026 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
33using namespace std;
34
35using namespace oomph;
36
37//=========================================================================
38/// Steady, straight 1D line in 2D space with stretch_ratio
39/// \f[ x = \zeta stretch_ratio \f]
40/// \f[ y = 0 \f]
41//=========================================================================
42class NewStraightLine: public GeomObject
43{
44public:
45
46 /// Constructor derives from GeomObject(1, 2)
47 /// Pass stretch_ratio to this class (defaults to 1 in which case
48 /// the Lagrangian coordinate is the arclength)
49 NewStraightLine(const double& stretch_ratio=1.0) : GeomObject(1, 2)
50 {
51 Stretch_ratio = stretch_ratio;
52 }
53
54 /// Broken copy constructor
55 NewStraightLine(const NewStraightLine& dummy) = delete;
56
57 /// Broken assignment operator
58 void operator=(const NewStraightLine&) = delete;
59
60 /// Position Vector at Lagrangian coordinate zeta
61 void position(const Vector<double>& zeta, Vector<double>& r) const
62 {
63 r[0] = zeta[0] * Stretch_ratio;
64 r[1] = 0.0;
65 }
66
67
68 /// Derivative of position Vector w.r.t. to coordinates:
69 /// \f$ \frac{dR_i}{d \zeta_\alpha}\f$ = drdzeta(alpha,i).
70 /// Evaluated at current time.
71 virtual void dposition(const Vector<double>& zeta,
72 DenseMatrix<double>& drdzeta) const
73 {
74 // Tangent vector
75 drdzeta(0, 0) = Stretch_ratio;
76 drdzeta(0, 1) = 0.0;
77 }
78
79
80 /// 2nd derivative of position Vector w.r.t. to coordinates:
81 /// \f$ \frac{d^2R_i}{d \zeta_\alpha d \zeta_\beta}\f$ =
82 /// ddrdzeta(alpha,beta,i). Evaluated at current time.
83 virtual void d2position(const Vector<double>& zeta,
84 RankThreeTensor<double>& ddrdzeta) const
85 {
86 // Derivative of tangent vector
87 ddrdzeta(0, 0, 0) = 0.0;
88 ddrdzeta(0, 0, 1) = 0.0;
89 }
90
91
92 /// Posn Vector and its 1st & 2nd derivatives
93 /// w.r.t. to coordinates:
94 /// \f$ \frac{dR_i}{d \zeta_\alpha}\f$ = drdzeta(alpha,i).
95 /// \f$ \frac{d^2R_i}{d \zeta_\alpha d \zeta_\beta}\f$ =
96 /// ddrdzeta(alpha,beta,i).
97 /// Evaluated at current time.
98 virtual void d2position(const Vector<double>& zeta,
99 Vector<double>& r,
100 DenseMatrix<double>& drdzeta,
101 RankThreeTensor<double>& ddrdzeta) const
102 {
103 // Position Vector
104 r[0] = zeta[0] * Stretch_ratio;
105 r[1] = 0.0;
106
107 // Tangent vector
108 drdzeta(0, 0) = Stretch_ratio;
109 drdzeta(0, 1) = 0.0;
110
111 // Derivative of tangent vector
112 ddrdzeta(0, 0, 0) = 0.0;
113 ddrdzeta(0, 0, 1) = 0.0;
114 }
115
116private:
117
118 /// Define the length of the beam
120
121};
122
123
124
125
126//========start_of_namespace========================
127/// Namespace for physical parameters
128//==================================================
130{
131 /// Non-dimensional thickness
132 double H;
133
134 /// 2nd Piola Kirchhoff pre-stress
135 double Sigma0;
136
137 /// Pressure load
138 double P_ext;
139
140 /// Load function: Apply a constant external pressure to the beam
141 void load(const Vector<double>& xi, const Vector<double> &x,
142 const Vector<double>& N, Vector<double>& load)
143 {
144 for(unsigned i=0;i<2;i++) {load[i] = -P_ext*N[i];}
145 }
146
147} // end of namespace
148
149//======start_of_problem_class==========================================
150/// Beam problem object
151//======================================================================
152class ElasticBeamProblem : public Problem
153{
154public:
155
156 /// Constructor: The arguments are the number of elements,
157 /// the length of domain
158 ElasticBeamProblem(const unsigned &n_elem, const double &length,
159 const double& stretch_ratio);
160
161 /// Conduct a parameter study
163
164 /// Return pointer to the mesh
165 OneDLagrangianMesh<HermiteBeamElement>* mesh_pt()
166 {return dynamic_cast<OneDLagrangianMesh<HermiteBeamElement>*>
167 (Problem::mesh_pt());}
168
169 /// No actions need to be performed after a solve
171
172 /// No actions need to be performed before a solve
174
175private:
176
177 /// Pointer to the node whose displacement is documented
178 Node* Doc_node_pt;
179
180 /// Length of domain (in terms of the Lagrangian coordinates)
181 double Length;
182
183 /// Pointer to geometric object that represents the beam's undeformed shape
184 GeomObject* Undef_beam_pt;
185
186}; // end of problem class
187
188
189//=============start_of_constructor=====================================
190/// Constructor for elastic beam problem
191//======================================================================
193 const double &length,
194 const double& stretch_ratio) : Length(length)
195{
196 // Set the undeformed beam to be a straight line at y=0
197 Undef_beam_pt=new NewStraightLine(stretch_ratio);
198
199 // Create the (Lagrangian!) mesh, using the geometric object
200 // Undef_beam_pt to specify the initial (Eulerian) position of the
201 // nodes.
202
203 // Compensate for stretched coordinate by adjusting nominal length (in terms
204 // of the coordinate)
205 double compensated_length=length/stretch_ratio;
206 Problem::mesh_pt() =
207 new OneDLagrangianMesh<HermiteBeamElement>(n_elem,compensated_length,Undef_beam_pt);
208
209 // Set the boundary conditions: Each end of the beam is fixed in space
210 // Loop over the boundaries (ends of the beam)
211 for(unsigned b=0;b<2;b++)
212 {
213 // Pin displacements in both x and y directions
214 // [Note: The mesh_pt() function has been overloaded
215 // to return a pointer to the actual mesh, rather than
216 // a pointer to the Mesh base class. The current mesh is derived
217 // from the SolidMesh class. In such meshes, all access functions
218 // to the nodes, such as boundary_node_pt(...), are overloaded
219 // to return pointers to SolidNodes (whose position can be
220 // pinned) rather than "normal" Nodes.]
221 mesh_pt()->boundary_node_pt(b,0)->pin_position(0);
222 mesh_pt()->boundary_node_pt(b,0)->pin_position(1);
223 }
224
225 //Find number of elements in the mesh
226 unsigned n_element = mesh_pt()->nelement();
227
228 //Loop over the elements to set physical parameters etc.
229 for(unsigned e=0;e<n_element;e++)
230 {
231 // Upcast to the specific element type
232 HermiteBeamElement *elem_pt =
233 dynamic_cast<HermiteBeamElement*>(mesh_pt()->element_pt(e));
234
235 // Set physical parameters for each element:
236 elem_pt->sigma0_pt() = &Global_Physical_Variables::Sigma0;
237 elem_pt->h_pt() = &Global_Physical_Variables::H;
238
239 // Set the load Vector for each element
240 elem_pt->load_vector_fct_pt() = &Global_Physical_Variables::load;
241
242 // Set the undeformed shape for each element
243 elem_pt->undeformed_beam_pt() = Undef_beam_pt;
244 } // end of loop over elements
245
246 // Choose node at which displacement is documented (halfway along -- provided
247 // we have an odd number of nodes; complain if this is not the
248 // case because the comparison with the exact solution will be wrong
249 // otherwise!)
250 unsigned n_nod=mesh_pt()->nnode();
251 if (n_nod%2!=1)
252 {
253 cout << "Warning: Even number of nodes " << n_nod << std::endl;
254 cout << "Comparison with exact solution will be misleading..." << std::endl;
255 }
256 Doc_node_pt=mesh_pt()->node_pt((n_nod+1)/2-1);
257
258 // Assign the global and local equation numbers
259 cout << "# of dofs " << assign_eqn_numbers() << std::endl;
260
261} // end of constructor
262
263
264//=======start_of_parameter_study==========================================
265/// Solver loop to perform parameter study
266//=========================================================================
268{
269 // Over-ride the default maximum value for the residuals
270 Problem::Max_residuals = 1.0e10;
271
272 // Set the increments in control parameters
273 double pext_increment = 0.001;
274
275 // Set initial values for control parameters
276 Global_Physical_Variables::P_ext = 0.0 - pext_increment;
277
278 // Create label for output
279 DocInfo doc_info;
280
281 // Set output directory -- this function checks if the output
282 // directory exists and issues a warning if it doesn't.
283 doc_info.set_directory("RESLT");
284
285 // Open a trace file
286 ofstream trace("RESLT/trace_beam.dat");
287
288 // Write a header for the trace file
289 trace <<
290 "VARIABLES=\"p_e_x_t\",\"d\"" <<
291 ", \"p_e_x_t_(_e_x_a_c_t_)\"" << std::endl;
292
293 // Output file stream used for writing results
294 ofstream file;
295 // String used for the filename
296 char filename[100];
297
298 // Loop over parameter increments
299 unsigned nstep=10;
300 for(unsigned i=1;i<=nstep;i++)
301 {
302 // Increment pressure
303 Global_Physical_Variables::P_ext += pext_increment;
304
305 // Solve the system
306 newton_solve();
307
308 // Calculate exact solution for `string under tension' (applicable for
309 // small wall thickness and pinned ends)
310
311 // The tangent of the angle beta
312 double tanbeta =-2.0*Doc_node_pt->x(1)/Length;
313
314 double exact_pressure = 0.0;
315 //If the beam has deformed, calculate the pressure required
316 if(tanbeta!=0)
317 {
318
319 //Calculate the opening angle alpha
320 double alpha = 2.0*atan(2.0*tanbeta/(1.0-tanbeta*tanbeta));
321
322 // Jump back onto the main branch if alpha>180 degrees
323 if (alpha<0) alpha+=2.0*MathematicalConstants::Pi;
324
325 // Green strain:
326 double gamma=0.5*(0.25*alpha*alpha/(sin(0.5*alpha)*sin(0.5*alpha))-1.0);
327
328 //Calculate the exact pressure
329 exact_pressure=Global_Physical_Variables::H*
331 }
332
333 // Document the solution
334 sprintf(filename,"RESLT/beam%i.dat",i);
335 file.open(filename);
336 mesh_pt()->output(file,5);
337 file.close();
338
339 // Write trace file: Pressure, displacement and exact solution
340 // (for string under tension)
341 trace << Global_Physical_Variables::P_ext << " "
342 << abs(Doc_node_pt->x(1))
343 << " " << exact_pressure
344 << std::endl;
345 }
346
347} // end of parameter study
348
349//========start_of_main================================================
350/// Driver for beam (string under tension) test problem.
351/// Redo with non-arclength parametrisation.
352//=====================================================================
353int main(int argc, char** argv)
354{
355
356 // Store command line arguments
357 CommandLineArgs::setup(argc, argv);
358
359 // Stretch ratio
360 double stretch_ratio=1.0;
361 CommandLineArgs::specify_command_line_flag(
362 "--stretch_ratio",&stretch_ratio);
363
364 // Set the 2nd Piola Kirchhoff prestress
366 CommandLineArgs::specify_command_line_flag(
368
369 // Parse command line
370 CommandLineArgs::parse_and_assign();
371
372 // Doc what has actually been specified on the command line
373 CommandLineArgs::doc_specified_flags();
374
375 // Set the non-dimensional thickness
377
378 // Set the length of domain
379 double L = 10.0;
380
381 // Number of elements (choose an even number if you want the control point
382 // to be located at the centre of the beam)
383 unsigned n_element = 10;
384
385 // Construst the problem
386 ElasticBeamProblem problem(n_element,L,stretch_ratio);
387
388 // Conduct parameter study
389 problem.parameter_study();
390
391} // end of main
392
Beam problem object.
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.
GeomObject * Undef_beam_pt
Pointer to geometric object that represents the beam's undeformed shape.
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.
OneDLagrangianMesh< HermiteBeamElement > * mesh_pt()
Return pointer to the mesh.
double Length
Length of domain (in terms of the Lagrangian coordinates)
Steady, straight 1D line in 2D space with stretch_ratio.
void operator=(const NewStraightLine &)=delete
Broken assignment operator.
NewStraightLine(const NewStraightLine &dummy)=delete
Broken copy constructor.
virtual void d2position(const Vector< double > &zeta, Vector< double > &r, DenseMatrix< double > &drdzeta, RankThreeTensor< double > &ddrdzeta) const
Posn Vector and its 1st & 2nd derivatives w.r.t. to coordinates: = drdzeta(alpha,...
virtual void d2position(const Vector< double > &zeta, RankThreeTensor< double > &ddrdzeta) const
2nd derivative of position Vector w.r.t. to coordinates: = ddrdzeta(alpha,beta,i)....
NewStraightLine(const double &stretch_ratio=1.0)
Constructor derives from GeomObject(1, 2) Pass stretch_ratio to this class (defaults to 1 in which ca...
virtual void dposition(const Vector< double > &zeta, DenseMatrix< double > &drdzeta) const
Derivative of position Vector w.r.t. to coordinates: = drdzeta(alpha,i). Evaluated at current time.
double Stretch_ratio
Define the length of the beam.
void position(const Vector< double > &zeta, Vector< double > &r) const
Position Vector at Lagrangian coordinate zeta.
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.