!-------------------------------------------------------------------------------

!                      Code_Saturne version 4.0.2
!                      --------------------------
! This file is part of Code_Saturne, a general-purpose CFD tool.
!
! Copyright (C) 1998-2015 EDF S.A.
!
! This program is free software; you can redistribute it and/or modify it under
! the terms of the GNU General Public License as published by the Free Software
! Foundation; either version 2 of the License, or (at your option) any later
! version.
!
! This program is distributed in the hope that it will be useful, but WITHOUT
! ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
! FOR A PARTICULAR PURPOSE.  See the GNU General Public License for more
! details.
!
! You should have received a copy of the GNU General Public License along with
! this program; if not, write to the Free Software Foundation, Inc., 51 Franklin
! Street, Fifth Floor, Boston, MA 02110-1301, USA.

!-------------------------------------------------------------------------------

!===============================================================================
! Purpose:
! -------

!> \file cs_user_initialization.f90
!>
!> \brief Initialize variables
!>
!> This subroutine is called at beginning of the computation
!> (restart or not) before the loop time step.
!>
!> This subroutine enables to initialize or modify (for restart)
!> unkown variables and time step values.
!>
!> \c rom and \c viscl values are equal to \c ro0 and \c viscl0 or initialize
!> by reading the restart file.
!> variable diffusivity and cp variables (when there are defined) have no value
!> excepted if they are read from a restart file.
!>
!> Modification of the behaviour law of physical quantities (rom, viscl,
!> viscls, cp) is not done here. It is the purpose of the user subroutine
!> \ref usphyv (in cs_user_physical_properties.f90)
!>
!> \section cs_user_initialization_cell_id Cells identification
!>
!> Cell value field ids
!>
!> - Density:                        \c iprpfl(irom)
!> - Dynamic molecular viscosity:    \c iprpfl(iviscl)
!> - Turbulent viscosity:            \c iprpfl(ivisct)
!> - Specific heat:                  \c iprpfl(icp)
!> - Diffusivity(lambda):            \c field_get_key_int(ivarfl(isca(iscal)),
!>                                      kivisl, ...)
!>
!> Cells may be identified using the \ref getcel subroutine.
!> The syntax of this subroutine is described in the
!> \ref cs_user_boundary_conditions subroutine,
!> but a more thorough description can be found in the user guide.
!
!-------------------------------------------------------------------------------

!-------------------------------------------------------------------------------
! Arguments
!______________________________________________________________________________.
!  mode           name          role                                           !
!______________________________________________________________________________!
!> \param[in]     nvar          total number of variables
!> \param[in]     nscal         total number of scalars
!> \param[in]     dt            time step (per cell)
!_______________________________________________________________________________

subroutine cs_user_initialization &
 ( nvar   , nscal  ,                                              &
   dt     )

!===============================================================================

!===============================================================================
! Module files
!===============================================================================

use paramx
use pointe
use numvar
use optcal
use cstphy
use cstnum
use entsor
use parall
use period
use ppppar
use ppthch
use coincl
use cpincl
use ppincl
use atincl
use ctincl
use elincl
use ppcpfu
use cs_coal_incl
use cs_fuel_incl
use mesh
use field
use turbomachinery

!===============================================================================

implicit none

! Arguments

integer          nvar   , nscal

double precision dt(ncelet)

! Local variables
integer          iel
double precision d2s3
double precision zent,xuent,xvent,xkent,xeent,tpent

! INSERT_VARIABLE_DEFINITIONS_HERE
double precision, dimension(:,:), pointer :: cvar_vel
double precision, dimension(:), pointer :: cvar_pr

integer, allocatable, dimension(:) :: lstelt

double precision, dimension(:), pointer :: cvar_k, cvar_ep, cvar_phi, cvar_fb
double precision, dimension(:), pointer :: cvar_omg, cvar_nusa
double precision, dimension(:), pointer :: cvar_r11, cvar_r22, cvar_r33
double precision, dimension(:), pointer :: cvar_r12, cvar_r13, cvar_r23
double precision, dimension(:), pointer :: cvar_scalt


!===============================================================================


!===============================================================================
! Initialization
!===============================================================================

allocate(lstelt(ncel)) ! temporary array for cells selection

! INSERT_MAIN_CODE_HERE

!< [init]
d2s3 = 2.d0/3.d0

if (isuite.eq.0) then
! Call the field of velocity and pressure
  call field_get_val_v(ivarfl(iu), cvar_vel)
  call field_get_val_s(ivarfl(ipr), cvar_pr)
! Call the field of turbulent scales.
  if (itytur.eq.2) then
    call field_get_val_s(ivarfl(ik), cvar_k)
    call field_get_val_s(ivarfl(iep), cvar_ep)
  elseif (itytur.eq.3) then
    call field_get_val_s(ivarfl(ir11), cvar_r11)
    call field_get_val_s(ivarfl(ir22), cvar_r22)
    call field_get_val_s(ivarfl(ir33), cvar_r33)
    call field_get_val_s(ivarfl(ir12), cvar_r12)
    call field_get_val_s(ivarfl(ir13), cvar_r13)
    call field_get_val_s(ivarfl(ir23), cvar_r23)
    call field_get_val_s(ivarfl(iep), cvar_ep)
  elseif (iturb.eq.50) then
    call field_get_val_s(ivarfl(ik), cvar_k)
    call field_get_val_s(ivarfl(iep), cvar_ep)
    call field_get_val_s(ivarfl(iphi), cvar_phi)
    call field_get_val_s(ivarfl(ifb), cvar_fb)
  elseif (iturb.eq.60) then
    call field_get_val_s(ivarfl(ik), cvar_k)
    call field_get_val_s(ivarfl(iomg), cvar_omg)
  elseif (iturb.eq.70) then
    call field_get_val_s(ivarfl(inusa), cvar_nusa)
  endif

! --- Velocity components
  do iel = 1, ncel
    cvar_vel(1,iel) = 0.0d0
    cvar_vel(2,iel) = 0.00d0
    cvar_vel(3,iel) = 0.0d0
  enddo

! --- 1. Pressure (Pa)
  if(.true.) then
!    ithvar = ithvar*2
    do iel = 1, ncel
      cvar_pr(iel) = p0
    enddo
  endif

! --- 2. Turbulence Components
 do iel = 1, ncel

    zent = xyzcen(3,iel)

!    call intprf                                                   &
    !==========
!   (nbmetd, nbmetm,                                               &
!    zdmet, tmmet, umet , zent  , ttcabs, xuent )

!    call intprf                                                   &
    !==========
!   (nbmetd, nbmetm,                                               &
!    zdmet, tmmet, vmet , zent  , ttcabs, xvent )

!    call intprf                                                   &
    !==========
!   (nbmetd, nbmetm,                                               &
!    zdmet, tmmet, ekmet, zent  , ttcabs, xkent )

!    call intprf                                                   &
    !==========
!   (nbmetd, nbmetm,                                               &
!    zdmet, tmmet, epmet, zent  , ttcabs, xeent )
! --- Velocity components
!    cvar_vel(1,iel) = xuent
!    cvar_vel(2,iel) = xvent
!    cvar_vel(3,iel) = 0.d0

!     ITYTUR est un indicateur qui vaut ITURB/10
    if    (itytur.eq.2) then

!      cvar_k(iel)  = xkent
!      cvar_ep(iel) = xeent

!    elseif (itytur.eq.3) then

!      cvar_r11(iel) = d2s3*xkent
!      cvar_r22(iel) = d2s3*xkent
!      cvar_r33(iel) = d2s3*xkent
!      cvar_r12(iel) = 0.d0
!      cvar_r13(iel) = 0.d0
!      cvar_r23(iel) = 0.d0
!      cvar_ep(iel)  = xeent

!    elseif (iturb.eq.50) then

!      cvar_k(iel)   = xkent
!      cvar_ep(iel)  = xeent
!      cvar_phi(iel) = d2s3
!      cvar_fb(iel)  = 0.d0

    elseif (iturb.eq.60) then

      cvar_k(iel)   = 0.75
      cvar_omg(iel) = 0.03

!    elseif (iturb.eq.70) then

!      cvar_nusa(iel) = cmu*xkent**2/xeent

    endif

    if (iscalt.ge.0) then
! On suppose que le scalaire est la temperature potentielle :
      call intprf                                                 &
      !==========
   (nbmett, nbmetm,                                               &
    ztmet, tmmet, tpmet, zent  , ttcabs, tpent )

      call field_get_val_s(ivarfl(isca(iscalt)), cvar_scalt)
      cvar_scalt(iel) = tpent

    endif
  enddo

endif
!< [init]

!--------
! Formats
!--------

!----
! End
!----

deallocate(lstelt) ! temporary array for cells selection

return
end subroutine cs_user_initialization
