
/*
  This file is part of code_saturne, a general-purpose CFD tool.

  Copyright (C) 1998-2023 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.
*/

/*----------------------------------------------------------------------------*/

#include "cs_defs.h"

/*----------------------------------------------------------------------------
 * Standard C library headers
 *----------------------------------------------------------------------------*/

#include <assert.h>
#include <math.h>
#include <string.h>

#if defined(HAVE_MPI)
#include <mpi.h>
#endif

/*----------------------------------------------------------------------------
 * PLE library headers
 *----------------------------------------------------------------------------*/

#include <ple_coupling.h>

/*----------------------------------------------------------------------------
 * Local headers
 *----------------------------------------------------------------------------*/

#include "cs_headers.h"

/*----------------------------------------------------------------------------*/

BEGIN_C_DECLS

/*----------------------------------------------------------------------------*/
/*!
 * \file cs_user_les_inflow-base.c
 *
 * \brief Generation of synthetic turbulence at LES inlets initialization.
 *
 * See \ref les_inflow for examples.
 */
/*----------------------------------------------------------------------------*/

/*============================================================================
 * Local (user defined) function definitions
 *============================================================================*/

/*----------------------------------------------------------------------------*/
/*!
 * \brief Define parameters of synthetic turbulence at LES inflow.
 */
/*----------------------------------------------------------------------------*/

void
cs_user_les_inflow_define(void)
{
  /*! [set_restart] */
  cs_les_inflow_set_restart(false,   /* allow_read */
                            false);  /* allow_write */
  /*! [set_restart] */

  /* First synthetic turbulence inlet: the Batten Method is used
     for boundary faces of zone "INLET_1" */

  /*! [init_1] */
  {
    /* Number of turbulence structures */
    int n_entities = 50;

    /* Velocity, turbulent kinetic energy and dissipation rate */
    cs_real_t vel_r[3] = {0.13, 0, 0};
    cs_real_t k_r =1;
    cs_real_t eps_r = 1;

    cs_les_inflow_add_inlet(CS_INFLOW_BATTEN,
                            false,
                            cs_boundary_zone_by_name("INLET_1"),
                            n_entities,
                            0,  /* verbosity */
                            vel_r,
                            k_r,
                            eps_r);
  }
}
/*----------------------------------------------------------------------------*/
/*!
 * \brief Update of the characteristics of a given synthetic turbulence inlet.
 *
 * \param[in]   zone       pointer to associated boundary zone
 * \param[out]  vel_r      reference mean velocity
 * \param[out]  k_r        reference turbulent kinetic energy
 * \param[out]  eps_r      reference turbulent dissipation
 */
/*----------------------------------------------------------------------------*/




END_C_DECLS
