PRISMS-PF Manual
Loading...
Searching...
No Matches
pde_operator_base.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: © 2025 PRISMS Center at the University of Michigan
2// SPDX-License-Identifier: GNU Lesser General Public Version 2.1
3
4#pragma once
5
6#include <deal.II/base/exceptions.h>
7
11#include <prismspf/core/types.h>
12
14
15#include <prismspf/config.h>
16
18
22template <unsigned int dim, unsigned int degree, typename number>
24{
25public:
26 using SizeType = dealii::VectorizedArray<number>;
27
31 explicit PDEOperatorBase(const UserInputParameters<dim> &_user_inputs,
32 const PhaseFieldTools<dim> &_pf_tools)
33 : user_inputs(&_user_inputs)
34 , pf_tools(&_pf_tools)
35 {}
36
40 virtual ~PDEOperatorBase() = default;
41
45 virtual void
46 set_initial_condition([[maybe_unused]] const unsigned int &index,
47 [[maybe_unused]] const unsigned int &component,
48 [[maybe_unused]] const dealii::Point<dim> &point,
49 [[maybe_unused]] number &scalar_value,
50 [[maybe_unused]] number &vector_component_value) const
51 {}
52
57 virtual void
58 set_dirichlet([[maybe_unused]] const unsigned int &index,
59 [[maybe_unused]] const unsigned int &boundary_id,
60 [[maybe_unused]] const unsigned int &component,
61 [[maybe_unused]] const dealii::Point<dim> &point,
62 [[maybe_unused]] number &scalar_value,
63 [[maybe_unused]] number &vector_component_value) const
64 {
65 this->set_initial_condition(index,
66 component,
67 point,
68 scalar_value,
69 vector_component_value);
70 }
71
75 virtual void
76 compute_rhs([[maybe_unused]] FieldContainer<dim, degree, number> &variable_list,
77 [[maybe_unused]] const SimulationTimer &sim_timer,
78 [[maybe_unused]] unsigned int solver_id) const
79 {}
80
84 virtual void
85 compute_lhs([[maybe_unused]] FieldContainer<dim, degree, number> &variable_list,
86 [[maybe_unused]] const SimulationTimer &sim_timer,
87 [[maybe_unused]] unsigned int solver_id) const
88 {}
89
93 [[nodiscard]] const UserInputParameters<dim> &
95 {
96 Assert(user_inputs != nullptr, dealii::StandardExceptions::ExcNotInitialized());
97 return *user_inputs;
98 }
99
103 [[nodiscard]] const PhaseFieldTools<dim> &
105 {
106 Assert(pf_tools != nullptr, dealii::StandardExceptions::ExcNotInitialized());
107 return *pf_tools;
108 }
109
110private:
115
120};
121
122PRISMS_PF_END_NAMESPACE
This class permits the access of a subset of indexed fields and gives an error if any non-allowed fie...
Definition field_container.h:47
PDEOperatorBase(const UserInputParameters< dim > &_user_inputs, const PhaseFieldTools< dim > &_pf_tools)
Constructor.
Definition pde_operator_base.h:31
dealii::VectorizedArray< number > SizeType
Definition pde_operator_base.h:26
virtual void set_dirichlet(const unsigned int &index, const unsigned int &boundary_id, const unsigned int &component, const dealii::Point< dim > &point, number &scalar_value, number &vector_component_value) const
User-implemented function for setting Dirichlet boundary conditions. Default behavior is to call init...
Definition pde_operator_base.h:58
const UserInputParameters< dim > * user_inputs
The user-inputs.
Definition pde_operator_base.h:114
virtual void compute_rhs(FieldContainer< dim, degree, number > &variable_list, const SimulationTimer &sim_timer, unsigned int solver_id) const
User-implemented class for the RHS of explicit equations.
Definition pde_operator_base.h:76
const PhaseFieldTools< dim > * pf_tools
Phase field tools.
Definition pde_operator_base.h:119
const UserInputParameters< dim > & get_user_inputs() const
Get the user inputs (constant reference).
Definition pde_operator_base.h:94
virtual ~PDEOperatorBase()=default
Destructor.
const PhaseFieldTools< dim > & get_pf_tools() const
Get the phase field tools (constant reference).
Definition pde_operator_base.h:104
virtual void set_initial_condition(const unsigned int &index, const unsigned int &component, const dealii::Point< dim > &point, number &scalar_value, number &vector_component_value) const
User-implemented class for the setting initial conditions.
Definition pde_operator_base.h:46
virtual void compute_lhs(FieldContainer< dim, degree, number > &variable_list, const SimulationTimer &sim_timer, unsigned int solver_id) const
User-implemented class for the RHS of nonexplicit equations.
Definition pde_operator_base.h:85
Definition simulation_timer.h:13
Definition user_input_parameters.h:32
Definition conditional_ostreams.cc:20
This class contains mutable utilities for phase field problems.
Definition phase_field_tools.h:21