1 /*! 2 * 3 * JSBeams : Numerical solution of beam structures in JavaScript 4 * version 0.6 (2012-04-20) 5 * 6 * Copyright (C) 2011 - 2012 Jan Stransky 7 * 8 * Czech Technical University, Faculty of Civil Engineering, 9 * Department of Structural Mechanics, 166 29 Prague, Czech Republic 10 * 11 * 12 * JSBeams is free software: you can redistribute it and/or modify 13 * it under the terms of the GNU General Public License as published by 14 * the Free Software Foundation, either version 3 of the License, or 15 * (at your option) any later version. 16 * 17 * JSBeams is distributed in the hope that it will be useful, 18 * but WITHOUT ANY WARRANTY; without even the implied warranty of 19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 20 * GNU General Public License for more details. 21 * 22 * You should have received a copy of the GNU General Public License 23 * along with this program. If not, see <http://www.gnu.org/licenses/>. 24 */ 25 26 /** 27 * @fileOverview JavaScript implementation of numerical solution (direct stiffness method) of beam structures. Modeled structures consist of nodes and beams. Nodes have only geometrical meaning and contains only coordinates, dofs and their values are stored in beams. Different static schemes and methods (clamped vs. pin ends, rigid arms, master x slave nodes etc.) can be modeled by special using of beam DOFs. Requires <a href="http://mech.fsv.cvut.cz/~stransky/software/jsmatrix/">jsmatrix.js</a> loaded before loading this file. For more information see <a href="http://mech.fsv.cvut.cz/~stransky/software/jsbeams/">project homepage</a>. 28 <br /><br/>JSBeams is a free software distributed under <a href='http://www.gnu.org/licenses/gpl.html'>GNU GPL license</a>. 29 * @author <a href="http://mech.fsv.cvut.cz/~stransky">Jan Stránský</a> 30 * @version 0.6 (2012-04-20) 31 */ 32 33 // constants for easier manipulation 34 /**constant*/ var JSB2D_CLAMPED_CLAMPED = 0; // should be 0 35 /**constant*/ var JSB2D_CC = 0; 36 /**constant*/ var JSB2D_CLAMPED_HINGE = 1; 37 /**constant*/ var JSB2D_CH = 1; 38 /**constant*/ var JSB2D_HINGE_CLAMPED = 2; 39 /**constant*/ var JSB2D_HC = 2; 40 /**constant*/ var JSB2D_HINGE_HINGE = 3; 41 /**constant*/ var JSB2D_HH = 3; 42 /**constant*/ var JSB2D_TRUSS = 3; 43 44 45 /** Returns elements of a that are not in b (Matlab-like function), ret = a - b 46 * @param {Array} a 47 * @param {Array} b 48 * @returns {Array} array of elements of a that are not in b 49 * @example a = [1,2,3,4,5,6]; 50 * b = [2,5,1,99,54]; 51 * c = setdiff(a,b) // c = [3,4,6] 52 */ 53 setdiff = function(A,B) { 54 var ret = []; 55 for (a=0; a<A.length; a++) { 56 var noMatch = 1; 57 for (var b=0; b<B.length; b++) { 58 if (A[a] == B[b]) { noMatch = 0; break; } 59 } 60 if (noMatch) { ret.push(A[a]); } 61 } 62 return ret; 63 } 64 65 66 /** Node implementation. The meaning of current implementation is only geometric, because DOFs itself are stored in beams (for both computational and graphic postprocessing part) 67 * @class Represents nodes 68 * @param {float} x x coordinate of node 69 * @param {float} y y coordinate of node 70 * @param {float} z z coordinate of node 71 * @property {float} x x coordinate 72 * @property {float} y y coordinate 73 * @property {float} z z coordinate 74 */ 75 Node = function(x,y,z) { 76 this.x = x || 0.; 77 this.y = y || 0.; 78 this.z = z || 0.; 79 } 80 /**#@+ 81 * @function 82 * @memberOf Node# 83 */ 84 85 /** String representation 86 * @returns {String} string representation 87 */ 88 Node.prototype.toString = function() { 89 return "Node ( "+this.x+", "+this.y+", "+this.z+" )"; 90 }, 91 92 /** Sets new position to receiver 93 * @param {float} x x coordinate of new position 94 * @param {float} y y coordinate of new position 95 * @param {float} z z coordinate of new position 96 */ 97 Node.prototype.setXYZ = function(x,y,z) { 98 this.x = x; 99 this.y = y; 100 this.z = z; 101 }, 102 103 /** Sets new position to receiver 104 * @param {float} x x coordinate of new position 105 * @param {float} z z coordinate of new position 106 */ 107 Node.prototype.setXZ = function(x,z) { 108 this.x = x; 109 this.z = z; 110 }, 111 112 /** Sets new position to receiver 113 * @param {float} y y coordinate of new position 114 * @param {float} z z coordinate of new position 115 */ 116 Node.prototype.setYZ = function(y,z) { 117 this.y = y; 118 this.z = z; 119 }, 120 121 /** Set new x coordinate to receiver 122 * @param x x coordinate of new position 123 */ 124 Node.prototype.setX = function(x) { 125 this.x = x; 126 }, 127 128 /** Set new y coordinate to receiver 129 * @param y y coordinate of new position 130 */ 131 Node.prototype.setY = function(y) { 132 this.y = y; 133 }, 134 135 /** Set new z coordinate to receiver 136 * @param z z coordinate of new position 137 */ 138 Node.prototype.setZ = function(z) { 139 this.z = z; 140 }, 141 142 /** Change position of receiver by given values 143 * @param {float} dx length of translation in x direction 144 * @param {float} dy length of translation in y direction 145 * @param {float} dz length of translation in z direction 146 */ 147 Node.prototype.translate = function(dx,dy,dz) { 148 this.x += dx; 149 this.y += dy; 150 this.z += dz; 151 }, 152 153 /** Change position of receiver by given values 154 * @param {float} dx length of translation in x direction 155 * @param {float} dz length of translation in z direction 156 */ 157 Node.prototype.translateXZ = function(dx,dz) { 158 this.translate(dx,0.,dz); 159 }, 160 161 /** Change position of receiver by given values 162 * @param {float} dy length of translation in y direction 163 * @param {float} dz length of translation in z direction 164 */ 165 Node.prototype.translateYZ = function(dy,dz) { 166 this.translate(0.,dy,dz); 167 }, 168 169 /** Change horizontal position of receiver by given values 170 * @param {float} dx length of translation in x direction 171 */ 172 Node.prototype.translateX = function(dx) { 173 this.translate(dx,0.,0.); 174 }, 175 176 /** Change horizontal position of receiver by given values 177 * @param {float} dy length of translation in y direction 178 */ 179 Node.prototype.translateY = function(dy) { 180 this.translate(0.,dy,0.); 181 }, 182 183 /** Change vertical position of receiver by given values 184 * @param {float} dz length of translation in z direction 185 */ 186 Node.prototype.translateZ = function(dz) { 187 this.translate(0.,0.,dz); 188 }, 189 190 /** Returns difference between receiver and given node 191 * @param {Node} n node to be substracted 192 * @returns {Node} difference this-node 193 */ 194 Node.prototype.sub = function(n) { 195 return new Node(this.x-n.x,this.y-n.y,this.z-n.z); 196 }, 197 198 /** Returns sum of receiver and given node 199 * @param {Node} n node to be substracted 200 * @returns {Node} difference this+node 201 */ 202 Node.prototype.add = function(n) { 203 return new Node(this.x+n.x,this.y+x.y,this.z+n.z); 204 }, 205 206 /** Returns receiver multiplied by given number 207 * @param {float} f multiplicator 208 * @returns {Node} product of multiplication f*this 209 */ 210 Node.prototype.mulf = function(f) { 211 return new Node(f*this.x,f*this.y,f*this.z); 212 }, 213 214 /** Returns cross product of receiver and given node 215 * @param {Node} n node to be crossed 216 * @returns {Node} cross product this x node 217 */ 218 Node.prototype.cross = function(n) { 219 return new Node(this.y*n.z-this.z*n.y, this.z*n.x-this.x*n.z, this.x*n.y-this.y*n.x); 220 }, 221 222 /** Returns dot product of receiver and given node 223 * @param {Node} n node to be crossed 224 * @returns {float} dot product this . node 225 */ 226 Node.prototype.dot = function(n) { 227 return this.x*n.x + this.y*n.y + this.z*n.z; 228 }, 229 230 /** Returns squared (Euclidean) norm of receiver this.dot(this) 231 * @returns {float} squared norm of receiver 232 */ 233 Node.prototype.squaredNorm = function() { 234 return this.dot(this); 235 }, 236 237 /** Returns (Euclidean) norm of receiver sqrt(this.dot(this)) 238 * @returns {float} norm of receiver 239 */ 240 Node.prototype.norm = function() { 241 return Math.sqrt(this.dot(this)); 242 }, 243 /**#@-*/ 244 245 /** Constructor, see {@link Node} for input parameters description 246 * @returns {Node} new Node object 247 */ 248 Node.create = function(x,y,z) { 249 var ret = new Node(x,y,z); 250 return ret; 251 } 252 253 254 255 256 /** 2d beam implementation 257 * @class Represents structural 2d beam with 3 degrees of freedom in each node 258 * @param {Node} node1 1st node 259 * @param {Node} node2 2nd node 260 * @param {Object} [params={}] object of parameters in format {param1:val1,param2:val2}, i.e. beam=Beam2d(n1,n2,{dofs:[1,2,3,4],ea:4}); 261 * @param {int} [params.type=JSB2D_CC] type of beam. JSB2D_CC for both ends clamped, JSB2D_HC for hinge-clamped, JSB2D_CH for clamped-hinge, JSB2D_HH for hinge-hinge (truss) 262 * @param {[ints]} [params.dofs=[0,1,2,3,4,5]] degrees of freedom on both ends. number (6,5,4) should correscpond to this.type 263 * @param {float} [params.ei=1.] bending stiffness [Nm2] 264 * @param {float} [params.ea=1.] normal stiffness [N] 265 * @param {float} [params.fzloc=0.] constant continuous load parallel to local z coordinate [N/m] 266 * @param {float} [params.fxloc=0.] constant continuous load parallel to local x coordinate [N/m] 267 * @param {float} [params.mu=1.] length mass [kg/m] 268 * @param {float} [params.dT=0.] positive temperature change (=warming) of beam center line [K] 269 * @param {float} [params.dTh=0.] (temperature of bottom fibers - temperature of top fibers)/height ((Tb-Tt)/h) [K/m] 270 * @param {float} [params.alpha=12e-6] thermal dillatation coefficient [K-1] 271 * @property {Node} node1 1st node 272 * @property {Node} node2 2nd node 273 * @property {[ints]} dofs degrees of freedom on both ends 274 * @property {float} ei bending stiffness [Nm2] 275 * @property {float} ea normal stiffness [N] 276 * @property {float} fzloc constant continuous load parallel to local z coordinate [N/m] 277 * @property {float} fxloc constant continuous load parallel to local x coordinate [N/m] 278 * @property {float} mu length mass [kg/m] 279 * @property {float} dT positive temperature change (=warming) of beam center line [K] 280 * @property {float} dTh (temperature of bottom fibers - temperature of top fibers)/height ((Tb-Tt)/h) [K/m] 281 * @property {float} alpha thermal dillatation coefficient [K-1] 282 * @property {Matrix} kloc stiffness matrix in local coordinate system, computed by setStiffAndTrsfMats method 283 * @property {Matrix} t transformation matrix (from global to local), computed by setStiffAndTrsfMats method 284 * @property {Matrix} tt transposed transformation matrix (from local to global), computed by setStiffAndTrsfMats method 285 * @property {Matrix} mloc mass matrix in local coordinate system, computed by setMassMats method 286 * @property {Matrix} m mass matrix in global coordinate system, computed by setMassMats method 287 * @property {Vector} load vector of end forces 288 * @property {Vector} dspl vector of end displacements 289 * @property {[Vectors]} eigShapes array of eigenshapes 290 * @property {float} len length of the beam 291 * @property {float} sin sin of node1->node2 vector 292 * @property {float} cos cos of node1->node2 vector 293 */ 294 // example? 295 Beam2d = function(node1,node2,params) { 296 var params = params==undefined? {} : params; 297 this.node1 = node1 || null; 298 this.node2 = node2 || null; 299 this.type = params.type || JSB2D_CC; 300 this.dofs = params.dofs || [0,1,2,3,4,5]; 301 this.ei = params.ei || 1.; 302 this.ea = params.ea || 1.; 303 this.fzloc = params.fzloc || 0.; 304 this.fxloc = params.fxloc || 0.; 305 this.mu = params.mu || 1.; 306 this.dT = params.dT || 0.; 307 this.dTh = params.dTh || 0.; 308 this.alpha = params.alpha || 12e-6; 309 this.kloc = null; 310 this.a = null; 311 this.b = null; 312 this.kaa = null; 313 this.kab = null; 314 this.kba = null; 315 this.kbb = null; 316 this.t = new Matrix(); 317 this.tt = new Matrix(); 318 this.mloc = new Matrix(); 319 this.m = new Matrix(); 320 this.load = Vector.Zeros(6); 321 this.loadCC = Vector.Zeros(6); 322 this.dspl = Vector.Zeros(6); 323 this.eigShapes = []; 324 this.len = null; 325 this.sin = null; 326 this.cos = null; 327 328 } 329 330 /**#@+ 331 * @function 332 * @memberOf Beam2d# 333 */ 334 335 /** Sets geometric parameters of received according to node positions 336 */ 337 Beam2d.prototype.setGeom = function() { 338 var dx = (this.node2.x - this.node1.x); 339 var dz = (this.node2.z - this.node1.z); 340 this.len = Math.sqrt(dx*dx + dz*dz); 341 this.cos = dx/this.len; 342 this.sin = dz/this.len; 343 }, 344 345 Beam2d.prototype.hasBodyLoad = function() { 346 return this.fzloc || this.fxloc || this.dT || this.dTh; 347 }, 348 349 /** Sets stiffness matrices in local and global coordinate system according to stored material and crossection parameters and transformation matrix and its transposition according to stored geometry 350 */ 351 Beam2d.prototype.setStiffAndTrsfMats = function() { 352 var l = this.len; 353 var eal = this.ea/l; 354 var eil = this.ei/l; 355 var eil2 = eil/l; 356 var eil3 = eil2/l; 357 var s = this.sin; 358 var c = this.cos; 359 // T - from gloal to local, TT - from local to global 360 switch (this.type) { 361 case JSB2D_HC: 362 this.t.fromArray([ 363 [ c, s, 0, 0, 0], 364 [-s, c, 0, 0, 0], 365 [ 0, 0, c, s, 0], 366 [ 0, 0,-s, c, 0], 367 [ 0, 0, 0, 0, 1] 368 ]); 369 this.tt.fromArray([ 370 [ c,-s, 0, 0, 0], 371 [ s, c, 0, 0, 0], 372 [ 0, 0, c,-s, 0], 373 [ 0, 0, s, c, 0], 374 [ 0, 0, 0, 0, 1] 375 ]); 376 break; 377 case JSB2D_CH: 378 this.t.fromArray([ 379 [ c, s, 0, 0, 0], 380 [-s, c, 0, 0, 0], 381 [ 0, 0, 1, 0, 0], 382 [ 0, 0, 0, c, s], 383 [ 0, 0, 0,-s, c] 384 ]); 385 this.tt.fromArray([ 386 [ c,-s, 0, 0, 0], 387 [ s, c, 0, 0, 0], 388 [ 0, 0, 1, 0, 0], 389 [ 0, 0, 0, c,-s], 390 [ 0, 0, 0, s, c] 391 ]); 392 break; 393 case JSB2D_HH: 394 this.t.fromArray([ 395 [ c, s, 0, 0], 396 [-s, c, 0, 0], 397 [ 0, 0, c, s], 398 [ 0, 0,-s, c] 399 ]); 400 this.tt.fromArray([ 401 [ c,-s, 0, 0], 402 [ s, c, 0, 0], 403 [ 0, 0, c,-s], 404 [ 0, 0, s, c] 405 ]); 406 break; 407 default: // c-c 408 this.t.fromArray([ 409 [ c, s, 0, 0, 0, 0], 410 [-s, c, 0, 0, 0, 0], 411 [ 0, 0, 1, 0, 0, 0], 412 [ 0, 0, 0, c, s, 0], 413 [ 0, 0, 0,-s, c, 0], 414 [ 0, 0, 0, 0, 0, 1] 415 ]); 416 this.tt.fromArray([ 417 [ c,-s, 0, 0, 0, 0], 418 [ s, c, 0, 0, 0, 0], 419 [ 0, 0, 1, 0, 0, 0], 420 [ 0, 0, 0, c,-s, 0], 421 [ 0, 0, 0, s, c, 0], 422 [ 0, 0, 0, 0, 0, 1] 423 ]); 424 } 425 var k = Matrix.create([ 426 [ eal, 0 , 0 , -eal, 0 , 0 ], 427 [ 0 , 12*eil3, -6*eil2, 0 , -12*eil3, -6*eil2], 428 [ 0 , - 6*eil2, 4*eil , 0 , 6*eil2, 2*eil ], 429 [-eal, 0 , 0 , eal, 0 , 0 ], 430 [ 0 , -12*eil3, 6*eil2, 0 , 12*eil3, 6*eil2], 431 [ 0 , - 6*eil2, 2*eil , 0 , 6*eil2, 4*eil ]]); 432 if (this.type) { // type is not c-c 433 switch (this.type) { 434 case JSB2D_HC: // hinged-clamped 435 this.a = [0,1,3,4,5]; 436 this.b = [2]; 437 break; 438 case JSB2D_CH: // clamped hinge 439 this.a = [0,1,2,3,4]; 440 this.b = [5]; 441 break; 442 case JSB2D_HH: // truss 443 this.a = [0,1,3,4]; 444 this.b = [2,5]; 445 break; 446 } 447 this.kaa = k.getSubMatrix(this.a,this.a); 448 this.kab = k.getSubMatrix(this.a,this.b); 449 this.kba = k.getSubMatrix(this.b,this.a); 450 this.kbbinv = k.getSubMatrix(this.b,this.b).inv(); 451 this.kloc = this.kaa.sub(this.kab.mulm(this.kbbinv).mulm(this.kba)) 452 } 453 else { this.kloc = k; } 454 }, 455 456 /** Computes stiffness matrix in global coordinates system 457 * @returns Matrix stiffness matrix and corresponding dofs 458 */ 459 Beam2d.prototype.computeGlobStiffMatrix = function() { 460 return this.tt.mulm(this.kloc).mulm(this.t); 461 }, 462 463 /** Sets mass matrices in local and gloal coordinate system from stored material and geometric parameters 464 */ 465 Beam2d.prototype.setMassMats = function() { 466 var l = this.len; 467 var mu = this.mu; 468 var a = mu*l/420; 469 var s = this.sin; 470 var c = this.cos; 471 this.Mloc.fromArray([ 472 [ 140*a, 0 , 0 , 70*a , 0 , 0 ], 473 [ 0 , 156*a ,-22*l*a , 0 , 54*a , 13*a*l ], 474 [ 0 ,-22*l*a, 4*a*l*l, 0 ,-13*a*l,-3*a*l*l], 475 [ 70*a , 0 , 0 , 140*a, 0 , 0 ], 476 [ 0 , 54*a ,-13*a*l , 0 , 156*a , 22*a*l ], 477 [ 0 , 13*a*l,-3*a*l*l, 0 , 22*a*l, 4*a*l*l]]); 478 /*this.M.fromArray([ 479 [ 140*a*c*c + 156*a*s*s, 16*a*c*s, -22*a*l*s, 70*a*c*c + 54*a*s*s, -16*a*c*s, 13*a*l*s], 480 [ 16*a*c*s, 156*a*c*c + 140*a*s*s, -22*a*c*l, -16*a*c*s, 54*a*c*c + 70*a*s*s, 13*a*c*l], 481 [ -22*a*l*s, -22*a*c*l, 4*a*l*l, -13*a*l*s, -13*a*c*l, -3*a*l*l], 482 [ 70*a*c*c + 54*a*s*s, -16*a*c*s, -13*a*l*s, 140*a*c*c + 156*a*s*s, 16*a*c*s, 22*a*l*s], 483 [ -16*a*c*s, 54*a*c*c + 70*a*s*s, -13*a*c*l, 16*a*c*s, 156*a*c*c + 140*a*s*s, 22*a*c*l], 484 [ 13*a*l*s, 13*a*c*l, -3*a*l*l, 22*a*l*s, 22*a*c*l, 4*a*l*l]]);*/ 485 }, 486 487 /** Sets geometry, stiffness matrices and "beam" load (end forces induced by continuous load and temperature load) 488 */ 489 Beam2d.prototype.linearStaticPreproc = function() { 490 this.setGeom(); 491 this.setStiffAndTrsfMats(); 492 this.load.zero(); 493 if (this.hasBodyLoad()) { 494 this.loadCC.fromArray([ 495 this.ea*this.alpha*this.dT - 0.5*this.fxloc*this.len , 496 -.5*this.fzloc*this.len , 497 this.fzloc*this.len*this.len/12. + this.ei*this.alpha*this.dTh , 498 -this.ea*this.alpha*this.dT - 0.5*this.fxloc*this.len , 499 -.5*this.fzloc*this.len , 500 -this.fzloc*this.len*this.len/12. - this.ei*this.alpha*this.dTh ]); 501 if (this.type) { 502 this.load.zero(); 503 this.load.setSubVector(this.a,this.loadCC.getSubVector(this.a).sub(this.kab.mulm(this.kbbinv).mulv(this.loadCC.getSubVector(this.b)))) 504 } 505 else { 506 this.load = this.loadCC; 507 } 508 } 509 }, 510 511 /** Computes end displacement in local coordinate system from given global components and respective end local end forces 512 * @param {Vector} dspl vector of beam displacement in globalComponents (if global param is true) or local (if globalComponents param is false) coordinate system 513 * @param {bool} [globalComponents=true] if true, dspl is considered in global coordinate system. If false, dspl is considered in local coordinate system 514 */ 515 Beam2d.prototype.linearStaticPostpro = function(dspl,globalComponents) { 516 var globalComponents = globalComponents || true; 517 var dspl = globalComponents? this.t.mulv(dspl) : dspl; 518 if (this.type) { 519 this.dspl.zero(); 520 this.dspl.setSubVector(this.a,dspl); 521 this.dspl.setSubVector(this.b,this.kbbinv.mulv(this.loadCC.getSubVector(this.b).add(this.kba.mulv(dspl))).negated()); 522 this.load.zero(); 523 this.load.setSubVector(this.a,this.kloc.mulv(dspl).add(this.loadCC.getSubVector(this.a).sub(this.kab.mulm(this.kbbinv).mulv(this.loadCC.getSubVector(this.b))))) 524 } 525 else { 526 this.dspl = dspl; 527 this.load = this.kloc.mulv(this.dspl).add(this.loadCC); 528 } 529 }, 530 531 /** Sets geometry and stiffness and mass matrices 532 */ 533 Beam2d.prototype.eigenvalueDynamicPreproc = function() { 534 this.setGeom(); 535 this.setStiffAndTrsfMats(); 536 this.setMassMats(); 537 }, 538 539 /** Compute local components of given eigenshape 540 * @param {int} eigMode number of assumed eigenshape - numbering from 0 541 * @param {Vector} eigShape eigMode-th eigenshape in global coordinate system 542 */ 543 Beam2d.prototype.eigenvalueDynamicPostpro = function(eigMode,eigShape) { 544 this.eigShapes[eigMode] = this.tt.mulv(eigShape); 545 }, 546 547 /** Returns middle deflection of receiver due to continuous load (clamped ends are assumed) 548 * @returns {float} middle deflection from continuous load 549 */ 550 Beam2d.prototype.giveMidDeflFromLoad = function() { // local z deflection in the middle of doble-clamped beam from fzloc 551 return this.fzloc*this.len*this.len*this.len*this.len / (384*this.ei); 552 }, 553 554 /** Returns middle deflection of receiver due to end displacement (no load is assumed) 555 * @returns {float} middle deflection from end displacement 556 */ 557 Beam2d.prototype.giveMidDeflFromDspl = function() { // local z deflection in the middle of the beam form nodal dspl and fzloc=0 558 return .5*(this.dspl[1] + this.dspl[4]) + this.len/8*( -this.dspl[2] + this.dspl[5]); 559 }, 560 561 /** Returns middle deflection of receiver 562 * @returns {float} middle deflection 563 */ 564 Beam2d.prototype.giveMidDefl = function() { 565 return this.giveMidDeflFromDspl() + this.giveMidDeflFromLoad(); 566 }, 567 568 /** Returns middle deflection of receiver of eigMode-th eigenshape 569 * @param {int} eigMode number of asumed eigenshape - numbering from 0 570 * @returns {float} middle deflection of eigMode-th eigenshape 571 */ 572 Beam2d.prototype.giveMidDeflOfEigShape = function(eigMode) { // local z deflection in the middle of the beam of eigMode-th eigenShape 573 var s = this.eigShapes[eigMode]; 574 return .5*(s[1] + s[4]) + this.len/8*( -s[2] + s[5]); 575 }, 576 577 /** Returns middle bending of receiver 578 * @returns {float} middle deflection 579 */ 580 Beam2d.prototype.giveMidMoment = function() { // bending moment in the middle of the beam 581 return -this.load[2] - this.load[1]*this.len*.5 - this.fzloc*this.len*this.len*.125; 582 }, 583 /**#@-*/ 584 585 /** Constructor, see {@link Beam2d} for input parameters description 586 * @returns {Beam2d} new Beam2d object 587 */ 588 Beam2d.create = function(node1,node2,params) { 589 var ret = new Beam2d(node1,node2,params); 590 ret.setGeom(); 591 return ret; 592 } 593 594 595 596 /** Linear static solver implementation 597 * @class Represents linear static solver 598 * @param {[Beams]} beams array of beams representing domain 599 * @param {[ints]} u array of unsupported dofs 600 * @param {[ints]} p array of supported dofs 601 * @param {Vector} [dspl=Vector.Zeros] vector of prescribed nodal displacements 602 * @param {Vector} [nodalLoad=Vector.Zeros] vector of external load in nodes (load from beams, i.e. linear load, temperature etc. will be added separately from beams) 603 * @param {bool} [withPrescribedDspl=true] if there is some nonzero prescribed displacement of support 604 * @param {bool} [computeReactions=false] compute global reactions or not 605 * @property {[Beams]} beams array of beams to include 606 * @property {[ints]} u array of unsupported dofs 607 * @property {[ints]} p array of supported dofs 608 * @property {int} nDofs total number of degrees of freedom of solved structure 609 * @property {Vector} nodalLoad vector of external load in nodes (load from beams, i.e. linear load, temperature etc. will be added separately from beams) 610 * @property {bool} [withPrescribedDspl=false] if there is some nonzero prescribed displacement of support 611 * @property {Matrix} K global stiffnesss matrix of whole structure 612 * @property {Vector} dspl displacement vector 613 * @property {Vector} beamLoad vector of "beam loads", i.e. end forces from beam continuous load and temperature load 614 * @property {Vector} load total load (sum od nodalLoad and beamLoad) 615 */ 616 LinearStaticSolver = function(beams,u,p,dspl,nodalLoad,withPrescribedDspl,computeReactions) { 617 this.beams = beams; 618 this.u = u; 619 this.p = p; 620 this.nDofs = this.u.length + this.p.length; 621 this.nodalLoad = nodalLoad==undefined? Vector.Zeros(this.nDofs) : nodalLoad; 622 this.withPrescribedDspl = withPrescribedDspl==undefined? false : withPrescribedDspl; 623 this.computeReactions = computeReactions==undefined? false : computeReactions; 624 this.k = Matrix.Zeros(this.nDofs,this.nDofs); 625 this.dspl = dspl==undefined? Vector.Zeros(this.nDofs) : dspl; 626 this.beamLoad = Vector.Zeros(this.nDofs); 627 this.load = Vector.Zeros(this.nDofs); 628 } 629 /**#@+ 630 * @function 631 * @memberOf LinearStaticSolver# 632 */ 633 634 /** Assembles global stiffness matrix of whole structure from stiffness matrices of particular beams 635 */ 636 LinearStaticSolver.prototype.assemble = function() { 637 this.k.zero(); 638 for (var beamNum=0; beamNum<this.beams.length; beamNum++) { 639 var beam = this.beams[beamNum]; 640 var dofs = beam.dofs; 641 beam.computeGlobStiffMatrix().assemble(this.k,dofs,dofs); 642 } 643 }, 644 645 /** Solves unknown displacements and forces 646 */ 647 LinearStaticSolver.prototype.solve = function() { 648 this.load = this.nodalLoad.add(this.beamLoad); 649 var kuu = this.k.getSubMatrix(this.u,this.u); 650 var kpu = this.k.getSubMatrix(this.p,this.u); 651 var fu = this.load.getSubVector(this.u); 652 if (this.withPrescribedDspl) { 653 var kup = this.k.getSubMatrix(this.u,this.p); 654 var kpp = this.k.getSubMatrix(this.p,this.p); 655 var rp = this.dspl.getSubVector(this.p); 656 this.dspl.setSubVector(this.u,kuu.linSolve(fu.sub(kup.mulv(rp)))); 657 if (this.computeReactions) { 658 this.load.setSubVector(this.p,kpu.mulv(this.dspl.getSubVector(this.u)).add(kpp.mulv(rp))); 659 this.load.isub(this.beamLoad); 660 } 661 } 662 else { 663 this.dspl.setSubVector(this.u,kuu.linSolve(fu)); 664 if (this.computeReactions) { 665 this.load.setSubVector(this.p,kpu.mulv(this.dspl.getSubVector(this.u))); 666 this.load.isub(this.beamLoad); 667 } 668 } 669 }, 670 671 /** Preprocesses all beams (sets their "beam load" and assemble it to global load vecor, set their geometric parameters and stiffness matrices 672 */ 673 LinearStaticSolver.prototype.preprocBeams = function() { 674 this.beamLoad.zero() 675 for (var beamNum=0; beamNum<this.beams.length; beamNum++) { 676 var beam = this.beams[beamNum]; 677 beam.linearStaticPreproc(); 678 if (beam.type) { 679 this.beamLoad.decrSubVector(beam.dofs,beam.tt.mulv(beam.load.getSubVector(beam.a))); 680 } else { 681 this.beamLoad.decrSubVector(beam.dofs,beam.tt.mulv(beam.load)); 682 } 683 } 684 }, 685 686 /** Postprocesses all beams (sets appropriate displacements to them) 687 */ 688 LinearStaticSolver.prototype.updateBeams = function() { 689 for (var beamNum=0; beamNum<this.beams.length; beamNum++) { 690 var beam = this.beams[beamNum]; 691 beam.linearStaticPostpro(this.dspl.getSubVector(beam.dofs)); 692 } 693 }, 694 695 /** Global function for whole simulation process 696 */ 697 LinearStaticSolver.prototype.update = function(preprocBeams,updateBeams) { 698 var preprocBeams = preprocBeams==undefined? true : preprocBeams; 699 var updateBeams = updateBeams==undefined? true : updateBeams; 700 if (preprocBeams) { this.preprocBeams(); } 701 this.assemble(); 702 this.solve(); 703 if (updateBeams) { this.updateBeams(); } 704 }, 705 /**#@-*/ 706 707 /** Constructor, see {@link LinearStaticSolver} for input parameters description 708 * @returns {LinearStaticSolver} new LinearStaticSolver object 709 */ 710 LinearStaticSolver.create = function(beams,u,p,dspl,nodalLoad,withPrescribedDspl,computeReactions) { 711 var ret = new LinearStaticSolver(beams,u,p,dspl,nodalLoad,withPrescribedDspl,computeReactions); 712 return ret; 713 } 714 715 716 717 718 /** Eigenvalue dynamic solver implementation 719 * @class represents eigenvalue dynamic solver 720 * @param {[Beams]} beams array of beams to include 721 * @param {[ints]} u array of unsupported dofs 722 * @param {[ints]} p array of supported dofs 723 * @param {int} [nEigModes=1] number of first eigenvalues and eigenvectors to be solved 724 * @param {Vector} [concMass=Vector.Zeros] Vector of concentrated mass for each Dof (default zero Vector) 725 * @property {[Beams]} beams array of beams to include 726 * @property {[ints]} u array of unsupported dofs 727 * @property {[ints]} p array of supported dofs 728 * @property {int} nDofs total number of degrees of freedom of solved structure 729 * @property {int} [nEigModes=1] number of first eigenvalues and eigenvectors to be solved 730 * @property {Vector} [concMass=Vector.Zeros] Vector of concentrated mass for each Dof (default zero Vector) 731 * @property {Matrix} k global stiffness matrix 732 * @property {Matrix} M global mass matrix 733 * @property {[floats]} eigFreqs computed eigenfrequencies 734 * @property {[Vectors]} eigShapes computed eigenshapes 735 */ 736 EigenValueDynamicSolver = function(beams,u,p,nEigModes,concMass) { 737 this.beams = beams; 738 this.u = u; 739 this.p = p; 740 this.nDofs = this.u.length + this.p.length; 741 this.nEigModes = nEigModes==undefined? 1 : nEigModes; 742 this.concMass = concMass==undefined? Vector.Zeros(this.nDofs) : concMass; 743 this.k = Matrix.Zeros(this.nDofs,this.nDofs); 744 this.m = Matrix.Zeros(this.nDofs,this.nDofs); 745 this.eigFreqs = []; 746 for (var eigMode=0; eigMode<this.nEigModes; eigMode++) { this.eigFreqs[eigMode] = 0.; } 747 this.eigShapes = []; 748 for (var eigMode=0; eigMode<this.nEigModes; eigMode++) { this.eigShapes[eigMode] = Vector.Zeros(this.nDofs); } 749 } 750 /**#@+ 751 * @function 752 * @memberOf EigenValueDynamicSolver# 753 */ 754 755 /** Assembles global stiffness and mass matrix of whole structure from stiffness and mass matrices of particular beams and concMass attribute 756 */ 757 EigenValueDynamicSolver.prototype.assemble = function() { 758 this.k.zero(); 759 this.m.zero() 760 for (var beamNum=0; beamNum<this.beams.length; beamNum++) { 761 var beam = this.beams[beamNum]; 762 var dofs = beam.dofs; 763 beam.computeGlobStiffMatrix().assemble(this.k,dofs,dofs); 764 this.m.incrSubMatrix(dofs,dofs,beam.m); 765 } 766 for (var i=0; i<this.nDofs; i++) { this.m[i][i] += this.concMass[i]; } 767 }, 768 769 /** Solves unknown eigenfrequencies and eigenshapes 770 */ 771 EigenValueDynamicSolver.prototype.solve = function() { 772 var kuu = this.k.getSubMatrix(this.u,this.u); 773 var kuu = this.m.getSubMatrix(this.u,this.u); 774 var results = eigSolve(Kuu,Muu,this.nEigModes); 775 for (var eigMode=0; eigMode<this.nEigModes; eigMode++) { 776 this.eigFreqs[eigMode] = results[0][eigMode]; 777 this.eigShapes[eigMode].setSubVector(this.u,results[1][eigMode]); 778 } 779 }, 780 781 /** Preprocesses all beams (sets their "beam load" and assemble it to global load vecor, set their geometric parameters and stiffness matrices 782 */ 783 EigenValueDynamicSolver.prototype.preprocBeams = function() { 784 for (var beamNum = 0; beamNum<this.beams.length; beamNum++) { this.beams[beamNum].eigenvalueDynamicPreproc(); } 785 }, 786 787 /** Postprocesses all beams (from computed eigenshapes sets appropriate displacements to them) 788 */ 789 EigenValueDynamicSolver.prototype.updateBeams = function() { 790 for (var eigMode = 0; eigMode<this.nEigModes; eigMode++) { 791 for (var beamNum=0; beamNum<this.beams.length; beamNum++) { 792 var beam = this.beams[beamNum]; 793 beam.eigenvalueDynamicPostpro(eigMode,this.eigShapes[eigMode].getSubVector(beam.dofs)); 794 } 795 } 796 }, 797 798 /** Global function for whole simulation process 799 */ 800 EigenValueDynamicSolver.prototype.update = function(preprocBeams,updateBeams) { 801 var preprocBeams = preprocBeams==undefined? true : preprocBeams; 802 var updateBeams = updateBeams==undefined? true : updateBeams; 803 if (preprocBeams) { this.preprocBeams(); } 804 this.assemble(); 805 this.solve(); 806 if (updateBeams) { this.updateBeams(); } 807 }, 808 /**#@-*/ 809 810 /** Constructor, see {@link EigenValueDynamicSolver} for input parameters description 811 * @returns {EigenValueDynamicSolver} new EigenValueDynamicSolver object 812 */ 813 EigenValueDynamicSolver.create = function(beams,u,p,nEigModes,concMass) { 814 var ret = new EigenValueDynamicSolver(beams,u,p,nEigModes,concMass); 815 return ret; 816 } 817 818 819 820 821 822 823 824 825 826 827 /** Cross section implementation 828 * @class represents (so far polygonal only) cross section 829 * @param {[Nodes]} nodes array of nodes to create yz (!!) polygon in anti-clockwise (!!) order 830 * @param {[[Nodes]]} [holes=[[]]] array of arrays of nodes representing holes in polygon created by nodes 831 * @property {[Nodes]} nodes array of polygon nodes in yz plane ( this.nodes[0]==this.nodes[this.nodes.length-1] ) 832 * @property {float} [a=null] area 833 * @property {float} [sy=null] static moment (first moment of mass) with respect to y coordinate (global) 834 * @property {float} [sz=null] static moment (first moment of mass) with respect to z coordinate (global) 835 * @property {float} [cy=null] y coordinate of center of mass 836 * @property {float} [cz=null] z coordinate of center of mass 837 * @property {float} [iy=null] moment of inertia (second moment of area) with respect to y coordinate (local) 838 * @property {float} [iz=null] moment of inertia (second moment of area) with respect to z coordinate (local) 839 * @property {float} [dyz=null] deviation (product) moment of inertia with respect to yz coordinate system 840 * @property {float} [i1=null] first principal value of moment of inertia (i1>i2) 841 * @property {float} [i2=null] second principal value of moment of inertia (i1>i2) 842 * @property {float} [beta=null] angle from global axes to principal inertia orientation 843 * @property {[Nodes]} [core=null] core of cross section 844 */ 845 CrossSection = function(nodes,holes) { 846 this.nodes = nodes; 847 if (this.nodes[this.nodes.length-1] != this.nodes[0]) { 848 this.nodes.push(this.nodes[0]); 849 } 850 this.holes = holes || [[]]; 851 this.a = null; 852 this.sy = null; 853 this.sz = null; 854 this.cy = null; 855 this.cz = null; 856 this.iy = null; 857 this.iz = null; 858 this.dyz = null; 859 this.i1 = null; 860 this.i2 = null; 861 this.beta = null; 862 this.core = null; 863 } 864 /**#@+ 865 * @function 866 * @memberOf CrossSection# 867 */ 868 869 /** Reset receiver (sets all values to null)*/ 870 CrossSection.prototype.reset = function() { 871 this.a = null; 872 this.sy = null; 873 this.sz = null; 874 this.cy = null; 875 this.cz = null; 876 this.iy = null; 877 this.iz = null; 878 this.dyz = null; 879 this.i1 = null; 880 this.i2 = null; 881 this.beta = null; 882 this.convexHull = null; 883 this.core = null; 884 }, 885 886 /** Computes area of receiver 887 * @returns {float} area of receiver*/ 888 CrossSection.prototype.computeA = function() { 889 var ret = 0.; 890 var n0, nL; 891 for (var i=0; i<this.nodes.length-1; i++) { 892 n0 = this.nodes[i]; 893 nL = this.nodes[i+1]; 894 ret += n0.y*nL.z - nL.y*n0.z; 895 } 896 return .5*ret; 897 }, 898 899 /** Computes static moment (first order moment of mass) of receiver with respect to y coordinate axis 900 * @param {bool} [local=false] if computed with respect to this.cy center of mass or global yz coordinates 901 * @returns {float} Sy of receiver*/ 902 CrossSection.prototype.computeSy = function(local) { 903 var local = local==undefined? false : local; 904 var zc = 0.; 905 if (local) { 906 if (this.cz==null) { return null; } 907 zc = this.cz; 908 } 909 var ret = 0.; 910 var n0, nL, y0, yL, z0, zL; 911 for (var i=0; i<this.nodes.length-1; i++) { 912 n0 = this.nodes[i]; 913 nL = this.nodes[i+1]; 914 y0 = n0.y; yL = nL.y; 915 z0 = n0.z-zc; zL = nL.z-zc; 916 ret += (zL+z0)*(y0*zL-yL*z0); 917 } 918 return 1/6.*ret; 919 }, 920 921 /** Computes static moment (first order moment of mass) of receiver with respect to z coordinate axis 922 * @param {bool} [local=false] if computed with respect to this.cz center of mass or global yz coordinates 923 * @returns {float} Sz of receiver*/ 924 CrossSection.prototype.computeSz = function(local) { 925 var local = local==undefined? false : local; 926 var yc = 0.; 927 if (local) { 928 if (this.cy==null) { return null; } 929 yc = this.cy; 930 } 931 var ret = 0.; 932 var n0, nL, y0, yL, z0, zL; 933 for (var i=0; i<this.nodes.length-1; i++) { 934 n0 = this.nodes[i]; 935 nL = this.nodes[i+1]; 936 y0 = n0.y-yc; yL = nL.y-yc; 937 z0 = n0.z; zL = nL.z; 938 ret += (y0+yL)*(y0*zL-yL*z0); 939 } 940 return 1/6.*ret; 941 }, 942 943 /** Computes y coordinate of center of mass of receiver 944 * @returns {float} Cy of receiver*/ 945 CrossSection.prototype.computeCy = function() { 946 var a = this.a==null? this.computeA() : this.a; 947 var sz = this.sz==null? this.computeSz() : this.sz; 948 return sz/a; 949 }, 950 951 /** Computes z coordinate of center of mass of receiver 952 * @returns {float} Cz of receiver*/ 953 CrossSection.prototype.computeCz = function() { 954 var a = this.a==null? this.computeA() : this.a; 955 var sy = this.sy==null? this.computeSy() : this.sy; 956 return sy/a; 957 }, 958 959 /** Computes moment of inertia (second moment of area) of receiver with respect to y coordinate axis 960 * @param {bool} [local=true] if computed with respect to this.cy center of mass or global yz coordinates 961 * @returns {float} moment of inertia (y) of receiver*/ 962 CrossSection.prototype.computeIy = function(local) { 963 var local = local==undefined? true : local; 964 var zc = 0.; 965 if (local) { 966 if (this.cz==null) { return null; } 967 zc = this.cz; 968 } 969 var ret = 0.; 970 var n0, nL, y0, yL, z0, zL; 971 for (var i=0; i<this.nodes.length-1; i++) { 972 n0 = this.nodes[i]; 973 nL = this.nodes[i+1]; 974 y0 = n0.y; yL = nL.y; 975 z0 = n0.z-zc; zL = nL.z-zc; 976 ret += (z0*z0+z0*zL+zL*zL)*(y0*zL-yL*z0); 977 } 978 return 1/12.*ret; 979 }, 980 981 /** Computes moment of inertia (second moment of area) of receiver with respect to z coordinate axis 982 * @param {bool} [local=true] if computed with respect to this.cz center of mass or global yz coordinates 983 * @returns {float} moment of inertia (z) of receiver*/ 984 CrossSection.prototype.computeIz = function(local) { 985 var local = local==undefined? true : local; 986 var yc = 0.; 987 if (local) { 988 if (this.cy==null) { return null; } 989 yc = this.cy; 990 } 991 var ret = 0.; 992 var n0, nL, y0, yL, z0, zL; 993 for (var i=0; i<this.nodes.length-1; i++) { 994 n0 = this.nodes[i]; 995 nL = this.nodes[i+1]; 996 y0 = n0.y-yc; yL = nL.y-yc; 997 z0 = n0.z; zL = nL.z; 998 ret += (y0*y0+y0*yL+yL*yL)*(y0*zL-yL*z0); 999 } 1000 return 1/12.*ret; 1001 }, 1002 1003 /** Computes deviation (prodict) moment of inertia of receiver with respect to yz coordinates 1004 * @param {bool} [local=true] if computed with respect to this.cy and this.cz (center of mass) or global yz coordinates 1005 * @returns {float} deviation (product) moment of inertia (yz) of receiver*/ 1006 CrossSection.prototype.computeDyz = function(local) { 1007 var local = local==undefined? true : local; 1008 var yc = 0.; 1009 var zc = 0.; 1010 if (local) { 1011 if (this.cy==null) { return null; } 1012 yc = this.cy; 1013 if (this.cz==null) { return null; } 1014 zc = this.cz; 1015 } 1016 var ret = 0.; 1017 var n0, nL, y0, yL, z0, zL; 1018 for (var i=0; i<this.nodes.length-1; i++) { 1019 n0 = this.nodes[i]; 1020 nL = this.nodes[i+1]; 1021 y0 = n0.y-yc; yL = nL.y-yc; 1022 z0 = n0.z-zc; zL = nL.z-zc; 1023 ret += (y0*zL+2*y0*z0+2*yL*zL+yL*z0)*(y0*zL-yL*z0); 1024 } 1025 return 1/24.*ret; 1026 }, 1027 1028 /** Computes polar moment of inertia of receiver with respect to yz coordinate axis 1029 * @returns {float polar} moment of inertia (yz) of receiver*/ 1030 CrossSection.prototype.computeI0 = function() { 1031 var iy = this.iy==null? this.computeIy() : this.iy; 1032 var iz = this.iz==null? this.computeIz() : this.iz; 1033 return iy+iz; 1034 }, 1035 1036 /** Computes moment of inertia tensor of receiver (i.e. [[Iy,-Dyz],[-Dyz,Iz]] ) 1037 * @returns {[[float,float],[float,float]]} moment of inertia tensor (yz) of receiver*/ 1038 CrossSection.prototype.computeInertiaTensor = function() { 1039 var iy = this.iy==null? this.computeIy() : this.iy; 1040 var iz = this.iz==null? this.computeIz() : this.iz; 1041 var dyz = this.dyz==null? this.computeDyz() : this.dyz; 1042 return [[iy,-dyz],[-dyz,iz]]; 1043 }, 1044 1045 /** Computes principal inertia values and angle from global axes to principal inertia orientation 1046 * @returns {[float,float,float]} principal values and angle of principal direction [I1,I2,angle], I1 > I2 */ 1047 CrossSection.prototype.computePrincipalValues = function() { 1048 var iy = this.iy==null? this.computeIy() : this.iy; 1049 var iz = this.iz==null? this.computeIz() : this.iz; 1050 var dyz = this.dyz==null? this.computeDyz() : this.dyz; 1051 if (iy==iz) { return [iy,iz,0.]; } 1052 if (dyz == 0.) { 1053 if (iy>iz) { return [iy,iz,0.]; } 1054 else { return [iz,iy,.5*Math.PI]; } 1055 } 1056 var beta = 0.5*Math.atan(2*dyz/(iz-iy)); 1057 var c = Math.cos(beta) 1058 var s = Math.sin(beta) 1059 var i1 = c*c*iy + s*s*iz - 2*c*s*dyz; 1060 var i2 = s*s*iy + c*c*iz + 2*c*s*dyz; 1061 if (i1>i2) { return [i1,i2,beta]; } 1062 else { return [i2,i1,beta+0.5*Math.PI]; } 1063 }, 1064 1065 /** Computes normal stress at certain point with coordinates x and y 1066 * @param {float} y y coordinate of desired point 1067 * @param {float} z z coordinate of desired point 1068 * @param {float} N normal force 1069 * @param {float} My moment y 1070 * @param {float} Mz moment z 1071 * @param {bool} [local=false] if y z are local coordinates or not 1072 * @returns {float} stress at [x,z] caused by N,My,Mz*/ 1073 CrossSection.prototype.computeStress = function(y,z,N,My,Mz,local) { 1074 var local = local || false; 1075 var a = this.a==null? this.computeA() : this.a; 1076 var iy = this.iy==null? this.computeIy() : this.iy; 1077 var iz = this.iz==null? this.computeIz() : this.iz; 1078 var dyz = this.dyz==null? this.computeDyz() : this.dyz; 1079 if (local) { var y=y; var z=z; } 1080 else { 1081 var y = y - (this.cy==null? this.computeCy() : this.cy); 1082 var z = z - (this.cz==null? this.computeCz() : this.cz); 1083 } 1084 return N/a + (My*iz+Mz*dyz)/(iy*iz-dyz*dyz)*z - (Mz*iy+My*dyz)/(iy*iz-dyz*dyz)*y; 1085 }, 1086 1087 /** Computes normal stress at certain polygon node 1088 * @param {Node} node polygon node 1089 * @param {float} N normal force 1090 * @param {float} My moment y 1091 * @param {float} Mz moment z 1092 * @returns {float} stress at node caused by N,My,Mz*/ 1093 CrossSection.prototype.computeStressInNode = function(node,N,My,Mz) { 1094 return this.computeStress(node.y,node.z,N,My,Mz); 1095 }, 1096 1097 /** Computes neutral axis 1098 * @param {float} N normal force 1099 * @param {float} My moment y 1100 * @param {float} Mz moment z 1101 * @param {bool} [local=false] if returned equation is in local coordinates or not 1102 * @returns {[float,float,float]} [a,b,c] as ay+bz+c=0;*/ 1103 CrossSection.prototype.computeNeutralAxis = function(N,My,Mz,local) { 1104 var local = local || false; 1105 var a = this.a==null? this.computeA() : this.a; 1106 var iy = this.iy==null? this.computeIy() : this.iy; 1107 var iz = this.iz==null? this.computeIz() : this.iz; 1108 var dyz = this.dyz==null? this.computeDyz() : this.dyz; 1109 var aa = -(Mz*iy+My*dyz)/(iy*iz-dyz*dyz); 1110 var bb = (My*iz+Mz*dyz)/(iy*iz-dyz*dyz); 1111 var cc = N/a; 1112 if (local) { 1113 var yc = this.cy==null? this.computeCy() : this.cy; 1114 var zc = this.cz==null? this.computeCz() : this.cz; 1115 cc += aa*yc+bb*zc; 1116 } 1117 return [aa,bb,cc]; 1118 }, 1119 1120 /** Returns those nodes of receiver that forms the convex hull 1121 * @returns {[Nodes]} those nodes that forms convex hull (ret[0]==ret[ret.length-1])*/ 1122 CrossSection.prototype.giveConvexHull = function() { 1123 if (this.nodes.length <= 4) { return this.nodes.slice(); } 1124 // choose two nodes (most lower-left and upper-right) that for sure forms a convex part 1125 var c1 = this.nodes[0]; 1126 var c2 = this.nodes[0]; 1127 var i1 = 0; 1128 var i2 = 0; 1129 for (var i=1; i<this.nodes.length; i++) { 1130 if (this.nodes[i].z < c1.z) { c1 = this.nodes[i]; i1 = i; continue; } 1131 if (this.nodes[i].z > c2.z) { c2 = this.nodes[i]; i2 = i; continue; } 1132 if (this.nodes[i].z == c1.z && this.nodes[i].y < c1.y) { c1 = this.nodes[i]; i1 = i; } 1133 if (this.nodes[i].z == c2.z && this.nodes[i].y > c2.y) { c2 = this.nodes[i]; i2 = i; } 1134 } 1135 1136 n = this.nodes.slice(0,this.nodes.length-1) 1137 var nodes12 = i1<i2? n.slice(i1,i2) : n.slice(i1).concat(n.slice(0,i2)); 1138 var nodes21 = i2<i1? n.slice(i2,i1) : n.slice(i2).concat(n.slice(0,i1)); 1139 nodes12.push(c2); nodes21.push(c1); 1140 1141 var v, c, s, cMin, iCur, cCur, vCur, l12 = nodes12.length, l21 = nodes21.length; 1142 var ret12 = [] 1143 iCur = 0; 1144 if (l12 > 2) { 1145 do { 1146 ret12.push(nodes12[iCur]); 1147 cCur = nodes12[iCur] 1148 v = c2.sub(cCur) 1149 cMin = 1e20; 1150 for (var i=iCur+1; i<l12; i++) { 1151 vCur = nodes12[i].sub(cCur) 1152 c = vCur.dot(v)/(vCur.norm()*v.norm()) 1153 if (c <= cMin) { if (vCur.cross(v).x >= 0.) { iCur = i; cMin=c; } } 1154 } 1155 } while (iCur<l12-1) 1156 } else { ret12 = [c1]; } 1157 1158 var ret21 = []; 1159 iCur = 0; 1160 if (l21 > 2) { 1161 do { 1162 ret21.push(nodes21[iCur]); 1163 cCur = nodes21[iCur] 1164 v = c1.sub(cCur) 1165 cMin = 1e20; 1166 for (var i=iCur+1; i<l21; i++) { 1167 vCur = nodes21[i].sub(cCur) 1168 c = vCur.dot(v)/(vCur.norm()*v.norm()) 1169 if (c <= cMin) { if (vCur.cross(v).x >= 0.) { iCur = i; cMin=c; } } 1170 } 1171 } while (iCur<l21-1) 1172 } else { ret21 = [c2]; } 1173 return ret12.concat(ret21).concat(c1) 1174 }, 1175 1176 /** Computes core (area, where applied force in x direction causes no tension) 1177 * @param {bool} [local=false] if computed with respect to this.cx center of mass or global yz coordinates 1178 * @returns {[Nodes]} array of nodes creating core polygon*/ 1179 CrossSection.prototype.computeCore = function(local) { 1180 var core = []; 1181 var nodes = this.convexHull==null? this.giveConvexHull() : this.convexHull; 1182 var local = local || false; 1183 var a = this.a==null? this.computeA() : this.a; 1184 var iy = this.iy==null? this.computeIy() : this.iy; 1185 var iz = this.iz==null? this.computeIz() : this.iz; 1186 var dyz = this.dyz==null? this.computeDyz() : this.dyz; 1187 var yc = this.cy==null? this.computeCy() : this.cy; 1188 var zc = this.cz==null? this.computeCz() : this.cz; 1189 var y,z,n1,n2,y1,y2,z1,z2,cc; 1190 for (var i=0; i<nodes.length-1; i++) { 1191 n1 = nodes[i]; 1192 n2 = nodes[i+1]; 1193 y1 = n1.y - yc; 1194 y2 = n2.y - yc; 1195 z1 = n1.z - zc; 1196 z2 = n2.z - zc; 1197 //y = iz2*(z1-z2)/(y1*z2-y2*z1); 1198 //z = iy2*(y2-y1)/(y1*z2-y2*z1); 1199 cc = a*(y1*z2-y2*z1); 1200 y = (dyz*(y2-y1)+iz*(z1-z2))/cc; 1201 z = (dyz*(z1-z2)+iy*(y2-y1))/cc; 1202 core[i] = local? new Node(0.,y,z) : new Node(0.,y+yc,z+zc); 1203 } 1204 return core; 1205 }, 1206 1207 /** computes all values (a,sy,sz,cy,cz,iy,iz,dyz,i1,i2,beta,core) within one function call*/ 1208 CrossSection.prototype.update = function() { 1209 this.a = this.computeA(); 1210 this.sy = this.computeSy(); 1211 this.sz = this.computeSz(); 1212 this.cy = this.computeCy(); 1213 this.cz = this.computeCz(); 1214 this.iy = this.computeIy(); 1215 this.iz = this.computeIz(); 1216 this.dyz = this.computeDyz(); 1217 var princ = this.computePrincipalValues(); 1218 this.i1 = princ[0]; 1219 this.i2 = princ[1]; 1220 this.beta = princ[2]; 1221 this.convexHull = this.giveConvexHull(); 1222 this.core = this.computeCore(); 1223 }, 1224 /**#@-*/ 1225 1226 /** Constructor, see {@link CrossSection} for input parameters description 1227 * @returns {CrossSection} new CrossSection object*/ 1228 CrossSection.create = function(nodes,holes) { 1229 var ret = new CrossSection(nodes,holes); 1230 return ret; 1231 } 1232