35   template<
unsigned DIM>
 
   48   template<
unsigned DIM>
 
   53     if ((strain.ncol() != DIM) || (strain.nrow() != DIM))
 
   55       std::ostringstream error_message;
 
   56       error_message << 
"Strain matrix is " << strain.ncol() << 
" x " 
   57                     << strain.nrow() << 
", but dimension of the equations is " 
   60         error_message.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
 
   64     unsigned n_position_type = this->nnodal_position_type();
 
   66     if (n_position_type != 1)
 
   68       std::ostringstream error_message;
 
   69       error_message << 
"PMLTimeHarmonicLinearElasticity is not yet " 
   70                     << 
"implemented for more than one position type" 
   73         error_message.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
 
   79     unsigned n_node = nnode();
 
   82     std::complex<unsigned> u_nodal_index[DIM];
 
   83     for (
unsigned i = 0; 
i < DIM; 
i++)
 
   85       u_nodal_index[
i] = u_index_time_harmonic_linear_elasticity(
i);
 
   90     DShape dpsidx(n_node, DIM);
 
   93     (void)dshape_eulerian(
s, psi, dpsidx);
 
   97       DIM, DIM, std::complex<double>(0.0, 0.0));
 
  100     for (
unsigned l = 0; l < n_node; l++)
 
  103       for (
unsigned i = 0; 
i < DIM; 
i++)
 
  106         const std::complex<double> u_value =
 
  107           std::complex<double>(this->nodal_value(l, u_nodal_index[
i].real()),
 
  108                                this->nodal_value(l, u_nodal_index[
i].imag()));
 
  111         for (
unsigned j = 0; j < DIM; j++)
 
  113           interpolated_dudx(
i, j) += u_value * dpsidx(l, j);
 
  119     for (
unsigned i = 0; 
i < DIM; 
i++)
 
  124       for (
int j = (DIM - 1); j >= 
static_cast<int>(
i); j--)
 
  127         if (
static_cast<int>(
i) != j)
 
  130             0.5 * (interpolated_dudx(
i, j) + interpolated_dudx(j, 
i));
 
  135           strain(
i, 
i) = interpolated_dudx(
i, 
i);
 
  139       for (
int j = (
i - 1); j >= 0; j--)
 
  141         strain(
i, j) = strain(j, 
i);
 
  151   template<
unsigned DIM>
 
  156     if ((stress.ncol() != DIM) || (stress.nrow() != DIM))
 
  158       std::ostringstream error_message;
 
  159       error_message << 
"Stress matrix is " << stress.ncol() << 
" x " 
  160                     << stress.nrow() << 
", but dimension of the equations is " 
  163         error_message.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
 
  169     this->get_strain(
s, strain);
 
  174     for (
unsigned i = 0; 
i < DIM; 
i++)
 
  176       for (
unsigned j = 0; j < DIM; j++)
 
  179         for (
unsigned k = 0; k < DIM; k++)
 
  181           for (
unsigned l = 0; l < DIM; l++)
 
  183             stress(
i, j) += this->
E(
i, j, k, l) * strain(k, l);
 
  195   template<
unsigned DIM>
 
  201     unsigned n_node = this->nnode();
 
  205     unsigned n_position_type = this->nnodal_position_type();
 
  207     if (n_position_type != 1)
 
  209       std::ostringstream error_message;
 
  210       error_message << 
"PMLTimeHarmonicLinearElasticity is not yet " 
  211                     << 
"implemented for more than one position type" 
  214         error_message.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
 
  218     if (this->Elasticity_tensor_pt == 0)
 
  221                           OOMPH_CURRENT_FUNCTION,
 
  222                           OOMPH_EXCEPTION_LOCATION);
 
  227     std::complex<unsigned> u_nodal_index[DIM];
 
  228     for (
unsigned i = 0; 
i < DIM; 
i++)
 
  230       u_nodal_index[
i] = this->u_index_time_harmonic_linear_elasticity(
i);
 
  235     const double omega_sq_local = this->omega_sq();
 
  239     DShape dpsidx(n_node, DIM);
 
  242     unsigned n_intpt = this->integral_pt()->nweight();
 
  248     int local_eqn_real = 0, local_unknown_real = 0;
 
  249     int local_eqn_imag = 0, local_unknown_imag = 0;
 
  252     for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
 
  255       for (
unsigned i = 0; 
i < DIM; ++
i)
 
  257         s[
i] = this->integral_pt()->knot(ipt, 
i);
 
  261       double w = this->integral_pt()->weight(ipt);
 
  264       double J = this->dshape_eulerian_at_knot(ipt, psi, dpsidx);
 
  271         DIM, std::complex<double>(0.0, 0.0));
 
  276         DIM, DIM, std::complex<double>(0.0, 0.0));
 
  279       for (
unsigned l = 0; l < n_node; l++)
 
  282         for (
unsigned i = 0; 
i < DIM; 
i++)
 
  285           interpolated_x[
i] += this->raw_nodal_position(l, 
i) * psi(l);
 
  289           const std::complex<double> u_value = std::complex<double>(
 
  290             this->raw_nodal_value(l, u_nodal_index[
i].real()),
 
  291             this->raw_nodal_value(l, u_nodal_index[
i].imag()));
 
  293           interpolated_u[
i] += u_value * psi(l);
 
  296           for (
unsigned j = 0; j < DIM; j++)
 
  298             interpolated_dudx(
i, j) += u_value * dpsidx(l, j);
 
  305       this->body_force(interpolated_x, body_force_vec);
 
  319       std::complex<double> pml_jacobian_det;
 
  320       this->compute_pml_coefficients(
 
  321         ipt, interpolated_x, pml_inverse_jacobian_diagonal, pml_jacobian_det);
 
  324       for (
unsigned l = 0; l < n_node; l++)
 
  327         for (
unsigned a = 0; a < DIM; a++)
 
  330           local_eqn_real = this->nodal_local_eqn(l, u_nodal_index[a].real());
 
  331           local_eqn_imag = this->nodal_local_eqn(l, u_nodal_index[a].imag());
 
  338           std::complex<double> mass_jacobian_contribution =
 
  339             -omega_sq_local * pml_jacobian_det;
 
  340           std::complex<double> mass_residual_contribution =
 
  341             (-omega_sq_local * interpolated_u[a] - body_force_vec[a]) *
 
  345           std::complex<double> stress_jacobian_contributions[DIM][DIM][DIM];
 
  346           std::complex<double> stress_residual_contributions[DIM];
 
  347           for (
unsigned b = 0; b < DIM; b++)
 
  349             stress_residual_contributions[b] = 0.0;
 
  350             for (
unsigned c = 0; c < DIM; c++)
 
  352               for (
unsigned d = 0; d < DIM; d++)
 
  354                 stress_jacobian_contributions[b][c][d] =
 
  355                   this->
E(a, b, c, d) * pml_jacobian_det *
 
  356                   pml_inverse_jacobian_diagonal[b] *
 
  357                   pml_inverse_jacobian_diagonal[d];
 
  359                 stress_residual_contributions[b] +=
 
  360                   stress_jacobian_contributions[b][c][d] *
 
  361                   interpolated_dudx(c, d);
 
  369           if (local_eqn_real >= 0)
 
  372             residuals[local_eqn_real] +=
 
  373               mass_residual_contribution.real() * psi(l) * 
W;
 
  376             for (
unsigned b = 0; b < DIM; b++)
 
  379               residuals[local_eqn_real] +=
 
  380                 stress_residual_contributions[b].real() * dpsidx(l, b) * 
W;
 
  387               for (
unsigned l2 = 0; l2 < n_node; l2++)
 
  390                 for (
unsigned c = 0; c < DIM; c++)
 
  393                     this->nodal_local_eqn(l2, u_nodal_index[c].real());
 
  395                     this->nodal_local_eqn(l2, u_nodal_index[c].imag());
 
  397                   if (local_unknown_real >= 0)
 
  402                       jacobian(local_eqn_real, local_unknown_real) +=
 
  403                         mass_jacobian_contribution.real() * psi(l) * psi(l2) *
 
  408                     for (
unsigned b = 0; b < DIM; b++)
 
  410                       for (
unsigned d = 0; d < DIM; d++)
 
  413                         jacobian(local_eqn_real, local_unknown_real) +=
 
  414                           stress_jacobian_contributions[b][c][d].real() *
 
  415                           dpsidx(l2, d) * dpsidx(l, b) * 
W;
 
  420                   if (local_unknown_imag >= 0)
 
  425                       jacobian(local_eqn_real, local_unknown_imag) -=
 
  426                         mass_jacobian_contribution.imag() * psi(l) * psi(l2) *
 
  431                     for (
unsigned b = 0; b < DIM; b++)
 
  433                       for (
unsigned d = 0; d < DIM; d++)
 
  436                         jacobian(local_eqn_real, local_unknown_imag) -=
 
  437                           stress_jacobian_contributions[b][c][d].imag() *
 
  438                           dpsidx(l2, d) * dpsidx(l, b) * 
W;
 
  450           if (local_eqn_imag >= 0)
 
  453             residuals[local_eqn_imag] +=
 
  454               mass_residual_contribution.imag() * psi(l) * 
W;
 
  457             for (
unsigned b = 0; b < DIM; b++)
 
  460               residuals[local_eqn_imag] +=
 
  461                 stress_residual_contributions[b].imag() * dpsidx(l, b) * 
W;
 
  468               for (
unsigned l2 = 0; l2 < n_node; l2++)
 
  471                 for (
unsigned c = 0; c < DIM; c++)
 
  474                     this->nodal_local_eqn(l2, u_nodal_index[c].imag());
 
  476                     this->nodal_local_eqn(l2, u_nodal_index[c].real());
 
  478                   if (local_unknown_imag >= 0)
 
  483                       jacobian(local_eqn_imag, local_unknown_imag) +=
 
  484                         mass_jacobian_contribution.real() * psi(l) * psi(l2) *
 
  489                     for (
unsigned b = 0; b < DIM; b++)
 
  491                       for (
unsigned d = 0; d < DIM; d++)
 
  494                         jacobian(local_eqn_imag, local_unknown_imag) +=
 
  495                           stress_jacobian_contributions[b][c][d].real() *
 
  496                           dpsidx(l2, d) * dpsidx(l, b) * 
W;
 
  501                   if (local_unknown_real >= 0)
 
  506                       jacobian(local_eqn_imag, local_unknown_real) +=
 
  507                         mass_jacobian_contribution.imag() * psi(l) * psi(l2) *
 
  512                     for (
unsigned b = 0; b < DIM; b++)
 
  514                       for (
unsigned d = 0; d < DIM; d++)
 
  517                         jacobian(local_eqn_imag, local_unknown_real) +=
 
  518                           stress_jacobian_contributions[b][c][d].imag() *
 
  519                           dpsidx(l2, d) * dpsidx(l, b) * 
W;
 
  537   template<
unsigned DIM>
 
  539     std::ostream& outfile,
 
  540     const unsigned& nplot,
 
  550     outfile << this->tecplot_zone_string(nplot);
 
  556     unsigned num_plot_points = this->nplot_points(nplot);
 
  557     for (
unsigned iplot = 0; iplot < num_plot_points; iplot++)
 
  560       this->get_s_plot(iplot, nplot, 
s);
 
  563       this->interpolated_x(
s, x);
 
  566       (*exact_soln_pt)(x, exact_soln);
 
  569       for (
unsigned i = 0; 
i < DIM; 
i++)
 
  571         outfile << x[
i] << 
" ";
 
  573       for (
unsigned i = 0; 
i < 2 * DIM; 
i++)
 
  575         outfile << exact_soln[
i] << 
" ";
 
  577       outfile << std::endl;
 
  581     this->write_tecplot_zone_footer(outfile, nplot);
 
  588   template<
unsigned DIM>
 
  590     std::ostream& outfile, 
const unsigned& nplot)
 
  598     outfile << this->tecplot_zone_string(nplot);
 
  601     unsigned num_plot_points = this->nplot_points(nplot);
 
  602     for (
unsigned iplot = 0; iplot < num_plot_points; iplot++)
 
  605       this->get_s_plot(iplot, nplot, 
s);
 
  608       this->interpolated_x(
s, x);
 
  609       this->interpolated_u_time_harmonic_linear_elasticity(
s, u);
 
  612       for (
unsigned i = 0; 
i < DIM; 
i++)
 
  614         outfile << x[
i] << 
" ";
 
  618       for (
unsigned i = 0; 
i < DIM; 
i++)
 
  620         outfile << u[
i].real() << 
" ";
 
  624       for (
unsigned i = 0; 
i < DIM; 
i++)
 
  626         outfile << u[
i].imag() << 
" ";
 
  629       outfile << std::endl;
 
  633     this->write_tecplot_zone_footer(outfile, nplot);
 
  650   template<
unsigned DIM>
 
  652     std::ostream& outfile,
 
  655     const unsigned& nplot)
 
  666     outfile << this->tecplot_zone_string(nplot);
 
  669     unsigned num_plot_points = this->nplot_points(nplot);
 
  670     for (
unsigned iplot = 0; iplot < num_plot_points; iplot++)
 
  673       this->get_s_plot(iplot, nplot, 
s);
 
  676       this->interpolated_x(
s, x);
 
  677       this->interpolated_u_time_harmonic_linear_elasticity(
s, u);
 
  680       (*incoming_wave_fct_pt)(x, incoming_soln);
 
  683       for (
unsigned i = 0; 
i < DIM; 
i++)
 
  685         outfile << x[
i] << 
" ";
 
  689       for (
unsigned i = 0; 
i < DIM; 
i++)
 
  691         outfile << (u[
i].real() + incoming_soln[2 * 
i]) * cos(phi) +
 
  692                      (u[
i].imag() + incoming_soln[2 * 
i + 1]) * sin(phi)
 
  696       outfile << std::endl;
 
  700     this->write_tecplot_zone_footer(outfile, nplot);
 
  714   template<
unsigned DIM>
 
  716     std::ostream& outfile, 
const double& phi, 
const unsigned& nplot)
 
  724     outfile << this->tecplot_zone_string(nplot);
 
  727     unsigned num_plot_points = this->nplot_points(nplot);
 
  728     for (
unsigned iplot = 0; iplot < num_plot_points; iplot++)
 
  731       this->get_s_plot(iplot, nplot, 
s);
 
  734       this->interpolated_x(
s, x);
 
  735       this->interpolated_u_time_harmonic_linear_elasticity(
s, u);
 
  738       for (
unsigned i = 0; 
i < DIM; 
i++)
 
  740         outfile << x[
i] << 
" ";
 
  744       for (
unsigned i = 0; 
i < DIM; 
i++)
 
  746         outfile << u[
i].real() * cos(phi) + u[
i].imag() * sin(phi) << 
" ";
 
  749       outfile << std::endl;
 
  753     this->write_tecplot_zone_footer(outfile, nplot);
 
  767   template<
unsigned DIM>
 
  769     std::ostream& outfile, 
const double& phi, 
const unsigned& nplot)
 
  777     outfile << this->tecplot_zone_string(nplot);
 
  780     unsigned num_plot_points = this->nplot_points(nplot);
 
  781     for (
unsigned iplot = 0; iplot < num_plot_points; iplot++)
 
  784       this->get_s_plot(iplot, nplot, 
s);
 
  787       this->interpolated_x(
s, x);
 
  788       this->interpolated_u_time_harmonic_linear_elasticity(
s, u);
 
  791       for (
unsigned i = 0; 
i < DIM; 
i++)
 
  793         outfile << x[
i] << 
" ";
 
  797       for (
unsigned i = 0; 
i < DIM; 
i++)
 
  799         outfile << u[
i].imag() * cos(phi) - u[
i].real() * sin(phi) << 
" ";
 
  802       outfile << std::endl;
 
  806     this->write_tecplot_zone_footer(outfile, nplot);
 
  813   template<
unsigned DIM>
 
  815     FILE* file_pt, 
const unsigned& nplot)
 
  821     fprintf(file_pt, 
"%s", this->tecplot_zone_string(nplot).c_str());
 
  824     unsigned num_plot_points = this->nplot_points(nplot);
 
  825     for (
unsigned iplot = 0; iplot < num_plot_points; iplot++)
 
  828       this->get_s_plot(iplot, nplot, 
s);
 
  831       for (
unsigned i = 0; 
i < DIM; 
i++)
 
  833         fprintf(file_pt, 
"%g ", this->interpolated_x(
s, 
i));
 
  837       for (
unsigned i = 0; 
i < DIM; 
i++)
 
  842           this->interpolated_u_time_harmonic_linear_elasticity(
s, 
i).real());
 
  844       for (
unsigned i = 0; 
i < DIM; 
i++)
 
  849           this->interpolated_u_time_harmonic_linear_elasticity(
s, 
i).imag());
 
  852     fprintf(file_pt, 
"\n");
 
  855     this->write_tecplot_zone_footer(file_pt, nplot);
 
  862   template<
unsigned DIM>
 
  878     unsigned n_node = this->nnode();
 
  883     unsigned n_intpt = this->integral_pt()->nweight();
 
  886     for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
 
  889       for (
unsigned i = 0; 
i < DIM; 
i++)
 
  891         s[
i] = this->integral_pt()->knot(ipt, 
i);
 
  895       double w = this->integral_pt()->weight(ipt);
 
  898       double J = this->J_eulerian(
s);
 
  904       this->interpolated_u_time_harmonic_linear_elasticity(
s, disp);
 
  907       for (
unsigned ii = 0; ii < DIM; ii++)
 
  909         norm += (disp[ii].real() * disp[ii].real() +
 
  910                  disp[ii].imag() * disp[ii].imag()) *
 
  922   template<
unsigned DIM>
 
  924     std::ostream& outfile,
 
  943     unsigned n_node = this->nnode();
 
  948     unsigned n_intpt = this->integral_pt()->nweight();
 
  954     for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
 
  957       for (
unsigned i = 0; 
i < DIM; 
i++)
 
  959         s[
i] = this->integral_pt()->knot(ipt, 
i);
 
  963       double w = this->integral_pt()->weight(ipt);
 
  966       double J = this->J_eulerian(
s);
 
  972       this->interpolated_x(
s, x);
 
  975       (*exact_soln_pt)(x, exact_soln);
 
  978       this->interpolated_u_time_harmonic_linear_elasticity(
s, disp);
 
  981       for (
unsigned ii = 0; ii < DIM; ii++)
 
  984         error += ((exact_soln[ii] - disp[ii].real()) *
 
  985                     (exact_soln[ii] - disp[ii].real()) +
 
  986                   (exact_soln[ii + DIM] - disp[ii].imag()) *
 
  987                     (exact_soln[ii + DIM] - disp[ii].imag())) *
 
  989         norm += (disp[ii].real() * disp[ii].real() +
 
  990                  disp[ii].imag() * disp[ii].imag()) *
 
 1005   template<
unsigned DIM>
 
A mapping function proposed by Bermudez et al, similar to the one above but is continuous across the ...
A Class for the derivatives of shape functions The class design is essentially the same as Shape,...
//////////////////////////////////////////////////////////////////////////// ////////////////////////...
void(* SteadyExactSolutionFctPt)(const Vector< double > &, Vector< double > &)
Function pointer for function that computes vector-valued steady "exact solution"  as .
An OomphLibError object which should be thrown when an run-time error is encountered....
A base class for elements that solve the equations of time-harmonic linear elasticity in Cartesian co...
static double Default_omega_sq_value
Static default value for square of frequency.
void get_strain(const Vector< double > &s, DenseMatrix< std::complex< double >> &strain) const
Return the strain tensor.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
void output_imag(std::ostream &outfile, const double &phi, const unsigned &n_plot)
Output function for imaginary part of full time-dependent solution u = Im( (u_r +i u_i) exp(-i omega ...
void get_stress(const Vector< double > &s, DenseMatrix< std::complex< double >> &sigma) const
Return the Cauchy stress tensor, as calculated from the elasticity tensor at specified local coordina...
void compute_norm(double &norm)
Compute norm of solution: square of the L2 norm.
void output(std::ostream &outfile)
Output: x,y,[z],u_r,v_r,[w_r],u_i,v_i,[w_i].
virtual void fill_in_generic_contribution_to_residuals_time_harmonic_linear_elasticity(Vector< double > &residuals, DenseMatrix< double > &jacobian, unsigned flag)
Private helper function to compute residuals and (if requested via flag) also the Jacobian matrix.
void output_total_real(std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt incoming_wave_fct_pt, const double &phi, const unsigned &nplot)
Output function for real part of full time-dependent solution constructed by adding the scattered fie...
void output_fct(std::ostream &outfile, const unsigned &nplot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Output exact solution x,y,[z],u_r,v_r,[w_r],u_i,v_i,[w_i].
void output_real(std::ostream &outfile, const double &phi, const unsigned &n_plot)
Output function for real part of full time-dependent solution u = Re( (u_r +i u_i) exp(-i omega t)) a...
void compute_error(std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
Get error against and norm of exact solution.
///////////////////////////////////////////////////////////////////// ///////////////////////////////...
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
//////////////////////////////////////////////////////////////////// ////////////////////////////////...