dune-fem  2.4.1-rc
hdivprojection.hh
Go to the documentation of this file.
1 #ifndef DUNE_FEM_HDIV_PROJECTION_HH
2 #define DUNE_FEM_HDIV_PROJECTION_HH
3 
4 //- Dune includes
5 #include <dune/geometry/referenceelements.hh>
6 
7 //- Dune-fem includes
14 
15 // make sure higher order Lagrange works (define USE_TWISTFREE_MAPPER)
18 
19 namespace Dune
20 {
21 
22  // External Forward Declarations
23  // -----------------------------
24 
25  template< int dim, class CoordCont >
26  class YaspGrid;
27 
28 #ifdef ENABLE_UG
29  template< int dim >
30  class UGGrid;
31 #endif // #ifdef ENABLE_UG
32 
33  namespace Fem
34  {
35 
57  template <class DiscreteFunctionType>
58  class HdivProjection : public SpaceOperatorInterface<DiscreteFunctionType>
59  {
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;
68 
70  //typedef ElementQuadrature <GridPartType , 0> VolumeQuadratureType;
71 
72  // face quadrature type
73  //typedef CachingQuadrature<GridPartType, 1> FaceQuadratureType;
75 
76  typedef typename GridPartType :: GridType GridType;
77 
78  enum { dimFaceRange = 1 };
79  enum { dimDomain = DiscreteFunctionSpaceType::dimDomain };
80  enum { dimFaceDomain = dimDomain - 1 };
81  enum { polOrdN = DiscreteFunctionSpaceType :: polynomialOrder };
82 
85 
86  enum { gradPolOrd = ((polOrdN - 1) < 0) ? 0 : (polOrdN - 1) };
87 
88  template <class Space> struct BubbleM;
89 
90  template <class FunctionSpaceImp,
91  class GridPartImp,
92  int polOrd,
93  template <class> class StorageImp,
94  template <class,class,int,template <class> class> class DiscreteFunctionSpaceImp>
95  struct BubbleM <DiscreteFunctionSpaceImp<FunctionSpaceImp,GridPartImp,polOrd,StorageImp> >
96  {
97  enum { bubblePModifier = dimDomain - 1 };
98  };
99 
100  template <class FunctionSpaceImp,
101  class GridPartImp,
102  int polOrd,
103  template <class> class StorageImp>
104  struct BubbleM < LegendreDiscontinuousGalerkinSpace<FunctionSpaceImp,GridPartImp,polOrd,StorageImp> >
105  {
106  enum { bubblePModifier = 1 };
107  };
108 
109  template <class DiscreteFunctionSpaceImp,
110  int N,
111  DofStoragePolicy policy>
112  struct BubbleM <CombinedSpace<DiscreteFunctionSpaceImp,N,policy> >
113  {
114  enum { bubblePModifier = BubbleM< DiscreteFunctionSpaceImp > :: bubblePModifier };
115  };
116 
117  typedef BubbleM <DiscreteFunctionSpaceType> BubbleMType;
118 
119 #ifdef USE_TWISTFREE_MAPPER
120  // modifier for bubble pol order
121  enum { bubblePModifier = BubbleMType :: bubblePModifier };
122 
123  // we need polOrd + 1 in 2d and polOrd + 2 in 3d
124  enum { bubblePolOrd = (polOrdN + bubblePModifier) };
125 #else
126 #ifndef NDEBUG
127 #warning "Hdiv-Projection only working for polOrd = 1 (enable higher order Lagrange with -DUSE_TWISTFREE_MAPPER)"
128 #endif
129  // modifier for bubble pol order
130  enum { bubblePModifier = 1 };
131 
132  // limit bubble polord otherwise no compilation possible
133  enum { bubblePolOrd = (polOrdN > 1) ? 2 : (polOrdN + 1) };
134 #endif
135 
136  template <class Space> struct Spaces;
137 
138  template <class FunctionSpaceImp,
139  class GridPartImp,
140  int polOrd,
141  template <class> class StorageImp,
142  template <class,class,int,template <class> class> class DiscreteFunctionSpaceImp>
143  struct Spaces<DiscreteFunctionSpaceImp<FunctionSpaceImp,GridPartImp,polOrd,StorageImp> >
144  {
145  //typedef LagrangeDiscreteFunctionSpace<ElSpaceType,GridPartImp,gradPolOrd,StorageImp> ElementGradientSpaceType;
146  //typedef LegendreDiscontinuousGalerkinSpace<ElSpaceType,GridPartImp,gradPolOrd,StorageImp> ElementGradientSpaceType;
147  typedef DiscreteFunctionSpaceImp<ElSpaceType,GridPartImp,gradPolOrd,StorageImp> ElementGradientSpaceType;
148  //typedef DiscontinuousGalerkinSpace<ElSpaceType,GridPartImp,gradPolOrd,StorageImp> ElementGradientSpaceType;
149  typedef DiscreteFunctionSpaceImp<FaceSpaceType,GridPartType,polOrdN,StorageImp> FaceDiscreteSpaceType;
150  //typedef LegendreDiscontinuousGalerkinSpace<FaceSpaceType,GridPartType,polOrdN,StorageImp> FaceDiscreteSpaceType;
151  //typedef DiscontinuousGalerkinSpace<FaceSpaceType,GridPartType,polOrdN,StorageImp> FaceDiscreteSpaceType;
153  };
154 
155  template <class DiscreteFunctionSpaceImp,
156  int N,
157  DofStoragePolicy policy>
158  struct Spaces<CombinedSpace<DiscreteFunctionSpaceImp,N,policy> >
159  {
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;
165  };
166 
167  typedef Spaces<DiscreteFunctionSpaceType> AllSpacesType;
168  typedef typename AllSpacesType :: ElementGradientSpaceType ElementGradientSpaceType;
169  typedef typename AllSpacesType :: ElementDiscreteSpaceType ElementDiscreteSpaceType;
170  typedef typename AllSpacesType :: FaceDiscreteSpaceType FaceDiscreteSpaceType;
171 
172  const DiscreteFunctionSpaceType& space_;
173  GridPartType & gridPart_;
174  const FaceDiscreteSpaceType faceSpace_;
175  const ElementDiscreteSpaceType elSpace_;
176  const ElementGradientSpaceType gradSpace_;
177 
178  public:
180  HdivProjection(const DiscreteFunctionSpaceType& space) :
181  space_(space),
182  gridPart_(space.gridPart()),
183  faceSpace_( gridPart_ ),
184  elSpace_( gridPart_ ),
185  gradSpace_( gridPart_ )
186  {}
187 
189  space_(org.space_),
190  gridPart_( org.gridPart_),
191  faceSpace_( gridPart_ ),
192  elSpace_( gridPart_ ),
193  gradSpace_( gridPart_ )
194  {}
195 
197  const DiscreteFunctionSpaceType& space() const
198  {
199  return space_;
200  }
201  void setTime(double) {
202  }
203  double timeStepEstimate() const {
204  return 0.;
205  }
206 
208  static double normalJump(const DiscreteFunctionType &discFunc, const int polyOrder = -1 )
209  {
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;
214 
215  enum { dim = GridType::dimension };
216 
217  const DiscreteFunctionSpaceType &space = discFunc.space();
218  const GridPartType &gridPart = space.gridPart();
219  const int polOrd = (polyOrder <0) ? (2 * space.order() + 2) : polyOrder;
220 
221  typedef typename DiscreteFunctionType::LocalFunctionType LocalFuncType;
222 
223  RangeType ret (0.0);
224  RangeType neighRet (0.0);
225 
226  // define type of quadrature
227  typedef ElementQuadrature <GridPartType , 1> FaceQuadratureType;
228 
229  double sum = 0.0;
230 
231  const LocalIdSetType &idSet = gridPart.grid().localIdSet();
232 
233  for(const auto & en : space)
234  {
235  const LocalFuncType lf = discFunc.localFunction(en);
236 
237  double localValue = 0.0;
238 
239  IntersectionIteratorType endnit = gridPart.iend(en);
240  for(IntersectionIteratorType nit = gridPart.ibegin(en);
241  nit != endnit; ++nit)
242  {
243  const IntersectionType& inter=*nit;
244  // only interior faces are considered
245  if(inter.neighbor() )
246  {
247  EntityType nb = make_entity( inter.outside() );
248 
249  if(idSet.id( en ) < idSet.id( nb ))
250  {
251  FaceQuadratureType faceQuadInner(gridPart, inter, polOrd, FaceQuadratureType::INSIDE);
252  FaceQuadratureType faceQuadOuter(gridPart, inter, polOrd, FaceQuadratureType::OUTSIDE);
253 
254  const LocalFuncType neighLf = discFunc.localFunction(nb);
255 
256  const int quadNop = faceQuadInner.nop();
257  for (int l = 0; l < quadNop ; ++l)
258  {
259  DomainType normal =
260  inter.unitOuterNormal(faceQuadInner.localPoint(l));
261 
262  lf.evaluate(faceQuadInner[l], ret);
263  neighLf.evaluate(faceQuadOuter[l], neighRet);
264 
265  ret -= neighRet;
266 
267  double val = ret * normal;
268  val *= faceQuadInner.weight(l);
269 
270  localValue += std::abs(val);
271  }
272  }
273  }
274  } // end of intersection iterator
275  sum += std::abs(localValue);
276  }
277 
278  return sum;
279  }
280 
281  private:
282 
283  // only works for 2d right now
284  void curl(const DomainType & arg, DomainType & dest, const int d ) const
285  {
286  if( DomainType :: dimension == 2 )
287  {
288  dest[0] = arg[1];
289  dest[1] = -arg[0];
290 
291  return ;
292  }
293  else if( DomainType :: dimension == 3 )
294  {
295  if( d == 0 )
296  {
297  dest[0] = arg[1];
298  dest[1] = -arg[2];
299  dest[2] = 0;
300  return ;
301  }
302  else if ( d == 1 )
303  {
304  dest[0] = 0;
305  dest[1] = arg[2];
306  dest[2] = -arg[0];
307  return ;
308  }
309  else
310  {
311  dest[0] = -arg[2];
312  dest[1] = 0;
313  dest[2] = arg[1];
314  return ;
315  }
316  }
317  else
318  {
319  assert( false );
320  abort();
321  }
322  }
323 
324  template <class EntityType,
325  class LocalFunctionType,
326  class ArrayType,
327  class MatrixType,
328  class VectorType>
329  void bubblePart(const ElementDiscreteSpaceType& space,
330  EntityType & en,
331  const LocalFunctionType & uLF, const int startRow ,
332  ArrayType& uRets,
333  MatrixType & matrix, VectorType& rhs ) const
334  {
335 
336  typedef typename ElementDiscreteSpaceType :: BaseFunctionSetType BaseFunctionSetType;
337  typedef typename ElementDiscreteSpaceType :: LagrangePointSetType LagrangePointSetType;
338  typedef typename ElementDiscreteSpaceType :: Traits::GridPartType GridPartType;
339 
340  typedef typename ElementDiscreteSpaceType :: JacobianRangeType JacobianRangeType;
341  typedef typename ElementDiscreteSpaceType :: DomainType DomainType;
342 
343  enum { dim = GridType::dimension };
344 
345  const LagrangePointSetType& lagrangePointSet = space.lagrangePointSet( en );
346  const BaseFunctionSetType bSet = space.baseFunctionSet( en );
347  const int polOrd = 2 * space.order(); // + 2;
348 
349  const int cols = uLF.numDofs();
350  assert( uRets.size() == (unsigned int)cols );
351 
352  VolumeQuadratureType quad (en,polOrd);
353  DomainType result;
354  std::vector< JacobianRangeType > valTmpVec( bSet.size() );
355  DomainType bVal;
356  DomainType aVal;
357 
358  // get geometry
359  typedef typename EntityType :: Geometry Geometry ;
360  const Geometry& geo = en.geometry();
361  // get geometry type
362  const GeometryType& type = geo.type();
363  const int bubbleOffset = (type.isSimplex()) ? 0 : baseFunctionOffset( 0 );
364 
365  // type of jacobian inverse
366  typedef typename Geometry :: ctype ctype;
367  enum { cdim = Geometry :: coorddimension };
368  enum { mydim = Geometry :: mydimension };
369  typedef typename Geometry::JacobianInverseTransposed JacobianInverseType;
370 
371  // get number of dofs for codim 0 (skip first for)
372  const int enDofs = numberOfBubbles( lagrangePointSet.numDofs( 0 ), type ,
373  cols, startRow );
374  const int bubbleMod = bubbleModifier( mydim );
375 
376  const int quadNop = quad.nop();
377  for (int l = 0; l < quadNop ; ++l)
378  {
379  // get jacobian inverse
380  const JacobianInverseType& inv = geo.jacobianInverseTransposed( quad.point(l) );
381 
382  // get integration element
383  const double intel = quad.weight(l) *
384  geo.integrationElement(quad.point(l));
385 
386  // evaluate u
387  uLF.evaluate(quad[l], result);
388 
389  // evaluate base functions of u
390  uLF.baseFunctionSet().evaluateAll( quad[l], uRets );
391 
392  // evaluate gradients
393  bSet.jacobianAll( quad[l], inv, valTmpVec );
394 
395  // for all bubble functions
396  for( int i = 0 ; i<enDofs; i += bubbleMod )
397  {
398  // we might have other row
399  int row = startRow + i;
400 
401  // map to lagrange base function number
402  const int localBaseFct = ((int) i/bubbleMod) + bubbleOffset;
403  // get dof number of 'localBaseFct' dof on codim 0 subentity 0
404  const int baseFct = lagrangePointSet.entityDofNumber( 0, 0, localBaseFct );
405 
406  // evaluate gradient
407  JacobianRangeType& val = valTmpVec[ baseFct ];
408 
409  //apply inverse
410  //inv.mv( valTmp[0], val[0] );
411 
412  for(int d = 0; d<bubbleMod; ++d )
413  {
414  // apply curl
415  curl(val[0], aVal, d );
416 
417  double r = aVal * result;
418  r *= intel;
419  rhs[row] += r;
420 
421  // for cols make matrix
422  for(int j=0; j<cols; ++j)
423  {
424  double t = aVal * uRets[ j ];
425  t *= intel;
426  matrix[ row ][ j ] += t;
427  }
428 
429  // increase row
430  ++row;
431  }
432  }
433  }
434  }
435 
436  template <class EntityType, class LocalFunctionType, class ArrayType,
437  class MatrixType, class VectorType>
438  void gradientPart(const ElementGradientSpaceType & space,
439  EntityType & en,
440  const LocalFunctionType & uLF, const int startRow ,
441  ArrayType& uRets, MatrixType& matrix, VectorType& rhs ) const
442  {
443  typedef typename ElementGradientSpaceType :: BaseFunctionSetType BaseFunctionSetType;
444  const BaseFunctionSetType bSet = space.baseFunctionSet( en );
445  int polOrd = 2 * space.order() + 1;
446 
447  const int localRows = gradientBaseFct( bSet );
448  const int cols = uLF.numDofs();
449 
450  VolumeQuadratureType quad (en,polOrd);
451 
452  RangeType result;
453  RangeType uPhi;
454 
455  typedef typename ElementGradientSpaceType :: JacobianRangeType GradJacobianRangeType;
456  std::vector< GradJacobianRangeType > gradTmpVec( bSet.size() );
457  GradJacobianRangeType gradPhi;
458 
459  typedef typename EntityType :: Geometry Geometry ;
460  const Geometry& geo = en.geometry();
461 
462  typedef typename Geometry :: ctype ctype;
463  enum { cdim = Geometry :: coorddimension };
464  enum { mydim = Geometry :: mydimension };
465  typedef typename Geometry::JacobianInverseTransposed JacobianInverseType;
466 
467  const int quadNop = quad.nop();
468  for (int l = 0; l < quadNop ; ++l)
469  {
470  // get jacobian inverse
471  const JacobianInverseType& inv = geo.jacobianInverseTransposed( quad.point( l ) );
472 
473  // get integration element
474  const double intel = quad.weight(l) *
475  geo.integrationElement( quad.point( l ) );
476 
477  // evaluate uLF
478  uLF.evaluate(quad[l], result);
479 
480  // evaluate base function on quadrature point
481  uLF.baseFunctionSet().evaluateAll( quad[l], uRets );
482 
483  // evaluate gradient (skip first function because this function
484  // is constant and the gradient therefore 0 )
485  bSet.jacobianAll( quad[l], inv, gradTmpVec );
486 
487  // apply jacobian Inverse
488  for(int i=0; i<localRows; ++i)
489  {
490  // we might have other row
491  const int row = startRow + i;
492 
493  // evaluate gradient (skip first function because this function
494  // is constant and the gradient therefore 0 )
495  GradJacobianRangeType& gradPhi = gradTmpVec[ baseFunctionOffset( i ) ];
496 
497  const double uDGVal = result * gradPhi[0];
498  rhs[row] += uDGVal * intel;
499 
500  // for cols make matrix
501  for(int j=0; j<cols; ++j)
502  {
503  //uLF.baseFunctionSet().evaluate(j, quad[l], uPhi);
504  //const double val = uPhi * gradPhi[0];
505  const double val = uRets[ j ] * gradPhi[0];
506  matrix[row][j] += val * intel;
507  }
508  }
509  }
510  }
511 
512  template <class FaceBSetType, class GridType>
513  struct GetSubBaseFunctionSet
514  {
515  template <class EntityType, class SpaceType>
516  static inline FaceBSetType faceBaseSet(const EntityType& en, const SpaceType& space)
517  {
518  return space.baseFunctionSet( (en.template subEntity<1> (0) )->type() );
519  }
520  };
521 
522  template <class FaceBSetType, int dim, class CoordCont>
523  struct GetSubBaseFunctionSet< FaceBSetType, YaspGrid< dim, CoordCont > >
524  {
525  template <class EntityType, class SpaceType>
526  static inline FaceBSetType faceBaseSet(const EntityType& en, const SpaceType& space)
527  {
528  const GeometryType geoType (GeometryType::cube,dim-1);
529  return space.baseFunctionSet( geoType );
530  }
531  };
532 
533 #ifdef ENABLE_UG
534  template <class FaceBSetType, int dim>
535  struct GetSubBaseFunctionSet<FaceBSetType, UGGrid<dim> >
536  {
537  template <class EntityType, class SpaceType>
538  static inline FaceBSetType faceBaseSet(const EntityType& en, const SpaceType& space)
539  {
540  const GeometryType geoType (en.geometry().type().basicType(),dim-1);
541  return space.baseFunctionSet( geoType );
542  }
543  };
544 #endif
545 
546  enum { gradFuncOffset = 1 };
547  template <class GradBaseFunctionSet>
548  int gradientBaseFct(const GradBaseFunctionSet& gradSet) const
549  {
550  return (gradPolOrd <= 0) ? 0 : gradSet.size() - gradFuncOffset;
551  }
552 
553  int baseFunctionOffset(const int i) const
554  {
555  return i + gradFuncOffset;
556  }
557 
558  int numberOfBubbles( const int bubbles , const GeometryType& type,
559  const int allDofs, const int allOther ) const
560  {
561  /*
562  // for hexahedrons this is different
563  if( type.isHexahedron() )
564  {
565  return (bubbleModifier( type.dim() )) * (bubbles - gradFuncOffset) - 1;
566  }
567  else
568  */
569  {
570  // the rest is padded with bubble functions
571  const int numBubble = allDofs - allOther ;
572  return (numBubble > 0) ? numBubble : 0;
573  }
574  }
575 
576  int bubbleModifier( const int dim ) const
577  {
578  // return 1 for 2d and 3 for 3d
579  return (dim - 2) * (dim - 1) + 1;
580  }
581 
583  void project(const DiscreteFunctionType &uDG,
584  DiscreteFunctionType & velo ) const
585  {
586  typedef typename DiscreteFunctionType::Traits::DiscreteFunctionSpaceType FunctionSpaceType;
587  enum { localBlockSize = FunctionSpaceType::localBlockSize };
588 
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;
596 
597  typedef typename FunctionSpaceType::RangeFieldType RangeFieldType;
598  typedef typename FunctionSpaceType::DomainFieldType DomainFieldType;
599  typedef typename FunctionSpaceType::RangeType RangeType;
600  typedef typename FunctionSpaceType :: DomainType DomainType;
601 
602  enum { dim = GridType::dimension };
603  typedef typename GridType :: ctype coordType;
604 
605  const FunctionSpaceType& space = uDG.space();
606 
607  // for polOrd 0 this is not working
608  // then just copy the function
609  if(space.order() < 1 )
610  {
611  velo.assign( uDG );
612  return ;
613  }
614 
615  const int polOrd = 2 * space.order() + 2;
616 
617  typedef typename FaceDiscreteSpaceType :: BaseFunctionSetType FaceBSetType ;
618  typedef typename FaceDiscreteSpaceType :: RangeType FaceRangeType;
619 
620  typedef typename ElementDiscreteSpaceType :: BaseFunctionSetType ElementBaseSetType ;
621  typedef typename ElementGradientSpaceType :: BaseFunctionSetType GradientBaseSetType ;
622 
623  typedef typename DiscreteFunctionType::LocalFunctionType LocalFuncType;
624 
625  Iterator start = space.begin();
626  // if empty grid, do nothing
627  if( start == space.end() ) return;
628 
629  // only working for spaces with one element type
630  if( space.multipleGeometryTypes() )
631  {
632  if( space.indexSet().geomTypes(0).size() > 1)
633  {
634  assert( space.indexSet().geomTypes(0).size() == 1 );
635  DUNE_THROW(NotImplemented,"H-div projection not implemented for hybrid grids");
636  }
637  }
638 
639  const GeometryType startType = start->type();
640 
641  // only working for spaces with one element type
642  if( startType.isHexahedron() && space.order() > 1 )
643  {
644  assert( ! startType.isHexahedron() || space.order() <= 1 );
645  DUNE_THROW(NotImplemented,"H-div projection not implemented for p > 1 on hexas! ");
646  }
647 
648  const int desiredOrder = space.order() + bubblePModifier;
649  // check Lagrange space present
650  if( bubblePolOrd != desiredOrder )
651  {
652  assert( bubblePolOrd == desiredOrder );
653  DUNE_THROW(NotImplemented,"H-div projection not working for "
654  << space.order() << " when LagrangeSpace of order "<< desiredOrder << " is missing");
655  }
656 
657  // colums are dofs of searched function
658  LocalFuncType lf = uDG.localFunction(*start);
659  const int numDofs = lf.numDofs();
660  //std::cout << numDofs << " numDofs \n";
661 
662  const FaceBSetType faceSet =
663  GetSubBaseFunctionSet<FaceBSetType,GridType>::faceBaseSet( *start , faceSpace_ );
664  // number of dofs on faces
665  const int numFaceDofs = faceSet.size();
666 
667 
668  const GradientBaseSetType gradSet = gradSpace_.baseFunctionSet(*start);
669  // in case of linear space the is zero
670  const int numGradDofs = gradientBaseFct( gradSet );
671 
672  const Dune::ReferenceElement< coordType, dim > & refElem =
673  Dune::ReferenceElements< coordType, dim >::general( startType );
674 
675  // get number of faces
676  const int overallFaceDofs = numFaceDofs * refElem.size(1);
677 
678  // get all element dofs from Lagrange space
679  const int numBubbleDofs = (space.order() <= 1) ? 0 :
680  numberOfBubbles( elSpace_.lagrangePointSet( *start ).numDofs ( 0 ) , startType,
681  numDofs, overallFaceDofs + numGradDofs );
682 
683  const int rows = (overallFaceDofs + numGradDofs + numBubbleDofs);
684 
685  // number of columns
686  const int cols = numDofs;
687 
688 #if 0
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";
693 #endif
694 
695  MutableArray< RangeFieldType > rets(numDofs);
696  MutableArray< RangeType > uRets(numDofs);
697 
698  typedef FieldMatrix<RangeFieldType,localBlockSize,localBlockSize> FieldMatrixType;
699 
700  // resulting system matrix
701  FieldMatrixType inv;
702 
703  typedef FieldVector<RangeFieldType,localBlockSize> VectorType;
704  VectorType fRhs(0.0);
705 
706  assert( numDofs == localBlockSize );
707  if( numDofs != localBlockSize )
708  {
709  DUNE_THROW(InvalidStateException,"wrong sizes ");
710  }
711 
712  // flag to say whether we have non-symetric or symetric situation
713  const bool nonSymetric = (cols != rows);
714 
715  // matrix type
716  typedef Fem::DenseMatrix<double> MatrixType;
717  MatrixType matrix(rows,cols);
718  typedef typename MatrixType :: RowType RowType;
719  // matrix type
720  MatrixType fakeMatrix(cols,cols);
721  RowType rhs(rows,0.0);
722  RowType fakeRhs(numDofs,0.0);
723 
724  // iterate over grid
725  for(const auto & en : space)
726  {
727  if( nonSymetric )
728  {
729  // reset values
730  matrix = 0.0;
731 
732  // reset rhs
733  for(int i=0; i<rows; ++i)
734  {
735  rhs[i] = 0.0;
736  }
737 
738  // fill non-symetric matrix
739  fillMatrix(gridPart_,en,uDG,
740  faceSpace_,
741  gradSpace_, overallFaceDofs,
742  elSpace_, rows - numBubbleDofs,
743  polOrd,numDofs,numFaceDofs,
744  rets,uRets, matrix,rhs);
745 
746  // apply least square
747  matrix.multTransposed(rhs, fakeRhs);
748  fakeMatrix.multiply_AT_A(matrix);
749 
750  // copy values
751  for(int i=0; i<numDofs; ++i)
752  {
753  fRhs[i] = fakeRhs[i];
754  for(int j=0; j<numDofs; ++j)
755  {
756  inv[i][j] = fakeMatrix[i][j];
757  }
758  }
759  }
760  else
761  {
762  // reset values
763  inv = 0.0;
764  // reset rhs
765  fRhs = 0.0;
766 
767  assert( cols == rows );
768  // fill inv and fRhs directly
769  fillMatrix(gridPart_,en,uDG,
770  faceSpace_,
771  gradSpace_, rows - numBubbleDofs - numGradDofs,
772  elSpace_, rows - numBubbleDofs,
773  polOrd,numDofs,numFaceDofs,
774  rets, uRets, inv,fRhs);
775  }
776 
777  // set new values to new velocity function
778  {
779  LocalFuncType veloLF = velo.localFunction( en );
780 #if 0
781  VectorType x( 0 );
782 
783  inv.solve(x, fRhs);
784 #else
785  // solve linear system
786  luSolve( inv, fRhs );
787  const VectorType& x = fRhs ;
788 #endif
789 
790  for(int i=0; i<localBlockSize; ++i)
791  {
792  veloLF[ i ] = x[ i ];
793  }
794  }
795  }
796  }
797 
798  template <class GridPartType,
799  class EntityType,
800  class ArrayType,
801  class Array2Type,
802  class MatrixType,
803  class VectorType>
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
813  {
814  typedef typename GridPartType :: IntersectionIteratorType IntersectionIteratorType;
815  typedef typename IntersectionIteratorType::Intersection IntersectionType;
816 
817  typedef typename FaceDiscreteSpaceType :: BaseFunctionSetType FaceBSetType ;
818  typedef typename FaceDiscreteSpaceType :: RangeType FaceRangeType;
819  FaceRangeType faceVal;
820 
821  typedef typename DiscreteFunctionSpaceType::RangeType RangeType;
822  RangeType ret (0.0);
823  RangeType neighRet (0.0);
824  RangeType uPhi (0.0);
825 
826  typedef typename DiscreteFunctionType :: LocalFunctionType LocalFuncType ;
827  typedef typename DiscreteFunctionType :: DiscreteFunctionSpaceType
828  :: BaseFunctionSetType BaseFunctionSetType;
829 
830  // get uDg local on entity
831  const LocalFuncType uLF = uDG.localFunction(en);
832 
833  // get base functions set
834  const BaseFunctionSetType & bSet = uLF.baseFunctionSet();
835 
836  // iterate over intersections
837  IntersectionIteratorType endnit = gridPart.iend(en);
838  for(IntersectionIteratorType nit = gridPart.ibegin(en);
839  nit != endnit; ++nit)
840  {
841  const IntersectionType& inter=*nit;
842  // get base function set of face
843  const FaceBSetType &faceSet = faceSpace.baseFunctionSet( inter.type() );
844 
845  const int firstRow = inter.indexInInside() * numFaceDofs;
846 
847  // only interior faces are considered
848  if(inter.neighbor())
849  {
850  // get neighbor entity
851  EntityType nb = make_entity( inter.outside() );
852 
853  // get local function of neighbor
854  const LocalFuncType uNeighLf = uDG.localFunction(nb);
855 
856  //typedef TwistUtility<GridType> TwistUtilityType;
857  // for conforming situations apply Quadrature given
858  //if( TwistUtilityType::conforming(gridPart.grid(),inter) )
859  if( inter.conforming() )
860  {
861  // create quadratures
862  FaceQuadratureType faceQuadInner(gridPart, inter, polOrd, FaceQuadratureType::INSIDE);
863  FaceQuadratureType faceQuadOuter(gridPart, inter, polOrd, FaceQuadratureType::OUTSIDE);
864 
865  applyLocalNeighbor(nit,faceQuadInner,faceQuadOuter,
866  bSet,faceSet, uLF, uNeighLf,
867  firstRow, numFaceDofs,
868  rets,
869  ret,neighRet,faceVal,
870  matrix, rhs);
871  }
872  else
873  {
874  // type of quadrature for non-conforming intersections
875  typedef typename FaceQuadratureType ::
876  NonConformingQuadratureType NonConformingQuadratureType;
877  // create quadratures
878  NonConformingQuadratureType faceQuadInner(gridPart, inter, polOrd, FaceQuadratureType::INSIDE);
879  NonConformingQuadratureType faceQuadOuter(gridPart, inter, polOrd, FaceQuadratureType::OUTSIDE);
880 
881  applyLocalNeighbor(nit,faceQuadInner,faceQuadOuter,
882  bSet,faceSet, uLF, uNeighLf,
883  firstRow, numFaceDofs,
884  rets,
885  ret,neighRet,faceVal,
886  matrix, rhs);
887 
888  }
889  }
890 
891  // only interior faces are considered
892  if(inter.boundary())
893  {
894  // create quadrature
895  FaceQuadratureType faceQuadInner(gridPart, inter, polOrd, FaceQuadratureType::INSIDE);
896  const int quadNop = faceQuadInner.nop();
897 
898  std::vector< RangeType > uPhiVec( numDofs );
899  std::vector< FaceRangeType > faceValVec( numFaceDofs );
900 
901  for (int l = 0; l < quadNop ; ++l)
902  {
903  DomainType unitNormal =
904  inter.integrationOuterNormal(faceQuadInner.localPoint(l));
905 
906  const double faceVol = unitNormal.two_norm();
907  unitNormal *= 1.0/faceVol;
908 
909  // get integration weight
910  const double intel = faceVol * faceQuadInner.weight(l);
911 
912  // evaluate function
913  uLF.evaluate(faceQuadInner[l], ret);
914 
915  RangeFieldType val = ret * unitNormal;
916  val *= intel;
917 
918  bSet.evaluateAll( faceQuadInner[l], uPhiVec );
919 
920  // evaluate base functions
921  for(int i=0; i<numDofs; ++i)
922  {
923  rets[i] = uPhiVec[ i ] * unitNormal;
924  rets[i] *= intel;
925  }
926 
927  faceSet.evaluateAll( faceQuadInner[ l ], faceValVec );
928 
929  int row = firstRow;
930  for(int j=0; j<numFaceDofs; ++j, ++row)
931  {
932  FaceRangeType& faceVal = faceValVec[ j ];
933  rhs[row] += val*faceVal[0];
934 
935  for(int i=0; i<numDofs; ++i)
936  {
937  matrix[row][i] += (faceVal[0] * rets[i]);
938  }
939  }
940  }
941  }
942  }
943 
944  // add gradient part
945  if( gradPolOrd > 0 )
946  {
947  gradientPart(gradSpace, en, uLF, startGradDofs, uRets, matrix, rhs );
948  }
949 
950  // add bubble part
951  if( bubblePolOrd - bubblePModifier > 1 )
952  {
953  bubblePart(elSpace, en, uLF, startBubbleDofs, uRets, matrix, rhs);
954  }
955 
956  // printMatrix( matrix );
957  }
958 
959  template <class MatrixType>
960  void printMatrix(const MatrixType& matrix) const
961  {
962  std::cout << "Print Matrix \n";
963  for(size_t row = 0; row < matrix.N(); ++row)
964  {
965  std::cout << row << ": ";
966  for(size_t col = 0; col< matrix.M(); ++col)
967  {
968  if( std::abs( matrix[row][col] ) < 1e-12 )
969  std::cout << "0 ";
970  else
971  std::cout << matrix[row][col] << " ";
972  }
973  std::cout << std::endl;
974  }
975  std::cout << "Finished print Matrix \n";
976  }
977 
978  template <class IntersectionIteratorType,
979  class QuadratureType,
980  class BaseFunctionSetType,
981  class FaceBaseFunctionSetType,
982  class LocalFunctionType,
983  class ArrayType,
984  class RangeType,
985  class FaceRangeType,
986  class MatrixType,
987  class RHSType>
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,
995  const int firstRow,
996  const int numFaceDofs,
997  ArrayType& rets,
998  RangeType& ret, RangeType& neighRet,
999  FaceRangeType& faceVal,
1000  MatrixType& matrix,
1001  RHSType& rhs)
1002  {
1003  const typename IntersectionIteratorType::Intersection& inter = *nit;
1004  const int quadNop = faceQuadInner.nop();
1005  const int numDofs = uLF.numDofs();
1006 
1007  std::vector< RangeType > phiVec( numDofs );
1008  std::vector< FaceRangeType > facePhiVec( numFaceDofs );
1009 
1010  for (int l = 0; l < quadNop ; ++l)
1011  {
1012  DomainType unitNormal =
1013  inter.integrationOuterNormal(faceQuadInner.localPoint(l));
1014 
1015  // get unit outer normal
1016  const double faceVol = unitNormal.two_norm();
1017  unitNormal *= 1.0/faceVol;
1018 
1019  // integration weight
1020  const double intel = faceVol * faceQuadInner.weight(l);
1021 
1022  // evaluate dg velocity
1023  uLF.evaluate(faceQuadInner[l], ret);
1024  uNeighLf.evaluate(faceQuadOuter[l], neighRet);
1025 
1026  // take mean value
1027  ret += neighRet;
1028  ret *= 0.5;
1029 
1030  RangeFieldType val = ret * unitNormal;
1031  val *= intel;
1032 
1033  bSet.evaluateAll( faceQuadInner[l], phiVec );
1034 
1035  // evaluate base functions
1036  for(int i=0; i<numDofs; ++i)
1037  {
1038  rets[i] = phiVec[ i ] * unitNormal;
1039  rets[i] *= intel;
1040  }
1041 
1042  // evaluate all basis functions
1043  faceSet.evaluateAll( faceQuadInner[ l ], facePhiVec );
1044 
1045  int row = firstRow;
1046  for(int j=0; j<numFaceDofs; ++j, ++row )
1047  {
1048  FaceRangeType& faceVal = facePhiVec[ j ];
1049 
1050  rhs[row] += val * faceVal[0];
1051 
1052  for(int i=0; i<numDofs; ++i)
1053  {
1054  matrix[row][i] += (faceVal[0] * rets[i]);
1055  }
1056  }
1057  }
1058  }
1059 
1060  public:
1062  virtual void operator () (const DiscreteFunctionType &arg,
1063  DiscreteFunctionType& dest) const
1064  {
1065  // apply H-div projection
1066  project(arg,dest);
1067  }
1068 
1069  template <class AdaptationType>
1070  static void estimator(const DiscreteFunctionType &velo,
1071  AdaptationType& adaptation)
1072  {
1073  typedef typename DiscreteFunctionType::LocalFunctionType LocalFuncType;
1074  typedef typename GridType :: Traits :: LocalIdSet LocalIdSetType;
1075 
1076  enum { dim = GridType :: dimension };
1077  typedef typename GridType :: template Codim<0> :: Entity EntityType;
1078 
1079  const DiscreteFunctionSpaceType& space = velo.space();
1080  GridPartType & gridPart = space.gridPart();
1081  int polOrd = space.order() + 1;
1082 
1083  // get local id set
1084  const LocalIdSetType& localIdSet = gridPart.grid().localIdSet();
1085 
1086  // define type of face quadrature
1087  typedef CachingQuadrature<GridPartType,1> FaceQuadratureType;
1088 
1089  for(const auto & en : space)
1090  {
1091  const LocalFuncType uLF = velo.localFunction(en);
1092  const double enVol = en.geometry().volume();
1093 
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)
1099  {
1100  const IntersectionType& inter=*nit;
1101  double enError = 0.0;
1102  // only interior faces are considered
1103  if(inter.neighbor())
1104  {
1105  EntityType nb = make_entity( inter.outside() );
1106  const double enVol_nbVol = 0.5 * (enVol + nb.geometry().volume());
1107 
1108 #if HAVE_MPI
1109  // get partition type
1110  const bool interiorEntity = (nb.partitionType() == InteriorEntity);
1111 #else
1112  const bool interiorEntity = true;
1113 #endif
1114  if( (localIdSet.id( en ) < localIdSet.id( nb ))
1115 #if HAVE_MPI
1116  || ( ! interiorEntity )
1117 #endif
1118  )
1119  {
1120  const LocalFuncType uNeighLf = velo.localFunction(nb);
1121 
1122  //typedef TwistUtility<GridType> TwistUtilityType;
1123  // for conforming situations apply Quadrature given
1124  //if( TwistUtilityType::conforming(gridPart.grid(),inter) )
1125  if( inter.conforming() )
1126  {
1127  // create quadratures
1128  FaceQuadratureType faceQuadInner(gridPart, inter, polOrd, FaceQuadratureType::INSIDE);
1129  FaceQuadratureType faceQuadOuter(gridPart, inter, polOrd, FaceQuadratureType::OUTSIDE);
1130 
1131  applyLocalNeighEstimator(nit,nb,faceQuadInner,faceQuadOuter,
1132  uLF, uNeighLf, enVol_nbVol, interiorEntity,
1133  enError, adaptation);
1134  }
1135  else
1136  {
1137  // type of quadrature for non-conforming intersections
1138  typedef typename FaceQuadratureType ::
1139  NonConformingQuadratureType NonConformingQuadratureType;
1140  // create quadratures
1141  NonConformingQuadratureType faceQuadInner(gridPart, inter, polOrd, FaceQuadratureType::INSIDE);
1142  NonConformingQuadratureType faceQuadOuter(gridPart, inter, polOrd, FaceQuadratureType::OUTSIDE);
1143 
1144  applyLocalNeighEstimator(nit,nb,faceQuadInner,faceQuadOuter,
1145  uLF, uNeighLf, enVol_nbVol, interiorEntity,
1146  enError, adaptation);
1147  }
1148 
1149  } // end enId < nbId
1150  } // end neighbor
1151 
1152  if(enError > 0.0)
1153  {
1154  adaptation.addToLocalIndicator( en , enError );
1155  }
1156  } // end intersection iterator
1157  }
1158  }
1159 
1160  private:
1161  template <class IntersectionIteratorType,
1162  class EntityType,
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,
1174  double& enError,
1175  AdaptationType& adaptation)
1176  {
1177  const typename IntersectionIteratorType::Intersection& inter=*nit;
1178  enum { dim = GridType :: dimension };
1179  RangeType jump;
1180  RangeType neighRet;
1181  double nbError = 0.0;
1182 
1183  const int quadNop = faceQuadInner.nop();
1184  for (int l = 0; l < quadNop ; ++l)
1185  {
1186  DomainType unitNormal =
1187  inter.integrationOuterNormal(faceQuadInner.localPoint(l));
1188 
1189  double faceVol = unitNormal.two_norm();
1190  unitNormal *= 1.0/faceVol;
1191 
1192  // in case of power != 0
1193  if(dim > 2)
1194  {
1195  const double power = (1.0/(dim-1));
1196  faceVol = pow(faceVol, power);
1197  }
1198 
1199  // evaluate | (u_l * n_l) + (u_r * n_r) | = | (u_l - u_r) * n_l |
1200  uLF.evaluate(faceQuadInner[l], jump);
1201  uNeighLf.evaluate(faceQuadOuter[l], neighRet);
1202 
1203  // get difference
1204  jump -= neighRet;
1205 
1206  const double weight = (enVol_nbVol) * faceQuadInner.weight(l);
1207 
1208  double error = weight * SQR(jump * unitNormal);
1209  error = std::abs( error );
1210 
1211  enError += error;
1212  nbError += error;
1213  }
1214 
1215  if( (nbError > 0.0)
1216 #if HAVE_MPI
1217  // only add neihgbor on interior entities
1218  && (interiorEntity)
1219 #endif
1220  )
1221  {
1222  adaptation.addToLocalIndicator( nb , nbError );
1223  }
1224  }
1225 
1226  // LU decomposition of matrix (matrix and b are overwritten)
1227  //
1228  // param[inout] a Matrix that LU decomposition is calculated for
1229  // param[in] b right hand side
1230  // param[out] solution solution of linear system
1231  template <class MatrixType, class VectorType>
1232  void luSolve(MatrixType& a,
1233  VectorType& x) const
1234  {
1235  typedef typename VectorType :: field_type ctype;
1236  enum { n = VectorType :: dimension };
1237 
1238  // make sure we got the right dimensions
1239  assert( a.N() == a.M() );
1240  assert( a.N() == n );
1241 
1242  // pivot storage
1243  int p[ n-1 ];
1244 
1245  for(int i=0; i<n-1; ++i)
1246  {
1247  // initialize
1248  p[i] = i;
1249 
1250  // Pivot search
1251  ctype max_abs = 0;
1252  for(int k=i; k<n; ++k)
1253  {
1254  if ( std::abs(a[k][i]) > max_abs )
1255  {
1256  max_abs = fabs(a[k][i]);
1257  p[i] = k;
1258  }
1259  }
1260 
1261  if( p[ i ] != i )
1262  {
1263  // toggle row i with row argmax=p[i]
1264  for(int j=0; j<n; ++j)
1265  {
1266  const ctype tmp = a[ p[i] ][j];
1267  a[ p[i] ][j] = a[i][j];
1268  a[i][j] = tmp;
1269  }
1270  }
1271 
1272  // elimination
1273  for(int k=i+1; k<n; ++k)
1274  {
1275  const ctype lambda = a[k][i] / a[i][i];
1276  a[k][i] = lambda;
1277  for(int j=i+1; j<n; ++j) a[k][j] -= lambda * a[i][j];
1278  }
1279  }
1280 
1281  // 1. x = Px_old, permutation with right hand side
1282  for(int i=0; i<n-1; ++i)
1283  {
1284  const ctype tmp = x[ i ];
1285  x[ i ] = x[ p[ i ] ];
1286  x[ p[ i ] ] = tmp;
1287  }
1288 
1289  // 1. Lx = x_old, forward loesen
1290  for(int i=0; i<n; ++i)
1291  {
1292  ctype dot = 0;
1293  for(int j=0; j<i; ++j) dot += a[i][j] * x[j];
1294  x[i] -= dot;
1295  }
1296 
1297  // 2. Ux = x_old, backward solve
1298  for(int i=n-1; i>=0; --i)
1299  {
1300  ctype dot = 0;
1301  for(int j=i+1; j<n; ++j) dot += a[i][j] * x[j];
1302  x[i] = (x[i] - dot) / a[i][i];
1303  }
1304  }
1305  };
1306 
1307  } // namespace Fem
1308 
1309 } // namespace Dune
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
Definition: arrays.hh:36
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