1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135
|
/*============================================================================
* Code_Saturne documentation page
*============================================================================*/
/*
This file is part of Code_Saturne, a general-purpose CFD tool.
Copyright (C) 1998-2018 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.
*/
/*-----------------------------------------------------------------------------*/
/*!
\page cs_head_losses Examples of data settings for head losses
(cs_user_zones.c and cs_user_head_losses.c)
\brief the \ref cs_user_head_losses function is used to compute
the values of the head loss term, and is called at each time step
for each previously defined head loss volume zone.
Volume zones may be defined using the GUI, or through
the \ref cs_user_zones (in cs_user_zones.c).
cku is the local head loss term.
It appears on the momentum as follows:
\f[ \rho \der{\vect{u}}{t} = - \grad p + \vect{headloss} \: (+\: \text{other terms})\f]
with \f[ \vect{headloss} = - \rho \tens{cku}\cdot \vect{u} \,\: (\text{in } kg\cdot m^{-2} \cdot s^{-1})\f]
For a distributed head loss, let \f${ \tens{\xi_l} = \dfrac{\tens{dh_l}}{(0.5 \rho u^2)}}\f$ given by the litterature
(\f$ \tens{dh_l} \f$ is the head loss per unit length)
the source term \c tspdc is equal to \f$\tens{dh_l} = - \tens{\xi_l}(0.5\rho\vect{u}^2)\f$
we have \f$ \tens{cku} = 0.5\tens{\xi_l}|\vect{U}| \f$
For a singular head loss, let \f$\tens{\xi_l} = \dfrac{\tens{dh_s}}{0.5\rho\vect{u}^2}\f$ given by the litterature
(\f$\tens{dh_s} \f$ is the singular head loss)
the source term \c tspdc is equal to \f[\frac{\tens{dh_s}}{L} = - \frac{\tens{\xi_l}}{L} (0.5 \rho\vect{u}^2)\f]. We have \f[\tens{cku} = 0.5\frac{\tens{\xi_s}}{L}|\vect{u}|\f]
where \f$ L \f$ is the length over which we have chosen to represent the
singular head loss.
\section cs_user_head_losses_examples Head loss setting examples
Here is the list of examples:
- \subpage base_head_losses_examples
*/
// _____________________________________________________________________________
/*!
\page base_head_losses_examples Basic examples
\section init_and_final Initialization and finalization
It is useful to map a field array to a local pointer for a clear and concise
access, such as done here fro the velocity:
\snippet cs_user_head_losses.c map_field_arrays
Otherwise, the zone entries (see \ref cs_volume_zone_t) should contain
the necessary information with no additional preparation.
\section body Body
\subsection beginning Defining a volume zone
A volume zone may be defined using the GUI, or in the \ref cs_user_zones
user function (in cs_user_zones.c), such as the following zone determined
by a geometric criterion:
\snippet cs_user_zones.c user_zones_head_loss_1
Note that if the \ref CS_VOLUME_ZONE_HEAD_LOSS flag is not set
(or the matching type set through the GUI), calls to \ref cs_user_head_losses
will ignore this zone.
\subsection head_loss_examples Head loss examples
Note that in the following examples, we checku the zone name, so we
know which zone we are dealing with using in case of multiple zones.
head loss tensor coefficients for each cell are organized as follows:
cku11, cku22, cku33, cku12, cku13, cku23.
Coefficients are set to zero (then computed based on definitions provided
through the GUI if this is the case) before calling this function, so
setting values to zero is usually not necessary, unless we want to fully
overwrite a GUI-based definition.
Note that diagonal coefficients must be positive; the calculation may
crash if this is not the case.
\subsection diagonal_tensor Example 1: head losses in direction \c x
Using the previously defined zone, we define head losses in direction \c x
\snippet cs_user_head_losses.c head_loss_1
\subsection alpha_tensor Example 2: alpha = 45 degres
3x3 tensor: Example of head losses at alpha = 45 degres x,y
direction \c x resists by \c cku1 and \c y by \c cku2 \n
<tt> cku2 = 0 </tt> represents vanes as follows:
in coordinate system \c x, \c y
\image html orthogonal_reference_frame_sketch.gif "Orthogonal reference frame sketch"
\snippet cs_user_head_losses.c head_loss_2
*/
// _____________________________________________________________________________
|