33 #include "meshes/one_d_mesh.h"
37 using namespace oomph;
51 bool operator()(
const complex<T> &x,
const complex<T> &y)
const
53 return x.real() < y.real();
70 class HarmonicEquations :
public virtual FiniteElement
75 HarmonicEquations() {}
80 virtual inline double u(
const unsigned& n)
const
81 {
return nodal_value(n,0);}
84 void output(ostream &outfile)
87 output(outfile,nplot);
92 void output(ostream &outfile,
const unsigned &nplot)
98 outfile << tecplot_zone_string(nplot);
101 unsigned num_plot_points=nplot_points(nplot);
102 for (
unsigned iplot=0;iplot<num_plot_points;iplot++)
105 get_s_plot(iplot,nplot,s);
107 outfile << interpolated_x(s,0) <<
" " << interpolated_u(s) << std::endl;
110 write_tecplot_zone_footer(outfile,nplot);
115 void fill_in_contribution_to_jacobian_and_mass_matrix(
116 Vector<double> &residuals,
117 DenseMatrix<double> &jacobian, DenseMatrix<double> &mass_matrix)
120 unsigned n_node = nnode();
124 DShape dpsidx(n_node,1);
127 unsigned n_intpt = integral_pt()->nweight();
130 int local_eqn=0, local_unknown=0;
133 for(
unsigned ipt=0;ipt<n_intpt;ipt++)
136 double w = integral_pt()->weight(ipt);
139 double J = dshape_eulerian_at_knot(ipt,psi,dpsidx);
146 for(
unsigned l=0;l<n_node;l++)
149 local_eqn = u_local_eqn(l);
154 for(
unsigned l2=0;l2<n_node;l2++)
156 local_unknown = u_local_eqn(l2);
158 if(local_unknown >= 0)
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;
170 inline double interpolated_u(
const Vector<double> &s)
const
172 unsigned n_node = nnode();
181 double interpolated_u = 0.0;
184 for(
unsigned l=0;l<n_node;l++) {interpolated_u+=u(l)*psi[l];}
187 return(interpolated_u);
194 virtual double dshape_eulerian(
const Vector<double> &s,
196 DShape &dpsidx)
const=0;
200 virtual double dshape_eulerian_at_knot(
const unsigned &ipt,
202 DShape &dpsidx)
const=0;
207 virtual inline int u_local_eqn(
const unsigned &n)
208 {
return nodal_local_eqn(n,0);}
221 template <
unsigned NNODE_1D>
223 public HarmonicEquations
238 {HarmonicEquations::output(outfile);}
241 void output(ostream &outfile,
const unsigned &Nplot)
242 {HarmonicEquations::output(outfile,Nplot);}
250 DShape &dpsidx)
const
251 {
return QElement<1,NNODE_1D>::dshape_eulerian(s,psi,dpsidx);}
258 DShape &dpsidx)
const
259 {
return QElement<1,NNODE_1D>::dshape_eulerian_at_knot(ipt,psi,dpsidx);}
267 template<
class ELEMENT,
class EIGEN_SOLVER>
279 void solve(
const unsigned &label);
283 void doc_solution(
const unsigned& label);
294 template<
class ELEMENT,
class EIGEN_SOLVER>
296 const unsigned& n_element)
299 this->eigen_solver_pt() =
new EIGEN_SOLVER;
305 Problem::mesh_pt() =
new OneDMesh<ELEMENT>(n_element,L);
313 mesh_pt()->boundary_node_pt(0,0)->pin(0);
317 mesh_pt()->boundary_node_pt(1,0)->pin(0);
320 assign_eqn_numbers();
329 template<
class ELEMENT,
class EIGEN_SOLVER>
341 sprintf(filename,
"soln%i.dat",label);
342 some_file.open(filename);
343 mesh_pt()->output(some_file,npts);
351 template<
class ELEMENT,
class EIGEN_SOLVER>
353 solve(
const unsigned& label)
356 Vector<complex<double> > eigenvalue;
358 Vector<DoubleVector> eigenvector_real;
359 Vector<DoubleVector> eigenvector_imag;
364 this->solve_eigenproblem(n_eval,eigenvalue,eigenvector_real,eigenvector_imag);
369 Vector<complex<double> > sorted_eigenvalue = eigenvalue;
370 sort(sorted_eigenvalue.begin(),sorted_eigenvalue.end(),
374 complex<double> temp_evalue = sorted_eigenvalue[1];
375 unsigned second_smallest_index=0;
378 for(
unsigned i=0;i<eigenvalue.size();i++)
382 if(eigenvalue[i] == temp_evalue) {second_smallest_index=i;
break;}
388 unsigned dim = eigenvector_real[second_smallest_index].nrow();
391 for(
unsigned i=0;i<dim;i++)
394 length += std::pow(eigenvector_real[second_smallest_index][i],2.0);
397 length = sqrt(length);
399 if(eigenvector_real[second_smallest_index][0] < 0) {length *= -1.0;}
401 for(
unsigned i=0;i<dim;i++)
403 eigenvector_real[second_smallest_index][i] /= length;
408 this->assign_eigenvector_to_dofs(eigenvector_real[second_smallest_index]);
410 this->doc_solution(label);
413 sprintf(filename,
"eigenvalues%i.dat",label);
416 ofstream evalues(filename);
417 for(
unsigned i=0;i<n_eval;i++)
420 cout << sorted_eigenvalue[i].real() <<
" "
421 << sorted_eigenvalue[i].imag() << std::endl;
423 evalues << sorted_eigenvalue[i].real() <<
" "
424 << sorted_eigenvalue[i].imag() << std::endl;
439 int main(
int argc,
char **argv)
444 MPI_Helpers::init(argc,argv);
448 unsigned n_element=100;
450 clock_t t_start1 = clock();
458 clock_t t_end1 = clock();
460 #ifdef OOMPH_HAS_TRILINOS
461 clock_t t_start2 = clock();
468 clock_t t_end2 = clock();
471 std::cout <<
"LAPACK TIME: " << (double)(t_end1 - t_start1)/CLOCKS_PER_SEC
474 #ifdef OOMPH_HAS_TRILINOS
475 std::cout <<
"ANASAZI TIME: " << (double)(t_end2 - t_start2)/CLOCKS_PER_SEC
482 MPI_Helpers::finalize();
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?
1D Harmonic problem in unit interval.
HarmonicProblem(const unsigned &n_element)
Constructor: Pass number of elements and pointer to source function.
~HarmonicProblem()
Destructor (empty)
void doc_solution(const unsigned &label)
Doc the solution, pass the number of the case considered, so that output files can be distinguished.
void solve(const unsigned &label)
Solve the problem.
QHarmonicElement<NNODE_1D> elements are 1D Elements with NNODE_1D nodal points that are used to solve...
QHarmonicElement()
Constructor: Call constructors for QElement and Poisson equations.
double dshape_eulerian(const Vector< double > &s, Shape &psi, DShape &dpsidx) const
Shape, test functions & derivs. w.r.t. to global coords. Return Jacobian.
void output(ostream &outfile, const unsigned &Nplot)
Output function overloaded from HarmonicEquations.
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....
unsigned required_nvalue(const unsigned &n) const
Required # of ‘values’ (pinned or dofs) at node n.
void output(ostream &outfile)
Output function overloaded from HarmonicEquations.
int main(int argc, char **argv)
///////////////////////////////////////////////////////////////////// ///////////////////////////////...