LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
dirac_operator_source_problem.h
1#ifndef HLDO_SPHERE_OPERATORS_DIRAC_OPERATOR_SOURCE_PROBLEM_H
2#define HLDO_SPHERE_OPERATORS_DIRAC_OPERATOR_SOURCE_PROBLEM_H
3
10#include <dirac_operator.h>
11#include <lf/assemble/assembler.h>
12#include <lf/assemble/coomatrix.h>
13#include <lf/assemble/dofhandler.h>
14#include <lf/mesh/entity.h>
15#include <lf/mesh/hybrid2d/mesh.h>
16#include <lf/mesh/hybrid2d/mesh_factory.h>
17#include <lf/mesh/mesh_interface.h>
18#include <lf/quad/quad.h>
19#include <load_vector_provider.h>
20#include <mass_matrix_provider.h>
21#include <sphere_triag_mesh_builder.h>
22#include <whitney_one_curl_curl_matrix_provider.h>
23#include <whitney_one_grad_matrix_provider.h>
24#include <whitney_one_mass_matrix_provider.h>
25#include <whitney_two_mass_matrix_provider.h>
26
27#include <Eigen/Dense>
28#include <cmath>
29#include <complex>
30#include <iomanip>
31#include <iostream>
32
33namespace projects::hldo_sphere {
34
35namespace operators {
36
37using complex = std::complex<double>;
38
57 public:
69
81 void Compute();
82
100 void Solve();
101
109 void SetMesh(std::shared_ptr<const lf::mesh::Mesh> mesh_p);
110
119 std::function<complex(const Eigen::Matrix<double, 3, 1> &)> f0,
120 std::function<
121 Eigen::Matrix<complex, 3, 1>(const Eigen::Matrix<double, 3, 1> &)>
122 f1,
123 std::function<complex(const Eigen::Matrix<double, 3, 1> &)> f2) {
124 f0_ = f0;
125 f1_ = f1;
126 f2_ = f2;
127 }
128
134 void SetK(double k) { k_ = std::complex<double>(k, 0); }
135
145 Eigen::Matrix<complex, Eigen::Dynamic, 1> GetLoadVector() { return phi_; }
146
157
184 Eigen::Matrix<double, Eigen::Dynamic, 1> GetMu(int index);
185
186 private:
187 std::shared_ptr<const lf::mesh::Mesh> mesh_p_;
188 std::function<complex(const Eigen::Matrix<double, 3, 1> &)> f0_;
189 std::function<Eigen::Matrix<complex, 3, 1>(
190 const Eigen::Matrix<double, 3, 1> &)>
192 std::function<complex(const Eigen::Matrix<double, 3, 1> &)> f2_;
194 Eigen::Matrix<complex, Eigen::Dynamic, 1> phi_;
195 Eigen::Matrix<complex, Eigen::Dynamic, 1> mu_;
197};
198
199} // namespace operators
200
201} // namespace projects::hldo_sphere
202
203#endif // HLDO_SPHERE_OPERATORS_DIRAC_OPERATOR_SOURCE_PROBLEM_H
Computes the Galerkin LSE for the Dirac Operator source problem.
void Compute()
Computes the Galerkin LSE for the dirac operator source problem.
std::function< complex(const Eigen::Matrix< double, 3, 1 > &)> f2_
Eigen::Matrix< double, Eigen::Dynamic, 1 > GetMu(int index)
retunrs the basis expansion coefficiants for the solution of the source problem
lf::assemble::COOMatrix< complex > GetGalerkinMatrix()
returns the Galerkin Matrix
Eigen::Matrix< complex, Eigen::Dynamic, 1 > GetLoadVector()
returns the Loadvector
DiracOperatorSourceProblem()
Constructor creates basic mesh (Octaeder with radius 1.0)
std::function< Eigen::Matrix< complex, 3, 1 >(const Eigen::Matrix< double, 3, 1 > &)> f1_
std::function< complex(const Eigen::Matrix< double, 3, 1 > &)> f0_
void SetLoadFunctions(std::function< complex(const Eigen::Matrix< double, 3, 1 > &)> f0, std::function< Eigen::Matrix< complex, 3, 1 >(const Eigen::Matrix< double, 3, 1 > &)> f1, std::function< complex(const Eigen::Matrix< double, 3, 1 > &)> f2)
Sets the load functions.
void SetMesh(std::shared_ptr< const lf::mesh::Mesh > mesh_p)
Sets the mesh.
std::complex< double > complex
Implementation of the thesis Hogde Laplacians and Dirac Operators on the surface of the 3-Sphere.
Definition: assemble.h:15