@@ -20,15 +20,20 @@ namespace impl
20
20
// / not collapsed)
21
21
// / @param[in] parent_space The parent space (The same space as V if not
22
22
// / collapsed)
23
+ // / @param[in] tol Tolerance for adding scaled basis values to MPC. Any
24
+ // / contribution that is less than this value is ignored. The tolerance is also
25
+ // / added as padding for the bounding box trees and corresponding collision
26
+ // / searches to determine periodic degrees of freedom.
23
27
// / @returns The multi point constraint
24
28
template <typename T, std::floating_point U>
25
29
dolfinx_mpc::mpc_data<T> _create_periodic_condition (
26
30
const dolfinx::fem::FunctionSpace<U>& V,
27
31
std::span<std::int32_t > slave_blocks,
28
32
const std::function<std::vector<U>(std::span<const U>)>& relation, T scale,
29
33
const std::function<const std::int32_t(const std::int32_t &)>& parent_map,
30
- const dolfinx::fem::FunctionSpace<U>& parent_space)
34
+ const dolfinx::fem::FunctionSpace<U>& parent_space, const U tol )
31
35
{
36
+
32
37
// Map a list of indices in collapsed space back to the parent space
33
38
auto sub_to_parent = [&parent_map](const std::vector<std::int32_t >& sub_dofs)
34
39
{
@@ -71,10 +76,6 @@ dolfinx_mpc::mpc_data<T> _create_periodic_condition(
71
76
" Periodic conditions for vector valued spaces are not "
72
77
" implemented" );
73
78
74
- // Tolerance for adding scaled basis values to MPC. Any scaled basis
75
- // value with lower absolute value than the tolerance is ignored
76
- const U tol = 500 * std::numeric_limits<U>::epsilon ();
77
-
78
79
auto mesh = V.mesh ();
79
80
auto dofmap = V.dofmap ();
80
81
auto imap = dofmap->index_map ;
@@ -138,7 +139,7 @@ dolfinx_mpc::mpc_data<T> _create_periodic_condition(
138
139
= dolfinx_mpc::find_local_collisions<U>(*mesh, tree, mapped_T_b, tol);
139
140
dolfinx::common::Timer t0 (" ~~Periodic: Local cell and eval basis" );
140
141
auto [basis_values, basis_shape] = dolfinx_mpc::evaluate_basis_functions<U>(
141
- V, mapped_T_b, local_cell_collisions);
142
+ V, mapped_T_b, local_cell_collisions, tol );
142
143
MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
143
144
const U, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t , 3 >>
144
145
tabulated_basis_values (basis_values.data (), basis_shape);
@@ -358,7 +359,7 @@ dolfinx_mpc::mpc_data<T> _create_periodic_condition(
358
359
= dolfinx_mpc::find_local_collisions<U>(*mesh, tree, coords_recvb, tol);
359
360
auto [remote_basis_valuesb, r_basis_shape]
360
361
= dolfinx_mpc::evaluate_basis_functions<U>(V, coords_recvb,
361
- remote_cell_collisions);
362
+ remote_cell_collisions, tol );
362
363
MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
363
364
const U, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t , 3 >>
364
365
remote_basis_values (remote_basis_valuesb.data (), r_basis_shape);
@@ -490,6 +491,10 @@ dolfinx_mpc::mpc_data<T> _create_periodic_condition(
490
491
// / @param[in] scale Scaling of the periodic condition
491
492
// / @param[in] collapse If true, the list of marked dofs is in the collapsed
492
493
// / input space
494
+ // / @param[in] tol Tolerance for adding scaled basis values to MPC. Any
495
+ // / contribution that is less than this value is ignored. The tolerance is also
496
+ // / added as padding for the bounding box trees and corresponding collision
497
+ // / searches to determine periodic degrees of freedom.
493
498
// / @returns The multi point constraint
494
499
template <typename T, std::floating_point U>
495
500
dolfinx_mpc::mpc_data<T> geometrical_condition (
@@ -502,7 +507,7 @@ dolfinx_mpc::mpc_data<T> geometrical_condition(
502
507
indicator,
503
508
const std::function<std::vector<U>(std::span<const U>)>& relation,
504
509
const std::vector<std::shared_ptr<const dolfinx::fem::DirichletBC<T>>>& bcs,
505
- T scale, bool collapse)
510
+ T scale, bool collapse, const U tol )
506
511
{
507
512
std::vector<std::int32_t > reduced_blocks;
508
513
if (collapse)
@@ -526,7 +531,7 @@ dolfinx_mpc::mpc_data<T> geometrical_condition(
526
531
auto sub_map
527
532
= [&parent_map](const std::int32_t & i) { return parent_map[i]; };
528
533
return _create_periodic_condition<T>(V_sub, std::span (reduced_blocks),
529
- relation , scale, sub_map, *V);
534
+ relation , scale, sub_map, *V, tol );
530
535
}
531
536
else
532
537
{
@@ -542,7 +547,7 @@ dolfinx_mpc::mpc_data<T> geometrical_condition(
542
547
reduced_blocks.push_back (slave_blocks[i]);
543
548
auto sub_map = [](const std::int32_t & dof) { return dof; };
544
549
return _create_periodic_condition<T>(*V, std::span (reduced_blocks),
545
- relation , scale, sub_map, *V);
550
+ relation , scale, sub_map, *V, tol );
546
551
}
547
552
}
548
553
@@ -558,6 +563,10 @@ dolfinx_mpc::mpc_data<T> geometrical_condition(
558
563
// / @param[in] scale Scaling of the periodic condition
559
564
// / @param[in] collapse If true, the list of marked dofs is in the collapsed
560
565
// / input space
566
+ // / @param[in] tol Tolerance for adding scaled basis values to MPC. Any
567
+ // / contribution that is less than this value is ignored. The tolerance is also
568
+ // / added as padding for the bounding box trees and corresponding collision
569
+ // / searches to determine periodic degrees of freedom.
561
570
// / @returns The multi point constraint
562
571
template <typename T, std::floating_point U>
563
572
dolfinx_mpc::mpc_data<T> topological_condition (
@@ -566,7 +575,7 @@ dolfinx_mpc::mpc_data<T> topological_condition(
566
575
const std::int32_t tag,
567
576
const std::function<std::vector<U>(std::span<const U>)>& relation,
568
577
const std::vector<std::shared_ptr<const dolfinx::fem::DirichletBC<T>>>& bcs,
569
- T scale, bool collapse)
578
+ T scale, bool collapse, const U tol )
570
579
{
571
580
std::vector<std::int32_t > entities = meshtag->find (tag);
572
581
V->mesh ()->topology_mutable ()->create_connectivity (
@@ -594,7 +603,7 @@ dolfinx_mpc::mpc_data<T> topological_condition(
594
603
= [&parent_map](const std::int32_t & i) { return parent_map[i]; };
595
604
// Create mpc on sub space
596
605
dolfinx_mpc::mpc_data<T> sub_data = _create_periodic_condition<T>(
597
- V_sub, std::span (reduced_blocks), relation , scale, sub_map, *V);
606
+ V_sub, std::span (reduced_blocks), relation , scale, sub_map, *V, tol );
598
607
return sub_data;
599
608
}
600
609
else
@@ -613,7 +622,7 @@ dolfinx_mpc::mpc_data<T> topological_condition(
613
622
const auto sub_map = [](const std::int32_t & dof) { return dof; };
614
623
615
624
return _create_periodic_condition<T, U>(*V, std::span (reduced_blocks),
616
- relation , scale, sub_map, *V);
625
+ relation , scale, sub_map, *V, tol );
617
626
}
618
627
};
619
628
@@ -633,10 +642,11 @@ mpc_data<double> create_periodic_condition_geometrical(
633
642
const std::function<std::vector<double>(std::span<const double >)>& relation,
634
643
const std::vector<std::shared_ptr<const dolfinx::fem::DirichletBC<double>>>&
635
644
bcs,
636
- double scale, bool collapse)
645
+ double scale, bool collapse,
646
+ const double tol = 500 * std::numeric_limits<double>::epsilon())
637
647
{
638
648
return impl::geometrical_condition<double , double >(V, indicator, relation ,
639
- bcs, scale, collapse);
649
+ bcs, scale, collapse, tol );
640
650
}
641
651
642
652
mpc_data<std::complex<double >> create_periodic_condition_geometrical (
@@ -651,10 +661,11 @@ mpc_data<std::complex<double>> create_periodic_condition_geometrical(
651
661
const std::vector<
652
662
std::shared_ptr<const dolfinx::fem::DirichletBC<std::complex<double>>>>&
653
663
bcs,
654
- std::complex<double> scale, bool collapse)
664
+ std::complex<double> scale, bool collapse,
665
+ const double tol = 500 * std::numeric_limits<double>::epsilon())
655
666
{
656
667
return impl::geometrical_condition<std::complex<double >, double >(
657
- V, indicator, relation , bcs, scale, collapse);
668
+ V, indicator, relation , bcs, scale, collapse, tol );
658
669
}
659
670
660
671
mpc_data<double > create_periodic_condition_topological (
@@ -664,10 +675,11 @@ mpc_data<double> create_periodic_condition_topological(
664
675
const std::function<std::vector<double >(std::span<const double >)>& relation,
665
676
const std::vector<std::shared_ptr<const dolfinx::fem::DirichletBC<double>>>&
666
677
bcs,
667
- double scale, bool collapse)
678
+ double scale, bool collapse,
679
+ const double tol = 500 * std::numeric_limits<double>::epsilon())
668
680
{
669
681
return impl::topological_condition<double , double >(V, meshtag, tag, relation ,
670
- bcs, scale, collapse);
682
+ bcs, scale, collapse, tol );
671
683
}
672
684
673
685
mpc_data<std::complex<double >> create_periodic_condition_topological (
@@ -678,10 +690,11 @@ mpc_data<std::complex<double>> create_periodic_condition_topological(
678
690
const std::vector<
679
691
std::shared_ptr<const dolfinx::fem::DirichletBC<std::complex<double>>>>&
680
692
bcs,
681
- std::complex<double> scale, bool collapse)
693
+ std::complex<double> scale, bool collapse,
694
+ const double tol = 500 * std::numeric_limits<double>::epsilon())
682
695
{
683
696
return impl::topological_condition<std::complex<double >, double >(
684
- V, meshtag, tag, relation , bcs, scale, collapse);
697
+ V, meshtag, tag, relation , bcs, scale, collapse, tol );
685
698
}
686
699
687
700
mpc_data<float > create_periodic_condition_geometrical (
@@ -695,10 +708,11 @@ mpc_data<float> create_periodic_condition_geometrical(
695
708
const std::function<std::vector<float>(std::span<const float >)>& relation,
696
709
const std::vector<std::shared_ptr<const dolfinx::fem::DirichletBC<float>>>&
697
710
bcs,
698
- float scale, bool collapse)
711
+ float scale, bool collapse,
712
+ const float tol = 500 * std::numeric_limits<float>::epsilon())
699
713
{
700
714
return impl::geometrical_condition<float , float >(V, indicator, relation , bcs,
701
- scale, collapse);
715
+ scale, collapse, tol );
702
716
}
703
717
704
718
mpc_data<std::complex<float >> create_periodic_condition_geometrical (
@@ -713,10 +727,11 @@ mpc_data<std::complex<float>> create_periodic_condition_geometrical(
713
727
const std::vector<
714
728
std::shared_ptr<const dolfinx::fem::DirichletBC<std::complex<float>>>>&
715
729
bcs,
716
- std::complex<float> scale, bool collapse)
730
+ std::complex<float> scale, bool collapse,
731
+ const float tol = 500 * std::numeric_limits<float>::epsilon())
717
732
{
718
733
return impl::geometrical_condition<std::complex<float >, float >(
719
- V, indicator, relation , bcs, scale, collapse);
734
+ V, indicator, relation , bcs, scale, collapse, tol );
720
735
}
721
736
722
737
mpc_data<float > create_periodic_condition_topological (
@@ -726,10 +741,11 @@ mpc_data<float> create_periodic_condition_topological(
726
741
const std::function<std::vector<float >(std::span<const float >)>& relation,
727
742
const std::vector<std::shared_ptr<const dolfinx::fem::DirichletBC<float>>>&
728
743
bcs,
729
- float scale, bool collapse)
744
+ float scale, bool collapse,
745
+ const float tol = 500 * std::numeric_limits<float>::epsilon())
730
746
{
731
747
return impl::topological_condition<float , float >(V, meshtag, tag, relation ,
732
- bcs, scale, collapse);
748
+ bcs, scale, collapse, tol );
733
749
}
734
750
735
751
mpc_data<std::complex<float >> create_periodic_condition_topological (
@@ -740,10 +756,11 @@ mpc_data<std::complex<float>> create_periodic_condition_topological(
740
756
const std::vector<
741
757
std::shared_ptr<const dolfinx::fem::DirichletBC<std::complex<float>>>>&
742
758
bcs,
743
- std::complex<float> scale, bool collapse)
759
+ std::complex<float> scale, bool collapse,
760
+ const float tol = 500 * std::numeric_limits<float>::epsilon())
744
761
{
745
762
return impl::topological_condition<std::complex<float >, float >(
746
- V, meshtag, tag, relation , bcs, scale, collapse);
763
+ V, meshtag, tag, relation , bcs, scale, collapse, tol );
747
764
}
748
765
749
766
} // namespace dolfinx_mpc
0 commit comments