update_variable_pressure_operator_cache Subroutine

public subroutine update_variable_pressure_operator_cache(mesh, flow, bc, transport)

Rebuild per-step variable-density pressure coefficients.

Arguments

Type IntentOptional Attributes Name
type(mesh_t), intent(in) :: mesh
type(flow_mpi_t), intent(in) :: flow
type(bc_set_t), intent(in) :: bc
type(transport_properties_t), intent(in) :: transport

Source Code

   subroutine update_variable_pressure_operator_cache(mesh, flow, bc, transport)
      type(mesh_t), intent(in) :: mesh
      type(flow_mpi_t), intent(in) :: flow
      type(bc_set_t), intent(in) :: bc
      type(transport_properties_t), intent(in) :: transport

      integer :: c, lf, f, nb
      real(rk) :: dist, coeff, beta_face

      call profiler_start('Pressure_VarCoeff_Update')
      call ensure_pressure_operator_cache(mesh, flow, bc)

      if (.not. allocated(pressure_cache%var_coeff) .or. .not. allocated(pressure_cache%var_diag)) then
         call fatal_error('numerics', 'variable pressure coefficient cache is not allocated')
      end if

      do c = flow%first_cell, flow%last_cell
         pressure_cache%var_diag(c) = zero
         do lf = 1, mesh%ncell_faces(c)
            pressure_cache%var_coeff(lf, c) = zero
            f = mesh%cell_faces(lf, c)
            nb = face_effective_neighbor(mesh, bc, f, c)

            if (nb > 0) then
               dist = face_normal_distance(mesh, bc, f, c, nb)
               beta_face = pressure_face_inverse_density(transport, c, nb)
               coeff = beta_face * mesh%faces(f)%area / dist
               pressure_cache%nb(lf, c) = nb
               pressure_cache%var_coeff(lf, c) = coeff
               pressure_cache%var_diag(c) = pressure_cache%var_diag(c) + coeff
            else if (boundary_pressure_type(mesh, bc, f) == bc_dirichlet) then
               dist = face_normal_distance(mesh, bc, f, c, 0)
               beta_face = pressure_face_inverse_density(transport, c, 0)
               coeff = beta_face * mesh%faces(f)%area / dist
               pressure_cache%var_diag(c) = pressure_cache%var_diag(c) + coeff
            end if
         end do

         if (pressure_cache%var_diag(c) <= tiny_safe) then
            call fatal_error('numerics', 'non-positive variable-density pressure diagonal')
         end if
      end do

      call profiler_stop('Pressure_VarCoeff_Update')
   end subroutine update_variable_pressure_operator_cache