LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
hodge_laplacian_experiment.h
1#ifndef THESIS_EXPERIMENTS_HODGE_LAPLACIAN_H
2#define THESIS_EXPERIMENTS_HODGE_LAPLACIAN_H
3
11#include <hodge_laplacians_source_problems.h>
12#include <lf/io/vtk_writer.h>
13#include <lf/mesh/hybrid2d/mesh_factory.h>
14#include <lf/mesh/utils/tp_triag_mesh_builder.h>
15#include <lf/mesh/utils/utils.h>
16#include <lf/quad/quad.h>
17#include <lf/refinement/refinement.h>
18#include <lf/uscalfe/uscalfe.h>
19#include <mesh_function_whitney_one.h>
20#include <mesh_function_whitney_two.h>
21#include <mesh_function_whitney_zero.h>
22#include <norms.h>
23#include <results_processing.h>
24#include <sphere_triag_mesh_builder.h>
25
26#include <cstdlib>
27#include <cstring>
28#include <iostream>
29#include <string>
30#include <vector>
31
33
43template <typename SCALAR>
45 public:
63 std::function<SCALAR(const Eigen::Matrix<double, 3, 1> &)> u_zero,
64 std::function<
65 Eigen::Matrix<SCALAR, 3, 1>(const Eigen::Matrix<double, 3, 1> &)>
66 u_one,
67 std::function<SCALAR(const Eigen::Matrix<double, 3, 1> &)> u_two,
68 std::function<SCALAR(const Eigen::Matrix<double, 3, 1> &)> f_zero,
69 std::function<
70 Eigen::Matrix<SCALAR, 3, 1>(const Eigen::Matrix<double, 3, 1> &)>
71 f_one,
72 std::function<SCALAR(const Eigen::Matrix<double, 3, 1> &)> f_two,
73 double &k, std::string name)
74 : u_zero_(u_zero),
75 u_one_(u_one),
76 u_two_(u_two),
77 f_zero_(f_zero),
78 f_one_(f_one),
79 f_two_(f_two),
80 k_(k),
81 name_(name) {}
82
92 void Compute(std::vector<unsigned> refinement_levels,
93 std::vector<double> ks) {
94 int nl = refinement_levels.size();
95 int nk = ks.size();
96
97 // Initialize solution wrapper
99 solutions;
100 solutions.k = ks;
101 solutions.levels = refinement_levels;
102 solutions.mesh = std::vector<std::shared_ptr<const lf::mesh::Mesh>>(nl);
103 solutions.solutions = std::vector<
105
107 lse_builder;
108
109 // functions are passed by reference hence changing the k still influences
110 // the functions
111 lse_builder.SetLoadFunctions(f_zero_, f_one_, f_two_);
112
113 // Read the mesh from the gmsh file
114 std::unique_ptr<lf::mesh::MeshFactory> factory =
115 std::make_unique<lf::mesh::hybrid2d::MeshFactory>(3);
116
119
120 // we always take the same radius
121 sphere.setRadius(1.);
122
123 // start timer for total time
124 std::chrono::steady_clock::time_point start_time_total =
125 std::chrono::steady_clock::now();
126
127 for (unsigned il = 0; il < nl; ++il) {
128 unsigned lvl = refinement_levels[il];
129 // Initialize solution vectors
130 solutions.solutions[il].mu_zero.resize(nk);
131 solutions.solutions[il].mu_one.resize(nk);
132 solutions.solutions[il].mu_two.resize(nk);
133
134 for (int ik = 0; ik < nk; ik++) {
135 k_ = ks[ik];
136 std::cout << "\nStart computation of refinement_level " << lvl
137 << " and k = " << k_ << std::flush;
138
139 // start timer
140 std::chrono::steady_clock::time_point start_time =
141 std::chrono::steady_clock::now();
142
143 // get mesh
144 sphere.setRefinementLevel(lvl);
145 const std::shared_ptr<lf::mesh::Mesh> mesh = sphere.Build();
146
147 // end timer
148 std::chrono::steady_clock::time_point end_time =
149 std::chrono::steady_clock::now();
150 double elapsed_sec =
151 std::chrono::duration_cast<std::chrono::milliseconds>(end_time -
152 start_time)
153 .count() /
154 1000.;
155
156 std::cout << " -> Built Mesh " << elapsed_sec << " [s] " << std::flush;
157
158 solutions.mesh[lvl] = mesh;
159
160 // setup the problem
161 lse_builder.SetMesh(mesh);
162 lse_builder.SetK(k_);
163
164 // start timer
165 start_time = std::chrono::steady_clock::now();
166
167 // compute the system
168 lse_builder.Compute();
169
170 // end timer
171 end_time = std::chrono::steady_clock::now();
172 elapsed_sec = std::chrono::duration_cast<std::chrono::milliseconds>(
173 end_time - start_time)
174 .count() /
175 1000.;
176
177 std::cout << " -> Computed LSE " << elapsed_sec << " [s] "
178 << std::flush;
179
180 // start timer
181 start_time = std::chrono::steady_clock::now();
182
183 // solve the system
184 lse_builder.Solve();
185
186 // end timer
187 end_time = std::chrono::steady_clock::now();
188
189 elapsed_sec = std::chrono::duration_cast<std::chrono::milliseconds>(
190 end_time - start_time)
191 .count() /
192 1000.;
193
194 std::cout << " -> Solved System " << elapsed_sec << " [s] "
195 << std::flush;
196
197 // store solutions
198 solutions.solutions[lvl].mu_zero[ik] = lse_builder.GetMuZero();
199 solutions.solutions[lvl].mu_one[ik] =
200 std::get<0>(lse_builder.GetMuOne());
201 solutions.solutions[lvl].mu_two[ik] =
202 std::get<1>(lse_builder.GetMuTwo());
203 } // end loop k
204 } // end loop level
205
206 // end timer total
207 std::chrono::steady_clock::time_point end_time_total =
208 std::chrono::steady_clock::now();
209 double elapsed_sec_total =
210 std::chrono::duration_cast<std::chrono::milliseconds>(end_time_total -
211 start_time_total)
212 .count() /
213 1000.;
214
215 std::cout << "\nTotal computation time for all levels " << elapsed_sec_total
216 << " [s]\n";
217
219 decltype(u_zero_), decltype(u_one_), decltype(u_two_), SCALAR>(
220 name_, solutions, u_zero_, u_one_, u_two_, k_);
221 }
222
223 private:
224 std::function<SCALAR(const Eigen::Matrix<double, 3, 1> &)> u_zero_;
225 std::function<Eigen::Matrix<SCALAR, 3, 1>(
226 const Eigen::Matrix<double, 3, 1> &)>
228 std::function<SCALAR(const Eigen::Matrix<double, 3, 1> &)> u_two_;
229 std::function<SCALAR(const Eigen::Matrix<double, 3, 1> &)> f_zero_;
230 std::function<Eigen::Matrix<SCALAR, 3, 1>(
231 const Eigen::Matrix<double, 3, 1> &)>
233 std::function<SCALAR(const Eigen::Matrix<double, 3, 1> &)> f_two_;
234 double &k_;
235 std::string name_;
236};
237
238} // namespace projects::hldo_sphere::experiments
239
240#endif // THESIS_EXPERIMENTS_HODGE_LAPLACIAN_H
Creates and solves the discretised Hodge Laplacian source problems for a given list of levels and val...
std::function< SCALAR(const Eigen::Matrix< double, 3, 1 > &)> f_zero_
HodgeLaplacianExperiment(std::function< SCALAR(const Eigen::Matrix< double, 3, 1 > &)> u_zero, std::function< Eigen::Matrix< SCALAR, 3, 1 >(const Eigen::Matrix< double, 3, 1 > &)> u_one, std::function< SCALAR(const Eigen::Matrix< double, 3, 1 > &)> u_two, std::function< SCALAR(const Eigen::Matrix< double, 3, 1 > &)> f_zero, std::function< Eigen::Matrix< SCALAR, 3, 1 >(const Eigen::Matrix< double, 3, 1 > &)> f_one, std::function< SCALAR(const Eigen::Matrix< double, 3, 1 > &)> f_two, double &k, std::string name)
Constructor setting all the functions and the reference k.
void Compute(std::vector< unsigned > refinement_levels, std::vector< double > ks)
Solves the hodge laplacian source problems for the tensor product of passed refinement levels and ks.
std::function< SCALAR(const Eigen::Matrix< double, 3, 1 > &)> u_two_
std::function< SCALAR(const Eigen::Matrix< double, 3, 1 > &)> u_zero_
std::function< Eigen::Matrix< SCALAR, 3, 1 >(const Eigen::Matrix< double, 3, 1 > &)> f_one_
std::function< SCALAR(const Eigen::Matrix< double, 3, 1 > &)> f_two_
std::function< Eigen::Matrix< SCALAR, 3, 1 >(const Eigen::Matrix< double, 3, 1 > &)> u_one_
void setRefinementLevel(lf::base::size_type n)
Set the refinement level.
std::shared_ptr< lf::mesh::Mesh > Build()
Build the mesh.
Creates and solves the discretised Hodge Laplacian source problems.
void SetLoadFunctions(std::function< SCALAR(const Eigen::Matrix< double, 3, 1 > &)> &f0, std::function< Eigen::Matrix< SCALAR, 3, 1 >(const Eigen::Matrix< double, 3, 1 > &)> &f1, std::function< SCALAR(const Eigen::Matrix< double, 3, 1 > &)> &f2)
Sets the load functions.
Wrapper classes for experiments.
void process_results(std::string name, ProblemSolutionWrapper< SCALAR > &results, U_ZERO &u_zero, U_ONE &u_one, U_TWO &u_two, double &k)
Generate vtk files containing meshdata and a csv table with results.
stores solutions for several k and one refinement level
stores solutions and setupinformation for a list of refinement levels and k values
std::vector< std::shared_ptr< const lf::mesh::Mesh > > mesh