Subroutine to solve a FE_FILM
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
type(FE_FILM), | intent(inout) | :: | fe_f | FE film |
||
real(kind=R8), | intent(out), | dimension(MAX_NNE, MAX_NNE) | :: | ke_ij | elementary matrix |
|
real(kind=R8), | intent(out), | dimension(MAX_NNE) | :: | be_i | elementary rhs member |
|
integer(kind=I4), | intent(out), | dimension(MAX_NNE) | :: | ind_e | elementary index member |
|
integer(kind=I4), | intent(in) | :: | e | element number |
||
integer(kind=I4), | intent(in) | :: | ass_c | assembly type |
subroutine elementary_assembly_FE_film_reynolds(fe_f, ke_ij, be_i, ind_e, e, ass_c)
implicit none
type(FE_FILM), intent(inout) :: fe_f !! *FE film*
real(kind=R8), intent(out ), dimension(MAX_NNE, MAX_NNE) :: ke_ij !! *elementary matrix*
real(kind=R8), intent(out ), dimension(MAX_NNE) :: be_i !! *elementary rhs member*
integer(kind=I4), intent(out ), dimension(MAX_NNE) :: ind_e !! *elementary index member*
integer(kind=I4), intent(in ) :: e !! *element number*
integer(kind=I4), intent(in ) :: ass_c !! *assembly type*
real(kind=R8), dimension(MAX_NNG, MAX_NNG) :: vni4x, vni4y
integer(kind=I4), dimension(MAX_NNE) :: con
integer(kind=I4) :: i, j, ii, jj, ng, ne
ng = fe_f%prc%ng
ne = fe_f%m%el_t(e)
! check of precomputed tables allocation
!~ if (.not.allocated(fe_f%prc%vcal)) call init_prc_tab(fe_f)
! assembly
call compute_prc_tables_reynolds_supg(fe_f, e)
con(1:ne) = fe_f%m%con(e, 1:ne)
ke_ij = 0._R8
ind_e = 0
be_i = 0._R8
do i = 1, ne
ii = con(i)
ind_e(i) = ii
! case of nodes on the boundary
if ((fe_f%bc(ii, REY) == 0) .and. (ass_c /= NO_BC)) then
ke_ij(i, i)=1.0e10
! case of nodes where the pressure is unknown
else
vni4x(1:ng, 1:ng) = fe_f%prc%vni4x(i, 1:ng, 1:ng)
vni4y(1:ng, 1:ng) = fe_f%prc%vni4y(i, 1:ng, 1:ng)
be_i(i) = be_i(i) - sum( vni4x(1:ng, 1:ng) * (fe_f%prc%vcal(6, 1:ng, 1:ng) - fe_f%prc%vcal(8, 1:ng, 1:ng)) ) & !
- sum( vni4y(1:ng, 1:ng) * (fe_f%prc%vcal(7, 1:ng, 1:ng) - fe_f%prc%vcal(9, 1:ng, 1:ng)) )
do j = 1, fe_f%m%el_t(e)
jj = con(j)
ke_ij(i, j) = sum( (vni4x(1:ng, 1:ng) * fe_f%prc%vni4x(j, 1:ng, 1:ng) &
+vni4y(1:ng, 1:ng) * fe_f%prc%vni4y(j, 1:ng, 1:ng))* fe_f%prc%vcal(2, 1:ng, 1:ng) ) &
+ sum( (vni4x(1:ng, 1:ng) * fe_f%prc%vcal(10, 1:ng, 1:ng) &
+ vni4y(1:ng, 1:ng) * fe_f%prc%vcal(11, 1:ng, 1:ng))* fe_f%prc%vni4(j, 1:ng, 1:ng))* fe_f%vn(jj, DRHODP_N) &
- sum( (vni4x(1:ng, 1:ng) * fe_f%prc%vcal( 3, 1:ng, 1:ng) &
+ vni4y(1:ng, 1:ng) * fe_f%prc%vcal( 4, 1:ng, 1:ng)) &
* fe_f%prc%vni4d(j, 1:ng, 1:ng) * fe_f%vn(jj, DRHODP_N) ) !better convergence
enddo
endif
enddo
return
endsubroutine elementary_assembly_FE_film_reynolds