Real Hex8::volume () const { // Make copies of our points. It makes the subsequent calculations a bit // shorter and avoids dereferencing the same pointer multiple times. Point x0 = point(0), x1 = point(1), x2 = point(2), x3 = point(3), x4 = point(4), x5 = point(5), x6 = point(6), x7 = point(7); // Construct constant data vectors. The notation is: // \vec{x}_{\xi} = \vec{a1}*eta*zeta + \vec{b1}*eta + \vec{c1}*zeta + \vec{d1} // \vec{x}_{\eta} = \vec{a2}*xi*zeta + \vec{b2}*xi + \vec{c2}*zeta + \vec{d2} // \vec{x}_{\zeta} = \vec{a3}*xi*eta + \vec{b3}*xi + \vec{c3}*eta + \vec{d3} // but it turns out that a1, a2, and a3 are not needed for the volume calculation. // Build up the 6 unique vectors which make up dx/dxi, dx/deta, and dx/dzeta. Point q[6] = { /*b1*/ x0 - x1 + x2 - x3 + x4 - x5 + x6 - x7, /*=b2*/ /*c1*/ x0 - x1 - x2 + x3 - x4 + x5 + x6 - x7, /*=b3*/ /*d1*/ -x0 + x1 + x2 - x3 - x4 + x5 + x6 - x7, /*c2*/ x0 + x1 - x2 - x3 - x4 - x5 + x6 + x7, /*=c3*/ /*d2*/ -x0 - x1 + x2 + x3 - x4 - x5 + x6 + x7, /*d3*/ -x0 - x1 - x2 - x3 + x4 + x5 + x6 + x7 }; // We could check for a linear element, but it's probably faster to // just compute the result... return (triple_product(q[0], q[4], q[3]) + triple_product(q[2], q[0], q[1]) + triple_product(q[1], q[3], q[5])) / 192. + triple_product(q[2], q[4], q[5]) / 64.; }
Real Pyramid5::volume () const { // The pyramid with a bilinear base has volume given by the // formula in: "Calculation of the Volume of a General Hexahedron // for Flow Predictions", AIAA Journal v.23, no.6, 1984, p.954- Point x0 = point(0), x1 = point(1), x2 = point(2), x3 = point(3), x4 = point(4); // Construct various edge and diagonal vectors. Point v40 = x0 - x4; Point v13 = x3 - x1; Point v02 = x2 - x0; Point v03 = x3 - x0; Point v01 = x1 - x0; // Finally, ready to return the volume! return triple_product(v40, v13, v02) / 6. + triple_product(v02, v01, v03) / 12.; }
Real Pyramid14::volume () const { // Make copies of our points. It makes the subsequent calculations a bit // shorter and avoids dereferencing the same pointer multiple times. Point x0 = point(0), x1 = point(1), x2 = point(2), x3 = point(3), x4 = point(4), x5 = point(5), x6 = point(6), x7 = point(7), x8 = point(8), x9 = point(9), x10 = point(10), x11 = point(11), x12 = point(12), x13 = point(13); // dx/dxi and dx/deta have 15 components while dx/dzeta has 20. // These are copied directly from the output of a Python script. Point dx_dxi[15] = { x6/2 - x8/2, x0/4 - x1/4 + x10 + x11 - x12 - x2/4 + x3/4 - 3*x6/2 + 3*x8/2 - x9, -x0/2 + x1/2 - 2*x10 - 2*x11 + 2*x12 + x2/2 - x3/2 + 3*x6/2 - 3*x8/2 + 2*x9, x0/4 - x1/4 + x10 + x11 - x12 - x2/4 + x3/4 - x6/2 + x8/2 - x9, x0/4 - x1/4 + x2/4 - x3/4, -3*x0/4 + 3*x1/4 - x10 + x11 - x12 - 3*x2/4 + 3*x3/4 + x9, x0/2 - x1/2 + x10 - x11 + x12 + x2/2 - x3/2 - x9, -x0/4 + x1/4 + x2/4 - x3/4 - x6/2 + x8/2, x0/4 - x1/4 - x2/4 + x3/4 + x6/2 - x8/2, -2*x13 + x6 + x8, 4*x13 - 2*x6 - 2*x8, -2*x13 + x6 + x8, -x0/2 - x1/2 + x2/2 + x3/2 + x5 - x7, x0/2 + x1/2 - x2/2 - x3/2 - x5 + x7, x0/2 + x1/2 + 2*x13 + x2/2 + x3/2 - x5 - x6 - x7 - x8 }; // dx/dxi and dx/deta have 15 components while dx/dzeta has 20. // These are copied directly from the output of a Python script. Point dx_deta[15] = { -x5/2 + x7/2, x0/4 + x1/4 - x10 + x11 + x12 - x2/4 - x3/4 + 3*x5/2 - 3*x7/2 - x9, -x0/2 - x1/2 + 2*x10 - 2*x11 - 2*x12 + x2/2 + x3/2 - 3*x5/2 + 3*x7/2 + 2*x9, x0/4 + x1/4 - x10 + x11 + x12 - x2/4 - x3/4 + x5/2 - x7/2 - x9, -2*x13 + x5 + x7, 4*x13 - 2*x5 - 2*x7, -2*x13 + x5 + x7, x0/4 - x1/4 + x2/4 - x3/4, -3*x0/4 + 3*x1/4 - x10 + x11 - x12 - 3*x2/4 + 3*x3/4 + x9, x0/2 - x1/2 + x10 - x11 + x12 + x2/2 - x3/2 - x9, -x0/2 + x1/2 + x2/2 - x3/2 - x6 + x8, x0/2 - x1/2 - x2/2 + x3/2 + x6 - x8, -x0/4 - x1/4 + x2/4 + x3/4 + x5/2 - x7/2, x0/4 + x1/4 - x2/4 - x3/4 - x5/2 + x7/2, x0/2 + x1/2 + 2*x13 + x2/2 + x3/2 - x5 - x6 - x7 - x8 }; // dx/dxi and dx/deta have 15 components while dx/dzeta has 20. // These are copied directly from the output of a Python script. Point dx_dzeta[20] = { -x0/4 - x1/4 + x10 + x11 + x12 - 2*x13 - x2/4 - x3/4 - x4 + x9, 5*x0/4 + 5*x1/4 - 5*x10 - 5*x11 - 5*x12 + 8*x13 + 5*x2/4 + 5*x3/4 + 7*x4 - 5*x9, -9*x0/4 - 9*x1/4 + 9*x10 + 9*x11 + 9*x12 - 12*x13 - 9*x2/4 - 9*x3/4 - 15*x4 + 9*x9, 7*x0/4 + 7*x1/4 - 7*x10 - 7*x11 - 7*x12 + 8*x13 + 7*x2/4 + 7*x3/4 + 13*x4 - 7*x9, -x0/2 - x1/2 + 2*x10 + 2*x11 + 2*x12 - 2*x13 - x2/2 - x3/2 - 4*x4 + 2*x9, x0/4 + x1/4 - x10 + x11 + x12 - x2/4 - x3/4 + x5/2 - x7/2 - x9, -3*x0/4 - 3*x1/4 + 3*x10 - 3*x11 - 3*x12 + 3*x2/4 + 3*x3/4 - 3*x5/2 + 3*x7/2 + 3*x9, 3*x0/4 + 3*x1/4 - 3*x10 + 3*x11 + 3*x12 - 3*x2/4 - 3*x3/4 + 3*x5/2 - 3*x7/2 - 3*x9, -x0/4 - x1/4 + x10 - x11 - x12 + x2/4 + x3/4 - x5/2 + x7/2 + x9, x0/4 - x1/4 + x10 + x11 - x12 - x2/4 + x3/4 - x6/2 + x8/2 - x9, -3*x0/4 + 3*x1/4 - 3*x10 - 3*x11 + 3*x12 + 3*x2/4 - 3*x3/4 + 3*x6/2 - 3*x8/2 + 3*x9, 3*x0/4 - 3*x1/4 + 3*x10 + 3*x11 - 3*x12 - 3*x2/4 + 3*x3/4 - 3*x6/2 + 3*x8/2 - 3*x9, -x0/4 + x1/4 - x10 - x11 + x12 + x2/4 - x3/4 + x6/2 - x8/2 + x9, -x0/4 + x1/4 - x10 + x11 - x12 - x2/4 + x3/4 + x9, x0/4 - x1/4 + x10 - x11 + x12 + x2/4 - x3/4 - x9, -x0/4 + x1/4 + x2/4 - x3/4 - x6/2 + x8/2, x0/4 - x1/4 - x2/4 + x3/4 + x6/2 - x8/2, -x0/4 - x1/4 + x2/4 + x3/4 + x5/2 - x7/2, x0/4 + x1/4 - x2/4 - x3/4 - x5/2 + x7/2, x0/2 + x1/2 + 2*x13 + x2/2 + x3/2 - x5 - x6 - x7 - x8 }; // The (xi, eta, zeta) exponents for each of the dx_dxi terms static const int dx_dxi_exponents[15][3] = { {0, 0, 0}, {0, 0, 1}, {0, 0, 2}, {0, 0, 3}, {0, 1, 0}, {0, 1, 1}, {0, 1, 2}, {0, 2, 0}, {0, 2, 1}, {1, 0, 0}, {1, 0, 1}, {1, 0, 2}, {1, 1, 0}, {1, 1, 1}, {1, 2, 0} }; // The (xi, eta, zeta) exponents for each of the dx_deta terms static const int dx_deta_exponents[15][3] = { {0, 0, 0}, {0, 0, 1}, {0, 0, 2}, {0, 0, 3}, {0, 1, 0}, {0, 1, 1}, {0, 1, 2}, {1, 0, 0}, {1, 0, 1}, {1, 0, 2}, {1, 1, 0}, {1, 1, 1}, {2, 0, 0}, {2, 0, 1}, {2, 1, 0} }; // The (xi, eta, zeta) exponents for each of the dx_dzeta terms static const int dx_dzeta_exponents[20][3] = { {0, 0, 0}, {0, 0, 1}, {0, 0, 2}, {0, 0, 3}, {0, 0, 4}, {0, 1, 0}, {0, 1, 1}, {0, 1, 2}, {0, 1, 3}, {1, 0, 0}, {1, 0, 1}, {1, 0, 2}, {1, 0, 3}, {1, 1, 0}, {1, 1, 1}, {1, 2, 0}, {1, 2, 1}, {2, 1, 0}, {2, 1, 1}, {2, 2, 0}, }; // Number of points in the quadrature rule const int N = 27; // Parameters of the quadrature rule static const Real // Parameters used for (xi, eta) quadrature points. a1 = Real(-7.1805574131988893873307823958101e-01L), a2 = Real(-5.0580870785392503961340276902425e-01L), a3 = Real(-2.2850430565396735359574351631314e-01L), // Parameters used for zeta quadrature points. b1 = Real( 7.2994024073149732155837979012003e-02L), b2 = Real( 3.4700376603835188472176354340395e-01L), b3 = Real( 7.0500220988849838312239847758405e-01L), // There are 9 unique weight values since there are three // for each of the three unique zeta values. w1 = Real( 4.8498876871878584357513834016440e-02L), w2 = Real( 4.5137737425884574692441981593901e-02L), w3 = Real( 9.2440441384508327195915094925393e-03L), w4 = Real( 7.7598202995005734972022134426305e-02L), w5 = Real( 7.2220379881415319507907170550242e-02L), w6 = Real( 1.4790470621521332351346415188063e-02L), w7 = Real( 1.2415712479200917595523541508209e-01L), w8 = Real( 1.1555260781026451121265147288039e-01L), w9 = Real( 2.3664752994434131762154264300901e-02L); // The points and weights of the 3x3x3 quadrature rule static const Real xi[N][3] = {// ^0 ^1 ^2 { 1., a1, a1*a1}, { 1., a2, a2*a2}, { 1., a3, a3*a3}, { 1., a1, a1*a1}, { 1., a2, a2*a2}, { 1., a3, a3*a3}, { 1., a1, a1*a1}, { 1., a2, a2*a2}, { 1., a3, a3*a3}, { 1., 0., 0. }, { 1., 0., 0. }, { 1., 0., 0. }, { 1., 0., 0. }, { 1., 0., 0. }, { 1., 0., 0. }, { 1., 0., 0. }, { 1., 0., 0. }, { 1., 0., 0. }, { 1., -a1, a1*a1}, { 1., -a2, a2*a2}, { 1., -a3, a3*a3}, { 1., -a1, a1*a1}, { 1., -a2, a2*a2}, { 1., -a3, a3*a3}, { 1., -a1, a1*a1}, { 1., -a2, a2*a2}, { 1., -a3, a3*a3} }; static const Real eta[N][3] = {// ^0 ^1 ^2 { 1., a1, a1*a1}, { 1., a2, a2*a2}, { 1., a3, a3*a3}, { 1., 0., 0. }, { 1., 0., 0. }, { 1., 0., 0. }, { 1., -a1, a1*a1}, { 1., -a2, a2*a2}, { 1., -a3, a3*a3}, { 1., a1, a1*a1}, { 1., a2, a2*a2}, { 1., a3, a3*a3}, { 1., 0., 0. }, { 1., 0., 0. }, { 1., 0., 0. }, { 1., -a1, a1*a1}, { 1., -a2, a2*a2}, { 1., -a3, a3*a3}, { 1., a1, a1*a1}, { 1., a2, a2*a2}, { 1., a3, a3*a3}, { 1., 0., 0. }, { 1., 0., 0. }, { 1., 0., 0. }, { 1., -a1, a1*a1}, { 1., -a2, a2*a2}, { 1., -a3, a3*a3} }; static const Real zeta[N][5] = {// ^0 ^1 ^2 ^3 ^4 { 1., b1, b1*b1, b1*b1*b1, b1*b1*b1*b1}, { 1., b2, b2*b2, b2*b2*b2, b2*b2*b2*b2}, { 1., b3, b3*b3, b3*b3*b3, b3*b3*b3*b3}, { 1., b1, b1*b1, b1*b1*b1, b1*b1*b1*b1}, { 1., b2, b2*b2, b2*b2*b2, b2*b2*b2*b2}, { 1., b3, b3*b3, b3*b3*b3, b3*b3*b3*b3}, { 1., b1, b1*b1, b1*b1*b1, b1*b1*b1*b1}, { 1., b2, b2*b2, b2*b2*b2, b2*b2*b2*b2}, { 1., b3, b3*b3, b3*b3*b3, b3*b3*b3*b3}, { 1., b1, b1*b1, b1*b1*b1, b1*b1*b1*b1}, { 1., b2, b2*b2, b2*b2*b2, b2*b2*b2*b2}, { 1., b3, b3*b3, b3*b3*b3, b3*b3*b3*b3}, { 1., b1, b1*b1, b1*b1*b1, b1*b1*b1*b1}, { 1., b2, b2*b2, b2*b2*b2, b2*b2*b2*b2}, { 1., b3, b3*b3, b3*b3*b3, b3*b3*b3*b3}, { 1., b1, b1*b1, b1*b1*b1, b1*b1*b1*b1}, { 1., b2, b2*b2, b2*b2*b2, b2*b2*b2*b2}, { 1., b3, b3*b3, b3*b3*b3, b3*b3*b3*b3}, { 1., b1, b1*b1, b1*b1*b1, b1*b1*b1*b1}, { 1., b2, b2*b2, b2*b2*b2, b2*b2*b2*b2}, { 1., b3, b3*b3, b3*b3*b3, b3*b3*b3*b3}, { 1., b1, b1*b1, b1*b1*b1, b1*b1*b1*b1}, { 1., b2, b2*b2, b2*b2*b2, b2*b2*b2*b2}, { 1., b3, b3*b3, b3*b3*b3, b3*b3*b3*b3}, { 1., b1, b1*b1, b1*b1*b1, b1*b1*b1*b1}, { 1., b2, b2*b2, b2*b2*b2, b2*b2*b2*b2}, { 1., b3, b3*b3, b3*b3*b3, b3*b3*b3*b3} }; static const Real w[N] = {w1, w2, w3, w4, w5, w6, // 0-5 w1, w2, w3, w4, w5, w6, // 6-11 w7, w8, w9, w4, w5, w6, // 12-17 w1, w2, w3, w4, w5, w6, // 18-23 w1, w2, w3}; // 24-26 Real vol = 0.; for (int q=0; q<N; ++q) { // Compute denominators for the current q. Real den2 = (1. - zeta[q][1])*(1. - zeta[q][1]), den3 = den2*(1. - zeta[q][1]); // Compute dx/dxi and dx/deta at the current q. Point dx_dxi_q, dx_deta_q; for (int c=0; c<15; ++c) { dx_dxi_q += xi[q][dx_dxi_exponents[c][0]]* eta[q][dx_dxi_exponents[c][1]]* zeta[q][dx_dxi_exponents[c][2]]*dx_dxi[c]; dx_deta_q += xi[q][dx_deta_exponents[c][0]]* eta[q][dx_deta_exponents[c][1]]* zeta[q][dx_deta_exponents[c][2]]*dx_deta[c]; } // Compute dx/dzeta at the current q. Point dx_dzeta_q; for (int c=0; c<20; ++c) { dx_dzeta_q += xi[q][dx_dzeta_exponents[c][0]]* eta[q][dx_dzeta_exponents[c][1]]* zeta[q][dx_dzeta_exponents[c][2]]*dx_dzeta[c]; } // Scale everything appropriately dx_dxi_q /= den2; dx_deta_q /= den2; dx_dzeta_q /= den3; // Compute scalar triple product, multiply by weight, and accumulate volume. vol += w[q] * triple_product(dx_dxi_q, dx_deta_q, dx_dzeta_q); } return vol; }
Real Tet10::volume () const { // Make copies of our points. It makes the subsequent calculations a bit // shorter and avoids dereferencing the same pointer multiple times. Point x0 = point(0), x1 = point(1), x2 = point(2), x3 = point(3), x4 = point(4), x5 = point(5), x6 = point(6), x7 = point(7), x8 = point(8), x9 = point(9); // The constant components of the dx/dxi vector, linear in xi, eta, zeta. // These were copied directly from the output of a Python script. Point dx_dxi[4] = { -3*x0 - x1 + 4*x4, // constant 4*x0 - 4*x4 - 4*x7 + 4*x8, // zeta 4*x0 - 4*x4 + 4*x5 - 4*x6, // eta 4*x0 + 4*x1 - 8*x4 // xi }; // The constant components of the dx/deta vector, linear in xi, eta, zeta. // These were copied directly from the output of a Python script. Point dx_deta[4] = { -3*x0 - x2 + 4*x6, // constant 4*x0 - 4*x6 - 4*x7 + 4*x9, // zeta 4*x0 + 4*x2 - 8*x6, // eta 4*x0 - 4*x4 + 4*x5 - 4*x6 // xi }; // The constant components of the dx/dzeta vector, linear in xi, eta, zeta. // These were copied directly from the output of a Python script. Point dx_dzeta[4] = { -3*x0 - x3 + 4*x7, // constant 4*x0 + 4*x3 - 8*x7, // zeta 4*x0 - 4*x6 - 4*x7 + 4*x9, // eta 4*x0 - 4*x4 - 4*x7 + 4*x8 // xi }; // 2x2x2 conical quadrature rule. Note: there is also a five point // rule for tets with a negative weight which would be cheaper, but // we'll use this one to preclude any possible issues with // cancellation error. const int N = 8; static const Real w[N] = { 3.6979856358852914509238091810505e-02L, 1.6027040598476613723156741868689e-02L, 2.1157006454524061178256145400082e-02L, 9.1694299214797439226823542540576e-03L, 3.6979856358852914509238091810505e-02L, 1.6027040598476613723156741868689e-02L, 2.1157006454524061178256145400082e-02L, 9.1694299214797439226823542540576e-03L }; static const Real xi[N] = { 1.2251482265544137786674043037115e-01L, 5.4415184401122528879992623629551e-01L, 1.2251482265544137786674043037115e-01L, 5.4415184401122528879992623629551e-01L, 1.2251482265544137786674043037115e-01L, 5.4415184401122528879992623629551e-01L, 1.2251482265544137786674043037115e-01L, 5.4415184401122528879992623629551e-01L }; static const Real eta[N] = { 1.3605497680284601717109468420738e-01L, 7.0679724159396903069267439165167e-02L, 5.6593316507280088053551297149570e-01L, 2.9399880063162286589079157179842e-01L, 1.3605497680284601717109468420738e-01L, 7.0679724159396903069267439165167e-02L, 5.6593316507280088053551297149570e-01L, 2.9399880063162286589079157179842e-01L }; static const Real zeta[N] = { 1.5668263733681830907933725249176e-01L, 8.1395667014670255076709592007207e-02L, 6.5838687060044409936029672711329e-02L, 3.4202793236766414300604458388142e-02L, 5.8474756320489429588282763292971e-01L, 3.0377276481470755305409673253211e-01L, 2.4571332521171333166171692542182e-01L, 1.2764656212038543100867773351792e-01L }; Real vol = 0.; for (int q=0; q<N; ++q) { // Compute dx_dxi, dx_deta, dx_dzeta at the current quadrature point. Point dx_dxi_q = dx_dxi[0] + zeta[q]*dx_dxi[1] + eta[q]*dx_dxi[2] + xi[q]*dx_dxi[3], dx_deta_q = dx_deta[0] + zeta[q]*dx_deta[1] + eta[q]*dx_deta[2] + xi[q]*dx_deta[3], dx_dzeta_q = dx_dzeta[0] + zeta[q]*dx_dzeta[1] + eta[q]*dx_dzeta[2] + xi[q]*dx_dzeta[3]; // Compute scalar triple product, multiply by weight, and accumulate volume. vol += w[q] * triple_product(dx_dxi_q, dx_deta_q, dx_dzeta_q); } return vol; }
bool TensorMechanicsPlasticMohrCoulombMulti::returnTip(const std::vector<Real> & eigvals, const std::vector<RealVectorValue> & n, std::vector<Real> & dpm, RankTwoTensor & returned_stress, Real intnl_old, Real & sinphi, Real & cohcos, Real initial_guess, bool & nr_converged, Real ep_plastic_tolerance, std::vector<Real> & yf) const { // This returns to the Mohr-Coulomb tip using the THREE directions // given in n, and yields the THREE dpm values. Note that you // must supply THREE suitable n vectors out of the total of SIX // flow directions, and then interpret the THREE dpm values appropriately. // // Eg1. You supply the flow directions n[0], n[1] and n[3] as // the "n" vectors. This is return-to-the-tip via 110100. // Then the three returned dpm values will be dpm[0], dpm[1] and dpm[3]. // Eg2. You supply the flow directions n[1], n[3] and n[5] as // the "n" vectors. This is return-to-the-tip via 010101. // Then the three returned dpm values will be dpm[1], dpm[3] and dpm[5]. // The returned point is defined by the three yield functions (corresonding // to the three supplied flow directions) all being zero. // that is, returned_stress = diag(cohcot, cohcot, cohcot), where // cohcot = cohesion*cosphi/sinphi // where intnl = intnl_old + dpm[0] + dpm[1] + dpm[2] // The 3 plastic multipliers, dpm, are defiend by the normality condition // eigvals - cohcot = dpm[0]*n[0] + dpm[1]*n[1] + dpm[2]*n[2] // (Kuhn-Tucker demands that all dpm are non-negative, but we leave // that checking for the end.) // (Remember that these "n" vectors and "dpm" values must be interpreted // differently depending on what you pass into this function.) // This is a vector equation with solution (A): // dpm[0] = triple(eigvals - cohcot, n[1], n[2])/trip; // dpm[1] = triple(eigvals - cohcot, n[2], n[0])/trip; // dpm[2] = triple(eigvals - cohcot, n[0], n[1])/trip; // where trip = triple(n[0], n[1], n[2]). // By adding the three components of that solution together // we can get an equation for x = dpm[0] + dpm[1] + dpm[2], // and then our Newton-Raphson only involves one variable (x). // In the following, i specialise to the isotropic situation. mooseAssert(n.size() == 3, "TensorMechanicsPlasticMohrCoulombMulti: Custom tip-return algorithm must be supplied with n of size 3, whereas yours is " << n.size()); mooseAssert(dpm.size() == 3, "TensorMechanicsPlasticMohrCoulombMulti: Custom tip-return algorithm must be supplied with dpm of size 3, whereas yours is " << dpm.size()); mooseAssert(yf.size() == 6, "TensorMechanicsPlasticMohrCoulombMulti: Custom tip-return algorithm must be supplied with yf of size 6, whereas yours is " << yf.size()); Real x = initial_guess; const Real trip = triple_product(n[0], n[1], n[2]); sinphi = std::sin(phi(intnl_old + x)); Real cosphi = std::cos(phi(intnl_old + x)); Real coh = cohesion(intnl_old + x); cohcos = coh*cosphi; Real cohcot = cohcos/sinphi; if (_cohesion.modelName().compare("Constant") != 0 || _phi.modelName().compare("Constant") != 0) { // Finding x is expensive. Therefore // although x!=0 for Constant Hardening, solution (A) // demonstrates that we don't // actually need to know x to find the dpm for // Constant Hardening. // // However, for nontrivial Hardening, the following // is necessary // cohcot_coeff = [1,1,1].(Cross[n[1], n[2]] + Cross[n[2], n[0]] + Cross[n[0], n[1]])/trip Real cohcot_coeff = (n[0](0)*(n[1](1) - n[1](2) - n[2](1)) + (n[1](2) - n[1](1))*n[2](0) + (n[1](0) - n[1](2))*n[2](1) + n[0](2)*(n[1](0) - n[1](1) - n[2](0) + n[2](1)) + n[0](1)*(n[1](2) - n[1](0) + n[2](0) - n[2](2)) + (n[0](0) - n[1](0) + n[1](1))*n[2](2))/trip; // eig_term = eigvals.(Cross[n[1], n[2]] + Cross[n[2], n[0]] + Cross[n[0], n[1]])/trip Real eig_term = eigvals[0]*(-n[0](2)*n[1](1) + n[0](1)*n[1](2) + n[0](2)*n[2](1) - n[1](2)*n[2](1) - n[0](1)*n[2](2) + n[1](1)*n[2](2))/trip; eig_term += eigvals[1]*(n[0](2)*n[1](0) - n[0](0)*n[1](2) - n[0](2)*n[2](0) + n[1](2)*n[2](0) + n[0](0)*n[2](2) - n[1](0)*n[2](2))/trip; eig_term += eigvals[2]*(n[0](0)*n[1](1) - n[1](1)*n[2](0) + n[0](1)*n[2](0) - n[0](1)*n[1](0) - n[0](0)*n[2](1) + n[1](0)*n[2](1))/trip; // and finally, the equation we want to solve is: // x - eig_term + cohcot*cohcot_coeff = 0 // but i divide by cohcot_coeff so the result has the units of // stress, so using _f_tol as a convergence check is reasonable eig_term /= cohcot_coeff; Real residual = x/cohcot_coeff - eig_term + cohcot; Real jacobian; Real deriv_phi; Real deriv_coh; unsigned int iter = 0; do { deriv_phi = dphi(intnl_old + x); deriv_coh = dcohesion(intnl_old + x); jacobian = 1.0/cohcot_coeff + deriv_coh*cosphi/sinphi - coh*deriv_phi/Utility::pow<2>(sinphi); x += -residual/jacobian; if (iter > _max_iters) // not converging { nr_converged = false; return false; } sinphi = std::sin(phi(intnl_old + x)); cosphi = std::cos(phi(intnl_old + x)); coh = cohesion(intnl_old + x); cohcos = coh*cosphi; cohcot = cohcos/sinphi; residual = x/cohcot_coeff - eig_term + cohcot; iter ++; } while (residual*residual > _f_tol*_f_tol/100); } // so the NR process converged, but we must // calculate the individual dpm values and // check Kuhn-Tucker nr_converged = true; if (x < -3*ep_plastic_tolerance) // obviously at least one of the dpm are < -ep_plastic_tolerance. No point in proceeding. This is a potential weak-point: if the user has set _f_tol quite large, and ep_plastic_tolerance quite small, the above NR process will quickly converge, but the solution may be wrong and violate Kuhn-Tucker. return false; // The following is the solution (A) written above // (dpm[0] = triple(eigvals - cohcot, n[1], n[2])/trip, etc) // in the isotropic situation RealVectorValue v; v(0) = eigvals[0] - cohcot; v(1) = eigvals[1] - cohcot; v(2) = eigvals[2] - cohcot; dpm[0] = triple_product(v, n[1], n[2])/trip; dpm[1] = triple_product(v, n[2], n[0])/trip; dpm[2] = triple_product(v, n[0], n[1])/trip; if (dpm[0] < -ep_plastic_tolerance || dpm[1] < -ep_plastic_tolerance || dpm[2] < -ep_plastic_tolerance) // Kuhn-Tucker failure. No point in proceeding return false; // Kuhn-Tucker has succeeded: just need returned_stress and yf values // I do not use the dpm to calculate returned_stress, because that // might add error (and computational effort), simply: returned_stress(0, 0) = returned_stress(1, 1) = returned_stress(2, 2) = cohcot; // So by construction the yield functions are all zero yf[0] = yf[1] = yf[2] = yf[3] = yf[4] = yf[5] = 0; return true; }
Real Prism6::volume () const { // Make copies of our points. It makes the subsequent calculations a bit // shorter and avoids dereferencing the same pointer multiple times. Point x0 = point(0), x1 = point(1), x2 = point(2), x3 = point(3), x4 = point(4), x5 = point(5); // constant and zeta terms only. These are copied directly from a // Python script. Point dx_dxi[2] = { -x0/2 + x1/2 - x3/2 + x4/2, // constant x0/2 - x1/2 - x3/2 + x4/2, // zeta }; // constant and zeta terms only. These are copied directly from a // Python script. Point dx_deta[2] = { -x0/2 + x2/2 - x3/2 + x5/2, // constant x0/2 - x2/2 - x3/2 + x5/2, // zeta }; // Constant, xi, and eta terms Point dx_dzeta[3] = { -x0/2 + x3/2, // constant x0/2 - x2/2 - x3/2 + x5/2, // eta x0/2 - x1/2 - x3/2 + x4/2 // xi }; // The quadrature rule the Prism6 is a tensor product between a // four-point TRI3 rule (in xi, eta) and a two-point EDGE2 rule (in // zeta) which is capable of integrating cubics exactly. // Number of points in the 2D quadrature rule. const int N2D = 4; static const Real w2D[N2D] = { 1.5902069087198858469718450103758e-01L, 9.0979309128011415302815498962418e-02L, 1.5902069087198858469718450103758e-01L, 9.0979309128011415302815498962418e-02L }; static const Real xi[N2D] = { 1.5505102572168219018027159252941e-01L, 6.4494897427831780981972840747059e-01L, 1.5505102572168219018027159252941e-01L, 6.4494897427831780981972840747059e-01L }; static const Real eta[N2D] = { 1.7855872826361642311703513337422e-01L, 7.5031110222608118177475598324603e-02L, 6.6639024601470138670269327409637e-01L, 2.8001991549907407200279599420481e-01L }; // Number of points in the 1D quadrature rule. The weights of the // 1D quadrature rule are equal to 1. const int N1D = 2; // Points of the 1D quadrature rule static const Real zeta[N1D] = { -std::sqrt(3.)/3, std::sqrt(3.)/3. }; Real vol = 0.; for (int i=0; i<N2D; ++i) { // dx_dzeta depends only on the 2D quadrature rule points. Point dx_dzeta_q = dx_dzeta[0] + eta[i]*dx_dzeta[1] + xi[i]*dx_dzeta[2]; for (int j=0; j<N1D; ++j) { // dx_dxi and dx_deta only depend on the 1D quadrature rule points. Point dx_dxi_q = dx_dxi[0] + zeta[j]*dx_dxi[1], dx_deta_q = dx_deta[0] + zeta[j]*dx_deta[1]; // Compute scalar triple product, multiply by weight, and accumulate volume. vol += w2D[i] * triple_product(dx_dxi_q, dx_deta_q, dx_dzeta_q); } } return vol; }
void cut_leaf (struct msscene *ms, struct surface *current_surface, double fine_pixel, struct leaf *lf) { int j, k, n, nx, nw, orn, sgn, i, used0, used1, nused; double t, xsmin, sangle, fsgn; double base[3], xpnts[2][3], vect[3]; char message[MAXLINE]; struct circle *cir; struct arc *al, *a; struct variety *vty; struct face *fac; struct leaf *lfcut; struct lax *xsptr, *xsend, *xsp; struct cycle *cyc; struct edge *edg; struct lax laxes[MAX_LAX]; struct cept *ex; fac = lf -> fac; /* count number of edges for face */ n = 2 * edges_in_face (fac); /* if toroidal, no cutting */ if (fac -> shape == SADDLE || n <= 0) { /* clip and render leaf */ clip_leaf (ms, current_surface, fine_pixel, lf); return; } /* determine sign (plus for convex, minus for concave) */ fsgn = ((lf -> shape == CONVEX) ? 1.0 : -1.0); /* create arc for leaf */ al = leaf_to_arc (lf); vty = fac -> vty; cir = lf -> cir; for (k = 0; k < 3; k++) base[k] = lf -> ends[0][k] - cir -> center[k]; normalize (base); /* determine intersection points between leaf and arcs of face */ xsptr = &(laxes[0]); for (cyc = fac -> first_cycle; cyc != NULL; cyc = cyc -> next) for (edg = cyc -> first_edge; edg != NULL; edg = edg -> next) { a = edg -> arcptr; orn = edg -> orn; sgn = 1 - 2 * orn; /* call arc-arc intersection function */ nx = arc_arc (al, a, xpnts); if (nx <= 0) continue; /* no intersections */ /* fill data into list */ for (i = 0; i < nx; i++) { if (xsptr - &(laxes[0]) >= n) { ex = new_cept (ARRAY_ERROR, MSOVERFLOW, FATAL_SEVERITY); add_function (ex, "cut_leaf"); add_source (ex, "msrender.c"); add_long (ex, "maximum number of leaf-arc intersections", n); return; } for (k = 0; k < 3; k++) { xsptr -> co[k] = xpnts[i][k]; vect[k] = fsgn * (xsptr -> co[k] - vty -> center[k]); } t = sgn * triple_product (al -> cir -> axis, a -> cir -> axis, vect); xsptr -> ent = (t < 0); xsptr -> used = 0; xsptr++; } } xsend = xsptr; /* mark end of list of intersection points */ nx = xsend - &(laxes[0]); /* compute number of intersection points */ /* intersection parity error checking */ if ((lf -> where[0] == lf -> where[1]) == (nx % 2)) { /* parity check fails, flag leaf to check accessibility of every pixel of leaf */ ms -> n_bdy_parity++; lf -> cep = 1; clip_leaf (ms, current_surface, fine_pixel, lf); lf -> cep = 0; frearc (al); return; } /* kludge to work around undiscovered bug for quartet cusps */ if (lf -> atmnum[3] != 0) { /* flag leaf to check accessibility of every pixel of leaf */ lf -> cep = 1; clip_leaf (ms, current_surface, fine_pixel, lf); lf -> cep = 0; frearc (al); return; } /* check for no intersection points */ if (nx == 0) { if (lf -> where[0] == INACCESSIBLE) { /* entire leaf inaccessible: nothing to render */ /* free temporary memory */ frearc (al); return; } /* entire leaf accessible: render entire leaf */ clip_leaf (ms, current_surface, fine_pixel, lf); if (error()) return; /* free temporary memory */ frearc (al); return; } /* calculate angles for transitions between accessible and inaccessible */ for (xsptr = &(laxes[0]); xsptr < xsend; xsptr++) { for (k = 0; k < 3; k++) vect[k] = xsptr -> co[k] - cir -> center[k]; normalize (vect); xsptr -> angle = positive_angle (base, vect, cir -> axis); } /* initialization */ used0 = 0; used1 = 0; nused = 0; nw = 0; /* count number of accessible endpoints of originial leaf */ for (j = 0; j < 2; j++) if (lf -> where[j]) nw++; /* create new leaves */ while (nused < nx + nw) { /* get starting point */ if (!used0 && lf -> where[0] == ACCESSIBLE) { /* start at original endpoint */ used0 = 1; /* first originial endpoint used */ nused++; /* increment number used */ lfcut = duplicate_leaf (lf); /* duplicate leaf */ sangle = 0.0; /* starting angle */ } else { /* look for unused cut vertex */ xsmin = 4 * PI; xsp = NULL; for (xsptr = &(laxes[0]); xsptr < xsend; xsptr++) { if (xsptr -> used) continue; if (xsptr -> angle < xsmin) { xsmin = xsptr -> angle; xsp = xsptr; } } if (xsp == NULL) { if (lf -> where[1] == INACCESSIBLE) break; ms -> n_missing_leaves++; return; } if (!xsp -> ent) { ms -> n_missing_leaves++; return; } xsp -> used = 1; /* mark as used */ nused++; /* increment number used */ sangle = xsp -> angle; /* starting angle */ lfcut = duplicate_leaf (lf); /* duplicate leaf */ for (k = 0; k < 3; k++) lfcut -> ends[0][k] = xsp -> co[k]; } /* get ending point */ /* initialization */ xsmin = 4 * PI; xsp = NULL; for (xsptr = &(laxes[0]); xsptr < xsend; xsptr++) { if (xsptr -> used) continue; /* already used */ if (xsptr -> angle < sangle) continue; /* end after start */ if (xsptr -> angle < xsmin) { /* best so far, save */ xsmin = xsptr -> angle; xsp = xsptr; } } if (xsp == NULL) { if (used1 || lf -> where[1] == INACCESSIBLE) { ms -> n_missing_leaves++; return; } /* use East */ used1 = 1; /* mark East original endpoint used */ nused++; /* increment number used */ for (k = 0; k < 3; k++) lfcut -> ends[1][k] = lf -> ends[1][k]; } else { /* use cut point */ for (k = 0; k < 3; k++) lfcut -> ends[1][k] = xsp -> co[k]; xsp -> used = 1; /* mark as used */ nused++; /* increment number used */ } /* we have a cut leaf; clip and render leaf */ clip_leaf (ms, current_surface, fine_pixel, lfcut); free_leaf (lfcut); /* free temporary memory */ } /* free temporary memory */ frearc (al); return; }
void clip_leaf (struct msscene *ms, struct surface *current_surface, double fine_pixel, struct leaf *lf) { int j, k, n, nx, nw, i, used0, used1, nused; double t, xsmin, sangle; double base[3], xpnts[2][3], vect[3]; double clip_center[3], clip_axis[3]; struct circle *cir; struct arc *al; struct leaf *lfcut; struct lax *xsptr, *xsend, *xsp; struct lax laxes[MAX_LAX]; /* count double number of clipping planes */ n = 0; if (ms -> clipping) n += 2; if (n <= 0 || !(current_surface -> clipping)) { lf -> clip_ep = 0; render_leaf (ms, current_surface, fine_pixel, lf); return; } /* determine accessibility of endpoints of leaf */ for (j = 0; j < 2; j++) lf -> when[j] = !(current_surface -> clipping && clipped (ms, lf -> ends[j])); /* create arc for leaf */ al = leaf_to_arc (lf); cir = lf -> cir; for (k = 0; k < 3; k++) base[k] = lf -> ends[0][k] - cir -> center[k]; normalize (base); /* determine intersection points between leaf and clipping planes */ xsptr = &(laxes[0]); if (ms -> clipping) { for (k = 0; k < 3; k++) { clip_center[k] = ms -> clip_center[k]; clip_axis[k] = ms -> clip_axis[k]; } /* call arc-plane intersection function */ nx = arc_plane (al, clip_center, clip_axis, xpnts); if (nx < 0) return; /* if (nx == 0) continue; no intersections */ /* fill data into list */ /* figure out enter or exit */ for (i = 0; i < nx; i++) { if (xsptr - &(laxes[0]) >= n) { ms -> n_missing_leaves++; return; } for (k = 0; k < 3; k++) xsptr -> co[k] = xpnts[i][k]; for (k = 0; k < 3; k++) vect[k] = xsptr -> co[k] - al -> cir -> center[k]; t = triple_product (al -> cir -> axis, clip_axis, vect); xsptr -> ent = (t > 0); /* or reverse? */ xsptr -> used = 0; xsptr++; } } xsend = xsptr; /* mark end of list of intersection points */ nx = xsend - &(laxes[0]); /* compute number of intersection points */ /* error checking */ if ((lf -> when[0] == lf -> when[1]) == (nx % 2)) { /* parity check fails, flag leaf to check clipping of every pixel of leaf */ ms -> n_clip_parity++; lf -> clip_ep = 1; render_leaf (ms, current_surface, fine_pixel, lf); if (error()) return; lf -> clip_ep = 0; frearc (al); return; } /* check for no intersection points */ if (nx == 0) { if (lf -> when[0] == CLIPPED) { /* entire leaf clipped: nothing to render */ /* free temporary memory */ frearc (al); return; } /* entire leaf unclipped: render entire leaf */ render_leaf (ms, current_surface, fine_pixel, lf); if (error()) return; /* free temporary memory */ frearc (al); return; } /* calculate angles for transitions between unclipped and clipped */ for (xsptr = &(laxes[0]); xsptr < xsend; xsptr++) { for (k = 0; k < 3; k++) vect[k] = xsptr -> co[k] - cir -> center[k]; normalize (vect); xsptr -> angle = positive_angle (base, vect, cir -> axis); } /* initialization */ used0 = 0; used1 = 0; nused = 0; nw = 0; /* count number of unclipped endpoints of originial leaf */ for (j = 0; j < 2; j++) if (lf -> when[j]) nw++; /* create new leaves */ while (nused < nx + nw) { /* get starting point */ if (!used0 && lf -> when[0] == UNCLIPPED) { /* start at original endpoint */ used0 = 1; /* first originial endpoint used */ nused++; /* increment number used */ lfcut = duplicate_leaf (lf); /* duplicate leaf */ sangle = 0.0; /* starting angle */ } else { /* look for unused cut vertex */ xsmin = 4 * PI; xsp = NULL; for (xsptr = &(laxes[0]); xsptr < xsend; xsptr++) { if (xsptr -> used) continue; if (xsptr -> angle < xsmin) { xsmin = xsptr -> angle; xsp = xsptr; } } if (xsp == NULL) { if (lf -> when[1] == CLIPPED) break; ms -> n_missing_leaves++; return; } if (!xsp -> ent) { ms -> n_missing_leaves++; return; } xsp -> used = 1; /* mark as used */ nused++; /* increment number used */ sangle = xsp -> angle; /* starting angle */ lfcut = duplicate_leaf (lf); /* duplicate leaf */ for (k = 0; k < 3; k++) lfcut -> ends[0][k] = xsp -> co[k]; } /* get ending point */ /* initialization */ xsmin = 4 * PI; xsp = NULL; for (xsptr = &(laxes[0]); xsptr < xsend; xsptr++) { if (xsptr -> used) continue; /* already used */ if (xsptr -> angle < sangle) continue; /* end after start */ if (xsptr -> angle < xsmin) { /* best so far, save */ xsmin = xsptr -> angle; xsp = xsptr; } } if (xsp == NULL) { if (used1 || lf -> when[1] == CLIPPED) { ms -> n_missing_leaves++; return; } /* use East */ used1 = 1; /* mark East original endpoint used */ nused++; /* increment number used */ for (k = 0; k < 3; k++) lfcut -> ends[1][k] = lf -> ends[1][k]; } else { /* use cut point */ for (k = 0; k < 3; k++) lfcut -> ends[1][k] = xsp -> co[k]; xsp -> used = 1; /* mark as used */ nused++; /* increment number used */ } /* we have a cut leaf */ /* clip and render leaf */ render_leaf (ms, current_surface, fine_pixel, lfcut); if (error()) return; free_leaf (lfcut); /* free temporary memory */ } /* free temporary memory */ frearc (al); return; }