FEAT 3
Finite Element Analysis Toolbox
Loading...
Searching...
No Matches
burgers_assembly_job.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#include <kernel/assembly/base.hpp>
9#include <kernel/assembly/asm_traits.hpp>
10#include <kernel/cubature/dynamic_factory.hpp>
11#include <kernel/lafem/dense_vector.hpp>
12#include <kernel/lafem/sparse_matrix_bcsr.hpp>
13#include <kernel/lafem/dense_vector_blocked.hpp>
14#include <kernel/global/vector.hpp>
15
16namespace FEAT
17{
18 namespace Assembly
19 {
100 template<typename DataType_, typename Space_, typename ConvVector_>
102 {
103 public:
105 typedef DataType_ DataType;
106
108 typedef Space_ SpaceType;
109
111 typedef ConvVector_ ConvVectorType;
112
115
117 DataType_ nu;
118
120 DataType_ theta;
121
123 DataType_ beta;
124
126 DataType_ frechet_beta;
127
129 DataType_ sd_delta;
130
132 DataType_ sd_nu;
133
135 DataType_ sd_v_norm;
136
138 const ConvVector_& convection_vector;
139
141 const Space_& space;
142
145
146 public:
159 explicit BurgersAssemblyJobBase(const ConvVector_& conv_vector, const Space_& space_, const String& cubature) :
160 deformation(false),
161 nu(DataType_(0)),
162 theta(DataType_(0)),
163 beta(DataType_(0)),
164 frechet_beta(DataType_(0)),
165 sd_delta(DataType_(0)),
166 sd_nu(DataType_(0)),
167 sd_v_norm(DataType_(0)),
168 convection_vector(conv_vector),
169 space(space_),
170 cubature_factory(cubature)
171 {
172 }
173
183 template<typename IndexType_, int conv_dim_>
185 {
186 const auto* vals = convect.elements();
187 DataType_ r = DataType_(0);
188 for(Index i(0); i < convect.size(); ++i)
189 r = Math::max(r, vals[i].norm_euclid());
190 return r;
191 }
192
205 template<typename LocalVector_, typename Mirror_>
207 {
208 const auto* gate = convect.get_gate();
209 if(gate != nullptr)
210 return gate->max(calc_sd_v_norm(convect.local()));
211 else
212 return calc_sd_v_norm(convect.local());
213 }
214
221 template<typename VectorType_>
222 void set_sd_v_norm(const VectorType_& convect)
223 {
224 this->sd_v_norm = calc_sd_v_norm(convect);
225 }
226 }; // class BurgersAssemblyJobBase<...>
227
253 template<typename Job_, typename DataType_>
255 {
256 public:
258 static constexpr bool need_scatter = true;
260 static constexpr bool need_combine = false;
261
263 typedef DataType_ DataType;
264
266 typedef typename Job_::SpaceType SpaceType;
267
269 typedef typename Job_::ConvVectorType ConvVectorType;
270
272 static constexpr int shape_dim = SpaceType::shape_dim;
273
275 static constexpr int conv_dim = ConvVectorType::BlockSize;
276
279
281 const Job_& job;
282
285
287 const bool deformation;
288
291
294
297
300
303
306
309
311 const bool need_diff, need_conv, need_conv_frechet, need_reac, need_streamdiff;
312
314 bool need_loc_v, need_loc_grad_v, need_mean_v, need_local_h;
315
333 typename ConvVectorType::GatherAxpy gather_conv;
334
336 static constexpr int max_local_dofs = AsmTraits::max_local_test_dofs;
337
340
343
346
349
352
355
358
361
362 public:
369 explicit BurgersAssemblyTaskBase(const Job_& job_) :
370 job(job_),
371 tol_eps(Math::sqrt(Math::eps<DataType>())),
373 nu(job_.nu),
374 theta(job_.theta),
375 beta(job_.beta),
377 sd_delta(job_.sd_delta),
378 sd_nu(job_.sd_nu),
379 sd_v_norm(job_.sd_v_norm),
380 need_diff(Math::abs(nu) > DataType(0)),
381 need_conv(Math::abs(beta) > DataType(0)),
382 need_conv_frechet(Math::abs(frechet_beta) > DataType(0)),
383 need_reac(Math::abs(theta) > DataType(0)),
384 need_streamdiff((Math::abs(sd_delta) > DataType(0)) && (sd_v_norm > tol_eps)),
385 need_loc_v(need_conv),
386 need_loc_grad_v(need_conv_frechet),
387 need_mean_v(need_streamdiff),
388 need_local_h(need_streamdiff),
389 space(job_.space),
390 trafo(space.get_trafo()),
394 cubature_rule(Cubature::ctor_factory, job_.cubature_factory),
395 trafo_data(),
396 space_data(),
397 gather_conv(job_.convection_vector),
400 loc_v(),
401 mean_v(),
402 loc_grad_v(),
404 barycenter(),
405 local_h(DataType(0)),
407 {
408 // compute reference element barycenter
409 for(int i(0); i < shape_dim; ++i)
411 }
412
419 void prepare(const Index cell)
420 {
421 // prepare dof mapping
422 dof_mapping.prepare(cell);
423
424 // prepare trafo evaluator
425 trafo_eval.prepare(cell);
426
427 // prepare space evaluator
428 space_eval.prepare(trafo_eval);
429
430 // fetch number of local dofs
431 num_local_dofs = space_eval.get_num_local_dofs();
432
433 // gather our local convection dofs
434 local_conv_dofs.format();
436
437 // compute mesh width if necessary
438 if(need_mean_v || need_local_h || need_streamdiff)
439 {
440 // reset local h and delta
442
443 // evaluate trafo and space at barycenter
446
447 // compute velocity at barycenter
448 mean_v.format();
449 for(int i(0); i < num_local_dofs; ++i)
450 mean_v.axpy(space_data.phi[i].value, local_conv_dofs[i]);
451
452 // compute norm of mean velocity
453 const DataType local_norm_v = mean_v.norm_euclid();
454
455 // do we have a non-zero velocity?
456 if(local_norm_v > tol_eps)
457 {
458 // compute local mesh width w.r.t. mean velocity
459 local_h = trafo_eval.width_directed(mean_v) * local_norm_v;
460
461 // compute local Re_T
462 const DataType local_re = (local_norm_v * local_h) / this->sd_nu;
463
464 // compute local delta
465 local_delta = this->sd_delta * (local_h / this->sd_v_norm) * (DataType(2)*local_re) / (DataType(1) + local_re);
466 }
467 }
468 }
469
480 void prepare_point(int cubature_point)
481 {
482 // compute trafo data
483 trafo_eval(trafo_data, cubature_rule.get_point(cubature_point));
484
485 // compute basis function data
487
488 // evaluate convection function and its gradient (if required)
489 if(need_loc_v || need_conv || need_streamdiff)
490 {
491 loc_v.format();
492 for(int i(0); i < num_local_dofs; ++i)
493 {
494 // update velocity value
495 loc_v.axpy(space_data.phi[i].value, local_conv_dofs[i]);
496 }
497 }
498 if(need_loc_grad_v || need_conv_frechet)
499 {
501 for(int i(0); i < num_local_dofs; ++i)
502 {
503 // update velocity gradient
505 }
506 }
507
508 // evaluate streamline diffusion coefficients
509 if(need_streamdiff)
510 {
511 for(int i(0); i < num_local_dofs; ++i)
512 {
513 streamdiff_coeffs[i] = Tiny::dot(loc_v, space_data.phi[i].grad);
514 }
515 }
516 }
517
521 void finish()
522 {
523 // finish evaluators
524 space_eval.finish();
525 trafo_eval.finish();
526
527 // finish dof mapping
528 dof_mapping.finish();
529 }
530
534 void combine()
535 {
536 // nothing to do here
537 }
538 }; // class BurgersAssemblyTaskBase<...>
539
559 template<typename Job_, typename DataType_, int block_size_>
561 public BurgersAssemblyTaskBase<Job_, DataType_>
562 {
563 public:
565
566 typedef DataType_ DataType;
567
570
571 protected:
575 LocalMatrixType local_matrix;
576
577 public:
584 explicit BurgersBlockedAssemblyTaskBase(const Job_& job_) :
585 BaseClass(job_),
586 local_matrix()
587 {
588 }
589
596 void assemble_burgers_point(const DataType weight)
597 {
598 // assemble diffusion matrix?
599 if(this->need_diff)
600 {
601 // assemble gradient-tensor diffusion
602
603 // test function loop
604 for(int i(0); i < this->num_local_dofs; ++i)
605 {
606 // trial function loop
607 for(int j(0); j < this->num_local_dofs; ++j)
608 {
609 // compute scalar value
610 const DataType value = this->nu * weight * Tiny::dot(this->space_data.phi[i].grad, this->space_data.phi[j].grad);
611
612 // update local matrix
613 local_matrix[i][j].add_scalar_main_diag(value);
614 }
615 }
616
617 // assemble deformation tensor?
618 if(this->deformation)
619 {
620 // test function loop
621 for(int i(0); i < this->num_local_dofs; ++i)
622 {
623 // trial function loop
624 for(int j(0); j < this->num_local_dofs; ++j)
625 {
626 // add outer product of grad(phi) and grad(psi)
627 local_matrix[i][j].add_outer_product(this->space_data.phi[j].grad, this->space_data.phi[i].grad, this->nu * weight);
628 }
629 }
630 }
631 }
632
633 // assemble convection?
634 if(this->need_conv)
635 {
636 // test function loop
637 for(int i(0); i < this->num_local_dofs; ++i)
638 {
639 // trial function loop
640 for(int j(0); j < this->num_local_dofs; ++j)
641 {
642 // compute scalar value
643 const DataType value = this->beta * weight * this->space_data.phi[i].value * Tiny::dot(this->loc_v, this->space_data.phi[j].grad);
644
645 // update local matrix
646 local_matrix[i][j].add_scalar_main_diag(value);
647 }
648 }
649 }
650
651 // assemble convection Frechet?
652 if(this->need_conv_frechet)
653 {
654 // test function loop
655 for(int i(0); i < this->num_local_dofs; ++i)
656 {
657 // trial function loop
658 for(int j(0); j < this->num_local_dofs; ++j)
659 {
660 // compute scalar value
661 const DataType value = this->frechet_beta * weight * this->space_data.phi[i].value * this->space_data.phi[j].value;
662
663 // update local matrix
664 local_matrix[i][j].axpy(value, this->loc_grad_v);
665 }
666 }
667 }
668
669 // assemble reaction?
670 if(this->need_reac)
671 {
672 // test function loop
673 for(int i(0); i < this->num_local_dofs; ++i)
674 {
675 // trial function loop
676 for(int j(0); j < this->num_local_dofs; ++j)
677 {
678 // compute scalar value
679 const DataType value = this->theta * weight * this->space_data.phi[i].value * this->space_data.phi[j].value;
680
681 // update local matrix
682 local_matrix[i][j].add_scalar_main_diag(value);
683 }
684 }
685 }
686 }
687
694 void assemble_streamline_diffusion(const DataType weight)
695 {
696 // assemble streamline diffusion?
697 if(this->need_streamdiff && (this->local_delta > this->tol_eps))
698 {
699 // test function loop
700 for(int i(0); i < this->num_local_dofs; ++i)
701 {
702 // trial function loop
703 for(int j(0); j < this->num_local_dofs; ++j)
704 {
705 // compute scalar value
706 const DataType value = this->local_delta * weight * this->streamdiff_coeffs[i] * this->streamdiff_coeffs[j];
707
708 // update local matrix
709 local_matrix[i][j].add_scalar_main_diag(value);
710 }
711 }
712 }
713 }
714
718 void assemble()
719 {
720 local_matrix.format();
721
722 // loop over all quadrature points and integrate
723 for(int point(0); point < this->cubature_rule.get_num_points(); ++point)
724 {
725 // prepare trafo and space for cubature point
726 this->prepare_point(point);
727
728 // precompute cubature weight
729 const DataType weight = this->trafo_data.jac_det * this->cubature_rule.get_weight(point);
730
731 // assemble the Burgers operator
732 this->assemble_burgers_point(weight);
733
734 // assemble the streamline diffusion
735 this->assemble_streamline_diffusion(weight);
736 } // continue with next cubature point
737 }
738 }; // class BurgersBlockedAssemblyTaskBase<...>
739
760 template<typename Matrix_, typename Space_, typename ConvVector_ = typename Matrix_::VectorTypeR>
762 public BurgersAssemblyJobBase<typename Matrix_::DataType, Space_, ConvVector_>
763 {
764 public:
766 typedef Matrix_ MatrixType;
768 typedef typename MatrixType::DataType DataType;
769
772
773 // no nonsense, please
774 static_assert(Matrix_::BlockHeight == Matrix_::BlockWidth, "only square matrix blocks are supported here");
775
777 static constexpr int block_size = Matrix_::BlockHeight;
778
781
782 public:
798 explicit BurgersBlockedMatrixAssemblyJob(Matrix_& matrix_, const ConvVector_& conv_vector,
799 const Space_& space_, const String& cubature_name) :
800 BaseClass(conv_vector, space_, cubature_name),
801 matrix(matrix_)
802 {
803 }
804
805 public:
807 class Task :
808 public BurgersBlockedAssemblyTaskBase<BurgersBlockedMatrixAssemblyJob, DataType, block_size>
809 {
810 public:
812
813 protected:
815 typename MatrixType::ScatterAxpy scatter_matrix;
816
817 public:
819 explicit Task(const BurgersBlockedMatrixAssemblyJob& job_) :
820 BaseClass(job_),
822 {
823 }
824
825 // prepare, assemble and finish are already implemented in the base classes
826
828 void scatter()
829 {
830 this->scatter_matrix(this->local_matrix, this->dof_mapping, this->dof_mapping);
831 }
832 }; // class BurgersBlockedMatrixAssemblyJob::Task
833 }; // class BurgersBlockedMatrixAssemblyJob
834
855 template<typename Vector_, typename Space_, typename ConvVector_ = Vector_>
857 public BurgersAssemblyJobBase<typename Vector_::DataType, Space_, ConvVector_>
858 {
859 public:
861 typedef Vector_ VectorType;
863 typedef typename VectorType::DataType DataType;
864
867
869 static constexpr int block_size = Vector_::BlockSize;
870
873
876
877 public:
898 explicit BurgersBlockedVectorAssemblyJob(Vector_& rhs_vector_, const Vector_& sol_vector_,
899 const ConvVector_& conv_vector, const Space_& space_, const String& cubature_name) :
900 BaseClass(conv_vector, space_, cubature_name),
901 vector_rhs(rhs_vector_),
902 vector_sol(sol_vector_)
903 {
904 XASSERTM(&rhs_vector_ != &sol_vector_, "rhs and solution vectors must not be the same object");
905 // compare void addresses to avoid warnings in case the classes are different
906 XASSERTM((void*)&rhs_vector_ != (void*)&conv_vector, "rhs and convection vectors must not be the same object");
907 }
908
909 public:
911 class Task :
912 public BurgersBlockedAssemblyTaskBase<BurgersBlockedVectorAssemblyJob, DataType, block_size>
913 {
914 public:
916
918
919 protected:
921 typename VectorType::GatherAxpy gather_vec_sol;
923 typename VectorType::ScatterAxpy scatter_vec_rhs;
925 const bool sol_equals_conv;
926
929
932
933 public:
935 explicit Task(const BurgersBlockedVectorAssemblyJob& job_) :
936 BaseClass(job_),
939 sol_equals_conv((void*)&job_.vector_sol == (void*)&job_.convection_vector),
942 {
943 }
944
945 void prepare(const Index cell)
946 {
947 BaseClass::prepare(cell);
948
949 // gather local solution dofs if required
950 if(!sol_equals_conv)
951 {
952 this->local_sol_dofs.format();
953 this->gather_vec_sol(this->local_sol_dofs, this->dof_mapping);
954 }
955 }
956
961 {
962 local_vector.format();
963
964 // compute local vector from local matrix and sol vector dofs
965 if(!sol_equals_conv)
966 {
967 for(int i(0); i < this->num_local_dofs; ++i)
968 for(int j(0); j < this->num_local_dofs; ++j)
969 this->local_vector[i].add_mat_vec_mult(this->local_matrix(i,j), this->local_sol_dofs[j]);
970 }
971 else
972 {
973 for(int i(0); i < this->num_local_dofs; ++i)
974 for(int j(0); j < this->num_local_dofs; ++j)
975 this->local_vector[i].add_mat_vec_mult(this->local_matrix(i,j), this->local_conv_dofs[j]);
976 }
977 }
978
979 void assemble()
980 {
981 // assemble the local matrix
983
984 // compute the local vector from the matrix
986 }
987
989 void scatter()
990 {
991 this->scatter_vec_rhs(this->local_vector, this->dof_mapping);
992 }
993 }; // class BurgersBlockedVectorAssemblyJob::Task
994 }; // class BurgersBlockedVectorAssemblyJob<...>
995
1012 template<typename Job_, typename DataType_>
1014 public BurgersAssemblyTaskBase<Job_, DataType_>
1015 {
1016 public:
1018
1019 typedef DataType_ DataType;
1020
1023
1024 protected:
1027 LocalMatrixType local_matrix;
1028
1029 public:
1036 explicit BurgersScalarAssemblyTaskBase(const Job_& job_) :
1037 BaseClass(job_),
1038 local_matrix()
1039 {
1040 }
1041
1048 void assemble_burgers_point(const DataType weight)
1049 {
1050 // assemble diffusion matrix?
1051 if(this->need_diff)
1052 {
1053 // assemble gradient-tensor diffusion
1054
1055 // test function loop
1056 for(int i(0); i < this->num_local_dofs; ++i)
1057 {
1058 // trial function loop
1059 for(int j(0); j < this->num_local_dofs; ++j)
1060 {
1061 local_matrix[i][j] += this->nu * weight
1062 * Tiny::dot(this->space_data.phi[i].grad, this->space_data.phi[j].grad);
1063 }
1064 }
1065 }
1066
1067 // assemble convection?
1068 if(this->need_conv)
1069 {
1070 // test function loop
1071 for(int i(0); i < this->num_local_dofs; ++i)
1072 {
1073 // trial function loop
1074 for(int j(0); j < this->num_local_dofs; ++j)
1075 {
1076 local_matrix[i][j] += this->beta * weight
1077 * this->space_data.phi[i].value * Tiny::dot(this->loc_v, this->space_data.phi[j].grad);
1078 }
1079 }
1080 }
1081
1082 // assemble reaction?
1083 if(this->need_reac)
1084 {
1085 // test function loop
1086 for(int i(0); i < this->num_local_dofs; ++i)
1087 {
1088 // trial function loop
1089 for(int j(0); j < this->num_local_dofs; ++j)
1090 {
1091 local_matrix[i][j] += this->theta * weight
1092 * this->space_data.phi[i].value * this->space_data.phi[j].value;
1093 }
1094 }
1095 }
1096 }
1097
1104 void assemble_streamline_diffusion(const DataType weight)
1105 {
1106 // assemble streamline diffusion?
1107 if(this->need_streamdiff && (this->local_delta > this->tol_eps))
1108 {
1109 // test function loop
1110 for(int i(0); i < this->num_local_dofs; ++i)
1111 {
1112 // trial function loop
1113 for(int j(0); j < this->num_local_dofs; ++j)
1114 {
1115 local_matrix[i][j] += this->local_delta * weight
1116 * this->streamdiff_coeffs[i] * this->streamdiff_coeffs[j];
1117 }
1118 }
1119 }
1120 }
1121
1126 {
1127 local_matrix.format();
1128
1129 // loop over all quadrature points and integrate
1130 for(int point(0); point < this->cubature_rule.get_num_points(); ++point)
1131 {
1132 // prepare trafo and space for cubature point
1133 this->prepare_point(point);
1134
1135 // precompute cubature weight
1136 const DataType weight = this->trafo_data.jac_det * this->cubature_rule.get_weight(point);
1137
1138 // assemble the Burgers operator
1139 this->assemble_burgers_point(weight);
1140
1141 // assembles the streamline diffusion stabilization
1142 this->assemble_streamline_diffusion(weight);
1143 } // continue with next cubature point
1144 }
1145 }; // class BurgersScalarAssemblyTaskBase<...>
1146
1167 template<typename Matrix_, typename Space_, typename ConvVector_>
1169 public BurgersAssemblyJobBase<typename Matrix_::DataType, Space_, ConvVector_>
1170 {
1171 public:
1173 typedef Matrix_ MatrixType;
1175 typedef typename MatrixType::DataType DataType;
1176
1179
1182
1183 public:
1199 explicit BurgersScalarMatrixAssemblyJob(Matrix_& matrix_, const ConvVector_& conv_vector,
1200 const Space_& space_, const String& cubature_name) :
1201 BaseClass(conv_vector, space_, cubature_name),
1202 matrix(matrix_)
1203 {
1204 }
1205
1206 public:
1208 class Task :
1209 public BurgersScalarAssemblyTaskBase<BurgersScalarMatrixAssemblyJob, DataType>
1210 {
1211 public:
1213
1214 protected:
1216 typename MatrixType::ScatterAxpy scatter_matrix;
1217
1218 public:
1220 explicit Task(const BurgersScalarMatrixAssemblyJob& job_) :
1221 BaseClass(job_),
1223 {
1224 }
1225
1226 // prepare, assemble and finish are already implemented in the base classes
1227
1229 void scatter()
1230 {
1231 this->scatter_matrix(this->local_matrix, this->dof_mapping, this->dof_mapping);
1232 }
1233 }; // class BurgersScalarMatrixAssemblyJob::Task
1234 }; // class BurgersScalarMatrixAssemblyJob
1235
1256 template<typename Vector_, typename Space_, typename ConvVector_>
1258 public BurgersAssemblyJobBase<typename Vector_::DataType, Space_, ConvVector_>
1259 {
1260 public:
1262 typedef Vector_ VectorType;
1264 typedef typename VectorType::DataType DataType;
1265
1268
1271
1274
1275 public:
1296 explicit BurgersScalarVectorAssemblyJob(Vector_& rhs_vector_, const Vector_& sol_vector_,
1297 const ConvVector_& conv_vector, const Space_& space_, const String& cubature_name) :
1298 BaseClass(conv_vector, space_, cubature_name),
1299 vector_rhs(rhs_vector_),
1300 vector_sol(sol_vector_)
1301 {
1302 XASSERTM(&rhs_vector_ != &sol_vector_, "rhs and solution vectors must not be the same object");
1303 }
1304
1305 public:
1307 class Task :
1308 public BurgersScalarAssemblyTaskBase<BurgersScalarVectorAssemblyJob, DataType>
1309 {
1310 public:
1312
1314
1315 protected:
1317 typename VectorType::GatherAxpy gather_vec_sol;
1319 typename VectorType::ScatterAxpy scatter_vec_rhs;
1320
1323
1326
1327 public:
1329 explicit Task(const BurgersScalarVectorAssemblyJob& job_) :
1330 BaseClass(job_),
1334 local_vector()
1335 {
1336 }
1337
1338 void prepare(const Index cell)
1339 {
1340 BaseClass::prepare(cell);
1341
1342 // gather local solution vector dofs
1343 this->local_sol_dofs.format();
1344 this->gather_vec_sol(this->local_sol_dofs, this->dof_mapping);
1345 }
1346
1351 {
1352 // compute local vector from local matrix and sol vector dofs
1353 this->local_vector.format();
1354 for(int i(0); i < this->num_local_dofs; ++i)
1355 for(int j(0); j < this->num_local_dofs; ++j)
1356 this->local_vector[i] += this->local_matrix(i,j) * this->local_sol_dofs[j];
1357 //this->local_vector.set_mat_vec_mult(this->local_matrix, this->local_sol_dofs);
1358 }
1359
1360 void assemble()
1361 {
1362 // assemble the local matrix
1364
1365 // compute local vector from the matrix
1367 }
1368
1370 void scatter()
1371 {
1372 this->scatter_vec_rhs(this->local_vector, this->dof_mapping);
1373 }
1374 }; // class BurgersScalarVectorAssemblyJob::Task
1375 }; // class BurgersScalarVectorAssemblyJob<...>
1376 } // namespace Assembly
1377} // namespace FEAT
#define XASSERTM(expr, msg)
Assertion macro definition with custom message.
Definition: assertion.hpp:263
Common single-space assembly traits class template.
Definition: asm_traits.hpp:83
TrafoType::template Evaluator< ShapeType, DataType >::Type TrafoEvaluator
trafo evaluator type
Definition: asm_traits.hpp:104
SpaceEvaluator::template ConfigTraits< space_config >::EvalDataType SpaceEvalData
space evaluation data types
Definition: asm_traits.hpp:142
Intern::CubatureTraits< TrafoEvaluator >::RuleType CubatureRuleType
cubature rule type
Definition: asm_traits.hpp:186
SpaceType::TrafoType TrafoType
trafo type
Definition: asm_traits.hpp:97
SpaceType::DofMappingType DofMapping
dof-mapping types
Definition: asm_traits.hpp:155
TrafoEvaluator::template ConfigTraits< trafo_config >::EvalDataType TrafoEvalData
trafo evaluation data type
Definition: asm_traits.hpp:138
SpaceType::template Evaluator< TrafoEvaluator >::Type SpaceEvaluator
space evaluator types
Definition: asm_traits.hpp:110
Base class for Burgers assembly jobs.
const ConvVector_ & convection_vector
the convection vector to be used
ConvVector_ ConvVectorType
the convection vector type
DataType_ sd_delta
scaling parameter for streamline diffusion stabilization operator S
static DataType calc_sd_v_norm(const Global::Vector< LocalVector_, Mirror_ > &convect)
Calculates the convection field norm for the streamline diffusion parameter delta_T.
void set_sd_v_norm(const VectorType_ &convect)
Sets the convection field norm for the local streamline diffusion parameter delta_T.
DataType_ beta
scaling parameter for convective operator K
DataType_ frechet_beta
scaling parameter for Frechet derivative of convective operator K'
bool deformation
specifies whether to use the deformation tensor
Cubature::DynamicFactory cubature_factory
the cubature factory to use
BurgersAssemblyJobBase(const ConvVector_ &conv_vector, const Space_ &space_, const String &cubature)
Constructor.
DataType_ nu
scaling parameter for diffusive operator L (aka viscosity)
static DataType calc_sd_v_norm(const LAFEM::DenseVectorBlocked< DataType_, IndexType_, conv_dim_ > &convect)
Calculates the convection field norm for the local streamline diffusion parameter delta_T.
DataType_ DataType
the datatype we use here
DataType_ sd_v_norm
velocity norm for streamline diffusion
const Space_ & space
the test-/trial-space to be used
DataType_ sd_nu
viscosity parameter nu_S for streamline diffusion (usually equal to nu)
DataType_ theta
scaling parameter for reactive operator M
Base class for Burgers assembly tasks.
static constexpr int max_local_dofs
maximum number of local dofs
AsmTraits::CubatureRuleType cubature_rule
the cubature rule used for integration
Tiny::Vector< DataType, shape_dim > barycenter
reference element barycenter
const AsmTraits::TrafoType & trafo
the cubature factory used for integration
const DataType tol_eps
tolerance for checks == 0: sqrt(eps)
void prepare(const Index cell)
Prepares the assembly task for a cell/element.
AsmTraits::TrafoEvaluator trafo_eval
the trafo evaluator
AsmTraits1< DataType, SpaceType, TrafoTags::jac_det, SpaceTags::value|SpaceTags::grad > AsmTraits
our assembly traits
ConvVectorType::GatherAxpy gather_conv
convection vector gather-axpy object
DataType local_delta
local delta for streamline diffusion
const DataType sd_nu
viscosity parameter nu_S for streamline diffusion (usually equal to nu)
const SpaceType & space
the test-/trial-space to be used
void finish()
Finishes the assembly on the current cell.
const DataType sd_delta
scaling parameter for streamline diffusion stabilization operator S
Tiny::Matrix< DataType, conv_dim, conv_dim > loc_grad_v
local convection field gradient
const DataType nu
scaling parameter for diffusive operator L (aka viscosity)
const Job_ & job
a reference to our job object
const DataType theta
scaling parameter for reactive operator M
const bool deformation
specifies whether to use the deformation tensor
Tiny::Vector< DataType, max_local_dofs > streamdiff_coeffs
local streamline diffusion coefficients
BurgersAssemblyTaskBase(const Job_ &job_)
Constructor.
DataType_ DataType
the data type to be used by the assembly
static constexpr bool need_combine
this task doesn't need to combine
static constexpr int conv_dim
the convection vector dimension
Tiny::Vector< Tiny::Vector< DataType, conv_dim >, max_local_dofs > local_conv_dofs
local convection field dofs
int num_local_dofs
actual number of local dofs on current element
AsmTraits::SpaceEvaluator space_eval
the space evaluator
Tiny::Vector< DataType, conv_dim > loc_v
local convection field value
static constexpr bool need_scatter
this task needs to scatter
bool need_loc_v
enable computation of certain quantities
AsmTraits::DofMapping dof_mapping
the space dof-mapping
AsmTraits::TrafoEvalData trafo_data
the trafo evaluation data
const DataType frechet_beta
scaling parameter for Frechet derivative of convective operator K'
DataType local_h
local mesh width for streamline diffusion
const DataType sd_v_norm
velocity norm for streamline diffusion
void prepare_point(int cubature_point)
Prepares the task for a cubature point.
const bool need_diff
keep track what we need to assemble
Job_::ConvVectorType ConvVectorType
the convection vector type
static constexpr int shape_dim
the shape dimension
AsmTraits::SpaceEvalData space_data
the space evaluation data
const DataType beta
scaling parameter for convective operator K
Base class for blocked Burgers assembly tasks.
BurgersBlockedAssemblyTaskBase(const Job_ &job_)
Constructor.
void assemble()
Performs the assembly of the local matrix.
void assemble_streamline_diffusion(const DataType weight)
Assembles the Streamline Diffusion stabilization operator in a single cubature point.
Tiny::Matrix< DataType_, block_size_, block_size_ > MatrixValue
our local matrix data
void assemble_burgers_point(const DataType weight)
Assembles the Burger operator in a single cubature point.
int num_local_dofs
actual number of local dofs on current element
MatrixType::ScatterAxpy scatter_matrix
matrix scatter-axpy object
Task(const BurgersBlockedMatrixAssemblyJob &job_)
constructor
Burgers assembly job for block matrix.
BurgersBlockedMatrixAssemblyJob(Matrix_ &matrix_, const ConvVector_ &conv_vector, const Space_ &space_, const String &cubature_name)
Constructor.
MatrixType & matrix
the matrix to be assembled
BurgersAssemblyJobBase< DataType, Space_, ConvVector_ > BaseClass
our base class
void compute_local_vector()
Computes the local vector by multiplying the local matrix by the solution vector.
const bool sol_equals_conv
specifies whether the solution vector and the convection vector are the same object
Task(const BurgersBlockedVectorAssemblyJob &job_)
constructor
Tiny::Vector< Tiny::Vector< DataType, block_size >, max_local_dofs > local_sol_dofs
local solution vector dofs
Tiny::Vector< Tiny::Vector< DataType, block_size >, max_local_dofs > local_vector
local rhs vector dofs
VectorType::ScatterAxpy scatter_vec_rhs
rhs vector scatter-axpy object
VectorType::GatherAxpy gather_vec_sol
sol vector gather-axpy object
Burgers assembly job for block right-hand-side vector.
BurgersBlockedVectorAssemblyJob(Vector_ &rhs_vector_, const Vector_ &sol_vector_, const ConvVector_ &conv_vector, const Space_ &space_, const String &cubature_name)
Constructor.
BurgersAssemblyJobBase< DataType, Space_, ConvVector_ > BaseClass
our base class
VectorType & vector_rhs
the RHS vector to be assembled
const VectorType & vector_sol
the solution vector
Base class for scalar Burgers assembly tasks.
void assemble_burgers_point(const DataType weight)
Assembles the Burger operator in a single cubature point.
void assemble_streamline_diffusion(const DataType weight)
Assembles the Streamline Diffusion stabilization operator in a single cubature point.
BurgersScalarAssemblyTaskBase(const Job_ &job_)
Constructor.
int num_local_dofs
actual number of local dofs on current element
void assemble()
Performs the assembly of the local matrix.
Tiny::Matrix< DataType, max_local_dofs, max_local_dofs > LocalMatrixType
our local matrix data
MatrixType::ScatterAxpy scatter_matrix
matrix scatter-axpy object
Task(const BurgersScalarMatrixAssemblyJob &job_)
constructor
Burgers assembly job for scalar matrix.
BurgersScalarMatrixAssemblyJob(Matrix_ &matrix_, const ConvVector_ &conv_vector, const Space_ &space_, const String &cubature_name)
Constructor.
BurgersAssemblyJobBase< DataType, Space_, ConvVector_ > BaseClass
our base class
MatrixType & matrix
the matrix to be assembled
Tiny::Vector< DataType, max_local_dofs > local_vector
local rhs vector dofs
VectorType::GatherAxpy gather_vec_sol
sol vector gather-axpy object
void compute_local_vector()
Computes the local vector by multiplying the local matrix by the solution vector.
Tiny::Vector< DataType, max_local_dofs > local_sol_dofs
local solution vector dofs
VectorType::ScatterAxpy scatter_vec_rhs
rhs vector scatter-axpy object
Task(const BurgersScalarVectorAssemblyJob &job_)
constructor
Burgers assembly job for scalar right-hand-side vector.
BurgersScalarVectorAssemblyJob(Vector_ &rhs_vector_, const Vector_ &sol_vector_, const ConvVector_ &conv_vector, const Space_ &space_, const String &cubature_name)
Constructor.
const VectorType & vector_sol
the solution vector
VectorType & vector_rhs
the RHS vector to be assembled
BurgersAssemblyJobBase< DataType, Space_, ConvVector_ > BaseClass
our base class
DataType max(DataType x) const
Computes the maximum of a scalar variable over all processes.
Definition: gate.hpp:619
Global vector wrapper class template.
Definition: vector.hpp:68
const GateType * get_gate() const
Returns a const pointer to the internal gate of the vector.
Definition: vector.hpp:148
LocalVector_ & local()
Returns a reference to the internal local LAFEM vector object.
Definition: vector.hpp:121
Blocked Dense data vector class template.
auto elements() const -> const typename Intern::DenseVectorBlockedPerspectiveHelper< DT_, BlockSize_, perspective_ >::Type *
Retrieve a pointer to the data array.
Index size() const
The number of elements.
String class implementation.
Definition: string.hpp:46
Tiny Matrix class template.
CUDA_HOST_DEVICE Matrix & add_scalar_main_diag(DataType alpha)
Adds a value onto the matrix's main diagonal.
CUDA_HOST_DEVICE Matrix & add_outer_product(const Vector< T_, m_, snx_ > &x, const Vector< T_, n_, sny_ > &y, const DataType alpha=DataType(1))
Adds the outer product of two vectors onto the matrix.
CUDA_HOST_DEVICE Matrix & axpy(DataType alpha, const Matrix< T_, m_, n_, sma_, sna_ > &a)
Adds another scaled matrix onto this matrix.
CUDA_HOST_DEVICE void format(DataType alpha=DataType(0))
Formats the matrix.
Tiny Vector class template.
CUDA_HOST_DEVICE void format(DataType alpha=DataType(0))
Formats the vector.
CUDA_HOST_DEVICE Vector & axpy(DataType alpha, const Vector< T_, n_, snx_ > &x)
Adds another scaled vector onto this vector.
CUDA_HOST_DEVICE DataType norm_euclid() const
Computes the euclid norm of this vector.
T_ max(T_ a, T_ b)
Returns the maximum of two values.
Definition: math.hpp:137
CUDA_HOST_DEVICE T_ dot(const T_ &a, const T_ &b)
Computes the dot-product of two scalars.
FEAT namespace.
Definition: adjactor.hpp:12
@ value
specifies whether the space should supply basis function values
std::uint64_t Index
Index data type.
Reference cell traits structure.
Definition: shape.hpp:230