Loading [MathJax]/extensions/tex2jax.js
Plugins
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
mech_output_writer_impl.h
Go to the documentation of this file.
1/*
2 * Copyright (c) 2014-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 MECH_OUTPUT_WRITER_IMPL_H_
34#define MECH_OUTPUT_WRITER_IMPL_H_
35
37
38// module intern headers
39#include "mech_output_writer.h"
40
41namespace ug {
42namespace SmallStrainMechanics{
43
44template <typename TGridFunction>
45void plast_ip(TGridFunction& plastIP,
47 TGridFunction& u)
48{
49 static const int dim = TGridFunction::dim;
50 typedef typename TGridFunction::template dim_traits<dim>::grid_base_object grid_base_object;
51 typedef typename TGridFunction::template dim_traits<dim>::const_iterator const_iterator;
52
53 // local indices and local algebra
54 LocalIndices ind; LocalVector locU, locPlastIP;
55
56 const_iterator iter = u.template begin<grid_base_object>();
57 const_iterator end = u.template end<grid_base_object>();
58
59 // loop over all elements
60 for(;iter != end; ++iter)
61 {
62 // get element
63 grid_base_object* elem = *iter;
64
65 // get global indices
66 u.indices(elem, ind);
67
68 // adapt local algebra
69 locU.resize(ind); locPlastIP.resize(ind);
70
71 // reset contribution of this element
72 locPlastIP = 0.0;
73
74 // local vector extract -> locU
75 GetLocalVector(locU, u); GetLocalVector(locPlastIP, plastIP);
76
77 mechOut.plastIP_elem(locPlastIP, elem, locU,
78 u.domain());
79
80 // send local to global plastIP
81 AddLocalVector(plastIP, locPlastIP);
82 }
83}
84
85template <typename TGridFunction>
86void equiv_plast_strain(TGridFunction& eqPlastStrain,
88 TGridFunction& u)
89{
90 static const int dim = TGridFunction::dim;
91 typedef typename TGridFunction::template dim_traits<dim>::grid_base_object grid_base_object;
92 typedef typename TGridFunction::template dim_traits<dim>::const_iterator const_iterator;
93
94 // local indices and local algebra
95 LocalIndices ind; LocalVector locU, locEqPlastStrain;
96
97 const_iterator iter = u.template begin<grid_base_object>();
98 const_iterator end = u.template end<grid_base_object>();
99
100 // loop over all elements
101 for(;iter != end; ++iter)
102 {
103 // get element
104 grid_base_object* elem = *iter;
105
106 // get global indices
107 u.indices(elem, ind);
108
109 // adapt local algebra
110 locU.resize(ind); locEqPlastStrain.resize(ind);
111
112 // reset contribution of this element
113 locEqPlastStrain = 0.0;
114
115 // local vector extract -> locU
116 GetLocalVector(locU, u); GetLocalVector(locEqPlastStrain, eqPlastStrain);
117
118 mechOut.eqPlastStrain_elem(locEqPlastStrain, elem, locU, u.domain());
119
120 // send local to global eqPlastStrain
121 AddLocalVector(eqPlastStrain, locEqPlastStrain);
122 }
123}
124
125template<typename TGridFunction>
126void
128 TGridFunction& sigma, TGridFunction& epsilon,
129 TGridFunction& stressFunc, TGridFunction& strainFunc, TGridFunction& u)
130{
131 static const int dim = TGridFunction::dim;
132 typedef typename TGridFunction::template dim_traits<dim>::grid_base_object grid_base_object;
133 typedef typename TGridFunction::template dim_traits<dim>::const_iterator const_iterator;
134
135 // local indices and local algebra
136 LocalIndices indU, indEps, indSig, indStressFunc, indStrainFunc;
137 LocalVector locU, locSig, locEps, locStressFunc, locStrainFunc;
138
139 const_iterator iter = u.template begin<grid_base_object>();
140 const_iterator end = u.template end<grid_base_object>();
141
142 // loop over all elements
143 for(;iter != end; ++iter)
144 {
145 // get element
146 grid_base_object* elem = *iter;
147
148 // get global indices
149 u.indices(elem, indU); epsilon.indices(elem, indEps);
150 sigma.indices(elem, indSig); stressFunc.indices(elem, indStressFunc);
151 strainFunc.indices(elem, indStrainFunc);
152
153 // adapt local algebra
154 locU.resize(indU); locEps.resize(indEps);
155 locSig.resize(indSig); locStressFunc.resize(indStressFunc);
156 locStrainFunc.resize(indStrainFunc);
157
158 // reset contribution of this element
159 locStressFunc = 0.0; locStrainFunc = 0.0; locSig = 0.0; locEps = 0.0;
160
161 // local vector extract -> locU
162 GetLocalVector(locU, u);
163
164 // call method of solid mechanics output writer
165 mechOut.normal_stress_strain_elem(locSig, locEps, locStressFunc, locStrainFunc, elem,
166 locU, u.domain());
167
168 // send local to global
169 // sigma (cauchy-stress tensor), epsilon (linearized strain tensor), devSigma (deviatoric part of Sigma)
170 AddLocalVector(stressFunc, locStressFunc);
171 AddLocalVector(strainFunc, locStrainFunc);
172 AddLocalVector(sigma, locSig);
173 AddLocalVector(epsilon, locEps);
174 }
175}
176
177template <typename TGridFunction>
178void invariants_kirchhoff_stress(TGridFunction& invarKirchhoffStress,
180 TGridFunction& u)
181{
182 static const int dim = TGridFunction::dim;
183 typedef typename TGridFunction::template dim_traits<dim>::grid_base_object grid_base_object;
184 typedef typename TGridFunction::template dim_traits<dim>::const_iterator const_iterator;
185
186 // local indices and local algebra
187 LocalIndices ind, indInvarKirchhoffStress;
188 LocalVector locU, locInvarKirchhoffStress;
189
190 const_iterator iter = u.template begin<grid_base_object>();
191 const_iterator end = u.template end<grid_base_object>();
192
193 // loop over all elements
194 for(;iter != end; ++iter)
195 {
196 // get element
197 grid_base_object* elem = *iter;
198
199 // get global indices
200 u.indices(elem, ind); invarKirchhoffStress.indices(elem, indInvarKirchhoffStress);
201
202 // adapt local algebra
203 locU.resize(ind); locInvarKirchhoffStress.resize(indInvarKirchhoffStress);
204
205 // local vector extract -> locU
206 GetLocalVector(locU, u); GetLocalVector(locInvarKirchhoffStress, invarKirchhoffStress);
207
208 // reset contribution of this element
209 locInvarKirchhoffStress = 0.0;
210
211 mechOut.invariants_kirchhoff_stress_elem(locInvarKirchhoffStress, elem, locU,
212 u.domain());
213
214 // send local to global tensor
215 AddLocalVector(invarKirchhoffStress, locInvarKirchhoffStress);
216 }
217}
218
220template<typename TDomain>
221void
223plastIP_elem(LocalVector& locPlastIP, TBaseElem* elem,
224 const LocalVector& locU, SmartPtr<TDomain> dom)
225{
226 // get vertices and extract corner coordinates
227 typedef typename TDomain::position_accessor_type position_accessor_type;
228 position_accessor_type& aaPos = dom->position_accessor();
230 const size_t numVertices = elem->num_vertices();
231 for (size_t i = 0; i < numVertices; ++i) {
232 coCoord[i] = aaPos[elem->vertex(i)];
233 };
234
235 // update the geometry for the Element elem
237 try{
238 geo.update(elem, &(coCoord[0]),
239 LFEID(LFEID::LAGRANGE, dim, 1), m_quadOrder);
240 }
241 UG_CATCH_THROW("SmallStrainMechOutput::plastIP_elem:"
242 " Cannot update Finite Element Geometry.");
243
244 // get all neighbor elems which share a vertex with the given element 'elem'
245 typename TDomain::grid_type& grid = *(dom->grid());
246 typedef typename vector<TBaseElem*>::iterator neighborElemIter;
247 vector<TBaseElem*> vNeighborElems;
248 CollectNeighbors(vNeighborElems, elem, grid, NHT_VERTEX_NEIGHBORS);
249
250 // pointer to internal variable of current elem
251 m_spMatLaw->internal_vars(elem);
252
253 /*MathMatrix<dim, dim> GradU;
254 number plasticIP = 0.0;
255 for (size_t ip = 0; ip < geo.num_ip(); ++ip)
256 {
257 // get displacementGradient (GradU)
258 m_spMatLaw->template DisplacementGradient<DimFEGeometry<dim> >(GradU, ip, geo, locU);
259
260 number gamma = m_spMatLaw->plastic_multiplier(ip, GradU);
261 if (gamma > 0.0)
262 plasticIP += 1.0;
263 if (gamma < 0.0)
264 UG_THROW("gamma: " << gamma << "in plastIP_elem \n");
265 }*/
266
267 MathMatrix<dim, dim> plastStrain;
268 number plasticIP = 0.0;
269 for (size_t ip = 0; ip < geo.num_ip(); ++ip)
270 {
271 plastStrain = *(m_spMatLaw->inelastic_strain_tensor(ip));
272
273 number plastStrainNorm = MatFrobeniusNorm(plastStrain);
274 if (plastStrainNorm > 0.0)
275 plasticIP += 1.0;
276 }
277
278 for (size_t co = 0; co < elem->num_vertices(); ++co)
279 {
280 // init 'elemsWithCo' with 1, because co lies in 'elem'
281 size_t elemsWithCo = 1;
282
283 // iterate neighbor elems and count how many elems contain the corner 'co'
284 for (neighborElemIter it = vNeighborElems.begin();
285 it != vNeighborElems.end(); ++it){
286 if (ContainsPoint(*it, coCoord[co], aaPos) == true)
287 elemsWithCo++;
288 }
289
290 // scaling factor for averaging values out of all ips
291 // of the associated elements of corner co
292 size_t scaleFac = elemsWithCo * geo.num_ip();
293
294 locPlastIP(0,co) = plasticIP / scaleFac;
295 }
296}
297
299template<typename TDomain>
300void
302eqPlastStrain_elem(LocalVector& locEqPlastStrain, TBaseElem* elem,
303 const LocalVector& locU, SmartPtr<TDomain> dom)
304{
305 // get vertices and extract corner coordinates
306 typedef typename TDomain::position_accessor_type position_accessor_type;
307 position_accessor_type& aaPos = dom->position_accessor();
309 const size_t numVertices = elem->num_vertices();
310 for (size_t i = 0; i < numVertices; ++i) {
311 coCoord[i] = aaPos[elem->vertex(i)];
312 };
313
314 // prepare geometry for type and order
316 try{
317 geo.update(elem, &(coCoord[0]),
318 LFEID(LFEID::LAGRANGE, dim, 1), m_quadOrder);
319 }
320 UG_CATCH_THROW("SmallStrainMechOutput::eqPlastStrain_elem:"
321 " Cannot update Finite Element Geometry.");
322
323 // get all neighbor elems which share a vertex with the given element 'elem'
324 typename TDomain::grid_type& grid = *(dom->grid());
325 typedef typename vector<TBaseElem*>::iterator neighborElemIter;
326 vector<TBaseElem*> vNeighborElems;
327 CollectNeighbors(vNeighborElems, elem, grid, NHT_VERTEX_NEIGHBORS);
328
329 // pointer to internal variable of current elem
330 m_spMatLaw->internal_vars(elem);
331
333 number SumEqPlastStrainIP = 0.0;
334 for (size_t ip = 0; ip < geo.num_ip(); ++ip)
335 {
336 // get displacementGradient (GradU)
337 m_spMatLaw->template DisplacementGradient<DimFEGeometry<dim> >(GradU, ip, geo, locU);
338
339 number gamma = m_spMatLaw->plastic_multiplier(ip, GradU);
340
341 // sum over all ip`s
342 SumEqPlastStrainIP += gamma * sqrt(2.0/3.0);
343 }
344
345 for (size_t co = 0; co < elem->num_vertices(); ++co)
346 {
347 // init 'elemsWithCo' with 1, because co lies in 'elem'
348 size_t elemsWithCo = 1;
349
350 // iterate neighbor elems and count how many elems contain the corner 'co'
351 for (neighborElemIter it = vNeighborElems.begin();
352 it != vNeighborElems.end(); ++it){
353 if (ContainsPoint(*it, coCoord[co], aaPos) == true)
354 elemsWithCo++;
355 }
356
357 // scaling factor for averaging values out of all ips
358 // of the associated elements of corner co
359 size_t scaleFac = elemsWithCo * geo.num_ip();
360
361 locEqPlastStrain(0, co) = SumEqPlastStrainIP / scaleFac;
362 }
363}
364template<typename TDomain>
365void
367normal_stress_strain_elem(LocalVector& locSigma, LocalVector& locEps, LocalVector& locStressFunc,
368 LocalVector& locStrainFunc, TBaseElem* elem, const LocalVector& locU,
370{
371// get all neighbor elems which share a vertex with the given element 'elem'
372 typename TDomain::grid_type& grid = *(dom->grid());
373 vector<TBaseElem*> vNeighborElems;
374 CollectNeighbors(vNeighborElems, elem, grid, NHT_VERTEX_NEIGHBORS);
375 typedef typename vector<TBaseElem*>::iterator neighborElemIter;
376
377 typedef typename TDomain::position_accessor_type position_accessor_type;
378 position_accessor_type& aaPos = dom->position_accessor();
379
380 // coord and vertex array
382
383 // get vertices and extract corner coordinates
384 const size_t numVertices = elem->num_vertices();
385 for (size_t i = 0; i < numVertices; ++i) {
386 coCoord[i] = aaPos[elem->vertex(i)];
387 };
388
389 // prepare geometry for type and order
391 try{
392 geo.update(elem, &(coCoord[0]),
393 LFEID(LFEID::LAGRANGE, dim, 1), m_quadOrder);
394 }
395 UG_CATCH_THROW("SmallStrainMechOutput::normal_stress_strain_elem:"
396 " Cannot update Finite Element Geometry.");
397
398 // pointer to internal variable of current elem
399 m_spMatLaw->internal_vars(elem);
400
401 MathMatrix<dim, dim> GradU, eps, sigma, devSigma, SumEpsIP, SumSigmaIP, inelasticStrainTens, eps_elast;
402 number normDevSig, normSig, normEps, SumNormDevSigIP, SumNormSigIP, SumNormEpsIP,
403 SumVolSigmaIP, normElastStrain, SumNormElastStrainIP, normInelastStrain, SumNormInelastStrainIP;
404
405 // init summation-values
406 SumEpsIP = 0.0; SumSigmaIP = 0.0; SumNormDevSigIP = 0.0; SumNormSigIP = 0.0;
407 SumNormEpsIP = 0.0; SumVolSigmaIP = 0.0; SumNormElastStrainIP = 0.0; SumNormInelastStrainIP = 0.0;
408
409 for(size_t ip = 0; ip < geo.num_ip(); ++ip)
410 {
411 // get displacementGradient (GradU)
412 m_spMatLaw->template DisplacementGradient<DimFEGeometry<dim> >(GradU, ip, geo, locU);
413 // get Cauchy stress-tensor (sigma) at ip
414 m_spMatLaw->stressTensor(sigma, ip, GradU);
415
416 // get inelastic strain tensor at ip
417 inelasticStrainTens = *(m_spMatLaw->inelastic_strain_tensor(ip));
418
419 // get linearized strain tensor (eps) at ip
420 for (size_t i = 0; i < (size_t) dim; ++i)
421 for (size_t J = 0; J < (size_t) dim; ++J)
422 {
423 eps[i][J] = 0.5 * (GradU[i][J] + GradU[J][i]);
424 eps_elast[i][J] = eps[i][J] - inelasticStrainTens[i][J];
425 }
426
427 // add values to local vector
428 for (size_t i = 0; i < (size_t) dim; ++i)
429 {
430 SumEpsIP[i][i] += eps[i][i];
431 SumSigmaIP[i][i] += sigma[i][i];
432 }
433
434 // stress functions
435 number trSigma = MatDeviatorTrace(sigma, devSigma);
436 normDevSig = MatFrobeniusNorm(devSigma);
437 normSig = MatFrobeniusNorm(sigma);
438 SumVolSigmaIP += 1.0/3.0 * trSigma;
439 SumNormDevSigIP += normDevSig;
440 SumNormSigIP += normSig;
441
442 // strain functions
443 normEps = MatFrobeniusNorm(eps);
444 SumNormEpsIP += normEps;
445
446 normElastStrain = MatFrobeniusNorm(eps_elast);
447 SumNormElastStrainIP += normElastStrain;
448
449 normInelastStrain = MatFrobeniusNorm(inelasticStrainTens);
450 SumNormInelastStrainIP += normInelastStrain;
451
452 } //end (ip)
453
454 // loop corner
455 for (size_t co = 0; co < numVertices; ++co)
456 {
457 // TODO: is this geometry update necessary here?
458
459 // update geometry
460 try{
461 geo.update(elem, &(coCoord[0]),
462 LFEID(LFEID::LAGRANGE, dim, 1), m_quadOrder);
463 }
464 UG_CATCH_THROW("SmallStrainMechOutput::normal_stress_strain_loc:"
465 " Cannot update Finite Element Geometry.");
466
467 // init 'elemsWithCo' with 1, because co lies in 'elem'
468 size_t elemsWithCo = 1;
469
470 // iterate neighbor elems and count how many elems contain the corner 'co'
471 for (neighborElemIter it = vNeighborElems.begin();
472 it != vNeighborElems.end(); ++it){
473 if (ContainsPoint(*it, coCoord[co], aaPos) == true)
474 elemsWithCo++;
475 }
476
477 // scaling factor for averaging values out of all ips
478 // of the associated elements of corner co
479 size_t scaleFac = elemsWithCo * geo.num_ip();
480
481 // add values to local vector
482 for (size_t i = 0; i < (size_t) dim; ++i)
483 {
484 locEps(i, co) += SumEpsIP[i][i] / scaleFac;
485 locSigma(i, co) += SumSigmaIP[i][i] / scaleFac;
486 }
487
488 locStressFunc(0, co) += SumNormSigIP / scaleFac;
489 locStressFunc(1, co) += SumNormDevSigIP / scaleFac;
490 locStressFunc(2, co) += SumVolSigmaIP / scaleFac;
491
492 locStrainFunc(0, co) += SumNormEpsIP / scaleFac;
493 locStrainFunc(1, co) += SumNormElastStrainIP / scaleFac;
494 locStrainFunc(2, co) += SumNormInelastStrainIP / scaleFac;
495 } //end(co)
496
497
498 // temporarily output for linear-elastic-benchmark-output:
499 for (size_t co = 0; co < numVertices; ++co)
500 {
501 if ((coCoord[co][0] == 90.0) && (coCoord[co][1] == 0.0))
502 {
503 UG_LOG("At 2 (90.0,0.0): \n");
504 for (size_t i = 0; i < (size_t) dim; ++i)
505 UG_LOG("u("<< i << "," << co << "): " << locU(i,co) << "\n");
506
507 UG_LOG("\n");
508
509 vector<size_t> vNextIP;
510 next_ips_to_point(vNextIP, coCoord[co], geo);
511
513 for (vector<size_t>::iterator it = vNextIP.begin();
514 it != vNextIP.end(); ++it)
515 {
516 // get displacementGradient (GradU) at ip
517 m_spMatLaw->template DisplacementGradient<DimFEGeometry<dim> >(GradU, *it, geo, locU);
518 // get the Cauchy Stress Tensor (sigma) at ip
519 m_spMatLaw->stressTensor(sigma, *it, GradU);
520
521 UG_LOG("At " << geo.global_ip(*it) << ": \n");
522 UG_LOG("sigma_yy: " << sigma[1][1] << "\n");
523 UG_LOG("sigma: " << sigma << "\n");
524 }
525 }
526
527 if ((coCoord[co][0] == 100.0) && (coCoord[co][1] == 100.0))
528 {
529 UG_LOG("At 4 (100.0,100.0): \n");
530 for (size_t i = 0; i < (size_t) dim; ++i)
531 UG_LOG("u("<< i << "," << co << "): " << locU(i,co) << "\n");
532
533 UG_LOG("\n");
534 }
535
536 if ((coCoord[co][0] == 0.0) && (coCoord[co][1] == 100.0))
537 {
538 UG_LOG("At 5 (0.0,100.0): \n");
539 for (size_t i = 0; i < (size_t) dim; ++i)
540 UG_LOG("u("<< i << "," << co << "): " << locU(i,co) << "\n");
541
542 UG_LOG("\n");
543 }
544
545 };
546}
547
549template<typename TDomain>
550void
552invariants_kirchhoff_stress_elem(LocalVector& locInvarKirchhoffStress,
553 TBaseElem* elem, const LocalVector& u, SmartPtr<TDomain> dom)
554{
555 // get vertices and extract corner coordinates
556 typedef typename TDomain::position_accessor_type position_accessor_type;
557 position_accessor_type& aaPos = dom->position_accessor();
559 const size_t numVertices = elem->num_vertices();
560 for (size_t i = 0; i < numVertices; ++i) {
561 coCoord[i] = aaPos[elem->vertex(i)];
562 };
563
564 // prepare geometry for type and order
566 try{
567 geo.update(elem, &(coCoord[0]),
568 LFEID(LFEID::LAGRANGE, dim, 1), m_quadOrder);
569 }
570 UG_CATCH_THROW("MechOutputWriterFinite::invariants_kirchhoff_stress_elem:"
571 " Cannot update Finite Element Geometry.");
572
573 // get all neighbor elems which share a vertex with the given element 'elem'
574 typename TDomain::grid_type& grid = *(dom->grid());
575 typedef typename vector<TBaseElem*>::iterator neighborElemIter;
576 vector<TBaseElem*> vNeighborElems;
577 CollectNeighbors(vNeighborElems, elem, grid, NHT_VERTEX_NEIGHBORS);
578
579 // pointer to internal variable of current elem
580 m_spMatLaw->internal_vars(elem);
581
582 MathMatrix<dim, dim> tauTens, tauTens2, GradU, sigma, F, Ident;
583 MatIdentity(Ident);
584 for (size_t co = 0; co < elem->num_vertices(); ++co)
585 {
586 number invar1_ip = 0.0;
587 number invar2_ip = 0.0;
588 number invar3_ip = 0.0;
589
590 // loop ip`s and determine the next ip`s to corner co!
591 vector<size_t> vNextIP;
592 next_ips_to_point(vNextIP, coCoord[co], geo);
593
594 for (vector<size_t>::iterator it = vNextIP.begin();
595 it != vNextIP.end(); ++it)
596 {
597 // get deformation gradient at ip
598 m_spMatLaw->template DisplacementGradient<DimFEGeometry<dim> >(GradU, *it, geo, u);
599
600 // get Cauchy stress-tensor (sigma) at ip
601 m_spMatLaw->stressTensor(sigma, *it, GradU);
602
603 // compute the kirchhoff stress-tensor at ip
604 MatAdd(F, GradU, Ident);
605 number detF = Determinant(F);
606 MatScale(tauTens, 1.0/detF, sigma);
607
608
609 MatMultiply(tauTens2, tauTens, tauTens);
610 number traceTauTens_ip = Trace(tauTens);
611 number traceTauTens2_ip = Trace(tauTens2);
612
613 invar1_ip += traceTauTens_ip;
614 invar2_ip += 0.5 * (traceTauTens_ip * traceTauTens_ip - traceTauTens2_ip);
615 invar3_ip += Determinant(tauTens);
616 }
617
618 // init 'elemsWithCo' with 1, because co lies in 'elem'
619 size_t elemsWithCo = 1;
620
621 // iterate neighbor elems and count how many elems contain the corner 'co'
622 for (neighborElemIter it = vNeighborElems.begin();
623 it != vNeighborElems.end(); ++it){
624 if (ContainsPoint(*it, coCoord[co], aaPos) == true)
625 elemsWithCo++;
626 }
627
628 // scaling factor for averaging values out of all ips
629 // of the associated elements of corner co
630 const size_t scaleFac = elemsWithCo * vNextIP.size();
631
632 locInvarKirchhoffStress(0, co) = invar1_ip/scaleFac;
633 locInvarKirchhoffStress(1, co) = invar2_ip/scaleFac;
634 locInvarKirchhoffStress(2, co) = invar3_ip/scaleFac;
635 }
636}
637
638template<typename TDomain>
639void
642{
643 /*if (m_stressEV)
644 {
645 #ifdef UG_PARALLEL
646 if (pcl::ProcRank() == 0){
647 m_fileStressEV = fopen("sig_eigen.dat", "w");
648 }
649 #else
650 m_fileStressEV = fopen("sig_eigen.dat", "w");
651 #endif
652 }*/
653}
654
655template<typename TDomain>
656template<typename TFEGeom>
657void
659post_timestep_elem(const number time, SmartPtr<TDomain> dom, TFEGeom& geo,
660 TBaseElem* elem, const LocalVector& u)
661{
662 typedef typename TDomain::position_accessor_type position_accessor_type;
663 position_accessor_type& aaPos = dom->position_accessor();
664
665 if (m_stressEV && (!m_bIP_values_written))
666 if (ContainsPoint(elem, m_evalPointEV, aaPos))
667 stress_eigenvalues_near_point(time, geo, u);
668
669 // COMPUTE NORMAL STRESSES OF THE STRESSTENSOR SIGMA AT IP
670 if (m_normalStress)
671 if (ContainsPoint(elem, m_evalPointNormStress, aaPos))
672 normal_stress_near_point(time, geo, u);
673}
674
675template<typename TDomain>
676void
678post_timestep(const number time)
679{
680 // write some material law data to console
681 m_spMatLaw->write_data_to_console(time);
682}
683
684template<typename TDomain>
685void
688{
689 /*if (m_stressEV)
690 {
691 #ifdef UG_PARALLEL
692 if (pcl::ProcRank() == 0){
693 fclose(m_fileStressEV);
694 }
695 #else
696 fclose(m_fileStressEV);
697 #endif
698 }*/
699}
700
701template<typename TDomain>
702template<typename TFEGeom>
703void
705stress_eigenvalues_near_point(const number time, TFEGeom& geo,
706 const LocalVector& u)
707{
708 vector<size_t> vNextIP;
709 next_ips_to_point(vNextIP, m_evalPointEV, geo);
710
711 MathMatrix<dim, dim> SYMMSig, Sig, GradU;
712 for (vector<size_t>::iterator it = vNextIP.begin();
713 it != vNextIP.end(); ++it)
714 {
715 m_spMatLaw->template DisplacementGradient<TFEGeom>(GradU, *it, geo, u);
716 m_spMatLaw->stressTensor(Sig, *it, GradU);
717
718 if (dim == 3)
719 {
720 MatMultiplyMTM(SYMMSig, Sig);
721
723 for (size_t a = 0; a < 3; ++a)
724 for (size_t b = 0; b < 3; ++b){
725 SigTens[a][b] = SYMMSig[a][b];
726 }
727
728 MathVector<3, number> evMin, evMed, evMax;
729 number SiglambdaMin, SiglambdaMed, SiglambdaMax;
730
731 CalculateEigenvalues(SigTens, SiglambdaMin,
732 SiglambdaMed, SiglambdaMax, evMin, evMed,
733 evMax);
734
735 number sqrSigLambdaMin = sqrt(SiglambdaMin);
736 number sqrSigLambdaMed = sqrt(SiglambdaMed);
737 number sqrSigLambdaMax = sqrt(SiglambdaMax);
738
739 UG_LOG("minimal EigenValueSigma: " << sqrSigLambdaMin << "\n");
740 UG_LOG("medium EigenValueSigma: " << sqrSigLambdaMed << "\n");
741 UG_LOG("maximal EigenValueSigma: " << sqrSigLambdaMax << "\n");
742
743 /* number absSiglambdamax = abs(sqrSigLambdaMax);
744
745 #ifdef UG_PARALLEL
746 if (pcl::ProcRank() == 0){
747 fprintf(m_fileStressEV, "%lg %lg \n ", time, absSiglambdamax);
748 }
749 #else
750 fprintf(m_fileStressEV, "%lg %lg \n ", time, absSiglambdamax);
751 #endif*/
752
753 }
754
755 // to be sure that data is only written once:
756 m_bIP_values_written = true;
757 } // end vNextIP-iteration
758}
759
760template<typename TDomain>
761template<typename TFEGeom>
762void
764normal_stress_near_point(const number time, TFEGeom& geo, const LocalVector& u)
765{
766 vector<size_t> vNextIP;
767 next_ips_to_point(vNextIP, m_evalPointNormStress, geo);
768
769 MathMatrix<dim, dim> Sig, GradU;
770 for (vector<size_t>::iterator it = vNextIP.begin();
771 it != vNextIP.end(); ++it)
772 {
773 m_spMatLaw->template DisplacementGradient<TFEGeom>(GradU, *it, geo, u);
774 m_spMatLaw->stressTensor(Sig, *it, GradU);
775
776 UG_LOG("At " << geo.global_ip(*it) << ": \n");
777 UG_LOG("sigma_xx: " << Sig[0][0] << "\n");
778 UG_LOG("sigma_yy: " << Sig[1][1] << "\n");
779 UG_LOG("sigma: " << Sig << "\n");
780 UG_LOG("######### \n");
781
782 // to be sure that data is only written once:
783 //m_bIP_values_written = true;
784 } // end vNextIP-iteration
785}
786
787template<typename TDomain>
788template<typename TFEGeom>
789void
791next_ips_to_point(vector<size_t>& vNextIP, const MathVector<dim>& point,
792 const TFEGeom& geo)
793{
794 number dist = std::numeric_limits<number>::max();
795
796 // determine shortest distance from an ip to a point
797 for (size_t ip = 0; ip < geo.num_ip(); ++ip)
798 {
799 number dist_ip = VecDistance(geo.global_ip(ip), point);
800 if (dist_ip < dist)
801 dist = dist_ip;
802 }
803
804 // determine all ip`s with shortest distance to a point
805 for (size_t ip = 0; ip < geo.num_ip(); ++ip)
806 {
807 number dist_ip = VecDistance(geo.global_ip(ip), point);
808 if (dist_ip == dist)
809 vNextIP.push_back(ip);
810 }
811}
812
813
814} // namespace SmallStrainMechanics
815} // namespace ug
816
817#endif /* MECH_OUTPUT_WRITER_IMPL_H_ */
void update(GridObject *pElem, const MathVector< worldDim > *vCorner)
size_t num_ip() const
void resize(const LocalIndices &ind)
Definition mech_output_writer.h:52
domain_traits< TDomain::dim >::grid_base_object TBaseElem
base element type of associated domain
Definition mech_output_writer.h:55
void post_timestep(const number time)
Definition mech_output_writer_impl.h:678
void preprocess()
Definition mech_output_writer_impl.h:641
void stress_eigenvalues_near_point(const number time, TFEGeom &geo, const LocalVector &u)
Definition mech_output_writer_impl.h:705
void postprocess()
Definition mech_output_writer_impl.h:687
void post_timestep_elem(const number time, SmartPtr< TDomain > dom, TFEGeom &geo, TBaseElem *elem, const LocalVector &u)
Definition mech_output_writer_impl.h:659
void invariants_kirchhoff_stress_elem(LocalVector &locInvarKirchhoffStress, TBaseElem *elem, const LocalVector &u, SmartPtr< TDomain > dom)
get the invariants of the kirchhoff stress tensor using the vector of unknowns u
Definition mech_output_writer_impl.h:552
void plastIP_elem(LocalVector &locPlastIP, TBaseElem *elem, const LocalVector &locU, SmartPtr< TDomain > dom)
marks the integration point, where a plastic evolution occurs
Definition mech_output_writer_impl.h:223
void next_ips_to_point(vector< size_t > &vNextIP, const MathVector< dim > &point, const TFEGeom &geo)
Definition mech_output_writer_impl.h:791
void normal_stress_strain_elem(LocalVector &locSigma, LocalVector &locEps, LocalVector &locStressFunc, LocalVector &locStrainFunc, TBaseElem *elem, const LocalVector &locU, SmartPtr< TDomain > dom)
get some normalized values of the Cauchy-stress and strain tensor
Definition mech_output_writer_impl.h:367
void eqPlastStrain_elem(LocalVector &locEqPlastStrain, TBaseElem *elem, const LocalVector &locU, SmartPtr< TDomain > dom)
get the equivalent plastic strain
Definition mech_output_writer_impl.h:302
void normal_stress_near_point(const number time, TFEGeom &geo, const LocalVector &u)
Definition mech_output_writer_impl.h:764
SmartPtr< TGrid > grid()
UG_API bool ContainsPoint(const EdgeVertices *e, const vector_t &p, TAAPos aaPos)
void CollectNeighbors(std::vector< Edge * > &vNeighborsOut, Edge *e, Grid &grid, NeighborhoodType nbhType=NHT_VERTEX_NEIGHBORS)
NHT_VERTEX_NEIGHBORS
matrix_t::value_type MatDeviatorTrace(const matrix_t &m, matrix_t &dev)
MathMatrix< 1, 1, T >::value_type Trace(const MathMatrix< 1, 1, T > &m)
void MatIdentity(matrix_t &mOut)
void MatMultiplyMTM(MathMatrix< N, N, T > &mOut, const MathMatrix< M, N, T > &m)
void MatAdd(matrix_t &mOut, const matrix_t &m, typename matrix_t::value_type s)
void MatScale(matrix_t &mOut, typename matrix_t::value_type s, const matrix_t &m)
matrix_t::value_type MatFrobeniusNorm(matrix_t &m)
void MatMultiply(MathMatrix< N, M, T > &mOut, const MathMatrix< N, L, T > &m1, const MathMatrix< L, M, T > &m2)
void equiv_plast_strain(TGridFunction &eqPlastStrain, MechOutputWriter< typename TGridFunction::domain_type > &mechOut, TGridFunction &u)
Definition mech_output_writer_impl.h:86
void normal_stresses_strains(MechOutputWriter< typename TGridFunction::domain_type > &mechOut, TGridFunction &sigma, TGridFunction &epsilon, TGridFunction &stressFunc, TGridFunction &strainFunc, TGridFunction &u)
Definition mech_output_writer_impl.h:127
void invariants_kirchhoff_stress(TGridFunction &invarKirchhoffStress, MechOutputWriter< typename TGridFunction::domain_type > &mechOut, TGridFunction &u)
Definition mech_output_writer_impl.h:178
void plast_ip(TGridFunction &plastIP, MechOutputWriter< typename TGridFunction::domain_type > &mechOut, TGridFunction &u)
Definition mech_output_writer_impl.h:45
#define UG_CATCH_THROW(msg)
#define UG_LOG(msg)
double number
bool CalculateEigenvalues(const ug::matrix33 &mat, number &lambdaMinOut, number &lambdaMedOut, number &lambdaMaxOut, ug::vector3 &evMinOut, ug::vector3 &evMedOut, ug::vector3 &evMaxOut)
vector_t::value_type VecDistance(const vector_t &v1, const vector_t &v2)
MathMatrix< 0, 0, T >::value_type Determinant(const MathMatrix< 0, 0, T > &m)
void GetLocalVector(LocalVector &lvec, const TVector &vec)
void AddLocalVector(TVector &vec, const LocalVector &lvec)