8#include <kernel/assembly/asm_traits.hpp>
9#include <kernel/lafem/dense_vector.hpp>
10#include <kernel/lafem/sparse_matrix_bcsr.hpp>
11#include <kernel/lafem/dense_vector_blocked.hpp>
12#include <kernel/global/vector.hpp>
61 template<
typename DataType_,
typename IndexType_,
int dim_>
123 template<
typename Space_,
typename CubatureFactory_>
128 const CubatureFactory_& cubature_factory,
129 const DataType_ scale = DataType_(1)
133 XASSERTM(matrix.
rows() == space.get_num_dofs(),
"invalid matrix dimensions");
134 XASSERTM(matrix.
columns() == space.get_num_dofs(),
"invalid matrix dimensions");
135 XASSERTM(convect.
size() == space.get_num_dofs(),
"invalid vector size");
154 const typename AsmTraits::TrafoType& trafo = space.get_trafo();
157 typename AsmTraits::TrafoEvaluator trafo_eval(trafo);
160 typename AsmTraits::SpaceEvaluator space_eval(space);
163 typename AsmTraits::DofMapping dof_mapping(space);
166 typename AsmTraits::TrafoEvalData trafo_data;
169 typename AsmTraits::SpaceEvalData space_data;
172 typename AsmTraits::CubatureRuleType cubature_rule(Cubature::ctor_factory, cubature_factory);
175 typename MatrixType::ScatterAxpy scatter_matrix(matrix);
178 typename VectorType::GatherAxpy gather_conv(convect);
181 static constexpr int max_local_dofs = AsmTraits::max_local_test_dofs;
186 LocalMatrixType local_matrix;
193 LocalVectorType local_conv_dofs;
207 for(
int i(0); i < dim_; ++i)
212 streamdiff_coeffs.
format();
218 for(
typename AsmTraits::CellIterator cell(trafo_eval.begin()); cell != trafo_eval.end(); ++cell)
221 trafo_eval.prepare(cell);
224 space_eval.prepare(trafo_eval);
227 dof_mapping.prepare(cell);
230 const int num_loc_dofs = space_eval.get_num_local_dofs();
233 local_conv_dofs.format();
234 gather_conv(local_conv_dofs, dof_mapping);
243 trafo_eval(trafo_data, barycentre);
244 space_eval(space_data, trafo_data);
248 for(
int i(0); i < num_loc_dofs; ++i)
249 mean_v.
axpy(space_data.phi[i].value, local_conv_dofs[i]);
255 if(local_norm_v > tol_eps)
258 const DataType local_h = trafo_eval.width_directed(mean_v) * local_norm_v;
261 const DataType local_re = (local_norm_v * local_h) / this->sd_nu;
269 local_matrix.format();
272 for(
int point(0); point < cubature_rule.get_num_points(); ++point)
275 trafo_eval(trafo_data, cubature_rule.get_point(point));
278 space_eval(space_data, trafo_data);
281 const DataType weight = trafo_data.jac_det * cubature_rule.get_weight(point);
284 if(need_conv || need_streamdiff)
287 for(
int i(0); i < num_loc_dofs; ++i)
290 loc_v.
axpy(space_data.phi[i].value, local_conv_dofs[i]);
293 if(need_conv_frechet)
296 for(
int i(0); i < num_loc_dofs; ++i)
306 for(
int i(0); i < num_loc_dofs; ++i)
308 streamdiff_coeffs[i] =
Tiny::dot(loc_v, space_data.phi[i].grad);
318 for(
int i(0); i < num_loc_dofs; ++i)
321 for(
int j(0); j < num_loc_dofs; ++j)
327 local_matrix[i][j].add_scalar_main_diag(
value);
335 for(
int i(0); i < num_loc_dofs; ++i)
338 for(
int j(0); j < num_loc_dofs; ++j)
346 local_matrix[i][j].add_scalar_main_diag(
value);
349 local_matrix[i][j].add_outer_product(space_data.phi[j].grad, space_data.phi[i].grad,
nu * weight);
358 for(
int i(0); i < num_loc_dofs; ++i)
361 for(
int j(0); j < num_loc_dofs; ++j)
367 local_matrix[i][j].add_scalar_main_diag(
value);
373 if(need_conv_frechet)
376 for(
int i(0); i < num_loc_dofs; ++i)
379 for(
int j(0); j < num_loc_dofs; ++j)
386 local_matrix[i][j].axpy(
value, loc_grad_v);
395 for(
int i(0); i < num_loc_dofs; ++i)
398 for(
int j(0); j < num_loc_dofs; ++j)
401 const DataType value =
theta * weight * space_data.phi[i].value * space_data.phi[j].value;
404 local_matrix[i][j].add_scalar_main_diag(
value);
410 if(need_streamdiff && (local_delta > tol_eps))
413 for(
int i(0); i < num_loc_dofs; ++i)
416 for(
int j(0); j < num_loc_dofs; ++j)
419 const DataType value = local_delta * weight * streamdiff_coeffs[i] * streamdiff_coeffs[j];
422 local_matrix[i][j].add_scalar_main_diag(
value);
430 scatter_matrix(local_matrix, dof_mapping, dof_mapping, scale);
433 dof_mapping.finish();
459 template<
typename Matrix_,
typename Space_,
typename CubatureFactory_>
464 const CubatureFactory_& cubature_factory,
465 const DataType_ scale = DataType_(1)
469 XASSERTM(matrix.rows() == space.get_num_dofs(),
"invalid matrix dimensions");
470 XASSERTM(matrix.columns() == space.get_num_dofs(),
"invalid matrix dimensions");
471 XASSERTM(convect.
size() == space.get_num_dofs(),
"invalid vector size");
474 typedef Matrix_ MatrixType;
493 const typename AsmTraits::TrafoType& trafo = space.get_trafo();
496 typename AsmTraits::TrafoEvaluator trafo_eval(trafo);
499 typename AsmTraits::SpaceEvaluator space_eval(space);
502 typename AsmTraits::DofMapping dof_mapping(space);
505 typename AsmTraits::TrafoEvalData trafo_data;
508 typename AsmTraits::SpaceEvalData space_data;
511 typename AsmTraits::CubatureRuleType cubature_rule(Cubature::ctor_factory, cubature_factory);
514 typename MatrixType::ScatterAxpy scatter_matrix(matrix);
517 typename VectorType::GatherAxpy gather_conv(convect);
520 static constexpr int max_local_dofs = AsmTraits::max_local_test_dofs;
523 typedef typename AsmTraits::template TLocalMatrix<DataType> LocalMatrixType;
524 LocalMatrixType local_matrix;
531 LocalVectorType local_conv_dofs;
540 for(
int i(0); i < dim_; ++i)
545 streamdiff_coeffs.
format();
551 for(
typename AsmTraits::CellIterator cell(trafo_eval.begin()); cell != trafo_eval.end(); ++cell)
554 trafo_eval.prepare(cell);
557 space_eval.prepare(trafo_eval);
560 dof_mapping.prepare(cell);
563 const int num_loc_dofs = space_eval.get_num_local_dofs();
566 local_conv_dofs.format();
567 gather_conv(local_conv_dofs, dof_mapping);
576 trafo_eval(trafo_data, barycentre);
577 space_eval(space_data, trafo_data);
581 for(
int i(0); i < num_loc_dofs; ++i)
582 mean_v.
axpy(space_data.phi[i].value, local_conv_dofs[i]);
588 if(local_norm_v > tol_eps)
591 const DataType local_h = trafo_eval.width_directed(mean_v) * local_norm_v;
594 const DataType local_re = (local_norm_v * local_h) / this->sd_nu;
602 local_matrix.format();
605 for(
int point(0); point < cubature_rule.get_num_points(); ++point)
608 trafo_eval(trafo_data, cubature_rule.get_point(point));
611 space_eval(space_data, trafo_data);
614 const DataType weight = trafo_data.jac_det * cubature_rule.get_weight(point);
617 if(need_conv || need_streamdiff)
620 for(
int i(0); i < num_loc_dofs; ++i)
623 loc_v.
axpy(space_data.phi[i].value, local_conv_dofs[i]);
630 for(
int i(0); i < num_loc_dofs; ++i)
632 streamdiff_coeffs[i] =
Tiny::dot(loc_v, space_data.phi[i].grad);
640 for(
int i(0); i < num_loc_dofs; ++i)
643 for(
int j(0); j < num_loc_dofs; ++j)
645 local_matrix[i][j] +=
nu * weight *
Tiny::dot(space_data.phi[i].grad, space_data.phi[j].grad);
654 for(
int i(0); i < num_loc_dofs; ++i)
657 for(
int j(0); j < num_loc_dofs; ++j)
659 local_matrix[i][j] +=
beta * weight * space_data.phi[i].value *
Tiny::dot(loc_v, space_data.phi[j].grad);
668 for(
int i(0); i < num_loc_dofs; ++i)
671 for(
int j(0); j < num_loc_dofs; ++j)
673 local_matrix[i][j] +=
theta * weight * space_data.phi[i].value * space_data.phi[j].value;
679 if(need_streamdiff && (local_delta > tol_eps))
682 for(
int i(0); i < num_loc_dofs; ++i)
685 for(
int j(0); j < num_loc_dofs; ++j)
687 local_matrix[i][j] += local_delta * weight * streamdiff_coeffs[i] * streamdiff_coeffs[j];
696 scatter_matrix(local_matrix, dof_mapping, dof_mapping, scale);
699 dof_mapping.finish();
728 template<
typename Space_,
typename CubatureFactory_>
734 const CubatureFactory_& cubature_factory,
735 const DataType_ scale = DataType_(1)
739 XASSERTM(vector.
size() == space.get_num_dofs(),
"invalid vector size");
740 XASSERTM(convect.
size() == space.get_num_dofs(),
"invalid vector size");
754 const typename AsmTraits::TrafoType& trafo = space.get_trafo();
757 typename AsmTraits::TrafoEvaluator trafo_eval(trafo);
760 typename AsmTraits::SpaceEvaluator space_eval(space);
763 typename AsmTraits::DofMapping dof_mapping(space);
766 typename AsmTraits::TrafoEvalData trafo_data;
769 typename AsmTraits::SpaceEvalData space_data;
772 typename AsmTraits::CubatureRuleType cubature_rule(Cubature::ctor_factory, cubature_factory);
775 typename VectorType::ScatterAxpy scatter_vector(vector);
778 typename VectorType::GatherAxpy gather_conv(convect);
781 typename VectorType::GatherAxpy gather_prim(primal);
784 static constexpr int max_local_dofs = AsmTraits::max_local_test_dofs;
789 LocalVectorType local_vector;
792 LocalVectorType local_conv_dofs;
795 LocalVectorType local_prim_dofs;
807 for(
typename AsmTraits::CellIterator cell(trafo_eval.begin()); cell != trafo_eval.end(); ++cell)
810 trafo_eval.prepare(cell);
813 space_eval.prepare(trafo_eval);
816 dof_mapping.prepare(cell);
819 const int num_loc_dofs = space_eval.get_num_local_dofs();
822 local_conv_dofs.format();
823 gather_conv(local_conv_dofs, dof_mapping);
826 local_prim_dofs.format();
827 gather_prim(local_prim_dofs, dof_mapping);
830 local_vector.format();
833 for(
int point(0); point < cubature_rule.get_num_points(); ++point)
836 trafo_eval(trafo_data, cubature_rule.get_point(point));
839 space_eval(space_data, trafo_data);
842 const DataType weight = trafo_data.jac_det * cubature_rule.get_weight(point);
848 for(
int i(0); i < num_loc_dofs; ++i)
851 loc_v.
axpy(space_data.phi[i].value, local_conv_dofs[i]);
870 for(
int i(0); i < num_loc_dofs; ++i)
873 for(
int j(0); j < num_loc_dofs; ++j)
879 local_vector[i].axpy(
value, local_prim_dofs[j]);
888 for(
int i(0); i < num_loc_dofs; ++i)
891 for(
int j(0); j < num_loc_dofs; ++j)
894 const DataType value1 =
nu * weight *
Tiny::dot(space_data.phi[i].grad, space_data.phi[j].grad);
897 const DataType value2 =
nu * weight *
Tiny::dot(local_prim_dofs[j], space_data.phi[i].grad);
900 local_vector[i].axpy(value1, local_prim_dofs[j]);
901 local_vector[i].axpy(value2, space_data.phi[j].grad);
910 for(
int i(0); i < num_loc_dofs; ++i)
913 for(
int j(0); j < num_loc_dofs; ++j)
919 local_vector[i].axpy(
value, local_prim_dofs[j]);
928 for(
int i(0); i < num_loc_dofs; ++i)
931 for(
int j(0); j < num_loc_dofs; ++j)
934 const DataType value =
theta * weight * space_data.phi[i].value * space_data.phi[j].value;
937 local_vector[i].axpy(
value, local_prim_dofs[j]);
946 scatter_vector(local_vector, dof_mapping, scale);
949 dof_mapping.finish();
965 const auto* vals = convect.
elements();
981 template<
typename LocalVector_,
typename Mirror_>
985 const auto* gate = convect.
get_gate();
987 this->sd_v_norm = gate->
max(this->sd_v_norm);
#define XASSERTM(expr, msg)
Assertion macro definition with custom message.
Common single-space assembly traits class template.
Burgers operator assembly class.
void assemble_scalar_matrix(Matrix_ &matrix, const LAFEM::DenseVectorBlocked< DataType_, IndexType_, dim_ > &convect, const Space_ &space, const CubatureFactory_ &cubature_factory, const DataType_ scale=DataType_(1)) const
Assembles the Burgers operator into a scalar matrix.
DataType_ sd_delta
scaling parameter for streamline diffusion stabilization operator S
DataType_ DataType
the datatype we use here
DataType_ frechet_beta
scaling parameter for Frechet derivative of convective operator K'
BurgersAssembler()
default constructor
DataType_ beta
scaling parameter for convective operator K
DataType_ nu
scaling parameter for diffusive operator L (aka viscosity)
void set_sd_v_norm(const LAFEM::DenseVectorBlocked< DataType_, IndexType_, dim_ > &convect)
Sets the convection field norm for the local streamline diffusion parameter delta_T.
void assemble_vector(LAFEM::DenseVectorBlocked< DataType_, IndexType_, dim_ > &vector, const LAFEM::DenseVectorBlocked< DataType_, IndexType_, dim_ > &convect, const LAFEM::DenseVectorBlocked< DataType_, IndexType_, dim_ > &primal, const Space_ &space, const CubatureFactory_ &cubature_factory, const DataType_ scale=DataType_(1)) const
Assembles the Burgers operator into a vector.
DataType_ sd_v_norm
velocity norm for streamline diffusion
void assemble_matrix(LAFEM::SparseMatrixBCSR< DataType_, IndexType_, dim_, dim_ > &matrix, const LAFEM::DenseVectorBlocked< DataType_, IndexType_, dim_ > &convect, const Space_ &space, const CubatureFactory_ &cubature_factory, const DataType_ scale=DataType_(1)) const
Assembles the Burgers operator into a matrix.
DataType_ theta
scaling parameter for reactive operator M
DataType_ sd_nu
viscosity parameter nu_S for streamline diffusion (usually equal to nu)
void set_sd_v_norm(const Global::Vector< LocalVector_, Mirror_ > &convect)
Sets the convection field norm for the streamline diffusion parameter delta_T.
bool deformation
specifies whether to use the deformation tensor
DataType max(DataType x) const
Computes the maximum of a scalar variable over all processes.
Global vector wrapper class template.
const GateType * get_gate() const
Returns a const pointer to the internal gate of the vector.
LocalVector_ & local()
Returns a reference to the internal local LAFEM vector object.
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.
CSR based blocked sparse matrix.
Index rows() const
Retrieve matrix row count.
Index columns() const
Retrieve matrix column count.
Tiny Matrix class template.
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 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_ sqrt(T_ x)
Returns the square-root of a value.
T_ abs(T_ x)
Returns the absolute value.
T_ max(T_ a, T_ b)
Returns the maximum of two values.
CUDA_HOST_DEVICE T_ dot(const T_ &a, const T_ &b)
Computes the dot-product of two scalars.
@ value
specifies whether the space should supply basis function values
std::uint64_t Index
Index data type.
Reference cell traits structure.