예제 #1
0
  void LowMachNavierStokes<Mu,SH,TC>::assemble_mass_time_deriv( bool /*compute_jacobian*/, 
								AssemblyContext& context,
								CachedValues& cache )
  {
    // The number of local degrees of freedom in each variable.
    const unsigned int n_p_dofs = context.get_dof_indices(this->_p_var).size();

    // Element Jacobian * quadrature weights for interior integration.
    const std::vector<libMesh::Real> &JxW =
      context.get_element_fe(this->_u_var)->get_JxW();

    // The pressure shape functions at interior quadrature points.
    const std::vector<std::vector<libMesh::Real> >& p_phi =
      context.get_element_fe(this->_p_var)->get_phi();

    libMesh::DenseSubVector<libMesh::Number> &Fp = context.get_elem_residual(this->_p_var); // R_{p}

    unsigned int n_qpoints = context.get_element_qrule().n_points();

    for (unsigned int qp=0; qp != n_qpoints; qp++)
      {
	libMesh::Number u, v, T;
	u = cache.get_cached_values(Cache::X_VELOCITY)[qp];
	v = cache.get_cached_values(Cache::Y_VELOCITY)[qp];

	T = cache.get_cached_values(Cache::TEMPERATURE)[qp];

	libMesh::Gradient grad_u = cache.get_cached_gradient_values(Cache::X_VELOCITY_GRAD)[qp];
	libMesh::Gradient grad_v = cache.get_cached_gradient_values(Cache::Y_VELOCITY_GRAD)[qp];

	libMesh::Gradient grad_T = cache.get_cached_gradient_values(Cache::TEMPERATURE_GRAD)[qp];

	libMesh::NumberVectorValue U(u,v);
	if (this->_dim == 3)
	  U(2) = cache.get_cached_values(Cache::Z_VELOCITY)[qp]; // w

	libMesh::Number divU = grad_u(0) + grad_v(1);
	if (this->_dim == 3)
          {
	    libMesh::Gradient grad_w = cache.get_cached_gradient_values(Cache::Z_VELOCITY_GRAD)[qp];
	    divU += grad_w(2);
          }

	// Now a loop over the pressure degrees of freedom.  This
	// computes the contributions of the continuity equation.
	for (unsigned int i=0; i != n_p_dofs; i++)
	  {
	    Fp(i) += (-U*grad_T/T + divU)*p_phi[i][qp]*JxW[qp];
	  }
      }

    return;
  }
예제 #2
0
  void LowMachNavierStokes<Mu,SH,TC>::assemble_energy_time_deriv( bool /*compute_jacobian*/,
								  AssemblyContext& context,
								  CachedValues& cache )
  {
    // The number of local degrees of freedom in each variable.
    const unsigned int n_T_dofs = context.get_dof_indices(this->_T_var).size();

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

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

    // The temperature shape functions gradients at interior quadrature points.
    const std::vector<std::vector<libMesh::RealGradient> >& T_gradphi =
      context.get_element_fe(this->_T_var)->get_dphi();

    libMesh::DenseSubVector<libMesh::Number> &FT = context.get_elem_residual(this->_T_var); // R_{T}

    unsigned int n_qpoints = context.get_element_qrule().n_points();
    for (unsigned int qp=0; qp != n_qpoints; qp++)
      {
	libMesh::Number u, v, T, p0;
	u = cache.get_cached_values(Cache::X_VELOCITY)[qp];
	v = cache.get_cached_values(Cache::Y_VELOCITY)[qp];
	T = cache.get_cached_values(Cache::TEMPERATURE)[qp];
	p0 = cache.get_cached_values(Cache::THERMO_PRESSURE)[qp];

	libMesh::Gradient grad_T = cache.get_cached_gradient_values(Cache::TEMPERATURE_GRAD)[qp];

	libMesh::NumberVectorValue U(u,v);
	if (this->_dim == 3)
	  U(2) = cache.get_cached_values(Cache::Z_VELOCITY)[qp]; // w

	libMesh::Number k = this->_k(T);
	libMesh::Number cp = this->_cp(T);

	libMesh::Number rho = this->rho( T, p0 );

	// Now a loop over the pressure degrees of freedom.  This
	// computes the contributions of the continuity equation.
	for (unsigned int i=0; i != n_T_dofs; i++)
	  {
	    FT(i) += ( -rho*cp*U*grad_T*T_phi[i][qp] // convection term
		       - k*grad_T*T_gradphi[i][qp]            // diffusion term
		       )*JxW[qp]; 
	  }
      }

    return;
  }
예제 #3
0
  void LowMachNavierStokes<Mu,SH,TC>::assemble_momentum_time_deriv( bool /*compute_jacobian*/, 
								    AssemblyContext& context,
								    CachedValues& cache )
  {
    // The number of local degrees of freedom in each variable.
    const unsigned int n_u_dofs = context.get_dof_indices(this->_u_var).size();

    // Check number of dofs is same for _u_var, v_var and w_var.
    libmesh_assert (n_u_dofs == context.get_dof_indices(this->_v_var).size());
    if (this->_dim == 3)
      libmesh_assert (n_u_dofs == context.get_dof_indices(this->_w_var).size());

    // Element Jacobian * quadrature weights for interior integration.
    const std::vector<libMesh::Real> &JxW =
      context.get_element_fe(this->_u_var)->get_JxW();

    // The pressure shape functions at interior quadrature points.
    const std::vector<std::vector<libMesh::Real> >& u_phi =
      context.get_element_fe(this->_u_var)->get_phi();

    // The velocity shape function gradients at interior quadrature points.
    const std::vector<std::vector<libMesh::RealGradient> >& u_gradphi =
      context.get_element_fe(this->_u_var)->get_dphi();

    libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(this->_u_var); // R_{u}
    libMesh::DenseSubVector<libMesh::Number> &Fv = context.get_elem_residual(this->_v_var); // R_{v}
    libMesh::DenseSubVector<libMesh::Number> &Fw = context.get_elem_residual(this->_w_var); // R_{w}

    unsigned int n_qpoints = context.get_element_qrule().n_points();
    for (unsigned int qp=0; qp != n_qpoints; qp++)
      {
	libMesh::Number u, v, p, p0, T;
	u = cache.get_cached_values(Cache::X_VELOCITY)[qp];
	v = cache.get_cached_values(Cache::Y_VELOCITY)[qp];

	T = cache.get_cached_values(Cache::TEMPERATURE)[qp];
	p = cache.get_cached_values(Cache::PRESSURE)[qp];
	p0 = cache.get_cached_values(Cache::THERMO_PRESSURE)[qp];

	libMesh::Gradient grad_u = cache.get_cached_gradient_values(Cache::X_VELOCITY_GRAD)[qp];
	libMesh::Gradient grad_v = cache.get_cached_gradient_values(Cache::Y_VELOCITY_GRAD)[qp];

	libMesh::Gradient grad_w;
	if (this->_dim == 3)
	  grad_w = cache.get_cached_gradient_values(Cache::Z_VELOCITY_GRAD)[qp];

	libMesh::NumberVectorValue grad_uT( grad_u(0), grad_v(0) ); 
	libMesh::NumberVectorValue grad_vT( grad_u(1), grad_v(1) );
	libMesh::NumberVectorValue grad_wT;
	if( this->_dim == 3 )
	  {
	    grad_uT(2) = grad_w(0);
	    grad_vT(2) = grad_w(1);
	    grad_wT = libMesh::NumberVectorValue( grad_u(2), grad_v(2), grad_w(2) );
	  }

	libMesh::NumberVectorValue U(u,v);
	if (this->_dim == 3)
	  U(2) = cache.get_cached_values(Cache::Z_VELOCITY)[qp]; // w

	libMesh::Number divU = grad_u(0) + grad_v(1);
	if (this->_dim == 3)
	  divU += grad_w(2);

	libMesh::Number rho = this->rho( T, p0 );
      
	// Now a loop over the pressure degrees of freedom.  This
	// computes the contributions of the continuity equation.
	for (unsigned int i=0; i != n_u_dofs; i++)
	  {
	    Fu(i) += ( -rho*U*grad_u*u_phi[i][qp]                 // convection term
		       + p*u_gradphi[i][qp](0)                           // pressure term
		       - this->_mu(T)*(u_gradphi[i][qp]*grad_u + u_gradphi[i][qp]*grad_uT
				       - 2.0/3.0*divU*u_gradphi[i][qp](0) )    // diffusion term
		       + rho*this->_g(0)*u_phi[i][qp]                 // hydrostatic term
		       )*JxW[qp]; 

	    Fv(i) += ( -rho*U*grad_v*u_phi[i][qp]                 // convection term
		       + p*u_gradphi[i][qp](1)                           // pressure term
		       - this->_mu(T)*(u_gradphi[i][qp]*grad_v + u_gradphi[i][qp]*grad_vT
				       - 2.0/3.0*divU*u_gradphi[i][qp](1) )    // diffusion term
		       + rho*this->_g(1)*u_phi[i][qp]                 // hydrostatic term
		       )*JxW[qp];
	    if (this->_dim == 3)
	      {
		Fw(i) += ( -rho*U*grad_w*u_phi[i][qp]                 // convection term
			   + p*u_gradphi[i][qp](2)                           // pressure term
			   - this->_mu(T)*(u_gradphi[i][qp]*grad_w + u_gradphi[i][qp]*grad_wT
					   - 2.0/3.0*divU*u_gradphi[i][qp](2) )    // diffusion term
			   + rho*this->_g(2)*u_phi[i][qp]                 // hydrostatic term
			   )*JxW[qp];
	      }

	    /*
	      if (compute_jacobian && context.get_elem_solution_derivative())
	      {
              libmesh_assert (context.get_elem_solution_derivative() == 1.0);

              for (unsigned int j=0; j != n_u_dofs; j++)
	      {
	      // TODO: precompute some terms like:
	      //   (Uvec*vel_gblgradphivec[j][qp]),
	      //   vel_phi[i][qp]*vel_phi[j][qp],
	      //   (vel_gblgradphivec[i][qp]*vel_gblgradphivec[j][qp])

	      Kuu(i,j) += JxW[qp] *
	      (-_rho*vel_phi[i][qp]*(Uvec*vel_gblgradphivec[j][qp])       // convection term
	      -_rho*vel_phi[i][qp]*graduvec_x*vel_phi[j][qp]             // convection term
	      -_mu*(vel_gblgradphivec[i][qp]*vel_gblgradphivec[j][qp])); // diffusion term
	      Kuv(i,j) += JxW[qp] *
	      (-_rho*vel_phi[i][qp]*graduvec_y*vel_phi[j][qp]);           // convection term

	      Kvv(i,j) += JxW[qp] *
	      (-_rho*vel_phi[i][qp]*(Uvec*vel_gblgradphivec[j][qp])       // convection term
	      -_rho*vel_phi[i][qp]*gradvvec_y*vel_phi[j][qp]             // convection term
	      -_mu*(vel_gblgradphivec[i][qp]*vel_gblgradphivec[j][qp])); // diffusion term
	      Kvu(i,j) += JxW[qp] *
	      (-_rho*vel_phi[i][qp]*gradvvec_x*vel_phi[j][qp]);           // convection term

	      if (_dim == 3)
	      {
	      Kuw(i,j) += JxW[qp] *
	      (-_rho*vel_phi[i][qp]*graduvec_z*vel_phi[j][qp]);           // convection term

	      Kvw(i,j) += JxW[qp] *
	      (-_rho*vel_phi[i][qp]*gradvvec_z*vel_phi[j][qp]);           // convection term

	      Kww(i,j) += JxW[qp] *
	      (-_rho*vel_phi[i][qp]*(Uvec*vel_gblgradphivec[j][qp])       // convection term
	      -_rho*vel_phi[i][qp]*gradwvec_z*vel_phi[j][qp]             // convection term
	      -_mu*(vel_gblgradphivec[i][qp]*vel_gblgradphivec[j][qp])); // diffusion term
	      Kwu(i,j) += JxW[qp] *
	      (-_rho*vel_phi[i][qp]*gradwvec_x*vel_phi[j][qp]);           // convection term
	      Kwv(i,j) += JxW[qp] *
	      (-_rho*vel_phi[i][qp]*gradwvec_y*vel_phi[j][qp]);           // convection term
	      }
	      } // end of the inner dof (j) loop

              // Matrix contributions for the up, vp and wp couplings
              for (unsigned int j=0; j != n_p_dofs; j++)
	      {
	      Kup(i,j) += JxW[qp]*vel_gblgradphivec[i][qp](0)*p_phi[j][qp];
	      Kvp(i,j) += JxW[qp]*vel_gblgradphivec[i][qp](1)*p_phi[j][qp];
	      if (_dim == 3)
	      Kwp(i,j) += JxW[qp]*vel_gblgradphivec[i][qp](2)*p_phi[j][qp];
	      } // end of the inner dof (j) loop

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

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

    return;
  }