save_fe_f_vtk Subroutine

public subroutine save_fe_f_vtk(fe_f, nom_fic)


Arguments

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

FE_FILM element to store

character(len=*), intent(in) :: nom_fic

filename


Calls

proc~~save_fe_f_vtk~~CallsGraph proc~save_fe_f_vtk save_fe_f_vtk proc~vtk_ini VTK_INI proc~save_fe_f_vtk->proc~vtk_ini proc~vtk_con VTK_CON proc~save_fe_f_vtk->proc~vtk_con interface~vtk_var VTK_VAR proc~save_fe_f_vtk->interface~vtk_var proc~vtk_end VTK_END proc~save_fe_f_vtk->proc~vtk_end interface~vtk_geo VTK_GEO proc~save_fe_f_vtk->interface~vtk_geo proc~vtk_dat VTK_DAT proc~save_fe_f_vtk->proc~vtk_dat proc~upper_case Upper_Case proc~vtk_ini->proc~upper_case proc~getunit GetUnit proc~vtk_ini->proc~getunit proc~vtk_var_vect_r4 VTK_VAR_VECT_R4 interface~vtk_var->proc~vtk_var_vect_r4 proc~vtk_var_vect_i4 VTK_VAR_VECT_I4 interface~vtk_var->proc~vtk_var_vect_i4 proc~vtk_var_vect_r8 VTK_VAR_VECT_R8 interface~vtk_var->proc~vtk_var_vect_r8 proc~vtk_var_text_r4 VTK_VAR_TEXT_R4 interface~vtk_var->proc~vtk_var_text_r4 proc~vtk_var_scal_r4 VTK_VAR_SCAL_R4 interface~vtk_var->proc~vtk_var_scal_r4 proc~vtk_var_text_r8 VTK_VAR_TEXT_R8 interface~vtk_var->proc~vtk_var_text_r8 proc~vtk_var_scal_i4 VTK_VAR_SCAL_I4 interface~vtk_var->proc~vtk_var_scal_i4 proc~vtk_var_scal_r8 VTK_VAR_SCAL_R8 interface~vtk_var->proc~vtk_var_scal_r8 proc~vtk_geo_unst_r8 VTK_GEO_UNST_R8 interface~vtk_geo->proc~vtk_geo_unst_r8 proc~vtk_geo_unst_r4 VTK_GEO_UNST_R4 interface~vtk_geo->proc~vtk_geo_unst_r4 proc~vtk_geo_rect_r4 VTK_GEO_RECT_R4 interface~vtk_geo->proc~vtk_geo_rect_r4 proc~vtk_geo_strp_r4 VTK_GEO_STRP_R4 interface~vtk_geo->proc~vtk_geo_strp_r4 proc~vtk_geo_strg_r4 VTK_GEO_STRG_R4 interface~vtk_geo->proc~vtk_geo_strg_r4 proc~vtk_geo_rect_r8 VTK_GEO_RECT_R8 interface~vtk_geo->proc~vtk_geo_rect_r8 proc~vtk_geo_strg_r8 VTK_GEO_STRG_R8 interface~vtk_geo->proc~vtk_geo_strg_r8 proc~vtk_geo_strp_r8 VTK_GEO_STRP_R8 interface~vtk_geo->proc~vtk_geo_strp_r8 proc~vtk_dat->proc~upper_case proc~vtk_var_vect_r4->proc~upper_case proc~vtk_var_vect_r8->proc~upper_case

Contents

Source Code


Source Code

   subroutine save_fe_f_vtk(fe_f, nom_fic)
   implicit none
   type(FE_FILM),    intent(in) :: fe_f      !! [[FE_FILM]] *element to store*
   character(len=*), intent(in) :: nom_fic   !! *filename*

      integer(kind=I4) :: E_IO   ! parametre identification vtk
      integer(kind=I4) :: i, k

      real(kind=R8),    dimension(:), allocatable :: z            ! tableau supplementaire pour la sortie
      integer(kind=I4), dimension(:), allocatable :: connec, tipo

      allocate (z(fe_f%m%n), connec(fe_f%m%ne*(1+4)), tipo(fe_f%m%ne))
      z    = 0._R8
      tipo = 9       ! elements a quatre noeuds: 9 dans le formalisme vtk

      k = 0
      do i = 1, fe_f%m%ne
         connec(k+1) = 4
         connec(k+2) = fe_f%m%con(i, 1) -1
         connec(k+3) = fe_f%m%con(i, 2) -1
         connec(k+4) = fe_f%m%con(i, 3) -1
         connec(k+5) = fe_f%m%con(i, 4) -1
         k=k+5
      enddo
      ! ouverture du fichier
      E_IO = VTK_INI(output_format = 'binary', &
                          filename = nom_fic, &
                             title = 'Resultats fe_f fluide', &
                     mesh_topology = 'UNSTRUCTURED_GRID')
      ! ecriture des coordonnees
      E_IO = VTK_GEO(NN = fe_f%m%n, &
                      X = fe_f%m%x, &
                      Y = fe_f%m%y, &
                      Z = Z)
      ! definition des cellules
      E_IO = VTK_CON(NC      = fe_f%m%ne, &
                     connect = connec,    &
                   cell_type = tipo)
      ! choix du type de variable (elementaire)
      E_IO = VTK_DAT(NC_NN   = fe_f%m%ne, &
                var_location = 'cell')
      ! passage des variables cellules
      do i = 1, fe_f%n_vc
          E_IO = VTK_VAR(NC_NN = fe_f%m%ne,       &
                       varname = fe_f%vc_name(i), &
                           var = fe_f%vc(:,i))
      enddo
      ! choix du type de variable (nodale)
      E_IO = VTK_DAT(NC_NN   = fe_f%m%n, &
                var_location = 'node')
      ! passage des variables nodales
      do i = 1, fe_f%n_vn
          E_IO = VTK_VAR(NC_NN = fe_f%m%n,        &
                       varname = fe_f%vn_name(i), &
                           var = fe_f%vn(:, i))
      enddo
      ! fermeture du fichier vtk
      E_IO = VTK_END()
      ! liberation des tableaux
      deallocate(z, connec, tipo)
   return
   endsubroutine save_fe_f_vtk