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

!                      Code_Saturne version 6.0.5
!                      --------------------------
! This file is part of Code_Saturne, a general-purpose CFD tool.
!
! Copyright (C) 1998-2019 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.

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

!===============================================================================
! Function:
! ---------

!> \file cs_user_boundary_conditions-mapped_inlet.f90
!>
!> Example of cs_f_user_boundary_conditions subroutine.f90 for inlet
!> with inlet profile mapped to profile inside the domain.

!> This example assumes the mesh is orthogonal at the inlet.
!
!-------------------------------------------------------------------------------

!-------------------------------------------------------------------------------
! Arguments
!______________________________________________________________________________.
!  mode           name          role                                           !
!______________________________________________________________________________!
!> \param[in]     nvar          total number of variables
!> \param[in]     nscal         total number of scalars
!> \param[out]    icodcl        boundary condition code:
!>                               - 1 Dirichlet
!>                               - 2 Radiative outlet
!>                               - 3 Neumann
!>                               - 4 sliding and
!>                                 \f$ \vect{u} \cdot \vect{n} = 0 \f$
!>                               - 5 smooth wall and
!>                                 \f$ \vect{u} \cdot \vect{n} = 0 \f$
!>                               - 6 rough wall and
!>                                 \f$ \vect{u} \cdot \vect{n} = 0 \f$
!>                               - 9 free inlet/outlet
!>                                 (input mass flux blocked to 0)
!>                               - 13 Dirichlet for the advection operator and
!>                                    Neumann for the diffusion operator
!> \param[in]     itrifb        indirection for boundary faces ordering
!> \param[in,out] itypfb        boundary face types
!> \param[out]    izfppp        boundary face zone number
!> \param[in]     dt            time step (per cell)
!> \param[in,out] rcodcl        boundary condition values:
!>                               - rcodcl(1) value of the dirichlet
!>                               - rcodcl(2) value of the exterior exchange
!>                                 coefficient (infinite if no exchange)
!>                               - rcodcl(3) value flux density
!>                                 (negative if gain) in w/m2 or roughness
!>                                 in m if icodcl=6
!>                                 -# for the velocity \f$ (\mu+\mu_T)
!>                                    \gradt \, \vect{u} \cdot \vect{n}  \f$
!>                                 -# for the pressure \f$ \Delta t
!>                                    \grad P \cdot \vect{n}  \f$
!>                                 -# for a scalar \f$ cp \left( K +
!>                                     \dfrac{K_T}{\sigma_T} \right)
!>                                     \grad T \cdot \vect{n} \f$
!_______________________________________________________________________________

subroutine cs_f_user_boundary_conditions &
 ( nvar   , nscal  ,                                              &
   icodcl , itrifb , itypfb , izfppp ,                            &
   dt     ,                                                       &
   rcodcl )

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

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

use paramx
use numvar
use optcal
use cstphy
use cstnum
use entsor
use parall
use period
use ihmpre
use ppppar
use ppthch
use coincl
use cpincl
use ppincl
use ppcpfu
use atincl
use atsoil
use ctincl
use cs_fuel_incl
use mesh
use field
use turbomachinery
use iso_c_binding
use cs_c_bindings

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

implicit none

! Arguments

integer          nvar   , nscal

integer          icodcl(nfabor,nvar)
integer          itrifb(nfabor), itypfb(nfabor)
integer          izfppp(nfabor)

double precision dt(ncelet)
double precision rcodcl(nfabor,nvar,3)

! Local variables

!< [loc_var_dec]
integer          ifac, iel, ii, ivar, iscal, ilelt, nlfac

integer          keyvar, keysca
integer          n_fields, f_id, normalize, interpolate
double precision xdh, rhomoy
double precision fmprsc, uref2

integer, allocatable, dimension(:) :: lstfac, lstcel
double precision, dimension(:), pointer :: brom

double precision, dimension(:,:), pointer :: vel
double precision, dimension(3,1) :: coord_shift
double precision, dimension(1) :: rvoid

type(c_ptr), save :: inlet_l1 = c_null_ptr
type(c_ptr), save :: inlet_l2 = c_null_ptr
type(c_ptr), save :: inlet_l3 = c_null_ptr
type(c_ptr), save :: inlet_l4 = c_null_ptr
!< [loc_var_dec]

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

!< [init]
allocate(lstfac(nfabor))  ! temporary array for boundary faces selection

! Map field arrays
call field_get_val_v(ivarfl(iu), vel)

call field_get_val_s(ibrom, brom)

fmprsc = 1.d0 ! mean prescribed velocity
!< [init]

!===============================================================================
! Assign a pseudo-periodic channel type inlet to a set of boundary faces,
! using boundary condition mapping.

! For each subset:
! - use selection criteria to filter boundary faces of a given subset
! - use boundary_conditions_map_set and boundary_conditions_mapped_set
!   to apply a profile from inside the domain to the inlet, renormalizing
!   for some variables.
!
! The impled feedback loop allows progressively reaching a state similar
! to that of a periodic channel at the inlet.
!===============================================================================

call getfbr('regular_fan_mapped_inlet_minus012', nlfac, lstfac)

if (ntcabs == ntpabs) then
  allocate(lstcel(ncel))
  do iel = 1, ncel
    lstcel(iel) = iel
  enddo
  coord_shift(1,1) = 0.d0
  coord_shift(2,1) = -0.12
  coord_shift(3,1) = 0.d0
  inlet_l1 = boundary_conditions_map(mesh_location_cells, ncel,                &
                                    nlfac, lstcel, lstfac,                    &
                                    coord_shift, 0,                           &
                                    0.01d0)
  deallocate(lstcel)
endif

  call field_get_n_fields(n_fields)
  call field_get_key_id("velocity", keyvar)
  call field_get_key_id("scalar1", keysca)

  interpolate = 0

 do f_id = 0, n_fields-1
    normalize = 0
    call boundary_conditions_mapped_set(f_id, inlet_l1, MESH_LOCATION_CELLS, &
                                        normalize, interpolate,             &
                                        nlfac, lstfac, rvoid,               &
                                        nvar, rcodcl)
  enddo

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

call getfbr('regular_fan_mapped_inlet_plus012', nlfac, lstfac)

if (ntcabs == ntpabs) then
  allocate(lstcel(ncel))
  do iel = 1, ncel
    lstcel(iel) = iel
  enddo
  coord_shift(1,1) = 0.d0
  coord_shift(2,1) = 0.12
  coord_shift(3,1) = 0.d0
  inlet_l2 = boundary_conditions_map(mesh_location_cells, ncel,                &
                                    nlfac, lstcel, lstfac,                    &
                                    coord_shift, 0,                           &
                                    0.01d0)
  deallocate(lstcel)
endif

  call field_get_n_fields(n_fields)
  call field_get_key_id("velocity", keyvar)
  call field_get_key_id("scalar1", keysca)

  interpolate = 0

 do f_id = 0, n_fields-1
    normalize = 0
    call boundary_conditions_mapped_set(f_id, inlet_l2, MESH_LOCATION_CELLS, &
                                        normalize, interpolate,             &
                                        nlfac, lstfac, rvoid,               &
                                        nvar, rcodcl)
  enddo

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

call getfbr('reduced_fan_mapped_inlet_minus012', nlfac, lstfac)

if (ntcabs == ntpabs) then
  allocate(lstcel(ncel))
  do iel = 1, ncel
    lstcel(iel) = iel
  enddo
  coord_shift(1,1) = 0.d0
  coord_shift(2,1) = -0.12
  coord_shift(3,1) = 0.d0
  inlet_l3 = boundary_conditions_map(mesh_location_cells, ncel,                &
                                    nlfac, lstcel, lstfac,                    &
                                    coord_shift, 0,                           &
                                    0.01d0)
  deallocate(lstcel)
endif

  call field_get_n_fields(n_fields)
  call field_get_key_id("velocity", keyvar)
  call field_get_key_id("scalar1", keysca)

  interpolate = 0

 do f_id = 0, n_fields-1
    normalize = 0
    call boundary_conditions_mapped_set(f_id, inlet_l3, MESH_LOCATION_CELLS, &
                                        normalize, interpolate,             &
                                        nlfac, lstfac, rvoid,               &
                                        nvar, rcodcl)
  enddo

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

call getfbr('extraction_mapped_inlet_minus001', nlfac, lstfac)

if (ntcabs == ntpabs) then
  allocate(lstcel(ncel))
  do iel = 1, ncel
    lstcel(iel) = iel
  enddo
  coord_shift(1,1) = 0.d0
  coord_shift(2,1) = 0.d0
  coord_shift(3,1) = -0.01
  inlet_l4 = boundary_conditions_map(mesh_location_cells, ncel,                &
                                    nlfac, lstcel, lstfac,                    &
                                    coord_shift, 0,                           &
                                    0.01d0)
  deallocate(lstcel)
endif

  call field_get_n_fields(n_fields)
  call field_get_key_id("velocity", keyvar)
  call field_get_key_id("scalar1", keysca)

  interpolate = 0

 do f_id = 0, n_fields-1
    normalize = 0
    call boundary_conditions_mapped_set(f_id, inlet_l4, MESH_LOCATION_CELLS, &
                                        normalize, interpolate,             &
                                        nlfac, lstfac, rvoid,               &
                                        nvar, rcodcl)
  enddo

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

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

deallocate(lstfac)  ! temporary array for boundary faces selection

return
end subroutine cs_f_user_boundary_conditions
