8 #ifndef DUNE_GRIDGLUE_MERGING_CONTACTMERGE_HH 9 #define DUNE_GRIDGLUE_MERGING_CONTACTMERGE_HH 21 #include <dune/common/fvector.hh> 22 #include <dune/common/function.hh> 23 #include <dune/common/exceptions.hh> 24 #include <dune/common/bitsetvector.hh> 25 #include <dune/common/deprecated.hh> 27 #include <dune/grid/common/grid.hh> 40 template<
int dimworld,
typename T =
double>
44 static constexpr
int dim = dimworld-1;
46 static_assert( dim==1 || dim==2,
47 "ContactMerge yet only handles the cases dim==1 and dim==2!");
73 std::function<
WorldCoords(WorldCoords)> domainDirections=
nullptr,
74 std::function<
WorldCoords(WorldCoords)> targetDirections=
nullptr,
76 : domainDirections_(domainDirections), targetDirections_(targetDirections),
77 overlap_(allowedOverlap), type_(type)
86 : overlap_(allowedOverlap),
100 std::function<
WorldCoords(WorldCoords)> targetDirections)
102 domainDirections_ = domainDirections;
103 targetDirections_ = targetDirections;
125 maxNormalProduct_ = cos(angle);
134 return acos(maxNormalProduct_);
144 std::function<WorldCoords(WorldCoords)> domainDirections_;
145 std::vector<WorldCoords> nodalDomainDirections_;
155 std::function<WorldCoords(WorldCoords)> targetDirections_;
156 std::vector<WorldCoords> nodalTargetDirections_;
167 T maxNormalProduct_ = T(-0.1);
173 void computeIntersections(
const Dune::GeometryType& grid1ElementType,
174 const std::vector<Dune::FieldVector<T,dimworld> >& grid1ElementCorners,
175 std::bitset<(1<<dim)>& neighborIntersects1,
176 unsigned int grid1Index,
177 const Dune::GeometryType& grid2ElementType,
178 const std::vector<Dune::FieldVector<T,dimworld> >& grid2ElementCorners,
179 std::bitset<(1<<dim)>& neighborIntersects2,
180 unsigned int grid2Index,
181 std::vector<SimplicialIntersection>& intersections)
override;
187 void build(
const std::vector<Dune::FieldVector<T,dimworld> >& grid1Coords,
188 const std::vector<unsigned int>& grid1Elements,
189 const std::vector<Dune::GeometryType>& grid1ElementTypes,
190 const std::vector<Dune::FieldVector<T,dimworld> >& grid2Coords,
191 const std::vector<unsigned int>& grid2Elements,
192 const std::vector<Dune::GeometryType>& grid2ElementTypes)
override 194 std::cout<<
"ContactMerge building grid!\n";
197 grid2Coords, grid2Elements, grid2ElementTypes);
199 Base::build(grid1Coords, grid1Elements, grid1ElementTypes,
200 grid2Coords, grid2Elements, grid2ElementTypes);
207 static LocalCoords localCornerCoords(
int i,
const Dune::GeometryType& gt)
209 const auto& ref = Dune::ReferenceElements<T,dim>::general(gt);
210 return ref.position(i,dim);
216 void computeCyclicOrder(
const std::vector<std::array<LocalCoords,2> >& polytopeCorners,
217 const LocalCoords& center, std::vector<int>& ordering)
const;
221 const std::vector<unsigned int>& elements1,
222 const std::vector<Dune::GeometryType>& elementTypes1,
223 const std::vector<WorldCoords>& coords2,
224 const std::vector<unsigned int>& elements2,
225 const std::vector<Dune::GeometryType>& elementTypes2);
229 const std::vector<unsigned int>& elements,
230 const std::vector<Dune::GeometryType>& elementTypes,
231 std::vector<WorldCoords>& normals);
234 void removeDoubles(std::vector<std::array<LocalCoords,2> >& polytopeCorners);
242 #endif // DUNE_GRIDGLUE_MERGING_CONTACTMERGE_HH void build(const std::vector< Dune::FieldVector< T, dimworld > > &grid1Coords, const std::vector< unsigned int > &grid1Elements, const std::vector< Dune::GeometryType > &grid1ElementTypes, const std::vector< Dune::FieldVector< T, dimworld > > &grid2Coords, const std::vector< unsigned int > &grid2Elements, const std::vector< Dune::GeometryType > &grid2ElementTypes) override
builds the merged grid
Definition: contactmerge.hh:187
Dune::FieldVector< T, dim > LocalCoords
the coordinate type used in this interface
Definition: contactmerge.hh:61
void computeCyclicOrder(const std::vector< std::array< LocalCoords, 2 > > &polytopeCorners, const LocalCoords ¢er, std::vector< int > &ordering) const
Order the corners of the intersection polytope in cyclic order.
Definition: contactmerge.cc:212
T ctype
the numeric type used in this interface
Definition: contactmerge.hh:55
Definition: contactmerge.hh:64
void setOverlap(T overlap)
Set the allowed overlap of the surfaces.
Definition: contactmerge.hh:108
void computeOuterNormalField(const std::vector< WorldCoords > &coords, const std::vector< unsigned int > &elements, const std::vector< Dune::GeometryType > &elementTypes, std::vector< WorldCoords > &normals)
If no direction field was specified compute the outer normal field.
Definition: contactmerge.cc:294
bool valid
Definition: standardmerge.hh:84
ContactMerge(const T allowedOverlap, ProjectionType type)
Construct merger given overlap and type of the projection.
Definition: contactmerge.hh:85
void build(const std::vector< Dune::FieldVector< T, dimworld > > &grid1_Coords, const std::vector< unsigned int > &grid1_elements, const std::vector< Dune::GeometryType > &grid1_element_types, const std::vector< Dune::FieldVector< T, dimworld > > &grid2_coords, const std::vector< unsigned int > &grid2_elements, const std::vector< Dune::GeometryType > &grid2_element_types) override
builds the merged grid
Definition: standardmerge.hh:388
Common base class for many merger implementations: produce pairs of entities that may intersect...
Definition: standardmerge.hh:54
void removeDoubles(std::vector< std::array< LocalCoords, 2 > > &polytopeCorners)
Remove all multiples.
Definition: contactmerge.cc:333
ProjectionType
Type of the projection, closest point or outer normal projection.
Definition: contactmerge.hh:64
Definition: contactmerge.hh:64
StandardMerge< T, dimworld-1, dimworld-1, dimworld >::SimplicialIntersection SimplicialIntersection
Definition: contactmerge.hh:138
Common base class for many merger implementations: produce pairs of entities that may intersect...
Merge two codimension-1 surfaces that may be a positive distance apart.
Definition: contactmerge.hh:41
ContactMerge(const T allowedOverlap=T(0), std::function< WorldCoords(WorldCoords)> domainDirections=nullptr, std::function< WorldCoords(WorldCoords)> targetDirections=nullptr, ProjectionType type=OUTER_NORMAL)
Construct merger given overlap and possible projection directions.
Definition: contactmerge.hh:72
T getOverlap() const
Get the allowed overlap of the surfaces.
Definition: contactmerge.hh:114
void minNormalAngle(T angle)
set minimum angle in radians between normals at x and Φ(x)
Definition: contactmerge.hh:122
Dune::FieldVector< T, dimworld > WorldCoords
the coordinate type used in this interface
Definition: contactmerge.hh:58
void setSurfaceDirections(std::function< WorldCoords(WorldCoords)> domainDirections, std::function< WorldCoords(WorldCoords)> targetDirections)
Set surface direction functions.
Definition: contactmerge.hh:99
void setupNodalDirections(const std::vector< WorldCoords > &coords1, const std::vector< unsigned int > &elements1, const std::vector< Dune::GeometryType > &elementTypes1, const std::vector< WorldCoords > &coords2, const std::vector< unsigned int > &elements2, const std::vector< Dune::GeometryType > &elementTypes2)
Setup the direction vectors containing the directions for each vertex.
Definition: contactmerge.cc:267
Definition: gridglue.hh:35
Central component of the module implementing the coupling of two grids.
T minNormalAngle() const
get minimum angle in radians between normals at x and Φ(x)
Definition: contactmerge.hh:131