void
MAST::MindlinBendingOperator::
initialize_bending_strain_operator_for_z(const MAST::FEBase& fe,
                                         const unsigned int qp,
                                         const Real z,
                                         MAST::FEMOperatorMatrix& Bmat_bend) {
    
    const std::vector<std::vector<libMesh::RealVectorValue> >& dphi = fe.get_dphi();
    const std::vector<std::vector<Real> >& phi = fe.get_phi();
    
    const unsigned int n_phi = (unsigned int)phi.size();
    
    RealVectorX phi_vec = RealVectorX::Zero(n_phi);
    for ( unsigned int i_nd=0; i_nd<n_phi; i_nd++ )
        phi_vec(i_nd) = dphi[i_nd][qp](0);  // dphi/dx
    
    phi_vec   *= z;
    Bmat_bend.set_shape_function(0, 4, phi_vec); // epsilon-x: thetay
    phi_vec   *= -1.0;
    Bmat_bend.set_shape_function(2, 3, phi_vec); // gamma-xy : thetax
    
    
    for ( unsigned int i_nd=0; i_nd<n_phi; i_nd++ )
        phi_vec(i_nd) = dphi[i_nd][qp](1);  // dphi/dy
    
    phi_vec   *= z;
    Bmat_bend.set_shape_function(2, 4, phi_vec); // gamma-xy : thetay
    //Bmat_trans.set_shape_function(1, 2, phi_vec); // gamma-yz : w
    phi_vec   *= -1.0;
    Bmat_bend.set_shape_function(1, 3, phi_vec); // epsilon-y: thetax
}
void
MAST::StructuralElement2D::
initialize_von_karman_strain_operator_sensitivity(const unsigned int qp,
                                                  RealMatrixX &vk_dwdxi_mat_sens) {
    
    const std::vector<std::vector<libMesh::RealVectorValue> >& dphi = _fe->get_dphi();
    const unsigned int n_phi = (unsigned int)dphi.size();
    
    libmesh_assert_equal_to(vk_dwdxi_mat_sens.rows(), 3);
    libmesh_assert_equal_to(vk_dwdxi_mat_sens.cols(), 2);
    
    Real dw=0.;
    vk_dwdxi_mat_sens.setConstant(0.);
    
    RealVectorX phi_vec  = RealVectorX::Zero(n_phi);
    
    dw = 0.;
    for ( unsigned int i_nd=0; i_nd<n_phi; i_nd++ ) {
        phi_vec(i_nd) = dphi[i_nd][qp](0);  // dphi/dx
        dw += phi_vec(i_nd)*_local_sol_sens(2*n_phi+i_nd); // dw/dx
    }
    vk_dwdxi_mat_sens(0, 0) = dw;  // epsilon-xx : dw/dx
    vk_dwdxi_mat_sens(2, 1) = dw;  // gamma-xy : dw/dx
    
    dw = 0.;
    for ( unsigned int i_nd=0; i_nd<n_phi; i_nd++ ) {
        phi_vec(i_nd) = dphi[i_nd][qp](1);  // dphi/dy
        dw += phi_vec(i_nd)*_local_sol_sens(2*n_phi+i_nd); // dw/dy
    }
    vk_dwdxi_mat_sens(1, 1) = dw;  // epsilon-yy : dw/dy
    vk_dwdxi_mat_sens(2, 0) = dw;  // gamma-xy : dw/dy
}
bool
MAST::HeatConductionElementBase::
surface_convection_residual(bool request_jacobian,
                            RealVectorX& f,
                            RealMatrixX& jac,
                            const unsigned int s,
                            MAST::BoundaryConditionBase& p) {
    
    // prepare the side finite element
    std::auto_ptr<libMesh::FEBase> fe;
    std::auto_ptr<libMesh::QBase> qrule;
    _get_side_fe_and_qrule(get_elem_for_quadrature(), s, fe, qrule);

    // get the function from this boundary condition
    const MAST::FieldFunction<Real>
    &coeff = p.get<MAST::FieldFunction<Real> >("convection_coeff"),
    &T_amb = p.get<MAST::FieldFunction<Real> >("ambient_temperature");
    
    const std::vector<Real> &JxW               = fe->get_JxW();
    const std::vector<libMesh::Point>& qpoint  = fe->get_xyz();
    const std::vector<std::vector<Real> >& phi = fe->get_phi();
    const unsigned int n_phi                   = (unsigned int)phi.size();
    
    
    RealVectorX  phi_vec  = RealVectorX::Zero(n_phi);
    RealMatrixX  mat      = RealMatrixX::Zero(n_phi, n_phi);
    Real temp, amb_temp, h_coeff;
    libMesh::Point pt;
    MAST::FEMOperatorMatrix Bmat;
    
    
    for (unsigned int qp=0; qp<qpoint.size(); qp++) {
        
        _local_elem->global_coordinates_location (qpoint[qp], pt);
        
        // 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];
        
        // value of flux
        coeff(pt, _time, h_coeff);
        T_amb(pt, _time, amb_temp);
        temp  = phi_vec.dot(_sol);
        
        // normal flux is given as:
        // qi_ni = h_coeff * (T - T_amb)
        //
        f   += JxW[qp] * phi_vec * h_coeff * (temp - amb_temp);
        
        if (request_jacobian) {
            
            Bmat.reinit(1, phi_vec);
            Bmat.right_multiply_transpose(mat, Bmat);
            jac += JxW[qp] * mat * h_coeff;
        }
    }
    
    return request_jacobian;
}
bool
MAST::HeatConductionElementBase::
surface_radiation_residual(bool request_jacobian,
                           RealVectorX& f,
                           RealMatrixX& jac,
                           MAST::BoundaryConditionBase& p) {
    
    // get the function from this boundary condition
    const MAST::FieldFunction<Real>
    &emissivity = p.get<MAST::FieldFunction<Real> >("emissivity");
    
    const MAST::Parameter
    &T_amb      = p.get<MAST::Parameter>("ambient_temperature"),
    &T_ref_zero = p.get<MAST::Parameter>("reference_zero_temperature"),
    &sb_const   = p.get<MAST::Parameter>("stefan_bolzmann_constant");
    
    
    const std::vector<Real> &JxW               = _fe->get_JxW();
    const std::vector<libMesh::Point>& qpoint  = _fe->get_xyz();
    const std::vector<std::vector<Real> >& phi = _fe->get_phi();
    const unsigned int n_phi                   = (unsigned int)phi.size();
    
    RealVectorX phi_vec  = RealVectorX::Zero(n_phi);
    RealMatrixX mat      = RealMatrixX::Zero(n_phi, n_phi);
    const Real
    sbc      = sb_const(),
    amb_temp = T_amb(),
    zero_ref = T_ref_zero();
    Real temp, emiss;
    libMesh::Point pt;
    MAST::FEMOperatorMatrix Bmat;
    
    for (unsigned int qp=0; qp<qpoint.size(); qp++) {
        
        _local_elem->global_coordinates_location (qpoint[qp], pt);
        
        // 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];
        
        // value of flux
        emissivity(pt, _time, emiss);
        temp  = phi_vec.dot(_sol);
        
        f   += JxW[qp] * phi_vec * sbc * emiss *
        (pow(temp-zero_ref, 4.) - pow(amb_temp-zero_ref, 4.));
        
        if (request_jacobian) {
            
            Bmat.reinit(1, phi_vec);
            Bmat.right_multiply_transpose(mat, Bmat);
            jac +=  JxW[qp] * mat * sbc * emiss * 4. * pow(temp-zero_ref, 3.);
        }
    }
    
    
    return request_jacobian;
}
void
MAST::StructuralElement2D::
initialize_von_karman_strain_operator(const unsigned int qp,
                                      RealVectorX& vk_strain,
                                      RealMatrixX& vk_dwdxi_mat,
                                      MAST::FEMOperatorMatrix& Bmat_vk) {
    
    const std::vector<std::vector<libMesh::RealVectorValue> >& dphi = _fe->get_dphi();
    const unsigned int n_phi = (unsigned int)dphi.size();
    
    libmesh_assert_equal_to(vk_strain.size(), 3);
    libmesh_assert_equal_to(vk_dwdxi_mat.rows(), 3);
    libmesh_assert_equal_to(vk_dwdxi_mat.cols(), 2);
    libmesh_assert_equal_to(Bmat_vk.m(), 2);
    libmesh_assert_equal_to(Bmat_vk.n(), 6*n_phi);
    
    Real dw=0.;
    vk_strain.setConstant(0.);
    vk_dwdxi_mat.setConstant(0.);
    
    RealVectorX phi_vec  = RealVectorX::Zero(n_phi);
    
    dw = 0.;
    for ( unsigned int i_nd=0; i_nd<n_phi; i_nd++ ) {
        phi_vec(i_nd) = dphi[i_nd][qp](0);  // dphi/dx
        dw += phi_vec(i_nd)*_local_sol(2*n_phi+i_nd); // dw/dx
    }
    Bmat_vk.set_shape_function(0, 2, phi_vec); // dw/dx
    vk_dwdxi_mat(0, 0) = dw;  // epsilon-xx : dw/dx
    vk_dwdxi_mat(2, 1) = dw;  // gamma-xy : dw/dx
    vk_strain(0) = 0.5*dw*dw; // 1/2 * (dw/dx)^2
    vk_strain(2) = dw;        // (dw/dx)*(dw/dy)  only dw/dx is provided here
    
    dw = 0.;
    for ( unsigned int i_nd=0; i_nd<n_phi; i_nd++ ) {
        phi_vec(i_nd) = dphi[i_nd][qp](1);  // dphi/dy
        dw += phi_vec(i_nd)*_local_sol(2*n_phi+i_nd); // dw/dy
    }
    Bmat_vk.set_shape_function(1, 2, phi_vec); // dw/dy
    vk_dwdxi_mat(1, 1) = dw;  // epsilon-yy : dw/dy
    vk_dwdxi_mat(2, 0) = dw;  // gamma-xy : dw/dy
    vk_strain(1) = 0.5*dw*dw; // 1/2 * (dw/dy)^2
    vk_strain(2) *= dw;       // (dw/dx)*(dw/dy)
}
bool
MAST::HeatConductionElementBase::
surface_flux_residual(bool request_jacobian,
                      RealVectorX& f,
                      RealMatrixX& jac,
                      const unsigned int s,
                      MAST::BoundaryConditionBase& p) {
    
    // prepare the side finite element
    std::auto_ptr<libMesh::FEBase> fe;
    std::auto_ptr<libMesh::QBase> qrule;
    _get_side_fe_and_qrule(get_elem_for_quadrature(), s, fe, qrule);
    
    
    // get the function from this boundary condition
    const MAST::FieldFunction<Real>& func =
    p.get<MAST::FieldFunction<Real> >("heat_flux");
    
    
    const std::vector<Real> &JxW                 = fe->get_JxW();
    const std::vector<libMesh::Point>& qpoint    = fe->get_xyz();
    const std::vector<std::vector<Real> >& phi   = fe->get_phi();
    const unsigned int n_phi                     = (unsigned int)phi.size();
    
    RealVectorX phi_vec  = RealVectorX::Zero(n_phi);
    libMesh::Point pt;
    Real  flux;
    
    for (unsigned int qp=0; qp<qpoint.size(); qp++) {
        
        _local_elem->global_coordinates_location (qpoint[qp], pt);
        
        // 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];
        
        // get the value of flux = q_i . n_i
        func(pt, _time, flux);
        
        f   +=  JxW[qp] * phi_vec * flux;
    }
    
    // calculation of the load vector is independent of solution
    return false;
}
bool
MAST::StructuralElementBase::
small_disturbance_surface_pressure_force(bool request_jacobian,
                                         DenseRealVector &f,
                                         DenseRealMatrix &jac,
                                         MAST::BoundaryCondition &p) {
    libmesh_assert(_elem.dim() < 3); // only applicable for lower dimensional elements.
    libmesh_assert(!follower_forces); // not implemented yet for follower forces
    libmesh_assert_equal_to(p.type(), MAST::SMALL_DISTURBANCE_MOTION);

    MAST::SmallDisturbanceMotion& sd_motion = dynamic_cast<MAST::SmallDisturbanceMotion&>(p);
    MAST::SurfaceMotionBase& surf_motion = sd_motion.get_deformation();
    MAST::SmallDisturbanceSurfacePressure& surf_press = sd_motion.get_pressure();

    FEMOperatorMatrix Bmat;
    
    // get the function from this boundary condition
    const std::vector<Real> &JxW = _fe->get_JxW();
    
    // Physical location of the quadrature points
    const std::vector<libMesh::Point>& qpoint = _fe->get_xyz();
    const std::vector<std::vector<Real> >& phi = _fe->get_phi();
    const unsigned int n_phi = (unsigned int)phi.size();
    const unsigned int n1=3, n2=6*n_phi;
    
    // normal for face integration
    libMesh::Point normal;
    // direction of pressure assumed to be normal (along local z-axis)
    // to the element face for 2D and along local y-axis for 1D element.
    normal(_elem.dim()) = -1.;
    
    Real press;
    ValType dpress;
    DenseRealVector phi_vec;
    libMesh::DenseVector<ValType> wtrans, utrans, dn_rot, force, local_f, vec_n2;
    wtrans.resize(3); utrans.resize(3); dn_rot.resize(3);
    phi_vec.resize(n_phi); force.resize(2*n1); local_f.resize(n2);
    vec_n2.resize(n2);
    
    libMesh::Point pt;
    
    for (unsigned int qp=0; qp<qpoint.size(); qp++)
    {
        this->global_coordinates(qpoint[qp], pt);
        
        // 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(2*n1, phi_vec);
        
        // get pressure and deformation information
        surf_press.surface_pressure<ValType>(_system.time, pt, press, dpress);
        surf_motion.surface_velocity(_system.time, pt, normal,
                                     wtrans, utrans, dn_rot);
//        libMesh::out << std::setw(15) << pt(0)
//        << std::setw(15) << std::real(press)
//        << std::setw(15) << std::imag(press)
//        << std::setw(15) << std::real(dpress)
//        << std::setw(15) << std::imag(dpress)
//        << std::setw(15) << std::real(utrans(1))
//        << std::setw(15) << std::imag(utrans(1))
//        << std::setw(15) << std::real(dn_rot(0))
//        << std::setw(15) << std::imag(dn_rot(0)) << std::endl;
        
        
        // calculate force
        for (unsigned int i_dim=0; i_dim<n1; i_dim++)
            force(i_dim) = ( press * dn_rot(i_dim) + // steady pressure
                            dpress * normal(i_dim)); // unsteady pressure
        
        Bmat.vector_mult_transpose(vec_n2, force);
        
        local_f.add(JxW[qp], vec_n2);
    }
    
    // now transform to the global system and add
    transform_to_global_system(local_f, vec_n2);
    MAST::add_to_assembled_vector(f, vec_n2);
    

    return (request_jacobian && follower_forces);
}
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::StructuralElementBase::small_disturbance_surface_pressure_force(bool request_jacobian,
                                                                      DenseRealVector &f,
                                                                      DenseRealMatrix &jac,
                                                                      const unsigned int side,
                                                                      MAST::BoundaryCondition &p) {
    libmesh_assert(!follower_forces); // not implemented yet for follower forces
    libmesh_assert_equal_to(p.type(), MAST::SMALL_DISTURBANCE_MOTION);
    
    MAST::SmallDisturbanceMotion& sd_motion = dynamic_cast<MAST::SmallDisturbanceMotion&>(p);
    MAST::SurfaceMotionBase& surf_motion = sd_motion.get_deformation();
    MAST::SmallDisturbanceSurfacePressure& surf_press = sd_motion.get_pressure();
    
    FEMOperatorMatrix Bmat;
    
    // get the function from this boundary condition
    std::auto_ptr<libMesh::FEBase> fe;
    std::auto_ptr<libMesh::QBase> qrule;
    _get_side_fe_and_qrule(this->local_elem().local_elem(), side, fe, qrule);
    
    const std::vector<Real> &JxW = fe->get_JxW();
    
    // Physical location of the quadrature points
    const std::vector<libMesh::Point>& qpoint = fe->get_xyz();
    const std::vector<std::vector<Real> >& phi = fe->get_phi();
    const unsigned int n_phi = (unsigned int)phi.size();
    const unsigned int n1=3, n2=6*n_phi;
    
    // boundary normals
    const std::vector<libMesh::Point>& face_normals = fe->get_normals();
    
    Real press;
    ValType dpress;
    DenseRealVector phi_vec;
    libMesh::DenseVector<ValType> wtrans, utrans, dn_rot, force, local_f, vec_n2;
    wtrans.resize(3); utrans.resize(3); dn_rot.resize(3);
    phi_vec.resize(n_phi); force.resize(2*n1); local_f.resize(n2);
    vec_n2.resize(n2);
    
    libMesh::Point pt;
    
    for (unsigned int qp=0; qp<qpoint.size(); qp++)
    {
        this->global_coordinates(qpoint[qp], pt);
        
        // 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(2*n1, phi_vec);

        // get pressure and deformation information
        surf_press.surface_pressure<ValType>(_system.time, pt, press, dpress);
        surf_motion.surface_velocity(_system.time, pt, face_normals[qp],
                                     wtrans, utrans, dn_rot);
        
        //            press = 0.;
        //            dpress = Complex(2./4.*std::real(dn_rot(0)),  2./4./.1*std::imag(utrans(1)));
        //            libMesh::out << q_point[qp](0)
        //            << "  " << std::real(utrans(1))
        //            << "  " << std::imag(utrans(1))
        //            << "  " << std::real(dn_rot(0))
        //            << "  " << std::imag(dn_rot(0))
        //            << "  " << std::real(press)
        //            << "  " << std::imag(press)
        //            << "  " << std::real(dpress)
        //            << "  " << std::imag(dpress) << std::endl;

        // calculate force
        for (unsigned int i_dim=0; i_dim<n1; i_dim++)
            force(i_dim) =  ( press * dn_rot(i_dim) + // steady pressure
                             dpress * face_normals[qp](i_dim)); // unsteady pressure

        
        Bmat.vector_mult_transpose(vec_n2, force);
        
        local_f.add(JxW[qp], vec_n2);
    }
    
    // now transform to the global system and add
    if (_elem.dim() < 3) {
        transform_to_global_system(local_f, vec_n2);
        MAST::add_to_assembled_vector(f, vec_n2);
    }
    else
        MAST::add_to_assembled_vector(f, local_f);
    

    return (request_jacobian && follower_forces);
}
bool
MAST::StructuralElementBase::surface_pressure_force(bool request_jacobian,
                                                    DenseRealVector &f,
                                                    DenseRealMatrix &jac,
                                                    MAST::BoundaryCondition &p) {

    libmesh_assert(_elem.dim() < 3); // only applicable for lower dimensional elements
    libmesh_assert(!follower_forces); // not implemented yet for follower forces
    
    FEMOperatorMatrix Bmat;
    
    // get the function from this boundary condition
    MAST::FieldFunction<Real>& func =
    dynamic_cast<MAST::FieldFunction<Real>&>(p.function());
    const std::vector<Real> &JxW = _fe->get_JxW();
    
    // Physical location of the quadrature points
    const std::vector<libMesh::Point>& qpoint = _fe->get_xyz();
    const std::vector<std::vector<Real> >& phi = _fe->get_phi();
    const unsigned int n_phi = (unsigned int)phi.size();
    const unsigned int n1=3, n2=6*n_phi;
    
    // normal for face integration
    libMesh::Point normal;
    // direction of pressure assumed to be normal (along local z-axis)
    // to the element face for 2D and along local y-axis for 1D element.
    normal(_elem.dim()) = -1.;
    
    Real press;
    
    DenseRealVector phi_vec, force, local_f, vec_n2;
    phi_vec.resize(n_phi); force.resize(2*n1); local_f.resize(n2);
    vec_n2.resize(n2);
    
    libMesh::Point pt;
    
    for (unsigned int qp=0; qp<qpoint.size(); qp++)
    {
        
        this->global_coordinates(qpoint[qp], pt);

        // 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(2*n1, phi_vec);
        
        // get pressure value
        func(pt, _system.time, press);
        
        // calculate force
        for (unsigned int i_dim=0; i_dim<n1; i_dim++)
            force(i_dim) = press * normal(i_dim);
        
        Bmat.vector_mult_transpose(vec_n2, force);
        
        local_f.add(JxW[qp], vec_n2);
    }
    
    // now transform to the global system and add
    transform_to_global_system(local_f, vec_n2);
    f.add(1., vec_n2);
    

    return (request_jacobian && follower_forces);
}
bool
MAST::StructuralElementBase::surface_pressure_force(bool request_jacobian,
                                                    DenseRealVector &f,
                                                    DenseRealMatrix &jac,
                                                    const unsigned int side,
                                                    MAST::BoundaryCondition &p) {

    libmesh_assert(!follower_forces); // not implemented yet for follower forces
    
    FEMOperatorMatrix Bmat;
    
    // get the function from this boundary condition
    MAST::FieldFunction<Real>& func =
    dynamic_cast<MAST::FieldFunction<Real>&>(p.function());
    std::auto_ptr<libMesh::FEBase> fe;
    std::auto_ptr<libMesh::QBase> qrule;
    _get_side_fe_and_qrule(this->get_elem_for_quadrature(), side, fe, qrule);
    
    const std::vector<Real> &JxW = fe->get_JxW();
    
    // Physical location of the quadrature points
    const std::vector<libMesh::Point>& qpoint = fe->get_xyz();
    const std::vector<std::vector<Real> >& phi = fe->get_phi();
    const unsigned int n_phi = (unsigned int)phi.size();
    const unsigned int n1=3, n2=6*n_phi;
    
    // boundary normals
    const std::vector<libMesh::Point>& face_normals = fe->get_normals();
    Real press;
    
    DenseRealVector phi_vec, force, local_f, vec_n2;
    phi_vec.resize(n_phi); force.resize(2*n1); local_f.resize(n2);
    vec_n2.resize(n2);
    
    libMesh::Point pt;
    
    for (unsigned int qp=0; qp<qpoint.size(); qp++)
    {
        this->global_coordinates(qpoint[qp], pt);
        
        // 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(2*n1, phi_vec);
        
        // get pressure value
        func(pt, _system.time, press);
        
        // calculate force
        for (unsigned int i_dim=0; i_dim<n1; i_dim++)
            force(i_dim) = press * face_normals[qp](i_dim);
        
        Bmat.vector_mult_transpose(vec_n2, force);
        
        local_f.add(JxW[qp], vec_n2);
    }
    
    // now transform to the global system and add
    if (_elem.dim() < 3) {
        transform_to_global_system(local_f, vec_n2);
        f.add(1., vec_n2);
    }
    else
        f.add(1., local_f);
    

    return (request_jacobian && follower_forces);
}
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;
}