FEAT 3
Finite Element Analysis Toolbox
Loading...
Searching...
No Matches
evaluator.hpp
1// FEAT3: Finite Element Analysis Toolbox, Version 3
2// Copyright (C) 2010 by Stefan Turek & the FEAT group
3// FEAT3 is released under the GNU General Public License version 3,
4// see the file 'copyright.txt' in the top level directory for details.
5
6#pragma once
7
8// includes, FEAT
9#include <kernel/space/parametric_evaluator.hpp>
10
11namespace FEAT
12{
13 namespace Space
14 {
15 namespace CroRavRanTur
16 {
23
29 template<
30 typename Space_,
31 typename TrafoEvaluator_,
32 typename SpaceEvalTraits_,
33 typename Shape_ = typename Space_::ShapeType>
34 class Evaluator DOXY({});
35
41 template<
42 typename Space_,
43 typename TrafoEvaluator_,
44 typename SpaceEvalTraits_>
45 class Evaluator<Space_, TrafoEvaluator_, SpaceEvalTraits_, Shape::Simplex<2> > :
47 Evaluator<
48 Space_,
49 TrafoEvaluator_,
50 SpaceEvalTraits_,
51 Shape::Simplex<2> >,
52 TrafoEvaluator_,
53 SpaceEvalTraits_,
54 ref_caps>
55 {
56 public:
59
61 typedef Space_ SpaceType;
62
64 typedef SpaceEvalTraits_ SpaceEvalTraits;
65
67 typedef typename SpaceEvalTraits::EvalPolicy EvalPolicy;
68
70 typedef typename EvalPolicy::DomainPointType DomainPointType;
71
73 typedef typename SpaceEvalTraits::DataType DataType;
74
75 public:
82 explicit Evaluator(const SpaceType& DOXY(space))
83 {
84 }
85
93 {
94 return 3;
95 }
96
106 template<typename EvalData_>
108 EvalData_& data,
109 const DomainPointType& point) const
110 {
111 data.phi[0].ref_value = DataType(2) * (point[0] + point[1]) - DataType(1);
112 data.phi[1].ref_value = DataType(1) - DataType(2)*point[0];
113 data.phi[2].ref_value = DataType(1) - DataType(2)*point[1];
114 }
115
125 template<typename EvalData_>
127 EvalData_& data,
128 const DomainPointType& DOXY(point)) const
129 {
130 data.phi[0].ref_grad[0] = DataType(2);
131 data.phi[0].ref_grad[1] = DataType(2);
132 data.phi[1].ref_grad[0] = -DataType(2);
133 data.phi[1].ref_grad[1] = DataType(0);
134 data.phi[2].ref_grad[0] = DataType(0);
135 data.phi[2].ref_grad[1] = -DataType(2);
136 }
137 }; // class Evaluator<...,Simplex<2>>
138
144 template<
145 typename Space_,
146 typename TrafoEvaluator_,
147 typename SpaceEvalTraits_>
148 class Evaluator<Space_, TrafoEvaluator_, SpaceEvalTraits_, Shape::Simplex<3> > :
149 public ParametricEvaluator<
150 Evaluator<
151 Space_,
152 TrafoEvaluator_,
153 SpaceEvalTraits_,
154 Shape::Simplex<3> >,
155 TrafoEvaluator_,
156 SpaceEvalTraits_,
157 ref_caps>
158 {
159 public:
162
164 typedef Space_ SpaceType;
165
167 typedef SpaceEvalTraits_ SpaceEvalTraits;
168
170 typedef typename SpaceEvalTraits::EvalPolicy EvalPolicy;
171
173 typedef typename EvalPolicy::DomainPointType DomainPointType;
174
176 typedef typename SpaceEvalTraits::DataType DataType;
177
178 public:
185 explicit Evaluator(const SpaceType& DOXY(space))
186 {
187 }
188
196 {
197 return 4;
198 }
199
209 template<typename EvalData_>
211 EvalData_& data,
212 const DomainPointType& point) const
213 {
214 data.phi[0].ref_value = DataType(3) * (point[0] + point[1] + point[2]) - DataType(2);
215 data.phi[1].ref_value = DataType(1) - DataType(3)*point[0];
216 data.phi[2].ref_value = DataType(1) - DataType(3)*point[1];
217 data.phi[3].ref_value = DataType(1) - DataType(3)*point[2];
218 }
219
229 template<typename EvalData_>
231 EvalData_& data,
232 const DomainPointType& DOXY(point)) const
233 {
234 data.phi[0].ref_grad[0] = DataType(3);
235 data.phi[0].ref_grad[1] = DataType(3);
236 data.phi[0].ref_grad[2] = DataType(3);
237 data.phi[1].ref_grad[0] = -DataType(3);
238 data.phi[1].ref_grad[1] = DataType(0);
239 data.phi[1].ref_grad[2] = DataType(0);
240 data.phi[2].ref_grad[0] = DataType(0);
241 data.phi[2].ref_grad[1] = -DataType(3);
242 data.phi[2].ref_grad[2] = DataType(0);
243 data.phi[3].ref_grad[0] = DataType(0);
244 data.phi[3].ref_grad[1] = DataType(0);
245 data.phi[3].ref_grad[2] = -DataType(3);
246 }
247 }; // class Evaluator<...,Simplex<3>>
248
254 template<
255 typename Space_,
256 typename TrafoEvaluator_,
257 typename SpaceEvalTraits_>
258 class Evaluator<Space_, TrafoEvaluator_, SpaceEvalTraits_, Shape::Quadrilateral> :
259 public EvaluatorBase<
260 Evaluator<Space_, TrafoEvaluator_, SpaceEvalTraits_, Shape::Quadrilateral>,
261 TrafoEvaluator_,
262 SpaceEvalTraits_>
263 {
264 public:
267
269 typedef Space_ SpaceType;
270
272 typedef SpaceEvalTraits_ SpaceEvalTraits;
273
275 typedef TrafoEvaluator_ TrafoEvaluator;
276
278 typedef typename TrafoEvaluator::EvalTraits TrafoEvalTraits;
279
281 typedef typename TrafoEvaluator::TrafoType TrafoType;
282
284 typedef typename TrafoType::MeshType MeshType;
285
287 typedef typename MeshType::template IndexSet<2,1>::Type FacetIndexSetType;
288
290 typedef typename SpaceEvalTraits::DataType DataType;
291
293 typedef typename SpaceEvalTraits::EvalPolicy EvalPolicy;
294
296 typedef typename EvalPolicy::DomainPointType DomainPointType;
298 typedef typename EvalPolicy::ImagePointType ImagePointType;
299
301 typedef typename EvalPolicy::JacobianMatrixType JacobianMatrixType;
303 typedef typename EvalPolicy::JacobianInverseType JacobianInverseType;
304
305 static constexpr SpaceTags eval_caps = SpaceTags::value | SpaceTags::grad;
306
307 template<SpaceTags cfg_>
308 struct ConfigTraits
309 {
311 static constexpr SpaceTags config = cfg_;
312
314 static constexpr TrafoTags trafo_config = TrafoTags::img_point;
315
318 };
319
320 protected:
322 static constexpr TrafoTags inv_lin_trafo_config = TrafoTags::dom_point | TrafoTags::img_point | TrafoTags::jac_inv;
323
325 typedef typename TrafoEvaluator::template ConfigTraits<inv_lin_trafo_config>::EvalDataType InvLinTrafoData;
326
328 typedef typename TrafoType::template Evaluator<Shape::Hypercube<1>, DataType>::Type FacetTrafoEvaluator;
330 typedef typename FacetTrafoEvaluator::EvalTraits FacetEvalTraits;
331
333 static constexpr TrafoTags facet_trafo_config = TrafoTags::img_point | TrafoTags::jac_det;
334
336 typedef typename FacetTrafoEvaluator::template ConfigTraits<facet_trafo_config>::EvalDataType FacetTrafoData;
337
340
343
344 // inverse linearized trafo vector
345 ImagePointType _inv_lin_vec;
346
349
350 void _build_inv_lin_trafo(const TrafoEvaluator& trafo_eval)
351 {
352 // create a domain point in the barycentre of the cell
354
355 // create the trafo data
356 InvLinTrafoData trafo_data;
357 trafo_eval(trafo_data, dom_point);
358
359 // store inverse trafo linearization
360 _inv_lin_mat = trafo_data.jac_inv;
361 _inv_lin_vec = trafo_data.img_point;
362 }
363
365 void _build_coeff_matrix(const TrafoEvaluator& trafo_eval)
366 {
367 // fetch the trafo
368 const TrafoType& trafo = trafo_eval.get_trafo();
369
370 // fetch the mesh
371 const MeshType& mesh = trafo.get_mesh();
372
373 // fetch the facet-index-set of the mesh
374 const FacetIndexSetType& facet_index_set = mesh.template get_index_set<2,1>();
375
376 // fetch the cell index of the currently active cell
377 const Index cell = trafo_eval.get_cell_index();
378
379 // create a nodal matrix
380 CoeffMatrixType _nodal_mat(DataType(0));
381
382 // create a facet trafo evaluator
383 FacetTrafoEvaluator facet_eval(trafo);
384
385 // create facet evaluation data
386 FacetTrafoData facet_data;
387
388 // define 2-point Gauss cubature point coordinate
389 const DataType g = Math::sqrt(DataType(1) / DataType(3));
390 const typename FacetEvalTraits::DomainPointType g1(-g), g2(+g);
391 DomainPointType q1, q2;
392
393 // loop over all 4 edges of the quad
394 for(int i(0); i < 4; ++i)
395 {
396 // initialize facet evaluator
397 facet_eval.prepare(Index(facet_index_set(cell, i)));
398
399 // map first cubature point
400 facet_eval(facet_data, g1);
401 DataType w1(facet_data.jac_det);
402 q1.set_mat_vec_mult(_inv_lin_mat, facet_data.img_point - _inv_lin_vec);
403
404 // map second cubature point
405 facet_eval(facet_data, g2);
406 DataType w2(facet_data.jac_det);
407 q2.set_mat_vec_mult(_inv_lin_mat, facet_data.img_point - _inv_lin_vec);
408
409 // release facet evaluator
410 facet_eval.finish();
411
412 // compute reciprocal edge length
413 DataType v = DataType(1) / (w1 + w2);
414
415 // evaluate node functionals
416 _nodal_mat(0,i) = DataType(1); // = v*(w1 + w2)
417 _nodal_mat(1,i) = v*(w1*q1[0] + w2*q2[0]);
418 _nodal_mat(2,i) = v*(w1*q1[1] + w2*q2[1]);
419 _nodal_mat(3,i) = v*(w1*((q1[0]+q1[1])*(q1[0]-q1[1])) + w2*((q2[0]+q2[1])*(q2[0]-q2[1])));
420 }
421
422 // invert nodal matrix to obtain coefficient matrix
423 _coeff_mat.set_inverse(_nodal_mat);
424 }
425
426 public:
433 explicit Evaluator(const SpaceType& DOXY(space))
434 {
435 }
436
438 virtual ~Evaluator()
439 {
440 }
441
449 {
450 return 4;
451 }
452
459 void prepare(const TrafoEvaluator& trafo_eval)
460 {
461 // prepare inverse linearized trafo
462 _build_inv_lin_trafo(trafo_eval);
463
464 // build coefficient matrix
465 _build_coeff_matrix(trafo_eval);
466 }
467
469 template<SpaceTags space_cfg_, TrafoTags trafo_cfg_>
472 const Trafo::EvalData<TrafoEvalTraits, trafo_cfg_>& trafo_data) const
473 {
474 // transform image point
476 pt.set_mat_vec_mult(_inv_lin_mat, trafo_data.img_point - _inv_lin_vec);
477
478 // "evaluate" monomials
479 DataType x(pt[0]);
480 DataType y(pt[1]);
481 DataType r((x+y)*(x-y));
482
483 // loop over all basis functions
484 for(int i(0); i < 4; ++i)
485 {
486 data.phi[i].value = _coeff_mat(i,0) + _coeff_mat(i,1)*x + _coeff_mat(i,2)*y + _coeff_mat(i,3)*r;
487 }
488 }
489
491 template<SpaceTags space_cfg_, TrafoTags trafo_cfg_>
494 const Trafo::EvalData<TrafoEvalTraits, trafo_cfg_>& trafo_data) const
495 {
496 // transform image point to local coordinate system
497 DomainPointType pt, loc_grad;
498 pt.set_mat_vec_mult(_inv_lin_mat, trafo_data.img_point - _inv_lin_vec);
499
500 // loop over all basis functions
501 for(int i(0); i < 4; ++i)
502 {
503 // calculate local gradients
504 loc_grad(0) = _coeff_mat(i,1) + DataType(2) * _coeff_mat(i,3) * pt[0];
505 loc_grad(1) = _coeff_mat(i,2) - DataType(2) * _coeff_mat(i,3) * pt[1];
506
507 // multiply by transpose inverse linearized trafo matrix for chain rule
508 data.phi[i].grad.set_vec_mat_mult(loc_grad, _inv_lin_mat);
509 }
510 }
511 }; // Evaluator<..., Shape::Quadrilateral>
512
513
519 template<
520 typename Space_,
521 typename TrafoEvaluator_,
522 typename SpaceEvalTraits_>
523 class Evaluator<Space_, TrafoEvaluator_, SpaceEvalTraits_, Shape::Hexahedron> :
524 public EvaluatorBase<
525 Evaluator<Space_, TrafoEvaluator_, SpaceEvalTraits_, Shape::Hexahedron>,
526 TrafoEvaluator_,
527 SpaceEvalTraits_>
528 {
529 public:
532
534 typedef Space_ SpaceType;
535
537 typedef SpaceEvalTraits_ SpaceEvalTraits;
538
540 typedef TrafoEvaluator_ TrafoEvaluator;
541
543 typedef typename TrafoEvaluator::EvalTraits TrafoEvalTraits;
544
546 typedef typename TrafoEvaluator::TrafoType TrafoType;
547
549 typedef typename TrafoType::MeshType MeshType;
550
552 typedef typename MeshType::template IndexSet<3,2>::Type FacetIndexSetType;
553
555 typedef typename SpaceEvalTraits::DataType DataType;
556
558 typedef typename SpaceEvalTraits::EvalPolicy EvalPolicy;
559
561 typedef typename EvalPolicy::DomainPointType DomainPointType;
563 typedef typename EvalPolicy::ImagePointType ImagePointType;
564
566 typedef typename EvalPolicy::JacobianMatrixType JacobianMatrixType;
568 typedef typename EvalPolicy::JacobianInverseType JacobianInverseType;
569
570 static constexpr SpaceTags eval_caps = SpaceTags::value | SpaceTags::grad;
571
572 template<SpaceTags cfg_>
573 struct ConfigTraits
574 {
576 static constexpr SpaceTags config = cfg_;
577
579 static constexpr TrafoTags trafo_config = TrafoTags::img_point;
580
583 };
584
585 protected:
587 static constexpr TrafoTags inv_lin_trafo_config = TrafoTags::dom_point | TrafoTags::img_point | TrafoTags::jac_inv;
588
590 typedef typename TrafoEvaluator::template ConfigTraits<inv_lin_trafo_config>::EvalDataType InvLinTrafoData;
591
593 typedef typename TrafoType::template Evaluator<Shape::Hypercube<2>, DataType>::Type FacetTrafoEvaluator;
595 typedef typename FacetTrafoEvaluator::EvalTraits FacetEvalTraits;
596
598 static constexpr TrafoTags facet_trafo_config = TrafoTags::img_point | TrafoTags::jac_det;
599
601 typedef typename FacetTrafoEvaluator::template ConfigTraits<facet_trafo_config>::EvalDataType FacetTrafoData;
602
605
608
609 // inverse linearized trafo vector
610 ImagePointType _inv_lin_vec;
611
614
615 void _build_inv_lin_trafo(const TrafoEvaluator& trafo_eval)
616 {
617 // create a domain point in the barycentre of the cell
619
620 // create the trafo data
621 InvLinTrafoData trafo_data;
622 trafo_eval(trafo_data, dom_point);
623
624 // store inverse trafo linearization
625 _inv_lin_mat = trafo_data.jac_inv;
626 _inv_lin_vec = trafo_data.img_point;
627 }
628
630 void _build_coeff_matrix(const TrafoEvaluator& trafo_eval)
631 {
632 // fetch the trafo
633 const TrafoType& trafo = trafo_eval.get_trafo();
634
635 // fetch the mesh
636 const MeshType& mesh = trafo.get_mesh();
637
638 // fetch the facet-index-set of the mesh
639 const FacetIndexSetType& facet_index_set = mesh.template get_index_set<3,2>();
640
641 // fetch the cell index of the currently active cell
642 const Index cell = trafo_eval.get_cell_index();
643
644 // create a nodal matrix
645 CoeffMatrixType _nodal_mat(DataType(0));
646
647 // create a facet trafo evaluator
648 FacetTrafoEvaluator facet_eval(trafo);
649
650 // create facet evaluation data
651 FacetTrafoData facet_data;
652
653 // define 2-point Gauss cubature point coordinate
654 const DataType g = Math::sqrt(DataType(1) / DataType(3));
655 typename FacetEvalTraits::DomainPointType g1, g2, g3, g4;
656 DomainPointType q1, q2, q3, q4;
657
658 // set up cubature points
659 g1[0] = g2[0] = -g;
660 g3[0] = g4[0] = +g;
661 g1[1] = g3[1] = -g;
662 g2[1] = g4[1] = +g;
663
664 // loop over all 6 faces of the hexa
665 for(int i(0); i < 6; ++i)
666 {
667 // initialize facet evaluator
668 facet_eval.prepare(Index(facet_index_set(cell, i)));
669
670 // map first cubature point
671 facet_eval(facet_data, g1);
672 DataType w1(facet_data.jac_det);
673 q1.set_mat_vec_mult(_inv_lin_mat, facet_data.img_point - _inv_lin_vec);
674
675 // map second cubature point
676 facet_eval(facet_data, g2);
677 DataType w2(facet_data.jac_det);
678 q2.set_mat_vec_mult(_inv_lin_mat, facet_data.img_point - _inv_lin_vec);
679
680 // map third cubature point
681 facet_eval(facet_data, g3);
682 DataType w3(facet_data.jac_det);
683 q3.set_mat_vec_mult(_inv_lin_mat, facet_data.img_point - _inv_lin_vec);
684
685 // map fourth cubature point
686 facet_eval(facet_data, g4);
687 DataType w4(facet_data.jac_det);
688 q4.set_mat_vec_mult(_inv_lin_mat, facet_data.img_point - _inv_lin_vec);
689
690 // release facet evaluator
691 facet_eval.finish();
692
693 // compute reciprocal quad area
694 DataType v = DataType(1) / (w1 + w2 + w3 + w4);
695
696 // evaluate node functionals
697 _nodal_mat(0,i) = DataType(1); // = v*(w1 + w2 + w3 + w4)
698 _nodal_mat(1,i) = v*(w1*q1[0] + w2*q2[0] + w3*q3[0] + w4*q4[0]);
699 _nodal_mat(2,i) = v*(w1*q1[1] + w2*q2[1] + w3*q3[1] + w4*q4[1]);
700 _nodal_mat(3,i) = v*(w1*q1[2] + w2*q2[2] + w3*q3[2] + w4*q4[2]);
701 _nodal_mat(4,i) = v*( // x^2 - y^2
702 w1*(q1[0] + q1[1])*(q1[0] - q1[1]) + w2*(q2[0] + q2[1])*(q2[0] - q2[1]) +
703 w3*(q3[0] + q3[1])*(q3[0] - q3[1]) + w4*(q4[0] + q4[1])*(q4[0] - q4[1]));
704 _nodal_mat(5,i) = v*( // y^2 - z^2
705 w1*(q1[1] + q1[2])*(q1[1] - q1[2]) + w2*(q2[1] + q2[2])*(q2[1] - q2[2]) +
706 w3*(q3[1] + q3[2])*(q3[1] - q3[2]) + w4*(q4[1] + q4[2])*(q4[1] - q4[2]));
707 }
708
709 // invert nodal matrix to obtain coefficient matrix
710 _coeff_mat.set_inverse(_nodal_mat);
711 }
712
713 public:
720 explicit Evaluator(const SpaceType& DOXY(space))
721 {
722 }
723
725 virtual ~Evaluator()
726 {
727 }
728
736 {
737 return 6;
738 }
739
746 void prepare(const TrafoEvaluator& trafo_eval)
747 {
748 // prepare inverse linearized trafo
749 _build_inv_lin_trafo(trafo_eval);
750
751 // build coefficient matrix
752 _build_coeff_matrix(trafo_eval);
753 }
754
756 template<SpaceTags space_cfg_, TrafoTags trafo_cfg_>
759 const Trafo::EvalData<TrafoEvalTraits, trafo_cfg_>& trafo_data) const
760 {
761 // transform image point
763 pt.set_mat_vec_mult(_inv_lin_mat, trafo_data.img_point - _inv_lin_vec);
764
765 // "evaluate" monomials
766 DataType x(pt[0]);
767 DataType y(pt[1]);
768 DataType z(pt[2]);
769 DataType rxy((x+y)*(x-y)); // = x^2 - y^2
770 DataType ryz((y+z)*(y-z)); // = y^2 - z^2
771
772 // loop over all basis functions
773 for(int i(0); i < 6; ++i)
774 {
775 data.phi[i].value = _coeff_mat(i,0) + _coeff_mat(i,1)*x + _coeff_mat(i,2)*y + _coeff_mat(i,3)*z
776 + _coeff_mat(i,4)*rxy + _coeff_mat(i,5)*ryz;
777 }
778 }
779
781 template<SpaceTags space_cfg_, TrafoTags trafo_cfg_>
784 const Trafo::EvalData<TrafoEvalTraits, trafo_cfg_>& trafo_data) const
785 {
786 // transform image point to local coordinate system
787 DomainPointType pt, loc_grad;
788 pt.set_mat_vec_mult(_inv_lin_mat, trafo_data.img_point - _inv_lin_vec);
789
790 // loop over all basis functions
791 for(int i(0); i < 6; ++i)
792 {
793 // calculate local gradients
794 loc_grad[0] = _coeff_mat(i,1) + DataType(2) * pt[0] * _coeff_mat(i,4);
795 loc_grad[1] = _coeff_mat(i,2) + DataType(2) * pt[1] * (_coeff_mat(i,5) - _coeff_mat(i,4));
796 loc_grad[2] = _coeff_mat(i,3) - DataType(2) * pt[2] * _coeff_mat(i,5);
797
798 // multiply by transpose inverse linearized trafo matrix for chain rule
799 data.phi[i].grad.set_vec_mat_mult(loc_grad, _inv_lin_mat);
800 }
801 }
802 }; // Evaluator<..., Shape::Quadrilateral>
803 } // namespace CroRavRanTur
804 } // namespace Space
805} // namespace FEAT
EvalTraits_::BasisValueType value
basis function value object
Definition: eval_data.hpp:47
EvalTraits_::BasisGradientType grad
basis gradient object
Definition: eval_data.hpp:45
void prepare(const TrafoEvaluator &trafo_eval)
Prepares the evaluator for a given cell.
Definition: evaluator.hpp:459
TrafoEvaluator::template ConfigTraits< inv_lin_trafo_config >::EvalDataType InvLinTrafoData
inverse linearized trafo data
Definition: evaluator.hpp:325
void _build_coeff_matrix(const TrafoEvaluator &trafo_eval)
computes the basis function coefficient matrix for the current cell
Definition: evaluator.hpp:365
void eval_gradients(EvalData< SpaceEvalTraits, space_cfg_ > &data, const Trafo::EvalData< TrafoEvalTraits, trafo_cfg_ > &trafo_data) const
Definition: evaluator.hpp:492
FacetTrafoEvaluator::template ConfigTraits< facet_trafo_config >::EvalDataType FacetTrafoData
facet trafo data
Definition: evaluator.hpp:336
TrafoType::template Evaluator< Shape::Hypercube< 1 >, DataType >::Type FacetTrafoEvaluator
trafo evaluator for facets (=edges)
Definition: evaluator.hpp:328
void eval_values(EvalData< SpaceEvalTraits, space_cfg_ > &data, const Trafo::EvalData< TrafoEvalTraits, trafo_cfg_ > &trafo_data) const
Evaluates the basis function values on the real cell.
Definition: evaluator.hpp:470
EvaluatorBase< Evaluator, TrafoEvaluator_, SpaceEvalTraits_ > BaseClass
base-class typedef
Definition: evaluator.hpp:266
void eval_ref_gradients(EvalData_ &data, const DomainPointType &point) const
Evaluates the basis function gradients on the reference cell.
Definition: evaluator.hpp:230
ParametricEvaluator< Evaluator, TrafoEvaluator_, SpaceEvalTraits_, ref_caps > BaseClass
base-class typedef
Definition: evaluator.hpp:161
void eval_ref_values(EvalData_ &data, const DomainPointType &point) const
Evaluates the basis function values on the reference cell.
Definition: evaluator.hpp:210
void prepare(const TrafoEvaluator &trafo_eval)
Prepares the evaluator for a given cell.
Definition: evaluator.hpp:746
TrafoEvaluator::template ConfigTraits< inv_lin_trafo_config >::EvalDataType InvLinTrafoData
inverse linearized trafo data
Definition: evaluator.hpp:590
FacetTrafoEvaluator::template ConfigTraits< facet_trafo_config >::EvalDataType FacetTrafoData
facet trafo data
Definition: evaluator.hpp:601
EvaluatorBase< Evaluator, TrafoEvaluator_, SpaceEvalTraits_ > BaseClass
base-class typedef
Definition: evaluator.hpp:531
void eval_values(EvalData< SpaceEvalTraits, space_cfg_ > &data, const Trafo::EvalData< TrafoEvalTraits, trafo_cfg_ > &trafo_data) const
Evaluates the basis function values on the real cell.
Definition: evaluator.hpp:757
void eval_gradients(EvalData< SpaceEvalTraits, space_cfg_ > &data, const Trafo::EvalData< TrafoEvalTraits, trafo_cfg_ > &trafo_data) const
Definition: evaluator.hpp:782
void _build_coeff_matrix(const TrafoEvaluator &trafo_eval)
computes the basis function coefficient matrix for the current cell
Definition: evaluator.hpp:630
TrafoType::template Evaluator< Shape::Hypercube< 2 >, DataType >::Type FacetTrafoEvaluator
trafo evaluator for facets
Definition: evaluator.hpp:593
Tiny::Matrix< DataType, 6, 6 > CoeffMatrixType
basis function coefficient matrix
Definition: evaluator.hpp:604
MeshType::template IndexSet< 3, 2 >::Type FacetIndexSetType
facet-index-set type
Definition: evaluator.hpp:552
void eval_ref_gradients(EvalData_ &data, const DomainPointType &point) const
Evaluates the basis function gradients on the reference cell.
Definition: evaluator.hpp:126
ParametricEvaluator< Evaluator, TrafoEvaluator_, SpaceEvalTraits_, ref_caps > BaseClass
base-class typedef
Definition: evaluator.hpp:58
void eval_ref_values(EvalData_ &data, const DomainPointType &point) const
Evaluates the basis function values on the reference cell.
Definition: evaluator.hpp:107
Crouzeix-Raviart / Rannacher-Turek Element Evaluator class template declaration.
Definition: evaluator.hpp:34
Space evaluation data structure.
Definition: eval_data.hpp:141
BasisDataType phi[max_local_dofs]
the basis function data vector
Definition: eval_data.hpp:150
Basic Space Evaluator CRTP base-class template.
Finite-Element Parametric Evaluator CRTP base-class template.
Tiny Matrix class template.
CUDA_HOST_DEVICE Matrix & set_inverse(const Matrix< T_, m_, n_, sma_, sna_ > &a)
Sets this matrix to the inverse of another matrix.
Trafo evaluation data structure.
Definition: eval_data.hpp:33
EvalTraits::ImagePointType img_point
image point
Definition: eval_data.hpp:53
T_ sqrt(T_ x)
Returns the square-root of a value.
Definition: math.hpp:300
static constexpr SpaceTags ref_caps
Crouzeix-Raviart Element Evaluator reference capabilities.
Definition: evaluator.hpp:22
FEAT namespace.
Definition: adjactor.hpp:12
SpaceTags
Space configuration tags enum.
Definition: eval_tags.hpp:97
@ value
specifies whether the space should supply basis function values
@ ref_value
specifies whether the space should supply reference basis function values
@ ref_grad
specifies whether the space should supply reference basis function gradients
@ grad
specifies whether the space should supply basis function gradients
std::uint64_t Index
Index data type.
TrafoTags
Trafo configuration tags enum.
Definition: eval_tags.hpp:22
@ img_point
specifies whether the trafo should supply image point coordinates
@ dom_point
specifies whether the trafo should supply domain point coordinates
@ jac_inv
specifies whether the trafo should supply inverse jacobian matrices
@ jac_det
specifies whether the trafo should supply jacobian determinants