NiHu  2.0
guiggiani_1992.hpp
Go to the documentation of this file.
1 // This file is a part of NiHu, a C++ BEM template library.
2 //
3 // Copyright (C) 2012-2014 Peter Fiala <fiala@hit.bme.hu>
4 // Copyright (C) 2012-2014 Peter Rucz <rucz@hit.bme.hu>
5 //
6 // This program is free software: you can redistribute it and/or modify
7 // it under the terms of the GNU General Public License as published by
8 // the Free Software Foundation, either version 3 of the License, or
9 // (at your option) any later version.
10 //
11 // This program is distributed in the hope that it will be useful,
12 // but WITHOUT ANY WARRANTY; without even the implied warranty of
13 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 // GNU General Public License for more details.
15 //
16 // You should have received a copy of the GNU General Public License
17 // along with this program. If not, see <http://www.gnu.org/licenses/>.
18 
28 #ifndef GUIGGIANI_1992_HPP_INCLUDED
29 #define GUIGGIANI_1992_HPP_INCLUDED
30 
31 #include <boost/math/constants/constants.hpp>
32 
33 #include "line_quad_store.hpp"
34 #include "location_normal.hpp"
35 #include "weighted_input.hpp"
36 #include "../core/field.hpp"
37 #include "../core/kernel.hpp"
38 #include "../core/shapeset.hpp"
39 #include "../util/block_product.hpp"
40 #include "../util/store_pattern.hpp"
41 
42 #include <cmath>
43 #if NIHU_MEX_DEBUGGING
44 #include <mex.h>
45 #endif
46 
47 
48 namespace NiHu
49 {
50 
63 template <class singularity_type>
65 
67 typedef std::integral_constant<int, -2> _m2;
69 typedef std::integral_constant<int, -1> _m1;
71 typedef std::integral_constant<int, 0> _0;
73 typedef std::integral_constant<int, 1> _1;
75 typedef std::integral_constant<int, 2> _2;
76 
83 template <class TrialField, class Kernel, unsigned RadialOrder, unsigned TangentialOrder = RadialOrder, class Enable = void>
84 class guiggiani;
85 
86 
88 template <class TrialField, class Kernel, unsigned RadialOrder, unsigned TangentialOrder>
89 class guiggiani<TrialField, Kernel, RadialOrder, TangentialOrder, typename std::enable_if<
90  element_traits::is_surface_element<typename TrialField::elem_t>::value
91 >::type>
92 {
93 public:
95  enum {
96  radial_order = RadialOrder,
97  tangential_order = TangentialOrder
98  };
99 
101  typedef TrialField trial_field_t;
103  typedef typename trial_field_t::elem_t elem_t;
104 
106  typedef typename elem_t::domain_t domain_t;
108  typedef typename domain_t::xi_t xi_t;
110  typedef typename elem_t::scalar_t scalar_t;
112  typedef Eigen::Matrix<scalar_t, domain_t::dimension, domain_t::dimension> trans_t;
113 
115  typedef typename elem_t::x_t x_t;
116 
118  typedef typename trial_field_t::nset_t trial_nset_t;
120  typedef typename trial_nset_t::shape_t trial_n_shape_t;
121 
123  typedef Kernel kernel_t;
124 
129 
132 
137 
139  typedef typename semi_block_product_result_type<
140  typename singular_core_t::result_t, trial_n_shape_t
142 
144  typedef typename semi_block_product_result_type<
145  typename kernel_t::result_t, trial_n_shape_t
146  >::type total_result_t;
147 
149  enum {
152  };
153 
158  guiggiani(elem_t const &elem, kernel_base<kernel_t> const &kernel)
159  : m_elem(elem), m_kernel(kernel)
160  {
161  }
162 
163 private:
168  void compute_xi0(xi_t const &xi0, x_t const &normal)
169  {
170  // store the singular point and the unit normal vector
171  m_xi0 = xi0;
172  m_n0 = normal.normalized();
173 
174  // get the physical location and its gradient
175  m_x0 = m_elem.get_x(m_xi0);
176  typename elem_t::dx_return_type dx = m_elem.get_dx(m_xi0);
177 
178  scalar_t cosgamma = dx.col(shape_derivative_index::dXI).dot(dx.col(shape_derivative_index::dETA)) /
179  dx.col(shape_derivative_index::dXI).norm() / dx.col(shape_derivative_index::dETA).norm();
180  scalar_t u = dx.col(shape_derivative_index::dETA).norm() / dx.col(shape_derivative_index::dXI).norm();
181  m_T <<
182  1.0, cosgamma * u,
183  0.0, std::sqrt(1.0-cosgamma*cosgamma) * u;
184  // this scaling ensures that A is always equal to 1.0
185  m_T *= dx.col(shape_derivative_index::dXI).norm();
186  m_Tinv = m_T.inverse();
187 
188  m_eta0 = m_T * m_xi0;
189 
190  m_Jvec_series[0] = m_elem.get_normal(m_xi0) / m_T.determinant();
191  m_N_series[0] = trial_nset_t::template eval_shape<0>(xi0);
192 
193  // geometrical parameters (planar helpers)
194  unsigned const N = domain_t::num_corners;
195  for (unsigned n = 0; n < N; ++n)
196  {
197  xi_t c1 = m_T * domain_t::get_corner(n); // corner
198  xi_t c2 = m_T * domain_t::get_corner((n + 1) % N); // next corner
199  xi_t l = xi_t::Zero();
200  l(0) = 1.0;
201  if ((c2-c1).norm() > 1e-3)
202  l = (c2 - c1).normalized(); // side unit vector
203 
204  xi_t d1 = c1 - m_eta0; // vector to corners
205  xi_t d0 = d1 - l*d1.dot(l); // perpendicular to side
206 
207  m_theta_lim[n] = std::atan2(d1(1), d1(0)); // corner angle
208  m_theta0[n] = std::atan2(d0(1), d0(0)); // mid angle
209  m_ref_distance[n] = d0.norm(); // distance to side
210  }
211  }
212 
216  void compute_theta(scalar_t theta)
217  {
218  // contains cos theta sin theta in xi system
219  xi_t xi = m_Tinv.col(shape_derivative_index::dXI) * std::cos(theta)
220  + m_Tinv.col(shape_derivative_index::dETA) * std::sin(theta);
221 
222  typename elem_t::dx_return_type dx = m_elem.get_dx(m_xi0);
223 
224  m_rvec_series[0] = dx.col(shape_derivative_index::dXI) * xi(shape_derivative_index::dXI, 0)
226  if (laurent_order > 1) // compile_time IF
227  {
228  typename elem_t::ddx_return_type ddx = m_elem.get_ddx(m_xi0);
229 
230  m_rvec_series[1] = ddx.col(shape_derivative_index::dXIXI) * xi(shape_derivative_index::dXI, 0)*xi(shape_derivative_index::dXI, 0) / 2.0 +
233 
234  m_Jvec_series[1] = (
236  (
239  )
240  +
242  (
245  )
246  ) / m_T.determinant();
247 
248  auto dN = trial_nset_t::template eval_shape<1>(m_xi0);
249  m_N_series[1] = dN.col(shape_derivative_index::dXI)*xi(shape_derivative_index::dXI, 0)
251  }
252  }
253 
254 public:
261  template <class result_t>
262  void integrate(result_t &&I, xi_t const &xi0, x_t const &normal)
263  {
264  using namespace boost::math::double_constants;
265 
266 #if NIHU_MEX_DEBUGGING
267  static bool printed = false;
268  if (!printed)
269  {
270  mexPrintf("Computing integral with Guiggiani.\nRadial_order: %d\nTangential order: %d\n", RadialOrder, TangentialOrder);
271  printed = true;
272  }
273 #endif
274 
275  // compute the xi0-related quantities
276  compute_xi0(xi0, normal);
277 
278  // bind the kernel at the test input
279  test_input_t test_input(m_elem, m_xi0);
280  auto bound = m_kernel.bind(test_input);
281 
282  // iterate through triangles
283  unsigned const N = domain_t::num_corners;
284  for (unsigned n = 0; n < N; ++n)
285  {
286  // get angular integration limits
287  scalar_t t1 = m_theta_lim[n];
288  scalar_t t2 = m_theta_lim[(n + 1) % N];
289 
291  if (std::abs(t2-t1) < 1e-3)
292  continue;
293 
294  // we assume that the domain's corners are listed in positive order */
295  if (std::abs(t2 - t1) > pi)
296  t2 += two_pi;
297 
298  // theta integration
299  for (auto const &q_theta : line_quad_store<tangential_order>::quadrature)
300  {
301  // compute theta base point and weight
302  scalar_t xx = q_theta.get_xi()(0);
303  scalar_t theta = ((1.0 - xx) * t1 + (1.0 + xx) * t2) / 2.0;
304  scalar_t w_theta = q_theta.get_w() * (t2 - t1) / 2.0;
305 
306  // compute theta-related members
307  compute_theta(theta);
308  laurent_t::eval(*this);
309 
310  // reference domain's limit
311  scalar_t rho_lim = m_ref_distance[n] / std::cos(theta - m_theta0[n]);
312 
313  laurent_coeff_t toadd = m_Fcoeffs[0] * std::log(rho_lim);
314  if (laurent_order > 1) // compile_time IF
315  toadd -= m_Fcoeffs[1] / rho_lim;
316  toadd *= w_theta;
317 
318  I += toadd;
319 
320 /* this was for old Eigen versions (scalar casting)
321  for (int r = 0; r < I.rows(); ++r) // loop needed for scalar casting
322  for (int c = 0; c < I.cols(); ++c)
323  I(r,c) += toadd(r,c);
324 */
325 
326  // radial part of surface integration
327  for (auto const &q_rho : line_quad_store<radial_order>::quadrature)
328  {
329  // quadrature base point and weight
330  scalar_t rho = (1.0 + q_rho.get_xi()(0)) * rho_lim / 2.0;
331  scalar_t w_rho = q_rho.get_w() * rho_lim / 2.0;
332 
333  // compute location in original reference domain
334  xi_t eta(rho*std::cos(theta), rho*std::sin(theta));
335  xi_t xi = m_Tinv * eta + m_xi0;
336 
337  // evaluate G * N * J * rho
338  w_trial_input_t trial_input(m_elem, xi);
339  typename kernel_t::result_t GJrho =
340  bound(trial_input) * (trial_input.get_jacobian() / m_T.determinant() * rho);
341  typename trial_nset_t::shape_t N =
342  trial_nset_t::template eval_shape<0>(xi);
343  total_result_t F = semi_block_product(GJrho, N);
344 
345  // subtract the analytical singularity
346  laurent_coeff_t singular_part = m_Fcoeffs[0];
347  if (laurent_order > 1) // compile_time IF
348  singular_part += m_Fcoeffs[1] / rho;
349  singular_part /= rho;
350 
351  F -= singular_part;
352 /* this loop was needed for casting from Matrix<double> to Matrix<complex>
353  for (int r = 0; r < F.rows(); ++r)
354  for (int c = 0; c < F.cols(); ++c)
355  F(r,c) -= singular_part(r,c);
356 */
357 
358  // surface integral accumulation
359  I += w_theta * w_rho * F;
360  } // rho loop
361  } // theta loop
362  } // element sides
363  } // end of function integrate
364 
366  template <class T, T order>
367  x_t const &get_rvec_series(std::integral_constant<T, order>) const
368  {
369  static_assert((order-1)>=0, "Required distance Taylor coefficient too low");
370  static_assert((order-1) < laurent_order, "Required distance Taylor coefficient too high");
371  return m_rvec_series[order-1];
372  }
373 
375  template <class T, T order>
376  x_t const &get_Jvec_series(std::integral_constant<T, order>) const
377  {
378  static_assert(order < laurent_order, "Required Jacobian Taylor coefficient too high");
379  static_assert(order >= 0, "Required Jacobian Taylor coefficient too low");
380  return m_Jvec_series[order];
381  }
382 
384  template <class T, T order>
385  trial_n_shape_t const &get_shape_series(std::integral_constant<T, order>) const
386  {
387  static_assert(order < laurent_order, "Required shape set Taylor coefficient too high");
388  static_assert(order >= 0, "Required Jacobian Taylor coefficient too low");
389  return m_N_series[order];
390  }
391 
393  template <class T, T order>
394  laurent_coeff_t &get_laurent_coeff(std::integral_constant<T, order>)
395  {
396  static_assert(-(order+1) < laurent_order, "Required Laurent coefficient too high");
397  static_assert(-(order+1) >= 0, "Required Laurent coefficient too low");
398  return m_Fcoeffs[-(order+1)];
399  }
400 
402  template <class T, T order>
403  void set_laurent_coeff(std::integral_constant<T, order>, laurent_coeff_t const &v)
404  {
405  static_assert(-(order+1) < laurent_order, "Required Laurent coefficient too high");
406  static_assert(-(order+1) >= 0, "Required Laurent coefficient too low");
407  m_Fcoeffs[-(order+1)] = v;
408  }
409 
411  x_t const &get_n0(void) const
412  {
413  return m_n0;
414  }
415 
416  kernel_t const &get_kernel(void) const
417  {
418  return m_kernel.derived();
419  }
420 
421 private:
422  elem_t const &m_elem;
423  kernel_base<kernel_t> const &m_kernel;
425  xi_t m_xi0;
426  x_t m_x0;
427  x_t m_n0;
429  trans_t m_T;
430  trans_t m_Tinv;
431  xi_t m_eta0;
434  scalar_t m_theta_lim[domain_t::num_corners];
436  scalar_t m_theta0[domain_t::num_corners];
438  scalar_t m_ref_distance[domain_t::num_corners];
439 
440  x_t m_rvec_series[laurent_order];
441  x_t m_Jvec_series[laurent_order];
442  trial_n_shape_t m_N_series[laurent_order];
443  laurent_coeff_t m_Fcoeffs[laurent_order];
444 };
445 
446 
447 
453 template <class TrialField, class Kernel, unsigned RadialOrder, unsigned TangentialOrder>
454 class guiggiani<TrialField, Kernel, RadialOrder, TangentialOrder, typename std::enable_if<
455  !element_traits::is_surface_element<typename TrialField::elem_t>::value
456 >::type>
457 {
458 public:
460  enum {
461  radial_order = RadialOrder,
462  tangential_order = TangentialOrder
463  };
464 
466  typedef TrialField trial_field_t;
468  typedef typename trial_field_t::elem_t elem_t;
469 
471  typedef typename elem_t::domain_t domain_t;
473  typedef typename domain_t::xi_t xi_t;
475  typedef typename elem_t::scalar_t scalar_t;
477  typedef Eigen::Matrix<scalar_t, domain_t::dimension, domain_t::dimension> trans_t;
478 
480  typedef typename elem_t::x_t x_t;
481 
483  typedef typename trial_field_t::nset_t trial_nset_t;
485  typedef typename trial_nset_t::shape_t trial_n_shape_t;
486 
488  typedef Kernel kernel_t;
489 
494 
497 
502 
504  typedef typename semi_block_product_result_type<
505  typename singular_core_t::result_t, trial_n_shape_t
507 
509  typedef typename semi_block_product_result_type<
510  typename kernel_t::result_t, trial_n_shape_t
511  >::type total_result_t;
512 
514  enum {
517  };
518 
523  guiggiani(elem_t const &elem, kernel_base<kernel_t> const &kernel)
524  : m_elem(elem), m_kernel(kernel)
525  {
526  }
527 
528 private:
532  void compute_xi0(xi_t const &xi0)
533  {
534  m_xi0 = xi0;
535  m_x0 = m_elem.get_x(m_xi0);
536  typename elem_t::dx_return_type dx = m_elem.get_dx(m_xi0);
537 
538  scalar_t cosgamma = dx.col(shape_derivative_index::dXI).dot(dx.col(shape_derivative_index::dETA)) /
539  dx.col(shape_derivative_index::dXI).norm() / dx.col(shape_derivative_index::dETA).norm();
540  scalar_t u = dx.col(shape_derivative_index::dETA).norm() / dx.col(shape_derivative_index::dXI).norm();
541  m_T <<
542  1.0, cosgamma * u,
543  0.0, std::sqrt(1.0-cosgamma*cosgamma) * u;
544  // this scaling ensures that A is always equal to 1.0
545  m_T *= dx.col(shape_derivative_index::dXI).norm();
546  m_Tinv = m_T.inverse();
547 
548  m_eta0 = m_T * m_xi0;
549 
550  m_J_series[0] = m_elem.get_dx(m_xi0).determinant() / m_T.determinant();
551  m_N_series[0] = trial_nset_t::template eval_shape<0>(xi0);
552 
553  // geometrical parameters (planar helpers)
554  unsigned const N = domain_t::num_corners;
555  for (unsigned n = 0; n < N; ++n)
556  {
557  xi_t c1 = m_T * domain_t::get_corner(n); // corner
558  xi_t c2 = m_T * domain_t::get_corner((n + 1) % N); // next corner
559  xi_t l = (c2 - c1).normalized(); // side vector
560 
561  xi_t d1 = c1 - m_eta0; // vector to corners
562  xi_t d0 = d1 - l*d1.dot(l); // perpendicular to side
563 
564  m_theta_lim[n] = std::atan2(d1(1), d1(0)); // corner angle
565  m_theta0[n] = std::atan2(d0(1), d0(0)); // mid angle
566  m_ref_distance[n] = d0.norm(); // distance to side
567  }
568  }
569 
570  template <class Input1, class Input2>
571  typename Input1::Scalar detvectors(Eigen::DenseBase<Input1> const &v1, Eigen::DenseBase<Input2> const &v2)
572  {
573  return v1(0)*v2(1) - v1(1)*v2(0);
574  }
575 
579  void compute_theta(scalar_t theta)
580  {
581  // contains cos theta sin theta in xi system
582  xi_t xi = m_Tinv.col(shape_derivative_index::dXI) * std::cos(theta)
583  + m_Tinv.col(shape_derivative_index::dETA) * std::sin(theta);
584 
585  typename elem_t::dx_return_type dx = m_elem.get_dx(m_xi0);
586 
587  m_rvec_series[0] = dx.col(shape_derivative_index::dXI) * xi(shape_derivative_index::dXI, 0)
589  if (laurent_order > 1) // compile_time IF
590  {
591  typename elem_t::ddx_return_type ddx = m_elem.get_ddx(m_xi0);
592 
593  m_rvec_series[1] =
596  +
599  +
602 
603  m_J_series[1] = (
604  xi(0, 0) *
605  (
606  detvectors(ddx.col(shape_derivative_index::dXIXI), dx.col(shape_derivative_index::dETA)) +
607  detvectors(dx.col(shape_derivative_index::dXI), ddx.col(shape_derivative_index::dXIETA))
608  )
609  +
610  xi(1, 0) *
611  (
612  detvectors(ddx.col(shape_derivative_index::dETAXI), dx.col(shape_derivative_index::dETA))+
614  )
615  ) / m_T.determinant();
616 
617  auto dN = trial_nset_t::template eval_shape<1>(m_xi0);
618  m_N_series[1] = dN.col(shape_derivative_index::dXI)*xi(shape_derivative_index::dXI, 0)
620  }
621  }
622 
623 public:
630  template <class result_t>
631  void integrate(result_t &&I, xi_t const &xi0)
632  {
633  using namespace boost::math::double_constants;
634 
635 #if NIHU_MEX_DEBUGGING
636  static bool printed = false;
637  if (!printed)
638  {
639  mexPrintf("Computing integral with Guiggiani.\nRadial_order: %d\nTangential order: %d\n", RadialOrder, TangentialOrder);
640  printed = true;
641  }
642 #endif
643 
644  // compute the xi0-related quantities
645  compute_xi0(xi0);
646 
647  // bind the kernel at the test input
648  test_input_t test_input(m_elem, m_xi0);
649  auto bound = m_kernel.bind(test_input);
650 
651  // iterate through triangles
652  unsigned const N = domain_t::num_corners;
653  for (unsigned n = 0; n < N; ++n)
654  {
655  // get angular integration limits
656  scalar_t t1 = m_theta_lim[n];
657  scalar_t t2 = m_theta_lim[(n + 1) % N];
658 
659  // we assume that the domain's corners are listed in positive order
660  if (std::abs(t2 - t1) > pi)
661  t2 += two_pi;
662 
663  // theta integration
664  for (auto const &q_theta : line_quad_store<tangential_order>::quadrature)
665  {
666  // compute theta base point and weight
667  scalar_t xx = q_theta.get_xi()(0);
668  scalar_t theta = ((1.0 - xx) * t1 + (1.0 + xx) * t2) / 2.0;
669  scalar_t w_theta = q_theta.get_w() * (t2 - t1) / 2.0;
670 
671  // compute theta-related members
672  compute_theta(theta);
673  laurent_t::eval(*this);
674 
675  // reference domain's limit
676  scalar_t rho_lim = m_ref_distance[n] / std::cos(theta - m_theta0[n]);
677 
678  laurent_coeff_t toadd = m_Fcoeffs[0] * std::log(rho_lim);
679  if (laurent_order > 1) // compile_time IF
680  toadd -= m_Fcoeffs[1] / rho_lim;
681  toadd *= w_theta;
682 
683  for (int r = 0; r < I.rows(); ++r) // loop needed for scalar casting
684  for (int c = 0; c < I.cols(); ++c)
685  I(r,c) += toadd(r,c);
686 
687  // radial part of surface integration
688  for (auto const &q_rho : line_quad_store<radial_order>::quadrature)
689  {
690  // quadrature base point and weight
691  scalar_t rho = (1.0 + q_rho.get_xi()(0)) * rho_lim / 2.0;
692  scalar_t w_rho = q_rho.get_w() * rho_lim / 2.0;
693 
694  // compute location in original reference domain
695  xi_t eta(rho*std::cos(theta), rho*std::sin(theta));
696  xi_t xi = m_Tinv * eta + m_xi0;
697 
698  // evaluate G * N * J * rho
699  w_trial_input_t trial_input(m_elem, xi);
700  typename kernel_t::result_t GJrho =
701  bound(trial_input) * (trial_input.get_jacobian() / m_T.determinant() * rho);
702  typename trial_nset_t::shape_t N =
703  trial_nset_t::template eval_shape<0>(xi);
704  total_result_t F = semi_block_product(GJrho, N);
705 
706  // subtract the analytical singularity
707  laurent_coeff_t singular_part = m_Fcoeffs[0];
708  if (laurent_order > 1) // compile_time IF
709  singular_part += m_Fcoeffs[1] / rho;
710  singular_part /= rho;
711 
712  for (int r = 0; r < F.rows(); ++r) // loop needed for scalar casting
713  for (int c = 0; c < F.cols(); ++c)
714  F(r,c) -= singular_part(r,c);
715 
716  // surface integral accumulation
717  I += w_theta * w_rho * F;
718  } // rho loop
719  } // theta loop
720  } // element sides
721  }
722 
724  template <class T, T order>
725  x_t const &get_rvec_series(std::integral_constant<T, order>) const
726  {
727  static_assert((order-1)>=0, "Required distance Taylor coefficient too low");
728  static_assert((order-1) < laurent_order, "Required distance Taylor coefficient too high");
729  return m_rvec_series[order-1];
730  }
731 
733  template <class T, T order>
734  double get_J_series(std::integral_constant<T, order>) const
735  {
736  static_assert(order < laurent_order, "Required Jacobian Taylor coefficient too high");
737  static_assert(order >= 0, "Required Jacobian Taylor coefficient too low");
738  return m_J_series[order];
739  }
740 
742  template <class T, T order>
743  trial_n_shape_t const &get_shape_series(std::integral_constant<T, order>) const
744  {
745  static_assert(order < laurent_order, "Required shape set Taylor coefficient too high");
746  static_assert(order >= 0, "Required Jacobian Taylor coefficient too low");
747  return m_N_series[order];
748  }
749 
751  template <class T, T order>
752  laurent_coeff_t &get_laurent_coeff(std::integral_constant<T, order>)
753  {
754  static_assert(-(order+1) < laurent_order, "Required Laurent coefficient too high");
755  static_assert(-(order+1) >= 0, "Required Laurent coefficient too low");
756  return m_Fcoeffs[-(order+1)];
757  }
758 
760  template <class T, T order>
761  void set_laurent_coeff(std::integral_constant<T, order>, laurent_coeff_t const &v)
762  {
763  static_assert(-(order+1) < laurent_order, "Required Laurent coefficient too high");
764  static_assert(-(order+1) >= 0, "Required Laurent coefficient too low");
765  m_Fcoeffs[-(order+1)] = v;
766  }
767 
768 
769  kernel_t const &get_kernel(void) const
770  {
771  return m_kernel.derived();
772  }
773 
774 
775 private:
776  elem_t const &m_elem;
777  kernel_base<kernel_t> const &m_kernel;
779  xi_t m_xi0;
780  x_t m_x0;
782  trans_t m_T;
783  trans_t m_Tinv;
784  xi_t m_eta0;
787  scalar_t m_theta_lim[domain_t::num_corners];
789  scalar_t m_theta0[domain_t::num_corners];
791  scalar_t m_ref_distance[domain_t::num_corners];
792 
793  x_t m_rvec_series[laurent_order];
794  scalar_t m_J_series[laurent_order];
795  trial_n_shape_t m_N_series[laurent_order];
796  laurent_coeff_t m_Fcoeffs[laurent_order];
797 };
798 
799 }
800 
801 
802 
803 #endif // GUIGGIANI_1992_HPP_INCLUDED
804 
NiHu::guiggiani< TrialField, Kernel, RadialOrder, TangentialOrder, typename std::enable_if< !element_traits::is_surface_element< typename TrialField::elem_t >::value >::type >::xi_t
domain_t::xi_t xi_t
the reference coordinate vector type
Definition: guiggiani_1992.hpp:473
NiHu::shape_derivative_index::dETAXI
@ dETAXI
index of eta_xi in 2nd derivative matrix
Definition: shapeset.hpp:53
NiHu::guiggiani< TrialField, Kernel, RadialOrder, TangentialOrder, typename std::enable_if< !element_traits::is_surface_element< typename TrialField::elem_t >::value >::type >::total_result_t
semi_block_product_result_type< typename kernel_t::result_t, trial_n_shape_t >::type total_result_t
value type of the integral result
Definition: guiggiani_1992.hpp:511
NiHu::guiggiani< TrialField, Kernel, RadialOrder, TangentialOrder, typename std::enable_if< !element_traits::is_surface_element< typename TrialField::elem_t >::value >::type >::laurent_coeff_t
semi_block_product_result_type< typename singular_core_t::result_t, trial_n_shape_t >::type laurent_coeff_t
value type of the Laurent coefficients
Definition: guiggiani_1992.hpp:506
NiHu::guiggiani< TrialField, Kernel, RadialOrder, TangentialOrder, typename std::enable_if< !element_traits::is_surface_element< typename TrialField::elem_t >::value >::type >::get_J_series
double get_J_series(std::integral_constant< T, order >) const
return Taylor coefficient of the Jacobian vector around the collocation point
Definition: guiggiani_1992.hpp:734
NiHu::guiggiani< TrialField, Kernel, RadialOrder, TangentialOrder, typename std::enable_if< element_traits::is_surface_element< typename TrialField::elem_t >::value >::type >::laurent_coeff_t
semi_block_product_result_type< typename singular_core_t::result_t, trial_n_shape_t >::type laurent_coeff_t
value type of the Laurent coefficients
Definition: guiggiani_1992.hpp:141
NiHu::volume_element
class describing a volume element that has no normal vector
Definition: element.hpp:455
NiHu::guiggiani< TrialField, Kernel, RadialOrder, TangentialOrder, typename std::enable_if< !element_traits::is_surface_element< typename TrialField::elem_t >::value >::type >::x_t
elem_t::x_t x_t
the physical coordinate vector type
Definition: guiggiani_1992.hpp:480
NiHu::singular_kernel_traits::singular_core_t
kernel_traits_ns::singular_core< Derived >::type singular_core_t
the kernel's singular core type
Definition: kernel.hpp:106
NiHu::guiggiani< TrialField, Kernel, RadialOrder, TangentialOrder, typename std::enable_if< !element_traits::is_surface_element< typename TrialField::elem_t >::value >::type >::trans_t
Eigen::Matrix< scalar_t, domain_t::dimension, domain_t::dimension > trans_t
the Rong's transformation matrix type
Definition: guiggiani_1992.hpp:477
NiHu::guiggiani< TrialField, Kernel, RadialOrder, TangentialOrder, typename std::enable_if< !element_traits::is_surface_element< typename TrialField::elem_t >::value >::type >::w_trial_input_t
weighted_input< trial_input_t, elem_t >::type w_trial_input_t
the kernel's weighted trial input type
Definition: guiggiani_1992.hpp:496
NiHu::guiggiani< TrialField, Kernel, RadialOrder, TangentialOrder, typename std::enable_if< !element_traits::is_surface_element< typename TrialField::elem_t >::value >::type >::elem_t
trial_field_t::elem_t elem_t
the element type
Definition: guiggiani_1992.hpp:468
NiHu::guiggiani< TrialField, Kernel, RadialOrder, TangentialOrder, typename std::enable_if< element_traits::is_surface_element< typename TrialField::elem_t >::value >::type >::set_laurent_coeff
void set_laurent_coeff(std::integral_constant< T, order >, laurent_coeff_t const &v)
set a Laurent coefficient
Definition: guiggiani_1992.hpp:403
NiHu::shape_derivative_index::dXIXI
@ dXIXI
index of xi_xi in 2nd derivative matrix
Definition: shapeset.hpp:51
NiHu::guiggiani< TrialField, Kernel, RadialOrder, TangentialOrder, typename std::enable_if< element_traits::is_surface_element< typename TrialField::elem_t >::value >::type >::get_shape_series
const trial_n_shape_t & get_shape_series(std::integral_constant< T, order >) const
return Taylor coefficient of the shape set around the collocation point
Definition: guiggiani_1992.hpp:385
NiHu::guiggiani< TrialField, Kernel, RadialOrder, TangentialOrder, typename std::enable_if< element_traits::is_surface_element< typename TrialField::elem_t >::value >::type >::xi_t
domain_t::xi_t xi_t
the reference coordinate vector type
Definition: guiggiani_1992.hpp:108
NiHu::guiggiani< TrialField, Kernel, RadialOrder, TangentialOrder, typename std::enable_if< element_traits::is_surface_element< typename TrialField::elem_t >::value >::type >::trial_field_t
TrialField trial_field_t
the trial field type
Definition: guiggiani_1992.hpp:101
NiHu::guiggiani< TrialField, Kernel, RadialOrder, TangentialOrder, typename std::enable_if< !element_traits::is_surface_element< typename TrialField::elem_t >::value >::type >::get_rvec_series
const x_t & get_rvec_series(std::integral_constant< T, order >) const
return Taylor coefficient of the distance measured from the collocation point
Definition: guiggiani_1992.hpp:725
NiHu::line_quad_store
store-wrapper of a statically stored line quadrature
Definition: line_quad_store.hpp:11
NiHu::volume_element::x_t
element_traits::location_value_type< Derived, 0 >::type x_t
type of the element's physical location variable
Definition: element.hpp:220
NiHu::_m2
std::integral_constant< int, -2 > _m2
shorthand for constant minus two
Definition: guiggiani_1992.hpp:64
NiHu::singular_kernel_traits
singular traits class of a kernel
Definition: kernel.hpp:101
NiHu::guiggiani< TrialField, Kernel, RadialOrder, TangentialOrder, typename std::enable_if< element_traits::is_surface_element< typename TrialField::elem_t >::value >::type >::trial_input_t
kernel_traits< kernel_t >::trial_input_t trial_input_t
the kernel's trial input type
Definition: guiggiani_1992.hpp:128
NiHu::_2
std::integral_constant< int, 2 > _2
shorthand for constant two
Definition: guiggiani_1992.hpp:75
NiHu::guiggiani< TrialField, Kernel, RadialOrder, TangentialOrder, typename std::enable_if< element_traits::is_surface_element< typename TrialField::elem_t >::value >::type >::laurent_t
polar_laurent_coeffs< singular_core_t > laurent_t
the Laurent coefficients computing class
Definition: guiggiani_1992.hpp:136
NiHu::semi_block_product
semi_block_product_result_type< mat, rightDerived >::type semi_block_product(mat const &m, Eigen::MatrixBase< rightDerived > const &r)
compute semi block product of a matrix and a vector m * v^T
Definition: block_product.hpp:194
NiHu::shape_derivative_index::dETAETA
@ dETAETA
index of eta_eta in 2nd derivative matrix
Definition: shapeset.hpp:54
NiHu::guiggiani< TrialField, Kernel, RadialOrder, TangentialOrder, typename std::enable_if< !element_traits::is_surface_element< typename TrialField::elem_t >::value >::type >::get_shape_series
const trial_n_shape_t & get_shape_series(std::integral_constant< T, order >) const
return Taylor coefficient of the shape set around the collocation point
Definition: guiggiani_1992.hpp:743
NiHu::guiggiani< TrialField, Kernel, RadialOrder, TangentialOrder, typename std::enable_if< element_traits::is_surface_element< typename TrialField::elem_t >::value >::type >::w_trial_input_t
weighted_input< trial_input_t, elem_t >::type w_trial_input_t
the kernel's weighted trial input type
Definition: guiggiani_1992.hpp:131
NiHu::exponential_covariance_kernel
Definition: covariance_kernel.hpp:42
NiHu::guiggiani< TrialField, Kernel, RadialOrder, TangentialOrder, typename std::enable_if< element_traits::is_surface_element< typename TrialField::elem_t >::value >::type >::kernel_t
Kernel kernel_t
the kernel type
Definition: guiggiani_1992.hpp:123
NiHu::guiggiani< TrialField, Kernel, RadialOrder, TangentialOrder, typename std::enable_if< !element_traits::is_surface_element< typename TrialField::elem_t >::value >::type >::guiggiani
guiggiani(elem_t const &elem, kernel_base< kernel_t > const &kernel)
constructor
Definition: guiggiani_1992.hpp:523
NiHu::guiggiani< TrialField, Kernel, RadialOrder, TangentialOrder, typename std::enable_if< !element_traits::is_surface_element< typename TrialField::elem_t >::value >::type >::set_laurent_coeff
void set_laurent_coeff(std::integral_constant< T, order >, laurent_coeff_t const &v)
set a Laurent coefficient
Definition: guiggiani_1992.hpp:761
NiHu::kernel_base< kernel_t >
NiHu::_m1
std::integral_constant< int, -1 > _m1
shorthand for constant minus one
Definition: guiggiani_1992.hpp:69
NiHu::kernel_traits
traits class of a kernel
Definition: kernel.hpp:71
NiHu::guiggiani< TrialField, Kernel, RadialOrder, TangentialOrder, typename std::enable_if< element_traits::is_surface_element< typename TrialField::elem_t >::value >::type >::scalar_t
elem_t::scalar_t scalar_t
the geometrical scalar type
Definition: guiggiani_1992.hpp:110
NiHu::guiggiani< TrialField, Kernel, RadialOrder, TangentialOrder, typename std::enable_if< element_traits::is_surface_element< typename TrialField::elem_t >::value >::type >::trial_n_shape_t
trial_nset_t::shape_t trial_n_shape_t
the shape function vector type
Definition: guiggiani_1992.hpp:120
NiHu::guiggiani< TrialField, Kernel, RadialOrder, TangentialOrder, typename std::enable_if< element_traits::is_surface_element< typename TrialField::elem_t >::value >::type >::trial_nset_t
trial_field_t::nset_t trial_nset_t
shape function set type
Definition: guiggiani_1992.hpp:118
NiHu::guiggiani< TrialField, Kernel, RadialOrder, TangentialOrder, typename std::enable_if< !element_traits::is_surface_element< typename TrialField::elem_t >::value >::type >::trial_n_shape_t
trial_nset_t::shape_t trial_n_shape_t
the shape function vector type
Definition: guiggiani_1992.hpp:485
NiHu::guiggiani< TrialField, Kernel, RadialOrder, TangentialOrder, typename std::enable_if< !element_traits::is_surface_element< typename TrialField::elem_t >::value >::type >::integrate
void integrate(result_t &&I, xi_t const &xi0)
the entry point of integration
Definition: guiggiani_1992.hpp:631
NiHu::guiggiani< TrialField, Kernel, RadialOrder, TangentialOrder, typename std::enable_if< !element_traits::is_surface_element< typename TrialField::elem_t >::value >::type >::kernel_t
Kernel kernel_t
the kernel type
Definition: guiggiani_1992.hpp:488
NiHu::_1
std::integral_constant< int, 1 > _1
shorthand for constant one
Definition: guiggiani_1992.hpp:73
NiHu::guiggiani< TrialField, Kernel, RadialOrder, TangentialOrder, typename std::enable_if< element_traits::is_surface_element< typename TrialField::elem_t >::value >::type >::get_laurent_coeff
laurent_coeff_t & get_laurent_coeff(std::integral_constant< T, order >)
return a Laurent coefficient
Definition: guiggiani_1992.hpp:394
NiHu::guiggiani< TrialField, Kernel, RadialOrder, TangentialOrder, typename std::enable_if< element_traits::is_surface_element< typename TrialField::elem_t >::value >::type >::guiggiani
guiggiani(elem_t const &elem, kernel_base< kernel_t > const &kernel)
constructor
Definition: guiggiani_1992.hpp:158
NiHu::semi_block_product_result_type
metafunction returning the value type of a semi block product
Definition: block_product.hpp:180
NiHu::guiggiani< TrialField, Kernel, RadialOrder, TangentialOrder, typename std::enable_if< element_traits::is_surface_element< typename TrialField::elem_t >::value >::type >::singular_core_t
singular_kernel_traits< kernel_t >::singular_core_t singular_core_t
the singular core type
Definition: guiggiani_1992.hpp:134
NiHu::guiggiani< TrialField, Kernel, RadialOrder, TangentialOrder, typename std::enable_if< !element_traits::is_surface_element< typename TrialField::elem_t >::value >::type >::domain_t
elem_t::domain_t domain_t
the original reference domain type
Definition: guiggiani_1992.hpp:471
NiHu::guiggiani< TrialField, Kernel, RadialOrder, TangentialOrder, typename std::enable_if< !element_traits::is_surface_element< typename TrialField::elem_t >::value >::type >::trial_input_t
kernel_traits< kernel_t >::trial_input_t trial_input_t
the kernel's trial input type
Definition: guiggiani_1992.hpp:493
NiHu::guiggiani< TrialField, Kernel, RadialOrder, TangentialOrder, typename std::enable_if< !element_traits::is_surface_element< typename TrialField::elem_t >::value >::type >::trial_nset_t
trial_field_t::nset_t trial_nset_t
shape function set type
Definition: guiggiani_1992.hpp:483
NiHu::guiggiani< TrialField, Kernel, RadialOrder, TangentialOrder, typename std::enable_if< element_traits::is_surface_element< typename TrialField::elem_t >::value >::type >::total_result_t
semi_block_product_result_type< typename kernel_t::result_t, trial_n_shape_t >::type total_result_t
value type of the integral result
Definition: guiggiani_1992.hpp:146
NiHu::shape_derivative_index::dXI
@ dXI
index of xi in 1st derivative matrix
Definition: shapeset.hpp:49
NiHu::guiggiani< TrialField, Kernel, RadialOrder, TangentialOrder, typename std::enable_if< !element_traits::is_surface_element< typename TrialField::elem_t >::value >::type >::laurent_t
polar_laurent_coeffs< singular_core_t > laurent_t
the Laurent coefficients computing class
Definition: guiggiani_1992.hpp:501
NiHu::guiggiani< TrialField, Kernel, RadialOrder, TangentialOrder, typename std::enable_if< element_traits::is_surface_element< typename TrialField::elem_t >::value >::type >::test_input_t
kernel_traits< kernel_t >::test_input_t test_input_t
the kernel's test input type
Definition: guiggiani_1992.hpp:126
NiHu::guiggiani< TrialField, Kernel, RadialOrder, TangentialOrder, typename std::enable_if< element_traits::is_surface_element< typename TrialField::elem_t >::value >::type >::get_rvec_series
const x_t & get_rvec_series(std::integral_constant< T, order >) const
return Taylor coefficient of the distance measured from the collocation point
Definition: guiggiani_1992.hpp:367
NiHu::volume_element::domain_t
typename base_t::domain_t domain_t
the domain type
Definition: element.hpp:637
location_normal.hpp
implementation of location and normal kernel inputs
NiHu::shape_derivative_index::dXIETA
@ dXIETA
index of xi_eta in 2nd derivative matrix
Definition: shapeset.hpp:52
NiHu::shape_derivative_index::dETA
@ dETA
index of eta in 1st derivative matrix
Definition: shapeset.hpp:50
NiHu::guiggiani< TrialField, Kernel, RadialOrder, TangentialOrder, typename std::enable_if< element_traits::is_surface_element< typename TrialField::elem_t >::value >::type >::elem_t
trial_field_t::elem_t elem_t
the element type
Definition: guiggiani_1992.hpp:103
NiHu::guiggiani< TrialField, Kernel, RadialOrder, TangentialOrder, typename std::enable_if< !element_traits::is_surface_element< typename TrialField::elem_t >::value >::type >::test_input_t
kernel_traits< kernel_t >::test_input_t test_input_t
the kernel's test input type
Definition: guiggiani_1992.hpp:491
NiHu::guiggiani< TrialField, Kernel, RadialOrder, TangentialOrder, typename std::enable_if< !element_traits::is_surface_element< typename TrialField::elem_t >::value >::type >::get_laurent_coeff
laurent_coeff_t & get_laurent_coeff(std::integral_constant< T, order >)
return a Laurent coefficient
Definition: guiggiani_1992.hpp:752
NiHu::guiggiani< TrialField, Kernel, RadialOrder, TangentialOrder, typename std::enable_if< !element_traits::is_surface_element< typename TrialField::elem_t >::value >::type >::scalar_t
elem_t::scalar_t scalar_t
the geometrical scalar type
Definition: guiggiani_1992.hpp:475
NiHu::guiggiani< TrialField, Kernel, RadialOrder, TangentialOrder, typename std::enable_if< !element_traits::is_surface_element< typename TrialField::elem_t >::value >::type >::trial_field_t
TrialField trial_field_t
the trial field type
Definition: guiggiani_1992.hpp:466
NiHu::polar_laurent_coeffs
definition of Laurent coefficients of singularities
Definition: guiggiani_1992.hpp:64
NiHu::_0
std::integral_constant< int, 0 > _0
shorthand for constant zero
Definition: guiggiani_1992.hpp:71
NiHu::guiggiani< TrialField, Kernel, RadialOrder, TangentialOrder, typename std::enable_if< element_traits::is_surface_element< typename TrialField::elem_t >::value >::type >::x_t
elem_t::x_t x_t
the physical coordinate vector type
Definition: guiggiani_1992.hpp:115
NiHu::guiggiani< TrialField, Kernel, RadialOrder, TangentialOrder, typename std::enable_if< element_traits::is_surface_element< typename TrialField::elem_t >::value >::type >::integrate
void integrate(result_t &&I, xi_t const &xi0, x_t const &normal)
the entry point of integration
Definition: guiggiani_1992.hpp:262
NiHu::guiggiani< TrialField, Kernel, RadialOrder, TangentialOrder, typename std::enable_if< element_traits::is_surface_element< typename TrialField::elem_t >::value >::type >::domain_t
elem_t::domain_t domain_t
the original reference domain type
Definition: guiggiani_1992.hpp:106
NiHu::guiggiani< TrialField, Kernel, RadialOrder, TangentialOrder, typename std::enable_if< element_traits::is_surface_element< typename TrialField::elem_t >::value >::type >::get_n0
const x_t & get_n0(void) const
return the unit normal at the collocation point
Definition: guiggiani_1992.hpp:411
weighted_input.hpp
implementation of metafunction NiHu::weighted_input
NiHu::guiggiani
Implementation of Guiggiani's method.
Definition: guiggiani_1992.hpp:84
NiHu::weighted_input
the weigthed input is an extended input that contains the jacobian a well
Definition: weighted_input.hpp:19
NiHu::guiggiani< TrialField, Kernel, RadialOrder, TangentialOrder, typename std::enable_if< !element_traits::is_surface_element< typename TrialField::elem_t >::value >::type >::singular_core_t
singular_kernel_traits< kernel_t >::singular_core_t singular_core_t
the singular core type
Definition: guiggiani_1992.hpp:499
NiHu::guiggiani< TrialField, Kernel, RadialOrder, TangentialOrder, typename std::enable_if< element_traits::is_surface_element< typename TrialField::elem_t >::value >::type >::get_Jvec_series
const x_t & get_Jvec_series(std::integral_constant< T, order >) const
return Taylor coefficient of the Jacobian vector around the collocation point
Definition: guiggiani_1992.hpp:376
NiHu::guiggiani< TrialField, Kernel, RadialOrder, TangentialOrder, typename std::enable_if< element_traits::is_surface_element< typename TrialField::elem_t >::value >::type >::trans_t
Eigen::Matrix< scalar_t, domain_t::dimension, domain_t::dimension > trans_t
the Rong's transformation matrix type
Definition: guiggiani_1992.hpp:112