*
* @author Matthias Maier, 2012 - 2015
*/
- template <typename FaceIterator>
+ template <typename FaceIterator, typename number>
void
make_periodicity_constraints(
const FaceIterator & face_1,
const typename identity<FaceIterator>::type &face_2,
- AffineConstraints<double> & constraints,
+ AffineConstraints<number> & constraints,
const ComponentMask & component_mask = ComponentMask(),
const bool face_orientation = true,
const bool face_flip = false,
*
* @author Daniel Arndt, Matthias Maier, 2013 - 2015
*/
- template <typename DoFHandlerType>
+ template <typename DoFHandlerType, typename number>
void
make_periodicity_constraints(
const std::vector<
GridTools::PeriodicFacePair<typename DoFHandlerType::cell_iterator>>
& periodic_faces,
- AffineConstraints<double> & constraints,
+ AffineConstraints<number> & constraints,
const ComponentMask & component_mask = ComponentMask(),
const std::vector<unsigned int> &first_vector_components =
std::vector<unsigned int>());
*
* @author Matthias Maier, 2012
*/
- template <typename DoFHandlerType>
+ template <typename DoFHandlerType, typename number>
void
make_periodicity_constraints(
const DoFHandlerType & dof_handler,
const types::boundary_id b_id1,
const types::boundary_id b_id2,
const int direction,
- AffineConstraints<double> &constraints,
+ AffineConstraints<number> &constraints,
const ComponentMask & component_mask = ComponentMask());
* @ref GlossPeriodicConstraints "Glossary entry on periodic boundary conditions"
* for further information.
*/
- template <typename DoFHandlerType>
+ template <typename DoFHandlerType, typename number>
void
make_periodicity_constraints(
const DoFHandlerType & dof_handler,
const types::boundary_id b_id,
const int direction,
- AffineConstraints<double> &constraints,
+ AffineConstraints<number> &constraints,
const ComponentMask & component_mask = ComponentMask());
/**
* such "identity constraints" if the opposite constraint already
* exists.
*/
- template <typename FaceIterator>
+ template <typename FaceIterator, typename number>
void
set_periodicity_constraints(
const FaceIterator & face_1,
const typename identity<FaceIterator>::type &face_2,
const FullMatrix<double> & transformation,
- AffineConstraints<double> & affine_constraints,
+ AffineConstraints<number> & affine_constraints,
const ComponentMask & component_mask,
const bool face_orientation,
const bool face_flip,
// a dependency cycle:
bool constraints_are_cyclic = true;
- double cycle_factor = factor;
+ number cycle_factor = factor;
for (auto test_dof = dof_right; test_dof != dof_left;)
{
// Low level interface:
- template <typename FaceIterator>
+ template <typename FaceIterator, typename number>
void
make_periodicity_constraints(
const FaceIterator & face_1,
const typename identity<FaceIterator>::type &face_2,
- AffineConstraints<double> & affine_constraints,
+ AffineConstraints<number> & affine_constraints,
const ComponentMask & component_mask,
const bool face_orientation,
const bool face_flip,
- template <typename DoFHandlerType>
+ template <typename DoFHandlerType, typename number>
void
make_periodicity_constraints(
const std::vector<
GridTools::PeriodicFacePair<typename DoFHandlerType::cell_iterator>>
& periodic_faces,
- AffineConstraints<double> & constraints,
+ AffineConstraints<number> & constraints,
const ComponentMask & component_mask,
const std::vector<unsigned int> &first_vector_components)
{
// High level interface variants:
- template <typename DoFHandlerType>
+ template <typename DoFHandlerType, typename number>
void
make_periodicity_constraints(const DoFHandlerType & dof_handler,
const types::boundary_id b_id1,
const types::boundary_id b_id2,
const int direction,
- dealii::AffineConstraints<double> &constraints,
+ dealii::AffineConstraints<number> &constraints,
const ComponentMask &component_mask)
{
static const int space_dim = DoFHandlerType::space_dimension;
- template <typename DoFHandlerType>
+ template <typename DoFHandlerType, typename number>
void
make_periodicity_constraints(const DoFHandlerType & dof_handler,
const types::boundary_id b_id,
const int direction,
- AffineConstraints<double> &constraints,
+ AffineConstraints<number> &constraints,
const ComponentMask & component_mask)
{
static const int dim = DoFHandlerType::dimension;
const FullMatrix<double> &,
const std::vector<unsigned int> &);
- template void DoFTools::make_periodicity_constraints<DH<deal_II_dimension>>(
+#ifdef DEAL_II_WITH_COMPLEX_VALUES
+ template void DoFTools::make_periodicity_constraints(
+ const DH<deal_II_dimension>::face_iterator &,
+ const DH<deal_II_dimension>::face_iterator &,
+ AffineConstraints<std::complex<double>> &,
+ const ComponentMask &,
+ bool,
+ bool,
+ bool,
+ const FullMatrix<double> &,
+ const std::vector<unsigned int> &);
+#endif
+
+ template void
+ DoFTools::make_periodicity_constraints<DH<deal_II_dimension>, double>(
const std::vector<
GridTools::PeriodicFacePair<DH<deal_II_dimension>::cell_iterator>> &,
AffineConstraints<double> &,
const ComponentMask &,
const std::vector<unsigned int> &);
+#ifdef DEAL_II_WITH_COMPLEX_VALUES
+ template void DoFTools::make_periodicity_constraints<DH<deal_II_dimension>,
+ std::complex<double>>(
+ const std::vector<
+ GridTools::PeriodicFacePair<DH<deal_II_dimension>::cell_iterator>> &,
+ AffineConstraints<std::complex<double>> &,
+ const ComponentMask &,
+ const std::vector<unsigned int> &);
+#endif
template void DoFTools::make_periodicity_constraints(
const DH<deal_II_dimension> &,
AffineConstraints<double> &,
const ComponentMask &);
+#ifdef DEAL_II_WITH_COMPLEX_VALUES
+ template void DoFTools::make_periodicity_constraints(
+ const DH<deal_II_dimension> &,
+ const types::boundary_id,
+ const types::boundary_id,
+ const int,
+ AffineConstraints<std::complex<double>> &,
+ const ComponentMask &);
+#endif
+
template void DoFTools::make_periodicity_constraints(
const DH<deal_II_dimension> &,
const types::boundary_id,
const int,
AffineConstraints<double> &,
const ComponentMask &);
+
+#ifdef DEAL_II_WITH_COMPLEX_VALUES
+ template void DoFTools::make_periodicity_constraints(
+ const DH<deal_II_dimension> &,
+ const types::boundary_id,
+ const int,
+ AffineConstraints<std::complex<double>> &,
+ const ComponentMask &);
+#endif
}
for (DH : DOFHANDLER_TEMPLATES; deal_II_dimension : DIMENSIONS;