void
MAST::MindlinBendingOperator::
_transverse_shear_operations(const std::vector<std::vector<Real> >& phi,
                             const std::vector<std::vector<libMesh::RealVectorValue> >& dphi,
                             const std::vector<Real>& JxW,
                             const unsigned int     qp,
                             const RealMatrixX&     material,
                             FEMOperatorMatrix&     Bmat,
                             RealVectorX&           phi_vec,
                             RealVectorX&           vec_n2,
                             RealVectorX&           vec_2,
                             RealMatrixX&           mat_n2n2,
                             RealMatrixX&           mat_2n2,
                             bool                   request_jacobian,
                             RealVectorX&           local_f,
                             RealMatrixX&           local_jac) {
    
    // initialize the strain operator
    for ( unsigned int i_nd=0; i_nd<phi.size(); i_nd++ )
        phi_vec(i_nd) = dphi[i_nd][qp](0);  // dphi/dx
    
    Bmat.set_shape_function(0, 2, phi_vec); // gamma-xz:  w
    
    for ( unsigned int i_nd=0; i_nd<phi.size(); i_nd++ )
        phi_vec(i_nd) = dphi[i_nd][qp](1);  // dphi/dy
    
    Bmat.set_shape_function(1, 2, phi_vec); // gamma-yz : w
    
    for ( unsigned int i_nd=0; i_nd<phi.size(); i_nd++ )
        phi_vec(i_nd) = phi[i_nd][qp];  // phi
    
    Bmat.set_shape_function(0, 4, phi_vec); // gamma-xz:  thetay
    phi_vec  *= -1.0;
    Bmat.set_shape_function(1, 3, phi_vec); // gamma-yz : thetax
    
    
    // now add the transverse shear component
    Bmat.vector_mult(vec_2, _structural_elem.local_solution());
    vec_2 = material * vec_2;
    Bmat.vector_mult_transpose(vec_n2, vec_2);
    local_f += JxW[qp] * vec_n2;
    
    if (request_jacobian) {
        
        // now add the transverse shear component
        Bmat.left_multiply(mat_2n2, material);
        Bmat.right_multiply_transpose(mat_n2n2, mat_2n2);
        local_jac += JxW[qp] * mat_n2n2;
    }
}
bool
MAST::StructuralElementBase::inertial_force (bool request_jacobian,
                                             DenseRealVector& f,
                                             DenseRealMatrix& jac)
{
    FEMOperatorMatrix Bmat;
    
    const std::vector<Real>& JxW = _fe->get_JxW();
    const std::vector<libMesh::Point>& xyz = _fe->get_xyz();
    const std::vector<std::vector<Real> >& phi = _fe->get_phi();
    const unsigned int n_phi = (unsigned int)phi.size(), n1=6, n2=6*n_phi;
    
    DenseRealMatrix material_mat, mat1_n1n2, mat2_n2n2, local_jac;
    DenseRealVector  phi_vec, vec1_n1, vec2_n2, local_f;
    mat1_n1n2.resize(n1, n2); mat2_n2n2.resize(n2, n2); local_jac.resize(n2, n2);
    phi_vec.resize(n_phi); vec1_n1.resize(n1); vec2_n2.resize(n2);
    local_f.resize(n2);
    
    std::auto_ptr<MAST::FieldFunction<DenseRealMatrix > > mat_inertia
    (_property.get_property(MAST::SECTION_INTEGRATED_MATERIAL_INERTIA_MATRIX,
                            *this).release());

    if (_property.if_diagonal_mass_matrix()) {
        // as an approximation, get matrix at the first quadrature point
        (*mat_inertia)(xyz[0], _system.time, material_mat);
        
        Real vol = 0.;
        const unsigned int nshp = _fe->n_shape_functions();
        for (unsigned int i=0; i<JxW.size(); i++)
            vol += JxW[i];
        vol /= (1.* nshp);
        for (unsigned int i_var=0; i_var<6; i_var++)
            for (unsigned int i=0; i<nshp; i++)
                local_jac(i_var*nshp+i, i_var*nshp+i) =
                vol*material_mat(i_var, i_var);
        local_jac.vector_mult(local_f, local_acceleration);
    }
    else {
        libMesh::Point p;
        
        for (unsigned int qp=0; qp<JxW.size(); qp++) {
            
            this->global_coordinates(xyz[qp], p);
            
            (*mat_inertia)(p, _system.time, material_mat);
            
            // now set the shape function values
            for ( unsigned int i_nd=0; i_nd<n_phi; i_nd++ )
                phi_vec(i_nd) = phi[i_nd][qp];
            
            Bmat.reinit(_system.n_vars(), phi_vec);
            
            Bmat.left_multiply(mat1_n1n2, material_mat);
            
            mat1_n1n2.vector_mult(vec1_n1, local_acceleration);
            Bmat.vector_mult_transpose(vec2_n2, vec1_n1);
            local_f.add(JxW[qp], vec2_n2);
            
            if (request_jacobian) {
                
                Bmat.right_multiply_transpose(mat2_n2n2,
                                              mat1_n1n2);
                local_jac.add(JxW[qp], mat2_n2n2);
            }
            
        }
    }
    
    // now transform to the global coorodinate system
    if (_elem.dim() < 3) {
        transform_to_global_system(local_f, vec2_n2);
        f.add(1., vec2_n2);
        if (request_jacobian) {
            transform_to_global_system(local_jac, mat2_n2n2);
            jac.add(1., mat2_n2n2);
        }
    }
    else {
        f.add(1., local_f);
        if (request_jacobian)
            jac.add(1., local_jac);
    }

    return request_jacobian;
}
void
MAST::StructuralElement2D::
_linearized_geometric_stiffness_sensitivity_with_static_solution
(const unsigned int n2,
 const unsigned int qp,
 const std::vector<Real>& JxW,
 RealMatrixX& local_jac,
 FEMOperatorMatrix& Bmat_mem,
 FEMOperatorMatrix& Bmat_bend,
 FEMOperatorMatrix& Bmat_vk,
 RealMatrixX& stress_l,
 RealMatrixX& vk_dwdxi_mat,
 RealMatrixX& material_A_mat,
 RealMatrixX& material_B_mat,
 RealVectorX& vec1_n1,
 RealVectorX& vec2_n1,
 RealMatrixX& mat1_n1n2,
 RealMatrixX& mat2_n2n2,
 RealMatrixX& mat3)
{
    this->initialize_direct_strain_operator(qp, Bmat_mem);
    _bending_operator->initialize_bending_strain_operator(qp, Bmat_bend);
    
    // first handle constant throught the thickness stresses: membrane and vonKarman
    Bmat_mem.vector_mult(vec1_n1, _local_sol_sens);
    vec2_n1 = material_A_mat * vec1_n1; // linear direct stress
    
    // copy the stress values to a matrix
    stress_l(0,0) = vec2_n1(0); // sigma_xx
    stress_l(0,1) = vec2_n1(2); // sigma_xy
    stress_l(1,0) = vec2_n1(2); // sigma_yx
    stress_l(1,1) = vec2_n1(1); // sigma_yy
    
    // get the von Karman operator matrix
    this->initialize_von_karman_strain_operator(qp,
                                                vec2_n1, // epsilon_vk
                                                vk_dwdxi_mat,
                                                Bmat_vk);
    
    // sensitivity of the vk_dwdxi matrix due to solution sensitivity
    this->initialize_von_karman_strain_operator_sensitivity(qp,
                                                            vk_dwdxi_mat);
    
    // membrane - vk
    mat3 = RealMatrixX::Zero(vk_dwdxi_mat.rows(), n2);
    Bmat_vk.left_multiply(mat3, vk_dwdxi_mat);
    mat3 = material_A_mat * mat3;
    Bmat_mem.right_multiply_transpose(mat2_n2n2, mat3);
    local_jac += JxW[qp] * mat2_n2n2;
    
    // vk - membrane
    Bmat_mem.left_multiply(mat1_n1n2, material_A_mat);
    mat3 = vk_dwdxi_mat.transpose() * mat1_n1n2;
    Bmat_vk.right_multiply_transpose(mat2_n2n2, mat3);
    local_jac += JxW[qp] * mat2_n2n2;
    
    // vk - vk: first order term
    mat3 = RealMatrixX::Zero(2, n2);
    Bmat_vk.left_multiply(mat3, stress_l);
    Bmat_vk.right_multiply_transpose(mat2_n2n2, mat3);
    local_jac += JxW[qp] * mat2_n2n2;
    
    // bending - vk
    mat3 = RealMatrixX::Zero(vk_dwdxi_mat.rows(), n2);
    Bmat_vk.left_multiply(mat3, vk_dwdxi_mat);
    mat3 = material_B_mat.transpose() * mat3;
    Bmat_bend.right_multiply_transpose(mat2_n2n2, mat3);
    local_jac += JxW[qp] * mat2_n2n2;
    
    // vk - bending
    Bmat_bend.left_multiply(mat1_n1n2, material_B_mat);
    mat3 =  vk_dwdxi_mat.transpose() * mat1_n1n2;
    Bmat_vk.right_multiply_transpose(mat2_n2n2, mat3);
    local_jac += JxW[qp] * mat2_n2n2;
}
bool
MAST::StructuralElementBase::inertial_force_sensitivity(bool request_jacobian,
                                                        DenseRealVector& f,
                                                        DenseRealMatrix& jac)
{
    // this should be true if the function is called
    libmesh_assert(this->sensitivity_param);
    libmesh_assert(!this->sensitivity_param->is_shape_parameter()); // this is not implemented for now
    
    
    // check if the material property or the provided exterior
    // values, like temperature, are functions of the sensitivity parameter
    bool calculate = false;
    calculate = calculate || _property.depends_on(*(this->sensitivity_param));
    
    // nothing to be calculated if the element does not depend on the
    // sensitivity parameter.
    if (!calculate)
        return false;
    
    FEMOperatorMatrix Bmat;
    
    const std::vector<Real>& JxW = _fe->get_JxW();
    const std::vector<libMesh::Point>& xyz = _fe->get_xyz();
    const std::vector<std::vector<Real> >& phi = _fe->get_phi();
    const unsigned int n_phi = (unsigned int)phi.size(), n1=6, n2=6*n_phi;
    
    DenseRealMatrix material_mat, mat1_n1n2, mat2_n2n2, local_jac;
    DenseRealVector  phi_vec, vec1_n1, vec2_n2, local_f;
    mat1_n1n2.resize(n1, n2); mat2_n2n2.resize(n2, n2); local_jac.resize(n2, n2);
    phi_vec.resize(n_phi); vec1_n1.resize(n1); vec2_n2.resize(n2);
    local_f.resize(n2);
    
    std::auto_ptr<MAST::FieldFunction<DenseRealMatrix > > mat_inertia
    (_property.get_property(MAST::SECTION_INTEGRATED_MATERIAL_INERTIA_MATRIX,
                            *this).release());
    
    if (_property.if_diagonal_mass_matrix()) {
        
        mat_inertia->total(*this->sensitivity_param,
                           xyz[0], _system.time, material_mat);
        
        Real vol = 0.;
        const unsigned int nshp = _fe->n_shape_functions();
        for (unsigned int i=0; i<JxW.size(); i++)
            vol += JxW[i];
        vol /= (1.* nshp);
        for (unsigned int i_var=0; i_var<6; i_var++)
            for (unsigned int i=0; i<nshp; i++)
                local_jac(i_var*nshp+i, i_var*nshp+i) =
                vol*material_mat(i_var, i_var);
        local_jac.vector_mult(local_f, local_acceleration);
    }
    else {
        libMesh::Point p;
        
        for (unsigned int qp=0; qp<JxW.size(); qp++) {
            
            this->global_coordinates(xyz[qp], p);
            
            mat_inertia->total(*this->sensitivity_param,
                               p, _system.time, material_mat);

            // now set the shape function values
            for ( unsigned int i_nd=0; i_nd<n_phi; i_nd++ )
                phi_vec(i_nd) = phi[i_nd][qp];
            
            Bmat.reinit(_system.n_vars(), phi_vec);
            
            Bmat.left_multiply(mat1_n1n2, material_mat);
            
            mat1_n1n2.vector_mult(vec1_n1, local_acceleration);
            Bmat.vector_mult_transpose(vec2_n2, vec1_n1);
            local_f.add(JxW[qp], vec2_n2);
            
            if (request_jacobian) {
                
                Bmat.right_multiply_transpose(mat2_n2n2,
                                              mat1_n1n2);
                local_jac.add(JxW[qp], mat2_n2n2);
            }
            
        }
    }
    
    // now transform to the global coorodinate system
    if (_elem.dim() < 3) {
        transform_to_global_system(local_f, vec2_n2);
        f.add(1., vec2_n2);
        if (request_jacobian) {
            transform_to_global_system(local_jac, mat2_n2n2);
            jac.add(1., mat2_n2n2);
        }
    }
    else {
        f.add(1., local_f);
        if (request_jacobian)
            jac.add(1., local_jac);
    }
    
    return request_jacobian;
}
void
MAST::StructuralElement2D::_internal_residual_operation
(bool if_bending,
 bool if_vk,
 const unsigned int n2,
 const unsigned int qp,
 const std::vector<Real>& JxW,
 bool request_jacobian,
 bool if_ignore_ho_jac,
 RealVectorX& local_f,
 RealMatrixX& local_jac,
 FEMOperatorMatrix& Bmat_mem,
 FEMOperatorMatrix& Bmat_bend,
 FEMOperatorMatrix& Bmat_vk,
 RealMatrixX& stress,
 RealMatrixX& stress_l,
 RealMatrixX& vk_dwdxi_mat,
 RealMatrixX& material_A_mat,
 RealMatrixX& material_B_mat,
 RealMatrixX& material_D_mat,
 RealVectorX& vec1_n1,
 RealVectorX& vec2_n1,
 RealVectorX& vec3_n2,
 RealVectorX& vec4_2,
 RealVectorX& vec5_2,
 RealMatrixX& mat1_n1n2,
 RealMatrixX& mat2_n2n2,
 RealMatrixX& mat3,
 RealMatrixX& mat4_2n2)
{
    this->initialize_direct_strain_operator(qp, Bmat_mem);
    
    // first handle constant throught the thickness stresses: membrane and vonKarman
    Bmat_mem.vector_mult(vec1_n1, _local_sol);
    vec2_n1 = material_A_mat * vec1_n1; // linear direct stress
    
    // copy the stress values to a matrix
    stress_l(0,0) = vec2_n1(0); // sigma_xx
    stress_l(0,1) = vec2_n1(2); // sigma_xy
    stress_l(1,0) = vec2_n1(2); // sigma_yx
    stress_l(1,1) = vec2_n1(1); // sigma_yy
    
    stress = stress_l;
    
    // get the bending strain operator
    vec2_n1.setConstant(0.); // used to store vk strain, if applicable
    if (if_bending) {
        _bending_operator->initialize_bending_strain_operator(qp, Bmat_bend);
        
        Bmat_bend.vector_mult(vec2_n1, _local_sol);
        vec1_n1 = material_B_mat * vec2_n1;
        stress_l(0,0) += vec2_n1(0); // sigma_xx
        stress_l(0,1) += vec2_n1(2); // sigma_xy
        stress_l(1,0) += vec2_n1(2); // sigma_yx
        stress_l(1,1) += vec2_n1(1); // sigma_yy
        
        stress(0,0) += vec2_n1(0); // sigma_xx
        stress(0,1) += vec2_n1(2); // sigma_xy
        stress(1,0) += vec2_n1(2); // sigma_yx
        stress(1,1) += vec2_n1(1); // sigma_yy
        
        
        if (if_vk)  { // get the vonKarman strain operator if needed
            this->initialize_von_karman_strain_operator(qp,
                                                        vec2_n1, // epsilon_vk
                                                        vk_dwdxi_mat,
                                                        Bmat_vk);
            
            vec1_n1 = material_A_mat * vec2_n1; // stress
            stress(0,0) += vec1_n1(0); // sigma_xx
            stress(0,1) += vec1_n1(2); // sigma_xy
            stress(1,0) += vec1_n1(2); // sigma_yx
            stress(1,1) += vec1_n1(1); // sigma_yy
            
        }
    }
    
    // add the linear and nonlinear direct strains
    Bmat_mem.vector_mult(vec1_n1, _local_sol);
    vec2_n1 += vec1_n1;  // epsilon_mem + epsilon_vk
    
    // copy the total integrated stress to the vector
    vec1_n1(0) = stress(0,0);
    vec1_n1(1) = stress(1,1);
    vec1_n1(2) = stress(0,1);
    
    // now the internal force vector
    // this includes the membrane strain operator with all A and B material operators
    Bmat_mem.vector_mult_transpose(vec3_n2, vec1_n1);
    local_f += JxW[qp] * vec3_n2;
    
    if (if_bending) {
        if (if_vk) {
            // von Karman strain
            vec4_2 = vk_dwdxi_mat.transpose() * vec1_n1;
            Bmat_vk.vector_mult_transpose(vec3_n2, vec4_2);
            local_f += JxW[qp] * vec3_n2;
        }
        
        // now coupling with the bending strain
        // B_bend^T [B] B_mem
        vec1_n1 = material_B_mat * vec2_n1;
        Bmat_bend.vector_mult_transpose(vec3_n2, vec1_n1);
        local_f += JxW[qp] * vec3_n2;
        
        // now bending stress
        Bmat_bend.vector_mult(vec2_n1, _local_sol);
        vec1_n1 = material_D_mat * vec2_n1;
        Bmat_bend.vector_mult_transpose(vec3_n2, vec1_n1);
        local_f += JxW[qp] * vec3_n2;
    }
    
    if (request_jacobian) {
        // membrane - membrane
        Bmat_mem.left_multiply(mat1_n1n2, material_A_mat);
        Bmat_mem.right_multiply_transpose(mat2_n2n2, mat1_n1n2);
        local_jac += JxW[qp] * mat2_n2n2;
        
        if (if_bending) {
            if (if_vk) {
                // membrane - vk
                mat3 = RealMatrixX::Zero(vk_dwdxi_mat.rows(), n2);
                Bmat_vk.left_multiply(mat3, vk_dwdxi_mat);
                mat3 = material_A_mat * mat3;
                Bmat_mem.right_multiply_transpose(mat2_n2n2, mat3);
                local_jac += JxW[qp] * mat2_n2n2;
                
                // vk - membrane
                Bmat_mem.left_multiply(mat1_n1n2, material_A_mat);
                mat3 = vk_dwdxi_mat.transpose() * mat1_n1n2;
                Bmat_vk.right_multiply_transpose(mat2_n2n2, mat3);
                local_jac += JxW[qp] * mat2_n2n2;
                
                // if only the first order term of the Jacobian is needed, for
                // example for linearized buckling analysis, then the linear
                // stress combined with the variation of the von Karman strain
                // is included. Otherwise, all terms are included
                if (if_ignore_ho_jac) {
                    // vk - vk: first order term
                    mat3 = RealMatrixX::Zero(2, n2);
                    Bmat_vk.left_multiply(mat3, stress_l);
                    Bmat_vk.right_multiply_transpose(mat2_n2n2, mat3);
                    local_jac += JxW[qp] * mat2_n2n2;
                }
                else {
                    // vk - vk
                    mat3 = RealMatrixX::Zero(2, n2);
                    Bmat_vk.left_multiply(mat3, stress);
                    Bmat_vk.right_multiply_transpose(mat2_n2n2, mat3);
                    local_jac += JxW[qp] * mat2_n2n2;
                    
                    mat3 = RealMatrixX::Zero(vk_dwdxi_mat.rows(), n2);
                    Bmat_vk.left_multiply(mat3, vk_dwdxi_mat);
                    mat3 = vk_dwdxi_mat.transpose() * material_A_mat * mat3;
                    Bmat_vk.right_multiply_transpose(mat2_n2n2, mat3);
                    local_jac += JxW[qp] * mat2_n2n2;
                }
                
                // bending - vk
                mat3 = RealMatrixX::Zero(vk_dwdxi_mat.rows(), n2);
                Bmat_vk.left_multiply(mat3, vk_dwdxi_mat);
                mat3 = material_B_mat.transpose() * mat3;
                Bmat_bend.right_multiply_transpose(mat2_n2n2, mat3);
                local_jac += JxW[qp] * mat2_n2n2;
                
                // vk - bending
                Bmat_bend.left_multiply(mat1_n1n2, material_B_mat);
                mat3 = vk_dwdxi_mat.transpose() * mat1_n1n2;
                Bmat_vk.right_multiply_transpose(mat2_n2n2, mat3);
                local_jac += JxW[qp] * mat2_n2n2;
            }
            
            // bending - membrane
            Bmat_mem.left_multiply(mat1_n1n2, material_B_mat);
            Bmat_bend.right_multiply_transpose(mat2_n2n2, mat1_n1n2);
            local_jac += JxW[qp] * mat2_n2n2;
            
            // membrane - bending
            Bmat_bend.left_multiply(mat1_n1n2, material_B_mat);
            Bmat_mem.right_multiply_transpose(mat2_n2n2, mat1_n1n2);
            local_jac += JxW[qp] * mat2_n2n2;
            
            // bending - bending
            Bmat_bend.left_multiply(mat1_n1n2, material_D_mat);
            Bmat_bend.right_multiply_transpose(mat2_n2n2, mat1_n1n2);
            local_jac += JxW[qp] * mat2_n2n2;
        }
    }
}