NiHu  2.0
laplace_3d_transparent_standalone.cpp
1 #include <boost/math/constants/constants.hpp>
2 
4 #include "library/laplace_3d.hpp"
5 #include "../library/lib_element.hpp"
6 #include "../interface/read_off_mesh.hpp"
7 
8 #include <iostream>
9 
10 typedef Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic> dMatrix;
11 typedef Eigen::Matrix<unsigned, Eigen::Dynamic, Eigen::Dynamic> uMatrix;
12 
13 template <class TestSpace, class TrialSpace>
14 static double tester(TestSpace const &test_space, TrialSpace const &trial_space)
15 {
16  using namespace boost::math::double_constants;
17 
18  // instantiate integral operators
22 
23  size_t N = test_space.get_num_dofs();
24 
25  // excitation and response
26  dMatrix q0(N, 1), p0(N,1);
27  p0.setZero();
28  q0.setZero();
29 
30  dMatrix x0(3, 2);
31  x0 <<
32  .2, .3,
33  .3, .2,
34  .1, .2;
35  dMatrix A(1, 2);
36  A << 1.0, -1.0;
37 
38  auto const &mesh = trial_space.get_mesh();
39  for (unsigned k = 0; k < N; ++k)
40  {
41  auto const &elem = mesh.template get_elem<NiHu::tria_1_elem>(k);
42  auto y = elem.get_center();
43  auto ny = elem.get_normal().normalized();
44  for (int s = 0; s < x0.cols(); ++s)
45  {
46  auto rvec = y - x0.col(s);
47  double r = rvec.norm();
48  double rdn = rvec.dot(ny) / r;
49  p0(k) += A(s) * 1./r/(4.0 * pi);
50  q0(k) += A(s) * -(1./r/r) * rdn /(4.0 * pi);
51  }
52  }
53 
54  // system matrices
55  dMatrix L(N, N), M(N, N), I(N,N);
56  L.setZero(); M.setZero(); I.setZero();
57 
58  I << test_space * I_op[trial_space];
59  L << test_space * L_op[trial_space];
60  M << test_space * M_op[trial_space];
61 
62  std::cout << I.block(0, 0, 5, 5) << std::endl << std::endl;
63  std::cout << L.block(0, 0, 5, 5) << std::endl << std::endl;
64  std::cout << M.block(0, 0, 5, 5) << std::endl << std::endl;
65 
66  dMatrix p = (M - .5 * I).colPivHouseholderQr().solve(L * q0);
67 
68  return (p-p0).norm() / p0.norm();
69 }
70 
71 int main(int argc, char *argv[])
72 {
73  if (argc < 2)
74  {
75  std::cerr << "Usage: " << argv[0] << " OFFname" << std::endl;
76  return 1;
77  }
78 
79  // read mesh
80  auto mesh = NiHu::read_off_mesh(argv[1], NiHu::tria_1_tag());
81  auto const &trial_space = NiHu::constant_view(mesh);
82 
83  std::cout << "Collocation" << std::endl;
84  double coll_error = tester(NiHu::dirac(trial_space), trial_space);
85  std::cout << coll_error << std::endl;
86 
87  std::cout << "Galerkin" << std::endl;
88  double gal_error = tester(trial_space, trial_space);
89  std::cout << gal_error << std::endl;
90 
91  return 0;
92 }
93 
NiHu::create_integral_operator
integral_operator< Kernel > create_integral_operator(Kernel &&kernel)
factory function of an integral operator
Definition: integral_operator.hpp:421
NiHu::read_off_mesh
mesh< tmp::vector< typename tag2type< Tags >::type... > > read_off_mesh(std::string const &fname, Tags...tags)
Read mesh from OFF format.
Definition: read_off_mesh.hpp:153
NiHu::type2tag
Metafunction assigning a tag to a type.
Definition: type2tag.hpp:17
NiHu::constant_view
const function_space_view< Mesh, field_option::constant, Dimension > & constant_view(mesh_base< Mesh > const &msh, Dimension dim=Dimension())
factory function to transform a mesh into a constant function space
Definition: function_space.hpp:304
weighted_residual.hpp
declaration of class NiHu::weighted_residual
NiHu::mex::real_matrix
Definition: mex_matrix_separate.hpp:133
NiHu::normal_derivative_kernel
Normal derivative of a distance dependent kernel.
Definition: normal_derivative_kernel.hpp:26
NiHu::identity_integral_operator
The identity integral operator .
Definition: integral_operator.hpp:310