solve_ms_prob Subroutine

private subroutine solve_ms_prob()


Arguments

None

Calls

proc~~solve_ms_prob~~CallsGraph proc~solve_ms_prob solve_ms_prob proc~save_ms_field save_ms_field proc~solve_ms_prob->proc~save_ms_field proc~multi_scale_solve_fe_film multi_scale_solve_fe_film proc~solve_ms_prob->proc~multi_scale_solve_fe_film proc~ms_fe_f_2_mat ms_fe_f_2_mat proc~solve_ms_prob->proc~ms_fe_f_2_mat proc~save_ms_field->proc~ms_fe_f_2_mat proc~write_surf write_surf proc~save_ms_field->proc~write_surf proc~empty empty proc~save_ms_field->proc~empty proc~init_scal init_scal proc~save_ms_field->proc~init_scal proc~convert_matrice_format convert_matrice_format proc~multi_scale_solve_fe_film->proc~convert_matrice_format proc~apply_bc_fe_film_simple apply_bc_FE_film_simple proc~multi_scale_solve_fe_film->proc~apply_bc_fe_film_simple omp_get_thread_num omp_get_thread_num proc~multi_scale_solve_fe_film->omp_get_thread_num proc~multi_scale_assembly_fe_film_reynolds multi_scale_assembly_fe_film_reynolds proc~multi_scale_solve_fe_film->proc~multi_scale_assembly_fe_film_reynolds proc~solve_syst solve_syst proc~multi_scale_solve_fe_film->proc~solve_syst proc~solve_fe_film solve_FE_film proc~multi_scale_solve_fe_film->proc~solve_fe_film proc~lower lower proc~write_surf->proc~lower proc~scal2surf scal2surf proc~write_surf->proc~scal2surf proc~get_unit get_unit proc~write_surf->proc~get_unit proc~surf2scal surf2scal proc~write_surf->proc~surf2scal proc~from_elemental_to_assembled from_elemental_to_assembled proc~convert_matrice_format->proc~from_elemental_to_assembled proc~multi_scale_assembly_fe_film_reynolds->omp_get_thread_num 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 proc~solve_fe_film->proc~convert_matrice_format proc~solve_fe_film->proc~apply_bc_fe_film_simple proc~solve_fe_film->omp_get_thread_num proc~solve_fe_film->proc~solve_syst 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~init_scal->proc~empty 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~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 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~c_f_string c_f_string proc~surf2scal->proc~c_f_string 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 proc~c_f_string->proc~empty interface~c_umfpack_di_report_control c_umfpack_di_report_control proc~umfpack_di_report_control->interface~c_umfpack_di_report_control interface~set_default_options set_default_options proc~init_superlu->interface~set_default_options 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_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_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~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 proc~ni4_up_1d ni4_up_1d proc~ni4_up_2d->proc~ni4_up_1d

Called by

proc~~solve_ms_prob~~CalledByGraph proc~solve_ms_prob solve_ms_prob proc~test_slider_ms test_slider_ms proc~test_slider_ms->proc~solve_ms_prob proc~test_rough_ms test_rough_ms proc~test_rough_ms->proc~solve_ms_prob proc~run_test run_test proc~run_test->proc~test_slider_ms proc~run_test->proc~test_rough_ms program~main main program~main->proc~run_test

Contents

Source Code


Source Code

   subroutine solve_ms_prob
   implicit none
      real(kind=R8), dimension(:, :), allocatable :: tab

      call system_clock(count=cinit)
      call cpu_time(t1)
         call multi_scale_solve_fe_film(ms_fe_f, ms_mat, bc)
      call cpu_time(t2)
      call system_clock(count=cend,count_rate=cr)

      write(OPU,*)          'MS cpu time (s):',          char(9), t2 - t1
      write(unit_num_res,*) 'MS_cpu_time_(s):',          char(9), t2 - t1
      write(OPU,*)          'MS real comp time (s):',    char(9), real(cend - cinit) / real(cr)
      write(unit_num_res,*) 'MS_real_comp_time_(s):',    char(9), real(cend - cinit) / real(cr)

      write(OPU,*)          'load MS_FE (N):',           char(9), ms_fe_f%ms_fz()
      write(unit_num_res,*) 'load_MS_FE_(N):',           char(9), ms_fe_f%ms_fz()
      write(OPU,*)          'fric/x MS FE (N):',         char(9), ms_fe_f%ms_fx()
      write(unit_num_res,*) 'fric/x_MS_FE_(N):',         char(9), ms_fe_f%ms_fx()
      write(OPU,*)          'fric/y MS FE (N):',         char(9), ms_fe_f%ms_fy()
      write(unit_num_res,*) 'fric/y_MS_FE_(N):',         char(9), ms_fe_f%ms_fy()

      call save_ms_fe_f_vtk(ms_fe_f, trim(ms_vtk))

      call save_profile_x_fe(ms_fe_f%ts_fe_f, trim(prof_ts), lx, ly/2)
      call save_profile_x_ms(ms_fe_f,         trim(prof_bs), lx, ly/2)

      call save_ms_field(ms_fe_f = ms_fe_f,                                   & !
                       file_name = "out/"//res_dir//"/"//"ms_pressure.sur",   & !
                            code = P_N,                                       & !
                           nodal = .true.                                     )

      if (compare_solution_file /= 0) then
         if ( allocated(tab_sol) ) deallocate(tab_sol)
         call read_surf(nom_fic = trim(pressure_solution_file),   & !
                             mu = -1._R8,                         & !
                             sq = -1._R8,                         & !
                          tab_s = tab_sol,                        & !
                           scal = scal_tmp)

         if ( allocated(tab_s) ) deallocate(tab_s)
         call ms_fe_f_2_mat(ms_fe_f = ms_fe_f,  & !
                               code = P_N,      & !
                              nodal = .true.,   & !
                                mat = tab_s)

         allocate( tab(1:nx, 1:ny) )
         tab(1:nx, 1:ny) = abs(tab_s(1:nx, 1:ny) -tab_sol(1:nx, 1:ny))

         scal_tmp%zlength_unit = 'Pa'
         scal_tmp%dz_unit      = 'Pa'
         call write_surf(nom_fic = "out/"//res_dir//"/"//"compare_pressure.sur",    & !
                           tab_s = tab(1:nx, 1:ny),                                 & !
                            scal = scal_tmp)

         where (tab_s < data_f%fl%p_0 .and. tab_sol > data_f%fl%p_0)
            tab = -1._R8
         elsewhere (tab_s > data_f%fl%p_0 .and. tab_sol < data_f%fl%p_0)
            tab = +1._R8
         elsewhere (tab_s < data_f%fl%p_0 .and. tab_sol < data_f%fl%p_0)
            tab =  0._R8
         elsewhere
            tab = -2._R8
         endwhere

         scal_tmp%zlength_unit = '  '
         scal_tmp%dz_unit      = '  '
         call write_surf(nom_fic = "out/"//res_dir//"/"//"compare_cavitation.sur",  & !
                           tab_s = tab(1:nx, 1:ny),                                 & !
                            scal = scal_tmp)

      endif

      if ( allocated(tab    ) ) deallocate(tab    )
      if ( allocated(tab_s  ) ) deallocate(tab_s  )
      if ( allocated(tab_sol) ) deallocate(tab_sol)
   return
   endsubroutine solve_ms_prob