4#ifndef nonuniform_dirichlet_h
5#define nonuniform_dirichlet_h
7#include <deal.II/base/function.h>
8#include <deal.II/base/point.h>
9#include <deal.II/lac/vector.h>
11#include <prismspf/config.h>
12#include <prismspf/core/type_enums.h>
13#include <prismspf/user_inputs/user_input_parameters.h>
15PRISMS_PF_BEGIN_NAMESPACE
26template <
int dim, fieldType field_type = fieldType::SCALAR>
34 const unsigned int &_boundary_id,
41 value(
const dealii::Point<dim> &p,
const unsigned int component = 0)
const override;
47 vector_value(
const dealii::Point<dim> &p, dealii::Vector<double> &
value)
const override;
50 const unsigned int index;
52 const unsigned int boundary_id;
59template <
int dim, fieldType field_type>
61 const unsigned int &_index,
62 const unsigned int &_boundary_id,
64 : dealii::Function<dim>((field_type == fieldType::VECTOR) ? dim : 1)
66 , boundary_id(_boundary_id)
67 , user_inputs(&_user_inputs)
70template <
int dim, fieldType field_type>
73 const dealii::Point<dim> &p,
74 [[maybe_unused]]
const unsigned int component)
const
77 double scalar_value = 0.0;
78 dealii::Vector<double> vector_value(dim);
81 custom_nonuniform_dirichlet.set_nonuniform_dirichlet(index,
92template <
int dim, fieldType field_type>
95 dealii::Vector<double> &value)
const
98 double scalar_value = 0.0;
99 dealii::Vector<double> vector_value(dim);
102 for (
unsigned int i = 0; i < dim; i++)
104 custom_nonuniform_dirichlet.set_nonuniform_dirichlet(index,
113 value = vector_value;
134 const unsigned int &boundary_id,
135 const unsigned int &component,
136 const dealii::Point<dim> &point,
137 double &scalar_value,
138 double &vector_component_value,
142PRISMS_PF_END_NAMESPACE