NumericType PostProcessedQuantities<NumericType>::component( const libMesh::FEMContext& context, 
							       unsigned int component,
							       const libMesh::Point& p,
							       libMesh::Real /*time*/ )
  {
    // Check if the Elem is the same between the incoming context and the cached one.
    // If not, reinit the cached MultiphysicsSystem context
    if( &(context.get_elem()) != &(_multiphysics_context->get_elem()) )
      {
	_multiphysics_context->pre_fe_reinit(*_multiphysics_sys,&context.get_elem());
	_multiphysics_context->elem_fe_reinit();
      }

    /* Optimization since we expect this function to be called many times with
       the same point. _prev_point initialized to something absurd so this should 
       always be false the first time. */
    if( _prev_point != p )
      {
	_prev_point = p;
	std::vector<libMesh::Point> point_vec(1,p);
	this->_cache.clear();
	_multiphysics_sys->compute_element_cache( *(this->_multiphysics_context), point_vec, this->_cache );
      }

    return this->compute_quantities( component );
  }
  void PostProcessedQuantities<NumericType>::init_context( const libMesh::FEMContext& context )
  {
    // Make sure we prepare shape functions for our output variables.
    /*! \todo I believe this is redundant because it's done in the project_vector call. Double check. */
    for( typename std::map<VariableIndex,unsigned int>::const_iterator it = _quantity_var_map.begin();
	 it != _quantity_var_map.end(); it++ )
      {
	libMesh::FEBase* elem_fe = NULL;
	context.get_element_fe( it->first, elem_fe );
	elem_fe->get_phi();
	elem_fe->get_dphi();
	elem_fe->get_JxW();
	elem_fe->get_xyz();

	libMesh::FEBase* side_fe = NULL;
	context.get_side_fe( it->first, side_fe );
	side_fe->get_phi();
	side_fe->get_dphi();
	side_fe->get_JxW();
	side_fe->get_xyz();
      }

    // Create the context we'll be using to compute MultiphysicsSystem quantities
    _multiphysics_context.reset( new AssemblyContext( *_multiphysics_sys ) );
    _multiphysics_sys->init_context(*_multiphysics_context);
    return;
  }
  void BoundaryConditions::pin_value( libMesh::FEMContext& context,
				      const CachedValues& /*cache*/,
				      const bool request_jacobian,
				      const VariableIndex var, 
				      const double pin_value,
				      const libMesh::Point& pin_location, 
				      const double penalty )
  {
    if (context.elem->contains_point(pin_location))
      {
	libMesh::DenseSubVector<libMesh::Number> &F_var = *context.elem_subresiduals[var]; // residual
	libMesh::DenseSubMatrix<libMesh::Number> &K_var = *context.elem_subjacobians[var][var]; // jacobian

	// The number of local degrees of freedom in p variable.
	const unsigned int n_var_dofs = context.dof_indices_var[var].size();

	libMesh::Number var_value = context.point_value(var, pin_location);

	libMesh::FEType fe_type = context.element_fe_var[var]->get_fe_type();
      
	libMesh::Point point_loc_in_masterelem = 
	  libMesh::FEInterface::inverse_map(context.dim, fe_type, context.elem, pin_location);

	std::vector<libMesh::Real> phi(n_var_dofs);

	for (unsigned int i=0; i != n_var_dofs; i++)
	  phi[i] = libMesh::FEInterface::shape( context.dim, fe_type, context.elem, i, 
						point_loc_in_masterelem );
      
	for (unsigned int i=0; i != n_var_dofs; i++)
	  {
	    F_var(i) += penalty*(var_value - pin_value)*phi[i];
	  
	    /** \todo What the hell is the context.elem_solution_derivative all about? */
	    if (request_jacobian && context.elem_solution_derivative)
	      {
		libmesh_assert (context.elem_solution_derivative == 1.0);
	      
		for (unsigned int j=0; j != n_var_dofs; j++)
		  K_var(i,j) += penalty*phi[i]*phi[j];

	      } // End if request_jacobian
	  } // End i loop
      } // End if pin_location

    return;
  }
Ejemplo n.º 4
0
  void HeatTransfer::side_time_derivative( bool compute_jacobian,
					   libMesh::FEMContext& context,
					   CachedValues& cache )
  {
#ifdef GRINS_USE_GRVY_TIMERS
    this->_timer->BeginTimer("HeatTransfer::side_time_derivative");
#endif

    std::vector<BoundaryID> ids = context.side_boundary_ids();

    for( std::vector<BoundaryID>::const_iterator it = ids.begin();
	 it != ids.end(); it++ )
      {
	libmesh_assert (*it != libMesh::BoundaryInfo::invalid_id);

	_bc_handler->apply_neumann_bcs( context, cache, compute_jacobian, *it );
      }

#ifdef GRINS_USE_GRVY_TIMERS
    this->_timer->EndTimer("HeatTransfer::side_time_derivative");
#endif

    return;
  }
Ejemplo n.º 5
0
  NumericType PostProcessedQuantities<NumericType>::component( const libMesh::FEMContext& context,
                                                               unsigned int component,
                                                               const libMesh::Point& p,
                                                               libMesh::Real /*time*/ )
  {
    // Check if the Elem is the same between the incoming context and the cached one.
    // If not, reinit the cached MultiphysicsSystem context
    if(context.has_elem() && _multiphysics_context->has_elem())
      {
        if( &(context.get_elem()) != &(_multiphysics_context->get_elem()) )
          {
            _multiphysics_context->pre_fe_reinit(*_multiphysics_sys,&context.get_elem());
            _multiphysics_context->elem_fe_reinit();
          }
      }
    else if( !context.has_elem() && _multiphysics_context->has_elem() )
      {
        // Incoming context has NULL elem ==> SCALAR variables
        _multiphysics_context->pre_fe_reinit(*_multiphysics_sys,NULL);
        _multiphysics_context->elem_fe_reinit();
      }
    else if( context.has_elem() && !_multiphysics_context->has_elem() )
      {
        _multiphysics_context->pre_fe_reinit(*_multiphysics_sys,&context.get_elem());
        _multiphysics_context->elem_fe_reinit();
      }
    //else
    /* If has_elem() is false for both contexts, we're still dealing with SCALAR variables
       and therefore don't need to reinit. */

    libMesh::Real value = 0.0;

    // Quantity we want had better be there.
    libmesh_assert(_quantity_index_var_map.find(component) != _quantity_index_var_map.end());
    unsigned int quantity_index = _quantity_index_var_map.find(component)->second;

    _multiphysics_sys->compute_postprocessed_quantity( quantity_index,
                                                       *(this->_multiphysics_context),
                                                       p, value );

    return value;
  }
Ejemplo n.º 6
0
  void HeatTransfer::element_time_derivative( bool compute_jacobian,
					      libMesh::FEMContext& context,
					      CachedValues& /*cache*/ )
  {
#ifdef GRINS_USE_GRVY_TIMERS
    this->_timer->BeginTimer("HeatTransfer::element_time_derivative");
#endif

    // The number of local degrees of freedom in each variable.
    const unsigned int n_T_dofs = context.dof_indices_var[_T_var].size();
    const unsigned int n_u_dofs = context.dof_indices_var[_u_var].size();

    //TODO: check n_T_dofs is same as n_u_dofs, n_v_dofs, n_w_dofs

    // We get some references to cell-specific data that
    // will be used to assemble the linear system.

    // Element Jacobian * quadrature weights for interior integration.
    const std::vector<libMesh::Real> &JxW =
      context.element_fe_var[_T_var]->get_JxW();

    // The temperature shape functions at interior quadrature points.
    const std::vector<std::vector<libMesh::Real> >& T_phi =
      context.element_fe_var[_T_var]->get_phi();

    // The velocity shape functions at interior quadrature points.
    const std::vector<std::vector<libMesh::Real> >& vel_phi =
      context.element_fe_var[_u_var]->get_phi();

    // The temperature shape function gradients (in global coords.)
    // at interior quadrature points.
    const std::vector<std::vector<libMesh::RealGradient> >& T_gradphi =
      context.element_fe_var[_T_var]->get_dphi();

    const std::vector<libMesh::Point>& u_qpoint = 
      context.element_fe_var[this->_u_var]->get_xyz();

    // We do this in the incompressible Navier-Stokes class and need to do it here too
    // since _w_var won't have been defined in the global map.
    if (_dim != 3)
      _w_var = _u_var; // for convenience

    libMesh::DenseSubMatrix<libMesh::Number> &KTT = *context.elem_subjacobians[_T_var][_T_var]; // R_{T},{T}

    libMesh::DenseSubMatrix<libMesh::Number> &KTu = *context.elem_subjacobians[_T_var][_u_var]; // R_{T},{u}
    libMesh::DenseSubMatrix<libMesh::Number> &KTv = *context.elem_subjacobians[_T_var][_v_var]; // R_{T},{v}
    libMesh::DenseSubMatrix<libMesh::Number> &KTw = *context.elem_subjacobians[_T_var][_w_var]; // R_{T},{w}

    libMesh::DenseSubVector<libMesh::Number> &FT = *context.elem_subresiduals[_T_var]; // R_{T}

    // Now we will build the element Jacobian and residual.
    // Constructing the residual requires the solution and its
    // gradient from the previous timestep.  This must be
    // calculated at each quadrature point by summing the
    // solution degree-of-freedom values by the appropriate
    // weight functions.
    unsigned int n_qpoints = context.element_qrule->n_points();

    for (unsigned int qp=0; qp != n_qpoints; qp++)
      {
	// Compute the solution & its gradient at the old Newton iterate.
	libMesh::Number u, v, w;
	u = context.interior_value(_u_var, qp);
	v = context.interior_value(_v_var, qp);
	if (_dim == 3)
	  w = context.interior_value(_w_var, qp);

	libMesh::Gradient grad_T;
	grad_T = context.interior_gradient(_T_var, qp);

	libMesh::NumberVectorValue U (u,v);
	if (_dim == 3)
	  U(2) = w;

        const libMesh::Number r = u_qpoint[qp](0);

        libMesh::Real jac = JxW[qp];

        if( _is_axisymmetric )
          {
            jac *= r;
          }

	// First, an i-loop over the  degrees of freedom.
	for (unsigned int i=0; i != n_T_dofs; i++)
	  {
	    FT(i) += jac *
	      (-_rho*_Cp*T_phi[i][qp]*(U*grad_T)    // convection term
	       -_k*(T_gradphi[i][qp]*grad_T) );  // diffusion term

	    if (compute_jacobian)
	      {
		for (unsigned int j=0; j != n_T_dofs; j++)
		  {
		    // TODO: precompute some terms like:
		    //   _rho*_Cp*T_phi[i][qp]*(vel_phi[j][qp]*T_grad_phi[j][qp])

		    KTT(i,j) += jac *
		      (-_rho*_Cp*T_phi[i][qp]*(U*T_gradphi[j][qp])  // convection term
		       -_k*(T_gradphi[i][qp]*T_gradphi[j][qp])); // diffusion term
		  } // end of the inner dof (j) loop

		// Matrix contributions for the Tu, Tv and Tw couplings (n_T_dofs same as n_u_dofs, n_v_dofs and n_w_dofs)
		for (unsigned int j=0; j != n_u_dofs; j++)
		  {
		    KTu(i,j) += jac*(-_rho*_Cp*T_phi[i][qp]*(vel_phi[j][qp]*grad_T(0)));
		    KTv(i,j) += jac*(-_rho*_Cp*T_phi[i][qp]*(vel_phi[j][qp]*grad_T(1)));
		    if (_dim == 3)
		      KTw(i,j) += jac*(-_rho*_Cp*T_phi[i][qp]*(vel_phi[j][qp]*grad_T(2)));
		  } // end of the inner dof (j) loop

	      } // end - if (compute_jacobian && context.elem_solution_derivative)

	  } // end of the outer dof (i) loop
      } // end of the quadrature point (qp) loop

#ifdef GRINS_USE_GRVY_TIMERS
    this->_timer->EndTimer("HeatTransfer::element_time_derivative");
#endif

    return;
  }
Ejemplo n.º 7
0
  void HeatTransfer::mass_residual( bool compute_jacobian,
				    libMesh::FEMContext& context,
				    CachedValues& /*cache*/ )
  {
#ifdef GRINS_USE_GRVY_TIMERS
    this->_timer->BeginTimer("HeatTransfer::mass_residual");
#endif

    // First we get some references to cell-specific data that
    // will be used to assemble the linear system.

    // Element Jacobian * quadrature weights for interior integration
    const std::vector<libMesh::Real> &JxW = 
      context.element_fe_var[_T_var]->get_JxW();

    // The shape functions at interior quadrature points.
    const std::vector<std::vector<libMesh::Real> >& phi = 
      context.element_fe_var[_T_var]->get_phi();

    const std::vector<libMesh::Point>& u_qpoint = 
      context.element_fe_var[this->_u_var]->get_xyz();

    // The number of local degrees of freedom in each variable
    const unsigned int n_T_dofs = context.dof_indices_var[_T_var].size();

    // The subvectors and submatrices we need to fill:
    libMesh::DenseSubVector<libMesh::Real> &F = *context.elem_subresiduals[_T_var];

    libMesh::DenseSubMatrix<libMesh::Real> &M = *context.elem_subjacobians[_T_var][_T_var];

    unsigned int n_qpoints = context.element_qrule->n_points();

    for (unsigned int qp = 0; qp != n_qpoints; ++qp)
      {
	// For the mass residual, we need to be a little careful.
	// The time integrator is handling the time-discretization
	// for us so we need to supply M(u_fixed)*u for the residual.
	// u_fixed will be given by the fixed_interior_* functions
	// while u will be given by the interior_* functions.
	libMesh::Real T_dot = context.interior_value(_T_var, qp);

        const libMesh::Number r = u_qpoint[qp](0);

        libMesh::Real jac = JxW[qp];

        if( _is_axisymmetric )
          {
            jac *= r;
          }

	for (unsigned int i = 0; i != n_T_dofs; ++i)
	  {
	    F(i) += _rho*_Cp*T_dot*phi[i][qp]*jac;

	    if( compute_jacobian )
              {
                for (unsigned int j=0; j != n_T_dofs; j++)
                  {
		    // We're assuming rho, cp are constant w.r.t. T here.
                    M(i,j) += _rho*_Cp*phi[j][qp]*phi[i][qp]*jac;
                  }
              }// End of check on Jacobian
          
	  } // End of element dof loop
      
      } // End of the quadrature point loop

#ifdef GRINS_USE_GRVY_TIMERS
    this->_timer->EndTimer("HeatTransfer::mass_residual");
#endif

    return;
  }