ug4
Loading...
Searching...
No Matches
obstacle_in_normal_dir_impl.h
Go to the documentation of this file.
1/*
2 * Copyright (c) 2013-2015: G-CSC, Goethe University Frankfurt
3 * Author: Raphael Prohl
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_ALGEBRA__OPERATOR__PRECONDITIONER__PROJECTED_GAUSS_SEIDEL__OBSTACLE_IN_NORMAL_DIR_IMPL__
34#define __H__UG__LIB_ALGEBRA__OPERATOR__PRECONDITIONER__PROJECTED_GAUSS_SEIDEL__OBSTACLE_IN_NORMAL_DIR_IMPL__
35
37
41
42namespace ug{
43
44template <typename TDomain, typename TAlgebra>
45void
47{
48 // for debugging
49 MathVector<dim> normal;
50 normal[0] = 2.0;
51 if (dim > 1) normal[1] = 1.5;
52 if (dim > 2) normal[2] = -3.0;
53
54 MathVector<dim> transformedONB[dim];
55 transform_eulerian_coord_sys(transformedONB, normal);
56
57}
58
59template <typename TDomain, typename TAlgebra>
60void
62 const MathVector<dim>& firstTransformedBaseVec)
63{
64 // create first unit eulerian base vector
65 MathVector<dim> unityX;
66 unityX[0] = 1.0;
67 if (dim > 1) unityX[1] = 0.0;
68 if (dim > 2) unityX[2] = 0.0;
69
70 // normalize first base vector of the transformed system
71 MathVector<dim> normalizedFirstBaseVec;
72 const number normOfFirstBaseVec = VecLength(firstTransformedBaseVec);
73 VecScale(normalizedFirstBaseVec, firstTransformedBaseVec, 1.0/normOfFirstBaseVec);
74
75 // compute vector, which is orthogonal to the householder hypersphere
76 MathVector<dim> orthoVec;
77 VecScaleAdd(orthoVec, 0.5, unityX, -0.5, normalizedFirstBaseVec);
78
79 // compute householder matrix
81 MatHouseholder(hMat, orthoVec);
82
83 // get and store transformed orthonormal base vectors
84 MathVector<dim> transBaseVec;
85 for(size_t i = 0; i < (size_t)dim; ++i)
86 {
87 for(size_t j = 0; j < (size_t)dim; ++j){
88 transBaseVec[j] = hMat(i,j);
89 }
90 transformedONB[i] = transBaseVec;
91 }
92}
93
94template <int dim> struct face_type_traits
95{
96 typedef void face_type0;
97 typedef void face_type1;
98 typedef void DimFEGeo;
99};
100
101template <> struct face_type_traits<1>
102{
106};
107
108template <> struct face_type_traits<2>
109{
113};
114
115template <> struct face_type_traits<3>
116{
120};
121
122template <typename TDomain, typename TAlgebra>
123template <typename TElem, typename TIterator>
124void
126adjust_sol_and_cor_elem(TIterator iterBegin,
127 TIterator iterEnd, value_type& sol_i, value_type& c_i, bool& dofIsActive,
128 const DoFIndex& dof)
129{
130 typedef typename face_type_traits<dim>::face_type0 face_type0;
131 typedef typename face_type_traits<dim>::face_type1 face_type1;
132
133 // storage for corner coordinates
134 vector<MathVector<dim> > vCorner;
135 vector<MathVector<dim> > vCoPos;
136
137 // outer normal vector
138 MathVector<dim> normal;
139
140 for(TIterator iter = iterBegin; iter != iterEnd; ++iter)
141 {
142 TElem* elem = *iter;
143
144 // reference object type
145 ReferenceObjectID roid = elem->reference_object_id();
146
147 const DimReferenceElement<dim-1>& rRefElem
148 = ReferenceElementProvider::get<dim-1>(roid);
149
150 // get corners of element
151 CollectCornerCoordinates(vCorner, *elem, *m_spDomain);
152
153 // here the ordering of the corners in the reference element is exploited
154 // in order to compute the outer normal later on
155 int nCorner = (int)vCorner.size();
156 for (int i = 0; i < nCorner; ++i)
157 vCoPos.push_back(vCorner[rRefElem.id(dim-1, 0, 0, i)]);
158
159 if ((int)vCorner.size() == dim)
160 ElementNormal<face_type0, dim>(normal, &vCoPos[0]);
161 else
162 ElementNormal<face_type1, dim>(normal, &vCoPos[0]);
163
164 /*for (int i = 0; i < (int)vCorner.size(); ++i)
165 UG_LOG("coPos: " << coPos[i] << "\n");
166 UG_LOG("normal: " << normal << "\n");
167
168 //const number normOfNormal = VecLength(normal);
169 //const value_type tmpVSol = sol_i + c_i;
170 //const number uTimesNormal = tmpSol / normOfNormal;
171 UG_LOG("\n");*/
172 }
173
174 const size_t comp = dof[1];
175
176 // tmpSol := u_{s-1/2} = u_{s-1} + c
177 const number tmpSol = BlockRef(sol_i, comp) + BlockRef(c_i, comp);
178
179 // get obstacle value corresponding to the dof
180 const number obsVal = m_mObstacleValues[dof];
181
182 // check, if dof is active (tmpSol >= obsVal)
183 // TODO: check if u * n > g, i.e. tmpSol * n > g!
184 if (!(tmpSol < obsVal))
185 {
186 // is active DoF
187 m_vActiveDofs.push_back(dof);
188
189 // adjust correction & set solution to obstacle-value
190 BlockRef(c_i, comp) = obsVal - BlockRef(sol_i, comp);
191 BlockRef(sol_i, comp) = obsVal;
192 dofIsActive = true;
193 }
194}
195
196template <typename TDomain, typename TAlgebra>
197void
199adjust_sol_and_cor(value_type& sol_i, value_type& c_i, bool& dofIsActive,
200 const DoFIndex& dof)
201{
202 // loop over all obstacle subsets
203 for (vector<int>::iterator obsSI = m_vObsSubsets.begin();
204 obsSI != m_vObsSubsets.end(); ++obsSI)
205 {
206 const int si = *obsSI;
207 UG_LOG("si: " <<si<<"\n");
208 const int subsetDim = DimensionOfSubset(*m_spDD->subset_handler(), si);
209
210 switch(subsetDim)
211 {
212 case 0:
213 break;
214 case 1:
215 adjust_sol_and_cor_elem<RegularEdge>
216 (m_spDD->template begin<RegularEdge>(si), m_spDD->template end<RegularEdge>(si),
217 sol_i, c_i, dofIsActive, dof);
218 break;
219 case 2:
220 adjust_sol_and_cor_elem<Triangle>
221 (m_spDD->template begin<Triangle>(si), m_spDD->template end<Triangle>(si),
222 sol_i, c_i, dofIsActive, dof);
223 adjust_sol_and_cor_elem<Quadrilateral>
224 (m_spDD->template begin<Quadrilateral>(si), m_spDD->template end<Quadrilateral>(si),
225 sol_i, c_i, dofIsActive, dof);
226 break;
227 default:
228 UG_THROW("ObstacleInNormalDir::adjust_sol_and_cor:"
229 "SubsetDimension "<< subsetDim <<" (subset="<< si <<") not supported.");
230 }
231 }
232}
233
234template <typename TDomain, typename TAlgebra>
235void
238{
239 for (vector<MultiIndex<2> >::iterator itActiveInd = m_vActiveDofs.begin();
240 itActiveInd < m_vActiveDofs.end(); ++itActiveInd)
241 {
242 // check, if Ax >= b. For that case the new defect is set to zero,
243 // since all equations/constraints are fulfilled
244 number defect = BlockRef(d[(*itActiveInd)[0]], (*itActiveInd)[1]);
245 if (defect > 0.0)
246 BlockRef(d[(*itActiveInd)[0]], (*itActiveInd)[1]) = 0.0;
247 }
248}
249
250template <typename TDomain, typename TAlgebra>
251void
255
256} // end namespace ug
257
258#endif /* __H__UG__LIB_ALGEBRA__OPERATOR__PRECONDITIONER__PROJECTED_GAUSS_SEIDEL__OBSTACLE_IN_NORMAL_DIR_IMPL__ */
259
Definition fe_geom.h:183
dimension dependent base class for reference elements
Definition reference_element.h:183
A class for fixed size, dense matrices.
Definition math_matrix.h:63
a mathematical Vector with N entries.
Definition math_vector.h:97
Definition multi_index.h:50
void adjust_sol_and_cor(value_type &sol_i, value_type &c_i, bool &dofIsActive, const DoFIndex &dof)
projects the i-th index of the solution onto the admissible set and adjusts the correction
Definition obstacle_in_normal_dir_impl.h:199
void transform_eulerian_coord_sys(MathVector< dim > transformedONB[], const MathVector< dim > &firstTransformedBaseVec)
Definition obstacle_in_normal_dir_impl.h:61
void adjust_defect_to_constraint(vector_type &d)
the defect needs to be adjusted for the active indices (those indices, which are in contact)
Definition obstacle_in_normal_dir_impl.h:237
algebra_type::vector_type vector_type
Vector type.
Definition obstacle_in_normal_dir.h:74
vector_type::value_type value_type
Value type.
Definition obstacle_in_normal_dir.h:77
void preprocess()
preprocess is useful to attach the normals for every obstacle DoF
Definition obstacle_in_normal_dir_impl.h:46
void restrict_obs_values()
restricts the obstacle values to a coarser grid in a multigrid hierarchy
Definition obstacle_in_normal_dir_impl.h:253
void adjust_sol_and_cor_elem(TIterator iterBegin, TIterator iterEnd, value_type &sol_i, value_type &c_i, bool &dofIsAdmissible, const DoFIndex &dof)
Definition obstacle_in_normal_dir_impl.h:126
Definition reference_element.h:345
static const DimReferenceElement< dim > & get(ReferenceObjectID roid)
returns a dimension dependent Reference Element
Definition reference_element.h:280
Definition reference_element.h:445
Definition reference_element.h:394
Definition reference_element.h:303
void CollectCornerCoordinates(std::vector< typename TAAPos::ValueType > &vCornerCoordsOut, const TElem &elem, const TAAPos &aaPos, bool clearContainer=true)
returns the corner coordinates of a geometric object
Definition domain_util_impl.h:75
void MatHouseholder(matrix_t &mOut, const vector_t &orthoVec)
Definition math_matrix_functions_common_impl.hpp:942
#define UG_THROW(msg)
Definition error.h:57
#define UG_LOG(msg)
Definition log.h:367
double number
Definition types.h:124
vector_t::value_type VecLength(const vector_t &v)
returns the length of v. Slower than VecLengthSq.
Definition math_vector_functions_common_impl.hpp:341
void VecScaleAdd(vector_t &vOut, typename vector_t::value_type s1, const vector_t &v1, typename vector_t::value_type s2, const vector_t &v2)
Scales two Vectors, adds them and returns the sum in a third vector.
Definition math_vector_functions_common_impl.hpp:265
void VecScale(vector_t &vOut, const vector_t &v, typename vector_t::value_type s)
scales a MathVector<N>
Definition math_vector_functions_common_impl.hpp:252
the ug namespace
ReferenceObjectID
these ids are used to identify the shape of a geometric object.
Definition grid_base_objects.h:74
double & BlockRef(T &vec, size_t i)
Definition blocks.h:66
int DimensionOfSubset(const ISubsetHandler &sh, int si)
returns the current dimension of the subset
Definition subset_dim_util.cpp:89
ReferenceVertex face_type1
Definition obstacle_in_normal_dir_impl.h:104
ReferenceVertex face_type0
Definition obstacle_in_normal_dir_impl.h:103
DimFEGeometry< 1, 1 > DimFEGeo
Definition obstacle_in_normal_dir_impl.h:105
ReferenceEdge face_type0
Definition obstacle_in_normal_dir_impl.h:110
ReferenceEdge face_type1
Definition obstacle_in_normal_dir_impl.h:111
DimFEGeometry< 2, 1 > DimFEGeo
Definition obstacle_in_normal_dir_impl.h:112
ReferenceTriangle face_type0
Definition obstacle_in_normal_dir_impl.h:117
ReferenceQuadrilateral face_type1
Definition obstacle_in_normal_dir_impl.h:118
DimFEGeometry< 3, 2 > DimFEGeo
Definition obstacle_in_normal_dir_impl.h:119
Definition active_set_impl.h:45
void face_type0
Definition obstacle_in_normal_dir_impl.h:96
void face_type1
Definition obstacle_in_normal_dir_impl.h:97
void DimFEGeo
Definition obstacle_in_normal_dir_impl.h:98