<html>
<head>
<meta content="text/html; charset=windows-1252"
http-equiv="Content-Type">
</head>
<body bgcolor="#FFFFFF" text="#000000">
<p>Hi Christian,<br>
<br>
I did not compile the module, but from the lines below I extract
that you are using a Dune::OneDGrid (i.e. dimension 1) with a 2d
first order QkDG local basis (check the order of the template
parameters), which is why the parameters of "evaluateFunction" are
not convertible.<br>
<br>
Best,<br>
Andreas<br>
<br>
/home/dune/Dune/dune-pdelab/dune/pdelab/newton/newton.hh:239:11:
required from 'void Dune::PDELab::NewtonSolver<GOS, S, TrlV,
TstV>::apply() [with GOS =
Dune::PDELab::GridOperator<Dune::PDELab::GridFunctionSpace<Dune::GridView<Dune::DefaultLeafGridViewTraits<const
Dune::OneDGrid, (Dune::PartitionIteratorType)4u> >,
Dune::PDELab::QkDGLocalFiniteElementMap<double, double, 1,
2>, Dune::PDELab::ConformingDirichletConstraints,
Dune::PDELab::istl::VectorBackend<> >,
Dune::PDELab::GridFunctionSpace<Dune::GridView<Dune::DefaultLeafGridViewTraits<const
Dune::OneDGrid, (Dune::PartitionIteratorType)4u> >,
Dune::PDELab::QkDGLocalFiniteElementMap<double, double, 1,
2>, Dune::PDELab::ConformingDirichletConstraints,
Dune::PDELab::istl::VectorBackend<> >,
Dune::PDELab::ConvectionDiffusionDG<Dune::SPL::ProblemA<Dune::GridView<Dune::DefaultLeafGridViewTraits<const
Dune::OneDGrid, (Dune::PartitionIteratorType)4u> >,
double>, Dune::PDELab::QkDGLocalFiniteElementMap<double,
double, 1, 2> >, Dune::PDELab::ISTLMatrixBackend, double,
double, double,
Dune::PDELab::ConstraintsTransformation<Dune::PDELab::DOFIndex<long
unsigned int, 1ul, 2ul>, Dune::PDELab::MultiIndex<long
unsigned int, 1ul>, double>,
Dune::PDELab::ConstraintsTransformation<Dune::PDELab::DOFIndex<long
unsigned int, 1ul, 2ul>, Dune::PDELab::MultiIndex<long
unsigned int, 1ul>, double> >; S =
Dune::PDELab::ISTLBackend_SEQ_BCGS_Jac; TrlV =
Dune::PDELab::istl::BlockVector<Dune::PDELab::GridFunctionSpace<Dune::GridView<Dune::DefaultLeafGridViewTraits<const
Dune::OneDGrid, (Dune::PartitionIteratorType)4u> >,
Dune::PDELab::QkDGLocalFiniteElementMap<double, double, 1,
2>, Dune::PDELab::ConformingDirichletConstraints,
Dune::PDELab::istl::VectorBackend<> >,
Dune::BlockVector<Dune::FieldVector<double, 1>,
std::allocator<Dune::FieldVector<double, 1> > >
>; TstV =
Dune::PDELab::istl::BlockVector<Dune::PDELab::GridFunctionSpace<Dune::GridView<Dune::DefaultLeafGridViewTraits<const
Dune::OneDGrid, (Dune::PartitionIteratorType)4u> >,
Dune::PDELab::QkDGLocalFiniteElementMap<double, double, 1,
2>, Dune::PDELab::ConformingDirichletConstraints,
Dune::PDELab::istl::VectorBackend<> >,
Dune::BlockVector<Dune::FieldVector<double, 1>,
std::allocator<Dune::FieldVector<double, 1> > >
>]'<br>
/home/dune/Dune/dune-spl/src/diffusion.cc:193:15: required from
here<br>
/home/dune/Dune/dune-pdelab/dune/pdelab/localoperator/convectiondiffusiondg.hh:939:102:
error: no matching function for call to
'Dune::PDELab::LocalBasisCache<Dune::QkStuff::QkLocalBasis<double,
double, 1, 2> >::evaluateFunction(const Vector&, const
LocalBasisType&) const'<br>
auto& phi =
cache[order].evaluateFunction(ip.position(),lfsv.finiteElement().localBasis());<br>
^<br>
/home/dune/Dune/dune-pdelab/dune/pdelab/localoperator/convectiondiffusiondg.hh:939:102:
note: candidate is:<br>
In file included from
/home/dune/Dune/dune-pdelab/dune/pdelab/localoperator/convectiondiffusiondg.hh:17:0,<br>
from
/home/dune/Dune/dune-spl/src/diffusion.cc:53:<br>
/home/dune/Dune/dune-pdelab/dune/pdelab/finiteelement/localbasiscache.hh:47:7:
note: const std::vector<typename LB::Traits::RangeType>&
Dune::PDELab::LocalBasisCache<LocalBasisType>::evaluateFunction(const
DomainType&, const LocalBasisType&) const [with
LocalBasisType = Dune::QkStuff::QkLocalBasis<double, double, 1,
2>; typename LB::Traits::RangeType =
Dune::FieldVector<double, 1>;
Dune::PDELab::LocalBasisCache<LocalBasisType>::DomainType =
Dune::FieldVector<double, 2>]<br>
evaluateFunction (const DomainType& position, const
LocalBasisType& localbasis) const<br>
^<br>
/home/dune/Dune/dune-pdelab/dune/pdelab/finiteelement/localbasiscache.hh:47:7:
note: no known conversion for argument 1 from 'const Vector {aka
const Dune::FieldVector<double, 1>}' to 'const
DomainType& {aka const Dune::FieldVector<double,
2>&}'<br>
</p>
<br>
<div class="moz-cite-prefix">On 08.06.2016 09:38, Christian
Kaltenecker wrote:<br>
</div>
<blockquote cite="mid:a2ce6565-9759-51f1-ec01-f954e6759f64@yahoo.de"
type="cite">
<pre wrap="">Hello Dune,
I am currently working with Alexander Grebhahn from the university of
Passau on creating an automatic analysis program to find several
variation points within a Dune application.
For the evaluation part, I want to generate several variants of the
diffusion-convection simulation by also using dune-testtools.
Above others, I am using ALUGrid (both dune-alugrid and the external
ALUGrid(v1.52) library), UGGrid(v3.11) and SuperLU(v4.3).
I have successfully compiled and linked those libraries.
Furthermore, I have checked out the "releases/2.4" branches of every
Dune module in order to have all Dune modules in the same version
(every module except my module compiles without errors).
Moreover, the whole Dune framework is running in a docker container.
I am executing dunecontrol with the following commands:
sudo dune-common/bin/dunecontrol configure "--with-ug=/home/dune/Libs/ug
--with-alugrid=/home/dune/Libs/ALUGrid
--with-superlu=/home/dune/Libs/SuperLU_4.3
--with-superlu-lib=libsuperlu_4.3.a"
sudo dune-common/bin/dunecontrol make all
The errors (as well as the output of configure) are included in the
attachment.
If you want to try it on your own, I have also included my Dune module
(dune-spl).
Note that I have commented out the variability that is not working in
the meta ini file, diffusion.mini. Also note, that you need dune-testtools.
What could I do to fix the problems that occurr in the log-files? Do I have to add some additional configuration parameters? I am writing to this mailing list because the errors are thrown in the dune-pdelab module.
Thanks in advance!
Regards,
Christian Kaltenecker
</pre>
<br>
<fieldset class="mimeAttachmentHeader"></fieldset>
<br>
<pre wrap="">_______________________________________________
dune-pdelab mailing list
<a class="moz-txt-link-abbreviated" href="mailto:dune-pdelab@dune-project.org">dune-pdelab@dune-project.org</a>
<a class="moz-txt-link-freetext" href="http://lists.dune-project.org/mailman/listinfo/dune-pdelab">http://lists.dune-project.org/mailman/listinfo/dune-pdelab</a>
</pre>
</blockquote>
<br>
</body>
</html>