1#ifndef HLDO_SPHERE_OPERATORS_DIRAC_OPERATOR_SOURCE_PROBLEM_H
2#define HLDO_SPHERE_OPERATORS_DIRAC_OPERATOR_SOURCE_PROBLEM_H
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>
37using complex = std::complex<double>;
109 void SetMesh(std::shared_ptr<const lf::mesh::Mesh> mesh_p);
119 std::function<
complex(
const Eigen::Matrix<double, 3, 1> &)> f0,
121 Eigen::Matrix<complex, 3, 1>(
const Eigen::Matrix<double, 3, 1> &)>
123 std::function<
complex(
const Eigen::Matrix<double, 3, 1> &)> f2) {
134 void SetK(
double k) {
k_ = std::complex<double>(k, 0); }
184 Eigen::Matrix<double, Eigen::Dynamic, 1>
GetMu(
int index);
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_;
Computes the Galerkin LSE for the Dirac Operator source problem.
lf::assemble::COOMatrix< complex > coo_matrix_
void Compute()
Computes the Galerkin LSE for the dirac operator source problem.
std::function< complex(const Eigen::Matrix< double, 3, 1 > &)> f2_
Eigen::Matrix< complex, Eigen::Dynamic, 1 > mu_
Eigen::Matrix< double, Eigen::Dynamic, 1 > GetMu(int index)
retunrs the basis expansion coefficiants for the solution of the source problem
Eigen::Matrix< complex, Eigen::Dynamic, 1 > phi_
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)
void Solve()
solves the linear systems build in Compute
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.
std::shared_ptr< const lf::mesh::Mesh > mesh_p_
void SetK(double k)
Sets k for the souce problems.
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.