Line data Source code
1 : //* This file is part of the MOOSE framework
2 : //* https://mooseframework.inl.gov
3 : //*
4 : //* All rights reserved, see COPYRIGHT for full restrictions
5 : //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 : //*
7 : //* Licensed under LGPL 2.1, please see LICENSE for details
8 : //* https://www.gnu.org/licenses/lgpl-2.1.html
9 :
10 : #include "ADComputeFiniteShellStrain.h"
11 : #include "DenseMatrix.h"
12 :
13 : #include "libmesh/quadrature.h"
14 : #include "libmesh/utility.h"
15 : #include "libmesh/enum_quadrature_type.h"
16 : #include "libmesh/fe_type.h"
17 : #include "libmesh/string_to_enum.h"
18 : #include "libmesh/quadrature_gauss.h"
19 :
20 : registerMooseObject("SolidMechanicsApp", ADComputeFiniteShellStrain);
21 :
22 : InputParameters
23 72 : ADComputeFiniteShellStrain::validParams()
24 : {
25 72 : InputParameters params = ADComputeIncrementalShellStrain::validParams();
26 72 : params.addClassDescription("Compute a large strain increment for the shell.");
27 72 : return params;
28 0 : }
29 :
30 54 : ADComputeFiniteShellStrain::ADComputeFiniteShellStrain(const InputParameters & parameters)
31 54 : : ADComputeIncrementalShellStrain(parameters), _B_nl()
32 : {
33 54 : _B_nl.resize(_t_points.size());
34 :
35 162 : for (unsigned int i = 0; i < _t_points.size(); ++i)
36 216 : _B_nl[i] = &declareADProperty<DenseMatrix<Real>>("B_nl_t_points_" + std::to_string(i));
37 54 : }
38 :
39 : void
40 352 : ADComputeFiniteShellStrain::initQpStatefulProperties()
41 : {
42 352 : computeGMatrix();
43 :
44 352 : ADDenseMatrix b(5, 20);
45 1056 : for (unsigned int t = 0; t < _t_points.size(); ++t)
46 704 : (*_B_nl[t])[_qp] = b;
47 352 : }
48 :
49 : void
50 1956 : ADComputeFiniteShellStrain::computeProperties()
51 : {
52 : // quadrature points in isoparametric space
53 1956 : _2d_points = _qrule->get_points(); // would be in 2D
54 :
55 1956 : unsigned int dim = _current_elem->dim();
56 :
57 : // derivatives of shape functions (dphidxi, dphideta and dphidzeta) evaluated at quadrature points
58 : // (in isoparametric space).
59 1956 : FEType fe_type(Utility::string_to_enum<Order>("First"),
60 1956 : Utility::string_to_enum<FEFamily>("LAGRANGE"));
61 1956 : auto & fe = _fe_problem.assembly(_tid, /*nl_sys_num=*/0).getFE(fe_type, dim);
62 1956 : _dphidxi_map = fe->get_fe_map().get_dphidxi_map();
63 1956 : _dphideta_map = fe->get_fe_map().get_dphideta_map();
64 1956 : _phi_map = fe->get_fe_map().get_phi_map();
65 :
66 9780 : for (unsigned int i = 0; i < 4; ++i)
67 7824 : _nodes[i] = _current_elem->node_ptr(i);
68 :
69 9780 : for (unsigned int i = 0; i < _2d_points.size(); ++i)
70 : {
71 23472 : for (unsigned int j = 0; j < _t_points.size(); ++j)
72 : {
73 15648 : (*_ge[j])[i] = (*_ge_old[j])[i];
74 15648 : (*_J_map[j])[i] = (*_J_map_old[j])[i];
75 : }
76 : }
77 :
78 1956 : computeSolnVector();
79 1956 : updatedxyz();
80 :
81 9780 : for (unsigned int k = 0; k < _nodes.size(); ++k)
82 7824 : _node_normal[k] = _node_normal_old[k];
83 :
84 1956 : computeBMatrix();
85 :
86 1956 : computeNodeNormal();
87 :
88 1956 : computeBNLMatrix();
89 :
90 9780 : for (unsigned int i = 0; i < _2d_points.size(); ++i)
91 : {
92 23472 : for (unsigned int j = 0; j < _t_points.size(); ++j)
93 : {
94 : // compute strain increment in covariant coordinate system using B and _soln_vector
95 93888 : for (unsigned int temp1 = 0; temp1 < 5; ++temp1)
96 : {
97 78240 : _strain_vector(temp1) = 0.0;
98 1643040 : for (unsigned int temp2 = 0; temp2 < 20; ++temp2)
99 3129600 : _strain_vector(temp1) += (*_B[j])[i](temp1, temp2) * _soln_vector(temp2);
100 : }
101 :
102 15648 : (*_strain_increment[j])[i](0, 0) = _strain_vector(0);
103 15648 : (*_strain_increment[j])[i](1, 1) = _strain_vector(1);
104 15648 : (*_strain_increment[j])[i](0, 1) = _strain_vector(2);
105 15648 : (*_strain_increment[j])[i](0, 2) = _strain_vector(3);
106 15648 : (*_strain_increment[j])[i](1, 2) = _strain_vector(4);
107 15648 : (*_strain_increment[j])[i](1, 0) = (*_strain_increment[j])[i](0, 1);
108 15648 : (*_strain_increment[j])[i](2, 0) = (*_strain_increment[j])[i](0, 2);
109 15648 : (*_strain_increment[j])[i](2, 1) = (*_strain_increment[j])[i](1, 2);
110 15648 : (*_total_strain[j])[i] = (*_total_strain_old[j])[i] + (*_strain_increment[j])[i];
111 62592 : for (unsigned int ii = 0; ii < 3; ++ii)
112 187776 : for (unsigned int jj = 0; jj < 3; ++jj)
113 140832 : _unrotated_total_strain(ii, jj) = MetaPhysicL::raw_value((*_total_strain[j])[i](ii, jj));
114 15648 : (*_total_global_strain[j])[i] = (*_contravariant_transformation_matrix[j])[i] *
115 15648 : _unrotated_total_strain *
116 15648 : (*_contravariant_transformation_matrix[j])[i].transpose();
117 : }
118 : }
119 1956 : }
120 :
121 : void
122 1956 : ADComputeFiniteShellStrain::computeNodeNormal()
123 : {
124 : // update _node_normal
125 9780 : for (unsigned int k = 0; k < _nodes.size(); ++k)
126 : {
127 7824 : _node_normal[k] =
128 7824 : -_v2[k] * _soln_vector(12 + k) + _v1[k] * _soln_vector(16 + k) + _node_normal_old[k];
129 15648 : _node_normal[k] /= _node_normal[k].norm();
130 : }
131 1956 : }
132 :
133 : void
134 1956 : ADComputeFiniteShellStrain::updatedxyz()
135 : {
136 9780 : for (unsigned int i = 0; i < _2d_points.size(); ++i)
137 : {
138 23472 : for (unsigned int j = 0; j < _t_points.size(); ++j)
139 : {
140 62592 : for (unsigned int component = 0; component < 3; ++component)
141 : {
142 46944 : (*_dxyz_dxi[j])[i](component) = 0.0;
143 46944 : (*_dxyz_deta[j])[i](component) = 0.0;
144 46944 : (*_dxyz_dzeta[j])[i](component) = 0.0;
145 234720 : for (unsigned int k = 0; k < _nodes.size(); ++k)
146 : {
147 187776 : (*_dxyz_dxi[j])[i](component) +=
148 375552 : _dphidxi_map[k][i] *
149 187776 : ((*_nodes[k])(component) + _sol_old(_soln_disp_index[k][component])) +
150 187776 : _t_points[j](0) / 2.0 * _thickness[i] * _dphidxi_map[k][i] *
151 187776 : _node_normal_old[k](component);
152 187776 : (*_dxyz_deta[j])[i](component) +=
153 375552 : _dphideta_map[k][i] *
154 187776 : ((*_nodes[k])(component) + _sol_old(_soln_disp_index[k][component])) +
155 187776 : _t_points[j](0) / 2.0 * _thickness[i] * _dphideta_map[k][i] *
156 187776 : _node_normal_old[k](component);
157 187776 : (*_dxyz_dzeta[j])[i](component) +=
158 187776 : _thickness[i] * _phi_map[k][i] * _node_normal_old[k](component) / 2.0;
159 : }
160 : }
161 : }
162 : }
163 1956 : }
164 :
165 : void
166 15648 : ADComputeFiniteShellStrain::updateGVectors()
167 : {
168 62592 : for (unsigned int component = 0; component < 3; ++component)
169 : {
170 : _g1_a(component) +=
171 46944 : 0.5 * (_sol_old(_soln_disp_index[2][component]) - _sol_old(_soln_disp_index[3][component]));
172 : _g1_c(component) +=
173 46944 : 0.5 * (_sol_old(_soln_disp_index[1][component]) - _sol_old(_soln_disp_index[0][component]));
174 : _g2_b(component) +=
175 46944 : 0.5 * (_sol_old(_soln_disp_index[3][component]) - _sol_old(_soln_disp_index[0][component]));
176 : _g2_d(component) +=
177 46944 : 0.5 * (_sol_old(_soln_disp_index[2][component]) - _sol_old(_soln_disp_index[1][component]));
178 : }
179 15648 : }
180 :
181 : void
182 1956 : ADComputeFiniteShellStrain::computeBNLMatrix()
183 : {
184 : // compute BNL matrix - rows correspond to [ux1, ux2, ux3, ux4, uy1, uy2, uy3, uy4, uz1, uz2, uz3,
185 : // uz4, a1, a2, a3, a4, b1, b2, b3, b4]
186 :
187 9780 : for (unsigned int i = 0; i < _2d_points.size(); ++i)
188 : {
189 23472 : for (unsigned int j = 0; j < _t_points.size(); ++j)
190 : {
191 15648 : (*_B_nl[j])[i].resize(5, 20);
192 15648 : (*_B_nl[j])[i].zero();
193 78240 : for (unsigned int k = 0; k < 4; ++k)
194 : {
195 312960 : for (unsigned int p = 0; p < 4; ++p) // loop over nodes
196 : {
197 : // corresponding to strain(0,0)
198 500736 : (*_B_nl[j])[i](0, k) += _dphidxi_map[k][i] * _dphidxi_map[p][i] *
199 250368 : (_soln_vector(p) + _t_points[j](0) / 2.0 * _thickness[i] *
200 250368 : (-_soln_vector(p + 12) * _v2[p](0) +
201 500736 : _soln_vector(p + 16) * _v1[p](0)));
202 250368 : (*_B_nl[j])[i](0, 4 + k) +=
203 250368 : _dphidxi_map[k][i] * _dphidxi_map[p][i] *
204 250368 : (_soln_vector(p + 4) +
205 250368 : _t_points[j](0) / 2.0 * _thickness[i] *
206 500736 : (-_soln_vector(p + 12) * _v2[p](1) + _soln_vector(p + 16) * _v1[p](1)));
207 250368 : (*_B_nl[j])[i](0, 8 + k) +=
208 250368 : _dphidxi_map[k][i] * _dphidxi_map[p][i] *
209 250368 : (_soln_vector(p + 8) +
210 250368 : _t_points[j](0) / 2.0 * _thickness[i] *
211 500736 : (-_soln_vector(p + 12) * _v2[p](2) + _soln_vector(p + 16) * _v1[p](2)));
212 250368 : (*_B_nl[j])[i](0, 12 + k) +=
213 250368 : _t_points[j](0) / 2.0 * _thickness[i] * _dphidxi_map[k][i] * _dphidxi_map[p][i] *
214 250368 : (-(_v2[p](0) * _soln_vector(p) + _v2[p](1) * _soln_vector(p + 4) +
215 250368 : _v2[p](2) * _soln_vector(p + 8)) +
216 250368 : _t_points[j](0) / 2.0 * _thickness[i] * _v2[k] *
217 500736 : (_v2[p] * _soln_vector(p + 12) - _v1[p] * _soln_vector(p + 16)));
218 250368 : (*_B_nl[j])[i](0, 16 + k) +=
219 250368 : _t_points[j](0) / 2.0 * _thickness[i] * _dphidxi_map[k][i] * _dphidxi_map[p][i] *
220 250368 : ((_v1[p](0) * _soln_vector(p) + _v1[p](1) * _soln_vector(p + 4) +
221 250368 : _v1[p](2) * _soln_vector(p + 8)) +
222 250368 : _t_points[j](0) / 2.0 * _thickness[i] * _v1[k] *
223 500736 : (-_v2[p] * _soln_vector(p + 12) + _v1[p] * _soln_vector(p + 16)));
224 :
225 : // corresponding to strain(1,1)
226 500736 : (*_B_nl[j])[i](1, k) += _dphideta_map[k][i] * _dphideta_map[p][i] *
227 250368 : (_soln_vector(p) + _t_points[j](0) / 2.0 * _thickness[i] *
228 250368 : (-_soln_vector(p + 12) * _v2[p](0) +
229 250368 : _soln_vector(p + 16) * _v1[p](0)));
230 250368 : (*_B_nl[j])[i](1, 4 + k) +=
231 250368 : _dphideta_map[k][i] * _dphideta_map[p][i] *
232 250368 : (_soln_vector(p + 4) +
233 250368 : _t_points[j](0) / 2.0 * _thickness[i] *
234 500736 : (-_soln_vector(p + 12) * _v2[p](1) + _soln_vector(p + 16) * _v1[p](1)));
235 250368 : (*_B_nl[j])[i](1, 8 + k) +=
236 250368 : _dphideta_map[k][i] * _dphideta_map[p][i] *
237 250368 : (_soln_vector(p + 8) +
238 250368 : _t_points[j](0) / 2.0 * _thickness[i] *
239 500736 : (-_soln_vector(p + 12) * _v2[p](2) + _soln_vector(p + 16) * _v1[p](2)));
240 250368 : (*_B_nl[j])[i](1, 12 + k) +=
241 250368 : _t_points[j](0) / 2.0 * _thickness[i] * _dphideta_map[k][i] * _dphideta_map[p][i] *
242 250368 : (-(_v2[p](0) * _soln_vector(p) + _v2[p](1) * _soln_vector(p + 4) +
243 250368 : _v2[p](2) * _soln_vector(p + 8)) +
244 250368 : _t_points[j](0) / 2.0 * _thickness[i] * _v2[k] *
245 500736 : (_v2[p] * _soln_vector(p + 12) - _v1[p] * _soln_vector(p + 16)));
246 250368 : (*_B_nl[j])[i](1, 16 + k) +=
247 250368 : _t_points[j](0) / 2.0 * _thickness[i] * _dphideta_map[k][i] * _dphideta_map[p][i] *
248 250368 : ((_v1[p](0) * _soln_vector(p) + _v1[p](1) * _soln_vector(p + 4) +
249 250368 : _v1[p](2) * _soln_vector(p + 8)) +
250 250368 : _t_points[j](0) / 2.0 * _thickness[i] * _v1[k] *
251 500736 : (-_v2[p] * _soln_vector(p + 12) + _v1[p] * _soln_vector(p + 16)));
252 :
253 : // terms corresponding to strain(2,2) are 0.
254 :
255 : // corresponding to strain(0,1)
256 250368 : (*_B_nl[j])[i](2, k) += 0.5 *
257 250368 : (_dphidxi_map[k][i] * _dphideta_map[p][i] +
258 250368 : _dphideta_map[k][i] * _dphidxi_map[p][i]) *
259 250368 : (_soln_vector(p) + _t_points[j](0) / 2.0 * _thickness[i] *
260 250368 : (-_soln_vector(p + 12) * _v2[p](0) +
261 250368 : _soln_vector(p + 16) * _v1[p](0)));
262 250368 : (*_B_nl[j])[i](2, 4 + k) +=
263 500736 : 0.5 *
264 250368 : (_dphidxi_map[k][i] * _dphideta_map[p][i] +
265 250368 : _dphideta_map[k][i] * _dphidxi_map[p][i]) *
266 250368 : (_soln_vector(p + 4) +
267 250368 : _t_points[j](0) / 2.0 * _thickness[i] *
268 500736 : (-_soln_vector(p + 12) * _v2[p](1) + _soln_vector(p + 16) * _v1[p](1)));
269 250368 : (*_B_nl[j])[i](2, 8 + k) +=
270 500736 : 0.5 *
271 250368 : (_dphidxi_map[k][i] * _dphideta_map[p][i] +
272 250368 : _dphideta_map[k][i] * _dphidxi_map[p][i]) *
273 250368 : (_soln_vector(p + 8) +
274 250368 : _t_points[j](0) / 2.0 * _thickness[i] *
275 500736 : (-_soln_vector(p + 12) * _v2[p](2) + _soln_vector(p + 16) * _v1[p](2)));
276 250368 : (*_B_nl[j])[i](2, 12 + k) +=
277 250368 : _t_points[j](0) * 0.25 *
278 250368 : (_dphidxi_map[k][i] * _dphideta_map[p][i] +
279 500736 : _dphideta_map[k][i] * _dphidxi_map[p][i]) *
280 250368 : _thickness[i] *
281 250368 : (-(_v2[k](0) * _soln_vector(p) + _v2[k](1) * _soln_vector(p + 4) +
282 250368 : _v2[k](2) * _soln_vector(p + 8)) +
283 250368 : _t_points[j](0) / 2.0 * _thickness[i] * _v2[k] *
284 500736 : (_v2[p] * _soln_vector(p + 12) - _v1[p] * _soln_vector(p + 16)));
285 250368 : (*_B_nl[j])[i](2, 16 + k) +=
286 250368 : _t_points[j](0) * 0.25 *
287 250368 : (_dphidxi_map[k][i] * _dphideta_map[p][i] +
288 500736 : _dphideta_map[k][i] * _dphidxi_map[p][i]) *
289 250368 : _thickness[i] *
290 250368 : ((_v1[k](0) * _soln_vector(p) + _v1[k](1) * _soln_vector(p + 4) +
291 250368 : _v1[k](2) * _soln_vector(p + 8)) +
292 250368 : _t_points[j](0) / 2.0 * _thickness[i] * _v1[k] *
293 500736 : (-_v2[p] * _soln_vector(p + 12) + _v1[p] * _soln_vector(p + 16)));
294 : }
295 : }
296 :
297 62592 : for (unsigned int component = 0; component < 3; ++component)
298 : {
299 : // corresponding to strain(0,2)
300 46944 : (*_B_nl[j])[i](3, 2 + component * 4) +=
301 46944 : 1.0 / 32.0 * (1.0 + _2d_points[i](1)) * _thickness[i] *
302 46944 : (-_soln_vector(12 + 2) * _v2[2](component) + _soln_vector(16 + 2) * _v1[2](component) -
303 46944 : _soln_vector(12 + 3) * _v2[3](component) + _soln_vector(16 + 3) * _v1[3](component));
304 46944 : (*_B_nl[j])[i](3, 3 + component * 4) += -(*_B_nl[j])[i](3, 2 + component * 4);
305 :
306 46944 : (*_B_nl[j])[i](3, 1 + component * 4) +=
307 46944 : 1.0 / 32.0 * (1.0 - _2d_points[i](1)) * _thickness[i] *
308 46944 : (-_soln_vector(12 + 1) * _v2[1](component) + _soln_vector(16 + 1) * _v1[1](component) -
309 46944 : _soln_vector(12 + 0) * _v2[0](component) + _soln_vector(16 + 0) * _v1[0](component));
310 46944 : (*_B_nl[j])[i](3, component * 4) += -(*_B_nl[j])[i](3, 1 + component * 4);
311 :
312 : // adding contributions corresponding to alpha 2 and 3 and beta 2 and 3
313 46944 : (*_B_nl[j])[i](3, 12 + 2) +=
314 46944 : -1.0 / 32.0 * (1.0 + _2d_points[i](1)) * _thickness[i] * _v2[2](component) *
315 46944 : (_soln_vector(2 + component * 4) - _soln_vector(3 + component * 4));
316 46944 : (*_B_nl[j])[i](3, 16 + 2) +=
317 46944 : 1.0 / 32.0 * (1.0 + _2d_points[i](1)) * _thickness[i] * _v1[2](component) *
318 46944 : (_soln_vector(2 + component * 4) - _soln_vector(3 + component * 4));
319 46944 : (*_B_nl[j])[i](3, 12 + 3) +=
320 46944 : -1.0 / 32.0 * (1.0 + _2d_points[i](1)) * _thickness[i] * _v2[3](component) *
321 46944 : (_soln_vector(2 + component * 4) - _soln_vector(3 + component * 4));
322 46944 : (*_B_nl[j])[i](3, 16 + 3) +=
323 46944 : 1.0 / 32.0 * (1.0 + _2d_points[i](1)) * _thickness[i] * _v1[3](component) *
324 46944 : (_soln_vector(2 + component * 4) - _soln_vector(3 + component * 4));
325 :
326 : // adding contributions corresponding to alpha 1 and 0 and beta 1 and 0
327 46944 : (*_B_nl[j])[i](3, 12 + 1) +=
328 46944 : -1.0 / 32.0 * (1.0 - _2d_points[i](1)) * _thickness[i] * _v2[1](component) *
329 46944 : (_soln_vector(1 + component * 4) - _soln_vector(component * 4));
330 46944 : (*_B_nl[j])[i](3, 16 + 1) +=
331 46944 : 1.0 / 32.0 * (1.0 - _2d_points[i](1)) * _thickness[i] * _v1[1](component) *
332 46944 : (_soln_vector(1 + component * 4) - _soln_vector(component * 4));
333 46944 : (*_B_nl[j])[i](3, 12 + 0) +=
334 46944 : -1.0 / 32.0 * (1.0 - _2d_points[i](1)) * _thickness[i] * _v2[0](component) *
335 46944 : (_soln_vector(1 + component * 4) - _soln_vector(component * 4));
336 46944 : (*_B_nl[j])[i](3, 16 + 0) +=
337 46944 : 1.0 / 32.0 * (1.0 - _2d_points[i](1)) * _thickness[i] * _v1[0](component) *
338 46944 : (_soln_vector(1 + component * 4) - _soln_vector(component * 4));
339 :
340 : // corresponding to strain(1,2)
341 46944 : (*_B_nl[j])[i](4, 2 + component * 4) +=
342 46944 : 1.0 / 32.0 * (1.0 + _2d_points[i](0)) * _thickness[i] *
343 46944 : (-_soln_vector(12 + 2) * _v2[2](component) + _soln_vector(16 + 2) * _v1[2](component) -
344 46944 : _soln_vector(12 + 1) * _v2[1](component) + _soln_vector(16 + 1) * _v1[1](component));
345 46944 : (*_B_nl[j])[i](4, 1 + component * 4) += -(*_B_nl[j])[i](3, 2 + component * 4);
346 :
347 46944 : (*_B_nl[j])[i](4, 3 + component * 4) +=
348 46944 : 1.0 / 32.0 * (1.0 - _2d_points[i](0)) * _thickness[i] *
349 46944 : (-_soln_vector(12 + 3) * _v2[3](component) + _soln_vector(16 + 3) * _v1[3](component) -
350 46944 : _soln_vector(12 + 0) * _v2[0](component) + _soln_vector(16 + 0) * _v1[0](component));
351 46944 : (*_B_nl[j])[i](4, component * 4) += -(*_B_nl[j])[i](3, 3 + component * 4);
352 :
353 : // adding contributions corresponding to alpha 2, 1 and beta 2 , 1
354 46944 : (*_B_nl[j])[i](4, 12 + 2) +=
355 46944 : -1.0 / 32.0 * (1.0 + _2d_points[i](0)) * _thickness[i] * _v2[2](component) *
356 46944 : (_soln_vector(2 + component * 4) - _soln_vector(1 + component * 4));
357 46944 : (*_B_nl[j])[i](4, 16 + 2) +=
358 46944 : 1.0 / 32.0 * (1.0 + _2d_points[i](0)) * _thickness[i] * _v1[2](component) *
359 46944 : (_soln_vector(2 + component * 4) - _soln_vector(1 + component * 4));
360 46944 : (*_B_nl[j])[i](4, 12 + 1) +=
361 46944 : -1.0 / 32.0 * (1.0 + _2d_points[i](0)) * _thickness[i] * _v2[1](component) *
362 46944 : (_soln_vector(2 + component * 4) - _soln_vector(1 + component * 4));
363 46944 : (*_B_nl[j])[i](4, 16 + 1) +=
364 46944 : 1.0 / 32.0 * (1.0 + _2d_points[i](0)) * _thickness[i] * _v1[1](component) *
365 46944 : (_soln_vector(2 + component * 4) - _soln_vector(1 + component * 4));
366 :
367 : // adding contributions corresponding to alpha 3, 0 and beta 3 , 0
368 46944 : (*_B_nl[j])[i](4, 12 + 3) +=
369 46944 : -1.0 / 32.0 * (1.0 - _2d_points[i](0)) * _thickness[i] * _v2[3](component) *
370 46944 : (_soln_vector(3 + component * 4) - _soln_vector(component * 4));
371 46944 : (*_B_nl[j])[i](4, 16 + 3) +=
372 46944 : 1.0 / 32.0 * (1.0 - _2d_points[i](0)) * _thickness[i] * _v1[3](component) *
373 46944 : (_soln_vector(3 + component * 4) - _soln_vector(component * 4));
374 46944 : (*_B_nl[j])[i](4, 12 + 0) +=
375 46944 : -1.0 / 32.0 * (1.0 - _2d_points[i](0)) * _thickness[i] * _v2[0](component) *
376 46944 : (_soln_vector(3 + component * 4) - _soln_vector(component * 4));
377 46944 : (*_B_nl[j])[i](4, 16 + 0) +=
378 46944 : 1.0 / 32.0 * (1.0 - _2d_points[i](0)) * _thickness[i] * _v1[0](component) *
379 46944 : (_soln_vector(3 + component * 4) - _soln_vector(component * 4));
380 : }
381 : }
382 : }
383 1956 : }
|