Dune Core Modules (2.8.0)
recipe-integration.cc
See explanation at Integrate a function on a grid
// always include the config file
#ifdef HAVE_CONFIG_H
#include "config.h"
#endif
// C++ includes
#include<math.h>
#include<iostream>
// dune-common includes
#include<dune/common/parallel/mpihelper.hh>
#include<dune/common/parametertreeparser.hh>
#include<dune/common/timer.hh>
// dune-geometry includes
#include<dune/geometry/referenceelements.hh>
#include<dune/geometry/quadraturerules.hh>
// dune-grid includes
#include <dune/grid/yaspgrid.hh>
int main(int argc, char** argv)
{
// Maybe initialize Mpi
[[maybe_unused]] Dune::MPIHelper& helper = Dune::MPIHelper::instance(argc, argv);
// [set up grid]
const int dim = 4;
std::array<int,dim> cells; for (auto& c : cells) c=5;
Grid grid(len,cells);
// [small vectors and matrices]
Dune::FieldVector<double,4> x({1,2,3,4}); // make a vector
auto y(x); // copy constructor
y *= 1.0/3.0; // scaling
[[maybe_unused]] auto s = x*y; // scalar product
A.mv(x,y); // matvec: y = Ax
A.usmv(0.5,x,y); // axpy: y += 0.5*Ax
// [a function to integrate]
auto u = [](const auto& x){return std::exp(x.two_norm());};
// [integration with midpoint rule]
double integral=0.0;
auto gv = grid.leafGridView(); // extract the grid view
integral += u(e.geometry().center())*e.geometry().volume();
std::cout << "integral = " << integral << std::endl;
// [integration with quadrature rule]
double integral2 = 0.0;
{
auto geo = e.geometry();
auto quadrature = QR::rule(geo.type(),5);
for (const auto& qp : quadrature)
integral2 += u(geo.global(qp.position()))
*geo.integrationElement(qp.position())*qp.weight();
}
std::cout << "integral2 = " << integral2 << std::endl;
// [integrating a flux]
auto f = [](const auto& x){return x;};
double divergence=0.0;
if (!I.neighbor())
{
auto geoI = I.geometry();
divergence += f(geoI.center())*I.centerUnitOuterNormal()*geoI.volume();
}
}
std::cout << "divergence = " << divergence << std::endl;
}
FieldTraits< value_type >::real_type two_norm() const
two norm sqrt(sum over squared values of entries)
Definition: densevector.hh:642
A container for all quadrature rules of dimension dim
Definition: quadraturerules.hh:198
IteratorRange<... > intersections(const GV &gv, const Entity &e)
Iterates over all Intersections of an Entity with respect to the given GridView.
IteratorRange<... > elements(const GV &gv)
Iterates over all elements / cells (entities with codimension 0) of a GridView.
Helpers for dealing with MPI.
Various parser methods to get data into a ParameterTree object.
A simple timing class.
|
Legal Statements / Impressum |
Hosted by TU Dresden |
generated with Hugo v0.111.3
(Nov 12, 23:30, 2024)