9#include <kernel/lafem/dense_vector.hpp>
10#include <kernel/lafem/dense_vector_blocked.hpp>
11#include <kernel/lafem/filter_chain.hpp>
12#include <kernel/lafem/filter_sequence.hpp>
13#include <kernel/lafem/sparse_matrix_csr.hpp>
14#include <kernel/lafem/sparse_matrix_bcsr.hpp>
15#include <kernel/lafem/sparse_matrix_bwrappedcsr.hpp>
16#include <kernel/lafem/saddle_point_matrix.hpp>
17#include <kernel/lafem/unit_filter.hpp>
18#include <kernel/lafem/unit_filter_blocked.hpp>
19#include <kernel/lafem/mean_filter.hpp>
20#include <kernel/lafem/none_filter.hpp>
21#include <kernel/lafem/slip_filter.hpp>
22#include <kernel/lafem/tuple_filter.hpp>
23#include <kernel/lafem/tuple_mirror.hpp>
24#include <kernel/lafem/tuple_diag_matrix.hpp>
25#include <kernel/lafem/transfer.hpp>
26#include <kernel/assembly/mirror_assembler.hpp>
27#include <kernel/assembly/symbolic_assembler.hpp>
28#include <kernel/assembly/bilinear_operator_assembler.hpp>
29#include <kernel/assembly/domain_assembler_basic_jobs.hpp>
30#include <kernel/assembly/common_operators.hpp>
31#include <kernel/assembly/gpdv_assembler.hpp>
32#include <kernel/assembly/grid_transfer.hpp>
33#include <kernel/assembly/mean_filter_assembler.hpp>
34#include <kernel/assembly/unit_filter_assembler.hpp>
35#include <kernel/global/gate.hpp>
36#include <kernel/global/muxer.hpp>
37#include <kernel/global/vector.hpp>
38#include <kernel/global/matrix.hpp>
39#include <kernel/global/filter.hpp>
40#include <kernel/global/mean_filter.hpp>
41#include <kernel/global/transfer.hpp>
42#include <kernel/global/splitter.hpp>
44#include <control/domain/domain_control.hpp>
45#include <control/asm/gate_asm.hpp>
46#include <control/asm/muxer_asm.hpp>
47#include <control/asm/splitter_asm.hpp>
48#include <control/asm/transfer_asm.hpp>
49#include <control/asm/transfer_voxel_asm.hpp>
50#include <control/asm/unit_filter_asm.hpp>
51#include <control/asm/mean_filter_asm.hpp>
52#include <control/asm/slip_filter_asm.hpp>
62 typename DataType_ =
Real,
63 typename IndexType_ =
Index,
64 typename MatrixBlockA_ = LAFEM::SparseMatrixBCSR<DataType_, IndexType_, dim_, dim_>,
65 typename MatrixBlockB_ = LAFEM::SparseMatrixBCSR<DataType_, IndexType_, dim_, 1>,
66 typename MatrixBlockD_ = LAFEM::SparseMatrixBCSR<DataType_, IndexType_, 1, dim_>,
67 typename ScalarMatrix_ = LAFEM::SparseMatrixCSR<DataType_, IndexType_>,
68 typename TransferMatrixV_ = LAFEM::SparseMatrixBWrappedCSR<DataType_, IndexType_, dim_>,
69 typename TransferMatrixP_ = LAFEM::SparseMatrixCSR<DataType_, IndexType_>,
70 typename TransferMatrixS_ = LAFEM::SparseMatrixCSR<DataType_, IndexType_>
76 typedef DataType_ DataType;
77 typedef IndexType_ IndexType;
78 static constexpr int dim = dim_;
81 typedef ScalarMatrix_ LocalScalarMatrix;
82 typedef typename LocalScalarMatrix::VectorTypeL LocalScalarVector;
85 typedef MatrixBlockA_ LocalMatrixBlockA;
86 typedef MatrixBlockB_ LocalMatrixBlockB;
87 typedef MatrixBlockD_ LocalMatrixBlockD;
88 typedef LocalScalarMatrix LocalSchurMatrix;
92 typedef typename LocalMatrixBlockB::VectorTypeL LocalVeloVector;
93 typedef typename LocalMatrixBlockD::VectorTypeL LocalPresVector;
97 typedef TransferMatrixV_ LocalVeloTransferMatrix;
98 typedef TransferMatrixP_ LocalPresTransferMatrix;
99 typedef TransferMatrixS_ LocalScalarTransferMatrix;
192 matrix_s(&gate_pres, &gate_pres),
194 transfer_pres(&coarse_muxer_pres),
195 transfer_sys(&coarse_muxer_sys),
196 transfer_scalar_velo(&coarse_muxer_scalar_velo)
211 return this->matrix_sys.
local().
bytes () + this->matrix_s.
local().bytes() + transfer_sys.
bytes();
214 void compile_system_matrix()
221 void compile_system_transfer()
229 transfer_sys.compile();
232 void compile_scalar_transfer()
237 transfer_scalar_velo.compile();
240 template<
typename D_,
typename I_,
typename SMA_,
typename SMB_,
typename SMD_,
typename SM_,
typename TV_,
typename TP_>
241 void convert(
const StokesBlockedSystemLevel<dim_, D_, I_, SMA_, SMB_, SMD_, SM_, TV_, TP_> & other)
245 Asm::build_gate_tuple(this->gate_sys, this->gate_velo, this->gate_pres);
246 this->gate_scalar_velo.
convert(this->gate_velo, LocalScalarVector(this->gate_velo.
get_freqs().template size<LAFEM::Perspective::native>()));
250 Asm::build_muxer_tuple(this->coarse_muxer_sys, this->gate_sys.
get_freqs(), this->coarse_muxer_velo, this->coarse_muxer_pres);
255 Asm::build_splitter_tuple(this->base_splitter_sys, this->gate_sys.
get_freqs(), this->base_splitter_velo, this->base_splitter_pres);
259 transfer_pres.
convert(&coarse_muxer_pres,
other.transfer_pres);
260 this->compile_system_transfer();
261 this->compile_scalar_transfer();
266 matrix_s.convert(&gate_pres, &gate_pres,
other.matrix_s);
267 this->compile_system_matrix();
270 GlobalSystemVector create_global_vector_sys()
const
275 GlobalVeloVector create_global_vector_velo()
const
280 GlobalPresVector create_global_vector_pres()
const
285 template<
typename DomainLevel_>
286 void assemble_gates(
const Domain::VirtualLevel<DomainLevel_>& virt_dom_lvl)
288 Asm::asm_gate(virt_dom_lvl, virt_dom_lvl->space_velo, this->gate_velo,
true);
289 Asm::asm_gate(virt_dom_lvl, virt_dom_lvl->space_pres, this->gate_pres,
true);
290 Asm::build_gate_tuple(this->gate_sys, this->gate_velo, this->gate_pres);
291 this->gate_scalar_velo.
convert(this->gate_velo, LocalScalarVector(virt_dom_lvl->space_velo.get_num_dofs()));
294 template<
typename DomainLevel_>
295 void assemble_coarse_muxers(
const Domain::VirtualLevel<DomainLevel_>& virt_lvl_coarse)
297 Asm::asm_muxer(virt_lvl_coarse, [](
const DomainLevel_& dl){
return &dl.space_velo;}, this->
coarse_muxer_velo);
298 Asm::asm_muxer(virt_lvl_coarse, [](
const DomainLevel_& dl){
return &dl.space_pres;}, this->coarse_muxer_pres);
299 Asm::build_muxer_tuple(this->coarse_muxer_sys, this->gate_sys.
get_freqs(), this->coarse_muxer_velo, this->coarse_muxer_pres);
303 template<
typename DomainLevel_>
304 void assemble_base_splitters(
const Domain::VirtualLevel<DomainLevel_>& virt_lvl)
306 Asm::asm_splitter(virt_lvl, [](
const DomainLevel_& dl){
return &dl.space_velo;}, this->
base_splitter_velo);
307 Asm::asm_splitter(virt_lvl, [](
const DomainLevel_& dl){
return &dl.space_pres;}, this->base_splitter_pres);
308 Asm::build_splitter_tuple(this->base_splitter_sys, this->gate_sys.
get_freqs(), this->base_splitter_velo, this->base_splitter_pres);
312 template<
typename DomainLevel_>
313 void assemble_transfers(
315 const Domain::VirtualLevel<DomainLevel_>& virt_lvl_fine,
316 const Domain::VirtualLevel<DomainLevel_>& virt_lvl_coarse,
317 const String& cubature,
bool trunc_v =
false,
bool trunc_p =
false,
bool shrink =
true)
319 Asm::asm_transfer_blocked(virt_lvl_fine, virt_lvl_coarse, cubature, trunc_v, shrink,
320 [](
const DomainLevel_& dl) {
return &dl.space_velo;},
321 this->transfer_velo.
local(), this->coarse_muxer_velo, this->gate_velo, sys_lvl_coarse.gate_velo);
322 Asm::asm_transfer_scalar(virt_lvl_fine, virt_lvl_coarse, cubature, trunc_p, shrink,
323 [](
const DomainLevel_& dl) {
return &dl.space_pres;},
324 this->transfer_pres.
local(), this->coarse_muxer_pres, this->gate_pres, sys_lvl_coarse.gate_pres);
326 this->transfer_velo.compile();
327 this->transfer_pres.compile();
328 this->compile_system_transfer();
329 this->compile_scalar_transfer();
332 template<
typename DomainLevel_>
333 void assemble_transfers(
334 const Domain::VirtualLevel<DomainLevel_>& virt_lvl_fine,
335 const Domain::VirtualLevel<DomainLevel_>& virt_lvl_coarse,
336 const String& cubature,
bool trunc_v =
false,
bool trunc_p =
false,
bool shrink =
true)
339 XASSERT(!virt_lvl_coarse.is_parent());
341 Asm::asm_transfer_blocked(virt_lvl_fine, virt_lvl_coarse, cubature, trunc_v, shrink,
342 [](
const DomainLevel_& dl) {
return &dl.space_velo;},
343 this->transfer_velo.
local(), this->coarse_muxer_velo, this->gate_velo, this->gate_velo);
344 Asm::asm_transfer_scalar(virt_lvl_fine, virt_lvl_coarse, cubature, trunc_p, shrink,
345 [](
const DomainLevel_& dl) {
return &dl.space_pres;},
346 this->transfer_pres.
local(), this->coarse_muxer_pres, this->gate_pres, this->gate_pres);
348 this->transfer_velo.compile();
349 this->transfer_pres.compile();
350 this->compile_system_transfer();
351 this->compile_scalar_transfer();
354 template<
typename DomainLevel_,
template<
typename>
class SlagDomainLevelWrapper_>
355 void assemble_transfers_voxel(
357 const Domain::VirtualLevel<SlagDomainLevelWrapper_<DomainLevel_>>& virt_lvl_fine,
358 const Domain::VirtualLevel<SlagDomainLevelWrapper_<DomainLevel_>>& virt_lvl_coarse,
359 const String& cubature,
bool trunc_v =
false,
bool trunc_p =
false,
bool shrink =
true)
361 Asm::asm_transfer_voxel_blocked(virt_lvl_fine, virt_lvl_coarse, cubature, trunc_v, shrink,
362 [](
const DomainLevel_& dl) {
return &dl.space_velo;},
363 this->transfer_velo.
local(), this->coarse_muxer_velo, this->gate_velo, sys_lvl_coarse.gate_velo);
364 Asm::asm_transfer_voxel_scalar(virt_lvl_fine, virt_lvl_coarse, cubature, trunc_p, shrink,
365 [](
const DomainLevel_& dl) {
return &dl.space_pres;},
366 this->transfer_pres.
local(), this->coarse_muxer_pres, this->gate_pres, sys_lvl_coarse.gate_pres);
368 this->transfer_velo.compile();
369 this->transfer_pres.compile();
370 this->compile_system_transfer();
371 this->compile_scalar_transfer();
374 template<
typename DomainLevel_,
template<
typename>
class SlagDomainLevelWrapper_>
375 void assemble_transfers_voxel(
376 const Domain::VirtualLevel<SlagDomainLevelWrapper_<DomainLevel_>>& virt_lvl_fine,
377 const Domain::VirtualLevel<SlagDomainLevelWrapper_<DomainLevel_>>& virt_lvl_coarse,
378 const String& cubature,
bool trunc_v =
false,
bool trunc_p =
false,
bool shrink =
true)
381 XASSERT(!virt_lvl_coarse.is_parent());
383 Asm::asm_transfer_voxel_blocked(virt_lvl_fine, virt_lvl_coarse, cubature, trunc_v, shrink,
384 [](
const DomainLevel_& dl) {
return &dl.space_velo;},
385 this->transfer_velo.
local(), this->coarse_muxer_velo, this->gate_velo, this->gate_velo);
386 Asm::asm_transfer_voxel_scalar(virt_lvl_fine, virt_lvl_coarse, cubature, trunc_p, shrink,
387 [](
const DomainLevel_& dl) {
return &dl.space_pres;},
388 this->transfer_pres.
local(), this->coarse_muxer_pres, this->gate_pres, this->gate_pres);
390 this->transfer_velo.compile();
391 this->transfer_pres.compile();
392 this->compile_system_transfer();
393 this->compile_scalar_transfer();
396 void clear_transfers()
398 this->transfer_velo.
local() = LocalVeloTransfer();
399 this->transfer_pres.
local() = LocalPresTransfer();
400 this->transfer_scalar_velo.
local() = LocalScalarTransfer();
401 this->compile_system_transfer();
402 this->compile_scalar_transfer();
405 template<
typename SpaceVelo_,
typename SpacePres_,
typename Cubature_>
406 void assemble_grad_div_matrices(
const SpaceVelo_& space_velo,
const SpacePres_& space_pres,
const Cubature_& cubature)
411 template<
typename Trafo_,
typename SpaceVelo_,
typename SpacePres_>
412 void assemble_grad_div_matrices(Assembly::DomainAssembler<Trafo_>& dom_asm,
413 const SpaceVelo_& space_velo,
const SpacePres_& space_pres,
const String& cubature)
416 if(this->matrix_b.
local().empty())
420 this->matrix_b.
local().format();
421 Assembly::Common::GradientTestOperatorBlocked<dim_> grad_op;
423 dom_asm, this->matrix_b.
local(), grad_op, space_velo, space_pres, cubature, -DataType(1));
426 this->matrix_d.
local().transpose(this->matrix_b.
local());
429 template<
typename Trafo_,
typename SpaceVelo_,
typename SpacePres_,
typename AsmDT_>
430 void assemble_grad_div_matrices_high_prec(Assembly::DomainAssembler<Trafo_>& dom_asm,
431 const SpaceVelo_& space_velo,
const SpacePres_& space_pres,
const String& cubature,
const AsmDT_ asm_weight)
434 if(this->matrix_b.
local().empty())
438 this->matrix_b.
local().format();
441 typename LocalMatrixBlockB::template ContainerTypeByDI<AsmDT_, IndexType_> tmp_mat_b;
442 tmp_mat_b.convert(this->matrix_b.
local());
445 Assembly::Common::GradientTestOperatorBlocked<dim_> grad_op;
447 dom_asm, tmp_mat_b, grad_op, space_velo, space_pres, cubature, -asm_weight);
450 this->matrix_b.
local().convert(tmp_mat_b);
454 this->matrix_d.
local().transpose(this->matrix_b.
local());
457 template<
typename SpaceVelo_>
458 void assemble_velo_struct(
const SpaceVelo_& space_velo)
462 this->matrix_a.
local().format();
465 template<
typename SpacePres_>
466 void assemble_pres_struct(
const SpacePres_& space_pres)
470 this->matrix_s.
local().format();
473 void compile_local_matrix_sys_type1()
484 typename DataType_ =
Real,
485 typename IndexType_ =
Index,
486 typename MatrixBlockA_ = LAFEM::SparseMatrixBCSR<DataType_, IndexType_, dim_, dim_>,
487 typename MatrixBlockB_ = LAFEM::SparseMatrixBCSR<DataType_, IndexType_, dim_, 1>,
488 typename MatrixBlockD_ = LAFEM::SparseMatrixBCSR<DataType_, IndexType_, 1, dim_>,
489 typename ScalarMatrix_ = LAFEM::SparseMatrixCSR<DataType_, IndexType_>,
490 typename TransferMatrixV_ = LAFEM::SparseMatrixBWrappedCSR<DataType_, IndexType_, dim_>,
491 typename TransferMatrixP_ = LAFEM::SparseMatrixCSR<DataType_, IndexType_>
495 MatrixBlockA_, MatrixBlockB_, MatrixBlockD_, ScalarMatrix_, TransferMatrixV_, TransferMatrixP_>
499 MatrixBlockA_, MatrixBlockB_, MatrixBlockD_, ScalarMatrix_, TransferMatrixV_, TransferMatrixP_>
BaseClass;
522 void compile_system_filter()
528 void compile_local_matrix_sys_type1()
530 BaseClass::compile_local_matrix_sys_type1();
541 typename DataType_ =
Real,
542 typename IndexType_ =
Index,
543 typename MatrixBlockA_ = LAFEM::SparseMatrixBCSR<DataType_, IndexType_, dim_, dim_>,
544 typename MatrixBlockB_ = LAFEM::SparseMatrixBCSR<DataType_, IndexType_, dim_, 1>,
545 typename MatrixBlockD_ = LAFEM::SparseMatrixBCSR<DataType_, IndexType_, 1, dim_>,
546 typename ScalarMatrix_ = LAFEM::SparseMatrixCSR<DataType_, IndexType_>,
547 typename TransferMatrixV_ = LAFEM::SparseMatrixBWrappedCSR<DataType_, IndexType_, dim_>,
548 typename TransferMatrixP_ = LAFEM::SparseMatrixCSR<DataType_, IndexType_>
552 MatrixBlockA_, MatrixBlockB_, MatrixBlockD_, ScalarMatrix_, TransferMatrixV_, TransferMatrixP_>
556 MatrixBlockA_, MatrixBlockB_, MatrixBlockD_, ScalarMatrix_, TransferMatrixV_, TransferMatrixP_>
BaseClass;
581 void compile_system_filter()
587 void compile_local_matrix_sys_type1()
589 BaseClass::compile_local_matrix_sys_type1();
592 this->filter_velo.local().template at<1>().
filter_mat(this->local_matrix_sys_type1.
block_a());
593 this->filter_velo.local().template at<1>().filter_offdiag_row_mat(this->local_matrix_sys_type1.
block_b());
605 typename DataType_ =
Real,
606 typename IndexType_ =
Index,
607 typename MatrixBlockA_ = LAFEM::SparseMatrixBCSR<DataType_, IndexType_, dim_, dim_>,
608 typename MatrixBlockB_ = LAFEM::SparseMatrixBCSR<DataType_, IndexType_, dim_, 1>,
609 typename MatrixBlockD_ = LAFEM::SparseMatrixBCSR<DataType_, IndexType_, 1, dim_>,
610 typename ScalarMatrix_ = LAFEM::SparseMatrixCSR<DataType_, IndexType_>,
611 typename TransferMatrixV_ = LAFEM::SparseMatrixBWrappedCSR<DataType_, IndexType_, dim_>,
612 typename TransferMatrixP_ = LAFEM::SparseMatrixCSR<DataType_, IndexType_>
616 MatrixBlockA_, MatrixBlockB_, MatrixBlockD_, ScalarMatrix_, TransferMatrixV_, TransferMatrixP_>
620 MatrixBlockA_, MatrixBlockB_, MatrixBlockD_, ScalarMatrix_, TransferMatrixV_, TransferMatrixP_>
BaseClass;
640 return this->matrix_sys.
local().
bytes () + this->matrix_s.
local().bytes() + filter_sys.local().
bytes();
643 void compile_system_filter()
649 template<
typename D_,
typename I_,
typename SM_>
650 void convert(
const StokesBlockedUnitVeloMeanPresSystemLevel<dim_, D_, I_, SM_> & other)
652 BaseClass::convert(other);
653 filter_velo.convert(other.filter_velo);
654 filter_pres.convert(other.filter_pres);
656 compile_system_filter();
659 template<
typename SpacePres_,
typename Cubature_>
660 void assemble_pressure_mean_filter(
const SpacePres_& space_pres,
const Cubature_& cubature)
663 LocalPresFilter& fil_loc_p = this->filter_pres.local();
666 typename BaseClass::GlobalPresVector vec_glob_v(&this->gate_pres), vec_glob_w(&this->gate_pres);
669 typename BaseClass::LocalPresVector& vec_loc_v = vec_glob_v.local();
670 typename BaseClass::LocalPresVector& vec_loc_w = vec_glob_w.local();
673 typename BaseClass::LocalPresVector& vec_loc_f = this->gate_pres.
_freqs;
683 fil_loc_p = LocalPresFilter(vec_loc_v.clone(), vec_loc_w.clone(), vec_loc_f.clone(), this->gate_pres.get_comm());
686 void compile_local_matrix_sys_type1()
688 BaseClass::compile_local_matrix_sys_type1();
699 typename DataType_ =
Real,
700 typename IndexType_ =
Index,
701 typename MatrixBlockA_ = LAFEM::SparseMatrixBCSR<DataType_, IndexType_, dim_, dim_>,
702 typename MatrixBlockB_ = LAFEM::SparseMatrixBCSR<DataType_, IndexType_, dim_, 1>,
703 typename MatrixBlockD_ = LAFEM::SparseMatrixBCSR<DataType_, IndexType_, 1, dim_>,
704 typename ScalarMatrix_ = LAFEM::SparseMatrixCSR<DataType_, IndexType_>,
705 typename TransferMatrixV_ = LAFEM::SparseMatrixBWrappedCSR<DataType_, IndexType_, dim_>,
706 typename TransferMatrixP_ = LAFEM::SparseMatrixCSR<DataType_, IndexType_>
710 MatrixBlockA_, MatrixBlockB_, MatrixBlockD_, ScalarMatrix_, TransferMatrixV_, TransferMatrixP_>
714 MatrixBlockA_, MatrixBlockB_, MatrixBlockD_, ScalarMatrix_, TransferMatrixV_, TransferMatrixP_>
BaseClass;
739 void compile_system_filter()
745 template<
typename D_,
typename I_,
typename SM_>
746 void convert(
const StokesBlockedUnitVeloMeanPresSystemLevel<dim_, D_, I_, SM_> & other)
748 BaseClass::convert(other);
749 filter_velo.convert(other.filter_velo);
750 filter_pres.convert(other.filter_pres);
752 compile_system_filter();
755 template<
typename SpacePres_,
typename Cubature_>
756 void assemble_global_filters(
const SpacePres_& space_pres,
const Cubature_& cubature)
759 LocalPresFilter& fil_loc_p = this->filter_pres.local();
762 typename BaseClass::GlobalPresVector vec_glob_v(&this->gate_pres), vec_glob_w(&this->gate_pres);
765 typename BaseClass::LocalPresVector& vec_loc_v = vec_glob_v.local();
766 typename BaseClass::LocalPresVector& vec_loc_w = vec_glob_w.local();
769 typename BaseClass::LocalPresVector& vec_loc_f = this->gate_pres.
_freqs;
779 fil_loc_p = LocalPresFilter(vec_loc_v.clone(), vec_loc_w.clone(), vec_loc_f.clone(), this->gate_pres.get_comm());
782 Asm::sync_slip_filter(this->gate_velo, filter_velo.local().template at<0>());
785 void compile_local_matrix_sys_type1()
787 BaseClass::compile_local_matrix_sys_type1();
790 this->filter_velo.local().template at<1>().
filter_mat(this->local_matrix_sys_type1.
block_a());
791 this->filter_velo.local().template at<1>().filter_offdiag_row_mat(this->local_matrix_sys_type1.
block_b());
803 typename DataType_ =
Real,
804 typename IndexType_ =
Index,
805 typename MatrixBlockA_ = LAFEM::SparseMatrixBCSR<DataType_, IndexType_, dim_, dim_>,
806 typename MatrixBlockB_ = LAFEM::SparseMatrixBCSR<DataType_, IndexType_, dim_, 1>,
807 typename MatrixBlockD_ = LAFEM::SparseMatrixBCSR<DataType_, IndexType_, 1, dim_>,
808 typename ScalarMatrix_ = LAFEM::SparseMatrixCSR<DataType_, IndexType_>,
809 typename TransferMatrixV_ = LAFEM::SparseMatrixBWrappedCSR<DataType_, IndexType_, dim_>,
810 typename TransferMatrixP_ = LAFEM::SparseMatrixCSR<DataType_, IndexType_>
814 MatrixBlockA_, MatrixBlockB_, MatrixBlockD_, ScalarMatrix_, TransferMatrixV_, TransferMatrixP_>
818 MatrixBlockA_, MatrixBlockB_, MatrixBlockD_, ScalarMatrix_, TransferMatrixV_, TransferMatrixP_>
BaseClass;
851 void compile_system_filter()
857 template<
typename D_,
typename I_,
typename SM_>
858 void convert(
const StokesBlockedCombinedSystemLevel<dim_, D_, I_, SM_> & other)
860 BaseClass::convert(other);
861 filter_velo.convert(other.filter_velo);
862 filter_pres.convert(other.filter_pres);
864 compile_system_filter();
867 LocalVeloSlipFilterSeq& get_local_velo_slip_filter_seq()
869 return this->filter_velo.local().template at<0>();
872 LocalVeloUnitFilterSeq& get_local_velo_unit_filter_seq()
874 return this->filter_velo.local().template at<1>();
877 LocalPresMeanFilter& get_local_pres_mean_filter()
879 return this->filter_pres.local().template at<0>();
882 LocalPresUnitFilter& get_local_pres_unit_filter()
884 return this->filter_pres.local().template at<1>();
887 template<
typename DomainLevel_,
typename Space_>
888 void assemble_velocity_unit_filter(
const DomainLevel_& dom_level,
const Space_& space,
const String& name,
const String& mesh_parts)
890 auto& loc_filter = get_local_velo_unit_filter_seq().find_or_add(name);
892 Asm::asm_unit_filter_blocked_homogeneous(loc_filter, dom_level, space, mesh_parts);
895 template<
typename DomainLevel_,
typename Space_,
typename Function_>
896 void assemble_velocity_unit_filter(
const DomainLevel_& dom_level,
const Space_& space,
const String& name,
const String& mesh_parts,
const Function_& function)
898 auto& loc_filter = get_local_velo_unit_filter_seq().find_or_add(name);
900 Asm::asm_unit_filter_blocked(loc_filter, dom_level, space, mesh_parts, function);
903 template<
typename DomainLevel_,
typename Space_>
904 void assemble_velocity_slip_filter(
const DomainLevel_& dom_level,
const Space_& space,
const String& name,
const String& mesh_parts)
906 auto& loc_filter = get_local_velo_slip_filter_seq().find_or_add(name);
908 Asm::asm_slip_filter(loc_filter, dom_level, space, mesh_parts);
909 Asm::sync_slip_filter(this->gate_velo, loc_filter);
912 template<
typename SpacePres_>
913 void assemble_pressure_mean_filter(
const SpacePres_& space_pres,
const String& cubature_name)
915 this->get_local_pres_mean_filter() = Asm::asm_mean_filter(this->gate_pres, space_pres, cubature_name);
918 void sync_velocity_slip_filters()
920 for(
auto& it : get_local_velo_slip_filter_seq())
921 Asm::sync_slip_filter(this->gate_velo, it.second);
924 void compile_local_matrix_sys_type1()
926 BaseClass::compile_local_matrix_sys_type1();
927 for(
const auto& loc_fil : get_local_velo_unit_filter_seq())
929 loc_fil.second.filter_mat(this->local_matrix_sys_type1.
block_a());
930 loc_fil.second.filter_offdiag_row_mat(this->local_matrix_sys_type1.
block_b());
#define XASSERT(expr)
Assertion macro definition.
static void assemble(LAFEM::SparseMatrixBCSR< DataType_, IndexType_, dim_, 1 > &matrix_b, LAFEM::SparseMatrixBCSR< DataType_, IndexType_, 1, dim_ > &matrix_d, const SpaceVelo_ &space_velo, const SpacePres_ &space_pres, const String &cubature_name, const DataType_ scale_b=-DataType_(1), const DataType_ scale_d=-DataType_(1))
Assembles the B and D matrices.
static void assemble(LAFEM::DenseVector< DataType_, IndexType_ > &vec_prim, LAFEM::DenseVector< DataType_, IndexType_ > &vec_dual, const Space_ &space, const String &cubature_name)
Assembles an integral mean filter.
static void assemble_matrix_std1(MatrixType_ &matrix, const Space_ &space)
Assembles a standard matrix structure from a single space.
static void assemble_matrix_std2(MatrixType_ &matrix, const TestSpace_ &test_space, const TrialSpace_ &trial_space)
Assembles a standard matrix structure from a test-trial-space pair.
Stokes blocked System level combining all supported types of boundary conditions.
std::size_t bytes() const
Returns the total amount of bytes allocated.
std::size_t bytes() const
Returns the total amount of bytes allocated.
std::size_t bytes() const
Returns the total amount of bytes allocated.
GlobalSystemMatrix matrix_sys
our global system matrix
VeloSplitter base_splitter_velo
our base-mesh multiplexer
LocalSystemMatrix local_matrix_sys_type1
local type-1 system matrix
std::size_t bytes() const
Returns the total amount of bytes allocated.
GlobalVeloTransfer transfer_velo
our global transfer operator
StokesBlockedSystemLevel()
CTOR.
VeloGate gate_velo
our system gate
VeloMuxer coarse_muxer_velo
our coarse-level system muxer
std::size_t bytes() const
Returns the total amount of bytes allocated.
Global Filter wrapper class template.
std::size_t bytes() const
Returns the total amount of bytes allocated.
Global gate implementation.
void convert(const Gate< LVT2_, MT2_ > &other)
Conversion function for same vector container type but with different MDI-Type.
const LocalVector_ & get_freqs() const
Returns a const reference to the frequencies vector.
LocalVector_ _freqs
frequency vector
Global Matrix wrapper class template.
LocalMatrix_ convert_to_1() const
Computes and returns the type-1 conversion of this matrix as a local matrix.
LocalMatrix_ & local()
Returns a reference to the internal local LAFEM matrix object.
Mean Filter class template.
MeanFilter clone(LAFEM::CloneMode clone_mode=LAFEM::CloneMode::Deep) const
Creates a clone of itself.
Global multiplexer/demultiplexer implementation.
void convert(const Muxer< LVT2_, MT2_ > &other)
Conversion function for same vector container type but with different MDI-Type.
Global base-mesh vector splitter (and joiner) implementation.
void convert(const Splitter< LVT2_, MT2_ > &other)
Conversion function for same vector container type but with different MDI-Type.
Global grid-transfer operator class template.
void convert(MuxerType *coarse_muxer, const Transfer< LocalTransfer2_, Mirror2_ > &other)
container conversion function
Transfer clone(LAFEM::CloneMode clone_mode=LAFEM::CloneMode::Weak) const
Creates a clone of this object.
std::size_t bytes() const
Global vector wrapper class template.
Filter Chainclass template.
FilterChain clone(CloneMode clone_mode=CloneMode::Deep) const
Creates a clone of itself.
void filter_mat(Matrix_ &matrix) const
Applies the filter onto a system matrix.
Sequence of filters of the same type.
None Filter class template.
NoneFilter clone(CloneMode=CloneMode::Deep) const
Creates a (empty) clone of itself.
Saddle-Point matrix meta class template.
MatrixTypeA & block_a()
Returns the sub-matrix block A.
std::size_t bytes() const
Returns the total amount of bytes allocated.
MatrixTypeB & block_b()
Returns the sub-matrix block B.
MatrixTypeD & block_d()
Returns the sub-matrix block D.
Slip Filter class template.
Grid-Transfer operator class template.
Tuple-Diag-Matrix meta class template.
TupleVector meta-filter class template.
std::size_t bytes() const
Returns the total amount of bytes allocated.
TupleVector meta-mirror class template.
Variadic TupleVector class template.
TupleVector clone(LAFEM::CloneMode mode=LAFEM::CloneMode::Weak) const
Creates and returns a copy of this vector.
Unit Filter Blocked class template.
void filter_offdiag_row_mat(MatrixType &matrix) const
Filter the non-diagonal row entries.
UnitFilterBlocked clone(CloneMode clone_mode=CloneMode::Deep) const
Creates a clone of itself.
void filter_mat(MatrixType &matrix) const
Applies the filter onto a system matrix.
Unit Filter class template.
Handles vector prolongation, restriction and serialization.
void assemble_bilinear_operator_matrix_2(DomainAssembler< Trafo_ > &dom_asm, Matrix_ &matrix, const BilOp_ &bilinear_operator, const TestSpace_ &test_space, const TrialSpace_ &trial_space, const String &cubature, const typename Matrix_::DataType alpha=typename Matrix_::DataType(1))
Assembles a bilinear operator into a matrix with different test- and trial-spaces.
@ other
generic/other permutation strategy
double Real
Real data type.
std::uint64_t Index
Index data type.
System level using a MeanFilter for the pressure.
std::size_t bytes() const
Returns the total amount of bytes allocated.