Newer
Older
#include <config.h>
// Includes for the ADOL-C automatic differentiation library
// Need to come before (almost) all others.
#include <adolc/adouble.h>

Oliver Sander
committed
#include <dune/fufem/utilities/adolcnamespaceinjections.hh>
#include <dune/common/bitsetvector.hh>
#include <dune/common/parametertree.hh>
#include <dune/common/parametertreeparser.hh>
#include <dune/grid/uggrid.hh>
#include <dune/grid/utility/structuredgridfactory.hh>
#include <dune/grid/io/file/gmshreader.hh>
#include <dune/functions/gridfunctions/discreteglobalbasisfunction.hh>

Oliver Sander
committed
#include <dune/functions/functionspacebases/pqknodalbasis.hh>
#include <dune/functions/functionspacebases/bsplinebasis.hh>

Oliver Sander
committed
#include <dune/fufem/boundarypatch.hh>
#include <dune/fufem/functiontools/basisinterpolator.hh>
#include <dune/fufem/functiontools/boundarydofs.hh>

Oliver Sander
committed
#include <dune/fufem/functionspacebases/dunefunctionsbasis.hh>

Oliver Sander
committed
#include <dune/fufem/discretizationerror.hh>
#include <dune/fufem/dunepython.hh>
#include <dune/solvers/solvers/iterativesolver.hh>
#include <dune/solvers/norms/energynorm.hh>
#include <dune/gfe/rotation.hh>
#include <dune/gfe/unitvector.hh>
#include <dune/gfe/realtuple.hh>
#include <dune/gfe/localgeodesicfeadolcstiffness.hh>
#include <dune/gfe/geodesicfeassembler.hh>
#include <dune/gfe/riemanniantrsolver.hh>
#include <dune/gfe/embeddedglobalgfefunction.hh>
// grid dimension
// Image space of the geodesic fe functions
// typedef Rotation<double,2> TargetSpace;
// typedef Rotation<double,3> TargetSpace;
// typedef UnitVector<double,2> TargetSpace;
typedef UnitVector<double,3> TargetSpace;
// typedef UnitVector<double,4> TargetSpace;
// typedef RealTuple<double,1> TargetSpace;
// Tangent vector of the image space
const int blocksize = TargetSpace::TangentVector::dimension;
#define LAGRANGE
using namespace Dune;
int main (int argc, char *argv[]) try
{
// Start Python interpreter
Python::start();
Python::Reference main = Python::import("__main__");
Python::run("import math");
//feenableexcept(FE_INVALID);
Python::runStream()
<< std::endl << "import sys"
<< std::endl << "sys.path.append('/home/sander/dune/dune-gfe/problems')"
<< std::endl;
typedef std::vector<TargetSpace> SolutionType;
// parse data file
ParameterTree parameterSet;
if (argc < 2)
DUNE_THROW(Exception, "Usage: ./harmonicmaps <parameter file>");
ParameterTreeParser::readINITree(argv[1], parameterSet);
ParameterTreeParser::readOptions(argc, argv, parameterSet);
const int numLevels = parameterSet.get<int>("numLevels");
#ifndef LAGRANGE
const double tolerance = parameterSet.get<double>("tolerance");
const int maxTrustRegionSteps = parameterSet.get<int>("maxTrustRegionSteps");
const double initialTrustRegionRadius = parameterSet.get<double>("initialTrustRegionRadius");
const int multigridIterations = parameterSet.get<int>("numIt");
const int nu1 = parameterSet.get<int>("nu1");
const int nu2 = parameterSet.get<int>("nu2");
const int mu = parameterSet.get<int>("mu");
const int baseIterations = parameterSet.get<int>("baseIt");
const double mgTolerance = parameterSet.get<double>("mgTolerance");
const double baseTolerance = parameterSet.get<double>("baseTolerance");
std::string resultPath = parameterSet.get("resultPath", "");
// ///////////////////////////////////////
// Create the grid
// ///////////////////////////////////////
typedef std::conditional<dim==1,OneDGrid,UGGrid<dim> >::type GridType;
shared_ptr<GridType> grid;
FieldVector<double,dim> lower(0), upper(1);
if (parameterSet.get<bool>("structuredGrid")) {
lower = parameterSet.get<FieldVector<double,dim> >("lower");
upper = parameterSet.get<FieldVector<double,dim> >("upper");
elements = parameterSet.get<array<unsigned int,dim> >("elements");
grid = StructuredGridFactory<GridType>::createCubeGrid(lower, upper, elements);
std::string path = parameterSet.get<std::string>("path");
std::string gridFile = parameterSet.get<std::string>("gridFile");
grid = shared_ptr<GridType>(GmshReader<GridType>::read(path + "/" + gridFile));
grid->globalRefine(numLevels-1);
//////////////////////////////////////////////////////////////////////////////////
// Construct the scalar function space basis corresponding to the GFE space
//////////////////////////////////////////////////////////////////////////////////
typedef Dune::Functions::PQkNodalBasis<typename GridType::LeafGridView, 1> FEBasis;
FEBasis feBasis(grid->leafGridView());
#else
typedef Dune::Functions::BSplineBasis<typename GridType::LeafGridView> FEBasis;
for (auto& i : elements)
i = i<<(numLevels-1);
FEBasis feBasis(grid->leafGridView(), lower, upper, elements, order);
#endif

Oliver Sander
committed
typedef DuneFunctionsBasis<FEBasis> FufemFEBasis;
FufemFEBasis fufemFeBasis(feBasis);
SolutionType x(fufemFeBasis.size());
// /////////////////////////////////////////
// Read Dirichlet values
// /////////////////////////////////////////
BitSetVector<1> allNodes(grid->size(dim));
BoundaryPatch<typename GridType::LeafGridView> dirichletBoundary(grid->leafGridView(), allNodes);
BitSetVector<blocksize> dirichletNodes;

Oliver Sander
committed
constructBoundaryDofs(dirichletBoundary,fufemFeBasis,dirichletNodes);
// //////////////////////////
// //////////////////////////
// Read initial iterate into a PythonFunction
typedef VirtualDifferentiableFunction<FieldVector<double, dim>, TargetSpace::CoordinateType> FBase;
Python::Module module = Python::import(parameterSet.get<std::string>("initialIterate"));
auto pythonInitialIterate = module.get("fdf").toC<std::shared_ptr<FBase>>();
std::vector<TargetSpace::CoordinateType> v;

Oliver Sander
committed
::Functions::interpolate(fufemFeBasis, v, *pythonInitialIterate);
#else
Dune::Functions::interpolate(feBasis, v, *pythonInitialIterate, lower, upper, elements, order);
#endif
for (size_t i=0; i<x.size(); i++)
x[i] = v[i];
// backup for error measurement later
SolutionType initialIterate = x;
// ////////////////////////////////////////////////////////////
// Create an assembler for the Harmonic Energy Functional
// ////////////////////////////////////////////////////////////
// Assembler using ADOL-C
typedef TargetSpace::rebind<adouble>::other ATargetSpace;

Oliver Sander
committed
std::shared_ptr<LocalGeodesicFEStiffness<FEBasis,ATargetSpace> > localEnergy;

Oliver Sander
committed
std::string energy = parameterSet.get<std::string>("energy");
if (energy == "harmonic")
{

Oliver Sander
committed
localEnergy.reset(new HarmonicEnergyLocalStiffness<FEBasis, ATargetSpace>);

Oliver Sander
committed
} else if (energy == "chiral_skyrmion")
{

Oliver Sander
committed
localEnergy.reset(new GFE::ChiralSkyrmionEnergy<FEBasis, adouble>(parameterSet.sub("energyParameters")));

Oliver Sander
committed
} else
DUNE_THROW(Exception, "Unknown energy type '" << energy << "'");

Oliver Sander
committed
LocalGeodesicFEADOLCStiffness<FEBasis,TargetSpace> localGFEADOLCStiffness(localEnergy.get());
GeodesicFEAssembler<FEBasis,TargetSpace> assembler(feBasis, &localGFEADOLCStiffness);
// /////////////////////////////////////////////////
// Create a Riemannian trust-region solver
// /////////////////////////////////////////////////
RiemannianTrustRegionSolver<FEBasis,TargetSpace> solver;
solver.setup(*grid,
&assembler,
x,
dirichletNodes,
tolerance,
maxTrustRegionSteps,
initialTrustRegionRadius,
multigridIterations,
mgTolerance,
mu, nu1, nu2,
baseIterations,
baseTolerance,
// /////////////////////////////////////////////////////
// Solve!
// /////////////////////////////////////////////////////
std::cout << "Energy: " << assembler.computeEnergy(x) << std::endl;
solver.solve();
x = solver.getSol();
// //////////////////////////////
// Output result
// //////////////////////////////
typedef BlockVector<TargetSpace::CoordinateType> EmbeddedVectorType;
EmbeddedVectorType xEmbedded(x.size());
for (size_t i=0; i<x.size(); i++)
xEmbedded[i] = x[i].globalCoordinates();
auto xFunction = Dune::Functions::makeDiscreteGlobalBasisFunction<TargetSpace::CoordinateType>(feBasis,TypeTree::hybridTreePath(),xEmbedded);

Oliver Sander
committed
auto localXFunction = localFunction(xFunction);

Oliver Sander
committed
VTKWriter<GridType::LeafGridView> vtkWriter(grid->leafGridView());
vtkWriter.addVertexData(localXFunction, VTK::FieldInfo("orientation", VTK::FieldInfo::Type::scalar, xEmbedded[0].size()));
vtkWriter.write(resultPath + "_" + energy + "_result");
// Write the corresponding coefficient vector: verbatim in binary, to be completely lossless
std::ofstream outFile("harmonicmaps-result-" + std::to_string(numLevels) + ".data", std::ios_base::binary);
GenericVector::writeBinary(outFile, xEmbedded);
outFile.close();
return 0;
} catch (Exception e) {
std::cout << e << std::endl;
}