solve_FE_film Subroutine

public subroutine solve_FE_film(fe_f, mat, bc, flag_ass)


Arguments

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

FE film data type

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

matrices for solving

real(kind=R8), intent(in), dimension(MAX_NNC):: bc

boundary conditions at the corners

logical(kind=I4), intent(in), optional :: flag_ass

optional flag for assembly


Calls

proc~~solve_fe_film~~CallsGraph proc~solve_fe_film solve_FE_film 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~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

Called by

proc~~solve_fe_film~~CalledByGraph proc~solve_fe_film solve_FE_film proc~elementary_full_domain_fe_film_reynolds elementary_full_domain_FE_film_reynolds proc~elementary_full_domain_fe_film_reynolds->proc~solve_fe_film proc~multi_scale_solve_fe_film multi_scale_solve_fe_film proc~multi_scale_solve_fe_film->proc~solve_fe_film proc~solve_fe_prob solve_fe_prob proc~solve_fe_prob->proc~solve_fe_film proc~solve_ms_prob solve_ms_prob proc~solve_ms_prob->proc~multi_scale_solve_fe_film proc~test_rough_fe test_rough_fe proc~test_rough_fe->proc~solve_fe_prob proc~test_bearing_x_fe test_bearing_x_fe proc~test_bearing_x_fe->proc~solve_fe_prob proc~test_pocket_fe test_pocket_fe proc~test_pocket_fe->proc~solve_fe_prob proc~test_bearing_y_fe test_bearing_y_fe proc~test_bearing_y_fe->proc~solve_fe_prob proc~test_slider_fe test_slider_fe proc~test_slider_fe->proc~solve_fe_prob proc~run_test run_test proc~run_test->proc~test_rough_fe proc~run_test->proc~test_bearing_x_fe proc~run_test->proc~test_pocket_fe proc~run_test->proc~test_bearing_y_fe proc~run_test->proc~test_slider_fe proc~test_slider_ms test_slider_ms proc~run_test->proc~test_slider_ms proc~test_rough_ms test_rough_ms proc~run_test->proc~test_rough_ms proc~test_slider_ms->proc~solve_ms_prob proc~test_rough_ms->proc~solve_ms_prob program~main main program~main->proc~run_test

Contents

Source Code


Source Code

   subroutine solve_FE_film(fe_f, mat, bc, flag_ass)
   implicit none
   type(FE_FILM),    intent(inout)                     :: fe_f     !! *FE film data type*
   type(MAT_SOLV),   intent(inout)                     :: mat      !! *matrices for solving*
   real(kind=R8),    intent(in   ), dimension(MAX_NNC) :: bc       !! *boundary conditions at the corners*
   logical(kind=I4), intent(in   ), optional           :: flag_ass !! *optional flag for assembly*

      logical(kind=I4) :: decomp
      integer(kind=I4) :: ass_c, i, it, e
      logical(kind=I4) :: conv
      real(kind=R8)    :: error
      integer(kind=I4), dimension(2) :: compt

      ! solution parameters
      ass_c = ASS
      if (present(flag_ass)) then
         if (.not.flag_ass) ass_c = NO_ASS
      endif
      decomp = (ass_c == ASS)

      ! apply boundary conditions
      if (BC_SPLINE) then
         call apply_bc_FE_film(fe_f = fe_f,    &
                                 bc = bc)
      else
         call apply_bc_FE_film_simple(fe_f = fe_f,    &
                                        bc = bc)
      endif

      if (VERBOSE >= 20) write(OPU,*) 'bc applied'

      ! update fluid properties
      do i = 1, fe_f%m%n
         fe_f%vn(i, RHO_N)    = fe_f%data_f%fl%rho(   fe_f%vn(i, P_N) , fe_f%vn(i, T_N))
         fe_f%vn(i, DRHODP_N) = fe_f%data_f%fl%drhodp(fe_f%vn(i, P_N) , fe_f%vn(i, T_N))
      enddo
      if (VERBOSE >= 20) write(OPU,*) 'fluid properties updated'

      if (mat%first) then
         mat%nn = fe_f%m%n
         mat%ne = fe_f%m%ne
         call solve_syst(mat = mat, &
                        step = 'ini')
         if (VERBOSE >= 10) write(OPU,*) 'Matrix initialized, thread ', omp_get_thread_num()

         ! matrices allocation
         compt(:) = 0
         do e = 1, fe_f%m%ne
            compt(1) = compt(1) +fe_f%m%el_n(e)       ! for each element, the number of lines is added
            do i = 1, fe_f%m%el_t(e)
               compt(2) = compt(2) + fe_f%m%el_t(e)   ! " ", for each node, the number of contributions is added
            enddo
         enddo
         mat%nvar = compt(1)
         mat%nt   = compt(2)
         if (.not.allocated(mat%eltvar)) allocate( mat%eltvar(mat%nvar  ) )
         if (.not.allocated(mat%a_elt )) allocate( mat%a_elt( mat%nt    ) )
         if (.not.allocated(mat%eltptr)) allocate( mat%eltptr(mat%ne +1 ) )
      endif

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

      ! convergence is false at zero iteration
      conv = .false.
      it   = 0
      ! solution loop
      do
         if (VERBOSE >= 20) write(OPU,*) '   Loop FE ', it
         if (conv) exit

         ! assembly of the system
         call assembly_FE_film_reynolds(fe_f, mat, ass_c)
         if (VERBOSE >= 30) write(OPU,*) '   System assembled, thread ', omp_get_thread_num()

         if (ass_c == ASS) then
            ! some stuff can be saved here, provided the reloading of jptr, irow, ... (instead of convert_matrice_format)
            call convert_matrice_format(mat)
            if (VERBOSE >= 30) write(OPU,*) '   Matrix converted, thread ', omp_get_thread_num()
         endif

         if (mat%first) then
            call solve_syst(mat = mat, &
                           step = 'ana')
            mat%first = .false.
            if (VERBOSE >= 10) write(OPU,*) '   Matrix analyzed, thread ', omp_get_thread_num()
         endif

         ! solution of the system
         if (ass_c == ASS) then
            call solve_syst(mat = mat, &
                           step = 'fac')
            if (VERBOSE >= 30) write(OPU,*) '   Matrix factorized, thread ', omp_get_thread_num()
         endif

         call solve_syst(mat = mat, &
                        step = 'sol')
         if (VERBOSE >= 30) write(OPU,*) '   System solved, thread ', omp_get_thread_num()
         if (ass_c == ASS) then
            call solve_syst(mat = mat, &
                           step = 'fre')
            if (VERBOSE >= 30) write(OPU,*) '   Matrix factors freed, thread ', omp_get_thread_num()
         endif

         ! error computation
         error = (sum(mat%x ** 2) / sum(fe_f%vn(:, P_N) ** 2)) ** (0.5_R8)
         it = it + 1
         if (VERBOSE >= 20) write(OPU,*) '   Error ', error

         ! convergence check
         if (error <= fe_f%num_p%eps) conv = .true.

         ! update of variables
         if (fe_f%data_f%fl%fluid_type == MIXT) then
            do i = 1, fe_f%m%n
               if (mat%x(i) < 0.) then
                  fe_f%vn(i, RHO_N) = fe_f%vn(i, RHO_N) + fe_f%vn(i, DRHODP_N) * mat%x(i) * fe_f%num_p%relax
                  if (fe_f%vn(i, RHO_N) < 0.) then
                      fe_f%vn(i,   P_N) = fe_f%data_f%fl%P_0/100
                  else
                      fe_f%vn(i,   P_N) = fe_f%data_f%fl%pres( fe_f%vn(i, RHO_N), & !
                                                               fe_f%vn(i,   T_N) )
                  endif
               else
                  fe_f%vn(i, P_N) = fe_f%vn(i, P_N) + mat%x(i) * fe_f%num_p%relax
               endif
            enddo
         else
            fe_f%vn(:, P_N) = fe_f%vn(:, P_N) + mat%x * fe_f%num_p%relax
         endif

         ! check pressure
         if ( fe_f%data_f%fl%fluid_type == GP ) then
            if (minval(fe_f%vn(:, P_N)) < 0._R8) then
               if (VERBOSE >= 30) write(OPU,*) 'P negative'
            endif
            where (fe_f%vn(:, P_N) < 0._R8) fe_f%vn(:, P_N) = fe_f%data_f%fl%p_0 / 1.e3_R8
         endif

         do i = 1, fe_f%m%n
            fe_f%vn(i, RHO_N)    = fe_f%data_f%fl%rho(   fe_f%vn(i, P_N) , fe_f%vn(i, T_N))
            fe_f%vn(i, DRHODP_N) = fe_f%data_f%fl%drhodp(fe_f%vn(i, P_N) , fe_f%vn(i, T_N))
         enddo

         if (it >= fe_f%num_p%it_max) then
            conv = .true.
            if (VERBOSE >= 20) write(OPU,*) 'maximum number of iteration reached before convergence'
         endif
      enddo

   return
   endsubroutine solve_FE_film