Example of solving Poisson equation using the DEC package.
#include "DGtal/io/readers/VolReader.h"
#include "DGtal/io/viewers/Viewer3D.h"
#include "DGtal/topology/SurfelAdjacency.h"
#include "DGtal/topology/LightImplicitDigitalSurface.h"
#include "DGtal/topology/DigitalSurface.h"
#include "DGtal/math/linalg/EigenSupport.h"
#include "DGtal/dec/DiscreteExteriorCalculus.h"
#include "DGtal/dec/DiscreteExteriorCalculusSolver.h"
#include "DGtal/dec/DiscreteExteriorCalculusFactory.h"
#include "ConfigExamples.h"
template <typename Predicate, typename Domain>
struct FalseOutsideDomain
{
myPredicate(&predicate), myDomain(&adomain)
{
}
bool
operator()(
const Point& point)
const
{
if (!myDomain->isInside(point)) return false;
return (*myPredicate)(point);
}
const Predicate* myPredicate;
};
void
alcapone_3d()
{
using std::endl;
const std::string filename = examplesPath + "samples/Al.100.vol";
typedef FalseOutsideDomain<Image, DGtal::Z3i::Domain> ImageExtended;
const ImageExtended image_extended(image,
domain);
const SurfelAdjacency surfel_adjacency(true);
const DigitalSurfaceContainer digital_surface_container(kspace, image_extended, surfel_adjacency, cell_bel);
const DigitalSurface digital_surface(digital_surface_container);
trace.
info() <<
"digital_surface_size=" << digital_surface.size() << endl;
{
viewer->setWindowTitle("alcapone surface");
{
(*viewer) << cell;
}
}
const Calculus
calculus = CalculusFactory::createFromNSCells<2>(digital_surface.begin(), digital_surface.end(),
false);
{
viewer->setWindowTitle("alcapone calculus");
}
for (int index=0; index<rho.length(); index++)
{
if (point[2]>1) continue;
rho.myContainer(index) = point[1]>100 ? 1 : -1;
}
trace.
info() <<
"rho_range=" << rho.myContainer.minCoeff() <<
"/" << rho.myContainer.maxCoeff() << endl;
{
viewer->setWindowTitle("alcapone poisson rho");
}
const Calculus::DualIdentity0 laplace =
calculus.laplace<
DUAL>();
PoissonSolver solver;
const Calculus::DualForm0
phi = solver.solve(rho);
trace.
info() <<
"phi_range=" <<
phi.myContainer.minCoeff() <<
"/" <<
phi.myContainer.maxCoeff() << endl;
{
viewer->setWindowTitle("alcapone poisson phi");
}
}
void
pyramid_3d()
{
using std::endl;
trace.
info() <<
"input_domain=" << input_domain << endl;
kspace.
init(input_domain.lowerBound(), input_domain.upperBound(),
true);
trace.
info() <<
"input_set_size=" << input_set.size() << endl;
{
viewer->setWindowTitle("input set");
(*viewer) << input_set;
}
const SurfelAdjacency surfel_adjacency(true);
const DigitalSurfaceContainer digital_surface_container(kspace, input_set, surfel_adjacency, cell_bel);
const DigitalSurface digital_surface(digital_surface_container);
for (int ii=0; ii<3; ii++)
{
std::stringstream ss;
for (int jj=0; jj<3; jj++)
ss << surfel_adjacency.getAdjacency(ii, jj) << " ";
}
trace.
info() <<
"digital_surface_container=" << digital_surface_container << endl;
trace.
info() <<
"digital_surface_size=" << digital_surface.size() << endl;
{
viewer->setWindowTitle("digital surface");
{
(*viewer) << cell;
}
}
const Calculus
calculus = CalculusFactory::createFromNSCells<2>(digital_surface.begin(), digital_surface.end());
{
viewer->setWindowTitle("discrete exterior calculus");
}
const Calculus::PrimalForm2 primal_surfel_area =
calculus.hodge<0,
DUAL>() *
Calculus::DualForm0(
calculus, Eigen::VectorXd::Ones(
calculus.kFormLength(0, DUAL)));
const Calculus::DualForm2 dual_surfel_area =
calculus.hodge<0,
PRIMAL>() *
Calculus::PrimalForm0(
calculus, Eigen::VectorXd::Ones(
calculus.kFormLength(0, PRIMAL)));
const double area_th =
calculus.kFormLength(2, PRIMAL);
const double area_primal = dual_surfel_area.myContainer.sum();
const double area_dual = primal_surfel_area.myContainer.sum();
trace.
info() <<
"area_th=" << area_th << endl;
trace.
info() <<
"area_primal=" << area_primal << endl;
trace.
info() <<
"area_dual=" << area_dual << endl;
FATAL_ERROR( area_th == area_primal );
FATAL_ERROR( area_th == area_dual );
const Calculus::DualForm1 dual_edge_length =
calculus.hodge<1,
PRIMAL>() *
Calculus::PrimalForm1(
calculus, Eigen::VectorXd::Ones(
calculus.kFormLength(1, PRIMAL)));
const bool dual_edge_length_match = dual_edge_length.myContainer == Eigen::VectorXd::Ones(
calculus.kFormLength(1, DUAL)) ;
trace.
info() <<
"dual_edge_length_match=" << dual_edge_length_match << endl;
FATAL_ERROR( dual_edge_length_match );
{
const Calculus::Scalar area = dual_surfel_area.myContainer(index);
trace.
info() << index <<
" " << cell <<
" " << area << endl;
}
}
int main(
int argc,
char* argv[])
{
QApplication app(argc, argv);
alcapone_3d();
return app.exec();
}
Aim: This class encapsulates its parameter class so that to indicate to the user that the object/poin...
Aim: A wrapper class around a STL associative container for storing sets of digital points within som...
Aim: Represents a set of n-1-cells in a nD space, together with adjacency relation between these cell...
Aim: This class provides static members to create DEC structures from various other DGtal structures.
Aim: This wraps a linear algebra solver around a discrete exterior calculus.
DiscreteExteriorCalculusSolver & compute(const Operator &linear_operator)
Aim: DiscreteExteriorCalculus represents a calculus in the dec package. This is the main structure in...
Aim: This class is a model of CCellularGridSpaceND. It represents the cubical grid as a cell complex,...
bool init(const Point &lower, const Point &upper, bool isClosed)
Specifies the upper and lower bounds for the maximal cells in this space.
const Point & sKCoords(const SCell &c) const
Return its Khalimsky coordinates.
Aim: A model of CDigitalSurfaceContainer which defines the digital surface as the boundary of an impl...
static Self diagonal(Component val=1)
static SCell findABel(const KSpace &K, const PointPredicate &pp, unsigned int nbtries=1000)
Aim: Represent adjacencies between surfel elements, telling if it follows an interior to exterior ord...
void beginBlock(const std::string &keyword="")
virtual void show()
Overload QWidget method in order to add a call to updateList() method (to ensure that the lists are w...
PolygonalCalculus< SH3::RealPoint, SH3::RealVector >::Vector phi(const Face f)
MyDigitalSurface::ConstIterator ConstIterator
Factory for GPL Display3D:
static void draw(Display3D< Space, KSpace > &display, const DGtal::DiscreteExteriorCalculus< dimEmbedded, dimAmbient, TLinearAlgebraBackend, TInteger > &calculus)
Eigen::SparseLU< SparseMatrix > SolverSparseLU
Represents a signed cell in a cellular grid space by its Khalimsky coordinates and a boolean value.
static ImageContainer importVol(const std::string &filename, const Functor &aFunctor=Functor())
Aim: This concept represents a digital domain, i.e. a non mutable subset of points of the given digit...
Aim: Defines a predicate on a point.
Aim: Define a simple functor using the static cast operator.
int main(int argc, char **argv)
InHalfPlaneBySimple3x3Matrix< Point, double > Functor
ImageContainerBySTLVector< Domain, Value > Image
HyperRectDomain< Space > Domain