multi_scale_assembly_fe_film_reynolds Subroutine

private subroutine multi_scale_assembly_fe_film_reynolds(ms_fe_f, ms_mat, ass_c)


Arguments

Type IntentOptional AttributesName
type(MS_FE_FILM), intent(inout) :: ms_fe_f
type(MS_MAT_SOLV), intent(inout) :: ms_mat

solver type matrices table

integer(kind=I4), intent(in) :: ass_c

assembly type


Calls

proc~~multi_scale_assembly_fe_film_reynolds~~CallsGraph proc~multi_scale_assembly_fe_film_reynolds multi_scale_assembly_fe_film_reynolds omp_get_thread_num omp_get_thread_num proc~multi_scale_assembly_fe_film_reynolds->omp_get_thread_num

Called by

proc~~multi_scale_assembly_fe_film_reynolds~~CalledByGraph proc~multi_scale_assembly_fe_film_reynolds multi_scale_assembly_fe_film_reynolds proc~multi_scale_solve_fe_film multi_scale_solve_fe_film proc~multi_scale_solve_fe_film->proc~multi_scale_assembly_fe_film_reynolds proc~solve_ms_prob solve_ms_prob proc~solve_ms_prob->proc~multi_scale_solve_fe_film 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

   subroutine multi_scale_assembly_fe_film_reynolds(ms_fe_f, ms_mat, ass_c)
   implicit none
   type(MS_FE_FILM),  intent(inout) :: ms_fe_f     !!
   type(MS_MAT_SOLV), intent(inout) :: ms_mat      !! *solver type matrices table*
   integer(kind=I4),  intent(in)    :: ass_c       !! *assembly type*

      real(kind=R8),    dimension(:,:,:), allocatable :: ke_ij
      real(kind=R8),    dimension(:,:)  , allocatable ::  be_i
      integer(kind=I4), dimension(:,:)  , allocatable :: ind_e
      integer(kind=I4), dimension(2)                  :: compt

      integer(kind=I4) :: e, i, ii, ne
      real(kind=R8)    :: val_t
      logical(kind=I4) :: first_assembly

      ne = ms_fe_f%ts_fe_f%m%ne

      allocate( ke_ij(ne, MAX_NNC, MAX_NNC) )
      allocate(  be_i(ne, MAX_NNC) )
      allocate( ind_e(ne, MAX_NNC) )

      ! assembly
      compt(:)        = 1
      ms_mat%ts_mat%b = 0._R8

      if (.not.allocated(ms_mat%ass_loc_in_mat)) then
         allocate(ms_mat%ass_loc_in_mat(ne))
         ms_mat%ass_loc_in_mat = -1
      endif

      ! check for first assembly: save the assembly location in ass_loc_in_mat
      if (ms_mat%ass_loc_in_mat(1) == -1) then
         first_assembly = .true.
         ms_mat%ass_loc_in_mat(1) = 1
      else
         first_assembly = .false.
      endif

      do e = 1, ne
         ! copy of the boundary conditions
         do i = 1, 4
            ii = ms_fe_f%ts_fe_f%m%con(e, i)
            be_i(e, i) = ms_fe_f%ts_fe_f%vn(ii, P_N)
         enddo
      enddo

      ! elementary matrices calculation
      !-------------------------------------------
      !open mp instructions (parallele computation)
      !$omp parallel
      !$omp do schedule(runtime)
      do e = 1, ms_fe_f%ts_fe_f%m%ne
         if (VERBOSE >= 3) write(OPU,*) "MS elements ", e, ' thread ', omp_get_thread_num()

         call elementary_full_domain_FE_film_reynolds(fe_f = ms_fe_f%bs_fe_f(e), &
                                                       mat = ms_mat%bs_mat(e),   &
                                                     ke_ij = ke_ij(e, :, :),     &
                                                      be_i = be_i( e, :),        &
                                                     ind_e = ind_e(e, :))
      enddo
      !$OMP end do
      !$OMP end parallel
      !end open mp instructions
      !-------------------------------------------

      do e = 1, ne
         ! consideration of the bc condition
         do i = 1, 4
            ii = ms_fe_f%ts_fe_f%m%con(e, i)
            ind_e(e, i) = ii
            if (ms_fe_f%ts_fe_f%bc(ii, REY) == 0) then
               val_t = ke_ij(e, i, i)
               ke_ij(e, i, :) = 0._R8
               ke_ij(e, i, i) = val_t
               be_i( e, i)    = 0._R8
            endif
         enddo
         ! copy of the rhs member
         do i = 1, 4
            ii = ms_fe_f%ts_fe_f%m%con(e, i)
            ms_mat%ts_mat%b(ii) = ms_mat%ts_mat%b(ii) + be_i(e, i)
         enddo
      enddo

      ! assembly of the elemental matrix in the solver matrix
      if (ass_c == ASS) then

         do e = 1, ne

            call assemble_in_mat_sol(mat = ms_mat%ts_mat,      &
                                     num = e,                  &
                                    nelt = 4,                  &
                                   nline = 4,                  &
                                    tind = ind_e(e, :),        &
                                   m_elt = ke_ij(e, :, :),     &
                                   compt = compt)
         enddo
         ms_mat%ts_mat%eltptr(1) = 1

      endif

      deallocate( ke_ij, be_i, ind_e )

   return
   endsubroutine multi_scale_assembly_fe_film_reynolds