Loading [MathJax]/extensions/tex2jax.js
Plugins
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
mech_output_writer.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_H_
34#define MECH_OUTPUT_WRITER_H_
35
36using namespace std;
37
38namespace ug{
39namespace SmallStrainMechanics{
40
44
50template<typename TDomain>
52{
53 private:
56
57 public:
59 static const int dim = TDomain::dim;
60
61 public:
66
68 void stress_eigenvalues_at(const number CoordX,
69 const number CoordY, const number CoordZ)
70 {
71 m_stressEV = true;
72 m_evalPointEV[0] = CoordX;
73 m_evalPointEV[1] = CoordY;
74 if (dim == 3)
75 m_evalPointEV[2] = CoordZ;
76 }
77
79 void normal_stresses_at(const number CoordX,
80 const number CoordY, const number CoordZ)
81 {
82 // TODO:check, if point is within the domain!
83 // if not, set m_normalStress to false!
84 m_normalStress = true;
85 m_evalPointNormStress[0] = CoordX;
86 m_evalPointNormStress[1] = CoordY;
87 if (dim == 3)
88 m_evalPointNormStress[2] = CoordZ;
89 }
90
93
94 public:
95 void preprocess();
97 { m_bIP_values_written = false;}
98
99 template<typename TFEGeom>
100 void post_timestep_elem(const number time, SmartPtr<TDomain> dom,
101 TFEGeom& geo, TBaseElem* elem, const LocalVector& u);
102
103 void post_timestep(const number time);
104 void postprocess();
105
107 {m_spMatLaw = spMatLaw;}
108 void quad_order(const int quadOrder)
109 {m_quadOrder = quadOrder;}
110
112 void plastIP_elem(LocalVector& locPlastIP, TBaseElem* elem,
113 const LocalVector& locU, SmartPtr<TDomain> dom);
114
116 void eqPlastStrain_elem(LocalVector& locEqPlastStrain, TBaseElem* elem, const LocalVector& locU,
118
120 void normal_stress_strain_elem(LocalVector& locSigma, LocalVector& locEps, LocalVector& locStressFunc,
121 LocalVector& locStrainFunc, TBaseElem* elem, const LocalVector& locU, SmartPtr<TDomain> dom);
122
124 void invariants_kirchhoff_stress_elem(LocalVector& locInvarKirchhoffStress,
125 TBaseElem* elem, const LocalVector& u, SmartPtr<TDomain> dom);
126
127 private:
128 template<typename TFEGeom>
129 void next_ips_to_point(vector<size_t>& vNextIP, const MathVector<dim>& point,
130 const TFEGeom& geo);
131
132 template<typename TFEGeom>
134 TFEGeom& geo, const LocalVector& u);
135
136 template<typename TFEGeom>
137 void normal_stress_near_point(const number time,
138 TFEGeom& geo, const LocalVector& u);
139
140 private:
143
146
149
153
157
160
161};
162
163// some global output functions
164
165template<typename TGridFunction>
167 TGridFunction& sigma, TGridFunction& epsilon, TGridFunction& stressFunc,
168 TGridFunction& strainFunc, TGridFunction& u);
169
170template <typename TGridFunction>
171void plast_ip(TGridFunction& plastIP,
173 TGridFunction& u);
174
175template <typename TGridFunction>
176void equiv_plast_strain(TGridFunction& eqPlastStrain,
178 TGridFunction& u);
179
180template <typename TGridFunction>
181void invariants_kirchhoff_stress(TGridFunction& invarKirchhoffStress,
183 TGridFunction& u);
184
185// end group small_strain_mechanics
187
188} //end of namespace SmallStrainMechanics
189} //end of namespace ug
190
192
193#endif /* MECH_OUTPUT_WRITER_H_ */
Definition mat_law_interface.h:44
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_at(const number CoordX, const number CoordY, const number CoordZ)
get stress eigenvalues at point
Definition mech_output_writer.h:68
void stress_eigenvalues_near_point(const number time, TFEGeom &geo, const LocalVector &u)
Definition mech_output_writer_impl.h:705
bool m_stressEV
flags for stress-eigenvalues / normalStress at specific corner
Definition mech_output_writer.h:151
bool m_bIP_values_written
flag indicating, if ip-values are already written
Definition mech_output_writer.h:148
int m_quadOrder
quadrature rule
Definition mech_output_writer.h:145
FILE * m_fileStressEV
output-file for stress eigenvalues
Definition mech_output_writer.h:159
void quad_order(const int quadOrder)
Definition mech_output_writer.h:108
MathVector< dim > m_evalPointEV
points, where to get stress-eigenvalues/normalStress
Definition mech_output_writer.h:155
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
SmartPtr< IMaterialLaw< TDomain > > m_spMatLaw
material law
Definition mech_output_writer.h:142
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
static const int dim
World dimension.
Definition mech_output_writer.h:59
~MechOutputWriter()
destructor
Definition mech_output_writer.h:92
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_stresses_at(const number CoordX, const number CoordY, const number CoordZ)
get normal stresses at point
Definition mech_output_writer.h:79
MathVector< dim > m_evalPointNormStress
Definition mech_output_writer.h:156
bool m_normalStress
Definition mech_output_writer.h:152
void pre_timestep()
Definition mech_output_writer.h:96
MechOutputWriter()
constructor
Definition mech_output_writer.h:63
void material_law(SmartPtr< IMaterialLaw< TDomain > > spMatLaw)
Definition mech_output_writer.h:106
void normal_stress_near_point(const number time, TFEGeom &geo, const LocalVector &u)
Definition mech_output_writer_impl.h:764
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
double number