10 #include "DECExamplesCommon.h" 11 #include "DGtal/math/linalg/EigenSupport.h" 12 #include "DGtal/dec/DiscreteExteriorCalculus.h" 13 #include "DGtal/dec/DiscreteExteriorCalculusSolver.h" 14 #include "DGtal/dec/DiscreteExteriorCalculusFactory.h" 15 #include "DGtal/io/boards/Board2D.h" 17 template <
typename Calculus>
18 typename Calculus::Scalar
19 sample_dual_0_form(
const typename Calculus::Properties& properties,
const typename Calculus::DualForm0& form,
const typename Calculus::Point& point)
21 const typename Calculus::Cell cell = form.myCalculus->myKSpace.uSpel(point);
22 const typename Calculus::Properties::const_iterator iter = properties.find(cell);
23 if (iter == properties.end())
return 0;
25 const typename Calculus::Index index = iter->second.index;
26 return form.myContainer(index);
29 template <
typename Calculus>
31 set_initial_phase_dir_yy(Calculus& calculus, Eigen::VectorXcd& initial_wave)
34 typedef typename Calculus::Index Index;
38 for (
int xx=9; xx<13; xx++)
39 for (
int yy=13; yy<17; yy++)
41 const Point point(xx,yy);
42 const Cell cell = calculus.myKSpace.uSpel(point);
43 const Index index = calculus.getCellIndex(cell);
44 initial_wave(index) = std::exp(std::complex<double>(0,2*M_PI*(yy-14.5)/8));
49 template <
typename Calculus>
51 set_initial_phase_null(Calculus& calculus, Eigen::VectorXcd& initial_wave)
54 typedef typename Calculus::Index Index;
58 for (
int xx=9; xx<13; xx++)
59 for (
int yy=13; yy<17; yy++)
61 const Point point(xx,yy);
62 const Cell cell = calculus.myKSpace.uSpel(point);
63 const Index index = calculus.getCellIndex(cell);
64 initial_wave(index) = 1;
69 using namespace DGtal;
82 const Calculus::Scalar cc = 8;
86 const Calculus calculus = CalculusFactory::createFromDigitalSet(generateDiskSet(
domain),
false);
89 const Calculus::DualIdentity0 laplace = calculus.laplace<
DUAL>() + 1e-8 * calculus.identity<0,
DUAL>();
91 trace.
info() <<
"laplace = " << laplace << endl;
95 typedef Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> EigenSolverMatrix;
96 const EigenSolverMatrix eigen_solver(laplace.myContainer);
98 const Eigen::VectorXd eigen_values = eigen_solver.eigenvalues();
99 const Eigen::MatrixXd eigen_vectors = eigen_solver.eigenvectors();
101 trace.
info() <<
"eigen_values_range = " << eigen_values.minCoeff() <<
"/" << eigen_values.maxCoeff() << endl;
105 const Eigen::VectorXd angular_frequencies = cc * eigen_values.array().sqrt();
108 Eigen::VectorXcd initial_wave = Eigen::VectorXcd::Zero(calculus.kFormLength(0,
DUAL));
110 set_initial_phase_dir_yy(calculus, initial_wave);
115 board <<
CustomStyle(
"KForm",
new KFormStyle2D(-1, 1));
116 board << Calculus::DualForm0(calculus, initial_wave.real());
117 board.
saveSVG(
"propagation_time_wave_initial_coarse.svg");
121 Eigen::VectorXcd initial_projections = eigen_vectors.transpose() * initial_wave;
126 const Calculus::Scalar lambda_cutoff = 4.5;
127 const Calculus::Scalar angular_frequency_cutoff = 2*M_PI * cc / lambda_cutoff;
129 for (
int kk=0; kk<initial_projections.rows(); kk++)
131 const Calculus::Scalar angular_frequency = angular_frequencies(kk);
132 if (angular_frequency < angular_frequency_cutoff)
continue;
133 initial_projections(kk) = 0;
137 trace.
info() <<
"cutted = " << cutted <<
"/" << initial_projections.rows() << endl;
140 const Eigen::VectorXcd wave = eigen_vectors * initial_projections;
143 board <<
CustomStyle(
"KForm",
new KFormStyle2D(-1, 1));
144 board << Calculus::DualForm0(calculus, wave.real());
145 board.
saveSVG(
"propagation_time_wave_initial_smoothed.svg");
148 const int kk_max = 60;
150 for (
int kk=0; kk<kk_max; kk++)
153 const Calculus::Scalar time = kk/20.;
154 const Eigen::VectorXcd current_projections = (angular_frequencies * std::complex<double>(0,time)).array().exp() * initial_projections.array();
155 const Eigen::VectorXcd current_wave = eigen_vectors * current_projections;
158 std::stringstream ss;
159 ss <<
"propagation_time_wave_solution_" << std::setfill(
'0') << std::setw(3) << kk <<
".svg";
163 board <<
CustomStyle(
"KForm",
new KFormStyle2D(-1, 1));
164 board << Calculus::DualForm0(calculus, current_wave.real());
165 board.
saveSVG(ss.str().c_str());
178 const Calculus calculus = CalculusFactory::createFromDigitalSet(generateDiskSet(domain),
false);
180 const Calculus::DualIdentity0 laplace = calculus.laplace<
DUAL>();
181 trace.
info() <<
"laplace = " << laplace << endl;
184 typedef Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> EigenSolverMatrix;
185 const EigenSolverMatrix eigen_solver(laplace.myContainer);
187 const Eigen::VectorXd laplace_eigen_values = eigen_solver.eigenvalues();
188 const Eigen::MatrixXd laplace_eigen_vectors = eigen_solver.eigenvectors();
189 trace.
info() <<
"eigen_values_range = " << laplace_eigen_values.minCoeff() <<
"/" << laplace_eigen_values.maxCoeff() << endl;
192 Calculus::DualForm0 concentration(calculus);
196 const Calculus::Index index = calculus.getCellIndex(cell);
197 concentration.myContainer(index) = 1;
203 board << concentration;
204 board.
saveSVG(
"propagation_forced_concentration.svg");
207 for (
int ll=0; ll<6; ll++)
210 const Calculus::Scalar lambda = 4*20.75/(1+2*ll);
212 trace.
info() <<
"lambda = " << lambda << endl;
215 const Eigen::VectorXd dalembert_eigen_values = (laplace_eigen_values.array() - (2*M_PI/lambda)*(2*M_PI/lambda)).array().inverse();
216 const Eigen::MatrixXd concentration_to_wave = laplace_eigen_vectors * dalembert_eigen_values.asDiagonal() * laplace_eigen_vectors.transpose();
220 Calculus::DualForm0 wave(calculus, concentration_to_wave * concentration.myContainer);
222 wave.myContainer /= wave.myContainer(calculus.getCellIndex(calculus.myKSpace.uSpel(
Z2i::Point(25,25))));
225 trace.
info() <<
"saving samples" << endl;
226 const Calculus::Properties properties = calculus.getProperties();
228 std::stringstream ss;
229 ss <<
"propagation_forced_samples_hv_" << ll <<
".dat";
230 std::ofstream handle(ss.str().c_str());
231 for (
int kk=0; kk<=50; kk++)
235 const Calculus::Scalar xx = 2 * (kk/50. - .5);
237 handle << sample_dual_0_form<Calculus>(properties, wave, point_horizontal) <<
" ";
238 handle << sample_dual_0_form<Calculus>(properties, wave, point_vertical) << endl;
243 std::stringstream ss;
244 ss <<
"propagation_forced_samples_diag_" << ll <<
".dat";
245 std::ofstream handle(ss.str().c_str());
246 for (
int kk=0; kk<=50; kk++)
250 const Calculus::Scalar xx = 2 * sqrt(2) * (kk/50. - .5);
252 handle << sample_dual_0_form<Calculus>(properties, wave, point_diag_pos) <<
" ";
253 handle << sample_dual_0_form<Calculus>(properties, wave, point_diag_neg) << endl;
259 std::stringstream ss;
260 ss <<
"propagation_forced_wave_" << ll <<
".svg";
263 board <<
CustomStyle(
"KForm",
new KFormStyle2D(-.5, 1));
265 board.
saveSVG(ss.str().c_str());
void beginBlock(const std::string &keyword="")
Aim: DiscreteExteriorCalculus represents a calculus in the dec package. This is the main structure in...
void progressBar(const double currentValue, const double maximalValue)
const Domain domain(Point(1, 2), Point(6, 5))
void saveSVG(const char *filename, PageSize size=Board::BoundingBox, double margin=10.0) const
Aim: This class provides static members to create DEC structures from various other DGtal structures...
int main(int argc, char **argv)
DGtal is the top-level namespace which contains all DGtal functions and types.
Aim: This class specializes a 'Board' class so as to display DGtal objects more naturally (with <<)...