fosanalysis
A framework to evaluate distributed fiber optic sensor data
Loading...
Searching...
No Matches
tensionstiffening.py
Go to the documentation of this file.
2r"""
3Contains class definitions for tension stiffening influences for concrete embedded and reinforcement attached sensors.
4\author Bertram Richter
5\date 2022
6"""
7
8from abc import abstractmethod
9
10import numpy as np
11
12from fosanalysis.utils import misc
13from . import compensator
14
16 r"""
17 Abstract base class for tension stiffening compensation approaches.
18 """
19 def __init__(self,
20 *args, **kwargs):
21 r"""
22 Constructs a TensionStiffeningCompensator object.
23 \param *args Additional positional arguments, will be passed to the superconstructor.
24 \param **kwargs Additional keyword arguments, will be passed to the superconstructor.
25 """
26 super().__init__(*args, **kwargs)
27 @abstractmethod
28 def run(self, x: np.array, strain: np.array, crack_list: list, *args, **kwargs) -> np.array:
29 r"""
30 Compensates for the strain, that does not contribute to a crack, but is located in the uncracked concrete.
31 An array with the compensation values for each measuring point is returned.
32 Every \ref TensionStiffeningCompensator has a \ref run() method, which takes the following parameters:
33 \param x Positional x values.
34 \param strain List of strain values.
35 \param crack_list \ref fosanalysis.crackmonitoring.cracks.CrackList with \ref fosanalysis.crackmonitoring.cracks.Crack objects, that already have assigned locations.
36 \param *args Additional positional arguments to customize the behavior.
37 \param **kwargs Additional keyword arguments to customize the behavior.
38 """
39 raise NotImplementedError()
40
42 r"""
43 Implements the tension stiffening approach according to the proposal by \cite Berrocal_2021_Crackmonitoringin.
44 The concrete strain \f$\varepsilon^{\mathrm{ts}}(x)\f$ is assumed to the difference between the real strain profile \f$\varepsilon^{\mathrm{DOFS}}(x)\f$
45 and the linear interpolation between the peaks \f$\hat{\varepsilon}(x)\f$ reduced by the reinforcement ratio \ref rho \f$\rho\f$ and Young's moduli ratio \ref alpha \f$\alpha\f$:
46 \f[
47 \varepsilon^{\mathrm{TS}}(x) = \rho \alpha \left(\hat{\varepsilon}(x) - \varepsilon^{\mathrm{DOFS}}(x)\right).
48 \f]
49 Outside of the outermost cracks (index \f$0\f$ and \f$n\f$), \f$\hat{\varepsilon}(x)\f$ assumed to be constant at the peak strain of the outermost crack:
50 \f[
51 \hat{\varepsilon}(x) =
52 \begin{cases}
53 \varepsilon_{\mathrm{cr},0} &\text{ if } x < x_{\mathrm{cr}, 0}\\
54 \varepsilon_{\mathrm{cr},n} &\text{ if } x_{\mathrm{cr}, n} < x.
55 \end{cases}
56 \f]
57 The interpolation is done using [`numpy.interp()`](https://numpy.org/doc/stable/reference/generated/numpy.interp.html).
58 """
59 def __init__(self,
60 alpha: float,
61 rho: float,
62 *args, **kwargs):
63 r"""
64 Constructs a TensionStiffeningCompensator object with according to the proposal by\cite Berrocal_2021_Crackmonitoringin.
65 \param alpha \copybrief alpha For more, see \ref alpha.
66 \param rho \copybrief rho For more, see \ref rho.
67 \param *args Additional positional arguments, will be passed to the superconstructor.
68 \param **kwargs Additional keyword arguments, will be passed to the superconstructor.
69 """
70 super().__init__(*args, **kwargs)
71
72 self.alpha = alpha
73
74 self.rho = rho
75 def run(self, x: np.array, strain: np.array, crack_list: list, *args, **kwargs) -> np.array:
76 r"""
77 \copydoc TensionStiffeningCompensator.run()
78 """
79 if not crack_list:
80 # crack_list is empty
81 return np.zeros_like(strain)
82 else:
83 tension_stiffening_values = np.interp(x=x, xp=crack_list.locations, fp=crack_list.max_strains)
84 # Difference of steel strain to the linear interpolation
85 tension_stiffening_values = tension_stiffening_values - strain
86 # Reduce by rho and alpha
87 tension_stiffening_values = tension_stiffening_values * self.alpha * self.rho
88 tension_stiffening_values = np.maximum(tension_stiffening_values, 0)
89 return tension_stiffening_values
90
92 r"""
93 Implements the tension stiffening approach based on \cite Fischer_2019_Distributedfiberoptic.
94 The calculative tension stiffening strain \f(\varepsilon^{\mathrm{ts}}_{\mathrm{concrete}}\f) is idealized to increase linearly from the crack's position
95 \f[
96 \varepsilon^{\mathrm{TS}}_{\mathrm{concrete}}(x) = \min{\left(\delta_{\varepsilon}(x) \times \varepsilon_{\mathrm{lim}}(x),\: \varepsilon^{\mathrm{DFOS}}(x)\right)}
97 \f]
98 with the normalized distance to the crack
99 \f[
100 \delta_{\varepsilon}(x) =
101 \begin{cases}
102 \frac{x_{\mathrm{cr}} - x}{l^{-}_{\mathrm{t}}} & \text{if } x \leq x_{\mathrm{cr}}, \\
103 \frac{x - x_{\mathrm{cr}}}{l^{+}_{\mathrm{t}}} & \text{if } x > x_{\mathrm{cr}}
104 \end{cases}
105 \f]
106 and the limit strain
107 \f[
108 \varepsilon_{\mathrm{lim}}(x) =
109 \begin{cases}
110 \min{\left(\varepsilon^{\mathrm{DFOS}}\left(x_{\mathrm{cr}} - l^{-}_{\mathrm{t}}\right),\: \varepsilon_{\mathrm{ctu}} \right)}& \text{if } x \leq x_{\mathrm{cr}}, \\
111 \min{\left(\varepsilon^{\mathrm{DFOS}}\left(x_{\mathrm{cr}} + l^{+}_{\mathrm{t}}\right),\: \varepsilon_{\mathrm{ctu}} \right)}& \text{if } x > x_{\mathrm{cr}}
112 \end{cases}
113 \f]
114 which is the minimum of the rupture strain \f(\varepsilon_{\mathrm{ctu}}\f) and the measured strain at the transfer length end.
115 The interpolation is done using [`numpy.interp()`](https://numpy.org/doc/stable/reference/generated/numpy.interp.html).
116 """
117 def __init__(self,
118 max_concrete_strain: int = 100,
119 *args, **kwargs):
120 r"""
121 Constructs a TensionStiffeningCompensator object with according to the proposal by \cite Fischer_2019_Distributedfiberoptic.
122 \param max_concrete_strain \copybrief max_concrete_strain For more, see \ref max_concrete_strain.
123 \param *args Additional positional arguments, will be passed to the superconstructor.
124 \param **kwargs Additional keyword arguments, will be passed to the superconstructor.
125 """
126 super().__init__(*args, **kwargs)
127
130 self.max_concrete_strain = max_concrete_strain
131 def run(self, x: np.array, strain: np.array, crack_list: list, *args, **kwargs) -> np.array:
132 r"""
133 \copydoc TensionStiffeningCompensator.run()
134 """
135 tension_stiffening_values = np.zeros_like(strain)
136 for crack in crack_list:
137 l_i, x_l = misc.find_closest_value(x, crack.x_l)
138 r_i, x_r = misc.find_closest_value(x, crack.x_r)
139 x_seg = x[l_i:r_i+1]
140 xp = [x_l, crack.location, x_r]
141 fp = np.minimum([strain[l_i], 0, strain[r_i]], self.max_concrete_strain)
142 tension_stiffening_values[l_i:r_i+1] = np.interp(x_seg, xp, fp)
143 tension_stiffening_values = np.minimum(tension_stiffening_values, strain)
144 tension_stiffening_values = np.maximum(tension_stiffening_values, 0)
145 return tension_stiffening_values
Implements the tension stiffening approach according to the proposal by berrocal_2021_crackmonitoring...
alpha
Ratio of Young's moduli of steel to concrete .
np.array run(self, np.array x, np.array strain, list crack_list, *args, **kwargs)
Compensates for the strain, that does not contribute to a crack, but is located in the uncracked conc...
__init__(self, float alpha, float rho, *args, **kwargs)
Constructs a TensionStiffeningCompensator object with according to the proposal byberrocal_2021_crack...
rho
Reinforcement ratio of steel to concrete .
Implements the tension stiffening approach based on fischer_2019_distributedfiberoptic.
np.array run(self, np.array x, np.array strain, list crack_list, *args, **kwargs)
Compensates for the strain, that does not contribute to a crack, but is located in the uncracked conc...
max_concrete_strain
Maximum strain in concrete that the concrete can bear, before a crack opens.
__init__(self, int max_concrete_strain=100, *args, **kwargs)
Constructs a TensionStiffeningCompensator object with according to the proposal by fischer_2019_distr...
Abstract base class for tension stiffening compensation approaches.
__init__(self, *args, **kwargs)
Constructs a TensionStiffeningCompensator object.
np.array run(self, np.array x, np.array strain, list crack_list, *args, **kwargs)
Compensates for the strain, that does not contribute to a crack, but is located in the uncracked conc...
Contains utility modules with general pupose.
Definition __init__.py:1