harmonic.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-2023 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 to solve the harmonic equation with homogeneous Dirichlet boundary
27 //conditions.
28 
29 // Generic oomph-lib routines
30 #include "generic.h"
31 
32 // Include the mesh
33 #include "meshes/one_d_mesh.h"
34 
35 using namespace std;
36 
37 using namespace oomph;
38 
39 
40 //===================================================================
41 /// Function-type-object to perform comparison of complex data types
42 /// Needed to sort the complex eigenvalues into order based on the
43 /// size of the real part.
44 //==================================================================
45 template <class T>
46 class ComplexLess
47 {
48 public:
49 
50  /// Comparison. Are the values identical or not?
51  bool operator()(const complex<T> &x, const complex<T> &y) const
52  {
53  return x.real() < y.real();
54  }
55 };
56 
57 
58 //=================================================================
59 /// A class for all elements that solve the simple one-dimensional
60 /// eigenvalue problem
61 /// \f[
62 /// \frac{\partial^2 u}{\partial x_i^2} + \lambda u = 0
63 /// \f]
64 /// These elements are very closely related to the Poisson
65 /// elements and could inherit from them. They are here developed
66 /// from scratch for pedagogical purposes.
67 /// This class contains the generic maths. Shape functions, geometric
68 /// mapping etc. must get implemented in derived class.
69 //================================================================
70 class HarmonicEquations : public virtual FiniteElement
71 {
72 
73 public:
74  /// Empty Constructor
75  HarmonicEquations() {}
76 
77  /// Access function: Eigenfunction value at local node n
78  /// Note that solving the eigenproblem does not assign values
79  /// to this storage space. It is used for output purposes only.
80  virtual inline double u(const unsigned& n) const
81  {return nodal_value(n,0);}
82 
83  /// Output the eigenfunction with default number of plot points
84  void output(ostream &outfile)
85  {
86  unsigned nplot=5;
87  output(outfile,nplot);
88  }
89 
90  /// Output FE representation of soln: x,y,u or x,y,z,u at
91  /// Nplot plot points
92  void output(ostream &outfile, const unsigned &nplot)
93  {
94  //Vector of local coordinates
95  Vector<double> s(1);
96 
97  // Tecplot header info
98  outfile << tecplot_zone_string(nplot);
99 
100  // Loop over plot points
101  unsigned num_plot_points=nplot_points(nplot);
102  for (unsigned iplot=0;iplot<num_plot_points;iplot++)
103  {
104  // Get local coordinates of plot point
105  get_s_plot(iplot,nplot,s);
106  //Output the coordinate and the eigenfunction
107  outfile << interpolated_x(s,0) << " " << interpolated_u(s) << std::endl;
108  }
109  // Write tecplot footer (e.g. FE connectivity lists)
110  write_tecplot_zone_footer(outfile,nplot);
111  }
112 
113  /// Assemble the contributions to the jacobian and mass
114  /// matrices
115  void fill_in_contribution_to_jacobian_and_mass_matrix(
116  Vector<double> &residuals,
117  DenseMatrix<double> &jacobian, DenseMatrix<double> &mass_matrix)
118  {
119  //Find out how many nodes there are
120  unsigned n_node = nnode();
121 
122  //Set up memory for the shape functions and their derivatives
123  Shape psi(n_node);
124  DShape dpsidx(n_node,1);
125 
126  //Set the number of integration points
127  unsigned n_intpt = integral_pt()->nweight();
128 
129  //Integers to store the local equation and unknown numbers
130  int local_eqn=0, local_unknown=0;
131 
132  //Loop over the integration points
133  for(unsigned ipt=0;ipt<n_intpt;ipt++)
134  {
135  //Get the integral weight
136  double w = integral_pt()->weight(ipt);
137 
138  //Call the derivatives of the shape and test functions
139  double J = dshape_eulerian_at_knot(ipt,psi,dpsidx);
140 
141  //Premultiply the weights and the Jacobian
142  double W = w*J;
143 
144  //Assemble the contributions to the mass matrix
145  //Loop over the test functions
146  for(unsigned l=0;l<n_node;l++)
147  {
148  //Get the local equation number
149  local_eqn = u_local_eqn(l);
150  /*IF it's not a boundary condition*/
151  if(local_eqn >= 0)
152  {
153  //Loop over the shape functions
154  for(unsigned l2=0;l2<n_node;l2++)
155  {
156  local_unknown = u_local_eqn(l2);
157  //If at a non-zero degree of freedom add in the entry
158  if(local_unknown >= 0)
159  {
160  jacobian(local_eqn,local_unknown) += dpsidx(l,0)*dpsidx(l2,0)*W;
161  mass_matrix(local_eqn, local_unknown) += psi(l)*psi(l2)*W;
162  }
163  }
164  }
165  }
166  }
167  } //end_of_fill_in_contribution_to_jacobian_and_mass_matrix
168 
169  /// Return FE representation of function value u(s) at local coordinate s
170  inline double interpolated_u(const Vector<double> &s) const
171  {
172  unsigned n_node = nnode();
173 
174  //Local shape functions
175  Shape psi(n_node);
176 
177  //Find values of basis function
178  this->shape(s,psi);
179 
180  //Initialise value of u
181  double interpolated_u = 0.0;
182 
183  //Loop over the local nodes and sum
184  for(unsigned l=0;l<n_node;l++) {interpolated_u+=u(l)*psi[l];}
185 
186  //Return the interpolated value of the eigenfunction
187  return(interpolated_u);
188  }
189 
190 protected:
191 
192  /// Shape/test functions and derivs w.r.t. to global coords at
193  /// local coord. s; return Jacobian of mapping
194  virtual double dshape_eulerian(const Vector<double> &s,
195  Shape &psi,
196  DShape &dpsidx) const=0;
197 
198  /// Shape/test functions and derivs w.r.t. to global coords at
199  /// integration point ipt; return Jacobian of mapping
200  virtual double dshape_eulerian_at_knot(const unsigned &ipt,
201  Shape &psi,
202  DShape &dpsidx) const=0;
203 
204  /// Access function that returns the local equation number
205  /// of the unknown in the problem. Default is to assume that it is the
206  /// first (only) value stored at the node.
207  virtual inline int u_local_eqn(const unsigned &n)
208  {return nodal_local_eqn(n,0);}
209 
210  private:
211 
212 };
213 
214 
215 
216 //======================================================================
217 /// QHarmonicElement<NNODE_1D> elements are 1D Elements with
218 /// NNODE_1D nodal points that are used to solve the Harmonic eigenvalue
219 /// Problem described by HarmonicEquations.
220 //======================================================================
221 template <unsigned NNODE_1D>
222 class QHarmonicElement : public virtual QElement<1,NNODE_1D>,
223  public HarmonicEquations
224 {
225 
226  public:
227 
228  /// Constructor: Call constructors for QElement and
229  /// Poisson equations
230  QHarmonicElement() : QElement<1,NNODE_1D>(), HarmonicEquations() {}
231 
232  /// Required # of `values' (pinned or dofs)
233  /// at node n
234  inline unsigned required_nvalue(const unsigned &n) const {return 1;}
235 
236  /// Output function overloaded from HarmonicEquations
237  void output(ostream &outfile)
238  {HarmonicEquations::output(outfile);}
239 
240  /// Output function overloaded from HarmonicEquations
241  void output(ostream &outfile, const unsigned &Nplot)
242  {HarmonicEquations::output(outfile,Nplot);}
243 
244 
245 protected:
246 
247 /// Shape, test functions & derivs. w.r.t. to global coords. Return Jacobian.
248  inline double dshape_eulerian(const Vector<double> &s,
249  Shape &psi,
250  DShape &dpsidx) const
251  {return QElement<1,NNODE_1D>::dshape_eulerian(s,psi,dpsidx);}
252 
253 
254  /// Shape, test functions & derivs. w.r.t. to global coords. at
255  /// integration point ipt. Return Jacobian.
256  inline double dshape_eulerian_at_knot(const unsigned& ipt,
257  Shape &psi,
258  DShape &dpsidx) const
259  {return QElement<1,NNODE_1D>::dshape_eulerian_at_knot(ipt,psi,dpsidx);}
260 
261 }; //end_of_QHarmonic_class_definition
262 
263 
264 //==start_of_problem_class============================================
265 /// 1D Harmonic problem in unit interval.
266 //====================================================================
267 template<class ELEMENT,class EIGEN_SOLVER>
268 class HarmonicProblem : public Problem
269 {
270 public:
271 
272  /// Constructor: Pass number of elements and pointer to source function
273  HarmonicProblem(const unsigned& n_element);
274 
275  /// Destructor (empty)
276  ~HarmonicProblem(){delete this->mesh_pt(); delete this->eigen_solver_pt();}
277 
278  /// Solve the problem
279  void solve(const unsigned &label);
280 
281  /// Doc the solution, pass the number of the case considered,
282  /// so that output files can be distinguished.
283  void doc_solution(const unsigned& label);
284 
285 }; // end of problem class
286 
287 
288 
289 //=====start_of_constructor===============================================
290 /// Constructor for 1D Harmonic problem in unit interval.
291 /// Discretise the 1D domain with n_element elements of type ELEMENT.
292 /// Specify function pointer to source function.
293 //========================================================================
294 template<class ELEMENT,class EIGEN_SOLVER>
296  const unsigned& n_element)
297 {
298  //Create the eigen solver
299  this->eigen_solver_pt() = new EIGEN_SOLVER;
300 
301  //Set domain length
302  double L=1.0;
303 
304  // Build mesh and store pointer in Problem
305  Problem::mesh_pt() = new OneDMesh<ELEMENT>(n_element,L);
306 
307  // Set the boundary conditions for this problem: By default, all nodal
308  // values are free -- we only need to pin the ones that have
309  // Dirichlet conditions.
310 
311  // Pin the single nodal value at the single node on mesh
312  // boundary 0 (= the left domain boundary at x=0)
313  mesh_pt()->boundary_node_pt(0,0)->pin(0);
314 
315  // Pin the single nodal value at the single node on mesh
316  // boundary 1 (= the right domain boundary at x=1)
317  mesh_pt()->boundary_node_pt(1,0)->pin(0);
318 
319  // Setup equation numbering scheme
320  assign_eqn_numbers();
321 
322 } // end of constructor
323 
324 
325 
326 //===start_of_doc=========================================================
327 /// Doc the solution in tecplot format. Label files with label.
328 //========================================================================
329 template<class ELEMENT,class EIGEN_SOLVER>
331 {
332 
333  ofstream some_file;
334  char filename[100];
335 
336  // Number of plot points
337  unsigned npts;
338  npts=5;
339 
340  // Output solution with specified number of plot points per element
341  sprintf(filename,"soln%i.dat",label);
342  some_file.open(filename);
343  mesh_pt()->output(some_file,npts);
344  some_file.close();
345 
346 } // end of doc
347 
348 //=======================start_of_solve==============================
349 /// Solve the eigenproblem
350 //===================================================================
351 template<class ELEMENT,class EIGEN_SOLVER>
353 solve(const unsigned& label)
354 {
355  //Set external storage for the eigenvalues
356  Vector<complex<double> > eigenvalue;
357  //Set external storage for the eigenvectors
358  Vector<DoubleVector> eigenvector_real;
359  Vector<DoubleVector> eigenvector_imag;
360  //Desired number eigenvalues
361  unsigned n_eval=4;
362 
363  //Solve the eigenproblem
364  this->solve_eigenproblem(n_eval,eigenvalue,eigenvector_real,eigenvector_imag);
365 
366  //We now need to sort the output based on the size of the real part
367  //of the eigenvalues.
368  //This is because the solver does not necessarily sort the eigenvalues
369  Vector<complex<double> > sorted_eigenvalue = eigenvalue;
370  sort(sorted_eigenvalue.begin(),sorted_eigenvalue.end(),
372 
373  //Read out the second smallest eigenvalue
374  complex<double> temp_evalue = sorted_eigenvalue[1];
375  unsigned second_smallest_index=0;
376  //Loop over the unsorted eigenvalues and find the entry that corresponds
377  //to our second smallest eigenvalue.
378  for(unsigned i=0;i<eigenvalue.size();i++)
379  {
380  //Note that equality tests for doubles are bad, but it was just
381  //sorted data, so should be fine
382  if(eigenvalue[i] == temp_evalue) {second_smallest_index=i; break;}
383  }
384 
385  //Normalise the eigenvector
386  {
387  //Get the dimension of the eigenvector
388  unsigned dim = eigenvector_real[second_smallest_index].nrow();
389  double length=0.0;
390  //Loop over all the entries
391  for(unsigned i=0;i<dim;i++)
392  {
393  //Add the contribution to the length
394  length += std::pow(eigenvector_real[second_smallest_index][i],2.0);
395  }
396  //Now take the magnitude
397  length = sqrt(length);
398  //Fix the sign
399  if(eigenvector_real[second_smallest_index][0] < 0) {length *= -1.0;}
400  //Finally normalise
401  for(unsigned i=0;i<dim;i++)
402  {
403  eigenvector_real[second_smallest_index][i] /= length;
404  }
405  }
406 
407  //Now assign the second eigenvector to the dofs of the problem
408  this->assign_eigenvector_to_dofs(eigenvector_real[second_smallest_index]);
409  //Output solution for this case (label output files with "1")
410  this->doc_solution(label);
411 
412  char filename[100];
413  sprintf(filename,"eigenvalues%i.dat",label);
414 
415  //Open an output file for the sorted eigenvalues
416  ofstream evalues(filename);
417  for(unsigned i=0;i<n_eval;i++)
418  {
419  //Print to screen
420  cout << sorted_eigenvalue[i].real() << " "
421  << sorted_eigenvalue[i].imag() << std::endl;
422  //Send to file
423  evalues << sorted_eigenvalue[i].real() << " "
424  << sorted_eigenvalue[i].imag() << std::endl;
425  }
426 
427  evalues.close();
428 } //end_of_solve
429 
430 
431 /// /////////////////////////////////////////////////////////////////////
432 /// /////////////////////////////////////////////////////////////////////
433 /// /////////////////////////////////////////////////////////////////////
434 
435 
436 //======start_of_main==================================================
437 /// Driver for 1D Poisson problem
438 //=====================================================================
439 int main(int argc, char **argv)
440 {
441 //Want to test Trilinos if we have it, so we must initialise MPI
442 //if we have compiled with it
443 #ifdef OOMPH_HAS_MPI
444  MPI_Helpers::init(argc,argv);
445 #endif
446 
447  // Set up the problem:
448  unsigned n_element=100; //Number of elements
449 
450  clock_t t_start1 = clock();
451  //Solve with LAPACK_QZ
452  {
454  problem(n_element);
455 
456  problem.solve(1);
457  }
458  clock_t t_end1 = clock();
459 
460 #ifdef OOMPH_HAS_TRILINOS
461  clock_t t_start2 = clock();
462 //Solve with Anasazi
463  {
464  // hierher Andrew: This doesn't seem to be included in the self tests
465  HarmonicProblem<QHarmonicElement<3>,ANASAZI> problem(n_element);
466  problem.solve(2);
467  }
468  clock_t t_end2 = clock();
469 #endif
470 
471  std::cout << "LAPACK TIME: " << (double)(t_end1 - t_start1)/CLOCKS_PER_SEC
472  << std::endl;
473 
474 #ifdef OOMPH_HAS_TRILINOS
475  std::cout << "ANASAZI TIME: " << (double)(t_end2 - t_start2)/CLOCKS_PER_SEC
476  << std::endl;
477 #else
478  //Need to skip test
479 #endif
480 
481 #ifdef OOMPH_HAS_MPI
482  MPI_Helpers::finalize();
483 #endif
484 
485 } // end of main
486 
487 
488 
489 
490 
491 
492 
493 
494 
495 
Function-type-object to perform comparison of complex data types Needed to sort the complex eigenvalu...
bool operator()(const complex< T > &x, const complex< T > &y) const
Comparison. Are the values identical or not?
Definition: harmonic.cc:51
1D Harmonic problem in unit interval.
Definition: harmonic.cc:269
HarmonicProblem(const unsigned &n_element)
Constructor: Pass number of elements and pointer to source function.
Definition: harmonic.cc:295
~HarmonicProblem()
Destructor (empty)
Definition: harmonic.cc:276
void doc_solution(const unsigned &label)
Doc the solution, pass the number of the case considered, so that output files can be distinguished.
Definition: harmonic.cc:330
void solve(const unsigned &label)
Solve the problem.
Definition: harmonic.cc:353
QHarmonicElement<NNODE_1D> elements are 1D Elements with NNODE_1D nodal points that are used to solve...
Definition: harmonic.cc:224
QHarmonicElement()
Constructor: Call constructors for QElement and Poisson equations.
Definition: harmonic.cc:230
double dshape_eulerian(const Vector< double > &s, Shape &psi, DShape &dpsidx) const
Shape, test functions & derivs. w.r.t. to global coords. Return Jacobian.
Definition: harmonic.cc:248
void output(ostream &outfile, const unsigned &Nplot)
Output function overloaded from HarmonicEquations.
Definition: harmonic.cc:241
double dshape_eulerian_at_knot(const unsigned &ipt, Shape &psi, DShape &dpsidx) const
Shape, test functions & derivs. w.r.t. to global coords. at integration point ipt....
Definition: harmonic.cc:256
unsigned required_nvalue(const unsigned &n) const
Required # of ‘values’ (pinned or dofs) at node n.
Definition: harmonic.cc:234
void output(ostream &outfile)
Output function overloaded from HarmonicEquations.
Definition: harmonic.cc:237
int main(int argc, char **argv)
///////////////////////////////////////////////////////////////////// ///////////////////////////////...
Definition: harmonic.cc:439