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;
}
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;
}
bool
MAST::StructuralElement2D::internal_residual (bool request_jacobian,
                                              RealVectorX& f,
                                              RealMatrixX& jac,
                                              bool if_ignore_ho_jac)
{
    FEMOperatorMatrix Bmat_mem, Bmat_bend, Bmat_vk;
    
    const std::vector<Real>& JxW = _fe->get_JxW();
    const std::vector<libMesh::Point>& xyz = _fe->get_xyz();
    const unsigned int n_phi = (unsigned int)_fe->get_phi().size();
    const unsigned int n1= this->n_direct_strain_components(), n2=6*n_phi,
    n3 = this->n_von_karman_strain_components();
    RealMatrixX
    material_A_mat,
    material_B_mat,
    material_D_mat,
    mat1_n1n2     = RealMatrixX::Zero(n1,n2),
    mat2_n2n2     = RealMatrixX::Zero(n2,n2),
    mat3,
    mat4_n3n2     = RealMatrixX::Zero(n3,n2),
    vk_dwdxi_mat  = RealMatrixX::Zero(n1,n3),
    stress        = RealMatrixX::Zero(2,2),
    stress_l      = RealMatrixX::Zero(2,2),
    local_jac     = RealMatrixX::Zero(n2,n2);
    RealVectorX
    vec1_n1    = RealVectorX::Zero(n1),
    vec2_n1    = RealVectorX::Zero(n1),
    vec3_n2    = RealVectorX::Zero(n2),
    vec4_n3    = RealVectorX::Zero(n3),
    vec5_n3    = RealVectorX::Zero(n3),
    local_f    = RealVectorX::Zero(n2);
    
    local_f.setZero();
    local_jac.setZero();
    
    Bmat_mem.reinit(n1, _system.n_vars(), n_phi); // three stress-strain components
    Bmat_bend.reinit(n1, _system.n_vars(), n_phi);
    Bmat_vk.reinit(n3, _system.n_vars(), n_phi); // only dw/dx and dw/dy
    
    bool if_vk = (_property.strain_type() == MAST::VON_KARMAN_STRAIN),
    if_bending = (_property.bending_model(_elem, _fe->get_fe_type()) != MAST::NO_BENDING);
    
    std::auto_ptr<MAST::FieldFunction<RealMatrixX > >
    mat_stiff_A  = _property.stiffness_A_matrix(*this),
    mat_stiff_B  = _property.stiffness_B_matrix(*this),
    mat_stiff_D  = _property.stiffness_D_matrix(*this);
    
    
    libMesh::Point p;
    
    for (unsigned int qp=0; qp<JxW.size(); qp++) {
        
        this->local_elem().global_coordinates_location(xyz[qp], p);
        
        // get the material matrix
        (*mat_stiff_A)(p, _time, material_A_mat);
        
        if (if_bending) {
            (*mat_stiff_B)(p, _time, material_B_mat);
            (*mat_stiff_D)(p, _time, material_D_mat);
        }
        
        // now calculte the quantity for these matrices
        _internal_residual_operation(if_bending, if_vk, n2, qp, JxW,
                                     request_jacobian, if_ignore_ho_jac,
                                     local_f, local_jac,
                                     Bmat_mem, Bmat_bend, Bmat_vk,
                                     stress, stress_l, vk_dwdxi_mat, material_A_mat,
                                     material_B_mat, material_D_mat, vec1_n1,
                                     vec2_n1, vec3_n2, vec4_n3,
                                     vec5_n3, mat1_n1n2, mat2_n2n2,
                                     mat3, mat4_n3n2);
        
    }
    
    
    // now calculate the transverse shear contribution if appropriate for the
    // element
    if (if_bending &&
        _bending_operator->include_transverse_shear_energy())
        _bending_operator->calculate_transverse_shear_residual(request_jacobian,
                                                               
                                                               local_f,
                                                               local_jac,
                                                               NULL);
    
    
    // now transform to the global coorodinate system
    transform_vector_to_global_system(local_f, vec3_n2);
    f += vec3_n2;
    if (request_jacobian) {
        // for 2D elements
        if (_elem.dim() == 2) {
            // add small values to the diagonal of the theta_z dofs
            for (unsigned int i=0; i<n_phi; i++)
                local_jac(5*n_phi+i, 5*n_phi+i) = -1.0e-8;
        }
        transform_matrix_to_global_system(local_jac, mat2_n2n2);
        jac += mat2_n2n2;
    }
    
    return request_jacobian;
}