Functions | |
template<int dim0, int dim1, int spacedim, typename Matrix > | |
void | assemble_coupling_mass_matrix_with_exact_intersections (const DoFHandler< dim0, spacedim > &, const DoFHandler< dim1, spacedim > &, const std::vector< std::tuple< typename Triangulation< dim0, spacedim >::cell_iterator, typename Triangulation< dim1, spacedim >::cell_iterator, Quadrature< spacedim >>> &, Matrix &matrix, const AffineConstraints< typename Matrix::value_type > &, const ComponentMask &, const ComponentMask &, const Mapping< dim0, spacedim > &, const Mapping< dim1, spacedim > &, const AffineConstraints< typename Matrix::value_type > &) |
Create the coupling mass matrix for non-matching, overlapping grids in an "exact" way, i.e. More... | |
template<int dim0, int dim1, int spacedim, typename Matrix > | |
void | assemble_nitsche_with_exact_intersections (const DoFHandler< dim0, spacedim > &, const std::vector< std::tuple< typename Triangulation< dim0, spacedim >::cell_iterator, typename Triangulation< dim1, spacedim >::cell_iterator, Quadrature< spacedim >>> &, Matrix &matrix, const AffineConstraints< typename Matrix::value_type > &, const ComponentMask &, const Mapping< dim0, spacedim > &, const Function< spacedim, typename Matrix::value_type > &nitsche_coefficient=Functions::ConstantFunction< spacedim >(1.0), const double penalty=1.) |
Given two non-matching, overlapping grids, this function computes the local contributions \(\)M_{ij}:= \int_B \gamma v_i v_j dx \(\) exactly, by integrating on the intersection of the two grids. More... | |
template<int dim0, int dim1, int spacedim> | |
Quadrature< spacedim > | compute_intersection (const typename Triangulation< dim0, spacedim >::cell_iterator &cell0, const typename Triangulation< dim1, spacedim >::cell_iterator &cell1, const unsigned int degree, const Mapping< dim0, spacedim > &mapping0=(ReferenceCells::get_hypercube< dim0 >() .template get_default_linear_mapping< dim0, spacedim >()), const Mapping< dim1, spacedim > &mapping1=(ReferenceCells::get_hypercube< dim1 >() .template get_default_linear_mapping< dim1, spacedim >())) |
Intersect cell0 and cell1 and construct a Quadrature<spacedim> of degree degree over the intersection, i.e. More... | |
template<int dim0, int dim1, int spacedim> | |
std::vector< std::tuple< typename Triangulation< dim0, spacedim >::cell_iterator, typename Triangulation< dim1, spacedim >::cell_iterator, Quadrature< spacedim > > > | compute_intersection (const GridTools::Cache< dim0, spacedim > &space_cache, const GridTools::Cache< dim1, spacedim > &immersed_cache, const unsigned int degree, const double tol=0.) |
Given two triangulations cached inside GridTools::Cache objects, compute all intersections between the two and return a vector where each entry is a tuple containing iterators to the respective cells and a Quadrature<spacedim> formula to integrate over the intersection. More... | |
template<int dim0, int dim1, int spacedim, typename Sparsity , typename number = double> | |
void | create_coupling_sparsity_pattern_with_exact_intersections (const std::vector< std::tuple< typename Triangulation< dim0, spacedim >::cell_iterator, typename Triangulation< dim1, spacedim >::cell_iterator, Quadrature< spacedim >>> &intersections_info, const DoFHandler< dim0, spacedim > &space_dh, const DoFHandler< dim1, spacedim > &immersed_dh, Sparsity &sparsity, const AffineConstraints< number > &constraints=AffineConstraints< number >(), const ComponentMask &space_comps=ComponentMask(), const ComponentMask &immersed_comps=ComponentMask(), const AffineConstraints< number > &immersed_constraints=AffineConstraints< number >()) |
Create a coupling sparsity pattern of two non-matching, overlapped grids. More... | |
template<int dim0, int dim1, int spacedim, typename VectorType > | |
void | create_nitsche_rhs_with_exact_intersections (const DoFHandler< dim0, spacedim > &, const std::vector< std::tuple< typename Triangulation< dim0, spacedim >::cell_iterator, typename Triangulation< dim1, spacedim >::cell_iterator, Quadrature< spacedim >>> &, VectorType &vector, const AffineConstraints< typename VectorType::value_type > &, const Mapping< dim0, spacedim > &, const Function< spacedim, typename VectorType::value_type > &, const Function< spacedim, typename VectorType::value_type > &, const double penalty=1.) |
Create the r.h.s. More... | |
template void | assemble_coupling_mass_matrix_with_exact_intersections (const DoFHandler< 3, 3 > &, const DoFHandler< 3, 3 > &, const std::vector< std::tuple< typename Triangulation< 3, 3 >::cell_iterator, typename Triangulation< 3, 3 >::cell_iterator, Quadrature< 3 >>> &, SparseMatrix< double > &, const AffineConstraints< typename SparseMatrix< double >::value_type > &, const ComponentMask &, const ComponentMask &, const Mapping< 3, 3 > &space_mapping, const Mapping< 3, 3 > &immersed_mapping, const AffineConstraints< typename SparseMatrix< double >::value_type > &) |
template void | NonMatching::assemble_coupling_mass_matrix_with_exact_intersections< 2, 1, 2, TrilinosWrappers::SparseMatrix > (DoFHandler< 2, 2 > const &, DoFHandler< 1, 2 > const &, std::vector< std::tuple< Triangulation< 2, 2 >::cell_iterator, Triangulation< 1, 2 >::cell_iterator, Quadrature< 2 >>, std::allocator< std::tuple< Triangulation< 2, 2 >::cell_iterator, Triangulation< 1, 2 >::cell_iterator, Quadrature< 2 >>>> const &, TrilinosWrappers::SparseMatrix &, AffineConstraints< TrilinosWrappers::SparseMatrix::value_type > const &, ComponentMask const &, ComponentMask const &, Mapping< 2, 2 > const &, Mapping< 1, 2 > const &, AffineConstraints< TrilinosWrappers::SparseMatrix::value_type > const &) |
template void | NonMatching::assemble_coupling_mass_matrix_with_exact_intersections< 2, 2, 2, TrilinosWrappers::SparseMatrix > (DoFHandler< 2, 2 > const &, DoFHandler< 2, 2 > const &, std::vector< std::tuple< Triangulation< 2, 2 >::cell_iterator, Triangulation< 2, 2 >::cell_iterator, Quadrature< 2 >>, std::allocator< std::tuple< Triangulation< 2, 2 >::cell_iterator, Triangulation< 2, 2 >::cell_iterator, Quadrature< 2 >>>> const &, TrilinosWrappers::SparseMatrix &, AffineConstraints< TrilinosWrappers::SparseMatrix::value_type > const &, ComponentMask const &, ComponentMask const &, Mapping< 2, 2 > const &, Mapping< 2, 2 > const &, AffineConstraints< TrilinosWrappers::SparseMatrix::value_type > const &) |
template void | NonMatching::assemble_coupling_mass_matrix_with_exact_intersections< 3, 2, 3, TrilinosWrappers::SparseMatrix > (DoFHandler< 3, 3 > const &, DoFHandler< 2, 3 > const &, std::vector< std::tuple< Triangulation< 3, 3 >::cell_iterator, Triangulation< 2, 3 >::cell_iterator, Quadrature< 3 >>, std::allocator< std::tuple< Triangulation< 3, 3 >::cell_iterator, Triangulation< 2, 3 >::cell_iterator, Quadrature< 3 >>>> const &, TrilinosWrappers::SparseMatrix &, AffineConstraints< TrilinosWrappers::SparseMatrix::value_type > const &, ComponentMask const &, ComponentMask const &, Mapping< 3, 3 > const &, Mapping< 2, 3 > const &, AffineConstraints< TrilinosWrappers::SparseMatrix::value_type > const &) |
template void | NonMatching::assemble_coupling_mass_matrix_with_exact_intersections< 3, 3, 3, TrilinosWrappers::SparseMatrix > (DoFHandler< 3, 3 > const &, DoFHandler< 3, 3 > const &, std::vector< std::tuple< Triangulation< 3, 3 >::cell_iterator, Triangulation< 3, 3 >::cell_iterator, Quadrature< 3 >>, std::allocator< std::tuple< Triangulation< 3, 3 >::cell_iterator, Triangulation< 3, 3 >::cell_iterator, Quadrature< 3 >>>> const &, TrilinosWrappers::SparseMatrix &, AffineConstraints< TrilinosWrappers::SparseMatrix::value_type > const &, ComponentMask const &, ComponentMask const &, Mapping< 3, 3 > const &, Mapping< 3, 3 > const &, AffineConstraints< TrilinosWrappers::SparseMatrix::value_type > const &) |
template void | NonMatching::assemble_coupling_mass_matrix_with_exact_intersections< 2, 1, 2, PETScWrappers::MPI::SparseMatrix > (DoFHandler< 2, 2 > const &, DoFHandler< 1, 2 > const &, std::vector< std::tuple< Triangulation< 2, 2 >::cell_iterator, Triangulation< 1, 2 >::cell_iterator, Quadrature< 2 >>, std::allocator< std::tuple< Triangulation< 2, 2 >::cell_iterator, Triangulation< 1, 2 >::cell_iterator, Quadrature< 2 >>>> const &, PETScWrappers::MPI::SparseMatrix &, AffineConstraints< PETScWrappers::MPI::SparseMatrix::value_type > const &, ComponentMask const &, ComponentMask const &, Mapping< 2, 2 > const &, Mapping< 1, 2 > const &, AffineConstraints< PETScWrappers::MPI::SparseMatrix::value_type > const &) |
template void | NonMatching::assemble_coupling_mass_matrix_with_exact_intersections< 2, 2, 2, PETScWrappers::MPI::SparseMatrix > (DoFHandler< 2, 2 > const &, DoFHandler< 2, 2 > const &, std::vector< std::tuple< Triangulation< 2, 2 >::cell_iterator, Triangulation< 2, 2 >::cell_iterator, Quadrature< 2 >>, std::allocator< std::tuple< Triangulation< 2, 2 >::cell_iterator, Triangulation< 2, 2 >::cell_iterator, Quadrature< 2 >>>> const &, PETScWrappers::MPI::SparseMatrix &, AffineConstraints< PETScWrappers::MPI::SparseMatrix::value_type > const &, ComponentMask const &, ComponentMask const &, Mapping< 2, 2 > const &, Mapping< 2, 2 > const &, AffineConstraints< PETScWrappers::MPI::SparseMatrix::value_type > const &) |
template void | NonMatching::assemble_coupling_mass_matrix_with_exact_intersections< 3, 2, 3, PETScWrappers::MPI::SparseMatrix > (DoFHandler< 3, 3 > const &, DoFHandler< 2, 3 > const &, std::vector< std::tuple< Triangulation< 3, 3 >::cell_iterator, Triangulation< 2, 3 >::cell_iterator, Quadrature< 3 >>, std::allocator< std::tuple< Triangulation< 3, 3 >::cell_iterator, Triangulation< 2, 3 >::cell_iterator, Quadrature< 3 >>>> const &, PETScWrappers::MPI::SparseMatrix &, AffineConstraints< PETScWrappers::MPI::SparseMatrix::value_type > const &, ComponentMask const &, ComponentMask const &, Mapping< 3, 3 > const &, Mapping< 2, 3 > const &, AffineConstraints< PETScWrappers::MPI::SparseMatrix::value_type > const &) |
template void | NonMatching::assemble_coupling_mass_matrix_with_exact_intersections< 3, 3, 3, PETScWrappers::MPI::SparseMatrix > (DoFHandler< 3, 3 > const &, DoFHandler< 3, 3 > const &, std::vector< std::tuple< Triangulation< 3, 3 >::cell_iterator, Triangulation< 3, 3 >::cell_iterator, Quadrature< 3 >>, std::allocator< std::tuple< Triangulation< 3, 3 >::cell_iterator, Triangulation< 3, 3 >::cell_iterator, Quadrature< 3 >>>> const &, PETScWrappers::MPI::SparseMatrix &, AffineConstraints< PETScWrappers::MPI::SparseMatrix::value_type > const &, ComponentMask const &, ComponentMask const &, Mapping< 3, 3 > const &, Mapping< 3, 3 > const &, AffineConstraints< PETScWrappers::MPI::SparseMatrix::value_type > const &) |
template Quadrature< 1 > | compute_cell_intersection (const Triangulation< 1, 1 >::cell_iterator &, const Triangulation< 1, 1 >::cell_iterator &, const unsigned int, const Mapping< 1, 1 > &, const Mapping< 1, 1 > &, const double) |
template Quadrature< 2 > | compute_cell_intersection (const Triangulation< 2, 2 >::cell_iterator &, const Triangulation< 1, 2 >::cell_iterator &, const unsigned int, const Mapping< 2, 2 > &, const Mapping< 1, 2 > &, const double) |
template Quadrature< 2 > | compute_cell_intersection (const Triangulation< 2, 2 >::cell_iterator &, const Triangulation< 2, 2 >::cell_iterator &, const unsigned int, const Mapping< 2, 2 > &, const Mapping< 2, 2 > &, const double) |
template Quadrature< 3 > | compute_cell_intersection (const Triangulation< 3, 3 >::cell_iterator &, const Triangulation< 2, 3 >::cell_iterator &, const unsigned int, const Mapping< 3, 3 > &, const Mapping< 2, 3 > &, const double) |
template Quadrature< 3 > | compute_cell_intersection (const Triangulation< 3, 3 >::cell_iterator &, const Triangulation< 3, 3 >::cell_iterator &, const unsigned int, const Mapping< 3, 3 > &, const Mapping< 3, 3 > &, const double) |
template void | create_coupling_sparsity_pattern_with_exact_intersections (const std::vector< std::tuple< typename Triangulation< 1, 1 >::cell_iterator, typename Triangulation< 1, 1 >::cell_iterator, Quadrature< 1 >>> &intersections_info, const DoFHandler< 1, 1 > &space_dh, const DoFHandler< 1, 1 > &immersed_dh, DynamicSparsityPattern &sparsity, const AffineConstraints< double > &constraints, const ComponentMask &space_comps, const ComponentMask &immersed_comps, const AffineConstraints< double > &immersed_constraint) |
template void | create_coupling_sparsity_pattern_with_exact_intersections (const std::vector< std::tuple< typename Triangulation< 2, 2 >::cell_iterator, typename Triangulation< 1, 2 >::cell_iterator, Quadrature< 2 >>> &intersections_info, const DoFHandler< 2, 2 > &space_dh, const DoFHandler< 1, 2 > &immersed_dh, DynamicSparsityPattern &sparsity, const AffineConstraints< double > &constraints, const ComponentMask &space_comps, const ComponentMask &immersed_comps, const AffineConstraints< double > &immersed_constraint) |
template void | create_coupling_sparsity_pattern_with_exact_intersections (const std::vector< std::tuple< typename Triangulation< 3, 3 >::cell_iterator, typename Triangulation< 1, 3 >::cell_iterator, Quadrature< 3 >>> &intersections_info, const DoFHandler< 3, 3 > &space_dh, const DoFHandler< 1, 3 > &immersed_dh, DynamicSparsityPattern &sparsity, const AffineConstraints< double > &constraints, const ComponentMask &space_comps, const ComponentMask &immersed_comps, const AffineConstraints< double > &immersed_constraint) |
template void | create_coupling_sparsity_pattern_with_exact_intersections (const std::vector< std::tuple< typename Triangulation< 2, 2 >::cell_iterator, typename Triangulation< 2, 2 >::cell_iterator, Quadrature< 2 >>> &intersections_info, const DoFHandler< 2, 2 > &space_dh, const DoFHandler< 2, 2 > &immersed_dh, DynamicSparsityPattern &sparsity, const AffineConstraints< double > &constraints, const ComponentMask &space_comps, const ComponentMask &immersed_comps, const AffineConstraints< double > &immersed_constraint) |
template void | create_coupling_sparsity_pattern_with_exact_intersections (const std::vector< std::tuple< typename Triangulation< 3, 3 >::cell_iterator, typename Triangulation< 2, 3 >::cell_iterator, Quadrature< 3 >>> &intersections_info, const DoFHandler< 3, 3 > &space_dh, const DoFHandler< 2, 3 > &immersed_dh, DynamicSparsityPattern &sparsity, const AffineConstraints< double > &constraints, const ComponentMask &space_comps, const ComponentMask &immersed_comps, const AffineConstraints< double > &immersed_constraint) |
template void | create_coupling_sparsity_pattern_with_exact_intersections (const std::vector< std::tuple< typename Triangulation< 3, 3 >::cell_iterator, typename Triangulation< 3, 3 >::cell_iterator, Quadrature< 3 >>> &intersections_info, const DoFHandler< 3, 3 > &space_dh, const DoFHandler< 3, 3 > &immersed_dh, DynamicSparsityPattern &sparsity, const AffineConstraints< double > &constraints, const ComponentMask &space_comps, const ComponentMask &immersed_comps, const AffineConstraints< double > &immersed_constraint) |
template void | create_coupling_sparsity_pattern_with_exact_intersections< 2, 1, 2, TrilinosWrappers::SparsityPattern, double > (std::vector< std::tuple< Triangulation< 2, 2 >::cell_iterator, Triangulation< 1, 2 >::cell_iterator, Quadrature< 2 >>, std::allocator< std::tuple< Triangulation< 2, 2 >::cell_iterator, Triangulation< 1, 2 >::cell_iterator, Quadrature< 2 >>>> const &, DoFHandler< 2, 2 > const &, DoFHandler< 1, 2 > const &, TrilinosWrappers::SparsityPattern &, AffineConstraints< double > const &, ComponentMask const &, ComponentMask const &, AffineConstraints< double > const &) |
template void | create_coupling_sparsity_pattern_with_exact_intersections< 2, 2, 2, TrilinosWrappers::SparsityPattern, double > (std::vector< std::tuple< Triangulation< 2, 2 >::cell_iterator, Triangulation< 2, 2 >::cell_iterator, Quadrature< 2 >>, std::allocator< std::tuple< Triangulation< 2, 2 >::cell_iterator, Triangulation< 2, 2 >::cell_iterator, Quadrature< 2 >>>> const &, DoFHandler< 2, 2 > const &, DoFHandler< 2, 2 > const &, TrilinosWrappers::SparsityPattern &, AffineConstraints< double > const &, ComponentMask const &, ComponentMask const &, AffineConstraints< double > const &) |
template void | create_coupling_sparsity_pattern_with_exact_intersections< 3, 2, 3, TrilinosWrappers::SparsityPattern, double > (std::vector< std::tuple< Triangulation< 3, 3 >::cell_iterator, Triangulation< 2, 3 >::cell_iterator, Quadrature< 3 >>, std::allocator< std::tuple< Triangulation< 3, 3 >::cell_iterator, Triangulation< 2, 3 >::cell_iterator, Quadrature< 3 >>>> const &, DoFHandler< 3, 3 > const &, DoFHandler< 2, 3 > const &, TrilinosWrappers::SparsityPattern &, AffineConstraints< double > const &, ComponentMask const &, ComponentMask const &, AffineConstraints< double > const &) |
template void | create_coupling_sparsity_pattern_with_exact_intersections< 3, 3, 3, TrilinosWrappers::SparsityPattern, double > (std::vector< std::tuple< Triangulation< 3, 3 >::cell_iterator, Triangulation< 3, 3 >::cell_iterator, Quadrature< 3 >>, std::allocator< std::tuple< Triangulation< 3, 3 >::cell_iterator, Triangulation< 3, 3 >::cell_iterator, Quadrature< 3 >>>> const &, DoFHandler< 3, 3 > const &, DoFHandler< 3, 3 > const &, TrilinosWrappers::SparsityPattern &, AffineConstraints< double > const &, ComponentMask const &, ComponentMask const &, AffineConstraints< double > const &) |
template void | create_coupling_sparsity_pattern_with_exact_intersections< 3, 2, 3, SparsityPattern, double > (std::vector< std::tuple< Triangulation< 3, 3 >::cell_iterator, Triangulation< 2, 3 >::cell_iterator, Quadrature< 3 >>, std::allocator< std::tuple< Triangulation< 3, 3 >::cell_iterator, Triangulation< 2, 3 >::cell_iterator, Quadrature< 3 >>>> const &, DoFHandler< 3, 3 > const &, DoFHandler< 2, 3 > const &, SparsityPattern &, AffineConstraints< double > const &, ComponentMask const &, ComponentMask const &, AffineConstraints< double > const &) |
template void | create_coupling_sparsity_pattern_with_exact_intersections< 2, 2, 2, SparsityPattern, double > (std::vector< std::tuple< Triangulation< 2, 2 >::cell_iterator, Triangulation< 2, 2 >::cell_iterator, Quadrature< 2 >>, std::allocator< std::tuple< Triangulation< 2, 2 >::cell_iterator, Triangulation< 2, 2 >::cell_iterator, Quadrature< 2 >>>> const &, DoFHandler< 2, 2 > const &, DoFHandler< 2, 2 > const &, SparsityPattern &, AffineConstraints< double > const &, ComponentMask const &, ComponentMask const &, AffineConstraints< double > const &) |
template void | create_coupling_sparsity_pattern_with_exact_intersections< 1, 1, 1, SparsityPattern, double > (std::vector< std::tuple< Triangulation< 1, 1 >::cell_iterator, Triangulation< 1, 1 >::cell_iterator, Quadrature< 1 >>, std::allocator< std::tuple< Triangulation< 1, 1 >::cell_iterator, Triangulation< 1, 1 >::cell_iterator, Quadrature< 1 >>>> const &, DoFHandler< 1, 1 > const &, DoFHandler< 1, 1 > const &, SparsityPattern &, AffineConstraints< double > const &, ComponentMask const &, ComponentMask const &, AffineConstraints< double > const &) |
template void | create_coupling_sparsity_pattern_with_exact_intersections< 3, 1, 3, SparsityPattern, double > (std::vector< std::tuple< Triangulation< 3, 3 >::cell_iterator, Triangulation< 1, 3 >::cell_iterator, Quadrature< 3 >>, std::allocator< std::tuple< Triangulation< 3, 3 >::cell_iterator, Triangulation< 1, 3 >::cell_iterator, Quadrature< 3 >>>> const &, DoFHandler< 3, 3 > const &, DoFHandler< 1, 3 > const &, SparsityPattern &, AffineConstraints< double > const &, ComponentMask const &, ComponentMask const &, AffineConstraints< double > const &) |
template void | create_coupling_sparsity_pattern_with_exact_intersections< 2, 1, 2, SparsityPattern, double > (std::vector< std::tuple< Triangulation< 2, 2 >::cell_iterator, Triangulation< 1, 2 >::cell_iterator, Quadrature< 2 >>, std::allocator< std::tuple< Triangulation< 2, 2 >::cell_iterator, Triangulation< 1, 2 >::cell_iterator, Quadrature< 2 >>>> const &, DoFHandler< 2, 2 > const &, DoFHandler< 1, 2 > const &, SparsityPattern &, AffineConstraints< double > const &, ComponentMask const &, ComponentMask const &, AffineConstraints< double > const &) |
template void | create_coupling_sparsity_pattern_with_exact_intersections< 3, 3, 3, SparsityPattern, double > (std::vector< std::tuple< Triangulation< 3, 3 >::cell_iterator, Triangulation< 3, 3 >::cell_iterator, Quadrature< 3 >>, std::allocator< std::tuple< Triangulation< 3, 3 >::cell_iterator, Triangulation< 3, 3 >::cell_iterator, Quadrature< 3 >>>> const &, DoFHandler< 3, 3 > const &, DoFHandler< 3, 3 > const &, SparsityPattern &, AffineConstraints< double > const &, ComponentMask const &, ComponentMask const &, AffineConstraints< double > const &) |
template void | create_nitsche_rhs_with_exact_intersections< 2, 2, 2 > (const DoFHandler< 2, 2 > &, const std::vector< std::tuple< typename Triangulation< 2, 2 >::cell_iterator, typename Triangulation< 2, 2 >::cell_iterator, Quadrature< 2 >>> &, Vector< double > &vector, const AffineConstraints< double > &, const Mapping< 2, 2 > &, const Function< 2, double > &, const Function< 2, double > &, const double) |
template void | create_nitsche_rhs_with_exact_intersections< 2, 1, 2 > (const DoFHandler< 2, 2 > &, const std::vector< std::tuple< typename Triangulation< 2, 2 >::cell_iterator, typename Triangulation< 1, 2 >::cell_iterator, Quadrature< 2 >>> &, Vector< double > &vector, const AffineConstraints< double > &, const Mapping< 2, 2 > &, const Function< 2, double > &, const Function< 2, double > &, const double) |
template void | create_nitsche_rhs_with_exact_intersections< 3, 3, 3 > (const DoFHandler< 3, 3 > &, const std::vector< std::tuple< typename Triangulation< 3, 3 >::cell_iterator, typename Triangulation< 3, 3 >::cell_iterator, Quadrature< 3 >>> &, Vector< double > &vector, const AffineConstraints< double > &, const Mapping< 3, 3 > &, const Function< 3, double > &, const Function< 3, double > &, const double) |
template void | create_nitsche_rhs_with_exact_intersections< 3, 2, 3 > (const DoFHandler< 3, 3 > &, const std::vector< std::tuple< typename Triangulation< 3, 3 >::cell_iterator, typename Triangulation< 2, 3 >::cell_iterator, Quadrature< 3 >>> &, Vector< double > &vector, const AffineConstraints< double > &, const Mapping< 3, 3 > &, const Function< 3, double > &, const Function< 3, double > &, const double) |
void dealii::NonMatching::assemble_coupling_mass_matrix_with_exact_intersections | ( | const DoFHandler< dim0, spacedim > & | , |
const DoFHandler< dim1, spacedim > & | , | ||
const std::vector< std::tuple< typename Triangulation< dim0, spacedim >::cell_iterator, typename Triangulation< dim1, spacedim >::cell_iterator, Quadrature< spacedim >>> & | , | ||
Matrix & | matrix, | ||
const AffineConstraints< typename Matrix::value_type > & | , | ||
const ComponentMask & | , | ||
const ComponentMask & | , | ||
const Mapping< dim0, spacedim > & | , | ||
const Mapping< dim1, spacedim > & | , | ||
const AffineConstraints< typename Matrix::value_type > & | |||
) |
Create the coupling mass matrix for non-matching, overlapping grids in an "exact" way, i.e.
by computing the local contributions \(\)M_{ij}:= \int_B v_i w_j dx \(\) as products of cellwise smooth functions on the intersection of the two grids. This information is described by a std::vector<std::tuple>>
where each tuple contains the two intersected cells and a Quadrature formula on their intersection.
dim0 | Intrinsic dimension of the first, space grid |
dim1 | Intrinsic dimension of the second, embedded space |
spacedim | Ambient space intrinsic dimension |
Matrix | Matrix type you wish to use |
Definition at line 251 of file assemble_coupling_mass_matrix_with_exact_intersections.cc.
References Assert, and StandardExceptions::ExcMessage().
void dealii::NonMatching::assemble_nitsche_with_exact_intersections | ( | const DoFHandler< dim0, spacedim > & | , |
const std::vector< std::tuple< typename Triangulation< dim0, spacedim >::cell_iterator, typename Triangulation< dim1, spacedim >::cell_iterator, Quadrature< spacedim >>> & | , | ||
Matrix & | matrix, | ||
const AffineConstraints< typename Matrix::value_type > & | , | ||
const ComponentMask & | , | ||
const Mapping< dim0, spacedim > & | , | ||
const Function< spacedim, typename Matrix::value_type > & | nitsche_coefficient = Functions::ConstantFunction< spacedim >(1.0) , |
||
const double | penalty = 1. |
||
) |
Given two non-matching, overlapping grids, this function computes the local contributions \(\)M_{ij}:= \int_B \gamma v_i v_j dx \(\) exactly, by integrating on the intersection of the two grids.
There's no need to change the sparsity pattern, as the DoFs are all living on the same cell.
dim0 | Intrinsic dimension of the first, space grid |
dim1 | Intrinsic dimension of the second, embedded space |
spacedim | Ambient space intrinsic dimension |
Matrix | Matrix type you wish to use |
Quadrature< spacedim > dealii::NonMatching::compute_intersection | ( | const typename Triangulation< dim0, spacedim >::cell_iterator & | cell0, |
const typename Triangulation< dim1, spacedim >::cell_iterator & | cell1, | ||
const unsigned int | degree, | ||
const Mapping< dim0, spacedim > & | mapping0 = (ReferenceCells::get_hypercube<dim0>() .template get_default_linear_mapping<dim0, spacedim>()) , |
||
const Mapping< dim1, spacedim > & | mapping1 = (ReferenceCells::get_hypercube<dim1>() .template get_default_linear_mapping<dim1, spacedim>()) |
||
) |
Intersect cell0
and cell1
and construct a Quadrature<spacedim>
of degree degree
over the intersection, i.e.
in the real space. Mappings for both cells are inmapping0 and
mapping1`, respectively.
dim0 | |
dim1 | |
spacedim |
cell0 | A cell_iteratator to the first cell |
cell1 | A cell_iteratator to the first cell |
degree | The degree of the Quadrature you want to build there |
mapping0 | The Mapping object describing the first cell |
mapping1 | The Mapping object describing the second cell |
Definition at line 119 of file compute_intersections.cc.
References Assert, and StandardExceptions::ExcMessage().
std::vector< std::tuple< typename Triangulation< dim0, spacedim >::cell_iterator, typename Triangulation< dim1, spacedim >::cell_iterator, Quadrature< spacedim > > > dealii::NonMatching::compute_intersection | ( | const GridTools::Cache< dim0, spacedim > & | space_cache, |
const GridTools::Cache< dim1, spacedim > & | immersed_cache, | ||
const unsigned int | degree, | ||
const double | tol = 0. |
||
) |
Given two triangulations cached inside GridTools::Cache
objects, compute all intersections between the two and return a vector where each entry is a tuple containing iterators to the respective cells and a Quadrature<spacedim>
formula to integrate over the intersection.
dim0 | Intrinsic dimension of the immersed grid |
dim1 | Intrinsic dimension of the ambient grid |
spacedim |
space_cache | |
immersed_cache | |
degree | Degree of the desired quadrature formula |
Definition at line 140 of file compute_intersections.cc.
References Assert, StandardExceptions::ExcMessage(), GridTools::Cache< int dim, int spacedim >::get_cell_bounding_boxes_rtree(), GridTools::Cache< int dim, int spacedim >::get_locally_owned_cell_bounding_boxes_rtree(), and GridTools::Cache< int dim, int spacedim >::get_mapping().
void dealii::NonMatching::create_coupling_sparsity_pattern_with_exact_intersections | ( | const std::vector< std::tuple< typename Triangulation< dim0, spacedim >::cell_iterator, typename Triangulation< dim1, spacedim >::cell_iterator, Quadrature< spacedim >>> & | intersections_info, |
const DoFHandler< dim0, spacedim > & | space_dh, | ||
const DoFHandler< dim1, spacedim > & | immersed_dh, | ||
Sparsity & | sparsity, | ||
const AffineConstraints< number > & | constraints = AffineConstraints<number>() , |
||
const ComponentMask & | space_comps = ComponentMask() , |
||
const ComponentMask & | immersed_comps = ComponentMask() , |
||
const AffineConstraints< number > & | immersed_constraints = AffineConstraints<number>() |
||
) |
Create a coupling sparsity pattern of two non-matching, overlapped grids.
As it relies on compute_intersection
, the "small" intersections do not enter in the sparsity pattern.
intersections_info | A vector of tuples where the i-th entry contains two active_cell_iterator s to the intersected cells |
space_dh | DoFHandler object for the space grid |
immersed_dh | DoFHandler object for the embedded grid |
sparsity | The sparsity pattern to be filled |
constraints | AffineConstraints for the space grid |
space_comps | Mask for the space space components of the finite element |
immersed_comps | Mask for the embedded components of the finite element |
immersed_constraints | AffineConstraints for the embedded grid |
Definition at line 194 of file create_coupling_sparsity_pattern_with_exact_intersections.cc.
References Assert, and StandardExceptions::ExcMessage().
void dealii::NonMatching::create_nitsche_rhs_with_exact_intersections | ( | const DoFHandler< dim0, spacedim > & | space_dh, |
const std::vector< std::tuple< typename Triangulation< dim0, spacedim >::cell_iterator, typename Triangulation< dim1, spacedim >::cell_iterator, Quadrature< spacedim >>> & | cells_and_quads, | ||
VectorType & | vector, | ||
const AffineConstraints< typename VectorType::value_type > & | space_constraints, | ||
const Mapping< dim0, spacedim > & | space_mapping, | ||
const Function< spacedim, typename VectorType::value_type > & | rhs_function, | ||
const Function< spacedim, typename VectorType::value_type > & | coefficient, | ||
const double | penalty = 1. |
||
) |
Create the r.h.s.
\(\int_{B} f v\), using exact quadratures on the intersections of the two grids.
dim0 | |
dim1 | |
spacedim | |
VectorType |
vector |
Definition at line 28 of file create_nitsche_rhs_with_exact_intersections.cc.
References Assert, AssertDimension, AffineConstraints< typename number >::distribute_local_to_global(), StandardExceptions::ExcMessage(), DoFHandler< int dim, int spacedim >::get_fe(), std::min(), DoFHandler< int dim, int spacedim >::n_dofs(), FiniteElement< int dim, int spacedim >::n_dofs_per_cell(), Mapping< int dim, int spacedim >::transform_points_real_to_unit_cell(), and Function< int dim, typename RangeNumberType >::value_list().
template void dealii::NonMatching::assemble_coupling_mass_matrix_with_exact_intersections | ( | const DoFHandler< 3, 3 > & | , |
const DoFHandler< 3, 3 > & | , | ||
const std::vector< std::tuple< typename Triangulation< 3, 3 >::cell_iterator, typename Triangulation< 3, 3 >::cell_iterator, Quadrature< 3 >>> & | , | ||
SparseMatrix< double > & | , | ||
const AffineConstraints< typename SparseMatrix< double >::value_type > & | , | ||
const ComponentMask & | , | ||
const ComponentMask & | , | ||
const Mapping< 3, 3 > & | space_mapping, | ||
const Mapping< 3, 3 > & | immersed_mapping, | ||
const AffineConstraints< typename SparseMatrix< double >::value_type > & | |||
) |
template void dealii::NonMatching::NonMatching::assemble_coupling_mass_matrix_with_exact_intersections< 2, 1, 2, TrilinosWrappers::SparseMatrix > | ( | DoFHandler< 2, 2 > const & | , |
DoFHandler< 1, 2 > const & | , | ||
std::vector< std::tuple< Triangulation< 2, 2 >::cell_iterator, Triangulation< 1, 2 >::cell_iterator, Quadrature< 2 >>, std::allocator< std::tuple< Triangulation< 2, 2 >::cell_iterator, Triangulation< 1, 2 >::cell_iterator, Quadrature< 2 >>>> const & | , | ||
TrilinosWrappers::SparseMatrix & | , | ||
AffineConstraints< TrilinosWrappers::SparseMatrix::value_type > const & | , | ||
ComponentMask const & | , | ||
ComponentMask const & | , | ||
Mapping< 2, 2 > const & | , | ||
Mapping< 1, 2 > const & | , | ||
AffineConstraints< TrilinosWrappers::SparseMatrix::value_type > const & | |||
) |
template void dealii::NonMatching::NonMatching::assemble_coupling_mass_matrix_with_exact_intersections< 2, 2, 2, TrilinosWrappers::SparseMatrix > | ( | DoFHandler< 2, 2 > const & | , |
DoFHandler< 2, 2 > const & | , | ||
std::vector< std::tuple< Triangulation< 2, 2 >::cell_iterator, Triangulation< 2, 2 >::cell_iterator, Quadrature< 2 >>, std::allocator< std::tuple< Triangulation< 2, 2 >::cell_iterator, Triangulation< 2, 2 >::cell_iterator, Quadrature< 2 >>>> const & | , | ||
TrilinosWrappers::SparseMatrix & | , | ||
AffineConstraints< TrilinosWrappers::SparseMatrix::value_type > const & | , | ||
ComponentMask const & | , | ||
ComponentMask const & | , | ||
Mapping< 2, 2 > const & | , | ||
Mapping< 2, 2 > const & | , | ||
AffineConstraints< TrilinosWrappers::SparseMatrix::value_type > const & | |||
) |
template void dealii::NonMatching::NonMatching::assemble_coupling_mass_matrix_with_exact_intersections< 3, 2, 3, TrilinosWrappers::SparseMatrix > | ( | DoFHandler< 3, 3 > const & | , |
DoFHandler< 2, 3 > const & | , | ||
std::vector< std::tuple< Triangulation< 3, 3 >::cell_iterator, Triangulation< 2, 3 >::cell_iterator, Quadrature< 3 >>, std::allocator< std::tuple< Triangulation< 3, 3 >::cell_iterator, Triangulation< 2, 3 >::cell_iterator, Quadrature< 3 >>>> const & | , | ||
TrilinosWrappers::SparseMatrix & | , | ||
AffineConstraints< TrilinosWrappers::SparseMatrix::value_type > const & | , | ||
ComponentMask const & | , | ||
ComponentMask const & | , | ||
Mapping< 3, 3 > const & | , | ||
Mapping< 2, 3 > const & | , | ||
AffineConstraints< TrilinosWrappers::SparseMatrix::value_type > const & | |||
) |
template void dealii::NonMatching::NonMatching::assemble_coupling_mass_matrix_with_exact_intersections< 3, 3, 3, TrilinosWrappers::SparseMatrix > | ( | DoFHandler< 3, 3 > const & | , |
DoFHandler< 3, 3 > const & | , | ||
std::vector< std::tuple< Triangulation< 3, 3 >::cell_iterator, Triangulation< 3, 3 >::cell_iterator, Quadrature< 3 >>, std::allocator< std::tuple< Triangulation< 3, 3 >::cell_iterator, Triangulation< 3, 3 >::cell_iterator, Quadrature< 3 >>>> const & | , | ||
TrilinosWrappers::SparseMatrix & | , | ||
AffineConstraints< TrilinosWrappers::SparseMatrix::value_type > const & | , | ||
ComponentMask const & | , | ||
ComponentMask const & | , | ||
Mapping< 3, 3 > const & | , | ||
Mapping< 3, 3 > const & | , | ||
AffineConstraints< TrilinosWrappers::SparseMatrix::value_type > const & | |||
) |
template void dealii::NonMatching::NonMatching::assemble_coupling_mass_matrix_with_exact_intersections< 2, 1, 2, PETScWrappers::MPI::SparseMatrix > | ( | DoFHandler< 2, 2 > const & | , |
DoFHandler< 1, 2 > const & | , | ||
std::vector< std::tuple< Triangulation< 2, 2 >::cell_iterator, Triangulation< 1, 2 >::cell_iterator, Quadrature< 2 >>, std::allocator< std::tuple< Triangulation< 2, 2 >::cell_iterator, Triangulation< 1, 2 >::cell_iterator, Quadrature< 2 >>>> const & | , | ||
PETScWrappers::MPI::SparseMatrix & | , | ||
AffineConstraints< PETScWrappers::MPI::SparseMatrix::value_type > const & | , | ||
ComponentMask const & | , | ||
ComponentMask const & | , | ||
Mapping< 2, 2 > const & | , | ||
Mapping< 1, 2 > const & | , | ||
AffineConstraints< PETScWrappers::MPI::SparseMatrix::value_type > const & | |||
) |
template void dealii::NonMatching::NonMatching::assemble_coupling_mass_matrix_with_exact_intersections< 2, 2, 2, PETScWrappers::MPI::SparseMatrix > | ( | DoFHandler< 2, 2 > const & | , |
DoFHandler< 2, 2 > const & | , | ||
std::vector< std::tuple< Triangulation< 2, 2 >::cell_iterator, Triangulation< 2, 2 >::cell_iterator, Quadrature< 2 >>, std::allocator< std::tuple< Triangulation< 2, 2 >::cell_iterator, Triangulation< 2, 2 >::cell_iterator, Quadrature< 2 >>>> const & | , | ||
PETScWrappers::MPI::SparseMatrix & | , | ||
AffineConstraints< PETScWrappers::MPI::SparseMatrix::value_type > const & | , | ||
ComponentMask const & | , | ||
ComponentMask const & | , | ||
Mapping< 2, 2 > const & | , | ||
Mapping< 2, 2 > const & | , | ||
AffineConstraints< PETScWrappers::MPI::SparseMatrix::value_type > const & | |||
) |
template void dealii::NonMatching::NonMatching::assemble_coupling_mass_matrix_with_exact_intersections< 3, 2, 3, PETScWrappers::MPI::SparseMatrix > | ( | DoFHandler< 3, 3 > const & | , |
DoFHandler< 2, 3 > const & | , | ||
std::vector< std::tuple< Triangulation< 3, 3 >::cell_iterator, Triangulation< 2, 3 >::cell_iterator, Quadrature< 3 >>, std::allocator< std::tuple< Triangulation< 3, 3 >::cell_iterator, Triangulation< 2, 3 >::cell_iterator, Quadrature< 3 >>>> const & | , | ||
PETScWrappers::MPI::SparseMatrix & | , | ||
AffineConstraints< PETScWrappers::MPI::SparseMatrix::value_type > const & | , | ||
ComponentMask const & | , | ||
ComponentMask const & | , | ||
Mapping< 3, 3 > const & | , | ||
Mapping< 2, 3 > const & | , | ||
AffineConstraints< PETScWrappers::MPI::SparseMatrix::value_type > const & | |||
) |
template void dealii::NonMatching::NonMatching::assemble_coupling_mass_matrix_with_exact_intersections< 3, 3, 3, PETScWrappers::MPI::SparseMatrix > | ( | DoFHandler< 3, 3 > const & | , |
DoFHandler< 3, 3 > const & | , | ||
std::vector< std::tuple< Triangulation< 3, 3 >::cell_iterator, Triangulation< 3, 3 >::cell_iterator, Quadrature< 3 >>, std::allocator< std::tuple< Triangulation< 3, 3 >::cell_iterator, Triangulation< 3, 3 >::cell_iterator, Quadrature< 3 >>>> const & | , | ||
PETScWrappers::MPI::SparseMatrix & | , | ||
AffineConstraints< PETScWrappers::MPI::SparseMatrix::value_type > const & | , | ||
ComponentMask const & | , | ||
ComponentMask const & | , | ||
Mapping< 3, 3 > const & | , | ||
Mapping< 3, 3 > const & | , | ||
AffineConstraints< PETScWrappers::MPI::SparseMatrix::value_type > const & | |||
) |
template Quadrature<1> dealii::NonMatching::compute_cell_intersection | ( | const Triangulation< 1, 1 >::cell_iterator & | , |
const Triangulation< 1, 1 >::cell_iterator & | , | ||
const unsigned int | , | ||
const Mapping< 1, 1 > & | , | ||
const Mapping< 1, 1 > & | , | ||
const double | |||
) |
template Quadrature<2> dealii::NonMatching::compute_cell_intersection | ( | const Triangulation< 2, 2 >::cell_iterator & | , |
const Triangulation< 1, 2 >::cell_iterator & | , | ||
const unsigned int | , | ||
const Mapping< 2, 2 > & | , | ||
const Mapping< 1, 2 > & | , | ||
const double | |||
) |
template Quadrature<2> dealii::NonMatching::compute_cell_intersection | ( | const Triangulation< 2, 2 >::cell_iterator & | , |
const Triangulation< 2, 2 >::cell_iterator & | , | ||
const unsigned int | , | ||
const Mapping< 2, 2 > & | , | ||
const Mapping< 2, 2 > & | , | ||
const double | |||
) |
template Quadrature<3> dealii::NonMatching::compute_cell_intersection | ( | const Triangulation< 3, 3 >::cell_iterator & | , |
const Triangulation< 2, 3 >::cell_iterator & | , | ||
const unsigned int | , | ||
const Mapping< 3, 3 > & | , | ||
const Mapping< 2, 3 > & | , | ||
const double | |||
) |
template Quadrature<3> dealii::NonMatching::compute_cell_intersection | ( | const Triangulation< 3, 3 >::cell_iterator & | , |
const Triangulation< 3, 3 >::cell_iterator & | , | ||
const unsigned int | , | ||
const Mapping< 3, 3 > & | , | ||
const Mapping< 3, 3 > & | , | ||
const double | |||
) |
template void dealii::NonMatching::create_coupling_sparsity_pattern_with_exact_intersections | ( | const std::vector< std::tuple< typename Triangulation< 1, 1 >::cell_iterator, typename Triangulation< 1, 1 >::cell_iterator, Quadrature< 1 >>> & | intersections_info, |
const DoFHandler< 1, 1 > & | space_dh, | ||
const DoFHandler< 1, 1 > & | immersed_dh, | ||
DynamicSparsityPattern & | sparsity, | ||
const AffineConstraints< double > & | constraints, | ||
const ComponentMask & | space_comps, | ||
const ComponentMask & | immersed_comps, | ||
const AffineConstraints< double > & | immersed_constraint | ||
) |
template void dealii::NonMatching::create_coupling_sparsity_pattern_with_exact_intersections | ( | const std::vector< std::tuple< typename Triangulation< 2, 2 >::cell_iterator, typename Triangulation< 1, 2 >::cell_iterator, Quadrature< 2 >>> & | intersections_info, |
const DoFHandler< 2, 2 > & | space_dh, | ||
const DoFHandler< 1, 2 > & | immersed_dh, | ||
DynamicSparsityPattern & | sparsity, | ||
const AffineConstraints< double > & | constraints, | ||
const ComponentMask & | space_comps, | ||
const ComponentMask & | immersed_comps, | ||
const AffineConstraints< double > & | immersed_constraint | ||
) |
template void dealii::NonMatching::create_coupling_sparsity_pattern_with_exact_intersections | ( | const std::vector< std::tuple< typename Triangulation< 3, 3 >::cell_iterator, typename Triangulation< 1, 3 >::cell_iterator, Quadrature< 3 >>> & | intersections_info, |
const DoFHandler< 3, 3 > & | space_dh, | ||
const DoFHandler< 1, 3 > & | immersed_dh, | ||
DynamicSparsityPattern & | sparsity, | ||
const AffineConstraints< double > & | constraints, | ||
const ComponentMask & | space_comps, | ||
const ComponentMask & | immersed_comps, | ||
const AffineConstraints< double > & | immersed_constraint | ||
) |
template void dealii::NonMatching::create_coupling_sparsity_pattern_with_exact_intersections | ( | const std::vector< std::tuple< typename Triangulation< 2, 2 >::cell_iterator, typename Triangulation< 2, 2 >::cell_iterator, Quadrature< 2 >>> & | intersections_info, |
const DoFHandler< 2, 2 > & | space_dh, | ||
const DoFHandler< 2, 2 > & | immersed_dh, | ||
DynamicSparsityPattern & | sparsity, | ||
const AffineConstraints< double > & | constraints, | ||
const ComponentMask & | space_comps, | ||
const ComponentMask & | immersed_comps, | ||
const AffineConstraints< double > & | immersed_constraint | ||
) |
template void dealii::NonMatching::create_coupling_sparsity_pattern_with_exact_intersections | ( | const std::vector< std::tuple< typename Triangulation< 3, 3 >::cell_iterator, typename Triangulation< 2, 3 >::cell_iterator, Quadrature< 3 >>> & | intersections_info, |
const DoFHandler< 3, 3 > & | space_dh, | ||
const DoFHandler< 2, 3 > & | immersed_dh, | ||
DynamicSparsityPattern & | sparsity, | ||
const AffineConstraints< double > & | constraints, | ||
const ComponentMask & | space_comps, | ||
const ComponentMask & | immersed_comps, | ||
const AffineConstraints< double > & | immersed_constraint | ||
) |
template void dealii::NonMatching::create_coupling_sparsity_pattern_with_exact_intersections | ( | const std::vector< std::tuple< typename Triangulation< 3, 3 >::cell_iterator, typename Triangulation< 3, 3 >::cell_iterator, Quadrature< 3 >>> & | intersections_info, |
const DoFHandler< 3, 3 > & | space_dh, | ||
const DoFHandler< 3, 3 > & | immersed_dh, | ||
DynamicSparsityPattern & | sparsity, | ||
const AffineConstraints< double > & | constraints, | ||
const ComponentMask & | space_comps, | ||
const ComponentMask & | immersed_comps, | ||
const AffineConstraints< double > & | immersed_constraint | ||
) |
template void dealii::NonMatching::create_coupling_sparsity_pattern_with_exact_intersections< 2, 1, 2, TrilinosWrappers::SparsityPattern, double > | ( | std::vector< std::tuple< Triangulation< 2, 2 >::cell_iterator, Triangulation< 1, 2 >::cell_iterator, Quadrature< 2 >>, std::allocator< std::tuple< Triangulation< 2, 2 >::cell_iterator, Triangulation< 1, 2 >::cell_iterator, Quadrature< 2 >>>> const & | , |
DoFHandler< 2, 2 > const & | , | ||
DoFHandler< 1, 2 > const & | , | ||
TrilinosWrappers::SparsityPattern & | , | ||
AffineConstraints< double > const & | , | ||
ComponentMask const & | , | ||
ComponentMask const & | , | ||
AffineConstraints< double > const & | |||
) |
template void dealii::NonMatching::create_coupling_sparsity_pattern_with_exact_intersections< 2, 2, 2, TrilinosWrappers::SparsityPattern, double > | ( | std::vector< std::tuple< Triangulation< 2, 2 >::cell_iterator, Triangulation< 2, 2 >::cell_iterator, Quadrature< 2 >>, std::allocator< std::tuple< Triangulation< 2, 2 >::cell_iterator, Triangulation< 2, 2 >::cell_iterator, Quadrature< 2 >>>> const & | , |
DoFHandler< 2, 2 > const & | , | ||
DoFHandler< 2, 2 > const & | , | ||
TrilinosWrappers::SparsityPattern & | , | ||
AffineConstraints< double > const & | , | ||
ComponentMask const & | , | ||
ComponentMask const & | , | ||
AffineConstraints< double > const & | |||
) |
template void dealii::NonMatching::create_coupling_sparsity_pattern_with_exact_intersections< 3, 2, 3, TrilinosWrappers::SparsityPattern, double > | ( | std::vector< std::tuple< Triangulation< 3, 3 >::cell_iterator, Triangulation< 2, 3 >::cell_iterator, Quadrature< 3 >>, std::allocator< std::tuple< Triangulation< 3, 3 >::cell_iterator, Triangulation< 2, 3 >::cell_iterator, Quadrature< 3 >>>> const & | , |
DoFHandler< 3, 3 > const & | , | ||
DoFHandler< 2, 3 > const & | , | ||
TrilinosWrappers::SparsityPattern & | , | ||
AffineConstraints< double > const & | , | ||
ComponentMask const & | , | ||
ComponentMask const & | , | ||
AffineConstraints< double > const & | |||
) |
template void dealii::NonMatching::create_coupling_sparsity_pattern_with_exact_intersections< 3, 3, 3, TrilinosWrappers::SparsityPattern, double > | ( | std::vector< std::tuple< Triangulation< 3, 3 >::cell_iterator, Triangulation< 3, 3 >::cell_iterator, Quadrature< 3 >>, std::allocator< std::tuple< Triangulation< 3, 3 >::cell_iterator, Triangulation< 3, 3 >::cell_iterator, Quadrature< 3 >>>> const & | , |
DoFHandler< 3, 3 > const & | , | ||
DoFHandler< 3, 3 > const & | , | ||
TrilinosWrappers::SparsityPattern & | , | ||
AffineConstraints< double > const & | , | ||
ComponentMask const & | , | ||
ComponentMask const & | , | ||
AffineConstraints< double > const & | |||
) |
template void dealii::NonMatching::create_coupling_sparsity_pattern_with_exact_intersections< 3, 2, 3, SparsityPattern, double > | ( | std::vector< std::tuple< Triangulation< 3, 3 >::cell_iterator, Triangulation< 2, 3 >::cell_iterator, Quadrature< 3 >>, std::allocator< std::tuple< Triangulation< 3, 3 >::cell_iterator, Triangulation< 2, 3 >::cell_iterator, Quadrature< 3 >>>> const & | , |
DoFHandler< 3, 3 > const & | , | ||
DoFHandler< 2, 3 > const & | , | ||
SparsityPattern & | , | ||
AffineConstraints< double > const & | , | ||
ComponentMask const & | , | ||
ComponentMask const & | , | ||
AffineConstraints< double > const & | |||
) |
template void dealii::NonMatching::create_coupling_sparsity_pattern_with_exact_intersections< 2, 2, 2, SparsityPattern, double > | ( | std::vector< std::tuple< Triangulation< 2, 2 >::cell_iterator, Triangulation< 2, 2 >::cell_iterator, Quadrature< 2 >>, std::allocator< std::tuple< Triangulation< 2, 2 >::cell_iterator, Triangulation< 2, 2 >::cell_iterator, Quadrature< 2 >>>> const & | , |
DoFHandler< 2, 2 > const & | , | ||
DoFHandler< 2, 2 > const & | , | ||
SparsityPattern & | , | ||
AffineConstraints< double > const & | , | ||
ComponentMask const & | , | ||
ComponentMask const & | , | ||
AffineConstraints< double > const & | |||
) |
template void dealii::NonMatching::create_coupling_sparsity_pattern_with_exact_intersections< 1, 1, 1, SparsityPattern, double > | ( | std::vector< std::tuple< Triangulation< 1, 1 >::cell_iterator, Triangulation< 1, 1 >::cell_iterator, Quadrature< 1 >>, std::allocator< std::tuple< Triangulation< 1, 1 >::cell_iterator, Triangulation< 1, 1 >::cell_iterator, Quadrature< 1 >>>> const & | , |
DoFHandler< 1, 1 > const & | , | ||
DoFHandler< 1, 1 > const & | , | ||
SparsityPattern & | , | ||
AffineConstraints< double > const & | , | ||
ComponentMask const & | , | ||
ComponentMask const & | , | ||
AffineConstraints< double > const & | |||
) |
template void dealii::NonMatching::create_coupling_sparsity_pattern_with_exact_intersections< 3, 1, 3, SparsityPattern, double > | ( | std::vector< std::tuple< Triangulation< 3, 3 >::cell_iterator, Triangulation< 1, 3 >::cell_iterator, Quadrature< 3 >>, std::allocator< std::tuple< Triangulation< 3, 3 >::cell_iterator, Triangulation< 1, 3 >::cell_iterator, Quadrature< 3 >>>> const & | , |
DoFHandler< 3, 3 > const & | , | ||
DoFHandler< 1, 3 > const & | , | ||
SparsityPattern & | , | ||
AffineConstraints< double > const & | , | ||
ComponentMask const & | , | ||
ComponentMask const & | , | ||
AffineConstraints< double > const & | |||
) |
template void dealii::NonMatching::create_coupling_sparsity_pattern_with_exact_intersections< 2, 1, 2, SparsityPattern, double > | ( | std::vector< std::tuple< Triangulation< 2, 2 >::cell_iterator, Triangulation< 1, 2 >::cell_iterator, Quadrature< 2 >>, std::allocator< std::tuple< Triangulation< 2, 2 >::cell_iterator, Triangulation< 1, 2 >::cell_iterator, Quadrature< 2 >>>> const & | , |
DoFHandler< 2, 2 > const & | , | ||
DoFHandler< 1, 2 > const & | , | ||
SparsityPattern & | , | ||
AffineConstraints< double > const & | , | ||
ComponentMask const & | , | ||
ComponentMask const & | , | ||
AffineConstraints< double > const & | |||
) |
template void dealii::NonMatching::create_coupling_sparsity_pattern_with_exact_intersections< 3, 3, 3, SparsityPattern, double > | ( | std::vector< std::tuple< Triangulation< 3, 3 >::cell_iterator, Triangulation< 3, 3 >::cell_iterator, Quadrature< 3 >>, std::allocator< std::tuple< Triangulation< 3, 3 >::cell_iterator, Triangulation< 3, 3 >::cell_iterator, Quadrature< 3 >>>> const & | , |
DoFHandler< 3, 3 > const & | , | ||
DoFHandler< 3, 3 > const & | , | ||
SparsityPattern & | , | ||
AffineConstraints< double > const & | , | ||
ComponentMask const & | , | ||
ComponentMask const & | , | ||
AffineConstraints< double > const & | |||
) |
template void dealii::NonMatching::create_nitsche_rhs_with_exact_intersections< 2, 2, 2 > | ( | const DoFHandler< 2, 2 > & | , |
const std::vector< std::tuple< typename Triangulation< 2, 2 >::cell_iterator, typename Triangulation< 2, 2 >::cell_iterator, Quadrature< 2 >>> & | , | ||
Vector< double > & | vector, | ||
const AffineConstraints< double > & | , | ||
const Mapping< 2, 2 > & | , | ||
const Function< 2, double > & | , | ||
const Function< 2, double > & | , | ||
const double | |||
) |
template void dealii::NonMatching::create_nitsche_rhs_with_exact_intersections< 2, 1, 2 > | ( | const DoFHandler< 2, 2 > & | , |
const std::vector< std::tuple< typename Triangulation< 2, 2 >::cell_iterator, typename Triangulation< 1, 2 >::cell_iterator, Quadrature< 2 >>> & | , | ||
Vector< double > & | vector, | ||
const AffineConstraints< double > & | , | ||
const Mapping< 2, 2 > & | , | ||
const Function< 2, double > & | , | ||
const Function< 2, double > & | , | ||
const double | |||
) |
template void dealii::NonMatching::create_nitsche_rhs_with_exact_intersections< 3, 3, 3 > | ( | const DoFHandler< 3, 3 > & | , |
const std::vector< std::tuple< typename Triangulation< 3, 3 >::cell_iterator, typename Triangulation< 3, 3 >::cell_iterator, Quadrature< 3 >>> & | , | ||
Vector< double > & | vector, | ||
const AffineConstraints< double > & | , | ||
const Mapping< 3, 3 > & | , | ||
const Function< 3, double > & | , | ||
const Function< 3, double > & | , | ||
const double | |||
) |
template void dealii::NonMatching::create_nitsche_rhs_with_exact_intersections< 3, 2, 3 > | ( | const DoFHandler< 3, 3 > & | , |
const std::vector< std::tuple< typename Triangulation< 3, 3 >::cell_iterator, typename Triangulation< 2, 3 >::cell_iterator, Quadrature< 3 >>> & | , | ||
Vector< double > & | vector, | ||
const AffineConstraints< double > & | , | ||
const Mapping< 3, 3 > & | , | ||
const Function< 3, double > & | , | ||
const Function< 3, double > & | , | ||
const double | |||
) |