Program Listing for File reference_cross_section.h¶
↰ Return to documentation for file (include/reference_cross_section.h)
// ---------------------------------------------------------------------
//
// Copyright (C) 2024 by Luca Heltai
//
// This file is part of the ImmersX application, based on
// the deal.II library.
//
// The ImmersX application is free software; you can use
// it, redistribute it, and/or modify it under the terms of the Apache-2.0
// License WITH LLVM-exception as published by the Free Software Foundation;
// either version 3.0 of the License, or (at your option) any later version. The
// full text of the license can be found in the file LICENSE.md at the top level
// of the ImmersX distribution.
//
// ---------------------------------------------------------------------
#ifndef rdlm_reference_inclusion
#define rdlm_reference_inclusion
#include <deal.II/base/parameter_acceptor.h>
#include <deal.II/base/patterns.h>
#include <deal.II/base/polynomials_p.h>
#include <deal.II/base/quadrature_lib.h>
#include <deal.II/dofs/dof_handler.h>
#include <deal.II/fe/fe.h>
#include <deal.II/fe/fe_system.h>
#include <deal.II/fe/mapping_q.h>
#include <deal.II/grid/tria.h>
#include <deal.II/lac/sparse_matrix.h>
#include <deal.II/lac/sparsity_pattern.h>
#include <deal.II/lac/vector.h>
#include <fstream>
using namespace dealii;
template <int dim, int spacedim = dim, int n_components = 1>
struct ReferenceCrossSectionParameters : public ParameterAcceptor
{
ReferenceCrossSectionParameters();
unsigned int refinement_level = 1;
std::string inclusion_type = "hyper_ball";
unsigned int inclusion_degree = 0;
mutable std::vector<unsigned int> selected_coefficients;
};
template <int dim, int spacedim = dim, int n_components = 1>
class ReferenceCrossSection
{
public:
ReferenceCrossSection(
const ReferenceCrossSectionParameters<dim, spacedim, n_components> &par);
const Quadrature<spacedim> &
get_global_quadrature() const;
const std::vector<Vector<double>> &
get_basis_functions() const;
const double &
shape_value(const unsigned int i,
const unsigned int q,
const unsigned int comp) const;
const SparseMatrix<double> &
get_mass_matrix() const;
unsigned int
n_selected_basis() const;
unsigned int
max_n_basis() const;
Quadrature<spacedim>
get_transformed_quadrature(const Point<spacedim> &new_origin,
const Tensor<1, spacedim> &new_vertical,
const double scale) const;
unsigned int
n_quadrature_points() const;
void
initialize();
double
measure(const double scale = 1.0) const;
private:
void
make_grid();
void
setup_dofs();
void
compute_basis();
const ReferenceCrossSectionParameters<dim, spacedim, n_components> ∥
PolynomialsP<spacedim> polynomials;
Triangulation<dim, spacedim> triangulation;
QGauss<dim> quadrature_formula;
FESystem<dim, spacedim> fe;
MappingQ<dim, spacedim> mapping;
DoFHandler<dim, spacedim> dof_handler;
std::vector<Vector<double>> basis_functions;
std::vector<Vector<double>> selected_basis_functions;
SparsityPattern sparsity_pattern;
SparseMatrix<double> mass_matrix;
Quadrature<spacedim> global_quadrature;
double reference_measure;
FullMatrix<double> basis_functions_on_qpoints;
};
#endif