Line data Source code
1 : /*
2 : * Copyright (c) 2013-2015: G-CSC, Goethe University Frankfurt
3 : * Author: Andreas Vogel
4 : *
5 : * This file is part of UG4.
6 : *
7 : * UG4 is free software: you can redistribute it and/or modify it under the
8 : * terms of the GNU Lesser General Public License version 3 (as published by the
9 : * Free Software Foundation) with the following additional attribution
10 : * requirements (according to LGPL/GPL v3 §7):
11 : *
12 : * (1) The following notice must be displayed in the Appropriate Legal Notices
13 : * of covered and combined works: "Based on UG4 (www.ug4.org/license)".
14 : *
15 : * (2) The following notice must be displayed at a prominent place in the
16 : * terminal output of covered works: "Based on UG4 (www.ug4.org/license)".
17 : *
18 : * (3) The following bibliography is recommended for citation and must be
19 : * preserved in all covered files:
20 : * "Reiter, S., Vogel, A., Heppner, I., Rupp, M., and Wittum, G. A massively
21 : * parallel geometric multigrid solver on hierarchically distributed grids.
22 : * Computing and visualization in science 16, 4 (2013), 151-164"
23 : * "Vogel, A., Reiter, S., Rupp, M., Nägel, A., and Wittum, G. UG4 -- a novel
24 : * flexible software system for simulating pde based models on high performance
25 : * computers. Computing and visualization in science 16, 4 (2013), 165-179"
26 : *
27 : * This program is distributed in the hope that it will be useful,
28 : * but WITHOUT ANY WARRANTY; without even the implied warranty of
29 : * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
30 : * GNU Lesser General Public License for more details.
31 : */
32 :
33 : #ifndef __H__UG__LIB_DISC__SPATIAL_DISC__ELEM_DISC__NEUMANN_BOUNDARY___NEUMANN_BOUNDARY_FV1__
34 : #define __H__UG__LIB_DISC__SPATIAL_DISC__ELEM_DISC__NEUMANN_BOUNDARY___NEUMANN_BOUNDARY_FV1__
35 :
36 : // other ug4 modules
37 : #include "common/common.h"
38 :
39 : // library intern headers
40 : #include "../neumann_boundary_base.h"
41 :
42 : namespace ug{
43 :
44 : template<typename TDomain>
45 : class NeumannBoundaryFV1
46 : : public NeumannBoundaryBase<TDomain>
47 : {
48 : private:
49 : /// Base class type
50 : typedef NeumannBoundaryBase<TDomain> base_type;
51 :
52 : /// Base class type
53 : typedef NeumannBoundaryFV1<TDomain> this_type;
54 :
55 : /// error estimator type
56 : typedef SideAndElemErrEstData<TDomain> err_est_type;
57 :
58 : public:
59 : /// World dimension
60 : static const int dim = base_type::dim;
61 :
62 : public:
63 : /// default constructor
64 : NeumannBoundaryFV1(const char* function);
65 :
66 : /// add a boundary value
67 : /// \{
68 : void add(SmartPtr<CplUserData<number, dim> > data, const char* BndSubsets, const char* InnerSubsets);
69 : void add(SmartPtr<CplUserData<number, dim, bool> > user, const char* BndSubsets, const char* InnerSubsets);
70 : void add(SmartPtr<CplUserData<MathVector<dim>, dim> > user, const char* BndSubsets, const char* InnerSubsets);
71 : /// \}
72 :
73 : protected:
74 : using typename base_type::Data;
75 :
76 : /// Unconditional scalar user data
77 : struct NumberData : public base_type::Data
78 : {
79 0 : NumberData(SmartPtr<CplUserData<number, dim> > data,
80 : std::string BndSubsets, std::string InnerSubsets)
81 0 : : base_type::Data(BndSubsets, InnerSubsets)
82 : {
83 0 : import.set_data(data);
84 0 : }
85 :
86 : template<typename TElem, typename TFVGeom>
87 : void extract_bip(const TFVGeom& geo);
88 :
89 : template <typename TElem, typename TFVGeom>
90 : void lin_def(const LocalVector& u,
91 : std::vector<std::vector<number> > vvvLinDef[],
92 : const size_t nip);
93 :
94 : template <int refDim>
95 : void set_local_ips(const MathVector<refDim>* ips, std::size_t nIPs)
96 0 : {import.template set_local_ips<refDim>(ips, nIPs);}
97 :
98 0 : void set_global_ips(const MathVector<dim>* ips, std::size_t nIPs)
99 0 : {import.set_global_ips(ips, nIPs);}
100 :
101 : template <int refDim>
102 : std::vector<MathVector<refDim> >* local_ips();
103 :
104 : DataImport<number, dim> import;
105 : std::vector<MathVector<3> > vLocIP_dim3;
106 : std::vector<MathVector<2> > vLocIP_dim2; // might have Neumann bnd for lower-dim elements!
107 : std::vector<MathVector<1> > vLocIP_dim1;
108 : std::vector<MathVector<dim> > vGloIP;
109 : };
110 :
111 : /// Conditional scalar user data
112 0 : struct BNDNumberData : public base_type::Data
113 : {
114 0 : BNDNumberData(SmartPtr<CplUserData<number, dim, bool> > functor_,
115 : std::string BndSubsets, std::string InnerSubsets)
116 0 : : base_type::Data(BndSubsets, InnerSubsets), functor(functor_) {}
117 :
118 : SmartPtr<CplUserData<number, dim, bool> > functor;
119 : };
120 :
121 : /// Unconditional vector user data
122 0 : struct VectorData : public base_type::Data
123 : {
124 0 : VectorData(SmartPtr<CplUserData<MathVector<dim>, dim> > functor_,
125 : std::string BndSubsets, std::string InnerSubsets)
126 0 : : base_type::Data(BndSubsets, InnerSubsets), functor(functor_) {}
127 :
128 : SmartPtr<CplUserData<MathVector<dim>, dim> > functor;
129 : };
130 :
131 : std::vector<NumberData> m_vNumberData;
132 : std::vector<BNDNumberData> m_vBNDNumberData;
133 : std::vector<VectorData> m_vVectorData;
134 :
135 : void update_subset_groups();
136 : using base_type::update_subset_groups;
137 :
138 : /// current inner subset
139 : int m_si;
140 :
141 : public:
142 : /// type of trial space for each function used
143 : virtual void prepare_setting(const std::vector<LFEID>& vLfeID, bool bNonRegularGrid);
144 :
145 : protected:
146 : /// assembling functions for fv1
147 : /// \{
148 : template<typename TElem, typename TFVGeom>
149 : void prep_elem_loop(const ReferenceObjectID roid, const int si);
150 : template<typename TElem, typename TFVGeom>
151 : void prep_elem(const LocalVector& u, GridObject* elem, const ReferenceObjectID roid, const MathVector<dim> vCornerCoords[]);
152 : template<typename TElem, typename TFVGeom>
153 : void fsh_elem_loop();
154 : template<typename TElem, typename TFVGeom>
155 : void add_rhs_elem(LocalVector& d, GridObject* elem, const MathVector<dim> vCornerCoords[]);
156 :
157 : /// prepares the loop over all elements of one type for the computation of the error estimator
158 : template <typename TElem, typename TFVGeom>
159 : void prep_err_est_elem_loop(const ReferenceObjectID roid, const int si);
160 :
161 : /// prepares the element for assembling the error estimator
162 : template <typename TElem, typename TFVGeom>
163 : void prep_err_est_elem(const LocalVector& u, GridObject* elem, const MathVector<dim> vCornerCoords[]);
164 :
165 : /// computes the error estimator contribution for one element
166 : template <typename TElem, typename TFVGeom>
167 : void compute_err_est_rhs_elem(GridObject* elem, const MathVector<dim> vCornerCoords[], const number& scale);
168 :
169 : /// postprocesses the loop over all elements of one type in the computation of the error estimator
170 : template <typename TElem, typename TFVGeom>
171 : void fsh_err_est_elem_loop();
172 :
173 :
174 : /// \}
175 :
176 : static const int _C_ = 0;
177 :
178 : private:
179 : void register_all_funcs(bool bHang);
180 : template <typename TElem, typename TFVGeom> void register_func();
181 : };
182 :
183 : } // end namespac ug
184 :
185 : #endif /*__H__UG__LIB_DISC__SPATIAL_DISC__ELEM_DISC__NEUMANN_BOUNDARY___NEUMANN_BOUNDARY_FV1__*/
|