DUNE-ACFEM (2.5.1)

neumannmodel.hh
1 #ifndef __DUNE_ACFEM_MODELS_MODULES_NEUMANNMODEL_HH__
2 #define __DUNE_ACFEM_MODELS_MODULES_NEUMANNMODEL_HH__
3 
4 #include "../modelexpressionbase.hh"
5 
6 namespace Dune {
7 
8  namespace ACFem {
9 
39  template<class GridFunction, class Indicator = EntireBoundaryIndicatorType>
41  : public ModelExpression<NeumannBoundaryModel<GridFunction, Indicator> >
42  {
45  public:
46  typedef Indicator OuterIndicatorType;
47  typedef typename TraitsType::NeumannBoundaryFunctionType NeumannBoundaryFunctionType;
48  typedef typename TraitsType::NeumannIndicatorType NeumannIndicatorType;
49  typedef typename TraitsType::FunctionSpaceType FunctionSpaceType;
50  typedef typename TraitsType::GridPartType GridPartType;
51  typedef typename GridPartType::IntersectionType IntersectionType;
52  typedef typename FunctionSpaceType::DomainType DomainType;
53  typedef typename FunctionSpaceType::RangeType RangeType;
54  typedef typename FunctionSpaceType::JacobianRangeType JacobianRangeType;
55 
56  NeumannBoundaryModel(const Fem::Function<FunctionSpaceType, GridFunction>& uExt,
57  const
59  indicator = NeumannIndicatorType(),
60  const std::string& name = "")
61  : bndryFunction_(asImp(indicator) * asImp(uExt)),
62  name_(name == "" ? "Neumann(" + bndryFunction_().name() + ")" : name)
63  {}
64 
66  std::string name() const
67  {
68  return name_;
69  }
70 
72  bool setIntersection(const IntersectionType& intersection) const
73  {
74  // Really, this is for Robin data, so we return indeed false.
75  return false;
76  }
77 
79  NeumannIndicatorType neumannIndicator() const
80  {
81  return bndryFunction_().indicator();
82  }
83 
89  NeumannBoundaryFunctionType neumannBoundaryFunction(const GridPartType& gridPart) const
90  {
91  return bndryFunction_();
92  }
93 
95  const GridPartType& gridPart() const
96  {
97  return bndryFunction_().gridPart();
98  }
99 
100  protected:
102  const std::string name_;
103  };
104 
105  template<class GridFunction, class Indicator>
106  struct ModelTraits<NeumannBoundaryModel<GridFunction, Indicator> >
107  : public DefaultModelTraits<typename GridFunction::DiscreteFunctionSpaceType::FunctionSpaceType,
108  typename GridFunction::GridPartType>
109  {
110  protected:
111  typedef
112  DefaultModelTraits<typename GridFunction::DiscreteFunctionSpaceType::FunctionSpaceType,
113  typename GridFunction::GridPartType>
114  BaseType;
115  typedef Indicator OuterIndicatorType;
116  public:
118  isLinear = true,
119  isSymmetric = true,
120  isSemiDefinite = true
121  };
122 
124  typedef
125  decltype(std::declval<OuterIndicatorType>() * std::declval<GridFunction>())
127 
129  typedef
130  typename NeumannBoundaryFunctionType::IndicatorType
132  };
133 
135  template<class GridFunction>
136  class NeumannBoundaryModel<GridFunction, EmptyBoundaryIndicatorType>
137  : public ZeroModel<typename GridFunction::FunctionSpaceType, typename GridFunction::GridPartType>
138  {
139  typedef
141  BaseType;
142  public:
143  typedef typename GridFunction::FunctionSpaceType FunctionSpaceType;
144  typedef EmptyBoundaryIndicatorType NeumannIndicatorType;
145 
146  NeumannBoundaryModel(const Fem::Function<FunctionSpaceType, GridFunction>& function,
148  const std::string& name = "")
149  : BaseType(name == "" ? ("Neumann(" + asImp(function).name() + ")") : name)
150  {}
151  };
152 
154 
167  template<class GridFunction, class Indicator = EntireBoundaryIndicatorType>
169  neumannBoundaryModel(const Fem::Function<typename GridFunction::FunctionSpaceType, GridFunction>& values,
170  const BoundaryIndicatorInterface<Indicator>& where = Indicator())
171  {
172  return NeumannBoundaryModel<GridFunction, Indicator>(values, where);
173  }
174 
199  template<class Object, class Indicator = EntireBoundaryIndicatorType>
200  NeumannBoundaryModel<ZeroGridFunction<typename Object::FunctionSpaceType,
201  typename Object::GridPartType>, Indicator>
202  neumannZeroModel(const Object& object,
203  const BoundaryIndicatorInterface<Indicator>& where = Indicator())
204  {
205  typedef
206  ZeroGridFunction<typename Object::FunctionSpaceType,
207  typename Object::GridPartType>
208  ZeroType;
209  return neumannBoundaryModel(ZeroType(object.gridPart()), where);
210  }
211 
213 
215 
217 
218  } // namespace ACFem
219 
220 } //Namespace Dune
221 
222 
223 #endif // __DUNE_ACFEM_MODELS_MODULES_NEUMANNMODEL_HH__
A simple interface class for a boundary-indicator.
Definition: boundaryindicator.hh:41
A function with potentially partial support on the boundary.
Definition: boundarysupportedfunction.hh:234
A Neumann-boundary model.
Definition: neumannmodel.hh:42
A grid-function always returning 0.
Definition: constantfunction.hh:352
Define a simple zero model to optimize expression templates.
Definition: zeromodel.hh:28
const Implementation & asImp(const Fem::BartonNackmanInterface< Interface, Implementation > &arg)
Up-cast to the implementation for any Fem::BartonNackmanInterface.
Definition: expressionoperations.hh:71
NeumannBoundaryModel< GridFunction, Indicator > neumannBoundaryModel(const Fem::Function< typename GridFunction::FunctionSpaceType, GridFunction > &values, const BoundaryIndicatorInterface< Indicator > &where=Indicator())
Generate a NeumannBoundaryModel from given grid-function and boundary indicator.
Definition: neumannmodel.hh:169
NeumannIndicatorType neumannIndicator() const
Generate an object to identify parts of the boundary subject to Neumann boundary conditions.
Definition: neumannmodel.hh:79
NeumannBoundaryFunctionType::IndicatorType NeumannIndicatorType
Something satisfying the BoundaryIndicatorInterface.
Definition: neumannmodel.hh:131
NeumannBoundaryModel< ZeroGridFunction< typename Object::FunctionSpaceType, typename Object::GridPartType >, Indicator > neumannZeroModel(const Object &object, const BoundaryIndicatorInterface< Indicator > &where=Indicator())
Generate homogeneous Neumann boundary conditions fitting the specified object.
Definition: neumannmodel.hh:202
const GridPartType & gridPart() const
Return the GridPart (non-interface method)
Definition: neumannmodel.hh:95
bool setIntersection(const IntersectionType &intersection) const
Definition: neumannmodel.hh:72
NeumannBoundaryFunctionType neumannBoundaryFunction(const GridPartType &gridPart) const
Generate an instance of a class defining Neumann boundary values as a Fem grid-function.
Definition: neumannmodel.hh:89
std::string name() const
Print a descriptive name for debugging and output.
Definition: neumannmodel.hh:66
Traits-template which has to be specialized for each individual model.
Definition: modelinterface.hh:48
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.80.0 (May 16, 22:29, 2024)