Loading [MathJax]/extensions/tex2jax.js
Plugins
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
contact_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 CONTACT_IMPL_H_
34#define CONTACT_IMPL_H_
35
36// other ug4 modules
37#include "common/common.h"
39
40#include "contact.h"
41
42namespace ug{
43namespace SmallStrainMechanics{
44
45template <int dim> struct face_type_traits
46{
47 typedef void face_type0;
48 typedef void face_type1;
49 typedef void DimFEGeo;
50};
51
58
65
72
73template <typename TDomain>
74template <typename TElem>
76 LocalVector& locForce, GridObject* side,
77 TElem* elem, const MathVector<dim> sideCoPos[], int numElemCorners,
78 const LocalVector& locU, vector<DoFIndex> vActiveSetLoc)
79{
80 typedef typename face_type_traits<dim>::face_type0 face_type0;
81 typedef typename face_type_traits<dim>::face_type1 face_type1;
82 typedef typename face_type_traits<dim>::DimFEGeo sideGeo;
83
84 ReferenceObjectID sideRoid = side->reference_object_id();
85 sideGeo geo(sideRoid, 3, LFEID(LFEID::LAGRANGE, dim, 1));
86
87 // prepare geometry for type and order
88 try{
89 geo.update_local(sideRoid, LFEID(LFEID::LAGRANGE, dim, 1), 3);
90 }
91 UG_CATCH_THROW("SmallStrainMechanicsElemDisc::contact_forces_elem_ips_avg:"
92 " Cannot update local values of finite element geometry.");
93
94 try{
95 geo.update(side, sideCoPos, LFEID(LFEID::LAGRANGE, dim, 1), 3);
96 }
97 UG_CATCH_THROW("SmallStrainMechanicsElemDisc::contact_forces_elem_ips_avg:"
98 " Cannot update finite element geometry.");
99
100 MathVector<dim> normal, normalizedNormal;
101 if (numElemCorners == dim)
102 ElementNormal<face_type0, dim>(normal, sideCoPos);
103 else
104 ElementNormal<face_type1, dim>(normal, sideCoPos);
105
106 number normOfNormal = VecLength(normal);
107 VecScale(normalizedNormal, normal, 1.0/normOfNormal);
108
109 MathMatrix<dim, dim> GradU, sigma;
110 for(size_t ip = 0; ip < geo.num_ip(); ++ip) // loop ip
111 {
112 // get displacementGradient (GradU), linearized strain tensor (eps)
113 // and Cauchy stress tensor sigma at ip
114
115 m_spMatLaw->template DisplacementGradient<sideGeo>(GradU, ip, geo, locU);
116 m_spMatLaw->stressTensor(sigma, ip, GradU);
117
118 /*for (size_t i = 0; i < (size_t) dim; ++i){
119 for (size_t J = 0; J < (size_t) dim; ++J)
120 {
121 GradU[i][J] = 0.0;
122 for (size_t a = 0; a < geo.num_sh(); ++a) // loop shapes
123 GradU[i][J] += geo.global_grad(ip, a)[J] * locU(i, a);
124 }
125 }
126
127 for (size_t i = 0; i < (size_t) dim; ++i)
128 for (size_t J = 0; J < (size_t) dim; ++J)
129 eps[i][J] = 0.5 * (GradU[i][J] + GradU[J][i]);
130
131 TensContract4(sigma, m_ElastTensorFunct, eps);*/
132
133 number sigma_n = 0.0;
134 for (size_t i = 0; i < (size_t) dim; ++i)
135 for (size_t j = 0; j < (size_t) dim; ++j)
136 sigma_n += sigma[i][j] * normalizedNormal[i] * normalizedNormal[j];
137
138 // compute normal component of sigma; sigma_n = sigma_ij * normal_i * normal_j
139 // 'vActiveSetLoc' is a vector, which stores all active local (dof,fct)-pairs
140 // of the SideElem 'side'
141 for (vector<DoFIndex>::iterator itActiveInd = vActiveSetLoc.begin();
142 itActiveInd < vActiveSetLoc.end(); ++itActiveInd)
143 locForce((*itActiveInd)[1], (*itActiveInd)[0]) += sigma_n * geo.weight(ip);
144
145 } // end(ip)
146
147 // FOR 3D-FULL-HEXAEDER-ELEMENTS instead of sideElems it is possible
148 // to consider plasticity by means of plastic ip's via
149 // 'm_pElemData = &m_aaElemData[static_cast<TElem*>(elem)];'
150}
151
152template <typename TDomain>
153template <typename TSide, typename TElem>
155 LocalVector& locForce, TSide* side,
156 TElem* elem, const MathVector<dim> sideCoPos[],
157 const LocalVector& locU, vector<DoFIndex> vActiveSetLoc)
158{
159 // TODO: inclusion of plastic variables, which are defined at the elements ip`s
160 typedef typename face_type_traits<dim>::face_type0 face_type0;
161 typedef typename face_type_traits<dim>::face_type1 face_type1;
162 SmartPtr<TDomain> dom = this->domain();
163 typename domain_type::position_accessor_type& aaPos = dom->position_accessor();
164
165 // reference object type
166 ReferenceObjectID sideRoid = side->reference_object_id();
167
168// some storage
169 MathMatrix<dim, dim-1> JTInv;
170 vector<MathVector<dim-1> > vLocalGrad;
171 vector<MathVector<dim> > vGlobalGrad;
172 vector<MathVector<dim> > vCorner;
173
174// get trial space
175 const LocalShapeFunctionSet<dim-1>& lsfs =
176 LocalFiniteElementProvider::get<dim-1>(sideRoid, LFEID(LFEID::LAGRANGE, dim, 1));
177
178// create a reference mapping
179 DimReferenceMapping<dim-1, dim>& refMap
180 = ReferenceMappingProvider::get<dim-1, dim>(sideRoid);
181
182// number of shape functions
183 const size_t numSH = lsfs.num_sh();
184 vLocalGrad.resize(numSH);
185 vGlobalGrad.resize(numSH);
186
187// get local Mid Point (for ROID_EDGE, ROID_QUADRILATERAL, ROID_TRIANGLE)
188 if (dim-1 < 1 || dim-1 > 2)
189 UG_THROW("SmallStrainMechanicsElemDisc::contact_forces_elem_midpoint:"
190 "side dimension " << dim-1 << " is not supported!");
191
192 MathVector<dim-1> localIP;
193 for (int i = 0; i < dim-1; ++i)
194 localIP[i] = 0.5;
195
196 if (sideRoid == ROID_TRIANGLE){
197 localIP[0] = 1.0/3.0; localIP[1] = 1.0/3.0;
198 }
199
200// evaluate reference gradient at local midpoint
201 lsfs.grads(&vLocalGrad[0], localIP);
202
203// get corners of element
204 CollectCornerCoordinates(vCorner, *side, aaPos);
205
206// update mapping
207 refMap.update(&vCorner[0]);
208
209// compute jacobian
210 refMap.jacobian_transposed_inverse(JTInv, localIP);
211
212// call 'ElementNormal' depending on elem-type
213// faces have dim corners in 1d, 2d
214// in 3d they have dim corners (triangle) or dim+1 corners (quadrilateral)
215 MathVector<dim> normal, normalizedNormal;
216
217 int numElemCorners = (int)vCorner.size();
218 if (numElemCorners == dim)
219 ElementNormal<face_type0, dim>(normal, sideCoPos);
220 else
221 ElementNormal<face_type1, dim>(normal, sideCoPos);
222
223 number normOfNormal = VecLength(normal);
224 VecScale(normalizedNormal, normal, 1.0/normOfNormal);
225
226// compute gradient at mid point by summing contributions of all shape fct
227 MathMatrix<dim, dim> GradU, eps, sigma;
228 GradU = 0.0;
229
230 sigma = 0.0;
231 /*for(size_t sh = 0; sh < numSH; ++sh)
232 {
233 // get global Gradient
234 MatVecMult(vGlobalGrad[sh], JTInv, vLocalGrad[sh]);
235 for(size_t fct = 0; fct < locU.num_fct(); ++fct)
236 for(size_t i = 0; i < (size_t) dim; ++i)
237 GradU[fct][i] += vGlobalGrad[sh][i] * locU(fct, sh);
238 }
239
240 // get strain tensor 'eps' at local midPoint-IP
241 for(size_t fct = 0; fct < locU.num_fct(); ++fct)
242 for (size_t i = 0; i < (size_t) dim; ++i)
243 eps[fct][i] = 0.5 * (GradU[fct][i] + GradU[i][fct]);
244
245 // get Cauchy stress tensor 'sigma'
246 TensContract4(sigma, m_ElastTensorFunct, eps);*/
247
248 // compute normal component -sigma_n
249 number sigma_n = 0.0;
250 for (size_t i = 0; i < (size_t) dim; ++i)
251 for (size_t j = 0; j < (size_t) dim; ++j)
252 sigma_n += sigma[i][j] * normalizedNormal[i] * normalizedNormal[j];
253
254/* UG_LOG("sigma_n: " << sigma_n << "\n");
255 MathVector<dim> sigma_n_normal;
256 VecScale(sigma_n_normal, normalizedNormal, sigma_n);
257 UG_LOG("sigma_n_normal: " << sigma_n_normal << "\n");*/
258
259 // add contributions to contactForces for active indices
260 for (vector<DoFIndex>::iterator itAcvtiveInd = vActiveSetLoc.begin();
261 itAcvtiveInd < vActiveSetLoc.end(); ++itAcvtiveInd)
262 locForce((*itAcvtiveInd)[1], (*itAcvtiveInd)[0]) = sigma_n; //sigma_n_normal[activeLocFct];
263
264}
265
266template <typename TDomain, typename TGridFunction>
272
273template <typename TDomain, typename TGridFunction>
274template <typename TSideElem, typename TIterator>
276 TIterator iterBegin,
277 TIterator iterEnd,
278 const TGridFunction& u,
279 TGridFunction& contactForce,
280 vector<DoFIndex> vActiveSetGlob)
281{
282// check if at least an element exists, else return
283 if(iterBegin == iterEnd) return;
284
285 typedef typename TGridFunction::domain_type domain_type;
286 typedef typename domain_type::grid_type grid_type;
287 typedef typename TGridFunction::element_type element_type;
288 typename grid_type::template traits<element_type>::secure_container associatedElems;
289 typename domain_type::position_accessor_type& aaPos
290 = contactForce.domain()->position_accessor();
291
292 grid_type& grid = *contactForce.domain()->grid();
293 vector<DoFIndex> vActiveSetLoc;
294
295 // local indices and local algebra
296 LocalIndices ind; LocalVector locU, locContactForce;
297
298// loop over all elements on activeSubset si
299 for(TIterator iter = iterBegin; iter != iterEnd; ++iter)
300 {
301 TSideElem* sideElem = *iter;
302
303 // get all associated elements to 'elem' of activeSubset
304 grid.associated_elements(associatedElems, sideElem);
305
306 if (associatedElems.size() > 1)
307 UG_THROW("Temporarily it is only possible to compute contact Forces "
308 "for an element of dimension " << dim-1 << " with exactly "
309 "one associated element of dimension " << dim << " in "
310 "ContactSmallStrainMechanics::contactForcesElem \n");
311
312 // get global indices
313 u.indices(*iter, ind);
314
315 // adapt local algebra
316 locU.resize(ind); locContactForce.resize(ind);
317
318 vActiveSetLoc.resize(0);
319
320 // create local activeSet-vector out of global activeSet-vector
321 for (vector<DoFIndex>::iterator itActiveInd = vActiveSetGlob.begin();
322 itActiveInd < vActiveSetGlob.end(); ++itActiveInd)
323 {
324 for(size_t fct=0; fct < locU.num_all_fct(); ++fct){
325 for(size_t dof=0; dof < locU.num_all_dof(fct); ++dof)
326 {
327 size_t globIndex = ind.index(fct, dof);
328 size_t globComp = ind.comp(fct, dof);
329
330 // check if local (fct,dof) corresponds to a global activeSetMultiIndex
331 if ((*itActiveInd)[0] == globIndex &&
332 (*itActiveInd)[1] == globComp)
333 {
334 UG_LOG("active multiIndex (ind,comp) in "
335 "ContactSmallStrainMechanics::contactForces: ("
336 << globIndex << "," << globComp << ") \n");
337
338 // store local (dof,fct)-pairs which are active
339 // in the current element
340 DoFIndex activeMultiIndexLoc(dof, fct);
341
342 UG_LOG("active locU (dof,fct) in "
343 "ContactSmallStrainMechanics::contactForces: ("
344 << dof << "," << fct << ") \n");
345
346 // create list of active MultiIndex-pairs (local)
347 vActiveSetLoc.push_back(activeMultiIndexLoc);
348 }
349 }
350 }
351 } // end (vActiveSet-loop)
352
353 if (vActiveSetLoc.size() != 0.0)
354 {
355 UG_LOG("activeDoFs in Elem" << *iter << "in "
356 "ContactSmallStrainMechanics::contactForces: "
357 << vActiveSetLoc.size() << "\n");
358
359 // reset contribution of this element
360 locContactForce = 0.0;
361
362 // local vector extract -> locU
363 GetLocalVector(locU, u);
364
365 // storage for corner coordinates
366 vector<MathVector<dim> > vCorner;
367 MathVector<dim> sideCoPos[dim+1];
368
369 // reference object type
370 ReferenceObjectID roid = sideElem->reference_object_id();
371
372 const DimReferenceElement<dim-1>& rRefElem
373 = ReferenceElementProvider::get<dim-1>(roid);
374
375 // get corner coordinates
376 CollectCornerCoordinates(vCorner, *sideElem, aaPos);
377
378 // here the ordering of the corners in the reference element is exploited
379 // in order to compute the outer normal later on
380 int numElemCorners = (int)vCorner.size();
381
382 for (int i = 0; i < numElemCorners; ++i)
383 sideCoPos[i] = vCorner[rRefElem.id(dim-1, 0, 0, i)];
384
385 // pass the local indices (dof,fct) which are active
386 // to elementwise computation of contactForces
387 //m_spMechElemDisc->contact_forces_elem_ips_avg(locForce, elem,
388 // assoElement[0], sideCoPos, numElemCorners, locU, vActiveSetLoc);
389
390 m_spMechElemDisc->contact_forces_elem_midpoint(locContactForce, sideElem,
391 associatedElems[0], sideCoPos, locU, vActiveSetLoc);
392
393 /*for(size_t fct=0; fct < locForce.num_all_fct(); ++fct){
394 for(size_t dof=0; dof < locForce.num_all_dof(fct); ++dof)
395 {
396 UG_LOG("locForce(fct " << fct << ", dof " << dof <<
397 "): " << locForce.value(fct,dof) << "\n");
398 //const size_t globIndex = indF.index(fct, dof);
399 //const size_t globComp = indF.comp(fct, dof);
400 //BlockRef(force[globIndex], globComp) = locForce.value(fct,dof);
401 }
402 }*/
403
404 // send local to global contact force
405 AddLocalVector(contactForce, locContactForce);
406 // TODO: not add forces, set forces?!
407 }
408 } // end (elements)
409}
410
411template <typename TDomain, typename TGridFunction>
413 TGridFunction& contactForce, const TGridFunction& u,
414 vector<DoFIndex> vActiveSetGlob,
415 vector<int> vActiveSubsets)
416{
417 if (m_spMechElemDisc.invalid())
418 UG_THROW("No element discretization set in "
419 "ContactSmallStrainMechanics:lagrange_multiplier \n");
420
421 UG_LOG("Active subsets: " << vActiveSubsets.size() << "\n");
422
423 SmartPtr<DoFDistribution> dd = contactForce.dof_distribution();
424
425 for (vector<int>::iterator siContact = vActiveSubsets.begin();
426 siContact != vActiveSubsets.end(); ++siContact)
427 {
428 UG_LOG("siContact: " << *siContact << "\n");
429 const int subsetDim = DimensionOfSubset(*dd->subset_handler(), *siContact);
430
431 switch(subsetDim)
432 {
433 case 1:
434 contact_forces_elem<RegularEdge>
435 (dd->template begin<RegularEdge>(*siContact), dd->template end<RegularEdge>(*siContact),
436 u, contactForce, vActiveSetGlob);
437 break;
438 case 2:
439 contact_forces_elem<Triangle>
440 (dd->template begin<Triangle>(*siContact), dd->template end<Triangle>(*siContact),
441 u, contactForce, vActiveSetGlob);
442 contact_forces_elem<Quadrilateral>
443 (dd->template begin<Quadrilateral>(*siContact), dd->template end<Quadrilateral>(*siContact),
444 u, contactForce, vActiveSetGlob);
445 break;
446 default:
447 UG_THROW("ContactSmallStrainMechanics::lagrange_multiplier:"
448 "SubsetDimension "<< subsetDim <<" (subset="<< *siContact <<") not supported.");
449 }
450 }
451
452}
453
454}// end of namespace SmallStrainMechanics
455}// end of namespace ug
456
457
458#endif /* CONTACT_IMPL_H_ */
virtual ReferenceObjectID reference_object_id() const=0
index_type & index(size_t fct, size_t dof)
comp_type & comp(size_t fct, size_t dof)
size_t num_all_fct() const
void resize(const LocalIndices &ind)
size_t num_all_dof(size_t fct) const
std::size_t size() const
static const DimReferenceElement< dim > & get(ReferenceObjectID roid)
static DimReferenceMapping< TDim, TWorldDim > & get(ReferenceObjectID roid)
virtual void lagrange_multiplier(TGridFunction &lagMult, const TGridFunction &u, vector< DoFIndex > vActiveSet, vector< int > vActiveSubsets)
Definition contact_impl.h:412
base_type::domain_type domain_type
Domain type.
Definition contact.h:77
void contact_forces_elem(TIterator iterBegin, TIterator iterEnd, const TGridFunction &u, TGridFunction &contactForce, vector< DoFIndex > vActiveSetGlob)
Definition contact_impl.h:275
ContactSmallStrainMechanics(SmartPtr< SmallStrainMechanicsElemDisc< TDomain > > spLinElastPlast)
Definition contact_impl.h:267
void contact_forces_elem_ips_avg(LocalVector &locForce, GridObject *side, TElem *elem, const MathVector< dim > sideCoPos[], int numElemCorners, const LocalVector &locU, std::vector< DoFIndex > vActiveSetLoc)
computing contact forces elementwise by averaging over all integration points
void contact_forces_elem_midpoint(LocalVector &locForce, TSide *side, TElem *elem, const MathVector< dim > sideCoPos[], const LocalVector &locU, std::vector< DoFIndex > vActiveSetLoc)
computing contact forces elementwise at every element midpoint
SmartPtr< TGrid > grid()
int element_type() const
void CollectCornerCoordinates(int base_object_id, std::vector< typename TDomain::position_type > &vCornerCoordsOut, GridObject &elem, const TDomain &domain, bool clearContainer)
#define UG_CATCH_THROW(msg)
#define UG_THROW(msg)
#define UG_LOG(msg)
double number
vector_t::value_type VecLength(const vector_t &v)
void VecScale(vector_t &vOut, const vector_t &v, typename vector_t::value_type s)
ReferenceObjectID
ROID_TRIANGLE
void GetLocalVector(LocalVector &lvec, const TVector &vec)
int DimensionOfSubset(const ISubsetHandler &sh, int si)
void AddLocalVector(TVector &vec, const LocalVector &lvec)
ReferenceVertex face_type0
Definition contact_impl.h:54
DimFEGeometry< 1, 1 > DimFEGeo
Definition contact_impl.h:56
ReferenceVertex face_type1
Definition contact_impl.h:55
DimFEGeometry< 2, 1 > DimFEGeo
Definition contact_impl.h:63
ReferenceEdge face_type0
Definition contact_impl.h:61
ReferenceEdge face_type1
Definition contact_impl.h:62
ReferenceQuadrilateral face_type1
Definition contact_impl.h:69
ReferenceTriangle face_type0
Definition contact_impl.h:68
DimFEGeometry< 3, 2 > DimFEGeo
Definition contact_impl.h:70
Definition contact_impl.h:46
void DimFEGeo
Definition contact_impl.h:49
void face_type1
Definition contact_impl.h:48
void face_type0
Definition contact_impl.h:47