elementary_full_domain_FE_film_reynolds Subroutine

public subroutine elementary_full_domain_FE_film_reynolds(fe_f, mat, ke_ij, be_i, ind_e)


Arguments

Type IntentOptional AttributesName
type(FE_FILM), intent(inout) :: fe_f

FE film

type(MAT_SOLV), intent(inout) :: mat

solver type matrices

real(kind=R8), intent(out), dimension(MAX_NNC, MAX_NNC):: ke_ij

elementary matrix

real(kind=R8), intent(inout), dimension(MAX_NNC ):: be_i

elementary rhs member

integer(kind=I4), intent(out), dimension(MAX_NNC ):: ind_e

elementary index member


Calls

proc~~elementary_full_domain_fe_film_reynolds~~CallsGraph proc~elementary_full_domain_fe_film_reynolds elementary_full_domain_FE_film_reynolds proc~solve_fe_film solve_FE_film proc~elementary_full_domain_fe_film_reynolds->proc~solve_fe_film proc~compute_corner_fluxes compute_corner_fluxes proc~elementary_full_domain_fe_film_reynolds->proc~compute_corner_fluxes proc~apply_bc_fe_film_simple apply_bc_FE_film_simple proc~solve_fe_film->proc~apply_bc_fe_film_simple proc~convert_matrice_format convert_matrice_format proc~solve_fe_film->proc~convert_matrice_format omp_get_thread_num omp_get_thread_num proc~solve_fe_film->omp_get_thread_num proc~apply_bc_fe_film apply_bc_FE_film proc~solve_fe_film->proc~apply_bc_fe_film proc~assembly_fe_film_reynolds assembly_FE_film_reynolds proc~solve_fe_film->proc~assembly_fe_film_reynolds proc~solve_syst solve_syst proc~solve_fe_film->proc~solve_syst proc~compute_corner_fluxes->proc~assembly_fe_film_reynolds proc~from_elemental_to_assembled from_elemental_to_assembled proc~convert_matrice_format->proc~from_elemental_to_assembled proc~elementary_assembly_fe_film_reynolds elementary_assembly_FE_film_reynolds proc~assembly_fe_film_reynolds->proc~elementary_assembly_fe_film_reynolds proc~assemble_in_mat_sol assemble_in_mat_sol proc~assembly_fe_film_reynolds->proc~assemble_in_mat_sol proc~analyse_solver analyse_solver proc~solve_syst->proc~analyse_solver proc~factorize_solver factorize_solver proc~solve_syst->proc~factorize_solver proc~close_solver close_solver proc~solve_syst->proc~close_solver proc~solution_solver solution_solver proc~solve_syst->proc~solution_solver proc~freefact_solver freefact_solver proc~solve_syst->proc~freefact_solver proc~init_solver init_solver proc~solve_syst->proc~init_solver ma48_analyse ma48_analyse proc~analyse_solver->ma48_analyse proc~prep_superlu prep_superlu proc~analyse_solver->proc~prep_superlu dmumps dmumps proc~analyse_solver->dmumps proc~s_umfpack_di_symbolic s_umfpack_di_symbolic proc~analyse_solver->proc~s_umfpack_di_symbolic proc~s_umfpack_di_numeric s_umfpack_di_numeric proc~factorize_solver->proc~s_umfpack_di_numeric proc~umfpack_di_free_numeric umfpack_di_free_numeric proc~factorize_solver->proc~umfpack_di_free_numeric proc~factorize_solver->dmumps ma48_factorize ma48_factorize proc~factorize_solver->ma48_factorize proc~compute_prc_tables_reynolds_supg compute_prc_tables_reynolds_supg proc~elementary_assembly_fe_film_reynolds->proc~compute_prc_tables_reynolds_supg proc~umfpack_di_free_symbolic umfpack_di_free_symbolic proc~close_solver->proc~umfpack_di_free_symbolic proc~close_solver->proc~umfpack_di_free_numeric proc~close_solver->dmumps mpi_finalize mpi_finalize proc~close_solver->mpi_finalize proc~close_superlu close_superlu proc~close_solver->proc~close_superlu proc~umfpack_di_report_info umfpack_di_report_info proc~close_solver->proc~umfpack_di_report_info ma48_finalize ma48_finalize proc~close_solver->ma48_finalize proc~solv_superlu solv_superlu proc~solution_solver->proc~solv_superlu ma48_solve ma48_solve proc~solution_solver->ma48_solve proc~solution_solver->dmumps proc~solution_solver->mpi_finalize proc~s_umfpack_di_solve s_umfpack_di_solve proc~solution_solver->proc~s_umfpack_di_solve proc~freefact_solver->proc~umfpack_di_free_numeric proc~free_superlu free_superlu proc~freefact_solver->proc~free_superlu mpi_init mpi_init proc~init_solver->mpi_init ma48_initialize ma48_initialize proc~init_solver->ma48_initialize proc~init_solver->dmumps proc~umfpack_di_report_control umfpack_di_report_control proc~init_solver->proc~umfpack_di_report_control proc~init_superlu init_superlu proc~init_solver->proc~init_superlu proc~init_solver->mpi_finalize proc~umfpack_di_defaults umfpack_di_defaults proc~init_solver->proc~umfpack_di_defaults proc~umfpack_di_numeric umfpack_di_numeric proc~s_umfpack_di_numeric->proc~umfpack_di_numeric proc~solv_superlu->proc~prep_superlu interface~statfree StatFree proc~solv_superlu->interface~statfree interface~destroy_supernode_matrix Destroy_SuperNode_Matrix proc~solv_superlu->interface~destroy_supernode_matrix interface~destroy_compcol_matrix Destroy_CompCol_Matrix proc~solv_superlu->interface~destroy_compcol_matrix interface~dgssvx dgssvx proc~solv_superlu->interface~dgssvx interface~statinit StatInit proc~solv_superlu->interface~statinit interface~dcreate_dense_matrix dCreate_Dense_Matrix proc~prep_superlu->interface~dcreate_dense_matrix interface~dcreate_compcol_matrix dCreate_CompCol_Matrix proc~prep_superlu->interface~dcreate_compcol_matrix interface~c_umfpack_di_free_symbolic c_umfpack_di_free_symbolic proc~umfpack_di_free_symbolic->interface~c_umfpack_di_free_symbolic interface~c_umfpack_di_free_numeric c_umfpack_di_free_numeric proc~umfpack_di_free_numeric->interface~c_umfpack_di_free_numeric interface~c_umfpack_di_report_control c_umfpack_di_report_control proc~umfpack_di_report_control->interface~c_umfpack_di_report_control proc~ni4_up_2d ni4_up_2d proc~compute_prc_tables_reynolds_supg->proc~ni4_up_2d proc~length_width_elem length_width_elem proc~compute_prc_tables_reynolds_supg->proc~length_width_elem proc~dj4 dj4 proc~compute_prc_tables_reynolds_supg->proc~dj4 interface~set_default_options set_default_options proc~init_superlu->interface~set_default_options proc~umfpack_di_solve umfpack_di_solve proc~s_umfpack_di_solve->proc~umfpack_di_solve interface~destroy_dense_matrix Destroy_Dense_Matrix proc~close_superlu->interface~destroy_dense_matrix proc~close_superlu->interface~destroy_compcol_matrix proc~umfpack_di_symbolic umfpack_di_symbolic proc~s_umfpack_di_symbolic->proc~umfpack_di_symbolic interface~c_umfpack_di_report_info c_umfpack_di_report_info proc~umfpack_di_report_info->interface~c_umfpack_di_report_info interface~c_umfpack_di_defaults c_umfpack_di_defaults proc~umfpack_di_defaults->interface~c_umfpack_di_defaults proc~ni4_up_1d ni4_up_1d proc~ni4_up_2d->proc~ni4_up_1d interface~c_umfpack_di_symbolic c_umfpack_di_symbolic proc~umfpack_di_symbolic->interface~c_umfpack_di_symbolic interface~c_umfpack_di_numeric c_umfpack_di_numeric proc~umfpack_di_numeric->interface~c_umfpack_di_numeric interface~c_umfpack_di_solve c_umfpack_di_solve proc~umfpack_di_solve->interface~c_umfpack_di_solve

Contents


Source Code

   subroutine elementary_full_domain_FE_film_reynolds(fe_f, mat, ke_ij, be_i, ind_e)
   implicit none
   type(FE_FILM),    intent(inout)                              :: fe_f   !! *FE film*
   type(MAT_SOLV),   intent(inout)                              :: mat    !! *solver type matrices*
   real(kind=R8),    intent(out  ), dimension(MAX_NNC, MAX_NNC) :: ke_ij  !! *elementary matrix*
   real(kind=R8),    intent(inout), dimension(MAX_NNC         ) :: be_i   !! *elementary rhs member*
   integer(kind=I4), intent(out  ), dimension(MAX_NNC         ) :: ind_e  !! *elementary index member*

      integer(kind=I4) :: e, ee, i, j, ii, nbc, nc
      logical(kind=I4) :: flag

      integer(kind=I4), dimension(MAX_NBS) :: ind_n
      real(kind=R8),    dimension(MAX_NBS) :: sav_vn
      real(kind=R8),    dimension(MAX_NNC) :: bc, dp

      nc = fe_f%m%nc

      ! check of precomputed tables allocation
      if (.not.allocated(fe_f%prc%vcal)) call init_prc_tab(fe_f)

      ! boundary condition initialization
      if (BC_SPLINE) then
         bc(1:nc)  = 0
      else
         bc(1:nc) = be_i(1:nc)
      endif

      j = 1
      do e = 1, fe_f%m%ned
      do ee = 1, fe_f%m%ed(e)%n
         ind_n(j)  = fe_f%m%ed(e)%nm(ee)
         sav_vn(j) = fe_f%vn(ind_n(j), P_N)
         j = j +1
      enddo
      enddo
      nbc = j -1

      ! matrices
      ke_ij = 0._R8
      ind_e = 0
      be_i  = 0._R8

      ! calculation of the elementary matrix
      do i = 1, nc ! number of corners
         ii = fe_f%m%cor(i)
         ind_e(i) = ii

         dp(i) = fe_f%vn(ii, P_N)/1000
         if (BC_SPLINE) then
            bc(i) = dp(i)
         else
            bc(i) = bc(i) + dp(i)
         endif

         if (i==1) flag = .true.

         call solve_FE_film(fe_f = fe_f,     &  !
                             mat = mat,      &  !
                              bc = bc,       &  ! in
                        flag_ass = flag)        !

         flag = .false.

         call compute_corner_fluxes(fe_f = fe_f,      &  !
                                     mat = mat,       &  !
                                      bf = be_i)         ! out

         ke_ij(1:nc, i) = - be_i(1:nc) / dp(i)
         if (BC_SPLINE) then
            bc(i) = 0
         else
            bc(i) = bc(i) -dp(i)
         endif

         do j = 1, nbc
            fe_f%vn(ind_n(j), P_N) = sav_vn(j)
         enddo

      enddo

      ! calculation of the RHS table
      call solve_FE_film(fe_f = fe_f,     &  !
                          mat = mat,      &  !
                           bc = bc,       &  !
                     flag_ass = .false.)     !

      call compute_corner_fluxes(fe_f = fe_f,      &  !
                                  mat = mat,       &  !
                                   bf = be_i)         ! out

      ! update of the matrix: the derivative is the difference
      ! between the residuals divided by the delta p
      do i = 1, nc
         ke_ij(1:nc, i) = ke_ij(1:nc, i) + be_i(1:nc)/dp(i)
      enddo
   return
   endsubroutine elementary_full_domain_FE_film_reynolds