OOFEM 3.0
Loading...
Searching...
No Matches
lattice3d_mt.C
Go to the documentation of this file.
1/*
2 *
3 * ##### ##### ###### ###### ### ###
4 * ## ## ## ## ## ## ## ### ##
5 * ## ## ## ## #### #### ## # ##
6 * ## ## ## ## ## ## ## ##
7 * ## ## ## ## ## ## ## ##
8 * ##### ##### ## ###### ## ##
9 *
10 *
11 * OOFEM : Object Oriented Finite Element Code
12 *
13 * Copyright (C) 1993 - 2025 Borek Patzak
14 *
15 *
16 *
17 * Czech Technical University, Faculty of Civil Engineering,
18 * Department of Structural Mechanics, 166 29 Prague, Czech Republic
19 *
20 * This program is free software; you can redistribute it and/or modify
21 * it under the terms of the GNU General Public License as published by
22 * the Free Software Foundation; either version 2 of the License, or
23 * (at your option) any later version.
24 *
25 * This program is distributed in the hope that it will be useful,
26 * but WITHOUT ANY WARRANTY; without even the implied warranty of
27 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
28 * GNU General Public License for more details.
29 *
30 * You should have received a copy of the GNU General Public License
31 * along with this program; if not, write to the Free Software
32 * Foundation, Inc., 675 Mass Ave, Cambridge, MA 02139, USA.
33 */
34
38#include "node.h"
39#include "material.h"
40#include "crosssection.h"
41#include "gausspoint.h"
43#include "floatmatrix.h"
44#include "floatarray.h"
45#include "intarray.h"
46#include "domain.h"
47#include "mathfem.h"
48#include "engngm.h"
49#include "load.h"
50#include "classfactory.h"
51#include "parametermanager.h"
52#include "paramkey.h"
53
54#ifdef __OOFEG
55 #include "oofeggraphiccontext.h"
56 #include "connectivitytable.h"
57#endif
58
59namespace oofem {
61
70
71Lattice3d_mt :: Lattice3d_mt(int n, Domain *aDomain, ElementMode em) :
72 LatticeTransportElement(n, aDomain, em)
73{
74 this->numberOfDofMans = 2;
76 minLength = 1.e-20;
77 dimension = 3.;
78 couplingFlag = 0;
79}
80
81double Lattice3d_mt :: giveLength()
82{
83 if ( this->geometryFlag == 0 ) {
85 }
86
87 return this->length;
88}
89
90void Lattice3d_mt :: giveCrackLengths(FloatArray &lengths)
91{
92 if ( this->geometryFlag == 0 ) {
94 }
95
96 lengths = this->crackLengths;
97}
98
99
100double
101Lattice3d_mt :: giveArea()
102{
103 if ( this->geometryFlag == 0 ) {
105 }
106
107 return this->area;
108}
109
110void
111Lattice3d_mt :: computeNmatrixAt(FloatMatrix &answer, const FloatArray &coords)
112{
113 this->computeNSubMatrixAt(answer, coords);
114}
115
116
117void
118Lattice3d_mt :: computeNSubMatrixAt(FloatMatrix &answer, const FloatArray &coords)
119{
120 double ksi = coords.at(1);
121 double n1 = ( 1. - ksi ) * 0.5;
122 double n2 = ( 1. + ksi ) * 0.5;
123
124 answer.resize(1, 2);
125 answer.zero();
126
127 answer.at(1, 1) = n1;
128 answer.at(1, 2) = n2;
129}
130
131
132void
133Lattice3d_mt :: computeGradientMatrixAt(FloatMatrix &answer, const FloatArray &lcoords)
134{
135 double l = giveLength();
136
137 answer.resize(1, 2);
138 answer.zero();
139 answer.at(1, 1) = -1.;
140 answer.at(1, 2) = 1.;
141
142 answer.times(1. / l);
143}
144
145
146void
147Lattice3d_mt :: updateInternalState(TimeStep *tStep)
148{
149 FloatArray f, r;
150 FloatMatrix n;
151 TransportMaterial *mat = static_cast< TransportMaterial * >( this->giveMaterial() );
152
153 // force updating ip values
154 for ( auto &iRule: integrationRulesArray ) {
155 for ( auto &gp: * iRule ) {
156 this->computeNmatrixAt(n, gp->giveNaturalCoordinates() );
157 this->computeVectorOf({ P_f }, VM_Total, tStep, r);
158 f.beProductOf(n, r);
159 mat->updateInternalState(f, gp, tStep);
160 }
161 }
162}
163
164
165void
166Lattice3d_mt :: computeGaussPoints()
167{
168 integrationRulesArray.resize(1);
169 integrationRulesArray [ 0 ] = std :: make_unique< GaussIntegrationRule >(1, this, 1, 2);
170 integrationRulesArray [ 0 ]->SetUpPointsOnLine(1, _3dMTLattice);
171}
172
173void
174Lattice3d_mt :: giveDofManDofIDMask(int inode, IntArray &answer) const
175{
176 answer = { P_f };
177}
178
179void
180Lattice3d_mt :: initializeFrom(InputRecord &ir, int priority)
181{
182 this->Element :: initializeFrom(ir, priority);
183 ParameterManager &ppm = this->giveDomain()->elementPPM;
185 PM_UPDATE_PARAMETER(dimension, ppm, ir, this->number, IPK_Lattice3d_mt_dim, priority) ;
186 PM_UPDATE_PARAMETER(area, ppm, ir, this->number, IPK_Lattice3d_mt_area, priority) ;
191
192}
193
194void Lattice3d_mt :: postInitialize()
195{
196 this->LatticeTransportElement :: postInitialize();
198 numberOfPolygonVertices = (int) (polygonCoords.giveSize() / 3.);
199 crackWidths.resizeWithValues(numberOfPolygonVertices);
200 this->computeGaussPoints();
201}
202
203
204void
205Lattice3d_mt :: computeFlow(FloatArray &answer, GaussPoint *gp, TimeStep *tStep)
206{
207 FloatArray r;
208 IntArray dofid;
209 double dV;
210 double length = giveLength();
211 double k;
212 answer.resize(1);
213
214 k = static_cast< TransportMaterial * >( this->giveMaterial() )->giveCharacteristicValue(Conductivity, gp, tStep);
215 dV = this->computeVolumeAround(gp);
216
217 this->giveElementDofIDMask(dofid);
218 this->computeVectorOf(dofid, VM_Total, tStep, r);
219 double flow;
220 double dP;
221
222 dP = r.at(2) - r.at(1);
223 flow = k * ( dP ) * dV / pow(length, 2.);
224 answer.at(1) = fabs(flow);
225
226 if ( !isActivated(tStep) ) {
227 answer.zero();
228 }
229}
230
231double
232Lattice3d_mt :: computeVolumeAround(GaussPoint *aGaussPoint)
233{
234 return this->giveArea() * this->giveLength();
235}
236
237
238void
239Lattice3d_mt :: computeGeometryProperties()
240{
241 //coordinates of the two nodes
242 Node *nodeA, *nodeB;
243 FloatArray coordsA(3), coordsB(3);
244
245 nodeA = this->giveNode(1);
246 nodeB = this->giveNode(2);
247
248 for ( int i = 0; i < 3; i++ ) {
249 coordsA.at(i + 1) = nodeA->giveCoordinate(i + 1);
250 coordsB.at(i + 1) = nodeB->giveCoordinate(i + 1);
251 }
252
253 //Construct an initial temporary local coordinate system
254 FloatArray s(3), t(3);
255
256 //Calculate normal vector
257 this->normal.resize(3);
258 for ( int i = 0; i < 3; i++ ) {
259 this->normal.at(i + 1) = coordsB.at(i + 1) - coordsA.at(i + 1);
260 }
261
262 // Compute midpoint
263 this->midPoint.resize(3);
264 for ( int i = 0; i < 3; i++ ) {
265 this->midPoint.at(i + 1) = 0.5 * ( coordsB.at(i + 1) + coordsA.at(i + 1) );
266 }
267
268 this->length = sqrt(pow(normal.at(1), 2.) + pow(normal.at(2), 2.) + pow(normal.at(3), 2.) );
269
270 if ( this->length < this->minLength ) {
271 this->length = this->minLength;
272 printf("Length could be zero. Might not be possible to create local coordinate system. Need to calculate area differently. Crack lengths need to approximated.\n");
273 //This will only work for cross-sections with three points at the moment.
275 this->geometryFlag = 1;
276 return;
277 }
278
279 for ( int i = 0; i < 3; i++ ) {
280 this->normal.at(i + 1) /= length;
281 }
282
284
285 this->geometryFlag = 1;
286
287 return;
288}
289
290void
291Lattice3d_mt :: computeSpecialCrossSectionProperties()
292{
293 if ( this->numberOfPolygonVertices != 3 ) {
294 return;
295 }
296
297 //Here, the area of the cross-section is calculated from the global vertex points without the generation of a local coordinate system.
298
299 FloatArray pointOne(3);
300 FloatArray pointTwo(3);
301 FloatArray crackPoint(3);
302
303 FloatArray l(3);
304 this->crackLengths.resize(3);
305 this->crackLengths.zero();
306
307 for ( int m = 0; m < numberOfPolygonVertices; m++ ) {
308 if ( m < numberOfPolygonVertices - 1 ) {
309 for ( int n = 0; n < 3; n++ ) {
310 pointOne.at(n + 1) = polygonCoords.at(3 * m + n + 1);
311 pointTwo.at(n + 1) = polygonCoords.at(3 * ( m + 1 ) + n + 1);
312 }
313 } else {
314 for ( int n = 0; n < 3; n++ ) {
315 pointOne.at(n + 1) = polygonCoords.at(3 * m + n + 1);
316 pointTwo.at(n + 1) = polygonCoords.at(n + 1);
317 }
318 }
319 l.at(m + 1) = sqrt(pow(pointOne.at(1) - pointTwo.at(1), 2.) + pow(pointOne.at(2) - pointTwo.at(2), 2.) + pow(pointOne.at(3) - pointTwo.at(3), 2.) );
320 }
321
322 double s = ( l.at(1) + l.at(2) + l.at(3) ) / 2.;
323 this->area = sqrt(s * ( s - l.at(1) ) * ( s - l.at(2) ) * ( s - l.at(3) ) );
324
325 printf("area = %e, length = %e\n", area, length);
326
327 if ( this->area < pow(this->minLength, 2.) ) {
328 this->area = pow(this->minLength, 2.);
329 printf("small area fixed in special calculation\n");
330 }
331
332 //Approximate crack lengths assuming
333 //an equilateral triangle of the same area as real triangle has.
334 for ( int m = 0; m < numberOfPolygonVertices; m++ ) {
335 this->crackLengths.at(m + 1) = pow(3, 0.25) * sqrt(this->area);
336 }
337
338 return;
339}
340
341
342
343// void
344// Lattice3d_mt :: computeHomogenisedInternalForcesVectorAt(FloatArray &answer, TimeStep *tStep, ValueModeType mode, FloatArray &unknowns)
345// {
346// IntArray dofids;
347// FloatArray tmp;
348
349// FloatMatrix s;
350// this->giveElementDofIDMask(dofids);
351
352// ///@todo Integrate and compute as a nonlinear problem instead of doing this tangent.
353// this->computeConductivityMatrix(s, Conductivity, tStep);
354// answer.beProductOf(s, unknowns);
355
356// this->computeInternalSourceRhsVectorAt(tmp, tStep, mode);
357// answer.subtract(tmp);
358
359// FloatMatrix bc_tangent;
360// this->computeBCMtrxAt(bc_tangent, tStep, VM_Total);
361// if ( bc_tangent.isNotEmpty() ) {
362// tmp.beProductOf(bc_tangent, unknowns);
363// answer.add(tmp);
364// }
365// }
366
367
368void
369Lattice3d_mt :: computeCrossSectionProperties() {
370 if ( this->numberOfPolygonVertices < 3 ) {
371 return;
372 }
373
374 //Construct two perpendicular axis so that n is normal to the plane which they create
375 //Check, if one of the components of the normal-direction is zero
376 FloatArray s(3), t(3);
377 if ( this->normal.at(1) == 0 ) {
378 s.at(1) = 0.;
379 s.at(2) = this->normal.at(3);
380 s.at(3) = -this->normal.at(2);
381 } else if ( this->normal.at(2) == 0 ) {
382 s.at(1) = this->normal.at(3);
383 s.at(2) = 0.;
384 s.at(3) = -this->normal.at(1);
385 } else {
386 s.at(1) = this->normal.at(2);
387 s.at(2) = -this->normal.at(1);
388 s.at(3) = 0.;
389 }
390
391 s.normalize();
392
393 t.beVectorProductOf(this->normal, s);
394 t.normalize();
395
396 //Set up rotation matrix
397 FloatMatrix lcs(3, 3);
398
399 for ( int i = 1; i <= 3; i++ ) {
400 lcs.at(1, i) = this->normal.at(i);
401 lcs.at(2, i) = s.at(i);
402 lcs.at(3, i) = t.at(i);
403 }
404
405
406 //Calculate the local coordinates of the polygon vertices
407 FloatArray help(3), test(3);
409 for ( int k = 0; k < numberOfPolygonVertices; k++ ) {
410 for ( int n = 0; n < 3; n++ ) {
411 help(n) = polygonCoords(3 * k + n);
412 }
413
414 test.beProductOf(lcs, help);
415 for ( int n = 0; n < 3; n++ ) {
416 lpc(3 * k + n) = test(n);
417 }
418 }
419
420 this->area = 0.;
421
422 for ( int k = 0; k < numberOfPolygonVertices; k++ ) {
423 if ( k < numberOfPolygonVertices - 1 ) {
424 this->area += lpc(3 * k + 1) * lpc(3 * ( k + 1 ) + 2) - lpc(3 * ( k + 1 ) + 1) * lpc(3 * k + 2);
425 } else { //Back to zero for n+1
426 this->area += lpc(3 * k + 1) * lpc(2) - lpc(1) * lpc(3 * k + 2);
427 }
428 }
429
430 this->area *= 0.5;
431
432 FloatArray tempCoords(3 * numberOfPolygonVertices);
433 if ( this->area < 0 ) { //Set area to a positive value and rearrange the coordinate entries
434 this->area *= -1.;
435
436 for ( int k = 0; k < numberOfPolygonVertices; k++ ) {
437 for ( int m = 0; m < 3; m++ ) {
438 tempCoords.at(3 * k + m + 1) = polygonCoords.at(3 * ( numberOfPolygonVertices - k - 1 ) + m + 1);
439 }
440 }
441
442 polygonCoords = tempCoords;
443
444 // Calculate again local co-ordinate system for different order
445 for ( int k = 0; k < numberOfPolygonVertices; k++ ) {
446 for ( int n = 0; n < 3; n++ ) {
447 help(n) = polygonCoords(3 * k + n);
448 }
449
450 test.beProductOf(lcs, help);
451 for ( int n = 0; n < 3; n++ ) {
452 lpc(3 * k + n) = test(n);
453 }
454 }
455 }
456
457 if ( this->area < pow(minLength, 2.) ) {
458 this->area = pow(minLength, 2.);
459 printf("Small area fixed.\n");
460 }
461
462 //Calculate centroids
463 centroid.resize(3);
464 for ( int k = 0; k < numberOfPolygonVertices; k++ ) {
465 if ( k < numberOfPolygonVertices - 1 ) {
466 centroid.at(2) += ( lpc(3 * k + 1) + lpc(3 * ( k + 1 ) + 1) ) * ( lpc(3 * k + 1) * lpc(3 * ( k + 1 ) + 2) - lpc(3 * ( k + 1 ) + 1) * lpc(3 * k + 2) );
467 centroid.at(3) += ( lpc(3 * k + 2) + lpc(3 * ( k + 1 ) + 2) ) * ( lpc(3 * k + 1) * lpc(3 * ( k + 1 ) + 2) - lpc(3 * ( k + 1 ) + 1) * lpc(3 * k + 2) );
468 } else { //Back to zero for n+1
469 centroid.at(2) += ( lpc(3 * k + 1) + lpc(1) ) * ( lpc(3 * k + 1) * lpc(2) - lpc(1) * lpc(3 * k + 2) );
470 centroid.at(3) += ( lpc(3 * k + 2) + lpc(2) ) * ( lpc(3 * k + 1) * lpc(2) - lpc(1) * lpc(3 * k + 2) );
471 }
472 }
473
474 centroid.times(1. / ( 6. * this->area ) );
475
476 centroid.at(1) = lpc.at(1); //The first component of all lpcs should be the same
477
478 //Shift coordinates to centroi
479 for ( int k = 0; k < numberOfPolygonVertices; k++ ) {
480 for ( int l = 0; l < 3; l++ ) {
481 lpc(3 * k + l) -= centroid(l);
482 }
483 }
484
485 //Compute second moments of area.
486
487 //This is for the temporary coordinate system
488 double Ixx = 0.;
489 double Iyy = 0.;
490 double Ixy = 0.;
491 double a;
492
493 for ( int k = 0; k < numberOfPolygonVertices; k++ ) {
494 if ( k < numberOfPolygonVertices - 1 ) {
495 a = lpc(3 * k + 1) * lpc(3 * ( k + 1 ) + 2) - lpc(3 * ( k + 1 ) + 1) * lpc(3 * k + 2);
496
497 Ixx += ( ( pow(lpc(3 * k + 2), 2.) + lpc(3 * k + 2) * lpc(3 * ( k + 1 ) + 2) + pow(lpc(3 * ( k + 1 ) + 2), 2.) ) * a ) / 12.;
498
499 Iyy += ( ( pow(lpc(3 * k + 1), 2.) + lpc(3 * k + 1) * lpc(3 * ( k + 1 ) + 1) + pow(lpc(3 * ( k + 1 ) + 1), 2.) ) * a ) / 12.;
500
501 Ixy += ( ( lpc(3 * k + 1) * lpc(3 * ( k + 1 ) + 2) + 2. * lpc(3 * k + 1) * lpc(3 * k + 2) +
502 2 * lpc(3 * ( k + 1 ) + 1) * lpc(3 * ( k + 1 ) + 2) + lpc(3 * ( k + 1 ) + 1) * lpc(3 * k + 2) ) * a ) / 24.;
503 } else { //Back to zero for n+1
504 a = lpc(3 * k + 1) * lpc(2) - lpc(1) * lpc(3 * k + 2);
505
506 Ixx += ( ( pow(lpc(3 * k + 2), 2.) + lpc(3 * k + 2) * lpc(2) + pow(lpc(2), 2.) ) * a ) / 12.;
507
508 Iyy += ( ( pow(lpc(3 * k + 1), 2.) + lpc(3 * k + 1) * lpc(1) + pow(lpc(1), 2.) ) * a ) / 12.;
509
510 Ixy += ( ( lpc(3 * k + 1) * lpc(2) + 2. * lpc(3 * k + 1) * lpc(3 * k + 2) +
511 2 * lpc(1) * lpc(2) + lpc(1) * lpc(3 * k + 2) ) * a ) / 24.;
512 }
513 }
514
515 //Compute main axis of the cross-section
516 double angleChange = 0.;
517 double sum = fabs(Ixx + Iyy);
518 double pi = 3.14159265;
519 if ( ( fabs(Ixx - Iyy) / sum > 1.e-6 ) && fabs(Ixy) / sum > 1.e-6 ) {
520 angleChange = 0.5 * atan(-2 * Ixy / ( Ixx - Iyy ) );
521 } else if ( ( fabs(Ixx - Iyy) / sum < 1.e-6 ) && fabs(Ixy) / sum > 1.e-6 ) {
522 angleChange = pi / 4.;
523 }
524
525 if ( Iyy > Ixx ) {
526 angleChange = angleChange + pi / 2.;
527 }
528
529 //Moment of inertias saved in the element
530
531 this->I1 = ( Ixx + Iyy ) / 2. + sqrt(pow( ( Ixx - Iyy ) / 2., 2. ) + pow(Ixy, 2.) );
532 this->I2 = ( Ixx + Iyy ) / 2. - sqrt(pow( ( Ixx - Iyy ) / 2., 2. ) + pow(Ixy, 2.) );
533
534 this->Ip = I1 + I2;
535
536 //Rotation around normal axis by angleChange
537
538 FloatMatrix rotationChange(3, 3);
539 rotationChange.zero();
540
541 rotationChange.at(1, 1) = 1.;
542 rotationChange.at(2, 2) = cos(angleChange);
543 rotationChange.at(2, 3) = -sin(angleChange);
544
545 rotationChange.at(3, 2) = sin(angleChange);
546 rotationChange.at(3, 3) = cos(angleChange);
547
548 this->localCoordinateSystem.beProductOf(rotationChange, lcs);
549
550 //Calculate the polygon vertices in the new coordinate system
551 for ( int k = 0; k < numberOfPolygonVertices; k++ ) {
552 for ( int n = 0; n < 3; n++ ) {
553 help(n) = polygonCoords(3 * k + n);
554 }
555
556 test.beProductOf(this->localCoordinateSystem, help);
557 for ( int n = 0; n < 3; n++ ) {
558 lpc(3 * k + n) = test(n);
559 }
560 }
561
562 //Calculate centroid again in local coordinate system
563 centroid.zero();
564 for ( int k = 0; k < numberOfPolygonVertices; k++ ) {
565 if ( k < numberOfPolygonVertices - 1 ) {
566 centroid.at(2) += ( lpc(3 * k + 1) + lpc(3 * ( k + 1 ) + 1) ) * ( lpc(3 * k + 1) * lpc(3 * ( k + 1 ) + 2) - lpc(3 * ( k + 1 ) + 1) * lpc(3 * k + 2) );
567 centroid.at(3) += ( lpc(3 * k + 2) + lpc(3 * ( k + 1 ) + 2) ) * ( lpc(3 * k + 1) * lpc(3 * ( k + 1 ) + 2) - lpc(3 * ( k + 1 ) + 1) * lpc(3 * k + 2) );
568 } else { //Back to zero for n+1
569 centroid.at(2) += ( lpc(3 * k + 1) + lpc(1) ) * ( lpc(3 * k + 1) * lpc(2) - lpc(1) * lpc(3 * k + 2) );
570 centroid.at(3) += ( lpc(3 * k + 2) + lpc(2) ) * ( lpc(3 * k + 1) * lpc(2) - lpc(1) * lpc(3 * k + 2) );
571 }
572 }
573
574 centroid.times(1. / ( 6. * area ) );
575
576 centroid.at(1) = lpc.at(1); //The first component of all lpcs should be the same
577
578 /*Express midpoint in local coordinate system.
579 * Note that this point does not have to lie in the same plane as the cross-section.
580 * However, we do not use the first componnent of the coordinate */
581 FloatArray midPointLocal(3);
582 midPointLocal.beProductOf(this->localCoordinateSystem, midPoint);
583
584 //eccentricities stored in the element
585 this->eccS = centroid.at(2) - midPointLocal.at(2);
586 this->eccT = centroid.at(3) - midPointLocal.at(3);
587
588 FloatMatrix transposeLCS;
589 transposeLCS.beTranspositionOf(this->localCoordinateSystem);
590
591 globalCentroid.beProductOf(transposeLCS, centroid);
592
594 double crackPointOne, crackPointTwo;
595 for ( int m = 0; m < numberOfPolygonVertices; m++ ) {
596 if ( m < numberOfPolygonVertices - 1 ) {
597 crackPointOne = ( lpc(3 * ( m + 1 ) + 1) + lpc(3 * ( m ) + 1) ) / 2.;
598 crackPointTwo = ( lpc(3 * ( m + 1 ) + 2) + lpc(3 * ( m ) + 2) ) / 2.;
599 } else {
600 crackPointOne = ( lpc(1) + lpc(3 * ( m ) + 1) ) / 2.;
601 crackPointTwo = ( lpc(2) + lpc(3 * ( m ) + 2) ) / 2.;
602 }
603
604 crackLengths.at(m + 1) = sqrt(pow(crackPointOne - midPointLocal.at(2), 2.) + pow(crackPointTwo - midPointLocal.at(3), 2.) );
605 }
606
607 return;
608}
609
610void
611Lattice3d_mt :: computeConductivityMatrix(FloatMatrix &answer, MatResponseMode rmode, TimeStep *tStep)
612{
613 GaussPoint *gp = integrationRulesArray [ 0 ]->getIntegrationPoint(0);
614
615 answer.resize(2, 2);
616 answer.zero();
617 answer.at(1, 1) = 1.;
618 answer.at(1, 2) = -1.;
619 answer.at(2, 1) = -1;
620 answer.at(2, 2) = 1.;
621
622 double length = giveLength();
623 double k = static_cast< TransportMaterial * >( this->giveMaterial() )->giveCharacteristicValue(Conductivity, gp, tStep);
624 double dV = this->computeVolumeAround(gp);
625 double temp = k * dV / pow(length, 2.);
626 answer.times(temp);
627
628 return;
629}
630
631
632void
633Lattice3d_mt :: computeCapacityMatrix(FloatMatrix &answer, TimeStep *tStep)
634{
635 GaussPoint *gp = integrationRulesArray [ 0 ]->getIntegrationPoint(0);
636 answer.resize(2, 2);
637 answer.zero();
638 answer.at(1, 1) = 2.;
639 answer.at(1, 2) = 1.;
640 answer.at(2, 1) = 1.;
641 answer.at(2, 2) = 2.;
642 double c = static_cast< TransportMaterial * >( this->giveMaterial() )->giveCharacteristicValue(Capacity, gp, tStep);
643 double dV = this->computeVolumeAround(gp) / ( 6.0 * this->dimension );
644 answer.times(c * dV);
645}
646
647
648
649void
650Lattice3d_mt :: computeInternalSourceRhsVectorAt(FloatArray &answer, TimeStep *atTime, ValueModeType mode)
651{
652 int i, j, n, nLoads;
653 double dV;
654 bcGeomType ltype;
655 Load *load;
656 GaussPoint *gp = integrationRulesArray [ 0 ]->getIntegrationPoint(0);
657
658 FloatArray deltaX(3), normalVector(3);
659 FloatArray val, helpLoadVector, globalIPcoords;
660 FloatMatrix nm;
661 answer.resize(0);
662
663 FloatArray gravityHelp(2);
664
665 nLoads = this->giveBodyLoadArray()->giveSize();
666 for ( i = 1; i <= nLoads; i++ ) {
667 n = bodyLoadArray.at(i);
668 load = ( Load * ) domain->giveLoad(n);
669 ltype = load->giveBCGeoType();
670
671 if ( ltype == GravityPressureBGT ) {
672 //Compute change of coordinates
673
674 if ( this->geometryFlag == 0 ) {
676 }
677
678 deltaX.at(1) = this->normal.at(1) * this->length;
679 deltaX.at(2) = this->normal.at(2) * this->length;
680 deltaX.at(3) = this->normal.at(3) * this->length;
681
682 gravityHelp.at(1) = 1.;
683 gravityHelp.at(2) = -1.;
684
685 dV = this->computeVolumeAround(gp);
686 load->computeValueAt(val, atTime, deltaX, mode);
687
688 double k = static_cast< TransportMaterial * >( this->giveMaterial() )->giveCharacteristicValue(Conductivity, gp, atTime);
689
690 double helpFactor = val.at(1) * k * dV;
691
692 helpFactor /= pow(this->giveLength(), 2.);
693 gravityHelp.times(helpFactor);
694
695 if ( helpLoadVector.isEmpty() ) {
696 helpLoadVector.resize(gravityHelp.giveSize() );
697 }
698
699 for ( j = 1; j <= gravityHelp.giveSize(); j++ ) {
700 helpLoadVector.at(j) += gravityHelp.at(j);
701 }
702 }
703
704
705
706 answer.add(helpLoadVector);
707 }
708
709 return;
710}
711
712int
713Lattice3d_mt :: computeGlobalCoordinates(FloatArray &answer, const FloatArray &lcoords)
714{
715 if ( geometryFlag == 0 ) {
717 }
718
719 answer.resize(3);
720 answer = midPoint;
721
722 return 1;
723}
724
725bool
726Lattice3d_mt :: computeLocalCoordinates(FloatArray &answer, const FloatArray &coords)
727{
728 answer.resize(1);
729 answer.at(1) = 0.;
730
731 return 1;
732}
733
734
735#ifdef __OOFEG
736
737void
738Lattice3d_mt :: drawYourself(oofegGraphicContext &gc, TimeStep *tStep)
739{
740 OGC_PlotModeType mode = gc.giveIntVarPlotMode();
741
742 if ( mode == OGC_rawGeometry ) {
743 this->drawRawGeometry(gc, tStep);
744 this->drawRawCrossSections(gc, tStep);
745 } else {
746 OOFEM_ERROR("drawYourself : unsupported mode");
747 }
748}
749
750
751void Lattice3d_mt :: drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep)
752{
753 GraphicObj *go;
754
755 WCRec p [ 2 ]; /* points */
756 if ( !gc.testElementGraphicActivity(this) ) {
757 return;
758 }
759
760 EASValsSetLineWidth(OOFEG_RAW_GEOMETRY_WIDTH);
761 EASValsSetColor(gc.getElementColor() );
762 EASValsSetLayer(OOFEG_RAW_GEOMETRY_LAYER);
763
764 p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveCoordinate(1);
765 p [ 0 ].y = ( FPNum ) this->giveNode(1)->giveCoordinate(2);
766 p [ 0 ].z = ( FPNum ) this->giveNode(1)->giveCoordinate(3);
767 p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveCoordinate(1);
768 p [ 1 ].y = ( FPNum ) this->giveNode(2)->giveCoordinate(2);
769 p [ 1 ].z = ( FPNum ) this->giveNode(2)->giveCoordinate(3);
770
771 go = CreateLine3D(p);
772 EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
773 EGAttachObject(go, ( EObjectP ) this);
774 EMAddGraphicsToModel(ESIModel(), go);
775}
776
777
778void Lattice3d_mt :: drawRawCrossSections(oofegGraphicContext &gc, TimeStep *tStep)
779{
780 GraphicObj *go;
781
782 // if (!go) { // create new one
783 //Create as many points as we have polygon vertices
784 WCRec p [ numberOfPolygonVertices ]; /* poin */
785
786 if ( !gc.testElementGraphicActivity(this) ) {
787 return;
788 }
789
790 EASValsSetLineWidth(OOFEG_RAW_GEOMETRY_WIDTH);
791 EASValsSetColor(gc.getElementColor() );
792 EASValsSetLayer(OOFEG_RAW_CROSSSECTION_LAYER);
793
794 for ( int i = 0; i < numberOfPolygonVertices; i++ ) {
795 p [ i ].x = ( FPNum ) polygonCoords(3 * i);
796 p [ i ].y = ( FPNum ) polygonCoords(3 * i + 1);
797 p [ i ].z = ( FPNum ) polygonCoords(3 * i + 2);
798 }
799
800 WCRec pTemp [ 2 ]; /* points */
801 for ( int i = 0; i < numberOfPolygonVertices; i++ ) {
802 if ( i < numberOfPolygonVertices - 1 ) {
803 pTemp [ 0 ] = p [ i ];
804 pTemp [ 1 ] = p [ i + 1 ];
805 } else {
806 pTemp [ 0 ] = p [ i ];
807 pTemp [ 1 ] = p [ 0 ];
808 }
809
810 go = CreateLine3D(pTemp);
811 EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
812 EGAttachObject(go, ( EObjectP ) this);
813 EMAddGraphicsToModel(ESIModel(), go);
814 }
815}
816#endif
817} // end namespace oofem
double length(const Vector &a)
Definition CSG.h:88
#define REGISTER_Element(class)
double giveCoordinate(int i) const
Definition dofmanager.h:383
virtual void giveElementDofIDMask(IntArray &answer) const
Definition element.h:510
Node * giveNode(int i) const
Definition element.h:629
virtual bool isActivated(TimeStep *tStep)
Definition element.C:838
IntArray * giveBodyLoadArray()
Returns array containing load numbers of loads acting on element.
Definition element.C:411
int numberOfDofMans
Number of dofmanagers.
Definition element.h:136
void computeVectorOf(ValueModeType u, TimeStep *tStep, FloatArray &answer)
Definition element.C:103
std::vector< std ::unique_ptr< IntegrationRule > > integrationRulesArray
Definition element.h:157
int numberOfGaussPoints
Definition element.h:175
virtual double giveCharacteristicValue(CharType type, TimeStep *tStep)
Definition element.C:678
IntArray bodyLoadArray
Definition element.h:147
Domain * giveDomain() const
Definition femcmpnn.h:97
Domain * domain
Link to domain object, useful for communicating with other FEM components.
Definition femcmpnn.h:79
int number
Component number.
Definition femcmpnn.h:77
void resize(Index s)
Definition floatarray.C:94
double & at(Index i)
Definition floatarray.h:202
Index giveSize() const
Returns the size of receiver.
Definition floatarray.h:261
void zero()
Zeroes all coefficients of receiver.
Definition floatarray.C:683
void beProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Definition floatarray.C:689
bool isEmpty() const
Returns true if receiver is empty.
Definition floatarray.h:265
void beVectorProductOf(const FloatArray &v1, const FloatArray &v2)
Definition floatarray.C:476
void add(const FloatArray &src)
Definition floatarray.C:218
void times(double s)
Definition floatarray.C:834
void times(double f)
void resize(Index rows, Index cols)
Definition floatmatrix.C:79
void beTranspositionOf(const FloatMatrix &src)
void zero()
Zeroes all coefficient of receiver.
double at(std::size_t i, std::size_t j) const
virtual bcGeomType giveBCGeoType() const
double computeVolumeAround(GaussPoint *) override
static ParamKey IPK_Lattice3d_mt_polycoords
double giveArea() override
void computeSpecialCrossSectionProperties()
void computeCrossSectionProperties()
FloatArray polygonCoords
static ParamKey IPK_Lattice3d_mt_couplingflag
double giveLength() override
FloatMatrix localCoordinateSystem
static ParamKey IPK_Lattice3d_mt_area
static ParamKey IPK_Lattice3d_mt_dim
FloatArray globalCentroid
void drawRawCrossSections(oofegGraphicContext &gc, TimeStep *tStep)
static ParamKey IPK_Lattice3d_mt_mlength
FloatArray crackWidths
void drawRawGeometry(oofegGraphicContext &, TimeStep *tStep) override
void computeNmatrixAt(FloatMatrix &n, const FloatArray &) override
FloatArray crackLengths
void computeGeometryProperties()
void computeNSubMatrixAt(FloatMatrix &n, const FloatArray &)
static ParamKey IPK_Lattice3d_mt_crackwidths
static ParamKey IPK_Lattice3d_mt_ranarea
void computeGaussPoints() override
static ParamKey IPK_Lattice3d_mt_couplingnumber
LatticeTransportElement(int, Domain *, ElementMode)
virtual void computeValueAt(FloatArray &answer, TimeStep *tStep, const FloatArray &coords, ValueModeType mode)=0
Material * giveMaterial() override
virtual void updateInternalState(const FloatArray &state, GaussPoint *gp, TimeStep *tStep)
#define OOFEM_ERROR(...)
Definition error.h:79
bcGeomType
Type representing the geometric character of loading.
Definition bcgeomtype.h:40
@ GravityPressureBGT
Pressure due to distributed body load.
Definition bcgeomtype.h:47
double sum(const FloatArray &x)
#define OOFEG_RAW_CROSSSECTION_LAYER
oofem::oofegGraphicContext gc[OOFEG_LAST_LAYER]
#define OOFEG_RAW_GEOMETRY_WIDTH
#define OOFEG_RAW_GEOMETRY_LAYER
#define PM_UPDATE_PARAMETER(_val, _pm, _ir, _componentnum, _paramkey, _prio)

This page is part of the OOFEM-3.0 documentation. Copyright Copyright (C) 1994-2025 Borek Patzak Bořek Patzák
Project e-mail: oofem@fsv.cvut.cz
Generated at for OOFEM by doxygen 1.15.0 written by Dimitri van Heesch, © 1997-2011