deal2lkit: A ToolKit library for Deal.II
parsed_dirichlet_bcs.h
Go to the documentation of this file.
1 //-----------------------------------------------------------
2 //
3 // Copyright (C) 2015 by the deal2lkit authors
4 //
5 // This file is part of the deal2lkit library.
6 //
7 // The deal2lkit library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE at
12 // the top level of the deal2lkit distribution.
13 //
14 //-----------------------------------------------------------
15 
16 #ifndef _d2k_parsed_dirichlet_bcs_h
17 #define _d2k_parsed_dirichlet_bcs_h
18 
19 #include <deal2lkit/config.h>
26 
27 using namespace dealii;
28 
29 D2K_NAMESPACE_OPEN
30 
91 template <int dim, int spacedim=dim>
92 class ParsedDirichletBCs : public ParsedMappedFunctions<spacedim>
93 {
94 public:
112  ParsedDirichletBCs (const std::string &name = "Dirichlet BCs",
113  const unsigned int &n_components=1,
114  const std::string &component_names = "",
115  const std::string &default_id_components = "0=ALL",
116  const std::string &default_id_functions = "",
117  const std::string &default_constants = "");
118 
122  virtual void declare_parameters (ParameterHandler &prm);
123 
127  virtual void parse_parameters_call_back ();
128 
135  void interpolate_boundary_values (const DoFHandler<dim,spacedim> &dof_handler,
136  ConstraintMatrix &constraints) const;
137 
145  const DoFHandler<dim,spacedim> &dof_handler,
146  ConstraintMatrix &constraints) const;
147 
154  void interpolate_boundary_values (const DoFHandler<dim,spacedim> &dof_handler,
155  std::map<types::global_dof_index,double> &d_dofs) const;
156 
164  const DoFHandler<dim,spacedim> &dof_handler,
165  std::map<types::global_dof_index,double> &d_dofs) const;
166 
173  void project_boundary_values (const DoFHandler<dim,spacedim> &dof_handler,
174  const Quadrature<dim-1> &quadrature,
175  ConstraintMatrix &constraints) const;
176 
183  void project_boundary_values (const Mapping<dim,spacedim> &mapping,
184  const DoFHandler<dim,spacedim> &dof_handler,
185  const Quadrature<dim-1> &quadrature,
186  ConstraintMatrix &constraints) const;
187 
194  void project_boundary_values (const DoFHandler<dim,spacedim> &dof_handler,
195  const Quadrature<dim-1> &quadrature,
196  std::map<types::global_dof_index,double> &projected_bv) const;
197 
204  void project_boundary_values (const Mapping<dim,spacedim> &mapping,
205  const DoFHandler<dim,spacedim> &dof_handler,
206  const Quadrature<dim-1> &quadrature,
207  std::map<types::global_dof_index,double> &projected_bv) const;
208 
209 
221  ConstraintMatrix &constraints) const;
233  const Mapping< dim, spacedim > &mapping,
234  ConstraintMatrix &constraints) const;
249  ConstraintMatrix &constraints) const;
250 
265  const Mapping< dim, spacedim > &mapping,
266  ConstraintMatrix &constraints) const;
267 private:
271  const unsigned int n_components;
272 
273 };
274 
275 D2K_NAMESPACE_CLOSE
276 
277 #endif
278 
void compute_no_normal_flux_constraints(const DoFHandlerType< dim, spacedim > &dof_handler, const unsigned int first_vector_component, const std::set< types::boundary_id > &boundary_ids, ConstraintMatrix &constraints, const Mapping< dim, spacedim > &mapping=StaticMappingQ1< dim >::mapping)
Parsed DirichletBCs.
void compute_nonzero_normal_flux_constraints(const DoFHandlerType< dim, spacedim > &dof_handler, const unsigned int first_vector_component, const std::set< types::boundary_id > &boundary_ids, typename FunctionMap< spacedim >::type &function_map, ConstraintMatrix &constraints, const Mapping< dim, spacedim > &mapping=StaticMappingQ1< dim >::mapping)
ParsedMappedFunctions object.
void project_boundary_values(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const std::map< types::boundary_id, const Function< spacedim, number > *> &boundary_functions, const Quadrature< dim-1 > &q, ConstraintMatrix &constraints, std::vector< unsigned int > component_mapping=std::vector< unsigned int >())
void interpolate_boundary_values(const Mapping< dim, spacedim > &mapping, const DoFHandlerType< dim, spacedim > &dof, const std::map< types::boundary_id, const Function< spacedim, number > *> &function_map, ConstraintMatrix &constraints, const ComponentMask &component_mask=ComponentMask())
const unsigned int n_components
Number of components of the underlying Function objects.