Loading [MathJax]/extensions/tex2jax.js
ug4
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
integrate.h
Go to the documentation of this file.
1/*
2 * Copyright (c) 2011-2015: G-CSC, Goethe University Frankfurt
3 * Author: Andreas Vogel
4 *
5 * This file is part of UG4.
6 *
7 * UG4 is free software: you can redistribute it and/or modify it under the
8 * terms of the GNU Lesser General Public License version 3 (as published by the
9 * Free Software Foundation) with the following additional attribution
10 * requirements (according to LGPL/GPL v3):
11 *
12 * (1) The following notice must be displayed in the Appropriate Legal Notices
13 * of covered and combined works: "Based on UG4 (www.ug4.org/license)".
14 *
15 * (2) The following notice must be displayed at a prominent place in the
16 * terminal output of covered works: "Based on UG4 (www.ug4.org/license)".
17 *
18 * (3) The following bibliography is recommended for citation and must be
19 * preserved in all covered files:
20 * "Reiter, S., Vogel, A., Heppner, I., Rupp, M., and Wittum, G. A massively
21 * parallel geometric multigrid solver on hierarchically distributed grids.
22 * Computing and visualization in science 16, 4 (2013), 151-164"
23 * "Vogel, A., Reiter, S., Rupp, M., Nägel, A., and Wittum, G. UG4 -- a novel
24 * flexible software system for simulating pde based models on high performance
25 * computers. Computing and visualization in science 16, 4 (2013), 165-179"
26 *
27 * This program is distributed in the hope that it will be useful,
28 * but WITHOUT ANY WARRANTY; without even the implied warranty of
29 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
30 * GNU Lesser General Public License for more details.
31 */
32
33#ifndef __H__UG__LIB_DISC__FUNCTION_SPACES__INTEGRATE__
34#define __H__UG__LIB_DISC__FUNCTION_SPACES__INTEGRATE__
35
36#include <cmath>
37
38#include <boost/function.hpp>
39
40#include "common/common.h"
41
43
46#include "lib_disc/domain_util.h" // for CollectCornerCoordinates
53
54#ifdef UG_FOR_LUA
56#endif
57
58namespace ug{
59
60
62// Object encapsulating (scalar) GridFunction related data
64template <typename TGridFunction>
66{
67public:
68 typedef typename TGridFunction::domain_type domain_type;
69
70 ScalarGridFunctionData(TGridFunction& gridFct, size_t cmp)
71 : m_gridFct(gridFct), m_fct(cmp),
72 m_id(m_gridFct.local_finite_element_id(m_fct))
73 {}
74
75 TGridFunction& grid_function()
76 {return m_gridFct;}
77
78 const TGridFunction& grid_function() const
79 {return m_gridFct;}
80
81 size_t fct()
82 {return m_fct;}
83
84 const LFEID &id() const
85 {return m_id;}
86
88 bool is_def_in_subset(int si) const
89 { return m_gridFct.is_def_in_subset(m_fct, si); }
90
93
95 ConstSmartPtr<domain_type> domain() const {return m_gridFct.domain();}
96
97 template <typename TElem>
98 size_t dof_indices(TElem* elem, std::vector<DoFIndex>& ind, bool bHang = false, bool bClear = true) const
99 {return m_gridFct.dof_indices(elem, m_fct, ind, bHang, bClear);}
100
101private:
103 TGridFunction& m_gridFct;
104
106 size_t m_fct;
107
110};
111
113// Integrand Interface
115
116
118
119template <typename TData, int TWorldDim>
121{
122 public:
124 static const int worldDim = TWorldDim;
125
127 typedef TData data_type;
128
130
141 virtual void values(TData vValue[],
142 const MathVector<worldDim> vGlobIP[],
143 GridObject* pElem,
144 const MathVector<worldDim> vCornerCoords[],
145 const MathVector<1> vLocIP[],
146 const MathMatrix<1, worldDim> vJT[],
147 const size_t numIP) = 0;
148 virtual void values(TData vValue[],
149 const MathVector<worldDim> vGlobIP[],
150 GridObject* pElem,
151 const MathVector<worldDim> vCornerCoords[],
152 const MathVector<2> vLocIP[],
153 const MathMatrix<2, worldDim> vJT[],
154 const size_t numIP) = 0;
155 virtual void values(TData vValue[],
156 const MathVector<worldDim> vGlobIP[],
157 GridObject* pElem,
158 const MathVector<worldDim> vCornerCoords[],
159 const MathVector<3> vLocIP[],
160 const MathMatrix<3, worldDim> vJT[],
161 const size_t numIP) = 0;
163
164 virtual ~IIntegrand() {}
165
166
168 virtual void set_subset(int si) {m_si = si;}
169
171 inline int subset() const {return m_si;}
172
173 protected:
175 int m_si;
176};
177
179template <typename TData, int TWorldDim, typename TImpl>
180class StdIntegrand : public IIntegrand<TData, TWorldDim>
181{
182 public:
184 static const int worldDim = TWorldDim;
185
187 typedef TData data_type;
188
191 virtual void values(TData vValue[],
192 const MathVector<worldDim> vGlobIP[],
193 GridObject* pElem,
194 const MathVector<worldDim> vCornerCoords[],
195 const MathVector<1> vLocIP[],
196 const MathMatrix<1, worldDim> vJT[],
197 const size_t numIP)
198 {
199 getImpl().template evaluate<1>(vValue,vGlobIP,pElem,vCornerCoords,vLocIP,vJT,numIP);
200 }
201 virtual void values(TData vValue[],
202 const MathVector<worldDim> vGlobIP[],
203 GridObject* pElem,
204 const MathVector<worldDim> vCornerCoords[],
205 const MathVector<2> vLocIP[],
206 const MathMatrix<2, worldDim> vJT[],
207 const size_t numIP)
208 {
209 getImpl().template evaluate<2>(vValue,vGlobIP,pElem,vCornerCoords,vLocIP,vJT,numIP);
210 }
211 virtual void values(TData vValue[],
212 const MathVector<worldDim> vGlobIP[],
213 GridObject* pElem,
214 const MathVector<worldDim> vCornerCoords[],
215 const MathVector<3> vLocIP[],
216 const MathMatrix<3, worldDim> vJT[],
217 const size_t numIP)
218 {
219 getImpl().template evaluate<3>(vValue,vGlobIP,pElem,vCornerCoords,vLocIP,vJT,numIP);
220 }
222
223 protected:
225 TImpl& getImpl() {return static_cast<TImpl&>(*this);}
226
228 const TImpl& getImpl() const {return static_cast<const TImpl&>(*this);}
229};
230
233// Generic Volume Integration Routine
236
238
259template <int WorldDim, int dim, typename TConstIterator>
260number Integrate(TConstIterator iterBegin,
261 TConstIterator iterEnd,
264 int quadOrder, std::string quadType,
267 *paaElemContribs = NULL
268 )
269{
270 PROFILE_FUNC();
271
272// reset the result
273 number integral = 0.0;
274
275// note: this iterator is for the base elements, e.g. Face and not
276// for the special type, e.g. Triangle, Quadrilateral
277 TConstIterator iter = iterBegin;
278
279// this is the base element type (e.g. Face). This is the type when the
280// iterators above are dereferenciated.
281 typedef typename domain_traits<dim>::grid_base_object grid_base_object;
282
283// get quad type
284 if(quadType.empty()) quadType = "best";
285 QuadType type = GetQuadratureType(quadType);
286
287// accessing without dereferencing a pointer first is simpler...
289 if(paaElemContribs)
290 aaElemContribs = *paaElemContribs;
291
292// We'll reuse containers to avoid reallocations
293 std::vector<MathVector<WorldDim> > vCorner;
294 std::vector<MathVector<WorldDim> > vGlobIP;
295 std::vector<MathMatrix<dim, WorldDim> > vJT;
296 std::vector<number> vValue;
297
298// iterate over all elements
299 for(; iter != iterEnd; ++iter)
300 {
301 // get element
302 grid_base_object* pElem = *iter;
303
304 // get reference object id (i.e. Triangle, Quadrilateral, Tetrahedron, ...)
305 ReferenceObjectID roid = (ReferenceObjectID) pElem->reference_object_id();
306
307 // get quadrature Rule for reference object id and order
308 try{
309 const QuadratureRule<dim>& rQuadRule
310 = QuadratureRuleProvider<dim>::get(roid, quadOrder, type);
311
312 // get reference element mapping by reference object id
314 = ReferenceMappingProvider::get<dim, WorldDim>(roid);
315
316 // number of integration points
317 const size_t numIP = rQuadRule.size();
318
319 // get all corner coordinates
320 CollectCornerCoordinates(vCorner, *pElem, aaPos, true);
321
322 // update the reference mapping for the corners
323 mapping.update(&vCorner[0]);
324
325 // compute global integration points
326 vGlobIP.resize(numIP);
327 mapping.local_to_global(&(vGlobIP[0]), rQuadRule.points(), numIP);
328
329 // compute transformation matrices
330 vJT.resize(numIP);
331 mapping.jacobian_transposed(&(vJT[0]), rQuadRule.points(), numIP);
332
333 // compute integrand values at integration points
334 vValue.resize(numIP);
335 try
336 {
337 integrand.values(&(vValue[0]), &(vGlobIP[0]),
338 pElem, &vCorner[0], rQuadRule.points(),
339 &(vJT[0]),
340 numIP);
341 }
342 UG_CATCH_THROW("Unable to compute values of integrand at integration point.");
343
344 // reset contribution of this element
345 number intValElem = 0;
346
347 // loop integration points
348 for(size_t ip = 0; ip < numIP; ++ip)
349 {
350 // get quadrature weight
351 const number weightIP = rQuadRule.weight(ip);
352
353 // get determinate of mapping
354 const number det = SqrtGramDeterminant(vJT[ip]);
355
356 // add contribution of integration point
357 intValElem += vValue[ip] * weightIP * det;
358 }
359
360 // add to global sum
361 integral += intValElem;
362 if(aaElemContribs.valid())
363 aaElemContribs[pElem] = intValElem;
364
365 }UG_CATCH_THROW("SumValuesOnElems failed.");
366 } // end elem
367
368// return the summed integral contributions of all elements
369 return integral;
370}
371
372template <typename TGridFunction, int dim>
374 TGridFunction& spGridFct,
375 int si, int quadOrder, std::string quadType)
376{
377// integrate elements of subset
378 typedef typename TGridFunction::template dim_traits<dim>::grid_base_object grid_base_object;
379 typedef typename TGridFunction::template dim_traits<dim>::const_iterator const_iterator;
380
381 spIntegrand.set_subset(si);
382
383 return Integrate<TGridFunction::dim,dim,const_iterator>
384 (spGridFct.template begin<grid_base_object>(si),
385 spGridFct.template end<grid_base_object>(si),
386 spGridFct.domain()->position_accessor(),
387 spIntegrand,
388 quadOrder, quadType);
389}
390
391
392template <typename TGridFunction>
394 TGridFunction& spGridFct,
395 const char* subsets, int quadOrder,
396 std::string quadType = std::string())
397{
398// world dimensions
399 static const int dim = TGridFunction::dim;
400
401// read subsets
402 SubsetGroup ssGrp(spGridFct.domain()->subset_handler());
403 if(subsets != NULL)
404 {
405 ssGrp.add(TokenizeString(subsets));
406 UG_COND_THROW(!SameDimensionsInAllSubsets(ssGrp), "IntegrateSubsets: Subsets '"<<subsets<<"' do not have same dimension."
407 "Cannot integrate on subsets of different dimensions.");
408 UG_LOG("IntegrateSubsets for subsets="<<subsets<<"\n");
409 }
410 else
411 {
412 // add all subsets and remove lower dim subsets afterwards
413 ssGrp.add_all();
415 }
416
417// reset value
418 number value = 0;
419
420// loop subsets
421 for(size_t i = 0; i < ssGrp.size(); ++i)
422 {
423 // get subset index
424 const int si = ssGrp[i];
425
426 // check dimension
427 UG_COND_THROW(ssGrp.dim(i) > dim, "IntegrateSubsets: Dimension of subset is "<<ssGrp.dim(i)<<", but "
428 " world dimension is "<<dim<<". Cannot integrate this.");
429
430 // integrate elements of subset
431 try{
432 switch(ssGrp.dim(i))
433 {
434 case DIM_SUBSET_EMPTY_GRID: break;
435 case 1: value += IntegrateSubset<TGridFunction, 1>(spIntegrand, spGridFct, si, quadOrder, quadType); break;
436 case 2: value += IntegrateSubset<TGridFunction, 2>(spIntegrand, spGridFct, si, quadOrder, quadType); break;
437 case 3: value += IntegrateSubset<TGridFunction, 3>(spIntegrand, spGridFct, si, quadOrder, quadType); break;
438 default: UG_THROW("IntegrateSubsets: Dimension "<<ssGrp.dim(i)<<" not supported. "
439 " World dimension is "<<dim<<".");
440 }
441 }
442 UG_CATCH_THROW("IntegrateSubsets: Integration failed on subset "<<si);
443 }
444
445#ifdef UG_PARALLEL
446 // sum over processes
447 if(pcl::NumProcs() > 1)
448 {
450 number local = value;
451 com.allreduce(&local, &value, 1, PCL_DT_DOUBLE, PCL_RO_SUM);
452 }
453#endif
454
455// return the result
456 return value;
457}
458
459
461// UserData Integrand
463
465template <typename TData, typename TGridFunction>
467 : public StdIntegrand<TData, TGridFunction::dim, UserDataIntegrand<TData, TGridFunction> >
468{
469 public:
470 // world dimension of grid function
471 static const int worldDim = TGridFunction::dim;
472
473 // data type
474 typedef TData data_type;
475
476 private:
477 // data to integrate
479
480 // grid function
481 TGridFunction* m_spGridFct;
482
483 // time
485
486 public:
489 TGridFunction* spGridFct,
490 number time)
491 : m_spData(spData), m_spGridFct(spGridFct), m_time(time)
492 {
493 m_spData->set_function_pattern(spGridFct->function_pattern());
494 };
495
498 number time)
499 : m_spData(spData), m_spGridFct(NULL), m_time(time)
500 {
501 UG_COND_THROW(m_spData->requires_grid_fct(),
502 "UserDataIntegrand: Missing GridFunction, but data requires grid function.");
503 };
504
506 template <int elemDim>
507 void evaluate(TData vValue[],
508 const MathVector<worldDim> vGlobIP[],
509 GridObject* pElem,
510 const MathVector<worldDim> vCornerCoords[],
511 const MathVector<elemDim> vLocIP[],
513 const size_t numIP)
514 {
515 // get local solution if needed
516 if(m_spData->requires_grid_fct())
517 {
518 // create storage
519 LocalIndices ind;
520 LocalVector u;
521
522 // get global indices
523 m_spGridFct->indices(pElem, ind);
524
525 // adapt local algebra
526 u.resize(ind);
527
528 // read local values of u
530
531 // compute data
532 try{
533 (*m_spData)(vValue, vGlobIP, m_time, this->m_si, pElem,
534 vCornerCoords, vLocIP, numIP, &u, &vJT[0]);
535 }
536 UG_CATCH_THROW("UserDataIntegrand: Cannot evaluate data.");
537 }
538 else
539 {
540 // compute data
541 try{
542 (*m_spData)(vValue, vGlobIP, m_time, this->m_si, numIP);
543 }
544 UG_CATCH_THROW("UserDataIntegrand: Cannot evaluate data.");
545 }
546
547 }
548};
549
550
551
552
554template <typename TData, typename TGridFunction>
556 : public StdIntegrand<number, TGridFunction::dim, UserDataIntegrandSq<TData, TGridFunction> >
557{
558 public:
559 // world dimension of grid function
560 static const int worldDim = TGridFunction::dim;
561
562 // data type
563 typedef TData data_type;
564
565 private:
566 // data to integrate
568
569 // grid function
570 const TGridFunction* m_pGridFct;
571
572 // time
574
575 public:
578 const TGridFunction* pGridFct,
579 number time)
580 : m_spData(spData), m_pGridFct(pGridFct), m_time(time)
581 {
582 m_spData->set_function_pattern(pGridFct->function_pattern());
583 };
584
587 number time)
588 : m_spData(spData), m_pGridFct(NULL), m_time(time)
589 {
590 if(m_spData->requires_grid_fct())
591 UG_THROW("UserDataIntegrand: Missing GridFunction, but "
592 " data requires grid function.")
593 };
594
595
597 template <int elemDim>
598 void evaluate(number vValue[],
599 const MathVector<worldDim> vGlobIP[],
600 GridObject* pElem,
601 const MathVector<worldDim> vCornerCoords[],
602 const MathVector<elemDim> vLocIP[],
604 const size_t numIP)
605 {
606
607 std::vector<TData> tmpValues(numIP);
608
609 // get local solution if needed
610 if(m_spData->requires_grid_fct())
611 {
612 UG_ASSERT(m_pGridFct!=NULL, "Error: Requires valid pointer!")
613 // create storage
614 LocalIndices ind;
615 LocalVector u;
616
617 // get global indices
618 m_pGridFct->indices(pElem, ind);
619
620 // adapt local algebra
621 u.resize(ind);
622
623 // read local values of u
625
626 // compute data
627 try{
628 (*m_spData)(&tmpValues.front(), vGlobIP, m_time, this->m_si, pElem,
629 vCornerCoords, vLocIP, numIP, &u, &vJT[0]);
630 }
631 UG_CATCH_THROW("UserDataIntegrand: Cannot evaluate data.");
632 }
633 else
634 {
635 // compute data
636 try{
637 (*m_spData)(&tmpValues.front(), vGlobIP, m_time, this->m_si, numIP);
638 }
639 UG_CATCH_THROW("UserDataIntegrand: Cannot evaluate data.");
640 }
641
642 for (size_t i=0; i<numIP; ++i)
643 {
644 vValue[i]=inner_prod(tmpValues[i], tmpValues[i]);
645 }
646
647 }
648 protected:
649
650 number inner_prod(const number &d1, const number &d2)
651 {return d1*d2;}
652
654 {return VecDot(d1, d2);}
655
656 template<typename T>
657 number inner_prod(const T &d1, const T &d2)
658 { UG_ASSERT(0, "NOT IMPLEMENTED"); return 0.0;}
659};
660
661
662
663
665template <typename TData, typename TGridFunction>
667 : public StdIntegrand<number, TGridFunction::dim, UserDataDistIntegrandSq<TData, TGridFunction> >
668{
669 public:
671 static const int worldDim = TGridFunction::dim;
672 // typedef typename L2Integrand<TGridFunction>::weight_type weight_type;
673
674 protected:
676 const int m_fineTopLevel;
677
680
683
684
685
687 // ConstSmartPtr<weight_type> m_spWeight;
688
689 double m_time;
690 public:
691
693 UserDataDistIntegrandSq(SmartPtr<UserData<TData, worldDim> > spData, TGridFunction& fineGridFct, size_t fineCmp,
694 TGridFunction& coarseGridFct, size_t coarseCmp)
695 : m_fineData(fineGridFct, fineCmp), m_fineTopLevel(fineGridFct.dof_distribution()->grid_level().level()),
696 m_coarseData(coarseGridFct, coarseCmp), m_coarseTopLevel(coarseGridFct.dof_distribution()->grid_level().level()),
697 m_spMG(m_fineData.domain()->grid()), m_spData(spData), m_time(0.0) /*, m_spWeight(make_sp(new ConstUserNumber<TGridFunction::dim>(1.0)))*/
698 {
700 UG_THROW("UserDataDistIntegrandSq: fine and top level inverted.");
701
702 if(m_fineData.domain().get() != m_coarseData.domain().get())
703 UG_THROW("UserDataDistIntegrandSq: grid functions defined on different domains.");
704 };
705
707
709 virtual void set_subset(int si)
710 {
711 if(!m_fineData.is_def_in_subset(si))
712 UG_THROW("UserDataDistIntegrandSq: Grid function component"
713 <<m_fineData.fct()<<" not defined on subset "<<si);
714 if(!m_coarseData.is_def_in_subset(si))
715 UG_THROW("UserDataDistIntegrandSq: Grid function component"
716 <<m_coarseData.fct()<<" not defined on subset "<<si);
718 }
719
721 template <int elemDim>
722 void evaluate(number vValue[],
723 const MathVector<worldDim> vGlobIP[],
724 GridObject* pFineElem,
725 const MathVector<worldDim> vCornerCoords[],
726 const MathVector<elemDim> vFineLocIP[],
728 const size_t numIP)
729 {
730
731
732 // must return 0.0, if m_spData is independent of grid functions
733 if(!m_spData->requires_grid_fct()) {
734 for (size_t i=0; i<numIP; ++i) { vValue[i]=0.0; }
735 return;
736 }
737
738 // get coarse element
739 GridObject* pCoarseElem = pFineElem;
741 int parentLevel = m_spMG->get_level(pCoarseElem);
742 while(parentLevel > m_coarseTopLevel){
743 pCoarseElem = m_spMG->get_parent(pCoarseElem);
744 parentLevel = m_spMG->get_level(pCoarseElem);
745 }
746 }
747
748 // get corner coordinates
749 typedef typename TGridFunction::template dim_traits<elemDim>::grid_base_object TElem;
750 std::vector<MathVector<worldDim> > vCornerCoarse;
751 CollectCornerCoordinates(vCornerCoarse, *static_cast<TElem*>(pCoarseElem), *m_coarseData.domain());
752
753 // get Reference Mapping
754 const ReferenceObjectID coarseROID = pCoarseElem->reference_object_id();
756 = ReferenceMappingProvider::get<elemDim, worldDim>(coarseROID, vCornerCoarse);
757
758 std::vector<MathVector<elemDim> > vCoarseLocIP;
759 vCoarseLocIP.resize(numIP);
760 for(size_t ip = 0; ip < vCoarseLocIP.size(); ++ip) VecSet(vCoarseLocIP[ip], 0.0);
761 map.global_to_local(&vCoarseLocIP[0], vGlobIP, numIP);
762
763 // element weights
764 /* typedef typename weight_type::data_type ipdata_type;
765 std::vector<ipdata_type> fineElemWeights(numIP, 1.0);
766 UG_ASSERT(m_spWeight.valid(), "L2DistIntegrand::evaluate requires valid weights! ");
767 (*m_spWeight)(&fineElemWeights[0], vFineGlobIP, 0.0, IIntegrand<number, worldDim>::subset(), numIP);*/
768
769
770 // get trial space
771 /* const LocalShapeFunctionSet<elemDim>& rFineLSFS =
772 LocalFiniteElementProvider::get<elemDim>(fineROID, m_fineData.id());
773 const LocalShapeFunctionSet<elemDim>& rCoarseLSFS =
774 LocalFiniteElementProvider::get<elemDim>(coarseROID, m_coarseData.id());*/
775
776 // get multiindices of element
777 /*std::vector<DoFIndex> vFineMI, vCoarseMI;
778 m_fineData.dof_indices(pFineElem, vFineMI);
779 m_coarseData.dof_indices(pCoarseElem, vCoarseMI);*/
780
781 TData fineValues[numIP];
782 TData coarseValues[numIP];
783
784
785
786 // compute coarse data
787 try{
788 LocalVector uCoarse;
789 LocalIndices indCoarse;
790
791 m_coarseData.grid_function().indices(pCoarseElem, indCoarse);
792 uCoarse.resize(indCoarse);
793 GetLocalVector(uCoarse, m_coarseData.grid_function());
794
795 (*m_spData)(coarseValues, vGlobIP, m_time, this->m_si, pCoarseElem,
796 &vCornerCoarse[0], &vCoarseLocIP[0], numIP, &uCoarse, &vJT[0]);
797 } UG_CATCH_THROW("UserDataDistIntegrandSq: Cannot evaluate coarse data.");
798
799
800 // compute fine data
801 try{
802 LocalVector uFine;
803 LocalIndices indFine;
804
805 m_fineData.grid_function().indices(pFineElem, indFine);
806 uFine.resize(indFine);
807 GetLocalVector(uFine, m_fineData.grid_function());
808
809 (*m_spData)(fineValues, vGlobIP, m_time, this->m_si, pFineElem,
810 vCornerCoords, vFineLocIP, numIP, &uFine, &vJT[0]);
811 } UG_CATCH_THROW("UserDataDistIntegrandSq: Cannot evaluate fine data.");
812
813 // loop all integration points
814 for(size_t ip = 0; ip < numIP; ++ip)
815 {
816 vValue[ip] = inner_dist2(fineValues[ip], coarseValues[ip]);
817 }
818
819 }
820
821
822 protected:
823
824 number inner_prod(const number &d1, const number &d2)
825 {return d1*d2;}
826
828 {return VecDot(d1, d2);}
829
830 template<typename T>
831 number inner_prod(const T &d1, const T &d2)
832 { UG_ASSERT(0, "NOT IMPLEMENTED"); return 0.0;}
833
834
835
836 number inner_dist2(const number &v1, const number &v2)
837 { return (v1-v2)*(v1-v2); }
838
840 { return VecDistanceSq(v1, v2); }
841
842 template<typename T>
843 number inner_dist2(const T &d1, const T &d2)
844 { UG_ASSERT(0, "NOT IMPLEMENTED"); return 0.0;}
845
846};
847
848
851// Volume Integrand implementations
854
855template <typename TGridFunction>
857 TGridFunction& spGridFct,
858 const char* subsets, number time,
859 int quadOrder, std::string quadType)
860{
861 UserDataIntegrand<number, TGridFunction> spIntegrand(spData, &spGridFct, time);
862 return IntegrateSubsets(spIntegrand, spGridFct, subsets, quadOrder, quadType);
863}
864
865template <typename TGridFunction>
867 const char* subsets, number time, int quadOrder, std::string quadType)
868{ return Integral(spData, *spGridFct, subsets, time, quadOrder, quadType); }
869
870template <typename TGridFunction>
871number Integral(SmartPtr<UserData<number, TGridFunction::dim> > spData, SmartPtr<TGridFunction> spGridFct,const char* subsets,number time, int order)
872{return Integral(spData, spGridFct, subsets, time, order, "best");}
873
874template <typename TGridFunction>
876{return Integral(spData, spGridFct, subsets, time, 1, "best");}
877
878template <typename TGridFunction>
880{return Integral(spData, spGridFct, NULL, time, 1, "best");}
881
882template <typename TGridFunction>
884{return Integral(spData, spGridFct, subsets, 0.0, 1, "best");}
885
886template <typename TGridFunction>
888{return Integral(spData, spGridFct, NULL, 0.0, 1, "best");}
889
891// const data
893
894template <typename TGridFunction>
896 const char* subsets,
897 number time, int quadOrder)
898{
899 static const int dim = TGridFunction::dim;
902 return Integral(sp, spGridFct, subsets, time, quadOrder);
903}
904
905template <typename TGridFunction>
906number Integral(number val, SmartPtr<TGridFunction> spGridFct,const char* subsets,number time)
907{return Integral(val, spGridFct, subsets, time, 1);}
908
909template <typename TGridFunction>
911{return Integral(val, spGridFct, NULL, time, 1);}
912
913template <typename TGridFunction>
914number Integral(number val, SmartPtr<TGridFunction> spGridFct,const char* subsets)
915{return Integral(val, spGridFct, subsets, 0.0, 1);}
916
917template <typename TGridFunction>
919{return Integral(val, spGridFct, NULL, 0.0, 1);}
920
922// lua data
924
925#ifdef UG_FOR_LUA
926template <typename TGridFunction>
927number Integral(const char* luaFct,
928 SmartPtr<TGridFunction> spGridFct,
929 const char* subsets, number time,
930 int quadOrder, std::string quadType)
931{
932 static const int dim = TGridFunction::dim;
935 return Integral(sp, *spGridFct, subsets, time, quadOrder, quadType);
936}
937
938template <typename TGridFunction>
939number Integral(const char* luaFct, SmartPtr<TGridFunction> spGridFct,const char* subsets,number time, int quadOrder)
940{return Integral(luaFct, spGridFct, subsets, time, quadOrder, "best");}
941
942template <typename TGridFunction>
943number Integral(const char* luaFct, SmartPtr<TGridFunction> spGridFct,const char* subsets,number time)
944{return Integral(luaFct, spGridFct, subsets, time, 1, "best");}
945
946template <typename TGridFunction>
947number Integral(const char* luaFct, SmartPtr<TGridFunction> spGridFct,number time)
948{return Integral(luaFct, spGridFct, NULL, time, 1, "best");}
949
950template <typename TGridFunction>
951number Integral(const char* luaFct, SmartPtr<TGridFunction> spGridFct,const char* subsets)
952{return Integral(luaFct, spGridFct, subsets, 0.0, 1, "best");}
953
954template <typename TGridFunction>
955number Integral(const char* luaFct, SmartPtr<TGridFunction> spGridFct)
956{return Integral(luaFct, spGridFct, NULL, 0.0, 1, "best");}
957
958#endif
959
960
961
962template <typename TGridFunction>
964 : public StdIntegrand<number, TGridFunction::dim, MaximumDistIntegrand<TGridFunction> >
965{
966 public:
968 static const int worldDim = TGridFunction::dim;
969
970 protected:
972
973 double m_max, m_min;
974
975 public:
977 MaximumDistIntegrand(TGridFunction& gridFct, size_t cmp)
978
979 : m_scalarData(gridFct, cmp),
980 m_max(-std::numeric_limits<float>::infinity()),
981 m_min(std::numeric_limits<float>::infinity())
982 {};
983
985
987 virtual void set_subset(int si)
988 {
989 UG_COND_THROW(!m_scalarData.is_def_in_subset(si), "L2ErrorIntegrand: Grid function component" <<
990 m_scalarData.fct()<<" not defined on subset "<<si);
992 }
993
994 void reset()
995 {
996 m_min = std::numeric_limits<double>::infinity();
997 m_max = -std::numeric_limits<double>::infinity();
998 }
999
1000 double min() { return m_min; }
1001 double max() { return m_max; }
1002
1004 template <int elemDim>
1005 void evaluate(number vValue[],
1006 const MathVector<worldDim> vGlobIP[],
1007 GridObject* pElem,
1008 const MathVector<worldDim> vCornerCoords[],
1009 const MathVector<elemDim> vLocIP[],
1011 const size_t numIP)
1012 {
1013 // get reference object id (i.e. Triangle, Quadrilateral, Tetrahedron, ...)
1014 const ReferenceObjectID roid = pElem->reference_object_id();
1015
1016 try{
1017 // get trial space
1018 const LocalShapeFunctionSet<elemDim>& rTrialSpace =
1019 LocalFiniteElementProvider::get<elemDim>(roid, m_scalarData.id());
1020
1021 // number of dofs on element
1022 const size_t num_sh = rTrialSpace.num_sh();
1023
1024 // get multiindices of element
1025 std::vector<DoFIndex> ind; // aux. index array
1026 m_scalarData.dof_indices(pElem, ind);
1027
1028 // check multi indices
1029 UG_COND_THROW(ind.size() != num_sh, "L2ErrorIntegrand::evaluate: Wrong number of multi indices.");
1030
1031 // evaluate function in corners.
1032 for(size_t sh = 0; sh < num_sh; ++sh)
1033 {
1034 // get value at shape point (e.g. corner for P1 fct)
1035 const number val = DoFRef(m_scalarData.grid_function(), ind[sh]);
1036
1037 m_max = std::max(m_max, val);
1038 m_min = std::min(m_min, val);
1039 }
1040
1041
1042 } UG_CATCH_THROW("L2ErrorIntegrand::evaluate: trial space missing.");
1043 }
1044};
1045
1046
1047template <typename TGridFunction>
1048number GetMinimum(TGridFunction &gridFct, const char* cmp, const char* subsets)
1049{
1050 // get function id of name & check that function exists
1051 const size_t fct = gridFct.fct_id_by_name(cmp);
1052 UG_COND_THROW(fct >= gridFct.num_fct(),
1053 "L2Error: Function space does not contain a function with name " << cmp << ".");
1054
1055 // Find minimum.
1056 MaximumDistIntegrand<TGridFunction> spIntegrand(gridFct, fct);
1057 IntegrateSubsets(spIntegrand, gridFct, subsets, 2);
1058 return spIntegrand.min();
1059}
1060
1061template <typename TGridFunction>
1062number Minimum(SmartPtr<TGridFunction> spGridFct, const char* cmp, const char* subsets)
1063{ return GetMinimum(*spGridFct, cmp, subsets); }
1064
1065
1067// L2 Error Integrand
1069
1070template <typename TGridFunction>
1072 : public StdIntegrand<number, TGridFunction::dim, L2ErrorIntegrand<TGridFunction> >
1073{
1074 public:
1076 static const int worldDim = TGridFunction::dim;
1077
1078 private:
1080
1083
1086
1087 public:
1090 TGridFunction& gridFct, size_t cmp,
1091 number time)
1092 : m_scalarData(gridFct, cmp),
1093 m_spExactSolution(spExactSol), m_time(time)
1094 {};
1095
1096 virtual ~L2ErrorIntegrand() {};
1097
1099 virtual void set_subset(int si)
1100 {
1101 if(!m_scalarData.is_def_in_subset(si))
1102 UG_THROW("L2ErrorIntegrand: Grid function component"
1103 <<m_scalarData.fct()<<" not defined on subset "<<si);
1105 }
1106
1108 template <int elemDim>
1109 void evaluate(number vValue[],
1110 const MathVector<worldDim> vGlobIP[],
1111 GridObject* pElem,
1112 const MathVector<worldDim> vCornerCoords[],
1113 const MathVector<elemDim> vLocIP[],
1115 const size_t numIP)
1116 {
1117 // get reference object id (i.e. Triangle, Quadrilateral, Tetrahedron, ...)
1118 const ReferenceObjectID roid = pElem->reference_object_id();
1119
1120 try{
1121 // get trial space
1122 const LocalShapeFunctionSet<elemDim>& rTrialSpace =
1123 LocalFiniteElementProvider::get<elemDim>(roid, m_scalarData.id());
1124
1125 // number of dofs on element
1126 const size_t num_sh = rTrialSpace.num_sh();
1127
1128 // get multiindices of element
1129 std::vector<DoFIndex> ind; // aux. index array
1130 m_scalarData.dof_indices(pElem, ind);
1131
1132 // check multi indices
1133 if(ind.size() != num_sh)
1134 UG_THROW("L2ErrorIntegrand::evaluate: Wrong number of"
1135 " multi indices.");
1136
1137 // loop all integration points
1138 for(size_t ip = 0; ip < numIP; ++ip)
1139 {
1140 // compute exact solution at integration point
1141 number exactSolIP;
1142 (*m_spExactSolution)(exactSolIP, vGlobIP[ip], m_time, this->subset());
1143
1144 // compute approximated solution at integration point
1145 number approxSolIP = 0.0;
1146 for(size_t sh = 0; sh < num_sh; ++sh)
1147 {
1148 // get value at shape point (e.g. corner for P1 fct)
1149 const number valSH = DoFRef(m_scalarData.grid_function(), ind[sh]);
1150
1151 // add shape fct at ip * value at shape
1152 approxSolIP += valSH * rTrialSpace.shape(sh, vLocIP[ip]);
1153 }
1154
1155 // get squared of difference
1156 vValue[ip] = (exactSolIP - approxSolIP);
1157 vValue[ip] *= vValue[ip];
1158 }
1159
1160 }
1161 UG_CATCH_THROW("L2ErrorIntegrand::evaluate: trial space missing.");
1162 }
1163};
1164
1165
1166
1168// L2 Error Integrand
1170
1172
1184template <typename TGridFunction>
1186 TGridFunction& gridFct, const char* cmp,
1187 number time, int quadOrder, const char* subsets)
1188{
1189// get function id of name
1190 const size_t fct = gridFct.fct_id_by_name(cmp);
1191
1192// check that function exists
1193 UG_COND_THROW(fct >= gridFct.num_fct(),
1194 "L2Error: Function space does not contain a function with name " << cmp << ".");
1195
1196 L2ErrorIntegrand<TGridFunction> spIntegrand(spExactSol, gridFct, fct, time);
1197 return sqrt(IntegrateSubsets(spIntegrand, gridFct, subsets, quadOrder));
1198}
1199
1200template <typename TGridFunction>
1202 SmartPtr<TGridFunction> spGridFct, const char* cmp,
1203 number time, int quadOrder, const char* subsets)
1204{ return L2Error(spExactSol, *spGridFct, cmp, time, quadOrder, subsets); }
1205
1206template <typename TGridFunction>
1208 SmartPtr<TGridFunction> spGridFct, const char* cmp,
1209 number time, int quadOrder)
1210{ return L2Error(spExactSol, *spGridFct, cmp, time, quadOrder, NULL); }
1211
1212
1213
1214
1215#ifdef UG_FOR_LUA
1216template <typename TGridFunction>
1217number L2Error(const char* ExactSol,
1218 SmartPtr<TGridFunction> spGridFct, const char* cmp,
1219 number time, int quadOrder, const char* subsets)
1220{
1222 = make_sp(new LuaUserData<number, TGridFunction::domain_type::dim>(ExactSol));
1223 return L2Error(spExactSol, spGridFct, cmp, time, quadOrder, subsets);
1224}
1225
1226template <typename TGridFunction>
1227number L2Error(const char* ExactSol,
1228 SmartPtr<TGridFunction> spGridFct, const char* cmp,
1229 number time, int quadOrder)
1230{
1231 return L2Error(ExactSol, spGridFct, cmp, time, quadOrder, NULL);
1232}
1233#endif //UG_FOR_LUA
1234
1235
1236
1237#ifdef UG_DEBUG
1238// True, if max norm is less than SMALL.
1239template <int worldDim, int elemDim>
1240void CheckGeneralizedInverse(const MathMatrix<elemDim, worldDim> &JT, const MathMatrix<worldDim, elemDim> &JTInv)
1241{
1242 if (worldDim==elemDim) return;
1243
1244 MathMatrix<worldDim, elemDim> myTmp;
1245 GeneralizedInverse(myTmp, JT); //old algorithms with Inverse, but Inverse is only defined for NxN Matrix
1246 MatSubtract(myTmp, JTInv, myTmp);
1247
1248 UG_ASSERT(MatMaxNorm(myTmp)<SMALL,
1249 "The inverse matrix is not identity to the map jacobian transposed inverse.");
1250
1251}
1252#endif
1253
1254
1256// H1 Error Integrand
1258
1259template <typename TGridFunction>
1261 : public StdIntegrand<number, TGridFunction::dim, H1ErrorIntegrand<TGridFunction> >
1262{
1263 public:
1265 static const int worldDim = TGridFunction::dim;
1266
1267 private:
1270
1273
1276
1279
1280 public:
1284 TGridFunction& gridFct, size_t cmp,
1285 number time)
1286 : m_scalarData (gridFct, cmp),
1287 m_spExactSolution(spExactSol),
1288 m_spExactGrad(spExactGrad),
1289 m_time(time)
1290 {}
1291
1293 virtual void set_subset(int si)
1294 {
1295 if(!m_scalarData.is_def_in_subset(si))
1296 UG_THROW("H1Error: Grid function component"
1297 <<m_scalarData.fct()<<" not defined on subset "<<si);
1299 }
1300
1302 template <int elemDim>
1303 void evaluate(number vValue[],
1304 const MathVector<worldDim> vGlobIP[],
1305 GridObject* pElem,
1306 const MathVector<worldDim> vCornerCoords[],
1307 const MathVector<elemDim> vLocIP[],
1309 const size_t numIP)
1310 {
1311 // get reference object id (i.e. Triangle, Quadrilateral, Tetrahedron, ...)
1312 const ReferenceObjectID roid = pElem->reference_object_id();
1313
1315 = ReferenceMappingProvider::get<elemDim, worldDim>(roid, vCornerCoords);
1316
1317 try{
1318 // get trial space
1319 const LocalShapeFunctionSet<elemDim>& rTrialSpace =
1320 LocalFiniteElementProvider::get<elemDim>(roid, m_scalarData.id());
1321
1322 // number of dofs on element
1323 const size_t num_sh = rTrialSpace.num_sh();
1324
1325 // get multiindices of element
1326 std::vector<DoFIndex> ind; // aux. index array
1327 m_scalarData.dof_indices(pElem, ind);
1328
1329 // check multi indices
1330 if(ind.size() != num_sh)
1331 UG_THROW("H1ErrorIntegrand::evaluate: Wrong number of"
1332 " multi indices.");
1333
1334 // loop all integration points
1335 std::vector<MathVector<elemDim> > vLocGradient(num_sh);
1336 for(size_t ip = 0; ip < numIP; ++ip)
1337 {
1338 // compute exact solution at integration point
1339 number exactSolIP;
1340 (*m_spExactSolution)(exactSolIP, vGlobIP[ip], m_time, this->subset());
1341
1342 // compute exact gradient at integration point
1343 MathVector<worldDim> exactGradIP;
1344 (*m_spExactGrad)(exactGradIP, vGlobIP[ip], m_time, this->subset());
1345
1346 // compute shape gradients at ip
1347 rTrialSpace.grads(&vLocGradient[0], vLocIP[ip]);
1348
1349 // compute approximated solution at integration point
1350 number approxSolIP = 0.0;
1351 MathVector<elemDim> locTmp; VecSet(locTmp, 0.0);
1352 for(size_t sh = 0; sh < num_sh; ++sh)
1353 {
1354 // get value at shape point (e.g. corner for P1 fct)
1355 const number valSH = DoFRef(m_scalarData.grid_function(), ind[sh]);
1356
1357 // add shape fct at ip * value at shape
1358 approxSolIP += valSH * rTrialSpace.shape(sh, vLocIP[ip]);
1359
1360 // add gradient at ip
1361 VecScaleAppend(locTmp, valSH, vLocGradient[sh]);
1362 }
1363
1364 // compute global gradient
1365 MathVector<worldDim> approxGradIP;
1367 map.jacobian_transposed_inverse(JTInv, vLocIP[ip]); //Inverse(JTInv, vJT[ip]);
1368
1369#ifdef UG_DEBUG
1370 CheckGeneralizedInverse<worldDim,elemDim>(vJT[ip], JTInv);
1371#endif //UG_DEBUG
1372
1373 MatVecMult(approxGradIP, JTInv, locTmp);
1374
1375 // get squared of difference
1376 vValue[ip] = (exactSolIP - approxSolIP) * (exactSolIP - approxSolIP);
1377 vValue[ip] += VecDistanceSq(approxGradIP, exactGradIP);
1378 }
1379
1380 }
1381 UG_CATCH_THROW("H1ErrorIntegrand::evaluate: trial space missing.");
1382 }
1383};
1384
1386
1399template <typename TGridFunction>
1401 SmartPtr<UserData<MathVector<TGridFunction::dim>, TGridFunction::dim> > spExactGrad,
1402 SmartPtr<TGridFunction> spGridFct, const char* cmp,
1403 number time, int quadOrder, const char* subsets)
1404{
1405// get function id of name
1406 const size_t fct = spGridFct->fct_id_by_name(cmp);
1407
1408// check that function exists
1409 if(fct >= spGridFct->num_fct())
1410 UG_THROW("H1Error: Function space does not contain"
1411 " a function with name " << cmp << ".");
1412
1413 H1ErrorIntegrand<TGridFunction> spIntegrand(spExactSol, spExactGrad, *spGridFct, fct, time);
1414 return sqrt(IntegrateSubsets(spIntegrand, *spGridFct, subsets, quadOrder));
1415}
1416
1417template <typename TGridFunction>
1419 SmartPtr<UserData<MathVector<TGridFunction::dim>, TGridFunction::dim> > spExactGrad,
1420 SmartPtr<TGridFunction> spGridFct, const char* cmp,
1421 number time, int quadOrder)
1422{
1423 return H1Error(spExactSol, spExactGrad, spGridFct, cmp, time, quadOrder, NULL);
1424}
1425
1426#ifdef UG_FOR_LUA
1427template <typename TGridFunction>
1428number H1Error(const char* ExactSol, const char* ExactGrad,
1429 SmartPtr<TGridFunction> spGridFct, const char* cmp,
1430 number time, int quadOrder, const char* subsets)
1431{
1432 static const int dim = TGridFunction::domain_type::dim;
1434 = make_sp(new LuaUserData<number, dim>(ExactSol));
1435 SmartPtr<UserData<MathVector<dim>, dim> > spExactGrad
1436 = make_sp(new LuaUserData<MathVector<dim>, dim>(ExactGrad));
1437 return H1Error(spExactSol, spExactGrad, spGridFct, cmp, time, quadOrder, subsets);
1438}
1439
1440template <typename TGridFunction>
1441number H1Error(const char* ExactSol, const char* ExactGrad,
1442 SmartPtr<TGridFunction> spGridFct, const char* cmp,
1443 number time, int quadOrder)
1444{
1445 return H1Error(ExactSol, ExactGrad, spGridFct, cmp, time, quadOrder, NULL);
1446}
1447#endif
1448
1449
1450
1452
1465template <typename TDistIntegrand, typename TGridFunction>
1466number GridFunctionDistance2(TGridFunction& spGridFct1, const char* cmp1,
1467 TGridFunction& spGridFct2, const char* cmp2,
1468 int quadOrder, const char* subsets)
1469{
1470// get function id of name
1471 const size_t fct1 = spGridFct1.fct_id_by_name(cmp1);
1472 const size_t fct2 = spGridFct2.fct_id_by_name(cmp2);
1473
1474// check that function exists
1475 if(fct1 >= spGridFct1.num_fct())
1476 UG_THROW("GridFunctionDistance: Function space does not contain"
1477 " a function with name " << cmp1 << ".");
1478 if(fct2 >= spGridFct2.num_fct())
1479 UG_THROW("GridFunctionDistance: Function space does not contain"
1480 " a function with name " << cmp2 << ".");
1481
1482// get top level of grid functions
1483 const int level1 = spGridFct1.dof_distribution()->grid_level().level();
1484 const int level2 = spGridFct2.dof_distribution()->grid_level().level();
1485
1486 if(level1 > level2){
1487 // level check
1488 TDistIntegrand spIntegrand(spGridFct1, fct1, spGridFct2, fct2);
1489 return IntegrateSubsets(spIntegrand, spGridFct1, subsets, quadOrder);
1490 }else{
1491 TDistIntegrand spIntegrand(spGridFct2, fct2, spGridFct1, fct1);
1492 return IntegrateSubsets(spIntegrand, spGridFct2, subsets, quadOrder);
1493 }
1494
1495}
1496
1498template <typename TDistIntegrand, typename TGridFunction>
1499number GridFunctionDistance2(TGridFunction& spGridFct1, const char* cmp1,
1500 TGridFunction& spGridFct2, const char* cmp2,
1501 int quadOrder, const char* subsets, ConstSmartPtr<typename TDistIntegrand::weight_type> spWeights)
1502{
1503// get function id of name
1504 const size_t fct1 = spGridFct1.fct_id_by_name(cmp1);
1505 const size_t fct2 = spGridFct2.fct_id_by_name(cmp2);
1506
1507// check that function exists
1508 if(fct1 >= spGridFct1.num_fct())
1509 UG_THROW("GridFunctionDistance: Function space does not contain"
1510 " a function with name " << cmp1 << ".");
1511 if(fct2 >= spGridFct2.num_fct())
1512 UG_THROW("GridFunctionDistance: Function space does not contain"
1513 " a function with name " << cmp2 << ".");
1514
1515// get top level of gridfunctions
1516 const int level1 = spGridFct1.dof_distribution()->grid_level().level();
1517 const int level2 = spGridFct2.dof_distribution()->grid_level().level();
1518
1519 // w/ weights
1520 if(level1 > level2){
1521 TDistIntegrand spIntegrand(spGridFct1, fct1, spGridFct2, fct2, spWeights);
1522 return IntegrateSubsets(spIntegrand, spGridFct1, subsets, quadOrder);
1523 }else{
1524 TDistIntegrand spIntegrand(spGridFct2, fct2, spGridFct1, fct1, spWeights);
1525 return IntegrateSubsets(spIntegrand, spGridFct2, subsets, quadOrder);
1526 }
1527}
1528
1529
1531template <typename TDistIntegrand, typename TGridFunction>
1532number GridFunctionDistance2(TGridFunction& spGridFct1, const char* cmp1,
1533 TGridFunction& spGridFct2, const char* cmp2,
1534 int quadOrder, const char* subsets,
1536 number distAvg12)
1537{
1538// get function id of name
1539 const size_t fct1 = spGridFct1.fct_id_by_name(cmp1);
1540 const size_t fct2 = spGridFct2.fct_id_by_name(cmp2);
1541
1542// check that function exists
1543 if(fct1 >= spGridFct1.num_fct())
1544 UG_THROW("GridFunctionDistance: Function space does not contain"
1545 " a function with name " << cmp1 << ".");
1546 if(fct2 >= spGridFct2.num_fct())
1547 UG_THROW("GridFunctionDistance: Function space does not contain"
1548 " a function with name " << cmp2 << ".");
1549
1550// get top level of gridfunctions
1551 const int level1 = spGridFct1.dof_distribution()->grid_level().level();
1552 const int level2 = spGridFct2.dof_distribution()->grid_level().level();
1553
1554 // w/ weights
1555 if(level1 > level2){
1556 TDistIntegrand spIntegrand(spGridFct1, fct1, spGridFct2, fct2, spWeights, distAvg12);
1557 return IntegrateSubsets(spIntegrand, spGridFct1, subsets, quadOrder);
1558 }else{
1559 TDistIntegrand spIntegrand(spGridFct2, fct2, spGridFct1, fct1, spWeights, -distAvg12);
1560 return IntegrateSubsets(spIntegrand, spGridFct2, subsets, quadOrder);
1561 }
1562}
1563
1564
1565
1567// L2 Integrand
1569
1571template <typename TGridFunction>
1573 : public StdIntegrand<number, TGridFunction::dim, L2Integrand<TGridFunction> >
1574{
1575 public:
1576 // world dimension of grid function
1577 static const int worldDim = TGridFunction::dim;
1579
1580 protected:
1581 // grid function data
1583
1586
1587 public:
1589 L2Integrand(TGridFunction& spGridFct, size_t cmp)
1590 : m_scalarData(spGridFct, cmp), m_spWeight(make_sp(new ConstUserNumber<worldDim>(1.0)))
1591 {};
1592
1593 L2Integrand(TGridFunction& spGridFct, size_t cmp, ConstSmartPtr<weight_type> spWeight)
1594 : m_scalarData(spGridFct, cmp), m_spWeight(spWeight)
1595 {};
1596
1598 virtual ~L2Integrand() {};
1599
1601 virtual void set_subset(int si)
1602 {
1603 if(!m_scalarData.is_def_in_subset(si))
1604 UG_THROW("L2ErrorIntegrand: Grid function component" <<m_scalarData.fct() <<" not defined on subset "<<si);
1606 }
1607
1609 template <int elemDim>
1610 void evaluate(number vValue[],
1611 const MathVector<worldDim> vGlobIP[],
1612 GridObject* pElem,
1613 const MathVector<worldDim> vCornerCoords[],
1614 const MathVector<elemDim> vLocIP[],
1616 const size_t numIP)
1617 {
1618 // get reference object id (i.e. Triangle, Quadrilateral, Tetrahedron, ...)
1620
1621 // element weights
1622 typedef typename weight_type::data_type ipdata_type;
1623 std::vector<ipdata_type> locElemWeights(numIP, 1.0);
1624 UG_ASSERT(m_spWeight.valid(), "L2Integrand::evaluate requires valid weights!");
1625 (*m_spWeight)(&locElemWeights[0], vGlobIP, 0.0, IIntegrand<number, worldDim>::subset(), numIP);
1626
1627 try{
1628 // get trial space
1629 const LocalShapeFunctionSet<elemDim>& rTrialSpace =
1630 LocalFiniteElementProvider::get<elemDim>(roid, m_scalarData.id());
1631
1632 // number of dofs on element
1633 const size_t num_sh = rTrialSpace.num_sh();
1634
1635 // get multiindices of element
1636 std::vector<DoFIndex> ind; // aux. index array
1637 m_scalarData.dof_indices(pElem, ind);
1638
1639 // check multi indices
1640 if(ind.size() != num_sh)
1641 UG_THROW("L2Integrand::evaluate: Wrong number of multi indices.");
1642
1643 // loop all integration points
1644 for(size_t ip = 0; ip < numIP; ++ip)
1645 {
1646
1647 // compute approximated solution at integration point
1648 number approxSolIP = 0.0;
1649 for(size_t sh = 0; sh < num_sh; ++sh)
1650 {
1651 // get value at shape point (e.g. corner for P1 fct)
1652 // and add shape fct at ip * value at shape
1653 const number valSH = DoFRef(m_scalarData.grid_function(), ind[sh]);
1654 approxSolIP += valSH * rTrialSpace.shape(sh, vLocIP[ip]);
1655 }
1656
1657 // get square
1658 vValue[ip] = locElemWeights[ip]*approxSolIP*approxSolIP;
1659
1660 }
1661
1662 }
1663 UG_CATCH_THROW("L2FuncIntegrand::values: trial space missing.");
1664 }
1665};
1666
1679template <typename TGridFunction>
1680number L2Norm2(TGridFunction& u, const char* cmp,
1681 int quadOrder, const char* subsets,
1683{
1684// get function id of name
1685 const size_t fct = u.fct_id_by_name(cmp);
1686
1687// check that function exists
1688 UG_COND_THROW(fct >= u.num_fct(), "L2Norm: Function space does not contain"
1689 " a function with name " << cmp << ".");
1690
1691 L2Integrand<TGridFunction> integrandL2(u, fct, spWeight);
1692 return IntegrateSubsets(integrandL2, u, subsets, quadOrder);
1693}
1694
1695template <typename TGridFunction>
1696number L2Norm2(TGridFunction& u, const char* cmp,
1697 int quadOrder, const char* subsets)
1698{
1699// get function id of name
1700 const size_t fct = u.fct_id_by_name(cmp);
1701
1702// check that function exists
1703 UG_COND_THROW(fct >= u.num_fct(), "L2Norm: Function space does not contain"
1704 " a function with name " << cmp << ".");
1705
1706 L2Integrand<TGridFunction> integrandL2(u, fct);
1707 return IntegrateSubsets(integrandL2, u, subsets, quadOrder);
1708}
1709
1710template <typename TGridFunction>
1711number L2Norm(TGridFunction& u, const char* cmp,
1712 int quadOrder, const char* subsets)
1713{
1714 return sqrt(L2Norm2(u, cmp, quadOrder, subsets));
1715}
1724template <typename TGridFunction>
1725number L2Norm(TGridFunction& gridFct, const char* cmp, int quadOrder)
1726{ return L2Norm(gridFct, cmp, quadOrder, NULL); }
1727
1728template <typename TGridFunction>
1729number L2Norm(SmartPtr<TGridFunction> spGridFct, const char* cmp, int quadOrder, const char* subsets)
1730{ return L2Norm(*spGridFct, cmp, quadOrder, subsets); }
1731
1732template <typename TGridFunction>
1733number L2Norm(SmartPtr<TGridFunction> spGridFct, const char* cmp, int quadOrder)
1734{ return L2Norm(spGridFct, cmp, quadOrder, NULL); }
1735
1736
1738template <typename TGridFunction>
1740 : public StdIntegrand<number, TGridFunction::dim, L2DistIntegrand<TGridFunction> >
1741{
1742 public:
1744 static const int worldDim = TGridFunction::dim;
1746
1747 protected:
1750
1753
1756
1758
1761
1762 public:
1763
1765 L2DistIntegrand(TGridFunction& fineGridFct, size_t fineCmp,
1766 TGridFunction& coarseGridFct, size_t coarseCmp)
1767 : m_fineData(fineGridFct, fineCmp), m_fineTopLevel(fineGridFct.dof_distribution()->grid_level().level()),
1768 m_coarseData(coarseGridFct, coarseCmp), m_coarseTopLevel(coarseGridFct.dof_distribution()->grid_level().level()),
1769 m_spMG(m_fineData.domain()->grid()), m_spWeight(make_sp(new ConstUserNumber<TGridFunction::dim>(1.0))),
1771 {
1773 "L2DiffIntegrand: fine and top level inverted.");
1774 UG_COND_THROW(m_fineData.domain().get() != m_coarseData.domain().get(),
1775 "L2DiffIntegrand: grid functions defined on different domains.");
1776 };
1777
1779 L2DistIntegrand(TGridFunction& fineGridFct, size_t fineCmp,
1780 TGridFunction& coarseGridFct, size_t coarseCmp, ConstSmartPtr<weight_type> spWeight)
1781 : m_fineData(fineGridFct, fineCmp), m_fineTopLevel(fineGridFct.dof_distribution()->grid_level().level()),
1782 m_coarseData(coarseGridFct, coarseCmp), m_coarseTopLevel(coarseGridFct.dof_distribution()->grid_level().level()),
1783 m_spMG(m_fineData.domain()->grid()), m_spWeight(spWeight),
1785 {
1787 "L2DiffIntegrand: fine and top level inverted.");
1788 UG_COND_THROW(m_fineData.domain().get() != m_coarseData.domain().get(),
1789 "L2DiffIntegrand: grid functions defined on different domains.");
1790 };
1791
1793 L2DistIntegrand(TGridFunction& fineGridFct, size_t fineCmp,
1794 TGridFunction& coarseGridFct, size_t coarseCmp, ConstSmartPtr<weight_type> spWeight, number dist12)
1795 : m_fineData(fineGridFct, fineCmp), m_fineTopLevel(fineGridFct.dof_distribution()->grid_level().level()),
1796 m_coarseData(coarseGridFct, coarseCmp), m_coarseTopLevel(coarseGridFct.dof_distribution()->grid_level().level()),
1797 m_spMG(m_fineData.domain()->grid()), m_spWeight(spWeight),
1798 m_deltaFineCoarse(dist12)
1799 {
1801 "L2DiffIntegrand: fine and top level inverted.");
1802 UG_COND_THROW(m_fineData.domain().get() != m_coarseData.domain().get(),
1803 "L2DiffIntegrand: grid functions defined on different domains.");
1804 };
1805
1806
1807 virtual ~L2DistIntegrand() {}
1808
1810 virtual void set_subset(int si)
1811 {
1812 UG_COND_THROW(!m_fineData.is_def_in_subset(si),
1813 "L2DiffIntegrand: Grid function component" <<m_fineData.fct()<<" not defined on subset "<<si);
1814 UG_COND_THROW(!m_coarseData.is_def_in_subset(si),
1815 "L2DiffIntegrand: Grid function component" <<m_coarseData.fct()<<" not defined on subset "<<si);
1817 }
1818
1820 template <int elemDim>
1821 void evaluate(number vValue[],
1822 const MathVector<worldDim> vFineGlobIP[],
1823 GridObject* pFineElem,
1824 const MathVector<worldDim> vCornerCoords[],
1825 const MathVector<elemDim> vFineLocIP[],
1827 const size_t numIP)
1828 {
1829 typedef typename TGridFunction::template dim_traits<elemDim>::grid_base_object Element;
1830
1831 // get coarse element
1832 GridObject* pCoarseElem = pFineElem;
1834 int parentLevel = m_spMG->get_level(pCoarseElem);
1835 while(parentLevel > m_coarseTopLevel){
1836 pCoarseElem = m_spMG->get_parent(pCoarseElem);
1837 parentLevel = m_spMG->get_level(pCoarseElem);
1838 }
1839 }
1840
1841 // get reference object id (i.e. Triangle, Quadrilateral, Tetrahedron, ...)
1842 const ReferenceObjectID fineROID = pFineElem->reference_object_id();
1843 const ReferenceObjectID coarseROID = pCoarseElem->reference_object_id();
1844
1845 // get corner coordinates
1846 std::vector<MathVector<worldDim> > vCornerCoarse;
1847 CollectCornerCoordinates(vCornerCoarse, *static_cast<Element*>(pCoarseElem), *m_coarseData.domain());
1848
1849 // get Reference Mapping
1851 = ReferenceMappingProvider::get<elemDim, worldDim>(coarseROID, vCornerCoarse);
1852
1853 std::vector<MathVector<elemDim> > vCoarseLocIP;
1854 vCoarseLocIP.resize(numIP);
1855 for(size_t ip = 0; ip < vCoarseLocIP.size(); ++ip) VecSet(vCoarseLocIP[ip], 0.0);
1856 map.global_to_local(&vCoarseLocIP[0], vFineGlobIP, numIP);
1857
1858 // element weights
1859 typedef typename weight_type::data_type ipdata_type;
1860 std::vector<ipdata_type> fineElemWeights(numIP, 1.0);
1861 UG_ASSERT(m_spWeight.valid(), "L2DistIntegrand::evaluate requires valid weights! ");
1862 (*m_spWeight)(&fineElemWeights[0], vFineGlobIP, 0.0, IIntegrand<number, worldDim>::subset(), numIP);
1863
1864 try{
1865 // get trial space
1866 const LocalShapeFunctionSet<elemDim>& rFineLSFS =
1867 LocalFiniteElementProvider::get<elemDim>(fineROID, m_fineData.id());
1868 const LocalShapeFunctionSet<elemDim>& rCoarseLSFS =
1869 LocalFiniteElementProvider::get<elemDim>(coarseROID, m_coarseData.id());
1870
1871 // get multiindices of element
1872 std::vector<DoFIndex> vFineMI, vCoarseMI;
1873 m_fineData.dof_indices(pFineElem, vFineMI);
1874 m_coarseData.dof_indices(pCoarseElem, vCoarseMI);
1875
1876 // loop all integration points
1877 for(size_t ip = 0; ip < numIP; ++ip)
1878 {
1879 // compute approximated solution at integration point
1880 number fineSolIP = 0.0;
1881 for(size_t sh = 0; sh < vFineMI.size(); ++sh)
1882 {
1883 // get value at shape point (e.g. corner for P1 fct)
1884 const number val = DoFRef(m_fineData.grid_function(), vFineMI[sh]);
1885
1886 // add shape fct at ip * value at shape
1887 fineSolIP += val * rFineLSFS.shape(sh, vFineLocIP[ip]);
1888 }
1889 number coarseSolIP = 0.0;
1890 for(size_t sh = 0; sh < vCoarseMI.size(); ++sh)
1891 {
1892 // get value at shape point (e.g. corner for P1 fct)
1893 const number val = DoFRef(m_coarseData.grid_function(), vCoarseMI[sh]);
1894
1895 // add shape fct at ip * value at shape
1896 coarseSolIP += val * rCoarseLSFS.shape(sh, vCoarseLocIP[ip]);
1897 }
1898
1899 // get squared of difference
1900 vValue[ip] = fineElemWeights[ip]*(fineSolIP - coarseSolIP -m_deltaFineCoarse)*(fineSolIP-coarseSolIP-m_deltaFineCoarse);
1901 }
1902
1903 }
1904 UG_CATCH_THROW("L2DistIntegrand::evaluate: trial space missing.");
1905 }
1906};
1907
1908
1910template <typename TGridFunction>
1911number L2Distance2(TGridFunction& spGridFct1, const char* cmp1,
1912 TGridFunction& spGridFct2, const char* cmp2,
1913 int quadOrder, const char* subsets,
1914 ConstSmartPtr<typename L2Integrand<TGridFunction>::weight_type> spWeight, number avgDist12=0.0)
1915{
1916 return GridFunctionDistance2<L2DistIntegrand<TGridFunction>, TGridFunction>
1917 (spGridFct1, cmp1, spGridFct2, cmp2, quadOrder, subsets, spWeight, avgDist12);
1918}
1919
1920
1922template <typename TGridFunction>
1923number L2Distance2(TGridFunction& spGridFct1, const char* cmp1,
1924 TGridFunction& spGridFct2, const char* cmp2,
1925 int quadOrder, const char* subsets)
1926{
1927 return GridFunctionDistance2<L2DistIntegrand<TGridFunction>, TGridFunction>
1928 (spGridFct1, cmp1, spGridFct2, cmp2, quadOrder, subsets);
1929}
1930
1932template <typename TGridFunction>
1933number L2Distance(TGridFunction& spGridFct1, const char* cmp1,
1934 TGridFunction& spGridFct2, const char* cmp2,
1935 int quadOrder, const char* subsets)
1936{
1937 return sqrt(L2Distance2(spGridFct1, cmp1, spGridFct2, cmp2, quadOrder, subsets));
1938}
1939
1940
1941template <typename TGridFunction>
1942number L2Error(SmartPtr<TGridFunction> spGridFct1, const char* cmp1,
1943 SmartPtr<TGridFunction> spGridFct2, const char* cmp2,
1944 int quadOrder)
1945{
1946 return L2Distance(*spGridFct1, cmp1, *spGridFct2, cmp2, quadOrder, NULL);
1947}
1948
1949template <typename TGridFunction>
1950number L2Error(SmartPtr<TGridFunction> spGridFct1, const char* cmp1,
1951 SmartPtr<TGridFunction> spGridFct2, const char* cmp2,
1952 int quadOrder, const char* subsets)
1953{
1954 return L2Distance(*spGridFct1, cmp1, *spGridFct2, cmp2, quadOrder, subsets);
1955}
1956
1958// H1 semi-norm Integrand
1960
1961
1962
1964template <typename TGridFunction>
1966 : public StdIntegrand<number, TGridFunction::dim, H1SemiIntegrand<TGridFunction> >
1967{
1968 public:
1970 static const int worldDim = TGridFunction::dim;
1972
1973 protected:
1976
1979
1980 public:
1982 H1SemiIntegrand(TGridFunction& gridFct, size_t cmp)
1983 : m_scalarData(gridFct, cmp), m_spWeight(make_sp(new ConstUserMatrix<worldDim>(1.0))) {}
1984
1986 H1SemiIntegrand(TGridFunction& gridFct, size_t cmp, ConstSmartPtr<weight_type> spWeight)
1987 : m_scalarData(gridFct, cmp), m_spWeight(spWeight) {}
1988
1990 virtual ~H1SemiIntegrand() {};
1991
1993 virtual void set_subset(int si)
1994 {
1995
1996 UG_COND_THROW(!m_scalarData.is_def_in_subset(si), "H1Error: Grid function component"
1997 <<m_scalarData.fct()<<" not defined on subset "<<si);
1999 }
2000
2002 template <int elemDim>
2003 void evaluate(number vValue[],
2004 const MathVector<worldDim> vGlobIP[],
2005 GridObject* pElem,
2006 const MathVector<worldDim> vCornerCoords[],
2007 const MathVector<elemDim> vLocIP[],
2009 const size_t numIP)
2010 {
2011 // get reference object id (i.e. Triangle, Quadrilateral, Tetrahedron, ...)
2012 const ReferenceObjectID roid = pElem->reference_object_id();
2013 const TGridFunction &gridFct= m_scalarData.grid_function();
2014
2016 = ReferenceMappingProvider::get<elemDim, worldDim>(roid, vCornerCoords);
2017
2018 typedef typename weight_type::data_type ipdata_type;
2019
2020 std::vector<ipdata_type> elemWeights(numIP, MathMatrix<worldDim, worldDim>());
2021 UG_ASSERT(m_spWeight.valid(), "H1SemiIntegrand::evaluate requires valid weights!");
2022
2023
2025 {
2026 // get local solution if needed
2027 LocalIndices ind;
2028 LocalVector uloc;
2029 gridFct.indices(pElem, ind); // get global indices
2030 uloc.resize(ind); // adapt local algebra
2031 GetLocalVector(uloc, gridFct); // read local values of u
2032
2033 // compute data
2034 try{
2035 (*m_spWeight)(&elemWeights[0], vGlobIP, 0.0, IIntegrand<number, worldDim>::subset(), pElem,
2036 vCornerCoords, vLocIP, numIP, &uloc, NULL);
2037 } UG_CATCH_THROW("H1SemiIntegrand: Cannot evaluate weight data.");
2038 }
2039 else
2040 {
2041 // compute data
2042 try{
2043 (*m_spWeight)(&elemWeights[0], vGlobIP, 0.0, IIntegrand<number, worldDim>::subset(), numIP);
2044 } UG_CATCH_THROW("H1SemiIntegrand: Cannot evaluate weight data.");
2045 }
2046
2047
2048 try{
2049
2050
2051
2052
2053 // get trial space
2054 const LocalShapeFunctionSet<elemDim>& rTrialSpace =
2055 LocalFiniteElementProvider::get<elemDim>(roid, m_scalarData.id());
2056
2057 // number of dofs on element
2058 const size_t num_sh = rTrialSpace.num_sh();
2059
2060 // get multiindices of element
2061 std::vector<DoFIndex> ind; // aux. index array
2062 gridFct.dof_indices(pElem, m_scalarData.fct(), ind);
2063
2064 // check multi indices
2065 UG_COND_THROW(ind.size() != num_sh, "H1SemiNormFuncIntegrand::evaluate: Wrong number of multi-)indices.");
2066
2067 // loop all integration points
2068 std::vector<MathVector<elemDim> > vLocGradient(num_sh);
2069 for(size_t ip = 0; ip < numIP; ++ip)
2070 {
2071 // compute shape gradients at ip
2072 rTrialSpace.grads(&vLocGradient[0], vLocIP[ip]);
2073
2074 // compute approximated solution at integration point
2075 number approxSolIP = 0.0;
2076 MathVector<elemDim> tmpVec(0.0);
2077 for(size_t sh = 0; sh < num_sh; ++sh)
2078 {
2079 // get value at shape point (e.g. corner for P1 fct)
2080 const number valSH = DoFRef(gridFct, ind[sh]);
2081
2082 // add shape fct at ip * value at shape
2083 approxSolIP += valSH * rTrialSpace.shape(sh, vLocIP[ip]);
2084
2085 // add gradient at ip
2086 VecScaleAppend(tmpVec, valSH, vLocGradient[sh]);
2087 }
2088
2089 // compute gradient
2090 MathVector<worldDim> approxGradIP;
2092 map.jacobian_transposed_inverse(JTInv, vLocIP[ip]);//Inverse(JTInv, vJT[ip]);
2093
2094#ifdef UG_DEBUG
2095 CheckGeneralizedInverse<worldDim,elemDim>(vJT[ip], JTInv);
2096#endif //UG_DEBUG
2097
2098 MatVecMult(approxGradIP, JTInv, tmpVec);
2099
2100 // get norm squared
2101 MathVector<worldDim> approxDGradIP;
2102 MatVecMult(approxDGradIP, elemWeights[ip], approxGradIP);
2103 vValue[ip] = VecDot(approxDGradIP, approxGradIP);
2104 }
2105
2106 }
2107 UG_CATCH_THROW("H1SemiIntegrand::evaluate: trial space missing.");
2108 }
2109};
2110
2111
2113
2124template <typename TGridFunction>
2125number H1SemiNorm2(TGridFunction& gridFct, const char* cmp, int quadOrder, const char* subsets=NULL,
2127{
2128// get function id of name
2129 const size_t fct = gridFct.fct_id_by_name(cmp);
2130
2131// check that function exists
2132 UG_COND_THROW(fct >= gridFct.num_fct(), "H1SemiNorm: Function space does not contain"
2133 " a function with name " << cmp << ".");
2134 if (weights.invalid()) {
2135 H1SemiIntegrand<TGridFunction> integrand(gridFct, fct);
2136 return IntegrateSubsets(integrand, gridFct, subsets, quadOrder);
2137 } else {
2138 H1SemiIntegrand<TGridFunction> integrand(gridFct, fct, weights);
2139 return IntegrateSubsets(integrand, gridFct, subsets, quadOrder);
2140 }
2141}
2142
2144template <typename TGridFunction>
2145number H1SemiNorm(TGridFunction& gridFct, const char* cmp, int quadOrder, const char* subsets=NULL,
2147{
2148 return (sqrt(H1SemiNorm2(gridFct, cmp, quadOrder, subsets, weights)));
2149}
2150
2151// Delegating to H1SemiNorm
2152template <typename TGridFunction>
2153number H1SemiNorm(SmartPtr<TGridFunction> spGridFct, const char* cmp, int quadOrder, const char* subsets)
2154{ return H1SemiNorm(*spGridFct, cmp, quadOrder, subsets); }
2155
2156template <typename TGridFunction>
2157number H1SemiNorm(SmartPtr<TGridFunction> spGridFct, const char* cmp, int quadOrder,
2158 const char* subsets, ConstSmartPtr<typename H1SemiIntegrand<TGridFunction>::weight_type> weights = SPNULL)
2159{ return H1SemiNorm(*spGridFct, cmp, quadOrder, subsets, weights); }
2160
2161template <typename TGridFunction>
2162number H1SemiNorm( SmartPtr<TGridFunction> spGridFct, const char* cmp, int quadOrder)
2163{ return H1SemiNorm(*spGridFct, cmp, quadOrder, NULL); }
2164
2165template <typename TGridFunction>
2166number H1SemiNorm( SmartPtr<TGridFunction> spGridFct, const char* cmp, int quadOrder,
2168{ return H1SemiNorm(*spGridFct, cmp, quadOrder, NULL, weights); }
2169
2171template <typename TGridFunction>
2172class H1SemiDistIntegrand : public StdIntegrand<number, TGridFunction::dim, H1SemiDistIntegrand<TGridFunction> >
2173{
2174 public:
2176 static const int worldDim = TGridFunction::dim;
2178
2179 private:
2182
2185
2188
2191
2192 public:
2194 H1SemiDistIntegrand(TGridFunction& fineGridFct, size_t fineCmp,
2195 TGridFunction& coarseGridFct, size_t coarseCmp)
2196 : m_fineData(fineGridFct, fineCmp),
2197 m_fineTopLevel(fineGridFct.dof_distribution()->grid_level().level()),
2198 m_coarseData(coarseGridFct, coarseCmp),
2199 m_coarseTopLevel(coarseGridFct.dof_distribution()->grid_level().level()),
2200 m_spMG(fineGridFct.domain()->grid()),
2201 m_spWeight(new ConstUserNumber<TGridFunction::dim>(1.0))
2202 {
2204 UG_THROW("H1SemiDiffIntegrand: fine and top level inverted.");
2205
2206 if(m_fineData.domain().get() != m_coarseData.domain().get())
2207 UG_THROW("H1SemiDiffIntegrand: grid functions defined on different domains.");
2208 };
2209
2211 H1SemiDistIntegrand(TGridFunction& fineGridFct, size_t fineCmp,
2212 TGridFunction& coarseGridFct, size_t coarseCmp,
2214 : m_fineData(fineGridFct, fineCmp),
2215 m_fineTopLevel(fineGridFct.dof_distribution()->grid_level().level()),
2216 m_coarseData(coarseGridFct, coarseCmp),
2217 m_coarseTopLevel(coarseGridFct.dof_distribution()->grid_level().level()),
2218 m_spMG(fineGridFct.domain()->grid()),
2219 m_spWeight(spWeight)
2220 {
2221 UG_COND_THROW(m_fineTopLevel < m_coarseTopLevel, "H1SemiDiffIntegrand: fine and top level inverted.");
2222 UG_COND_THROW(m_fineData.domain().get() != m_coarseData.domain().get(), "H1SemiDiffIntegrand: grid functions defined on different domains.");
2223 }
2225
2227 virtual void set_subset(int si)
2228 {
2229
2230 UG_COND_THROW(!m_fineData.is_def_in_subset(si), "H1SemiDiffIntegrand: Grid function component"
2231 <<m_fineData.fct()<<" not defined on subset "<<si);
2232 UG_COND_THROW(!m_coarseData.is_def_in_subset(si), "H1SemiDiffIntegrand: Grid function component"
2233 <<m_coarseData.fct()<<" not defined on subset "<<si);
2235 }
2236
2238 template <int elemDim>
2239 void evaluate(number vValue[],
2240 const MathVector<worldDim> vFineGlobIP[],
2241 GridObject* pFineElem,
2242 const MathVector<worldDim> vCornerCoords[],
2243 const MathVector<elemDim> vFineLocIP[],
2245 const size_t numIP)
2246 {
2247 typedef typename TGridFunction::template dim_traits<elemDim>::grid_base_object Element;
2248
2249 const TGridFunction &fineGridFct = m_fineData.grid_function();
2250 const TGridFunction &coarseGridFct = m_coarseData.grid_function();
2251
2252 // get coarse element
2253 GridObject* pCoarseElem = pFineElem;
2255 int parentLevel = m_spMG->get_level(pCoarseElem);
2256 while(parentLevel > m_coarseTopLevel){
2257 pCoarseElem = m_spMG->get_parent(pCoarseElem);
2258 parentLevel = m_spMG->get_level(pCoarseElem);
2259 }
2260 }
2261
2262 // get reference object id (i.e. Triangle, Quadrilateral, Tetrahedron, ...)
2263 const ReferenceObjectID fineROID = pFineElem->reference_object_id();
2264 const ReferenceObjectID coarseROID = pCoarseElem->reference_object_id();
2265
2266 // get corner coordinates
2267 std::vector<MathVector<worldDim> > vCornerCoarse;
2268 CollectCornerCoordinates(vCornerCoarse, *static_cast<Element*>(pCoarseElem), *m_coarseData.domain());
2269
2270 // get reference Mapping
2272 = ReferenceMappingProvider::get<elemDim, worldDim>(coarseROID, vCornerCoarse);
2274 = ReferenceMappingProvider::get<elemDim, worldDim>(fineROID, vCornerCoords);
2275
2276 std::vector<MathVector<elemDim> > vCoarseLocIP;
2277 vCoarseLocIP.resize(numIP);
2278 for(size_t ip = 0; ip < vCoarseLocIP.size(); ++ip) VecSet(vCoarseLocIP[ip], 0.0);
2279 map.global_to_local(&vCoarseLocIP[0], vFineGlobIP, numIP);
2280
2281
2282 // determine weights
2283 std::vector<typename weight_type::data_type> elemWeights(numIP, MathMatrix<worldDim, worldDim>());
2284 UG_ASSERT(m_spWeight.valid(), "H1SemiDistIntegrand::evaluate requires valid weights!");
2285
2286 // get local solution (if required)
2288 {
2289
2290 LocalIndices ind;
2291 LocalVector u;
2292 fineGridFct.indices(pFineElem, ind); // get global indices
2293 u.resize(ind); // adapt local algebra
2294 GetLocalVector(u, fineGridFct); // read local values of u
2295
2296 // compute data
2297 try{
2298 (*m_spWeight)(&elemWeights[0], vFineGlobIP, 0.0, IIntegrand<number, worldDim>::subset(), pFineElem,
2299 vCornerCoords, vFineLocIP, numIP, &u, NULL);
2300 } UG_CATCH_THROW("H1SemiDistIntegrand: Cannot evaluate data.");
2301 }
2302 else
2303 {
2304 // compute data
2305 try{
2306 (*m_spWeight)(&elemWeights[0], vFineGlobIP, 0.0, IIntegrand<number, worldDim>::subset(), numIP);
2307 } UG_CATCH_THROW("H1SemiDistIntegrand: Cannot evaluate data.");
2308 }
2309
2310 try{
2311
2312
2313 // get trial space
2314 const LocalShapeFunctionSet<elemDim>& rFineLSFS =
2315 LocalFiniteElementProvider::get<elemDim>(fineROID, m_fineData.id());
2316 const LocalShapeFunctionSet<elemDim>& rCoarseLSFS =
2317 LocalFiniteElementProvider::get<elemDim>(coarseROID, m_coarseData.id());
2318
2319 // get multiindices of element
2320 std::vector<DoFIndex> vFineMI, vCoarseMI;
2321 m_fineData.dof_indices(pFineElem, vFineMI);
2322 m_coarseData.dof_indices(pCoarseElem, vCoarseMI);
2323
2324 std::vector<MathVector<elemDim> > vFineLocGradient(vFineMI.size());
2325 std::vector<MathVector<elemDim> > vCoarseLocGradient(vCoarseMI.size());
2326
2327 // loop all integration points
2328 for(size_t ip = 0; ip < numIP; ++ip)
2329 {
2330 // compute shape gradients at ip
2331 rFineLSFS.grads(&vFineLocGradient[0], vFineLocIP[ip]);
2332 rCoarseLSFS.grads(&vCoarseLocGradient[0], vCoarseLocIP[ip]);
2333
2334 // compute approximated solutions at integration point
2335 number fineSolIP = 0.0;
2336 MathVector<elemDim> fineLocTmp(0.0);
2337 for(size_t sh = 0; sh < vFineMI.size(); ++sh)
2338 {
2339 const number val = DoFRef(fineGridFct, vFineMI[sh]);
2340 fineSolIP += val * rFineLSFS.shape(sh, vFineLocIP[ip]);
2341 VecScaleAppend(fineLocTmp, val, vFineLocGradient[sh]);
2342 }
2343
2344 number coarseSolIP = 0.0;
2345 MathVector<elemDim> coarseLocTmp(0.0);
2346 for(size_t sh = 0; sh < vCoarseMI.size(); ++sh)
2347 {
2348 const number val = DoFRef(coarseGridFct, vCoarseMI[sh]);
2349 coarseSolIP += val * rCoarseLSFS.shape(sh, vCoarseLocIP[ip]);
2350 VecScaleAppend(coarseLocTmp, val, vCoarseLocGradient[sh]);
2351 }
2352
2353 // compute global gradient
2354 MathVector<worldDim> fineGradIP;
2356 mapF.jacobian_transposed_inverse(fineJTInv, vFineLocIP[ip]);//Inverse(fineJTInv, vJT[ip]);
2357
2358#ifdef UG_DEBUG
2359 CheckGeneralizedInverse<worldDim,elemDim>(vJT[ip], fineJTInv);
2360#endif //UG_DEBUG
2361
2362 MatVecMult(fineGradIP, fineJTInv, fineLocTmp);
2363
2364 // compute global gradient
2365 MathVector<worldDim> coarseGradIP;
2367 map.jacobian_transposed_inverse(coarseJTInv, vCoarseLocIP[ip]);
2368 MatVecMult(coarseGradIP, coarseJTInv, coarseLocTmp);
2369
2370 // get squared of difference
2371 /*vValue[ip] = (coarseSolIP - fineSolIP);
2372 vValue[ip] *= vValue[ip]; */
2373 vValue[ip] = VecDistanceSq(coarseGradIP, fineGradIP, elemWeights[ip]);
2374 }
2375
2376 }
2377 UG_CATCH_THROW("H1SemiDiffIntegrand::evaluate: trial space missing.");
2378 }
2379};
2380
2381
2383template <typename TGridFunction>
2384number H1SemiError2(SmartPtr<TGridFunction> spGridFct1, const char* cmp1,
2385 SmartPtr<TGridFunction> spGridFct2, const char* cmp2,
2386 int quadOrder, const char* subsets)
2387{
2388 return GridFunctionDistance2<H1SemiDistIntegrand<TGridFunction>, TGridFunction>
2389 (*spGridFct1, cmp1, *spGridFct2, cmp2, quadOrder, subsets);
2390}
2391
2393template <typename TGridFunction>
2394number H1SemiError(SmartPtr<TGridFunction> spGridFct1, const char* cmp1,
2395 SmartPtr<TGridFunction> spGridFct2, const char* cmp2,
2396 int quadOrder, const char* subsets)
2397{
2398 return sqrt(H1SemiError2(*spGridFct1, cmp1, *spGridFct2, cmp2, quadOrder, subsets));
2399}
2400
2402template <typename TGridFunction>
2403number H1SemiError(SmartPtr<TGridFunction> spGridFct1, const char* cmp1,
2404 SmartPtr<TGridFunction> spGridFct2, const char* cmp2,
2405 int quadOrder)
2406{
2407 return H1SemiError(*spGridFct1, cmp1, *spGridFct2, cmp2, quadOrder, NULL);
2408}
2409
2411template <typename TGridFunction>
2412number H1SemiDistance2(TGridFunction& spGridFct1, const char* cmp1,
2413 TGridFunction& spGridFct2, const char* cmp2,
2414 int quadOrder, const char* subsets,
2416{
2417 return GridFunctionDistance2<H1SemiDistIntegrand<TGridFunction>, TGridFunction>
2418 (spGridFct1, cmp1, spGridFct2, cmp2, quadOrder, subsets, weights);
2419}
2420
2422template <typename TGridFunction>
2423number H1SemiDistance(TGridFunction& spGridFct1, const char* cmp1,
2424 TGridFunction& spGridFct2, const char* cmp2,
2425 int quadOrder, const char* subsets,
2427{ return sqrt(H1SemiDistance2(spGridFct1, cmp1, spGridFct2, cmp2, quadOrder, subsets, weights)); }
2428
2430template <typename TGridFunction>
2431number H1SemiDistance2(TGridFunction& spGridFct1, const char* cmp1,
2432 TGridFunction& spGridFct2, const char* cmp2,
2433 int quadOrder, ConstSmartPtr<typename H1SemiDistIntegrand<TGridFunction>::weight_type > weights)
2434{ return H1SemiDistance2(spGridFct1, cmp1, spGridFct2, cmp2, quadOrder, NULL, weights); }
2435
2437template <typename TGridFunction>
2438number H1SemiDistance(TGridFunction& spGridFct1, const char* cmp1,
2439 TGridFunction& spGridFct2, const char* cmp2,
2440 int quadOrder, ConstSmartPtr<typename H1SemiDistIntegrand<TGridFunction>::weight_type > weights)
2441{ return sqrt(H1SemiDistance2(spGridFct1, cmp1, spGridFct2, cmp2, quadOrder, NULL, weights)); }
2442
2443
2444
2445
2447// H1 semi-norm Integrand
2449
2451template <typename TGridFunction>
2453 : public StdIntegrand<number, TGridFunction::dim, H1EnergyIntegrand<TGridFunction> >
2454{
2455 public:
2457 static const int worldDim = TGridFunction::dim;
2459
2460 protected:
2463
2466
2467 public:
2469 H1EnergyIntegrand(TGridFunction& gridFct, size_t cmp)
2470 : m_scalarData(gridFct, cmp), m_spWeight(make_sp(new ConstUserMatrix<worldDim>(1.0))) {}
2471
2473 H1EnergyIntegrand(TGridFunction& gridFct, size_t cmp, ConstSmartPtr<weight_type> spWeight)
2474 : m_scalarData(gridFct, cmp), m_spWeight(spWeight) {}
2475
2478
2480 virtual void set_subset(int si)
2481 {
2482
2483 UG_COND_THROW(!m_scalarData.is_def_in_subset(si), "H1EnergyIntegrand: Grid function component"
2484 <<m_scalarData.fct()<<" not defined on subset "<<si);
2486 }
2487
2489 template <int elemDim>
2490 void evaluate(number vValue[],
2491 const MathVector<worldDim> vGlobIP[],
2492 GridObject* pElem,
2493 const MathVector<worldDim> vCornerCoords[],
2494 const MathVector<elemDim> vLocIP[],
2496 const size_t numIP)
2497 {
2498 // get reference object id (i.e. Triangle, Quadrilateral, Tetrahedron, ...)
2499 const ReferenceObjectID roid = pElem->reference_object_id();
2500 const TGridFunction &gridFct= m_scalarData.grid_function();
2501
2503 = ReferenceMappingProvider::get<elemDim, worldDim>(roid, vCornerCoords);
2504
2505 typedef typename weight_type::data_type ipdata_type;
2506
2507 std::vector<ipdata_type> elemWeights(numIP, MathMatrix<worldDim, worldDim>());
2508 UG_ASSERT(m_spWeight.valid(), "H1EnergyIntegrand::evaluate requires valid weights!");
2509
2510
2511 if(m_spWeight->requires_grid_fct())
2512 {
2513 // get local solution if needed
2514 LocalIndices ind;
2515 LocalVector uloc;
2516 gridFct.indices(pElem, ind); // get global indices
2517 uloc.resize(ind); // adapt local algebra
2518 GetLocalVector(uloc, gridFct); // read local values of u
2519
2520 // compute data
2521 try{
2522 (*m_spWeight)(&elemWeights[0], vGlobIP, 0.0, IIntegrand<number, worldDim>::subset(), pElem,
2523 vCornerCoords, vLocIP, numIP, &uloc, NULL);
2524 } UG_CATCH_THROW("H1EnergyIntegrand: Cannot evaluate weight data.");
2525 }
2526 else
2527 {
2528 // compute data
2529 try{
2530 (*m_spWeight)(&elemWeights[0], vGlobIP, 0.0, IIntegrand<number, worldDim>::subset(), numIP);
2531 } UG_CATCH_THROW("H1EnergyIntegrand: Cannot evaluate weight data.");
2532 }
2533
2534
2535 try{
2536
2537
2538
2539
2540 // get trial space
2541 const LocalShapeFunctionSet<elemDim>& rTrialSpace =
2542 LocalFiniteElementProvider::get<elemDim>(roid, m_scalarData.id());
2543
2544 // number of dofs on element
2545 const size_t num_sh = rTrialSpace.num_sh();
2546
2547 // get multiindices of element
2548 std::vector<DoFIndex> ind; // aux. index array
2549 gridFct.dof_indices(pElem, m_scalarData.fct(), ind);
2550
2551 // check multi indices
2552 UG_COND_THROW(ind.size() != num_sh, "H1EnergyIntegrand::evaluate: Wrong number of multi-)indices.");
2553
2554 // loop all integration points
2555 std::vector<MathVector<elemDim> > vLocGradient(num_sh);
2556 for(size_t ip = 0; ip < numIP; ++ip)
2557 {
2558 // compute shape gradients at ip
2559 rTrialSpace.grads(&vLocGradient[0], vLocIP[ip]);
2560
2561 // compute approximated solution at integration point
2562 number approxSolIP = 0.0;
2563 MathVector<elemDim> tmpVec(0.0);
2564 for(size_t sh = 0; sh < num_sh; ++sh)
2565 {
2566 // get value at shape point (e.g. corner for P1 fct)
2567 const number valSH = DoFRef(gridFct, ind[sh]);
2568
2569 // add shape fct at ip * value at shape
2570 approxSolIP += valSH * rTrialSpace.shape(sh, vLocIP[ip]);
2571
2572 // add gradient at ip
2573 VecScaleAppend(tmpVec, valSH, vLocGradient[sh]);
2574 }
2575
2576 // compute gradient
2577 MathVector<worldDim> approxGradIP;
2579 map.jacobian_transposed_inverse(JTInv, vLocIP[ip]);//Inverse(JTInv, vJT[ip]);
2580
2581#ifdef UG_DEBUG
2582 CheckGeneralizedInverse<worldDim,elemDim>(vJT[ip], JTInv);
2583#endif //UG_DEBUG
2584
2585
2586 MatVecMult(approxGradIP, JTInv, tmpVec);
2587
2588 // get norm squared
2589 MathVector<worldDim> approxDGradIP;
2590 MatVecMult(approxDGradIP, elemWeights[ip], approxGradIP);
2591 vValue[ip] = VecTwoNormSq(approxDGradIP);
2592 }
2593
2594 }
2595 UG_CATCH_THROW("H1EnergyIntegrand::evaluate: trial space missing.");
2596 }
2597};
2598
2599
2601
2612template <typename TGridFunction>
2613number H1EnergyNorm2(TGridFunction& gridFct, const char* cmp, int quadOrder, const char* subsets=NULL,
2615{
2616// get function id of name
2617 const size_t fct = gridFct.fct_id_by_name(cmp);
2618
2619// check that function exists
2620 UG_COND_THROW(fct >= gridFct.num_fct(), "H1SemiNorm: Function space does not contain"
2621 " a function with name " << cmp << ".");
2622 if (weights.invalid()) {
2623 H1EnergyIntegrand<TGridFunction> integrand(gridFct, fct);
2624 return IntegrateSubsets(integrand, gridFct, subsets, quadOrder);
2625 } else {
2626 H1EnergyIntegrand<TGridFunction> integrand(gridFct, fct, weights);
2627 return IntegrateSubsets(integrand, gridFct, subsets, quadOrder);
2628 }
2629}
2630template <typename TGridFunction>
2631number H1EnergyNorm(TGridFunction& gridFct, const char* cmp, int quadOrder, const char* subsets=NULL,
2633{
2634 return (sqrt(H1EnergyNorm2(gridFct, cmp, quadOrder, subsets, weights)));
2635}
2636
2637template <typename TGridFunction>
2638number H1EnergyNorm(SmartPtr<TGridFunction> spGridFct, const char* cmp, int quadOrder,
2639 const char* subsets, ConstSmartPtr<typename H1SemiIntegrand<TGridFunction>::weight_type> weights = SPNULL)
2640{ return H1EnergyNorm(*spGridFct, cmp, quadOrder, subsets, weights); }
2641
2642template <typename TGridFunction>
2643number H1EnergyNorm( SmartPtr<TGridFunction> spGridFct, const char* cmp, int quadOrder)
2644{ return H1EnergyNorm(spGridFct, cmp, quadOrder, NULL); }
2645
2646template <typename TGridFunction>
2647number H1EnergyNorm( SmartPtr<TGridFunction> spGridFct, const char* cmp, int quadOrder,
2649{ return H1EnergyNorm(spGridFct, cmp, quadOrder, NULL, weights); }
2650
2651
2652
2653
2655template <typename TGridFunction>
2657 : public StdIntegrand<number, TGridFunction::dim, H1EnergyDistIntegrand<TGridFunction> >
2658{
2659 public:
2661 static const int worldDim = TGridFunction::dim;
2663
2664 private:
2667
2670
2673
2676
2677 public:
2679 H1EnergyDistIntegrand(TGridFunction& fineGridFct, size_t fineCmp,
2680 TGridFunction& coarseGridFct, size_t coarseCmp)
2681 : m_fineData(fineGridFct, fineCmp),
2682 m_fineTopLevel(fineGridFct.dof_distribution()->grid_level().level()),
2683 m_coarseData(coarseGridFct, coarseCmp),
2684 m_coarseTopLevel(coarseGridFct.dof_distribution()->grid_level().level()),
2685 m_spMG(fineGridFct.domain()->grid()),
2686 m_spWeight(new ConstUserNumber<TGridFunction::dim>(1.0))
2687 {
2689 UG_THROW("H1EnergyDistIntegrand: fine and top level inverted.");
2690
2691 if(m_fineData.domain().get() != m_coarseData.domain().get())
2692 UG_THROW("H1EnergyDistIntegrand: grid functions defined on different domains.");
2693 };
2694
2696 H1EnergyDistIntegrand(TGridFunction& fineGridFct, size_t fineCmp,
2697 TGridFunction& coarseGridFct, size_t coarseCmp,
2699 : m_fineData(fineGridFct, fineCmp),
2700 m_fineTopLevel(fineGridFct.dof_distribution()->grid_level().level()),
2701 m_coarseData(coarseGridFct, coarseCmp),
2702 m_coarseTopLevel(coarseGridFct.dof_distribution()->grid_level().level()),
2703 m_spMG(fineGridFct.domain()->grid()),
2704 m_spWeight(spWeight)
2705 {
2707 UG_THROW("H1EnergyDistIntegrand: fine and top level inverted.");
2708
2709 if(m_fineData.domain().get() != m_coarseData.domain().get())
2710 UG_THROW("H1EnergyDistIntegrand: grid functions defined on different domains.");
2711 }
2713
2715 virtual void set_subset(int si)
2716 {
2717 if(!m_fineData.is_def_in_subset(si))
2718 UG_THROW("H1EnergyDistIntegrand: Grid function component"
2719 <<m_fineData.fct()<<" not defined on subset "<<si);
2720 if(!m_coarseData.is_def_in_subset(si))
2721 UG_THROW("H1EnergyDistIntegrand: Grid function component"
2722 <<m_coarseData.fct()<<" not defined on subset "<<si);
2724 }
2725
2727 template <int elemDim>
2728 void evaluate(number vValue[],
2729 const MathVector<worldDim> vFineGlobIP[],
2730 GridObject* pFineElem,
2731 const MathVector<worldDim> vCornerCoords[],
2732 const MathVector<elemDim> vFineLocIP[],
2734 const size_t numIP)
2735 {
2736 typedef typename TGridFunction::template dim_traits<elemDim>::grid_base_object Element;
2737
2738 const TGridFunction &fineGridFct = m_fineData.grid_function();
2739 const TGridFunction &coarseGridFct = m_coarseData.grid_function();
2740
2741 // get coarse element
2742 GridObject* pCoarseElem = pFineElem;
2744 int parentLevel = m_spMG->get_level(pCoarseElem);
2745 while(parentLevel > m_coarseTopLevel){
2746 pCoarseElem = m_spMG->get_parent(pCoarseElem);
2747 parentLevel = m_spMG->get_level(pCoarseElem);
2748 }
2749 }
2750
2751 // get reference object id (i.e. Triangle, Quadrilateral, Tetrahedron, ...)
2752 const ReferenceObjectID fineROID = pFineElem->reference_object_id();
2753 const ReferenceObjectID coarseROID = pCoarseElem->reference_object_id();
2754
2755 // get corner coordinates
2756 std::vector<MathVector<worldDim> > vCornerCoarse;
2757 CollectCornerCoordinates(vCornerCoarse, *static_cast<Element*>(pCoarseElem), *m_coarseData.domain());
2758
2759 // get reference Mapping
2761 = ReferenceMappingProvider::get<elemDim, worldDim>(coarseROID, vCornerCoarse);
2763 = ReferenceMappingProvider::get<elemDim, worldDim>(fineROID, vCornerCoords);
2764
2765 std::vector<MathVector<elemDim> > vCoarseLocIP;
2766 vCoarseLocIP.resize(numIP);
2767 for(size_t ip = 0; ip < vCoarseLocIP.size(); ++ip) VecSet(vCoarseLocIP[ip], 0.0);
2768 map.global_to_local(&vCoarseLocIP[0], vFineGlobIP, numIP);
2769
2770
2771 // determine weights
2772 std::vector<typename weight_type::data_type> elemWeights(numIP, MathMatrix<worldDim, worldDim>());
2773 UG_ASSERT(m_spWeight.valid(), "H1SemiDistIntegrand::evaluate requires valid weights!");
2774
2775 // get local solution if needed
2777 {
2778
2779 LocalIndices ind;
2780 LocalVector u;
2781 fineGridFct.indices(pFineElem, ind); // get global indices
2782 u.resize(ind); // adapt local algebra
2783 GetLocalVector(u, fineGridFct); // read local values of u
2784
2785 // compute data
2786 try{
2787 (*m_spWeight)(&elemWeights[0], vFineGlobIP, 0.0, IIntegrand<number, worldDim>::subset(), pFineElem,
2788 vCornerCoords, vFineLocIP, numIP, &u, NULL);
2789 } UG_CATCH_THROW("H1SemiDistIntegrand: Cannot evaluate data.");
2790 }
2791 else
2792 {
2793 // compute data
2794 try{
2795 (*m_spWeight)(&elemWeights[0], vFineGlobIP, 0.0, IIntegrand<number, worldDim>::subset(), numIP);
2796 } UG_CATCH_THROW("H1SemiDistIntegrand: Cannot evaluate data.");
2797 }
2798
2799 try{
2800
2801
2802 // get trial space
2803 const LocalShapeFunctionSet<elemDim>& rFineLSFS =
2804 LocalFiniteElementProvider::get<elemDim>(fineROID, m_fineData.id());
2805 const LocalShapeFunctionSet<elemDim>& rCoarseLSFS =
2806 LocalFiniteElementProvider::get<elemDim>(coarseROID, m_coarseData.id());
2807
2808 // get multiindices of element
2809 std::vector<DoFIndex> vFineMI, vCoarseMI;
2810 m_fineData.dof_indices(pFineElem, vFineMI);
2811 m_coarseData.dof_indices(pCoarseElem, vCoarseMI);
2812
2813 std::vector<MathVector<elemDim> > vFineLocGradient(vFineMI.size());
2814 std::vector<MathVector<elemDim> > vCoarseLocGradient(vCoarseMI.size());
2815
2816 // loop all integration points
2817 for(size_t ip = 0; ip < numIP; ++ip)
2818 {
2819 // compute shape gradients at ip
2820 rFineLSFS.grads(&vFineLocGradient[0], vFineLocIP[ip]);
2821 rCoarseLSFS.grads(&vCoarseLocGradient[0], vCoarseLocIP[ip]);
2822
2823 // compute approximated solutions at integration point
2824 number fineSolIP = 0.0;
2825 MathVector<elemDim> fineLocTmp(0.0);
2826 for(size_t sh = 0; sh < vFineMI.size(); ++sh)
2827 {
2828 const number val = DoFRef(fineGridFct, vFineMI[sh]);
2829 fineSolIP += val * rFineLSFS.shape(sh, vFineLocIP[ip]);
2830 VecScaleAppend(fineLocTmp, val, vFineLocGradient[sh]);
2831 }
2832
2833 number coarseSolIP = 0.0;
2834 MathVector<elemDim> coarseLocTmp(0.0);
2835 for(size_t sh = 0; sh < vCoarseMI.size(); ++sh)
2836 {
2837 const number val = DoFRef(coarseGridFct, vCoarseMI[sh]);
2838 coarseSolIP += val * rCoarseLSFS.shape(sh, vCoarseLocIP[ip]);
2839 VecScaleAppend(coarseLocTmp, val, vCoarseLocGradient[sh]);
2840 }
2841
2842 // compute global D*gradient
2843 MathVector<worldDim> fineGradIP;
2845 mapF.jacobian_transposed_inverse(fineJTInv, vFineLocIP[ip]);//Inverse(fineJTInv, vJT[ip]);
2846
2847#ifdef UG_DEBUG
2848 CheckGeneralizedInverse<worldDim,elemDim>(vJT[ip], fineJTInv);
2849#endif //UG_DEBUG
2850
2851
2852 MatVecMult(fineGradIP, fineJTInv, fineLocTmp);
2853 MathVector<worldDim> fineWorldLocTmp(0.0);
2854 MatVecMult(fineWorldLocTmp, elemWeights[ip], fineGradIP);
2855
2856 // compute global D*gradient
2857 MathVector<worldDim> coarseGradIP;
2859 map.jacobian_transposed_inverse(coarseJTInv, vCoarseLocIP[ip]);
2860 MatVecMult(coarseGradIP, coarseJTInv, coarseLocTmp);
2861 MathVector<worldDim> coarseWorldLocTmp(0.0);
2862 MatVecMult(coarseWorldLocTmp, elemWeights[ip], coarseGradIP);
2863
2864 // get squared difference
2865 vValue[ip] = VecDistanceSq(fineWorldLocTmp, coarseWorldLocTmp);
2866 //vValue[ip] = VecDistanceSq(fineGradIP, coarseGradIP, elemWeights[ip]);
2867 }
2868
2869 }
2870 UG_CATCH_THROW("H1EnergyDiffIntegrand::evaluate: trial space missing.");
2871 }
2872};
2873
2874
2876template <typename TGridFunction>
2877number H1EnergyDistance2(TGridFunction& spGridFct1, const char* cmp1,
2878 TGridFunction& spGridFct2, const char* cmp2,
2879 int quadOrder, const char* subsets,
2881{
2882 return GridFunctionDistance2<H1EnergyDistIntegrand<TGridFunction>, TGridFunction>
2883 (spGridFct1, cmp1, spGridFct2, cmp2, quadOrder, subsets, weights);
2884}
2885
2887template <typename TGridFunction>
2888number H1EnergyDistance(TGridFunction& spGridFct1, const char* cmp1,
2889 TGridFunction& spGridFct2, const char* cmp2,
2890 int quadOrder, const char* subsets,
2892{ return sqrt(H1EnergyDistance2(spGridFct1, cmp1, spGridFct2, cmp2, quadOrder, subsets, weights)); }
2893
2895template <typename TGridFunction>
2896number H1EnergyDistance2(TGridFunction& spGridFct1, const char* cmp1,
2897 TGridFunction& spGridFct2, const char* cmp2,
2898 int quadOrder, ConstSmartPtr<typename H1SemiDistIntegrand<TGridFunction>::weight_type > weights)
2899{ return H1EnergyDistance2(spGridFct1, cmp1, spGridFct2, cmp2, quadOrder, NULL, weights); }
2900
2902template <typename TGridFunction>
2903number H1EnergyDistance(TGridFunction& spGridFct1, const char* cmp1,
2904 TGridFunction& spGridFct2, const char* cmp2,
2905 int quadOrder, ConstSmartPtr<typename H1SemiDistIntegrand<TGridFunction>::weight_type > weights)
2906{ return sqrt(H1EnergyDistance2(spGridFct1, cmp1, spGridFct2, cmp2, quadOrder, NULL, weights)); }
2907
2908
2909
2910
2912// H1 norm integrand
2914template <typename TGridFunction>
2916 : public StdIntegrand<number, TGridFunction::dim, H1NormIntegrand<TGridFunction> >
2917{
2918 public:
2920 static const int worldDim = TGridFunction::dim;
2921
2922 private:
2924
2925 public:
2927 H1NormIntegrand(TGridFunction& gridFct, size_t cmp)
2928 : m_scalarData(gridFct, cmp) {}
2929
2931 virtual ~H1NormIntegrand() {}
2932
2934 virtual void set_subset(int si)
2935 {
2936 if(!m_scalarData.is_def_in_subset(si))
2937 UG_THROW("H1Norm: Grid function component"
2938 <<m_scalarData.fct()<<" not defined on subset "<<si);
2940 }
2941
2943 template <int elemDim>
2944 void evaluate(number vValue[],
2945 const MathVector<worldDim> vGlobIP[],
2946 GridObject* pElem,
2947 const MathVector<worldDim> vCornerCoords[],
2948 const MathVector<elemDim> vLocIP[],
2950 const size_t numIP)
2951 {
2952 // get reference object id (i.e. Triangle, Quadrilateral, Tetrahedron, ...)
2953 const ReferenceObjectID roid = pElem->reference_object_id();
2954
2956 = ReferenceMappingProvider::get<elemDim, worldDim>(roid, vCornerCoords);
2957
2958 try{
2959 // get trial space
2960 const LocalShapeFunctionSet<elemDim>& rTrialSpace =
2961 LocalFiniteElementProvider::get<elemDim>(roid, m_scalarData.id());
2962
2963 // number of dofs on element
2964 const size_t num_sh = rTrialSpace.num_sh();
2965
2966 // get multiindices of element
2967 std::vector<DoFIndex> ind; // aux. index array
2968 m_scalarData.dof_indices(pElem, ind);
2969
2970 // check multi indices
2971 if(ind.size() != num_sh)
2972 UG_THROW("H1ErrorIntegrand::evaluate: Wrong number of"
2973 " multi indices.");
2974
2975 // loop all integration points
2976 std::vector<MathVector<elemDim> > vLocGradient(num_sh);
2977 for(size_t ip = 0; ip < numIP; ++ip)
2978 {
2979
2980 // compute shape gradients at ip
2981 rTrialSpace.grads(&vLocGradient[0], vLocIP[ip]);
2982
2983 // compute approximated solution at integration point
2984 number approxSolIP = 0.0;
2985 MathVector<elemDim> locTmp; VecSet(locTmp, 0.0);
2986 for(size_t sh = 0; sh < num_sh; ++sh)
2987 {
2988 // get value at shape point (e.g. corner for P1 fct)
2989 const number valSH = DoFRef(m_scalarData.grid_function(), ind[sh]);
2990
2991 // add shape fct at ip * value at shape
2992 approxSolIP += valSH * rTrialSpace.shape(sh, vLocIP[ip]);
2993
2994 // add gradient at ip
2995 VecScaleAppend(locTmp, valSH, vLocGradient[sh]);
2996 }
2997
2998 // compute global gradient
2999 MathVector<worldDim> approxGradIP;
3001 map.jacobian_transposed_inverse(JTInv, vLocIP[ip]);//RightInverse(JTInv, vJT[ip]);
3002
3003#ifdef UG_DEBUG
3004 CheckGeneralizedInverse<worldDim,elemDim>(vJT[ip], JTInv);
3005#endif //UG_DEBUG
3006
3007 MatVecMult(approxGradIP, JTInv, locTmp);
3008
3009 // get squared of difference
3010 vValue[ip] = approxSolIP * approxSolIP;
3011 vValue[ip] += VecDot(approxGradIP, approxGradIP);
3012 }
3013
3014 }
3015 UG_CATCH_THROW("H1SemiNormFuncIntegrand::evaluate: trial space missing.");
3016 }
3017};
3018
3019
3020
3021
3022template <typename TGridFunction>
3023number H1Norm2(TGridFunction& u, const char* cmp,
3024 int quadOrder, const char* subsets=NULL)
3025{
3026// get function id of name
3027 const size_t fct = u.fct_id_by_name(cmp);
3028
3029// check that function exists
3030 if(fct >= u.num_fct())
3031 UG_THROW("H1Norm: Function space does not contain"
3032 " a function with name " << cmp << ".");
3033
3034 H1NormIntegrand<TGridFunction> spIntegrand(u, fct);
3035 return IntegrateSubsets(spIntegrand, u, subsets, quadOrder);
3036}
3037
3038
3039template <typename TGridFunction>
3040number H1Norm(TGridFunction& u, const char* cmp,
3041 int quadOrder, const char* subsets=NULL)
3042{
3043 return sqrt(H1Norm2(u, cmp, quadOrder,subsets));
3044}
3045
3046template <typename TGridFunction>
3047number H1Norm(SmartPtr<TGridFunction> spGridFct, const char* cmp,
3048 int quadOrder, const char* subsets)
3049{
3050 return H1Norm(*spGridFct, cmp, quadOrder, subsets);
3051}
3052
3053template <typename TGridFunction>
3054number H1Norm(SmartPtr<TGridFunction> spGridFct, const char* cmp, int quadOrder)
3055{
3056 return H1Norm(*spGridFct, cmp, quadOrder, NULL);
3057}
3058
3059
3060
3062template <typename TGridFunction>
3064 : public StdIntegrand<number, TGridFunction::dim, H1DistIntegrand<TGridFunction> >
3065{
3066 public:
3068 static const int worldDim = TGridFunction::dim;
3069
3070 private:
3073
3076
3079
3080 public:
3082 H1DistIntegrand(TGridFunction& fineGridFct, size_t fineCmp,
3083 TGridFunction& coarseGridFct, size_t coarseCmp)
3084 : m_fineData(fineGridFct, fineCmp),
3085 m_fineTopLevel(fineGridFct.dof_distribution()->grid_level().level()),
3086 m_coarseData(coarseGridFct, coarseCmp),
3087 m_coarseTopLevel(coarseGridFct.dof_distribution()->grid_level().level()),
3088 m_spMG(fineGridFct.domain()->grid())
3089 {
3091 UG_THROW("H1DiffIntegrand: fine and top level inverted.");
3092
3093 if(fineGridFct.domain().get() !=
3094 coarseGridFct.domain().get())
3095 UG_THROW("H1DiffIntegrand: grid functions defined on different domains.");
3096 };
3097
3099 virtual ~H1DistIntegrand() {};
3100
3102 virtual void set_subset(int si)
3103 {
3104 if(!m_fineData.is_def_in_subset(si))
3105 UG_THROW("H1DiffIntegrand: Grid function component"
3106 <<m_fineData.fct()<<" not defined on subset "<<si);
3107 if(!m_coarseData.is_def_in_subset(si))
3108 UG_THROW("H1DiffIntegrand: Grid function component"
3109 <<m_coarseData.fct()<<" not defined on subset "<<si);
3111 }
3112
3114 template <int elemDim>
3115 void evaluate(number vValue[],
3116 const MathVector<worldDim> vFineGlobIP[],
3117 GridObject* pFineElem,
3118 const MathVector<worldDim> vCornerCoords[],
3119 const MathVector<elemDim> vFineLocIP[],
3121 const size_t numIP)
3122 {
3123 typedef typename TGridFunction::template dim_traits<elemDim>::grid_base_object Element;
3124
3125 // get coarse element
3126 GridObject* pCoarseElem = pFineElem;
3128 int parentLevel = m_spMG->get_level(pCoarseElem);
3129 while(parentLevel > m_coarseTopLevel){
3130 pCoarseElem = m_spMG->get_parent(pCoarseElem);
3131 parentLevel = m_spMG->get_level(pCoarseElem);
3132 }
3133 }
3134
3135 // get reference object id (i.e. Triangle, Quadrilateral, Tetrahedron, ...)
3136 const ReferenceObjectID fineROID = pFineElem->reference_object_id();
3137 const ReferenceObjectID coarseROID = pCoarseElem->reference_object_id();
3138
3139 // get corner coordinates
3140 std::vector<MathVector<worldDim> > vCornerCoarse;
3141 CollectCornerCoordinates(vCornerCoarse, *static_cast<Element*>(pCoarseElem), *m_coarseData.domain());
3142
3143 // get Reference Mapping
3145 = ReferenceMappingProvider::get<elemDim, worldDim>(coarseROID, vCornerCoarse);
3147 = ReferenceMappingProvider::get<elemDim, worldDim>(fineROID, vCornerCoords);
3148
3149 std::vector<MathVector<elemDim> > vCoarseLocIP;
3150 vCoarseLocIP.resize(numIP);
3151 for(size_t ip = 0; ip < vCoarseLocIP.size(); ++ip) VecSet(vCoarseLocIP[ip], 0.0);
3152 map.global_to_local(&vCoarseLocIP[0], vFineGlobIP, numIP);
3153
3154 try{
3155 // get trial space
3156 const LocalShapeFunctionSet<elemDim>& rFineLSFS =
3157 LocalFiniteElementProvider::get<elemDim>(fineROID, m_fineData.id());
3158 const LocalShapeFunctionSet<elemDim>& rCoarseLSFS =
3159 LocalFiniteElementProvider::get<elemDim>(coarseROID, m_coarseData.id());
3160
3161 // get multiindices of element
3162 std::vector<DoFIndex> vFineMI, vCoarseMI;
3163 m_fineData.dof_indices(pFineElem, vFineMI);
3164 m_coarseData.dof_indices(pCoarseElem, vCoarseMI);
3165
3166 std::vector<MathVector<elemDim> > vFineLocGradient(vFineMI.size());
3167 std::vector<MathVector<elemDim> > vCoarseLocGradient(vCoarseMI.size());
3168
3169 // loop all integration points
3170 for(size_t ip = 0; ip < numIP; ++ip)
3171 {
3172 // compute shape gradients at ip
3173 rFineLSFS.grads(&vFineLocGradient[0], vFineLocIP[ip]);
3174 rCoarseLSFS.grads(&vCoarseLocGradient[0], vCoarseLocIP[ip]);
3175
3176 // compute approximated solution at integration point
3177 number fineSolIP = 0.0;
3178 MathVector<elemDim> fineLocTmp; VecSet(fineLocTmp, 0.0);
3179 for(size_t sh = 0; sh < vFineMI.size(); ++sh)
3180 {
3181 const number val = DoFRef(m_fineData.grid_function(), vFineMI[sh]);
3182 fineSolIP += val * rFineLSFS.shape(sh, vFineLocIP[ip]);
3183 VecScaleAppend(fineLocTmp, val, vFineLocGradient[sh]);
3184 }
3185 number coarseSolIP = 0.0;
3186 MathVector<elemDim> coarseLocTmp; VecSet(coarseLocTmp, 0.0);
3187 for(size_t sh = 0; sh < vCoarseMI.size(); ++sh)
3188 {
3189 const number val = DoFRef(m_coarseData.grid_function(), vCoarseMI[sh]);
3190 coarseSolIP += val * rCoarseLSFS.shape(sh, vCoarseLocIP[ip]);
3191 VecScaleAppend(coarseLocTmp, val, vCoarseLocGradient[sh]);
3192 }
3193
3194 // compute global gradient
3195 MathVector<worldDim> fineGradIP;
3197 mapF.jacobian_transposed_inverse(fineJTInv, vFineLocIP[ip]);//RightInverse(fineJTInv, vJT[ip]);
3198
3199#ifdef UG_DEBUG
3200 CheckGeneralizedInverse<worldDim,elemDim>(vJT[ip], fineJTInv);
3201#endif //UG_DEBUG
3202
3203
3204 MatVecMult(fineGradIP, fineJTInv, fineLocTmp);
3205
3206 // compute global gradient
3207 MathVector<worldDim> coarseGradIP;
3209 map.jacobian_transposed_inverse(coarseJTInv, vCoarseLocIP[ip]);
3210 MatVecMult(coarseGradIP, coarseJTInv, coarseLocTmp);
3211
3212 // get squared of difference
3213 vValue[ip] = (coarseSolIP - fineSolIP);
3214 vValue[ip] *= vValue[ip];
3215 vValue[ip] += VecDistanceSq(coarseGradIP, fineGradIP);
3216 }
3217
3218 }
3219 UG_CATCH_THROW("H1DiffIntegrand::evaluate: trial space missing.");
3220 }
3221};
3222
3223
3224
3225template <typename TGridFunction>
3226number H1Distance2(TGridFunction& spGridFct1, const char* cmp1,
3227 TGridFunction& spGridFct2, const char* cmp2,
3228 int quadOrder, const char* subsets=NULL)
3229{
3230 return GridFunctionDistance2<H1DistIntegrand<TGridFunction>, TGridFunction>
3231 (spGridFct1, cmp1, spGridFct2, cmp2, quadOrder, subsets);
3232}
3233
3234
3235template <typename TGridFunction>
3236number H1Distance(TGridFunction& spGridFct1, const char* cmp1,
3237 TGridFunction& spGridFct2, const char* cmp2,
3238 int quadOrder, const char* subsets=NULL)
3239{
3240 return sqrt(H1Distance2(spGridFct1, cmp1, spGridFct2, cmp2, quadOrder, subsets));
3241}
3243template <typename TGridFunction>
3244number H1Error(SmartPtr<TGridFunction> spGridFct1, const char* cmp1,
3245 SmartPtr<TGridFunction> spGridFct2, const char* cmp2,
3246 int quadOrder, const char* subsets)
3247{
3248 return H1Distance(*spGridFct1, cmp1, *spGridFct2, cmp2, quadOrder, subsets);
3249}
3250
3252template <typename TGridFunction>
3253number H1Error(SmartPtr<TGridFunction> spGridFct1, const char* cmp1,
3254 SmartPtr<TGridFunction> spGridFct2, const char* cmp2,
3255 int quadOrder)
3256{
3257 return H1Distance(*spGridFct1, cmp1, *spGridFct2, cmp2, quadOrder, NULL);
3258}
3259
3260
3262// Standard Integrand
3264
3265template <typename TGridFunction>
3267 : public StdIntegrand<number, TGridFunction::dim, StdFuncIntegrand<TGridFunction> >
3268{
3269 public:
3270 // world dimension of grid function
3271 static const int worldDim = TGridFunction::dim;
3272
3273 private:
3274 // grid function
3275 TGridFunction* m_pGridFct;
3276
3277 // component of function
3278 const size_t m_fct;
3279
3280 public:
3282 StdFuncIntegrand(TGridFunction* pGridFct, size_t cmp)
3283 : m_pGridFct(pGridFct), m_fct(cmp)
3284 {};
3285
3287
3289 virtual void set_subset(int si)
3290 {
3291 if(!m_pGridFct->is_def_in_subset(m_fct, si))
3292 UG_THROW("L2ErrorIntegrand: Grid function component"
3293 <<m_fct<<" not defined on subset "<<si);
3295 }
3296
3298 template <int elemDim>
3299 void evaluate(number vValue[],
3300 const MathVector<worldDim> vGlobIP[],
3301 GridObject* pElem,
3302 const MathVector<worldDim> vCornerCoords[],
3303 const MathVector<elemDim> vLocIP[],
3305 const size_t numIP)
3306 {
3307 // get reference object id (i.e. Triangle, Quadrilateral, Tetrahedron, ...)
3309
3310 const LFEID m_id = m_pGridFct->local_finite_element_id(m_fct);
3311
3312 try{
3313 // get trial space
3314 const LocalShapeFunctionSet<elemDim>& rTrialSpace =
3315 LocalFiniteElementProvider::get<elemDim>(roid, m_id);
3316
3317 // number of dofs on element
3318 const size_t num_sh = rTrialSpace.num_sh();
3319
3320 // get multiindices of element
3321
3322 std::vector<DoFIndex> ind; // aux. index array
3323 m_pGridFct->dof_indices(pElem, m_fct, ind);
3324
3325 // check multi indices
3326 UG_COND_THROW(ind.size() != num_sh,
3327 "StdFuncIntegrand::evaluate: Wrong number of multi indices.");
3328
3329 // loop all integration points
3330 for(size_t ip = 0; ip < numIP; ++ip)
3331 {
3332
3333 // compute approximated solution at integration point
3334 number approxSolIP = 0.0;
3335 for(size_t sh = 0; sh < num_sh; ++sh)
3336 {
3337 // get value at shape point (e.g. corner for P1 fct)
3338 // and add shape fct at ip * value at shape
3339 const number valSH = DoFRef((*m_pGridFct), ind[sh]);
3340 approxSolIP += valSH * rTrialSpace.shape(sh, vLocIP[ip]);
3341 }
3342
3343 // get function value at ip
3344 vValue[ip] = approxSolIP;
3345
3346 }
3347
3348 }
3349 UG_CATCH_THROW("StdFuncIntegrand::evaluate: trial space missing.");
3350 }
3351};
3352
3353
3354template <typename TGridFunction>
3355number StdFuncIntegralOnVertex(TGridFunction& gridFct,
3356 size_t fct,
3357 int si)
3358{
3359// integrate elements of subset
3360 typedef typename TGridFunction::template dim_traits<0>::grid_base_object grid_base_object;
3361 typedef typename TGridFunction::template dim_traits<0>::const_iterator const_iterator;
3362
3363 // reset the result
3364 number integral = 0;
3365
3366// note: this iterator is for the base elements, e.g. Face and not
3367// for the special type, e.g. Triangle, Quadrilateral
3368 const_iterator iter = gridFct.template begin<grid_base_object>(si);
3369 const_iterator iterEnd = gridFct.template end<grid_base_object>(si);
3370
3371// iterate over all elements
3372 for(; iter != iterEnd; ++iter)
3373 {
3374 // get element
3375 grid_base_object* pElem = *iter;
3376
3377 std::vector<DoFIndex> ind; // aux. index array
3378 gridFct.dof_indices(pElem, fct, ind);
3379
3380 // compute approximated solution at integration point
3381 number value = 0.0;
3382 for(size_t sh = 0; sh < ind.size(); ++sh)
3383 {
3384 value += DoFRef(gridFct, ind[sh]);
3385 }
3386
3387 // add to global sum
3388 integral += value;
3389
3390 } // end elem
3391
3392// return the summed integral contributions of all elements
3393 return integral;
3394}
3395
3396template <typename TGridFunction>
3398{ return StdFuncIntegralOnVertex(*spGridFct, fct, si); }
3399
3400
3401template <typename TGridFunction>
3402number Integral(TGridFunction& gridFct, const char* cmp,
3403 const char* subsets, int quadOrder)
3404{
3405// get function id of name
3406 const size_t fct = gridFct.fct_id_by_name(cmp);
3407
3408// check that function exists
3409 UG_COND_THROW(fct >= gridFct.num_fct(), "L2Norm: Function space does not contain"
3410 " a function with name " << cmp << ".");
3411
3412// read subsets
3413 SubsetGroup ssGrp(gridFct.domain()->subset_handler());
3414 if(subsets != NULL)
3415 {
3416 ssGrp.add(TokenizeString(subsets));
3417 if(!SameDimensionsInAllSubsets(ssGrp))
3418 UG_THROW("IntegrateSubsets: Subsets '"<<subsets<<"' do not have same dimension."
3419 "Can not integrate on subsets of different dimensions.");
3420 }
3421 else
3422 {
3423 // add all subsets and remove lower dim subsets afterwards
3424 ssGrp.add_all();
3425 RemoveLowerDimSubsets(ssGrp);
3426 }
3427
3428 // \TODO: This should be generalite in IntegrateSubset instead of doing it directly here
3429 bool bOnlyVertex = true;
3430 for(size_t s = 0; s < ssGrp.size(); ++s)
3431 if(ssGrp.dim(s) != 0) bOnlyVertex = false;
3432
3433 if(bOnlyVertex)
3434 {
3435 number value = 0;
3436 for(size_t s = 0; s < ssGrp.size(); ++s)
3437 value += StdFuncIntegralOnVertex(gridFct, fct, ssGrp[s]);
3438
3439#ifdef UG_PARALLEL
3440 // sum over processes
3441 if(pcl::NumProcs() > 1)
3442 {
3444 number local = value;
3445 com.allreduce(&local, &value, 1, PCL_DT_DOUBLE, PCL_RO_SUM);
3446 }
3447#endif
3448 return value;
3449 }
3450
3451 StdFuncIntegrand<TGridFunction> integrand(&gridFct, fct);
3452 return IntegrateSubsets(integrand, gridFct, subsets, quadOrder);
3453}
3454
3455template <typename TGridFunction>
3456number Integral(SmartPtr<TGridFunction> spGridFct, const char* cmp,
3457 const char* subsets, int quadOrder)
3458{ return Integral(*spGridFct, cmp, subsets, quadOrder); }
3459
3460
3461template <typename TGridFunction>
3462number Integral(SmartPtr<TGridFunction> spGridFct, const char* cmp,
3463 const char* subsets)
3464{
3465 return Integral(spGridFct, cmp, subsets, 1);
3466}
3467template <typename TGridFunction>
3468number Integral(SmartPtr<TGridFunction> spGridFct, const char* cmp)
3469{
3470 return Integral(spGridFct, cmp, NULL, 1);
3471}
3472
3473
3476// Generic Boundary Integration Routine
3479
3480template <int WorldDim, int dim, typename TConstIterator>
3482 TConstIterator iterEnd,
3484 const ISubsetHandler* ish,
3485 IIntegrand<MathVector<WorldDim>, WorldDim>& integrand,
3486 const SubsetGroup& bndSSGrp)
3487{
3488// reset the result
3489 number integral = 0;
3490
3491// note: this iterator is for the base elements, e.g. Face and not
3492// for the special type, e.g. Triangle, Quadrilateral
3493 TConstIterator iter = iterBegin;
3494
3495// this is the base element type (e.g. Face). This is the type when the
3496// iterators above are dereferenciated.
3497 typedef typename domain_traits<dim>::grid_base_object grid_base_object;
3498
3499// create a FV1 Geometry
3501
3502// specify, which subsets are boundary
3503 for(size_t s = 0; s < bndSSGrp.size(); ++s)
3504 {
3505 // get subset index
3506 const int bndSubset = bndSSGrp[s];
3507
3508 // request this subset index as boundary subset. This will force the
3509 // creation of boundary subsets when calling geo.update
3510 geo.add_boundary_subset(bndSubset);
3511 }
3512
3513// vector of corner coordinates of element corners (to be filled for each elem)
3514 std::vector<MathVector<WorldDim> > vCorner;
3515
3516// iterate over all elements
3517 for(; iter != iterEnd; ++iter)
3518 {
3519 // get element
3520 grid_base_object* pElem = *iter;
3521
3522 // get all corner coordinates
3523 CollectCornerCoordinates(vCorner, *pElem, aaPos, true);
3524
3525 // compute bf and grads at bip for element
3526 try{
3527 geo.update(pElem, &vCorner[0], ish);
3528 }
3529 UG_CATCH_THROW("IntegralNormalComponentOnManifold: "
3530 "Cannot update Finite Volume Geometry.");
3531
3532 // specify, which subsets are boundary
3533 for(size_t s = 0; s < bndSSGrp.size(); ++s)
3534 {
3535 // get subset index
3536 const int bndSubset = bndSSGrp[s];
3537
3538 // get all bf of this subset
3539 typedef typename DimFV1Geometry<dim>::BF BF;
3540 const std::vector<BF>& vBF = geo.bf(bndSubset);
3541
3542 // loop boundary faces
3543 for(size_t b = 0; b < vBF.size(); ++b)
3544 {
3545 // get bf
3546 const BF& bf = vBF[b];
3547
3548 // value
3550
3551 // jacobian
3553 Inverse(JT, bf.JTInv());
3554
3555 // compute integrand values at integration points
3556 try
3557 {
3558 integrand.values(&value, &(bf.global_ip()),
3559 pElem, &vCorner[0], &(bf.local_ip()),
3560 &(JT),
3561 1);
3562 }
3563 UG_CATCH_THROW("IntegralNormalComponentOnManifold: Unable to compute values of "
3564 "integrand at integration point.");
3565
3566 // sum up contribution (normal includes area)
3567 integral += VecDot(value, bf.normal());
3568
3569 } // end bf
3570 } // end bnd subsets
3571 } // end elem
3572
3573// return the summed integral contributions of all elements
3574 return integral;
3575}
3576
3577
3578template <int WorldDim, int dim, typename TConstIterator>
3580 TConstIterator iterBegin,
3581 TConstIterator iterEnd,
3583 const ISubsetHandler* ish,
3584 IIntegrand<MathVector<WorldDim>, WorldDim>& integrand,
3585 const SubsetGroup& bndSSGrp,
3586 int quadOrder,
3587 Grid& grid)
3588{
3589// reset the result
3590 number integral = 0;
3591
3592// note: this iterator is for the base elements, e.g. Face and not
3593// for the special type, e.g. Triangle, Quadrilateral
3594 TConstIterator iter = iterBegin;
3595
3596// this is the base element type (e.g. Face). This is the type when the
3597// iterators above are dereferenciated.
3598 typedef typename domain_traits<dim>::element_type Element;
3599 typedef typename domain_traits<dim>::side_type Side;
3600
3601// vector of corner coordinates of element corners (to be filled for each elem)
3602 std::vector<MathVector<WorldDim> > vCorner;
3603 std::vector<int> vSubsetIndex;
3604
3605// iterate over all elements
3606 for(; iter != iterEnd; ++iter)
3607 {
3608 // get element
3609 Element* pElem = *iter;
3610
3611 // get all corner coordinates
3612 CollectCornerCoordinates(vCorner, *pElem, aaPos, true);
3613
3614 // get reference object id
3615 const ReferenceObjectID elemRoid = pElem->reference_object_id();
3616
3617 // get sides
3619 grid.associated_elements_sorted(vSide, pElem);
3620 vSubsetIndex.resize(vSide.size());
3621 for(size_t i = 0; i < vSide.size(); ++i)
3622 vSubsetIndex[i] = ish->get_subset_index(vSide[i]);
3623
3625 = ReferenceMappingProvider::get<dim, WorldDim>(elemRoid, vCorner);
3626
3627 const DimReferenceElement<dim>& rRefElem
3628 = ReferenceElementProvider::get<dim>(elemRoid);
3629
3630 // loop sub elements
3631 for(size_t side = 0; side < vSide.size(); ++side)
3632 {
3633 // check if side used
3634 if(!bndSSGrp.contains(vSubsetIndex[side])) continue;
3635
3636 // get side
3637 Side* pSide = vSide[side];
3638
3639 std::vector<MathVector<WorldDim> > vSideCorner(rRefElem.num(dim-1, side, 0));
3640 std::vector<MathVector<dim> > vLocalSideCorner(rRefElem.num(dim-1, side, 0));
3641 for(size_t co = 0; co < vSideCorner.size(); ++co){
3642 vSideCorner[co] = vCorner[rRefElem.id(dim-1, side, 0, co)];
3643 vLocalSideCorner[co] = rRefElem.corner(rRefElem.id(dim-1, side, 0, co));
3644 }
3645
3646 // side quad rule
3647 const ReferenceObjectID sideRoid = pSide->reference_object_id();
3648 const QuadratureRule<dim-1>& rSideQuadRule
3649 = QuadratureRuleProvider<dim-1>::get(sideRoid, quadOrder);
3650
3651 // normal
3652 MathVector<WorldDim> Normal;
3653 ElementNormal<WorldDim>(sideRoid, Normal, &vSideCorner[0]);
3654
3655 // quadrature points
3656 const number* vWeight = rSideQuadRule.weights();
3657 const size_t nip = rSideQuadRule.size();
3658 std::vector<MathVector<dim> > vLocalIP(nip);
3659 std::vector<MathVector<dim> > vGlobalIP(nip);
3660
3661 DimReferenceMapping<dim-1, dim>& map
3662 = ReferenceMappingProvider::get<dim-1, dim>(sideRoid, vLocalSideCorner);
3663
3664 for(size_t ip = 0; ip < nip; ++ip)
3665 map.local_to_global(vLocalIP[ip], rSideQuadRule.point(ip));
3666
3667 for(size_t ip = 0; ip < nip; ++ip)
3668 rMapping.local_to_global(vGlobalIP[ip], vLocalIP[ip]);
3669
3670 // compute transformation matrices
3671 std::vector<MathMatrix<dim-1, WorldDim> > vJT(nip);
3672 map.jacobian_transposed(&(vJT[0]), rSideQuadRule.points(), nip);
3673
3674 std::vector<MathMatrix<dim, WorldDim> > vElemJT(nip);
3675 rMapping.jacobian_transposed(&(vElemJT[0]), &vLocalIP[0], nip);
3676
3677 std::vector<MathVector<WorldDim> > vValue(nip);
3678
3679 // compute integrand values at integration points
3680 try
3681 {
3682 integrand.values(&vValue[0], &vGlobalIP[0],
3683 pElem, &vCorner[0], &vLocalIP[0],
3684 &(vElemJT[0]),
3685 nip);
3686 }
3687 UG_CATCH_THROW("IntegralNormalComponentOnManifold: Unable to compute values of "
3688 "integrand at integration point.");
3689
3690 // loop integration points
3691 for(size_t ip = 0; ip < nip; ++ip)
3692 {
3693 // get quadrature weight
3694 const number weightIP = vWeight[ip];
3695
3696 // get determinate of mapping
3697 const number det = SqrtGramDeterminant(vJT[ip]);
3698
3699 // add contribution of integration point
3700 integral += VecDot(vValue[ip], Normal) * weightIP * det;
3701 }
3702 } // end bf
3703 } // end elem
3704
3705// return the summed integral contributions of all elements
3706 return integral;
3707}
3708
3709template <typename TGridFunction, int dim>
3711 SmartPtr<TGridFunction> spGridFct,
3712 int si, const SubsetGroup& bndSSGrp, int quadOrder)
3713{
3714// integrate elements of subset
3715 typedef typename TGridFunction::template dim_traits<dim>::grid_base_object grid_base_object;
3716 typedef typename TGridFunction::template dim_traits<dim>::const_iterator const_iterator;
3717 static const int WorldDim = TGridFunction::dim;
3718
3719 spIntegrand->set_subset(si);
3720
3721 if(quadOrder == 1)
3722 return IntegralNormalComponentOnManifoldUsingFV1Geom<WorldDim,dim,const_iterator>
3723 (spGridFct->template begin<grid_base_object>(si),
3724 spGridFct->template end<grid_base_object>(si),
3725 spGridFct->domain()->position_accessor(),
3726 spGridFct->domain()->subset_handler().get(),
3727 *spIntegrand, bndSSGrp);
3728 else{
3729 UG_LOG(" #### IntegralNormalComponentOnManifoldSubset ####:\n")
3730 return IntegralNormalComponentOnManifoldGeneral<WorldDim,dim,const_iterator>
3731 (spGridFct->template begin<grid_base_object>(si),
3732 spGridFct->template end<grid_base_object>(si),
3733 spGridFct->domain()->position_accessor(),
3734 spGridFct->domain()->subset_handler().get(),
3735 *spIntegrand, bndSSGrp, quadOrder, *spGridFct->domain()->grid());
3736 }
3737}
3738
3739template <typename TGridFunction>
3741 SmartPtr<IIntegrand<MathVector<TGridFunction::dim>, TGridFunction::dim> > spIntegrand,
3742 SmartPtr<TGridFunction> spGridFct,
3743 const char* BndSubsets, const char* InnerSubsets,
3744 int quadOrder)
3745{
3746// world dimensions
3747 static const int dim = TGridFunction::dim;
3748
3749// read subsets
3750 SubsetGroup innerSSGrp(spGridFct->domain()->subset_handler());
3751 if(InnerSubsets != NULL)
3752 {
3753 innerSSGrp.add(TokenizeString(InnerSubsets));
3754 if(!SameDimensionsInAllSubsets(innerSSGrp))
3755 UG_THROW("IntegralNormalComponentOnManifold: Subsets '"<<InnerSubsets<<"' do not have same dimension."
3756 "Can not integrate on subsets of different dimensions.");
3757 }
3758 else
3759 {
3760 // add all subsets and remove lower dim subsets afterwards
3761 innerSSGrp.add_all();
3762 RemoveLowerDimSubsets(innerSSGrp);
3763 }
3764
3765// read subsets
3766 SubsetGroup bndSSGrp(spGridFct->domain()->subset_handler());
3767 if(BndSubsets != NULL)
3768 bndSSGrp.add(TokenizeString(BndSubsets));
3769 else
3770 UG_THROW("IntegralNormalComponentOnManifold: No boundary subsets passed.");
3771
3772// reset value
3773 number value = 0;
3774
3775// loop subsets
3776 for(size_t i = 0; i < innerSSGrp.size(); ++i)
3777 {
3778 // get subset index
3779 const int si = innerSSGrp[i];
3780
3781 // check dimension
3782 if(innerSSGrp.dim(i) != dim && innerSSGrp.dim(i) != DIM_SUBSET_EMPTY_GRID)
3783 UG_THROW("IntegralNormalComponentOnManifold: Dimension of inner subset is "<<
3784 innerSSGrp.dim(i)<<", but only World Dimension "<<dim<<
3785 " subsets can be used for inner subsets.");
3786
3787 // integrate elements of subset
3788 switch(innerSSGrp.dim(i))
3789 {
3790 case DIM_SUBSET_EMPTY_GRID: break;
3791 case dim: value += IntegralNormalComponentOnManifoldSubset<TGridFunction, dim>(spIntegrand, spGridFct, si, bndSSGrp, quadOrder); break;
3792 default: UG_THROW("IntegralNormalComponentOnManifold: Dimension "<<innerSSGrp.dim(i)<<" not supported. "
3793 " World dimension is "<<dim<<".");
3794 }
3795 }
3796
3797#ifdef UG_PARALLEL
3798 // sum over processes
3799 if(pcl::NumProcs() > 1)
3800 {
3802 number local = value;
3803 com.allreduce(&local, &value, 1, PCL_DT_DOUBLE, PCL_RO_SUM);
3804 }
3805#endif
3806
3807// return the result
3808 return value;
3809}
3810
3811template <typename TGridFunction>
3813 SmartPtr<UserData<MathVector<TGridFunction::dim>, TGridFunction::dim> > spData,
3814 SmartPtr<TGridFunction> spGridFct,
3815 const char* BndSubset, const char* InnerSubset,
3816 number time,
3817 int quadOrder)
3818{
3819 SmartPtr<IIntegrand<MathVector<TGridFunction::dim>, TGridFunction::dim> > spIntegrand
3820 = make_sp(new UserDataIntegrand<MathVector<TGridFunction::dim>, TGridFunction>(spData, &(*spGridFct), time));
3821
3822 return IntegralNormalComponentOnManifoldSubsets(spIntegrand, spGridFct, BndSubset, InnerSubset, quadOrder);
3823}
3824
3825template <typename TGridFunction>
3827 SmartPtr<UserData<MathVector<TGridFunction::dim>, TGridFunction::dim> > spData,
3828 SmartPtr<TGridFunction> spGridFct,
3829 const char* BndSubset, const char* InnerSubset,
3830 number time)
3831{return IntegralNormalComponentOnManifold(spData, spGridFct, BndSubset, InnerSubset, time, 1);}
3832
3833template <typename TGridFunction>
3835 SmartPtr<UserData<MathVector<TGridFunction::dim>, TGridFunction::dim> > spData,
3836 SmartPtr<TGridFunction> spGridFct,
3837 const char* BndSubset,
3838 number time)
3839{return IntegralNormalComponentOnManifold(spData, spGridFct, BndSubset, NULL, time, 1);}
3840
3841template <typename TGridFunction>
3843 SmartPtr<UserData<MathVector<TGridFunction::dim>, TGridFunction::dim> > spData,
3844 SmartPtr<TGridFunction> spGridFct,
3845 const char* BndSubset, const char* InnerSubset)
3846{return IntegralNormalComponentOnManifold(spData, spGridFct, BndSubset, InnerSubset, 0.0, 1);}
3847
3848template <typename TGridFunction>
3850 SmartPtr<UserData<MathVector<TGridFunction::dim>, TGridFunction::dim> > spData,
3851 SmartPtr<TGridFunction> spGridFct,
3852 const char* BndSubset)
3853{return IntegralNormalComponentOnManifold(spData, spGridFct, BndSubset, NULL, 0.0, 1);}
3854
3856// lua data
3858
3859#ifdef UG_FOR_LUA
3860template <typename TGridFunction>
3862 const char* luaFct, SmartPtr<TGridFunction> spGridFct,
3863 const char* BndSubset, const char* InnerSubset,
3864 number time, int quadOrder)
3865{
3866 static const int dim = TGridFunction::dim;
3868 LuaUserDataFactory<MathVector<dim>, dim>::create(luaFct);
3869 return IntegralNormalComponentOnManifold(sp, spGridFct, BndSubset, InnerSubset, time, quadOrder);
3870}
3871
3872template <typename TGridFunction>
3874 const char* luaFct, SmartPtr<TGridFunction> spGridFct,
3875 const char* BndSubset, const char* InnerSubset,
3876 number time)
3877{return IntegralNormalComponentOnManifold(luaFct, spGridFct, BndSubset, InnerSubset, time, 1);}
3878
3879template <typename TGridFunction>
3881 const char* luaFct, SmartPtr<TGridFunction> spGridFct,
3882 const char* BndSubset,
3883 number time)
3884{return IntegralNormalComponentOnManifold(luaFct, spGridFct, BndSubset, NULL, time, 1);}
3885
3886template <typename TGridFunction>
3888 const char* luaFct, SmartPtr<TGridFunction> spGridFct,
3889 const char* BndSubset, const char* InnerSubset)
3890{return IntegralNormalComponentOnManifold(luaFct, spGridFct, BndSubset, InnerSubset, 0.0, 1);}
3891
3892template <typename TGridFunction>
3894 const char* luaFct, SmartPtr<TGridFunction> spGridFct,
3895 const char* BndSubset)
3896{return IntegralNormalComponentOnManifold(luaFct, spGridFct, BndSubset, NULL, 0.0, 1);}
3897#endif
3898
3900// Boundary Integral
3903
3914template <typename TGridFunction>
3915number IntegrateNormalComponentOnManifold(TGridFunction& u, const char* cmp,
3916 const char* BndSubset)
3917{
3918// get function id of name
3919 const size_t fct = u.fct_id_by_name(cmp);
3920
3921// check that function exists
3922 if(fct >= u.num_fct())
3923 UG_THROW("IntegrateNormalComponentOnManifold: Function space does not contain"
3924 " a function with name " << cmp << ".");
3925
3926 if(u.local_finite_element_id(fct) != LFEID(LFEID::LAGRANGE, TGridFunction::dim, 1))
3927 UG_THROW("IntegrateNormalComponentOnManifold:"
3928 "Only implemented for Lagrange P1 functions.");
3929
3930// read bnd subsets
3931 SubsetGroup bndSSGrp(u.domain()->subset_handler());
3932 bndSSGrp.add(TokenizeString(BndSubset));
3933
3934// reset value
3935 number value = 0;
3936
3937// loop subsets
3938 for(size_t i = 0; i < bndSSGrp.size(); ++i)
3939 {
3940 // get subset index
3941 const int si = bndSSGrp[i];
3942
3943 // skip if function is not defined in subset
3944 if(!u.is_def_in_subset(fct, si)) continue;
3945
3946 // create integration kernel
3947 static const int dim = TGridFunction::dim;
3948
3949 // integrate elements of subset
3950 typedef typename domain_traits<dim>::grid_base_object grid_base_object;
3951
3952 // get iterators for all elems on subset
3953 typedef typename TGridFunction::template dim_traits<dim>::const_iterator const_iterator;
3954 const_iterator iter = u.template begin<grid_base_object>();
3955 const_iterator iterEnd = u.template end<grid_base_object>();
3956
3957 // create a FV1 Geometry
3959
3960 // specify, which subsets are boundary
3961 for(size_t s = 0; s < bndSSGrp.size(); ++s)
3962 {
3963 // get subset index
3964 const int bndSubset = bndSSGrp[s];
3965
3966 // request this subset index as boundary subset. This will force the
3967 // creation of boundary subsets when calling geo.update
3968 geo.add_boundary_subset(bndSubset);
3969 }
3970
3971 // vector of corner coordinates of element corners (to be filled for each elem)
3972 std::vector<MathVector<dim> > vCorner;
3973
3974 // loop elements of subset
3975 for( ; iter != iterEnd; ++iter)
3976 {
3977 // get element
3978 grid_base_object* elem = *iter;
3979
3980 // get all corner coordinates
3981 CollectCornerCoordinates(vCorner, *elem, u.domain()->position_accessor(), true);
3982
3983 // compute bf and grads at bip for element
3984 try{
3985 geo.update(elem, &vCorner[0], u.domain()->subset_handler().get());
3986 }
3987 UG_CATCH_THROW("IntegrateNormalComponenOnManifold: "
3988 "Cannot update Finite Volume Geometry.");
3989
3990 // get fct multi-indices of element
3991 std::vector<DoFIndex> ind;
3992 u.dof_indices(elem, fct, ind);
3993
3994 // specify, which subsets are boundary
3995 for(size_t s = 0; s < bndSSGrp.size(); ++s)
3996 {
3997 // get subset index
3998 const int bndSubset = bndSSGrp[s];
3999
4000 // get all bf of this subset
4001 typedef typename DimFV1Geometry<dim>::BF BF;
4002 const std::vector<BF>& vBF = geo.bf(bndSubset);
4003
4004 // loop boundary faces
4005 for(size_t b = 0; b < vBF.size(); ++b)
4006 {
4007 // get bf
4008 const BF& bf = vBF[b];
4009
4010 // get normal on bf
4011 const MathVector<dim>& normal = bf.normal();
4012
4013 // check multi indices
4014 UG_ASSERT(ind.size() == bf.num_sh(),
4015 "IntegrateNormalComponentOnManifold::values: Wrong number of"
4016 " multi indices, ind: "<<ind.size() << ", bf.num_sh: "
4017 << bf.num_sh());
4018
4019 MathVector<dim> val; VecSet(val, 0.0);
4020 for(size_t sh = 0; sh < bf.num_sh(); ++sh)
4021 {
4022 const number fctVal = DoFRef(u, ind[sh]);
4023 const MathVector<dim>& ip = bf.global_ip();
4024 VecScaleAdd(val, 1.0, val, fctVal, ip);
4025 }
4026
4027 // sum up contributions
4028 value -= VecDot(val, normal);
4029 }
4030 }
4031 }
4032 }
4033
4034#ifdef UG_PARALLEL
4035 // sum over processes
4036 if(pcl::NumProcs() > 1)
4037 {
4039 number local = value;
4040 com.allreduce(&local, &value, 1, PCL_DT_DOUBLE, PCL_RO_SUM);
4041 }
4042#endif
4043
4044// return the result
4045 return value;
4046}
4047
4049
4063template <typename TGridFunction>
4064number IntegrateNormalGradientOnManifold(TGridFunction& u, const char* cmp,
4065 const char* BndSubset, const char* InnerSubset = NULL)
4066{
4067// get function id of name
4068 const size_t fct = u.fct_id_by_name(cmp);
4069
4070// check that function exists
4071 if(fct >= u.num_fct())
4072 UG_THROW("IntegrateNormalGradientOnManifold: Function space does not contain"
4073 " a function with name " << cmp << ".");
4074
4075 if(u.local_finite_element_id(fct) != LFEID(LFEID::LAGRANGE, TGridFunction::dim, 1))
4076 UG_THROW("IntegrateNormalGradientOnManifold:"
4077 "Only implemented for Lagrange P1 functions.");
4078
4079// read subsets
4080 SubsetGroup innerSSGrp(u.domain()->subset_handler());
4081 if(InnerSubset != NULL)
4082 innerSSGrp.add(TokenizeString(InnerSubset));
4083 else // add all if no subset specified
4084 innerSSGrp.add_all();
4085
4086// read bnd subsets
4087 SubsetGroup bndSSGrp(u.domain()->subset_handler());
4088 if(InnerSubset != NULL){
4089 bndSSGrp.add(TokenizeString(BndSubset));
4090 }
4091 else{
4092 UG_THROW("IntegrateNormalGradientOnManifold: No boundary subsets specified. Aborting.");
4093 }
4094
4095// reset value
4096 number value = 0;
4097
4098// loop subsets
4099 for(size_t i = 0; i < innerSSGrp.size(); ++i)
4100 {
4101 // get subset index
4102 const int si = innerSSGrp[i];
4103
4104 // skip if function is not defined in subset
4105 if(!u.is_def_in_subset(fct, si)) continue;
4106
4107
4108 if (innerSSGrp.dim(i) != TGridFunction::dim)
4109 UG_THROW("IntegrateNormalGradientOnManifold: Element dimension does not match world dimension!");
4110
4111 // create integration kernel
4112 static const int dim = TGridFunction::dim;
4113
4114 // integrate elements of subset
4115 typedef typename domain_traits<dim>::grid_base_object grid_base_object;
4116
4117 // get iterators for all elems on subset
4118 typedef typename TGridFunction::template dim_traits<dim>::const_iterator const_iterator;
4119 const_iterator iter = u.template begin<grid_base_object>();
4120 const_iterator iterEnd = u.template end<grid_base_object>();
4121
4122 // create a FV1 Geometry
4124
4125 // specify, which subsets are boundary
4126 for(size_t s = 0; s < bndSSGrp.size(); ++s)
4127 {
4128 // get subset index
4129 const int bndSubset = bndSSGrp[s];
4130
4131 // request this subset index as boundary subset. This will force the
4132 // creation of boundary subsets when calling geo.update
4133 geo.add_boundary_subset(bndSubset);
4134 }
4135
4136 // vector of corner coordinates of element corners (to be filled for each elem)
4137 std::vector<MathVector<dim> > vCorner;
4138
4139 // loop elements of subset
4140 for( ; iter != iterEnd; ++iter)
4141 {
4142 // get element
4143 grid_base_object* elem = *iter;
4144
4145 // get all corner coordinates
4146 CollectCornerCoordinates(vCorner, *elem, u.domain()->position_accessor(), true);
4147
4148 // compute bf and grads at bip for element
4149 try{
4150 geo.update(elem, &vCorner[0], u.domain()->subset_handler().get());
4151 }
4152 UG_CATCH_THROW("IntegrateNormalGradientOnManifold: "
4153 "Cannot update Finite Volume Geometry.");
4154
4155 // get fct multi-indeces of element
4156 std::vector<DoFIndex> ind;
4157 u.dof_indices(elem, fct, ind);
4158
4159 // specify, which subsets are boundary
4160 for(size_t s = 0; s < bndSSGrp.size(); ++s)
4161 {
4162 // get subset index
4163 const int bndSubset = bndSSGrp[s];
4164
4165 // get all bf of this subset
4166 typedef typename DimFV1Geometry<dim>::BF BF;
4167 const std::vector<BF>& vBF = geo.bf(bndSubset);
4168
4169 // loop boundary faces
4170 for(size_t b = 0; b < vBF.size(); ++b)
4171 {
4172 // get bf
4173 const BF& bf = vBF[b];
4174
4175 // get normal on bf
4176 const MathVector<dim>& normal = bf.normal();
4177
4178 // check multi indices
4179 UG_ASSERT(ind.size() == bf.num_sh(),
4180 "IntegrateNormalGradientOnManifold::values: Wrong number of"
4181 " multi indices, ind: "<<ind.size() << ", bf.num_sh: "
4182 << bf.num_sh());
4183
4184 // compute gradient of solution at bip
4185 MathVector<dim> grad; VecSet(grad, 0.0);
4186 for(size_t sh = 0; sh < bf.num_sh(); ++sh)
4187 {
4188 const number fctVal = DoFRef(u, ind[sh]);
4189
4190 VecScaleAdd(grad, 1.0, grad, fctVal, bf.global_grad(sh));
4191 }
4192
4193 // sum up contributions
4194 value -= VecDot(grad, normal);
4195 }
4196 }
4197 }
4198 }
4199
4200#ifdef UG_PARALLEL
4201 // sum over processes
4202 if(pcl::NumProcs() > 1)
4203 {
4205 number local = value;
4206 com.allreduce(&local, &value, 1, PCL_DT_DOUBLE, PCL_RO_SUM);
4207 }
4208#endif
4209
4210// return the result
4211 return value;
4212}
4213
4214} // namespace ug
4215
4216#endif /*__H__UG__LIB_DISC__FUNCTION_SPACES__INTEGRATE__*/
parameterString s
Definition smart_pointer.h:296
bool valid() const
returns true if the pointer is valid, false if not.
Definition smart_pointer.h:414
Definition smart_pointer.h:108
T * get()
returns encapsulated pointer
Definition smart_pointer.h:197
Definition pcl_process_communicator.h:70
void allreduce(const void *sendBuf, void *recBuf, int count, DataType type, ReduceOperation op) const
performs MPI_Allreduce on the processes of the communicator.
Definition pcl_process_communicator.cpp:318
bool valid() const
Definition attachment_pipe.h:553
constant matrix user data
Definition const_user_data.h:232
constant scalar user data
Definition const_user_data.h:153
boundary face
Definition fv1_geom.h:903
Geometry and shape functions for 1st order Vertex-Centered Finite Volume.
Definition fv1_geom.h:665
const BF & bf(int si, size_t i) const
returns the boundary face i for subsetIndex
Definition fv1_geom.h:1135
void update(GridObject *elem, const MathVector< worldDim > *vCornerCoords, const ISubsetHandler *ish=NULL)
update data for given element
Definition fv1_geom.cpp:669
void add_boundary_subset(int subsetIndex)
add subset that is interpreted as boundary subset.
Definition fv1_geom.h:1104
dimension dependent base class for reference elements
Definition reference_element.h:183
const MathVector< dim > & corner(size_t i) const
coordinates of reference corner (i = 0 ... num(0))
Definition reference_element.h:192
virtual base class for reference mappings
Definition reference_mapping_provider.h:53
virtual void jacobian_transposed(MathMatrix< dim, worldDim > &JT, const MathVector< dim > &locPos) const =0
returns transposed of jacobian
virtual number jacobian_transposed_inverse(MathMatrix< worldDim, dim > &JTInv, const MathVector< dim > &locPos) const =0
returns transposed of the inverse of the jacobian and returns sqrt of gram determinante
virtual void local_to_global(MathVector< worldDim > &globPos, const MathVector< dim > &locPos) const =0
map local coordinate to global coordinate
virtual void update(const MathVector< worldDim > *vCornerCoord)=0
refresh mapping for new set of corners
virtual void global_to_local(MathVector< dim > &locPos, const MathVector< worldDim > &globPos, const size_t maxIter=1000, const number tol=1e-10) const =0
map global coordinate to local coordinate
the generic attachment-accessor for access to grids attachment pipes.
Definition grid.h:182
Manages the elements of a grid and their interconnection.
Definition grid.h:132
void associated_elements_sorted(traits< Vertex >::secure_container &elemsOut, TElem *e)
Puts all elements of type TAss which are contained in 'e' into elemsOut in the reference elements ord...
Definition grid_impl.hpp:503
The base class for all geometric objects, such as vertices, edges, faces, volumes,...
Definition grid_base_objects.h:157
virtual ReferenceObjectID reference_object_id() const =0
Integrand for the distance of two grid functions - evaluated in the H1 norm.
Definition integrate.h:3065
const int m_fineTopLevel
Definition integrate.h:3072
const int m_coarseTopLevel
Definition integrate.h:3075
ScalarGridFunctionData< TGridFunction > m_fineData
Definition integrate.h:3071
SmartPtr< MultiGrid > m_spMG
multigrid
Definition integrate.h:3078
virtual ~H1DistIntegrand()
DTOR.
Definition integrate.h:3099
ScalarGridFunctionData< TGridFunction > m_coarseData
Definition integrate.h:3074
static const int worldDim
world dimension of grid function
Definition integrate.h:3068
virtual void set_subset(int si)
sets subset
Definition integrate.h:3102
H1DistIntegrand(TGridFunction &fineGridFct, size_t fineCmp, TGridFunction &coarseGridFct, size_t coarseCmp)
constructor (1 is fine grid function)
Definition integrate.h:3082
void evaluate(number vValue[], const MathVector< worldDim > vFineGlobIP[], GridObject *pFineElem, const MathVector< worldDim > vCornerCoords[], const MathVector< elemDim > vFineLocIP[], const MathMatrix< elemDim, worldDim > vJT[], const size_t numIP)
returns the values of the integrand for a bunch of ips
Definition integrate.h:3115
Integrand for the distance of two grid functions - evaluated in the norm .
Definition integrate.h:2658
ConstSmartPtr< weight_type > m_spWeight
scalar weight (optional)
Definition integrate.h:2675
static const int worldDim
world dimension of grid function
Definition integrate.h:2661
SmartPtr< MultiGrid > m_spMG
multigrid
Definition integrate.h:2672
H1SemiIntegrand< TGridFunction >::weight_type weight_type
Definition integrate.h:2662
virtual void set_subset(int si)
sets subset
Definition integrate.h:2715
ScalarGridFunctionData< TGridFunction > m_fineData
Definition integrate.h:2665
virtual ~H1EnergyDistIntegrand()
Definition integrate.h:2712
void evaluate(number vValue[], const MathVector< worldDim > vFineGlobIP[], GridObject *pFineElem, const MathVector< worldDim > vCornerCoords[], const MathVector< elemDim > vFineLocIP[], const MathMatrix< elemDim, worldDim > vJT[], const size_t numIP)
returns the values of the integrand for a bunch of ips
Definition integrate.h:2728
const int m_fineTopLevel
Definition integrate.h:2666
ScalarGridFunctionData< TGridFunction > m_coarseData
Definition integrate.h:2668
H1EnergyDistIntegrand(TGridFunction &fineGridFct, size_t fineCmp, TGridFunction &coarseGridFct, size_t coarseCmp, ConstSmartPtr< weight_type > spWeight)
constructor
Definition integrate.h:2696
H1EnergyDistIntegrand(TGridFunction &fineGridFct, size_t fineCmp, TGridFunction &coarseGridFct, size_t coarseCmp)
constructor
Definition integrate.h:2679
const int m_coarseTopLevel
Definition integrate.h:2669
Norm of a grid function, evaluated in (weighted) H1-semi norm.
Definition integrate.h:2454
ConstSmartPtr< weight_type > m_spWeight
scalar weight (optional)
Definition integrate.h:2465
virtual void set_subset(int si)
sets subset
Definition integrate.h:2480
static const int worldDim
world dimension of grid function
Definition integrate.h:2457
UserData< MathMatrix< worldDim, worldDim >, worldDim > weight_type
Definition integrate.h:2458
ScalarGridFunctionData< TGridFunction > m_scalarData
grid function data
Definition integrate.h:2462
virtual ~H1EnergyIntegrand()
DTOR.
Definition integrate.h:2477
H1EnergyIntegrand(TGridFunction &gridFct, size_t cmp, ConstSmartPtr< weight_type > spWeight)
constructor
Definition integrate.h:2473
H1EnergyIntegrand(TGridFunction &gridFct, size_t cmp)
constructor
Definition integrate.h:2469
void evaluate(number vValue[], const MathVector< worldDim > vGlobIP[], GridObject *pElem, const MathVector< worldDim > vCornerCoords[], const MathVector< elemDim > vLocIP[], const MathMatrix< elemDim, worldDim > vJT[], const size_t numIP)
returns the values of the integrand for a bunch of ips
Definition integrate.h:2490
Definition integrate.h:1262
virtual void set_subset(int si)
sets subset
Definition integrate.h:1293
SmartPtr< UserData< MathVector< worldDim >, worldDim > > m_spExactGrad
exact gradient
Definition integrate.h:1275
H1ErrorIntegrand(SmartPtr< UserData< number, worldDim > > spExactSol, SmartPtr< UserData< MathVector< worldDim >, worldDim > > spExactGrad, TGridFunction &gridFct, size_t cmp, number time)
constructor
Definition integrate.h:1282
ScalarGridFunctionData< TGridFunction > m_scalarData
grid function
Definition integrate.h:1269
SmartPtr< UserData< number, worldDim > > m_spExactSolution
exact solution
Definition integrate.h:1272
static const int worldDim
world dimension of grid function
Definition integrate.h:1265
void evaluate(number vValue[], const MathVector< worldDim > vGlobIP[], GridObject *pElem, const MathVector< worldDim > vCornerCoords[], const MathVector< elemDim > vLocIP[], const MathMatrix< elemDim, worldDim > vJT[], const size_t numIP)
returns the values of the integrand for a bunch of ips
Definition integrate.h:1303
number m_time
time
Definition integrate.h:1278
Definition integrate.h:2917
void evaluate(number vValue[], const MathVector< worldDim > vGlobIP[], GridObject *pElem, const MathVector< worldDim > vCornerCoords[], const MathVector< elemDim > vLocIP[], const MathMatrix< elemDim, worldDim > vJT[], const size_t numIP)
returns the values of the integrand for a bunch of ips
Definition integrate.h:2944
ScalarGridFunctionData< TGridFunction > m_scalarData
Definition integrate.h:2923
virtual ~H1NormIntegrand()
DTOR.
Definition integrate.h:2931
virtual void set_subset(int si)
sets subset
Definition integrate.h:2934
static const int worldDim
world dimension of grid function
Definition integrate.h:2920
H1NormIntegrand(TGridFunction &gridFct, size_t cmp)
CTOR.
Definition integrate.h:2927
Integrand for the distance of two grid functions - evaluated in the (weighted) H1-semi norm.
Definition integrate.h:2173
void evaluate(number vValue[], const MathVector< worldDim > vFineGlobIP[], GridObject *pFineElem, const MathVector< worldDim > vCornerCoords[], const MathVector< elemDim > vFineLocIP[], const MathMatrix< elemDim, worldDim > vJT[], const size_t numIP)
returns the values of the integrand for a bunch of ips
Definition integrate.h:2239
H1SemiDistIntegrand(TGridFunction &fineGridFct, size_t fineCmp, TGridFunction &coarseGridFct, size_t coarseCmp, ConstSmartPtr< weight_type > spWeight)
constructor
Definition integrate.h:2211
const int m_fineTopLevel
Definition integrate.h:2181
H1SemiDistIntegrand(TGridFunction &fineGridFct, size_t fineCmp, TGridFunction &coarseGridFct, size_t coarseCmp)
constructor
Definition integrate.h:2194
ScalarGridFunctionData< TGridFunction > m_fineData
Definition integrate.h:2180
virtual ~H1SemiDistIntegrand()
Definition integrate.h:2224
ConstSmartPtr< weight_type > m_spWeight
scalar weight (optional)
Definition integrate.h:2190
H1SemiIntegrand< TGridFunction >::weight_type weight_type
Definition integrate.h:2177
virtual void set_subset(int si)
sets subset
Definition integrate.h:2227
static const int worldDim
world dimension of grid function
Definition integrate.h:2176
SmartPtr< MultiGrid > m_spMG
multigrid
Definition integrate.h:2187
ScalarGridFunctionData< TGridFunction > m_coarseData
Definition integrate.h:2183
const int m_coarseTopLevel
Definition integrate.h:2184
Norm of a grid function, evaluated in (weighted) H1-semi norm.
Definition integrate.h:1967
UserData< MathMatrix< worldDim, worldDim >, worldDim > weight_type
Definition integrate.h:1971
void evaluate(number vValue[], const MathVector< worldDim > vGlobIP[], GridObject *pElem, const MathVector< worldDim > vCornerCoords[], const MathVector< elemDim > vLocIP[], const MathMatrix< elemDim, worldDim > vJT[], const size_t numIP)
returns the values of the integrand for a bunch of ips
Definition integrate.h:2003
static const int worldDim
world dimension of grid function
Definition integrate.h:1970
H1SemiIntegrand(TGridFunction &gridFct, size_t cmp, ConstSmartPtr< weight_type > spWeight)
constructor
Definition integrate.h:1986
ScalarGridFunctionData< TGridFunction > m_scalarData
grid function data
Definition integrate.h:1975
ConstSmartPtr< weight_type > m_spWeight
scalar weight (optional)
Definition integrate.h:1978
virtual void set_subset(int si)
sets subset
Definition integrate.h:1993
H1SemiIntegrand(TGridFunction &gridFct, size_t cmp)
constructor
Definition integrate.h:1982
virtual ~H1SemiIntegrand()
DTOR.
Definition integrate.h:1990
Abstract integrand interface.
Definition integrate.h:121
virtual void values(TData vValue[], const MathVector< worldDim > vGlobIP[], GridObject *pElem, const MathVector< worldDim > vCornerCoords[], const MathVector< 1 > vLocIP[], const MathMatrix< 1, worldDim > vJT[], const size_t numIP)=0
returns the values of the integrand for a bunch of ips
virtual void set_subset(int si)
sets the subset
Definition integrate.h:168
static const int worldDim
world dimension
Definition integrate.h:124
int subset() const
returns the subset
Definition integrate.h:171
virtual void values(TData vValue[], const MathVector< worldDim > vGlobIP[], GridObject *pElem, const MathVector< worldDim > vCornerCoords[], const MathVector< 3 > vLocIP[], const MathMatrix< 3, worldDim > vJT[], const size_t numIP)=0
returns the values of the integrand for a bunch of ips
int m_si
subset
Definition integrate.h:175
TData data_type
data type
Definition integrate.h:127
virtual void values(TData vValue[], const MathVector< worldDim > vGlobIP[], GridObject *pElem, const MathVector< worldDim > vCornerCoords[], const MathVector< 2 > vLocIP[], const MathMatrix< 2, worldDim > vJT[], const size_t numIP)=0
returns the values of the integrand for a bunch of ips
virtual ~IIntegrand()
Definition integrate.h:164
Definition subset_handler_interface.h:223
int get_subset_index(GridObject *elem) const
Definition subset_handler_interface.cpp:560
Integrand for the distance of two grid functions - evaluated in the (weighted) H1-semi norm.
Definition integrate.h:1741
const int m_coarseTopLevel
Definition integrate.h:1752
void evaluate(number vValue[], const MathVector< worldDim > vFineGlobIP[], GridObject *pFineElem, const MathVector< worldDim > vCornerCoords[], const MathVector< elemDim > vFineLocIP[], const MathMatrix< elemDim, worldDim > vJT[], const size_t numIP)
returns the values of the integrand for a bunch of ips
Definition integrate.h:1821
virtual ~L2DistIntegrand()
Definition integrate.h:1807
double m_deltaFineCoarse
shift
Definition integrate.h:1760
static const int worldDim
world dimension of grid function
Definition integrate.h:1744
L2DistIntegrand(TGridFunction &fineGridFct, size_t fineCmp, TGridFunction &coarseGridFct, size_t coarseCmp, ConstSmartPtr< weight_type > spWeight)
constructor (1st is fine grid function)
Definition integrate.h:1779
ConstSmartPtr< weight_type > m_spWeight
Definition integrate.h:1757
const int m_fineTopLevel
Definition integrate.h:1749
L2DistIntegrand(TGridFunction &fineGridFct, size_t fineCmp, TGridFunction &coarseGridFct, size_t coarseCmp)
constructor (1st is fine grid function)
Definition integrate.h:1765
ScalarGridFunctionData< TGridFunction > m_coarseData
Definition integrate.h:1751
L2DistIntegrand(TGridFunction &fineGridFct, size_t fineCmp, TGridFunction &coarseGridFct, size_t coarseCmp, ConstSmartPtr< weight_type > spWeight, number dist12)
constructor (1st is fine grid function)
Definition integrate.h:1793
virtual void set_subset(int si)
sets subset
Definition integrate.h:1810
L2Integrand< TGridFunction >::weight_type weight_type
Definition integrate.h:1745
SmartPtr< MultiGrid > m_spMG
multigrid
Definition integrate.h:1755
ScalarGridFunctionData< TGridFunction > m_fineData
Definition integrate.h:1748
Definition integrate.h:1073
L2ErrorIntegrand(SmartPtr< UserData< number, worldDim > > spExactSol, TGridFunction &gridFct, size_t cmp, number time)
constructor
Definition integrate.h:1089
static const int worldDim
world dimension of grid function
Definition integrate.h:1076
number m_time
time
Definition integrate.h:1085
ScalarGridFunctionData< TGridFunction > m_scalarData
Definition integrate.h:1079
SmartPtr< UserData< number, worldDim > > m_spExactSolution
exact solution
Definition integrate.h:1082
virtual void set_subset(int si)
sets subset
Definition integrate.h:1099
virtual ~L2ErrorIntegrand()
Definition integrate.h:1096
void evaluate(number vValue[], const MathVector< worldDim > vGlobIP[], GridObject *pElem, const MathVector< worldDim > vCornerCoords[], const MathVector< elemDim > vLocIP[], const MathMatrix< elemDim, worldDim > vJT[], const size_t numIP)
returns the values of the integrand for a bunch of ips
Definition integrate.h:1109
Grid function as L2 integrand.
Definition integrate.h:1574
ScalarGridFunctionData< TGridFunction > m_scalarData
Definition integrate.h:1582
virtual void set_subset(int si)
sets subset
Definition integrate.h:1601
UserData< number, worldDim > weight_type
Definition integrate.h:1578
L2Integrand(TGridFunction &spGridFct, size_t cmp, ConstSmartPtr< weight_type > spWeight)
Definition integrate.h:1593
virtual ~L2Integrand()
DTOR.
Definition integrate.h:1598
void evaluate(number vValue[], const MathVector< worldDim > vGlobIP[], GridObject *pElem, const MathVector< worldDim > vCornerCoords[], const MathVector< elemDim > vLocIP[], const MathMatrix< elemDim, worldDim > vJT[], const size_t numIP)
returns the values of the integrand for a bunch of ips
Definition integrate.h:1610
static const int worldDim
Definition integrate.h:1577
L2Integrand(TGridFunction &spGridFct, size_t cmp)
CTOR.
Definition integrate.h:1589
ConstSmartPtr< weight_type > m_spWeight
scalar weight (optional, default is 1.0)
Definition integrate.h:1585
Identifier for Local Finite Elements.
Definition local_finite_element_id.h:98
@ LAGRANGE
Definition local_finite_element_id.h:104
virtual size_t num_sh() const
Definition local_dof_set.cpp:46
Definition local_algebra.h:50
virtual base class for local shape function sets
Definition local_shape_function_set.h:70
virtual void grads(grad_type *vGrad, const MathVector< dim > &x) const =0
returns all gradients evaluated at a point
virtual shape_type shape(size_t i, const MathVector< dim > &x) const =0
evaluates the shape function
Definition local_algebra.h:198
void resize(const LocalIndices &ind)
resize for current local indices
Definition local_algebra.h:214
static SmartPtr< LuaUserData< TData, dim, TRet > > create(const std::string &name)
Definition lua_user_data.h:223
A class for fixed size, dense matrices.
Definition math_matrix.h:63
a mathematical Vector with N entries.
Definition math_vector.h:97
Definition integrate.h:965
ScalarGridFunctionData< TGridFunction > m_scalarData
Definition integrate.h:971
double max()
Definition integrate.h:1001
virtual void set_subset(int si)
sets subset
Definition integrate.h:987
static const int worldDim
world dimension of grid function
Definition integrate.h:968
MaximumDistIntegrand(TGridFunction &gridFct, size_t cmp)
constructor
Definition integrate.h:977
void reset()
Definition integrate.h:994
double m_min
Definition integrate.h:973
double min()
Definition integrate.h:1000
void evaluate(number vValue[], const MathVector< worldDim > vGlobIP[], GridObject *pElem, const MathVector< worldDim > vCornerCoords[], const MathVector< elemDim > vLocIP[], const MathMatrix< elemDim, worldDim > vJT[], const size_t numIP)
returns the values of the integrand for a bunch of ips
Definition integrate.h:1005
double m_max
Definition integrate.h:973
virtual ~MaximumDistIntegrand()
Definition integrate.h:984
provides quadrature rule for a Reference Dimension
Definition quadrature.h:70
number weight(size_t i) const
return the i'th weight
Definition quadrature.h:105
size_t size() const
number of integration points
Definition quadrature.h:92
const MathVector< dim > * points() const
returns all positions in an array of size()
Definition quadrature.h:102
static const QuadratureRule< TDim > & get(size_t order, QuadType type=BEST)
gets quadrature rule of requested order
Definition quadrature_provider_impl.h:40
size_t num(int dim) const
returns the number of geometric objects of dim
Definition reference_element.h:95
int id(int dim_i, size_t i, int dim_j, size_t j) const
id of object j in dimension dim_j of obj i in dimension dim_i
Definition reference_element.h:127
static DimReferenceMapping< TDim, TWorldDim > & get(ReferenceObjectID roid)
returns a reference to a DimReferenceMapping
Definition reference_mapping_provider.h:225
Definition integrate.h:66
TGridFunction & grid_function()
Definition integrate.h:75
TGridFunction::domain_type domain_type
Definition integrate.h:68
const TGridFunction & grid_function() const
Definition integrate.h:78
size_t m_fct
component of function
Definition integrate.h:106
size_t fct()
Definition integrate.h:81
TGridFunction & m_gridFct
grid function
Definition integrate.h:103
ConstSmartPtr< domain_type > domain() const
returns const domain (forward)
Definition integrate.h:95
ScalarGridFunctionData(TGridFunction &gridFct, size_t cmp)
Definition integrate.h:70
SmartPtr< domain_type > domain()
returns domain (forward)
Definition integrate.h:92
size_t dof_indices(TElem *elem, std::vector< DoFIndex > &ind, bool bHang=false, bool bClear=true) const
Definition integrate.h:98
LFEID m_id
local finite element id
Definition integrate.h:109
const LFEID & id() const
Definition integrate.h:84
bool is_def_in_subset(int si) const
returns true, iff scalar function is defined in subset si
Definition integrate.h:88
Definition integrate.h:3268
const size_t m_fct
Definition integrate.h:3278
virtual void set_subset(int si)
sets subset
Definition integrate.h:3289
TGridFunction * m_pGridFct
Definition integrate.h:3275
static const int worldDim
Definition integrate.h:3271
virtual ~StdFuncIntegrand()
Definition integrate.h:3286
StdFuncIntegrand(TGridFunction *pGridFct, size_t cmp)
constructor
Definition integrate.h:3282
void evaluate(number vValue[], const MathVector< worldDim > vGlobIP[], GridObject *pElem, const MathVector< worldDim > vCornerCoords[], const MathVector< elemDim > vLocIP[], const MathMatrix< elemDim, worldDim > vJT[], const size_t numIP)
returns the values of the integrand for a bunch of ips
Definition integrate.h:3299
Abstract integrand interface (using CRTP)
Definition integrate.h:181
virtual void values(TData vValue[], const MathVector< worldDim > vGlobIP[], GridObject *pElem, const MathVector< worldDim > vCornerCoords[], const MathVector< 1 > vLocIP[], const MathMatrix< 1, worldDim > vJT[], const size_t numIP)
returns the values of the integrand for a bunch of ips
Definition integrate.h:191
virtual void values(TData vValue[], const MathVector< worldDim > vGlobIP[], GridObject *pElem, const MathVector< worldDim > vCornerCoords[], const MathVector< 2 > vLocIP[], const MathMatrix< 2, worldDim > vJT[], const size_t numIP)
returns the values of the integrand for a bunch of ips
Definition integrate.h:201
TImpl & getImpl()
access to implementation
Definition integrate.h:225
const TImpl & getImpl() const
const access to implementation
Definition integrate.h:228
static const int worldDim
world dimension
Definition integrate.h:184
TData data_type
data type
Definition integrate.h:187
virtual void values(TData vValue[], const MathVector< worldDim > vGlobIP[], GridObject *pElem, const MathVector< worldDim > vCornerCoords[], const MathVector< 3 > vLocIP[], const MathMatrix< 3, worldDim > vJT[], const size_t numIP)
returns the values of the integrand for a bunch of ips
Definition integrate.h:211
Group of subsets.
Definition subset_group.h:51
size_t size() const
number of subsets in this group
Definition subset_group.h:122
void add_all()
select all subsets of underlying subset handler
Definition subset_group.cpp:133
int dim(size_t i) const
dimension of subset
Definition subset_group.cpp:237
void add(int si)
adds a subset by number to this group
Definition subset_group.cpp:64
bool contains(int si) const
returns true if subset is contained in this group
Definition subset_group.cpp:272
For arbitrary UserData and grid functions and , this class (should) define the integrand .
Definition integrate.h:668
ScalarGridFunctionData< TGridFunction > m_coarseData
Definition integrate.h:678
number inner_prod(const MathVector< worldDim > &d1, const MathVector< worldDim > &d2)
Definition integrate.h:827
number inner_dist2(const number &v1, const number &v2)
Definition integrate.h:836
number inner_prod(const number &d1, const number &d2)
Definition integrate.h:824
virtual ~UserDataDistIntegrandSq()
Definition integrate.h:706
ScalarGridFunctionData< TGridFunction > m_fineData
Definition integrate.h:675
const int m_fineTopLevel
Definition integrate.h:676
static const int worldDim
world dimension of grid function
Definition integrate.h:671
SmartPtr< MultiGrid > m_spMG
multigrid
Definition integrate.h:682
number inner_prod(const T &d1, const T &d2)
Definition integrate.h:831
number inner_dist2(const MathVector< worldDim > &v1, const MathVector< worldDim > &v2)
Definition integrate.h:839
UserDataDistIntegrandSq(SmartPtr< UserData< TData, worldDim > > spData, TGridFunction &fineGridFct, size_t fineCmp, TGridFunction &coarseGridFct, size_t coarseCmp)
constructor (1st is fine grid function)
Definition integrate.h:693
double m_time
Definition integrate.h:689
const int m_coarseTopLevel
Definition integrate.h:679
void evaluate(number vValue[], const MathVector< worldDim > vGlobIP[], GridObject *pFineElem, const MathVector< worldDim > vCornerCoords[], const MathVector< elemDim > vFineLocIP[], const MathMatrix< elemDim, worldDim > vJT[], const size_t numIP)
returns the values of the integrand for a bunch of ips
Definition integrate.h:722
SmartPtr< UserData< TData, worldDim > > m_spData
Definition integrate.h:686
virtual void set_subset(int si)
sets subset
Definition integrate.h:709
number inner_dist2(const T &d1, const T &d2)
Definition integrate.h:843
Type based UserData.
Definition user_data.h:143
virtual bool requires_grid_fct() const =0
returns if grid function is needed for evaluation
number data_type
Definition user_data.h:145
For arbitrary UserData , this class defines the integrand .
Definition integrate.h:468
UserDataIntegrand(SmartPtr< UserData< TData, worldDim > > spData, number time)
constructor
Definition integrate.h:497
static const int worldDim
Definition integrate.h:471
TData data_type
Definition integrate.h:474
TGridFunction * m_spGridFct
Definition integrate.h:481
UserDataIntegrand(SmartPtr< UserData< TData, worldDim > > spData, TGridFunction *spGridFct, number time)
constructor
Definition integrate.h:488
number m_time
Definition integrate.h:484
SmartPtr< UserData< TData, worldDim > > m_spData
Definition integrate.h:478
void evaluate(TData vValue[], const MathVector< worldDim > vGlobIP[], GridObject *pElem, const MathVector< worldDim > vCornerCoords[], const MathVector< elemDim > vLocIP[], const MathMatrix< elemDim, worldDim > vJT[], const size_t numIP)
returns the values of the integrand for a bunch of ips
Definition integrate.h:507
For arbitrary UserData (of type TData), this class defines the integrand .
Definition integrate.h:557
number inner_prod(const MathVector< worldDim > &d1, const MathVector< worldDim > &d2)
Definition integrate.h:653
UserDataIntegrandSq(SmartPtr< UserData< TData, worldDim > > spData, const TGridFunction *pGridFct, number time)
constructor
Definition integrate.h:577
SmartPtr< UserData< TData, worldDim > > m_spData
Definition integrate.h:567
void evaluate(number vValue[], const MathVector< worldDim > vGlobIP[], GridObject *pElem, const MathVector< worldDim > vCornerCoords[], const MathVector< elemDim > vLocIP[], const MathMatrix< elemDim, worldDim > vJT[], const size_t numIP)
returns the values of the integrand for a bunch of ips
Definition integrate.h:598
static const int worldDim
Definition integrate.h:560
TData data_type
Definition integrate.h:563
number m_time
Definition integrate.h:573
const TGridFunction * m_pGridFct
Definition integrate.h:570
number inner_prod(const number &d1, const number &d2)
Definition integrate.h:650
number inner_prod(const T &d1, const T &d2)
Definition integrate.h:657
UserDataIntegrandSq(SmartPtr< UserData< TData, worldDim > > spData, number time)
constructor
Definition integrate.h:586
void CollectCornerCoordinates(std::vector< typename TAAPos::ValueType > &vCornerCoordsOut, const TElem &elem, const TAAPos &aaPos, bool clearContainer=true)
returns the corner coordinates of a geometric object
Definition domain_util_impl.h:75
MathMatrix< N, M, T >::value_type SqrtGramDeterminant(const MathMatrix< N, M, T > &m)
Square root of Gram Determinant of a matrix.
Definition math_matrix_functions_common_impl.hpp:508
MathMatrix< N, M, T >::value_type GeneralizedInverse(MathMatrix< N, M, T > &mOut, const MathMatrix< M, N, T > &m)
Definition math_matrix_functions_common_impl.hpp:751
matrix_t::value_type MatMaxNorm(matrix_t &m)
Definition math_matrix_functions_common_impl.hpp:1024
MathMatrix< N, M, T >::value_type Inverse(MathMatrix< N, M, T > &mOut, const MathMatrix< M, N, T > &m)
Inverse of a matrix.
Definition math_matrix_functions_common_impl.hpp:560
void MatSubtract(matrix_t &mOut, const matrix_t &m1, const matrix_t &m2)
subtracts m2 from m1 and stores the result in a mOut
Definition math_matrix_functions_common_impl.hpp:68
#define PCL_RO_SUM
Definition pcl_methods.h:63
#define PCL_DT_DOUBLE
Definition pcl_datatype.h:57
int NumProcs()
returns the number of processes
Definition pcl_base.cpp:91
const NullSmartPtr SPNULL
The equivalent to NULL for smart pointers.
Definition smart_pointer.h:90
#define UG_ASSERT(expr, msg)
Definition assert.h:70
#define UG_CATCH_THROW(msg)
Definition error.h:64
#define UG_THROW(msg)
Definition error.h:57
#define UG_LOG(msg)
Definition log.h:367
#define UG_COND_THROW(cond, msg)
UG_COND_THROW(cond, msg) : performs a UG_THROW(msg) if cond == true.
Definition error.h:61
double number
Definition types.h:124
void MatVecMult(vector_t_out &vOut, const matrix_t &m, const vector_t_in &v)
Matrix - Vector Multiplication.
Definition math_matrix_vector_functions_common_impl.hpp:49
void VecScaleAppend(vector_t &vOut, typename vector_t::value_type s1, const vector_t &v1)
Scales a Vector and adds it to a second vector.
Definition math_vector_functions_common_impl.hpp:126
void VecSet(vector_t &vInOut, typename vector_t::value_type s)
Set each vector component to scalar (componentwise)
Definition math_vector_functions_common_impl.hpp:539
void VecScaleAdd(vector_t &vOut, typename vector_t::value_type s1, const vector_t &v1, typename vector_t::value_type s2, const vector_t &v2)
Scales two Vectors, adds them and returns the sum in a third vector.
Definition math_vector_functions_common_impl.hpp:265
vector_t::value_type VecTwoNormSq(const vector_t &v)
Definition math_vector_functions_common_impl.hpp:585
vector_t::value_type VecDistanceSq(const vector_t &v1, const vector_t &v2)
returns the squared distance of two vector_ts.
Definition math_vector_functions_common_impl.hpp:351
vector_t::value_type VecDot(const vector_t &v1, const vector_t &v2)
returns the dot-product of two vector_ts
Definition math_vector_functions_common_impl.hpp:385
function util LuaCallbackHelper create(func)
Definition smart_pointer.h:814
the ug namespace
number H1SemiDistance2(TGridFunction &spGridFct1, const char *cmp1, TGridFunction &spGridFct2, const char *cmp2, int quadOrder, const char *subsets, ConstSmartPtr< typename H1SemiDistIntegrand< TGridFunction >::weight_type > weights)
Squared distance in H1 semi norm (with select subsets & weights)
Definition integrate.h:2412
number H1SemiError2(SmartPtr< TGridFunction > spGridFct1, const char *cmp1, SmartPtr< TGridFunction > spGridFct2, const char *cmp2, int quadOrder, const char *subsets)
Distance in H1 semi norm (with subset selection)
Definition integrate.h:2384
number IntegrateNormalGradientOnManifold(TGridFunction &u, const char *cmp, const char *BndSubset, const char *InnerSubset=NULL)
Integrates the Flux of a component over some boundary subsets.
Definition integrate.h:4064
number Integrate(TConstIterator iterBegin, TConstIterator iterEnd, typename domain_traits< WorldDim >::position_accessor_type &aaPos, IIntegrand< number, WorldDim > &integrand, int quadOrder, std::string quadType, Grid::AttachmentAccessor< typename domain_traits< dim >::grid_base_object, ANumber > *paaElemContribs=NULL)
integrates on the whole domain
Definition integrate.h:260
number L2Distance2(TGridFunction &spGridFct1, const char *cmp1, TGridFunction &spGridFct2, const char *cmp2, int quadOrder, const char *subsets, ConstSmartPtr< typename L2Integrand< TGridFunction >::weight_type > spWeight, number avgDist12=0.0)
computes the squared l2 distance between two functions
Definition integrate.h:1911
number H1Norm(TGridFunction &u, const char *cmp, int quadOrder, const char *subsets=NULL)
Definition integrate.h:3040
number H1SemiNorm2(TGridFunction &gridFct, const char *cmp, int quadOrder, const char *subsets=NULL, ConstSmartPtr< typename H1SemiIntegrand< TGridFunction >::weight_type > weights=SPNULL)
compute H1 semi-norm of a function on the whole domain (or on selected subsets)
Definition integrate.h:2125
QuadType GetQuadratureType(const std::string &name)
returns Identifier from string
Definition quadrature_provider.cpp:299
ReferenceObjectID
these ids are used to identify the shape of a geometric object.
Definition grid_base_objects.h:74
number IntegralNormalComponentOnManifoldGeneral(TConstIterator iterBegin, TConstIterator iterEnd, typename domain_traits< WorldDim >::position_accessor_type &aaPos, const ISubsetHandler *ish, IIntegrand< MathVector< WorldDim >, WorldDim > &integrand, const SubsetGroup &bndSSGrp, int quadOrder, Grid &grid)
Definition integrate.h:3579
number IntegrateNormalComponentOnManifold(TGridFunction &u, const char *cmp, const char *BndSubset)
Integrates a component over some boundary subsets.
Definition integrate.h:3915
number L2Norm2(TGridFunction &u, const char *cmp, int quadOrder, const char *subsets, ConstSmartPtr< typename L2Integrand< TGridFunction >::weight_type > spWeight)
Definition integrate.h:1680
const number SMALL
Definition math_constants.h:41
number L2Distance(TGridFunction &spGridFct1, const char *cmp1, TGridFunction &spGridFct2, const char *cmp2, int quadOrder, const char *subsets)
computes the l2 distance between two functions
Definition integrate.h:1933
number H1Error(SmartPtr< UserData< number, TGridFunction::dim > > spExactSol, SmartPtr< UserData< MathVector< TGridFunction::dim >, TGridFunction::dim > > spExactGrad, SmartPtr< TGridFunction > spGridFct, const char *cmp, number time, int quadOrder, const char *subsets)
compute H1 error of a function on the whole domain or on some subsets
Definition integrate.h:1400
number IntegrateSubsets(IIntegrand< number, TGridFunction::dim > &spIntegrand, TGridFunction &spGridFct, const char *subsets, int quadOrder, std::string quadType=std::string())
Definition integrate.h:393
number IntegralNormalComponentOnManifoldUsingFV1Geom(TConstIterator iterBegin, TConstIterator iterEnd, typename domain_traits< WorldDim >::position_accessor_type &aaPos, const ISubsetHandler *ish, IIntegrand< MathVector< WorldDim >, WorldDim > &integrand, const SubsetGroup &bndSSGrp)
Definition integrate.h:3481
number H1SemiDistance(TGridFunction &spGridFct1, const char *cmp1, TGridFunction &spGridFct2, const char *cmp2, int quadOrder, const char *subsets, ConstSmartPtr< typename H1SemiDistIntegrand< TGridFunction >::weight_type > weights)
Distance in H1 semi norm (with select subsets & weights)
Definition integrate.h:2423
void TokenizeString(const string &str, vector< string > &vToken, const char delimiter)
Definition string_util.cpp:56
number IntegrateSubset(IIntegrand< number, TGridFunction::dim > &spIntegrand, TGridFunction &spGridFct, int si, int quadOrder, std::string quadType)
Definition integrate.h:373
number H1EnergyNorm(TGridFunction &gridFct, const char *cmp, int quadOrder, const char *subsets=NULL, ConstSmartPtr< typename H1SemiIntegrand< TGridFunction >::weight_type > weights=SPNULL)
Definition integrate.h:2631
number H1Distance2(TGridFunction &spGridFct1, const char *cmp1, TGridFunction &spGridFct2, const char *cmp2, int quadOrder, const char *subsets=NULL)
Definition integrate.h:3226
number IntegralNormalComponentOnManifold(SmartPtr< UserData< MathVector< TGridFunction::dim >, TGridFunction::dim > > spData, SmartPtr< TGridFunction > spGridFct, const char *BndSubset, const char *InnerSubset, number time, int quadOrder)
Definition integrate.h:3812
bool SameDimensionsInAllSubsets(const SubsetGroup &subsetGroup)
Definition subset_group.cpp:299
number Minimum(SmartPtr< TGridFunction > spGridFct, const char *cmp, const char *subsets)
Definition integrate.h:1062
number StdFuncIntegralOnVertex(TGridFunction &gridFct, size_t fct, int si)
Definition integrate.h:3355
number IntegralNormalComponentOnManifoldSubsets(SmartPtr< IIntegrand< MathVector< TGridFunction::dim >, TGridFunction::dim > > spIntegrand, SmartPtr< TGridFunction > spGridFct, const char *BndSubsets, const char *InnerSubsets, int quadOrder)
Definition integrate.h:3740
number H1SemiNorm(TGridFunction &gridFct, const char *cmp, int quadOrder, const char *subsets=NULL, ConstSmartPtr< typename H1SemiIntegrand< TGridFunction >::weight_type > weights=SPNULL)
Computes the H1SemiNorm.
Definition integrate.h:2145
number IntegralNormalComponentOnManifoldSubset(SmartPtr< IIntegrand< MathVector< TGridFunction::dim >, TGridFunction::dim > > spIntegrand, SmartPtr< TGridFunction > spGridFct, int si, const SubsetGroup &bndSSGrp, int quadOrder)
Definition integrate.h:3710
void GetLocalVector(LocalVector &lvec, const TVector &vec)
Definition local_algebra.h:739
void RemoveLowerDimSubsets(SubsetGroup &subsetGroup)
Definition subset_group.cpp:315
number H1Norm2(TGridFunction &u, const char *cmp, int quadOrder, const char *subsets=NULL)
Definition integrate.h:3023
number L2Error(SmartPtr< UserData< number, TGridFunction::dim > > spExactSol, TGridFunction &gridFct, const char *cmp, number time, int quadOrder, const char *subsets)
computes the l2 error function on the whole domain or on some subsets
Definition integrate.h:1185
number H1Distance(TGridFunction &spGridFct1, const char *cmp1, TGridFunction &spGridFct2, const char *cmp2, int quadOrder, const char *subsets=NULL)
Definition integrate.h:3236
number GridFunctionDistance2(TGridFunction &spGridFct1, const char *cmp1, TGridFunction &spGridFct2, const char *cmp2, int quadOrder, const char *subsets)
computes an (abstract) distance between two functions
Definition integrate.h:1466
number H1SemiError(SmartPtr< TGridFunction > spGridFct1, const char *cmp1, SmartPtr< TGridFunction > spGridFct2, const char *cmp2, int quadOrder, const char *subsets)
Distance in H1 semi norm (with subset selection)
Definition integrate.h:2394
number H1EnergyDistance(TGridFunction &spGridFct1, const char *cmp1, TGridFunction &spGridFct2, const char *cmp2, int quadOrder, const char *subsets, ConstSmartPtr< typename H1SemiDistIntegrand< TGridFunction >::weight_type > weights)
Distance in H1 semi norm (with select subsets & weights)
Definition integrate.h:2888
number L2Norm(TGridFunction &u, const char *cmp, int quadOrder, const char *subsets)
Definition integrate.h:1711
number Integral(SmartPtr< UserData< number, TGridFunction::dim > > spData, TGridFunction &spGridFct, const char *subsets, number time, int quadOrder, std::string quadType)
Definition integrate.h:856
number H1EnergyDistance2(TGridFunction &spGridFct1, const char *cmp1, TGridFunction &spGridFct2, const char *cmp2, int quadOrder, const char *subsets, ConstSmartPtr< typename H1SemiDistIntegrand< TGridFunction >::weight_type > weights)
Squared distance in H1 semi norm (with select subsets & weights)
Definition integrate.h:2877
number & DoFRef(TMatrix &mat, const DoFIndex &iInd, const DoFIndex &jInd)
Definition multi_index.h:276
number GetMinimum(TGridFunction &gridFct, const char *cmp, const char *subsets)
Definition integrate.h:1048
number H1EnergyNorm2(TGridFunction &gridFct, const char *cmp, int quadOrder, const char *subsets=NULL, ConstSmartPtr< typename H1SemiIntegrand< TGridFunction >::weight_type > weights=SPNULL)
compute energy -norm of a function on the whole domain (or on selected subsets)
Definition integrate.h:2613
@ DIM_SUBSET_EMPTY_GRID
Definition subset_dim_util.h:92
#define PROFILE_FUNC()
Definition profiler.h:257
QuadType
types of quadratures
Definition quadrature_provider.h:41
SmartPtr< T, FreePolicy > make_sp(T *inst)
returns a SmartPtr for the passed raw pointer
Definition smart_pointer.h:836
PointerConstArray< TElem * > secure_container
Definition grid.h:146
Definition domain_traits.h:53
function ProblemDisc new(problemDesc, dom)