Skip to content
Snippets Groups Projects
Commit 3670b93d authored by Dr. Felix Tobias Schindler's avatar Dr. Felix Tobias Schindler
Browse files

f func

parent 131ef52e
No related branches found
No related tags found
No related merge requests found
...@@ -67,10 +67,8 @@ public: ...@@ -67,10 +67,8 @@ public:
const ElementFilterType& filter = ApplyOnAllElements()) const ElementFilterType& filter = ApplyOnAllElements())
: BaseType(assembly_grid_view, source_spc, vec) : BaseType(assembly_grid_view, source_spc, vec)
{ {
this->append(LocalFunctionalType(local_binary_to_unary_element_integrand(inducing_function, LocalIntegrandType(1)), this->append(
over_integrate), LocalFunctionalType(LocalIntegrandType(1).with_ansatz(inducing_function), over_integrate), param, filter);
param,
filter);
} }
L2VolumeVectorFunctional(AssemblyGridView assembly_grid_view, L2VolumeVectorFunctional(AssemblyGridView assembly_grid_view,
...@@ -81,10 +79,8 @@ public: ...@@ -81,10 +79,8 @@ public:
const ElementFilterType& filter = ApplyOnAllElements()) const ElementFilterType& filter = ApplyOnAllElements())
: BaseType(assembly_grid_view, source_spc) : BaseType(assembly_grid_view, source_spc)
{ {
this->append(LocalFunctionalType(local_binary_to_unary_element_integrand(inducing_function, LocalIntegrandType(1)), this->append(
over_integrate), LocalFunctionalType(LocalIntegrandType(1).with_ansatz(inducing_function), over_integrate), param, filter);
param,
filter);
} }
}; // class L2VolumeVectorFunctional }; // class L2VolumeVectorFunctional
......
...@@ -160,10 +160,11 @@ public: ...@@ -160,10 +160,11 @@ public:
using typename FunctionalBaseType::SourceVectorType; using typename FunctionalBaseType::SourceVectorType;
using typename WalkerBaseType::ElementType; using typename WalkerBaseType::ElementType;
using typename WalkerBaseType::IntersectionType;
using ElementFilterType = XT::Grid::ElementFilter<AssemblyGridViewType>; using ElementFilterType = XT::Grid::ElementFilter<AssemblyGridViewType>;
using IntersectionFilterType = XT::Grid::IntersectionFilter<AssemblyGridViewType>; using IntersectionFilterType = XT::Grid::IntersectionFilter<AssemblyGridViewType>;
using ApplyOnAllElements = XT::Grid::ApplyOn::AllElements<AssemblyGridViewType>; using ApplyOnAllElements = XT::Grid::ApplyOn::AllElements<AssemblyGridViewType>;
using ApplyOnAllIntersection = XT::Grid::ApplyOn::AllIntersections<AssemblyGridViewType>; using ApplyOnAllIntersections = XT::Grid::ApplyOn::AllIntersections<AssemblyGridViewType>;
/** /**
* \name Ctors which accept an existing vector into which to assemble. * \name Ctors which accept an existing vector into which to assemble.
...@@ -232,7 +233,23 @@ public: ...@@ -232,7 +233,23 @@ public:
local_functional, param, XT::Grid::ApplyOn::GenericFilteredElements<AssemblyGridViewType>(filter_lambda)); local_functional, param, XT::Grid::ApplyOn::GenericFilteredElements<AssemblyGridViewType>(filter_lambda));
} }
// similar append for LocalIntersectionFunctionalInterface ... ThisType& append(const LocalIntersectionFunctionalInterface<IntersectionType, r, rC, F, DofFieldType>& local_functional,
const XT::Common::Parameter& param = {},
const IntersectionFilterType& filter = ApplyOnAllIntersections())
{
LocalIntersectionFunctionalAssembler<V, AssemblyGridViewType, r, rC, F, GV> tmp(
this->source_space(), local_functional, VectorStorage::access(), param);
this->append(tmp, filter);
return *this;
}
ThisType& append(const LocalIntersectionFunctionalInterface<IntersectionType, r, rC, F, DofFieldType>& local_functional,
const XT::Common::Parameter& param,
std::function<bool(const AssemblyGridViewType&, const IntersectionType&)> filter_lambda)
{
return append(
local_functional, param, XT::Grid::ApplyOn::GenericFilteredIntersections<AssemblyGridViewType>(filter_lambda));
}
void assemble(const bool use_tbb = false) override final void assemble(const bool use_tbb = false) override final
{ {
......
...@@ -96,6 +96,85 @@ private: ...@@ -96,6 +96,85 @@ private:
}; // class LocalElementFunctionalAssembler }; // class LocalElementFunctionalAssembler
template <class Vector, class GridView, size_t r, size_t rC, class R = double, class SpaceGridView = GridView>
class LocalIntersectionFunctionalAssembler : public XT::Grid::IntersectionFunctor<GridView>
{
static_assert(XT::LA::is_vector<Vector>::value, "");
static_assert(XT::Grid::is_view<GridView>::value, "");
static_assert(XT::Grid::is_view<SpaceGridView>::value, "");
using ThisType = LocalIntersectionFunctionalAssembler<Vector, GridView, r, rC, R, SpaceGridView>;
using BaseType = XT::Grid::IntersectionFunctor<GridView>;
public:
using typename BaseType::ElementType;
using typename BaseType::IntersectionType;
using VectorType = Vector;
using FieldType = typename VectorType::ScalarType;
using SpaceType = SpaceInterface<SpaceGridView, r, rC, R>;
using LocalFunctionalType = LocalIntersectionFunctionalInterface<IntersectionType, r, rC, R, FieldType>;
LocalIntersectionFunctionalAssembler(const SpaceType& space,
const LocalFunctionalType& local_functional,
VectorType& global_vector,
const XT::Common::Parameter& param = {})
: space_(space)
, local_functional_(local_functional.copy())
, global_vector_(global_vector)
, param_(param)
, local_vector_(space_.mapper().max_local_size())
, global_indices_(space_.mapper().max_local_size())
, basis_(space_.basis().localize())
{
DUNE_THROW_IF(global_vector_.size() != space_.mapper().size(),
XT::Common::Exceptions::shapes_do_not_match,
"global_vector.size() = " << global_vector_.size() << "\n "
<< "space.mapper().size()" << space_.mapper().size());
}
LocalIntersectionFunctionalAssembler(const ThisType& other)
: BaseType()
, space_(other.space_)
, local_functional_(other.local_functional_->copy())
, global_vector_(other.global_vector_)
, param_(other.param_)
, local_vector_(other.local_vector_)
, global_indices_(other.global_indices_)
, basis_(space_.basis().localize())
{}
LocalIntersectionFunctionalAssembler(ThisType&& source) = default;
BaseType* copy() override final
{
return new ThisType(*this);
}
void apply_local(const IntersectionType& intersection,
const ElementType& inside_element,
const ElementType& outside_element) override final
{
const auto& element = local_functional_->inside() ? inside_element : outside_element;
// apply functional
basis_->bind(element);
local_functional_->apply(intersection, *basis_, local_vector_, param_);
// copy local vector to global
space_.mapper().global_indices(element, global_indices_);
for (size_t jj = 0; jj < basis_->size(param_); ++jj)
global_vector_.add_to_entry(global_indices_[jj], local_vector_[jj]);
}
private:
const SpaceType& space_;
std::unique_ptr<LocalFunctionalType> local_functional_;
VectorType& global_vector_;
XT::Common::Parameter param_;
DynamicVector<FieldType> local_vector_;
DynamicVector<size_t> global_indices_;
mutable std::unique_ptr<typename SpaceType::GlobalBasisType::LocalizedType> basis_;
}; // class LocalIntersectionFunctionalAssembler
} // namespace GDT } // namespace GDT
} // namespace Dune } // namespace Dune
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment