1 #ifndef DUNE_FEM_HDIV_PROJECTION_HH 2 #define DUNE_FEM_HDIV_PROJECTION_HH 5 #include <dune/geometry/referenceelements.hh> 25 template<
int dim,
class CoordCont >
31 #endif // #ifdef ENABLE_UG 57 template <
class DiscreteFunctionType>
60 typedef typename DiscreteFunctionType::DiscreteFunctionSpaceType DiscreteFunctionSpaceType;
61 typedef typename DiscreteFunctionType :: LocalFunctionType LocalFunctionType;
62 typedef typename DiscreteFunctionSpaceType::RangeType RangeType;
63 typedef typename DiscreteFunctionSpaceType::DomainType DomainType;
64 typedef typename DiscreteFunctionSpaceType::DomainFieldType DomainFieldType;
65 typedef typename DiscreteFunctionSpaceType::RangeFieldType RangeFieldType;
66 typedef typename DiscreteFunctionSpaceType::JacobianRangeType JacobianRangeType;
67 typedef typename DiscreteFunctionSpaceType::GridPartType GridPartType;
76 typedef typename GridPartType :: GridType GridType;
78 enum { dimFaceRange = 1 };
79 enum { dimDomain = DiscreteFunctionSpaceType::dimDomain };
80 enum { dimFaceDomain = dimDomain - 1 };
81 enum { polOrdN = DiscreteFunctionSpaceType :: polynomialOrder };
86 enum { gradPolOrd = ((polOrdN - 1) < 0) ? 0 : (polOrdN - 1) };
88 template <
class Space>
struct BubbleM;
90 template <
class FunctionSpaceImp,
93 template <
class>
class StorageImp,
94 template <
class,
class,
int,
template <
class>
class>
class DiscreteFunctionSpaceImp>
95 struct BubbleM <DiscreteFunctionSpaceImp<FunctionSpaceImp,GridPartImp,polOrd,StorageImp> >
97 enum { bubblePModifier = dimDomain - 1 };
100 template <
class FunctionSpaceImp,
103 template <
class>
class StorageImp>
104 struct BubbleM < LegendreDiscontinuousGalerkinSpace<FunctionSpaceImp,GridPartImp,polOrd,StorageImp> >
106 enum { bubblePModifier = 1 };
109 template <
class DiscreteFunctionSpaceImp,
112 struct BubbleM <CombinedSpace<DiscreteFunctionSpaceImp,N,policy> >
114 enum { bubblePModifier = BubbleM< DiscreteFunctionSpaceImp > :: bubblePModifier };
117 typedef BubbleM <DiscreteFunctionSpaceType> BubbleMType;
119 #ifdef USE_TWISTFREE_MAPPER 121 enum { bubblePModifier = BubbleMType :: bubblePModifier };
124 enum { bubblePolOrd = (polOrdN + bubblePModifier) };
127 #warning "Hdiv-Projection only working for polOrd = 1 (enable higher order Lagrange with -DUSE_TWISTFREE_MAPPER)" 130 enum { bubblePModifier = 1 };
133 enum { bubblePolOrd = (polOrdN > 1) ? 2 : (polOrdN + 1) };
136 template <
class Space>
struct Spaces;
138 template <
class FunctionSpaceImp,
141 template <
class>
class StorageImp,
142 template <
class,
class,
int,
template <
class>
class>
class DiscreteFunctionSpaceImp>
143 struct Spaces<DiscreteFunctionSpaceImp<FunctionSpaceImp,GridPartImp,polOrd,StorageImp> >
147 typedef DiscreteFunctionSpaceImp<ElSpaceType,GridPartImp,gradPolOrd,StorageImp> ElementGradientSpaceType;
149 typedef DiscreteFunctionSpaceImp<FaceSpaceType,GridPartType,polOrdN,StorageImp> FaceDiscreteSpaceType;
155 template <
class DiscreteFunctionSpaceImp,
158 struct Spaces<CombinedSpace<DiscreteFunctionSpaceImp,N,policy> >
160 typedef typename DiscreteFunctionSpaceImp :: GridPartType GridPartImp;
161 typedef Spaces<DiscreteFunctionSpaceImp> AllSpacesType;
162 typedef typename AllSpacesType :: ElementGradientSpaceType ElementGradientSpaceType;
163 typedef typename AllSpacesType :: ElementDiscreteSpaceType ElementDiscreteSpaceType;
164 typedef typename AllSpacesType :: FaceDiscreteSpaceType FaceDiscreteSpaceType;
167 typedef Spaces<DiscreteFunctionSpaceType> AllSpacesType;
168 typedef typename AllSpacesType :: ElementGradientSpaceType ElementGradientSpaceType;
169 typedef typename AllSpacesType :: ElementDiscreteSpaceType ElementDiscreteSpaceType;
170 typedef typename AllSpacesType :: FaceDiscreteSpaceType FaceDiscreteSpaceType;
172 const DiscreteFunctionSpaceType& space_;
173 GridPartType & gridPart_;
174 const FaceDiscreteSpaceType faceSpace_;
175 const ElementDiscreteSpaceType elSpace_;
176 const ElementGradientSpaceType gradSpace_;
182 gridPart_(space.gridPart()),
183 faceSpace_( gridPart_ ),
184 elSpace_( gridPart_ ),
185 gradSpace_( gridPart_ )
190 gridPart_( org.gridPart_),
191 faceSpace_( gridPart_ ),
192 elSpace_( gridPart_ ),
193 gradSpace_( gridPart_ )
197 const DiscreteFunctionSpaceType&
space()
const 208 static double normalJump(
const DiscreteFunctionType &discFunc,
const int polyOrder = -1 )
210 typedef typename GridPartType :: IntersectionIteratorType IntersectionIteratorType;
211 typedef typename IntersectionIteratorType :: Intersection IntersectionType;
212 typedef typename GridType :: template Codim<0> :: Entity EntityType;
213 typedef typename GridType :: Traits :: LocalIdSet LocalIdSetType;
215 enum { dim = GridType::dimension };
217 const DiscreteFunctionSpaceType &
space = discFunc.space();
218 const GridPartType &gridPart = space.gridPart();
219 const int polOrd = (polyOrder <0) ? (2 * space.order() + 2) : polyOrder;
221 typedef typename DiscreteFunctionType::LocalFunctionType LocalFuncType;
224 RangeType neighRet (0.0);
231 const LocalIdSetType &idSet = gridPart.grid().localIdSet();
233 for(
const auto & en : space)
235 const LocalFuncType lf = discFunc.localFunction(en);
237 double localValue = 0.0;
239 IntersectionIteratorType endnit = gridPart.iend(en);
240 for(IntersectionIteratorType nit = gridPart.ibegin(en);
241 nit != endnit; ++nit)
243 const IntersectionType& inter=*nit;
245 if(inter.neighbor() )
249 if(idSet.id( en ) < idSet.id( nb ))
251 FaceQuadratureType faceQuadInner(gridPart, inter, polOrd, FaceQuadratureType::INSIDE);
252 FaceQuadratureType faceQuadOuter(gridPart, inter, polOrd, FaceQuadratureType::OUTSIDE);
254 const LocalFuncType neighLf = discFunc.localFunction(nb);
256 const int quadNop = faceQuadInner.nop();
257 for (
int l = 0; l < quadNop ; ++l)
260 inter.unitOuterNormal(faceQuadInner.localPoint(l));
262 lf.evaluate(faceQuadInner[l], ret);
263 neighLf.evaluate(faceQuadOuter[l], neighRet);
267 double val = ret * normal;
268 val *= faceQuadInner.weight(l);
284 void curl(
const DomainType & arg, DomainType & dest,
const int d )
const 286 if( DomainType :: dimension == 2 )
293 else if( DomainType :: dimension == 3 )
324 template <
class EntityType,
325 class LocalFunctionType,
329 void bubblePart(
const ElementDiscreteSpaceType&
space,
331 const LocalFunctionType & uLF,
const int startRow ,
333 MatrixType & matrix, VectorType& rhs )
const 336 typedef typename ElementDiscreteSpaceType :: BaseFunctionSetType BaseFunctionSetType;
337 typedef typename ElementDiscreteSpaceType :: LagrangePointSetType LagrangePointSetType;
338 typedef typename ElementDiscreteSpaceType :: Traits::GridPartType GridPartType;
340 typedef typename ElementDiscreteSpaceType :: JacobianRangeType JacobianRangeType;
341 typedef typename ElementDiscreteSpaceType :: DomainType DomainType;
343 enum { dim = GridType::dimension };
345 const LagrangePointSetType& lagrangePointSet = space.lagrangePointSet( en );
346 const BaseFunctionSetType bSet = space.baseFunctionSet( en );
347 const int polOrd = 2 * space.order();
349 const int cols = uLF.numDofs();
350 assert( uRets.size() == (
unsigned int)cols );
352 VolumeQuadratureType quad (en,polOrd);
354 std::vector< JacobianRangeType > valTmpVec( bSet.size() );
359 typedef typename EntityType :: Geometry Geometry ;
360 const Geometry& geo = en.geometry();
362 const GeometryType& type = geo.type();
363 const int bubbleOffset = (type.isSimplex()) ? 0 : baseFunctionOffset( 0 );
366 typedef typename Geometry :: ctype ctype;
367 enum { cdim = Geometry :: coorddimension };
368 enum { mydim = Geometry :: mydimension };
369 typedef typename Geometry::JacobianInverseTransposed JacobianInverseType;
372 const int enDofs = numberOfBubbles( lagrangePointSet.numDofs( 0 ), type ,
374 const int bubbleMod = bubbleModifier( mydim );
376 const int quadNop = quad.nop();
377 for (
int l = 0; l < quadNop ; ++l)
380 const JacobianInverseType& inv = geo.jacobianInverseTransposed( quad.point(l) );
383 const double intel = quad.weight(l) *
384 geo.integrationElement(quad.point(l));
387 uLF.evaluate(quad[l], result);
390 uLF.baseFunctionSet().evaluateAll( quad[l], uRets );
393 bSet.jacobianAll( quad[l], inv, valTmpVec );
396 for(
int i = 0 ; i<enDofs; i += bubbleMod )
399 int row = startRow + i;
402 const int localBaseFct = ((int) i/bubbleMod) + bubbleOffset;
404 const int baseFct = lagrangePointSet.entityDofNumber( 0, 0, localBaseFct );
407 JacobianRangeType& val = valTmpVec[ baseFct ];
412 for(
int d = 0; d<bubbleMod; ++d )
415 curl(val[0], aVal, d );
417 double r = aVal * result;
422 for(
int j=0; j<cols; ++j)
424 double t = aVal * uRets[ j ];
426 matrix[ row ][ j ] += t;
436 template <
class EntityType,
class LocalFunctionType,
class ArrayType,
437 class MatrixType,
class VectorType>
438 void gradientPart(
const ElementGradientSpaceType & space,
440 const LocalFunctionType & uLF,
const int startRow ,
441 ArrayType& uRets, MatrixType& matrix, VectorType& rhs )
const 443 typedef typename ElementGradientSpaceType :: BaseFunctionSetType BaseFunctionSetType;
444 const BaseFunctionSetType bSet = space.baseFunctionSet( en );
445 int polOrd = 2 * space.order() + 1;
447 const int localRows = gradientBaseFct( bSet );
448 const int cols = uLF.numDofs();
450 VolumeQuadratureType quad (en,polOrd);
455 typedef typename ElementGradientSpaceType :: JacobianRangeType GradJacobianRangeType;
456 std::vector< GradJacobianRangeType > gradTmpVec( bSet.size() );
457 GradJacobianRangeType gradPhi;
459 typedef typename EntityType :: Geometry Geometry ;
460 const Geometry& geo = en.geometry();
462 typedef typename Geometry :: ctype ctype;
463 enum { cdim = Geometry :: coorddimension };
464 enum { mydim = Geometry :: mydimension };
465 typedef typename Geometry::JacobianInverseTransposed JacobianInverseType;
467 const int quadNop = quad.nop();
468 for (
int l = 0; l < quadNop ; ++l)
471 const JacobianInverseType& inv = geo.jacobianInverseTransposed( quad.point( l ) );
474 const double intel = quad.weight(l) *
475 geo.integrationElement( quad.point( l ) );
478 uLF.evaluate(quad[l], result);
481 uLF.baseFunctionSet().evaluateAll( quad[l], uRets );
485 bSet.jacobianAll( quad[l], inv, gradTmpVec );
488 for(
int i=0; i<localRows; ++i)
491 const int row = startRow + i;
495 GradJacobianRangeType& gradPhi = gradTmpVec[ baseFunctionOffset( i ) ];
497 const double uDGVal = result * gradPhi[0];
498 rhs[row] += uDGVal * intel;
501 for(
int j=0; j<cols; ++j)
505 const double val = uRets[ j ] * gradPhi[0];
506 matrix[row][j] += val * intel;
512 template <
class FaceBSetType,
class Gr
idType>
513 struct GetSubBaseFunctionSet
515 template <
class EntityType,
class SpaceType>
516 static inline FaceBSetType faceBaseSet(
const EntityType& en,
const SpaceType& space)
518 return space.baseFunctionSet( (en.template subEntity<1> (0) )->type() );
522 template <
class FaceBSetType,
int dim,
class CoordCont>
523 struct GetSubBaseFunctionSet< FaceBSetType, YaspGrid< dim, CoordCont > >
525 template <
class EntityType,
class SpaceType>
526 static inline FaceBSetType faceBaseSet(
const EntityType& en,
const SpaceType& space)
528 const GeometryType geoType (GeometryType::cube,dim-1);
529 return space.baseFunctionSet( geoType );
534 template <
class FaceBSetType,
int dim>
535 struct GetSubBaseFunctionSet<FaceBSetType, UGGrid<dim> >
537 template <
class EntityType,
class SpaceType>
538 static inline FaceBSetType faceBaseSet(
const EntityType& en,
const SpaceType& space)
540 const GeometryType geoType (en.geometry().type().basicType(),dim-1);
541 return space.baseFunctionSet( geoType );
546 enum { gradFuncOffset = 1 };
547 template <
class GradBaseFunctionSet>
548 int gradientBaseFct(
const GradBaseFunctionSet& gradSet)
const 550 return (gradPolOrd <= 0) ? 0 : gradSet.size() - gradFuncOffset;
553 int baseFunctionOffset(
const int i)
const 555 return i + gradFuncOffset;
558 int numberOfBubbles(
const int bubbles ,
const GeometryType& type,
559 const int allDofs,
const int allOther )
const 571 const int numBubble = allDofs - allOther ;
572 return (numBubble > 0) ? numBubble : 0;
576 int bubbleModifier(
const int dim )
const 579 return (dim - 2) * (dim - 1) + 1;
583 void project(
const DiscreteFunctionType &uDG,
584 DiscreteFunctionType & velo )
const 586 typedef typename DiscreteFunctionType::Traits::DiscreteFunctionSpaceType FunctionSpaceType;
587 enum { localBlockSize = FunctionSpaceType::localBlockSize };
589 typedef typename FunctionSpaceType::BaseFunctionSetType BaseFunctionSetType;
590 typedef typename FunctionSpaceType::GridType GridType;
591 typedef typename FunctionSpaceType::GridPartType GridPartType;
592 typedef typename FunctionSpaceType::IteratorType Iterator;
593 typedef typename GridPartType :: IntersectionIteratorType IntersectionIteratorType;
594 typedef typename GridType :: template Codim<0> :: Entity EntityType;
595 typedef typename GridType :: Traits :: LocalIdSet LocalIdSetType;
597 typedef typename FunctionSpaceType::RangeFieldType RangeFieldType;
598 typedef typename FunctionSpaceType::DomainFieldType DomainFieldType;
599 typedef typename FunctionSpaceType::RangeType RangeType;
600 typedef typename FunctionSpaceType :: DomainType DomainType;
602 enum { dim = GridType::dimension };
603 typedef typename GridType :: ctype coordType;
605 const FunctionSpaceType& space = uDG.space();
609 if(space.order() < 1 )
615 const int polOrd = 2 * space.order() + 2;
617 typedef typename FaceDiscreteSpaceType :: BaseFunctionSetType FaceBSetType ;
618 typedef typename FaceDiscreteSpaceType :: RangeType FaceRangeType;
620 typedef typename ElementDiscreteSpaceType :: BaseFunctionSetType ElementBaseSetType ;
621 typedef typename ElementGradientSpaceType :: BaseFunctionSetType GradientBaseSetType ;
623 typedef typename DiscreteFunctionType::LocalFunctionType LocalFuncType;
625 Iterator start = space.begin();
627 if( start == space.end() )
return;
630 if( space.multipleGeometryTypes() )
632 if( space.indexSet().geomTypes(0).size() > 1)
634 assert( space.indexSet().geomTypes(0).size() == 1 );
635 DUNE_THROW(NotImplemented,
"H-div projection not implemented for hybrid grids");
639 const GeometryType startType = start->type();
642 if( startType.isHexahedron() && space.order() > 1 )
644 assert( ! startType.isHexahedron() || space.order() <= 1 );
645 DUNE_THROW(NotImplemented,
"H-div projection not implemented for p > 1 on hexas! ");
648 const int desiredOrder = space.order() + bubblePModifier;
650 if( bubblePolOrd != desiredOrder )
652 assert( bubblePolOrd == desiredOrder );
653 DUNE_THROW(NotImplemented,
"H-div projection not working for " 654 << space.order() <<
" when LagrangeSpace of order "<< desiredOrder <<
" is missing");
658 LocalFuncType lf = uDG.localFunction(*start);
659 const int numDofs = lf.numDofs();
662 const FaceBSetType faceSet =
663 GetSubBaseFunctionSet<FaceBSetType,GridType>::faceBaseSet( *start , faceSpace_ );
665 const int numFaceDofs = faceSet.size();
668 const GradientBaseSetType gradSet = gradSpace_.baseFunctionSet(*start);
670 const int numGradDofs = gradientBaseFct( gradSet );
672 const Dune::ReferenceElement< coordType, dim > & refElem =
673 Dune::ReferenceElements< coordType, dim >::general( startType );
676 const int overallFaceDofs = numFaceDofs * refElem.size(1);
679 const int numBubbleDofs = (space.order() <= 1) ? 0 :
680 numberOfBubbles( elSpace_.lagrangePointSet( *start ).numDofs ( 0 ) , startType,
681 numDofs, overallFaceDofs + numGradDofs );
683 const int rows = (overallFaceDofs + numGradDofs + numBubbleDofs);
686 const int cols = numDofs;
689 std::cout << numBubbleDofs <<
" bubbleDofs | bubbleP = " << bubblePolOrd <<
"\n";
690 std::cout << numGradDofs <<
" numGradDofs \n";
691 std::cout << overallFaceDofs <<
" allFAceDofs \n";
692 std::cout <<
"cols " << cols <<
" | rows " << rows <<
"\n";
698 typedef FieldMatrix<RangeFieldType,localBlockSize,localBlockSize> FieldMatrixType;
703 typedef FieldVector<RangeFieldType,localBlockSize> VectorType;
704 VectorType fRhs(0.0);
706 assert( numDofs == localBlockSize );
707 if( numDofs != localBlockSize )
709 DUNE_THROW(InvalidStateException,
"wrong sizes ");
713 const bool nonSymetric = (cols != rows);
717 MatrixType matrix(rows,cols);
718 typedef typename MatrixType :: RowType
RowType;
720 MatrixType fakeMatrix(cols,cols);
721 RowType rhs(rows,0.0);
722 RowType fakeRhs(numDofs,0.0);
725 for(
const auto & en : space)
733 for(
int i=0; i<rows; ++i)
739 fillMatrix(gridPart_,en,uDG,
741 gradSpace_, overallFaceDofs,
742 elSpace_, rows - numBubbleDofs,
743 polOrd,numDofs,numFaceDofs,
744 rets,uRets, matrix,rhs);
747 matrix.multTransposed(rhs, fakeRhs);
748 fakeMatrix.multiply_AT_A(matrix);
751 for(
int i=0; i<numDofs; ++i)
753 fRhs[i] = fakeRhs[i];
754 for(
int j=0; j<numDofs; ++j)
756 inv[i][j] = fakeMatrix[i][j];
767 assert( cols == rows );
769 fillMatrix(gridPart_,en,uDG,
771 gradSpace_, rows - numBubbleDofs - numGradDofs,
772 elSpace_, rows - numBubbleDofs,
773 polOrd,numDofs,numFaceDofs,
774 rets, uRets, inv,fRhs);
779 LocalFuncType veloLF = velo.localFunction( en );
786 luSolve( inv, fRhs );
787 const VectorType& x = fRhs ;
790 for(
int i=0; i<localBlockSize; ++i)
792 veloLF[ i ] = x[ i ];
798 template <
class GridPartType,
804 void fillMatrix(
const GridPartType& gridPart,
805 const EntityType& en,
806 const DiscreteFunctionType& uDG,
807 const FaceDiscreteSpaceType& faceSpace,
808 const ElementGradientSpaceType& gradSpace,
const int startGradDofs,
809 const ElementDiscreteSpaceType& elSpace,
const int startBubbleDofs,
810 const int polOrd,
const int numDofs,
const int numFaceDofs,
811 ArrayType& rets, Array2Type& uRets,
812 MatrixType& matrix, VectorType& rhs)
const 814 typedef typename GridPartType :: IntersectionIteratorType IntersectionIteratorType;
815 typedef typename IntersectionIteratorType::Intersection IntersectionType;
817 typedef typename FaceDiscreteSpaceType :: BaseFunctionSetType FaceBSetType ;
818 typedef typename FaceDiscreteSpaceType :: RangeType FaceRangeType;
819 FaceRangeType faceVal;
821 typedef typename DiscreteFunctionSpaceType::RangeType RangeType;
823 RangeType neighRet (0.0);
824 RangeType uPhi (0.0);
826 typedef typename DiscreteFunctionType :: LocalFunctionType LocalFuncType ;
827 typedef typename DiscreteFunctionType :: DiscreteFunctionSpaceType
828 :: BaseFunctionSetType BaseFunctionSetType;
831 const LocalFuncType uLF = uDG.localFunction(en);
834 const BaseFunctionSetType & bSet = uLF.baseFunctionSet();
837 IntersectionIteratorType endnit = gridPart.iend(en);
838 for(IntersectionIteratorType nit = gridPart.ibegin(en);
839 nit != endnit; ++nit)
841 const IntersectionType& inter=*nit;
843 const FaceBSetType &faceSet = faceSpace.baseFunctionSet( inter.type() );
845 const int firstRow = inter.indexInInside() * numFaceDofs;
854 const LocalFuncType uNeighLf = uDG.localFunction(nb);
859 if( inter.conforming() )
862 FaceQuadratureType faceQuadInner(gridPart, inter, polOrd, FaceQuadratureType::INSIDE);
863 FaceQuadratureType faceQuadOuter(gridPart, inter, polOrd, FaceQuadratureType::OUTSIDE);
865 applyLocalNeighbor(nit,faceQuadInner,faceQuadOuter,
866 bSet,faceSet, uLF, uNeighLf,
867 firstRow, numFaceDofs,
869 ret,neighRet,faceVal,
875 typedef typename FaceQuadratureType ::
876 NonConformingQuadratureType NonConformingQuadratureType;
878 NonConformingQuadratureType faceQuadInner(gridPart, inter, polOrd, FaceQuadratureType::INSIDE);
879 NonConformingQuadratureType faceQuadOuter(gridPart, inter, polOrd, FaceQuadratureType::OUTSIDE);
881 applyLocalNeighbor(nit,faceQuadInner,faceQuadOuter,
882 bSet,faceSet, uLF, uNeighLf,
883 firstRow, numFaceDofs,
885 ret,neighRet,faceVal,
895 FaceQuadratureType faceQuadInner(gridPart, inter, polOrd, FaceQuadratureType::INSIDE);
896 const int quadNop = faceQuadInner.nop();
898 std::vector< RangeType > uPhiVec( numDofs );
899 std::vector< FaceRangeType > faceValVec( numFaceDofs );
901 for (
int l = 0; l < quadNop ; ++l)
903 DomainType unitNormal =
904 inter.integrationOuterNormal(faceQuadInner.localPoint(l));
906 const double faceVol = unitNormal.two_norm();
907 unitNormal *= 1.0/faceVol;
910 const double intel = faceVol * faceQuadInner.weight(l);
913 uLF.evaluate(faceQuadInner[l], ret);
915 RangeFieldType val = ret * unitNormal;
918 bSet.evaluateAll( faceQuadInner[l], uPhiVec );
921 for(
int i=0; i<numDofs; ++i)
923 rets[i] = uPhiVec[ i ] * unitNormal;
927 faceSet.evaluateAll( faceQuadInner[ l ], faceValVec );
930 for(
int j=0; j<numFaceDofs; ++j, ++row)
932 FaceRangeType& faceVal = faceValVec[ j ];
933 rhs[row] += val*faceVal[0];
935 for(
int i=0; i<numDofs; ++i)
937 matrix[row][i] += (faceVal[0] * rets[i]);
947 gradientPart(gradSpace, en, uLF, startGradDofs, uRets, matrix, rhs );
951 if( bubblePolOrd - bubblePModifier > 1 )
953 bubblePart(elSpace, en, uLF, startBubbleDofs, uRets, matrix, rhs);
959 template <
class MatrixType>
960 void printMatrix(
const MatrixType& matrix)
const 962 std::cout <<
"Print Matrix \n";
963 for(
size_t row = 0; row < matrix.N(); ++row)
965 std::cout << row <<
": ";
966 for(
size_t col = 0; col< matrix.M(); ++col)
968 if(
std::abs( matrix[row][col] ) < 1e-12 )
971 std::cout << matrix[row][col] <<
" ";
973 std::cout << std::endl;
975 std::cout <<
"Finished print Matrix \n";
978 template <
class IntersectionIteratorType,
979 class QuadratureType,
980 class BaseFunctionSetType,
981 class FaceBaseFunctionSetType,
982 class LocalFunctionType,
988 static void applyLocalNeighbor(
const IntersectionIteratorType& nit,
989 const QuadratureType& faceQuadInner,
990 const QuadratureType& faceQuadOuter,
991 const BaseFunctionSetType& bSet,
992 const FaceBaseFunctionSetType& faceSet,
993 const LocalFunctionType& uLF,
994 const LocalFunctionType& uNeighLf,
996 const int numFaceDofs,
998 RangeType& ret, RangeType& neighRet,
999 FaceRangeType& faceVal,
1003 const typename IntersectionIteratorType::Intersection& inter = *nit;
1004 const int quadNop = faceQuadInner.nop();
1005 const int numDofs = uLF.numDofs();
1007 std::vector< RangeType > phiVec( numDofs );
1008 std::vector< FaceRangeType > facePhiVec( numFaceDofs );
1010 for (
int l = 0; l < quadNop ; ++l)
1012 DomainType unitNormal =
1013 inter.integrationOuterNormal(faceQuadInner.localPoint(l));
1016 const double faceVol = unitNormal.two_norm();
1017 unitNormal *= 1.0/faceVol;
1020 const double intel = faceVol * faceQuadInner.weight(l);
1023 uLF.evaluate(faceQuadInner[l], ret);
1024 uNeighLf.evaluate(faceQuadOuter[l], neighRet);
1030 RangeFieldType val = ret * unitNormal;
1033 bSet.evaluateAll( faceQuadInner[l], phiVec );
1036 for(
int i=0; i<numDofs; ++i)
1038 rets[i] = phiVec[ i ] * unitNormal;
1043 faceSet.evaluateAll( faceQuadInner[ l ], facePhiVec );
1046 for(
int j=0; j<numFaceDofs; ++j, ++row )
1048 FaceRangeType& faceVal = facePhiVec[ j ];
1050 rhs[row] += val * faceVal[0];
1052 for(
int i=0; i<numDofs; ++i)
1054 matrix[row][i] += (faceVal[0] * rets[i]);
1063 DiscreteFunctionType& dest)
const 1069 template <
class AdaptationType>
1071 AdaptationType& adaptation)
1073 typedef typename DiscreteFunctionType::LocalFunctionType LocalFuncType;
1074 typedef typename GridType :: Traits :: LocalIdSet LocalIdSetType;
1076 enum { dim = GridType :: dimension };
1077 typedef typename GridType :: template Codim<0> :: Entity EntityType;
1079 const DiscreteFunctionSpaceType& space = velo.space();
1080 GridPartType & gridPart = space.gridPart();
1081 int polOrd = space.order() + 1;
1084 const LocalIdSetType& localIdSet = gridPart.grid().localIdSet();
1089 for(
const auto & en : space)
1091 const LocalFuncType uLF = velo.localFunction(en);
1092 const double enVol = en.geometry().volume();
1094 typedef typename GridPartType :: IntersectionIteratorType IntersectionIteratorType;
1095 typedef typename IntersectionIteratorType :: Intersection IntersectionType;
1096 IntersectionIteratorType endnit = gridPart.iend(en);
1097 for(IntersectionIteratorType nit = gridPart.ibegin(en);
1098 nit != endnit; ++nit)
1100 const IntersectionType& inter=*nit;
1101 double enError = 0.0;
1103 if(inter.neighbor())
1106 const double enVol_nbVol = 0.5 * (enVol + nb.geometry().volume());
1110 const bool interiorEntity = (nb.partitionType() == InteriorEntity);
1112 const bool interiorEntity =
true;
1114 if( (localIdSet.id( en ) < localIdSet.id( nb ))
1116 || ( ! interiorEntity )
1120 const LocalFuncType uNeighLf = velo.localFunction(nb);
1125 if( inter.conforming() )
1128 FaceQuadratureType faceQuadInner(gridPart, inter, polOrd, FaceQuadratureType::INSIDE);
1129 FaceQuadratureType faceQuadOuter(gridPart, inter, polOrd, FaceQuadratureType::OUTSIDE);
1131 applyLocalNeighEstimator(nit,nb,faceQuadInner,faceQuadOuter,
1132 uLF, uNeighLf, enVol_nbVol, interiorEntity,
1133 enError, adaptation);
1138 typedef typename FaceQuadratureType ::
1139 NonConformingQuadratureType NonConformingQuadratureType;
1141 NonConformingQuadratureType faceQuadInner(gridPart, inter, polOrd, FaceQuadratureType::INSIDE);
1142 NonConformingQuadratureType faceQuadOuter(gridPart, inter, polOrd, FaceQuadratureType::OUTSIDE);
1144 applyLocalNeighEstimator(nit,nb,faceQuadInner,faceQuadOuter,
1145 uLF, uNeighLf, enVol_nbVol, interiorEntity,
1146 enError, adaptation);
1154 adaptation.addToLocalIndicator( en , enError );
1161 template <
class IntersectionIteratorType,
1163 class QuadratureType,
1164 class LocalFunctionType,
1165 class AdaptationType>
1166 static inline void applyLocalNeighEstimator(
const IntersectionIteratorType& nit,
1167 const EntityType& nb,
1168 const QuadratureType& faceQuadInner,
1169 const QuadratureType& faceQuadOuter,
1170 const LocalFunctionType& uLF,
1171 const LocalFunctionType& uNeighLf,
1172 const double enVol_nbVol,
1173 const bool interiorEntity,
1175 AdaptationType& adaptation)
1177 const typename IntersectionIteratorType::Intersection& inter=*nit;
1178 enum { dim = GridType :: dimension };
1181 double nbError = 0.0;
1183 const int quadNop = faceQuadInner.nop();
1184 for (
int l = 0; l < quadNop ; ++l)
1186 DomainType unitNormal =
1187 inter.integrationOuterNormal(faceQuadInner.localPoint(l));
1189 double faceVol = unitNormal.two_norm();
1190 unitNormal *= 1.0/faceVol;
1195 const double power = (1.0/(dim-1));
1196 faceVol = pow(faceVol, power);
1200 uLF.evaluate(faceQuadInner[l], jump);
1201 uNeighLf.evaluate(faceQuadOuter[l], neighRet);
1206 const double weight = (enVol_nbVol) * faceQuadInner.weight(l);
1208 double error = weight * SQR(jump * unitNormal);
1222 adaptation.addToLocalIndicator( nb , nbError );
1231 template <
class MatrixType,
class VectorType>
1232 void luSolve(MatrixType& a,
1233 VectorType& x)
const 1235 typedef typename VectorType :: field_type ctype;
1236 enum { n = VectorType :: dimension };
1239 assert( a.N() == a.M() );
1240 assert( a.N() == n );
1245 for(
int i=0; i<n-1; ++i)
1252 for(
int k=i; k<n; ++k)
1256 max_abs = fabs(a[k][i]);
1264 for(
int j=0; j<n; ++j)
1266 const ctype tmp = a[ p[i] ][j];
1267 a[ p[i] ][j] = a[i][j];
1273 for(
int k=i+1; k<n; ++k)
1275 const ctype lambda = a[k][i] / a[i][i];
1277 for(
int j=i+1; j<n; ++j) a[k][j] -= lambda * a[i][j];
1282 for(
int i=0; i<n-1; ++i)
1284 const ctype tmp = x[ i ];
1285 x[ i ] = x[ p[ i ] ];
1290 for(
int i=0; i<n; ++i)
1293 for(
int j=0; j<i; ++j) dot += a[i][j] * x[j];
1298 for(
int i=n-1; i>=0; --i)
1301 for(
int j=i+1; j<n; ++j) dot += a[i][j] * x[j];
1302 x[i] = (x[i] - dot) / a[i][i];
1310 #endif // #ifndef DUNE_FEM_HDIV_PROJECTION_HH
Definition: subobjects.hh:57
static double normalJump(const DiscreteFunctionType &discFunc, const int polyOrder=-1)
return sum of jumps of discrete function normal to intersection
Definition: hdivprojection.hh:208
HdivProjection(const DiscreteFunctionSpaceType &space)
constructor taking space
Definition: hdivprojection.hh:180
Lagrange discrete function space.
Definition: cornerstorage.hh:23
DestinationType::DiscreteFunctionSpaceType SpaceType
type of discrete function space
Definition: spaceoperatorif.hh:113
A vector valued function space.
Definition: functionspace.hh:16
const DiscreteFunctionSpaceType & space() const
return reference to space
Definition: hdivprojection.hh:197
void setTime(double)
set time for operators
Definition: hdivprojection.hh:201
double timeStepEstimate() const
estimate maximum time step
Definition: hdivprojection.hh:203
DenseMatrix based on std::vector< std::vector< T > >
Definition: blockmatrix.hh:23
static constexpr T sum(T a)
Definition: utility.hh:33
Double abs(const Double &a)
Definition: double.hh:860
Dune::EntityPointer< Grid, Implementation >::Entity make_entity(const Dune::EntityPointer< Grid, Implementation > &entityPointer)
Definition: compatibility.hh:23
interface for time evolution operators
Definition: spaceoperatorif.hh:101
HdivProjection(const HdivProjection &org)
Definition: hdivprojection.hh:188
Definition: coordinate.hh:4
static void estimator(const DiscreteFunctionType &velo, AdaptationType &adaptation)
Definition: hdivprojection.hh:1070
quadrature on the codim-0 reference element
Definition: elementquadrature.hh:47
H-div Projection for discontinuous discrete functions. The projection is described in detail in: ...
Definition: hdivprojection.hh:58
virtual void operator()(const DiscreteFunctionType &arg, DiscreteFunctionType &dest) const
application operator projection arg to H-div space
Definition: hdivprojection.hh:1062
DofStoragePolicy
Definition: dofstorage.hh:16