The DEFINE
macros presented in this
section are used for multiphase applications, only.
Table 2.10: Quick Reference Guide for Multiphase DEFINE Macros provides
a quick reference guide to the multiphase-specific DEFINE
macros, the functions they are used to define, and the dialog boxes
in which they are activated in Ansys Fluent. Definitions of each DEFINE
macro are listed in the udf.h
header file (see Appendix C: Quick Reference Guide for Multiphase DEFINE
Macros).
Appendix B: DEFINE
Macro Definitions contains a list of general
purpose DEFINE
macros that can also be used
to define UDFs for multiphase cases. For example, the general purpose DEFINE_PROPERTY
macro is used to define a surface tension
coefficient UDF for the multiphase VOF model. See
DEFINE_PROPERTY
UDFs for details.
Table 2.10: Quick Reference Guide for Multiphase DEFINE Macros
Model | Function |
DEFINE Macro | Dialog Box Activated |
---|---|---|---|
VOF | linearized mass transfer |
DEFINE_LINEARIZED_MASS_TRANSFER
| Phase Interaction |
mass transfer |
DEFINE_MASS_TRANSFER
| Phase Interaction | |
heterogeneous reaction rate |
DEFINE_HET_RXN_RATE
| Phase Interaction | |
Mixture | linearized mass transfer |
DEFINE_LINEARIZED_MASS_TRANSFER
| Phase Interaction |
mass transfer |
DEFINE_MASS_TRANSFER
| Phase Interaction | |
drag exchange coefficient |
DEFINE_EXCHANGE_PROPERTY
| Phase Interaction | |
drag modification factor |
DEFINE_EXCHANGE_PROPERTY
| Phase Interaction | |
slip velocity |
DEFINE_VECTOR_EXCHANGE_PROPERTY
| Phase Interaction | |
cavitation rate |
DEFINE_CAVITATION_RATE
| User-Defined Function Hooks | |
heterogeneous reaction rate |
DEFINE_HET_RXN_RATE
| Phase Interaction | |
Eulerian | linearized mass transfer |
DEFINE_LINEARIZED_MASS_TRANSFER
| Phase Interaction |
mass transfer |
DEFINE_MASS_TRANSFER
| Phase Interaction | |
heat transfer |
DEFINE_EXCHANGE_PROPERTY
| Phase Interaction | |
drag exchange coefficient |
DEFINE_EXCHANGE_PROPERTY
| Phase Interaction | |
drag modification factor |
DEFINE_EXCHANGE_PROPERTY
| Phase Interaction | |
lift coefficient |
DEFINE_EXCHANGE_PROPERTY
| Phase Interaction | |
wall lubrication coefficient |
DEFINE_EXCHANGE_PROPERTY
| Phase Interaction | |
wall lubrication model parameters |
DEFINE_EXCHANGE_PROPERTY
| Phase Interaction | |
turbulent dispersion force |
DEFINE_VECTOR_EXCHANGE_PROPERTY
| Phase Interaction | |
turbulent dispersion model parameters |
DEFINE_EXCHANGE_PROPERTY
| Phase Interaction | |
turbulence interaction model parameters |
DEFINE_EXCHANGE_PROPERTY
| Phase Interaction | |
virtual mass coefficient |
DEFINE_EXCHANGE_PROPERTY
| Phase Interaction | |
interfacial area |
DEFINE_EXCHANGE_PROPERTY
| Phase Interaction | |
boiling model and quenching |
DEFINE_BOILING_PROPERTY
| Boiling Model | |
heterogeneous reaction rate |
DEFINE_HET_RXN_RATE
| Phase Interaction |
You can use DEFINE_BOILING_PROPERTY
to model the boiling model parameters and the quenching model correction.
The parameters include the Bubble Departure Diameter, Frequency of Bubble Departure, Nucleation
Site Density, Area Influence Coeff., and Liquid Reference Temperature for quenching
correction.
DEFINE_BOILING_PROPERTY
(name
, f
, t
, c0
, t0
, from_phase_index
, from_species_index
, to_phase_index
, to_species_index
)
Important: As with all the other user-defined functions, all of the arguments
to this DEFINE
macro must be placed on the
same line in your source code. Splitting the DEFINE
statement onto several lines will result in a compilation error.
Argument Type |
Description |
---|---|
|
UDF name. |
|
Index that identifies a wall face. |
|
Pointer to the wall face thread. |
|
Cell index that identifies the cell next to the wall. |
|
Pointer to mixture-level cell thread. |
|
Liquid phase in boiling models. |
|
ID of liquid species (ID=-1 if no mixture material). |
|
Vapor phase in boiling models. |
|
ID of vapor species (ID=-1 if no mixture material). |
Function returns
real
There are nine arguments to DEFINE_BOILING_PROPERTY
: name
, f
, t
, c0
, t0
, from_phase_index
, from_species_index
, to_phase_index
, and to_species_index
. You supply name
, the name of the UDF.
The remaining eight variables are passed by the Ansys Fluent solver to
your UDF. The defined UDF will return the desired real
value for a specific model parameter.
Important: Note that the arguments from_species_index
and to_species_index
are relevant for multiphase
species transport problems only, and only if the respective phase
has a mixture material associated with it.
The following UDF named bubble_depart_dia
, demonstrates how the bubble diameter is computed. All other boiling
parameters can use this example and can be modified accordingly.
/*********************************************************************** UDF that demonstrates how to compute the bubble diameter based on tolubinski-kostanchuk. Can be interpreted or compiled. ************************************************************************/ #include "udf.h" #define d_bw_max 0.0014 #define d_bw_coef 0.0006 #define subcool_ref 45.0 DEFINE_BOILING_PROPERTY(bubble_depart_dia,f,t,c0,t0,from_index,from_species_index,to_index,to_species_index) { real diam_b, subcool; int liq_phase = from_index; Thread **pt0 = THREAD_SUB_THREADS(t0); real T_SAT = C_STORAGE_R(c0,t0,SV_SAT_TEMPERATURE); real T_l = C_T(c0, pt0[liq_phase]); subcool = T_SAT - T_l; diam_b = MIN(d_bw_max,d_bw_coef*exp(-subcool/subcool_ref)); return diam_b; }
After the UDF that you have defined using DEFINE_BOILING_PROPERTY
is interpreted (Interpreting UDFs)
or compiled (Compiling UDFs), the name
of the argument that you supplied as the first DEFINE
macro argument (for example, bubble_depart_dia
) will become visible and selectable in the Boiling Models dialog box in Ansys Fluent. See Hooking DEFINE_BOILING_PROPERTY
UDFs for details.
You can use DEFINE_CAVITATION_RATE
to
model the cavitation source terms and in the vapor
mass fraction transport equation used in the Singhal et al model (see Equation 14–603 in the Theory Guide). Assuming denotes the mass-transfer rate between liquid and
vapor phases, we have
where and are the mass-fraction of the liquid and vapor phase, respectively.
DEFINE_CAVITATION_RATE
is used to calculate only. The values of and are computed
by the solver, accordingly.
Note: For the DEFINE_CAVITATION_RATE
UDF, the vaporization
pressure is set to 1.15 times the saturation pressure. To change this value, enter the
following scheme command in the console:
(rpsetvar 'mp/cvt/vapor-pressure-factor
)user-input
where
is your specified vaporization pressure value.user-input
DEFINE_CAVITATION_RATE
(name
, c
, t
, p
, rhoV
, rhoL
, mafV
, p_v
, cigma
, f_gas
, m_dot
)
Argument Type |
Description |
---|---|
|
UDF name. |
|
Cell index. |
|
Pointer to the mixture-level thread. |
|
Pointer to shared pressure. |
|
Pointer to vapor density. |
|
Pointer to liquid density. |
|
Pointer to vapor mass fraction. |
|
Pointer to vaporization pressure. |
|
Pointer to liquid surface tension coefficient. |
|
Pointer to the prescribed mass fraction of non condensable gases. |
|
Pointer to cavitation mass transfer rate. |
Function returns
void
There are eleven arguments to DEFINE_CAVITATION_RATE
: name
, c
, t
, p
, rhoV
, rhoL
, mafV
, p_v
, cigma
, f_gas
, and m_dot
. You supply name
, the name of the UDF. c
, t
, p
, rhoV
, rhoL
, mafV
, p_v
, cigma
, f_gas
,
and m_dot
are variables that are passed by
the Ansys Fluent solver to your UDF. Your UDF will need to set the value
referenced by the real
pointer m_dot
to the cavitation rate.
The following UDF named c_rate
, is
an example of a cavitation model for a multiphase mixture that is
different from the default model in Ansys Fluent. This cavitation model
calculates the cavitation mass transfer rates between the liquid and
vapor phase depending on fluid pressure (*p
), turbulence kinetic energy (C_K(c,t)
),
and the liquid vaporization pressure (*p_v
).
In general, the existence of turbulence enhances cavitation.
In this example, the turbulence effect is taken into account by increasing
the cavitation pressure by 0.195* C_R(c,t) * C_K(c,t)
. The pressure p_vapor
that determines whether
cavitation occurs increases from p_v
to
p_v + 0.195 * C_R(c,t) * C_K(c,t)
When the absolute fluid pressure (ABS_P
) is lower than p_vapor
, then liquid evaporates
to vapor (). When it is greater than p_vapor
, vapor condenses to liquid ().
The evaporation rate is calculated by
If ABS_P < p_vapor, then c_evap * rhoV[c] * sqrt(2.0/3.0*rhoL[c]) * ABS(p_vapor - ABS_P(p[c]))
The condensation rate is
If ABS_P > p_vapor, then -c_con*rhoL[c] * sqrt(2.0/3.0*rhoL[c]) * ABS(p_vapor - ABS_P(p[c]))
where c_evap
and c_con
are model coefficients.
/*********************************************************************** UDF that is an example of a cavitation model different from default. Can be interpreted or compiled. ************************************************************************/ #include "udf.h" #define c_evap 1.0 #define c_con 0.1 DEFINE_CAVITATION_RATE(c_rate,c,t,p,rhoV,rhoL,mafV,p_v,cigma,f_gas, m_dot) { real p_vapor = *p_v; real dp, dp0, source; p_vapor += MIN(0.195*C_R(c,t)*C_K(c,t), 5.0*p_vapor); dp = p_vapor - ABS_P(p[c], op_pres); dp0 = MAX(0.1, ABS(dp)); source = sqrt(2.0/3.0*rhoL[c])*dp0; if(dp > 0.0) *m_dot = c_evap*rhoV[c]*source; else *m_dot = -c_con*rhoL[c]*source; }
After the UDF that you have defined using DEFINE_CAVITATION_RATE
is interpreted (Interpreting UDFs)
or compiled (Compiling UDFs), the name
of the argument that you supplied as the first DEFINE
macro argument (for example, c_rate
) will
become visible and selectable in the User-Defined Function
Hooks dialog box in Ansys Fluent. See Hooking DEFINE_CAVITATION_RATE
UDFs for details.
You can use DEFINE_EXCHANGE_PROPERTY
to specify UDFs for some
phase interaction variables in multiphase models. These include net heat transfer rate
between phases, virtual mass coefficient, drag exchange coefficient, lift coefficient,
wall lubrication coefficient, and interfacial area (for the Eulerian multiphase boiling
model). Drag exchange coefficient may also be specified for the Mixture model. Below is a
list of user-defined functions that can be specified using
DEFINE_EXCHANGE_PROPERTY
for the multiphase models in Ansys Fluent.
Note that there are some phase interaction variables such as vaporization pressure and
surface tension coefficient (cavitation parameters) that are defined using
DEFINE_PROPERTY
. See
DEFINE_PROPERTY
UDFs for details.
Table 2.11:
DEFINE_EXCHANGE_PROPERTY
Variables
Mixture Model |
Eulerian Model |
---|---|
drag exchange coefficient drag modification factor, in Equation 14–293 of the Theory Guide |
volumetric heat transfer coefficient, in Equation 14–382 of the Theory Guide virtual mass coefficient, in Equation 14–343 in the Fluent Theory Guide drag exchange coefficient, in Equation 14–212 of the Theory Guide drag modification factor, in Equation 14–293 of the Theory Guide lift coefficient, in Equation 14–298 of the Theory Guide wall lubrication coefficient, in Equation 14–319 of the Theory Guide interfacial area, in Interfacial Area Concentration in the Theory Guide model coefficients for wall lubrication, turbulent dispersion, and turbulence interaction models. |
DEFINE_EXCHANGE_PROPERTY
(name
,
c
, mixture_thread
,
phase_index_1
, phase_index_2
)
Important: Note that all of the arguments to a DEFINE
macro must be
placed on the same line in your source code. Splitting the
DEFINE
statement onto several lines will result in a
compilation error.
Argument Type |
Description |
---|---|
|
UDF name. |
|
Cell index. |
|
Pointer to the mixture-level thread. |
|
Identifier of the first phase from the pair of phases for which you want to
customize a phase interaction property. An index of |
|
Identifier of the second phase from the pair of phases for which you want to
customize a phase interaction property. See |
Function returns
real
There are five arguments to DEFINE_EXCHANGE_PROPERTY
:
name
, c
,
mixture_thread
, phase_index_1
, and
phase_index_2
. You supply name
, the
name of the UDF. c
, mixture_thread
,
phase_index_1
, and phase_index_2
are
variables that are passed by the Ansys Fluent solver to your UDF. Your UDF will need to return
the real
value of the volumetric heat transfer coefficient,
virtual mass coefficient, lift coefficient, drag exchange coefficient, wall lubrication
coefficient, or interfacial area to the solver.
Note:
ANSYS Fluent passes values of phase identifiers in the same order as they appear in the Phase Pairs selection list in the Multiphase Model dialog box (Phase Interaction > Forces tab). In the example shown in Defining Forces in the Forces Tab in the Fluent User's Guide, if you want to customize a phase interaction property for the phase pair (liquid vapor), where liquid is the primary phase and vapor is the only secondary phase, ANSYS Fluent will pass
0
and1
forphase_index_1
andphase_index_2
, respectively.Your UDF should return the effective coefficient. So, if the quantity being modeled by your UDF depends on the interfacial area, your UDF should return the coefficient multiplied by the interfacial area.
The following UDF, named custom_drag
, can be used to customize the default Syamlal drag law in Ansys Fluent.
The default drag law uses 0.8 (for void <0.85) and 2.65 (void>0.85)
for bfac
. This results in a minimum fluid
velocity of 25 cm/s. The UDF modifies the drag law to result in a
minimum fluid velocity of 8 cm/s, using 0.28 and 9.07 for the bfac
parameters.
/****************************************************************** UDF for customizing the default Syamlal drag law in ANSYS Fluent *******************************************************************/ #include "udf.h" #define pi 4.*atan(1.) #define diam2 3.e-4 DEFINE_EXCHANGE_PROPERTY(custom_drag,cell,mix_thread,s_col,f_col) { Thread *thread_g, *thread_s; real x_vel_g, x_vel_s, y_vel_g, y_vel_s, abs_v, slip_x, slip_y, rho_g, rho_s, mu_g, reyp, afac, bfac, void_g, vfac, fdrgs, taup, k_g_s; /* find the threads for the gas (primary) */ /* and solids (secondary phases) */ thread_g = THREAD_SUB_THREAD(mix_thread, s_col);/* gas phase */ thread_s = THREAD_SUB_THREAD(mix_thread, f_col);/* solid phase*/ /* find phase velocities and properties*/ x_vel_g = C_U(cell, thread_g); y_vel_g = C_V(cell, thread_g); x_vel_s = C_U(cell, thread_s); y_vel_s = C_V(cell, thread_s); slip_x = x_vel_g - x_vel_s; slip_y = y_vel_g - y_vel_s; rho_g = C_R(cell, thread_g); rho_s = C_R(cell, thread_s); mu_g = C_MU_L(cell, thread_g); /*compute slip*/ abs_v = sqrt(slip_x*slip_x + slip_y*slip_y); /*compute Reynolds number*/ reyp = rho_g*abs_v*diam2/mu_g; /* compute particle relaxation time */ taup = rho_s*diam2*diam2/18./mu_g; void_g = C_VOF(cell, thread_g);/* gas vol frac*/ /*compute drag and return drag coeff, k_g_s*/ afac = pow(void_g,4.14); if(void_g <= 0.85) bfac = 0.281632*pow(void_g, 1.28); else bfac = pow(void_g, 9.076960); vfac = 0.5*(afac-0.06*reyp+sqrt(0.0036*reyp*reyp+0.12*reyp*(2.*bfac- afac)+afac*afac)); fdrgs = void_g*(pow((0.63*sqrt(reyp)/ vfac+4.8*sqrt(vfac)/vfac),2))/24.0; k_g_s = (1.-void_g)*rho_s*fdrgs/taup; return k_g_s; }
The following UDF, named custom_lift, computes the coefficient for the lift force using a formulation developed by Tomiyama in 2002:
/* this example uses "user-defined" to implement a lift model by Tomiyama et al in 2002 */ #include "udf.h" #include "flow.h" #define Eo_l1 4.0 #define Eo_l2 10.0 DEFINE_EXCHANGE_PROPERTY(custom_lift,c,t,i,j) { /* i -- liquid-phase; j -- vapor-phase */ Thread **pt = THREAD_SUB_THREADS(t); real v_x=0., v_y=0., v_z=0.; real vel, Rev, Eo, d2, T_sfc, sigma; real lift_coeff, lift_co, wk_co;; real diam = C_PHASE_DIAMETER(c,pt[j]); real rho_v = C_R(c,pt[j]); real rho_l = C_R(c,pt[i]); real mu_l = C_MU_L(c,pt[i]); real gravity = NV_MAG(M_gravity); Property *prop = DOMAIN_COMPLEX_PROP_PROPERTY(DOMAIN_INTERACTION(root_domain), COMPLEX_PROP_sfc_tension_coeff, i,j); T_sfc = (sg_temperature && NNULLP(THREAD_STORAGE(pt[i],SV_T)))? C_T(c,pt[i]) : T_REF; if(prop == NULL || PROPERTY_METHOD(prop,0) == PROP_METHOD_NONE) Error("Lift-Tomiyama: Please set value for surface tension !"); sigma = generic_property(c,t,prop,(Property_ID)0,T_sfc); if(sigma <= 0.) Error("Lift-Tomiyama: Please set nonzero value for surface tension !"); /* calculate bubble Reynolds Number */ v_x = C_U(c,pt[j]) - C_U(c,pt[i]); v_y = C_V(c,pt[j]) - C_V(c,pt[i]); #if RP_3D v_z = C_W(c,pt[j]) - C_W(c,pt[i]); #endif vel = sqrt(v_x*v_x + v_y*v_y + v_z*v_z); Rev = RE_NUMBER(rho_l,vel,diam,mu_l); d2 = diam*diam; Eo = gravity*(rho_l-rho_v)*d2/sigma; if (Eo <= Eo_l1) wk_co = 0.0; else if (Eo < Eo_l2) wk_co = -0.096*Eo + 0.384; else wk_co = -0.576; lift_co = 0.288*tanh(0.121*Rev); lift_coeff = lift_co + wk_co; return lift_coeff; }
The following UDF, named heat_udf
,
specifies a coefficient that when multiplied by the temperature difference
between the dispersed and continuous phases, is equal to the net rate
of heat transfer per unit volume.
#include "udf.h" #define PR_NUMBER(cp,mu,k) ((cp)*(mu)/(k)) #define IP_HEAT_COEFF(vof,k,nu,d) ((vof)*6.*(k)*(Nu)/(d)/(d)) static real heat_ranz_marshall(cell_t c, Thread *ti, Thread *tj) { real h; real d = C_PHASE_DIAMETER(c,tj); real k = C_K_L(c,ti); real NV_VEC(v), vel, Re, Pr, Nu; NV_DD(v,=,C_U(c,tj),C_V(c,tj),C_W(c,tj),-,C_U(c,ti),C_V(c,ti),C_W(c,ti)); vel = NV_MAG(v); Re = RE_NUMBER(C_R(c,ti),vel,d,C_MU_L(c,ti)); Pr = PR_NUMBER (C_CP(c,ti),C_MU_L(c,ti),k); Nu = 2. + 0.6*sqrt(Re)*pow(Pr,1./3.); h = IP_HEAT_COEFF(C_VOF(c,tj),k,Nu,d); return h; } DEFINE_EXCHANGE_PROPERTY(heat_udf, c, t, i, j) { Thread *ti = THREAD_SUB_THREAD(t,i); Thread *tj = THREAD_SUB_THREAD(t,j); real val; val = heat_ranz_marshall(c,ti, tj); return val; }
The following UDF named custom_ia
,
computes the interfacial area, while including the symmetric model.
#include "udf.h" DEFINE_EXCHANGE_PROPERTY(custom_ia,c,t,i,j) { /* i -- liquid-phase; j -- vapor-phase */ Thread **pt = THREAD_SUB_THREADS(t); real diam = C_PHASE_DIAMETER(c, pt[j]); real vof_i = C_VOF(c,pt[i]); real vof_j = C_VOF(c,pt[j]); real area_intf; area_intf = 6.*vof_i*vof_j/diam; return area_intf; }
After the UDF that you have defined using DEFINE_EXCHANGE_PROPERTY
is
interpreted (Interpreting UDFs) or compiled (Compiling UDFs), the name of the argument that you supplied as
the first DEFINE
macro argument (for example,
heat_udf
) will become visible and selectable in the
Phase Interaction tab of the Multiphase Model
dialog box. See Hooking DEFINE_EXCHANGE_PROPERTY
UDFs for details.
You need to use DEFINE_HET_RXN_RATE
to specify reaction rates for heterogeneous reactions. A heterogeneous
reaction is one that involves reactants and products from more than
one phase. Unlike DEFINE_VR_RATE
, a DEFINE_HET_RXN_RATE
UDF can be specified differently
for different heterogeneous reactions.
During Ansys Fluent execution, the DEFINE_HET_RXN_RATE
UDF for each heterogeneous reaction that is defined is called in
every fluid cell. Ansys Fluent will use the reaction rate specified by
the UDF to compute production/destruction of the species participating
in the reaction, as well as heat and momentum transfer across phases
due to the reaction.
A heterogeneous reaction is typically used to define reactions involving species of different phases. Heterogeneous reactions are defined in the Phase Interaction > Heat, Mass, Reactions > Reactions tab of the Multiphase Model dialog box.
DEFINE_HET_RXN_RATE
(name,c,t,r,mw,yi,rr,rr_t
)
Argument Type |
Description |
---|---|
|
UDF name. |
|
Cell index. |
|
Cell thread (mixture level) on which heterogeneous reaction rate is to be applied. |
|
Pointer to data structure
that represents the current heterogeneous reaction (see |
|
Matrix of species molecular weights. |
|
Matrix of species mass fractions. |
|
Pointer to laminar reaction rate. |
|
Currently not used. Provided for future use. |
Function returns
void
There are eight arguments to DEFINE_HET_RXN_RATE
: name
, c
, t
, r
, mw
, yi
, rr
, and rr_t
. You supply name
, the
name of the UDF. c
, t
, r
, mw
, yi
, rr
, and rr_t
are variables that are passed by the Ansys Fluent solver to your UDF.
Your UDF will need to set the values referenced by the real
pointer rr
. The values
must be specified in (where
the volume is the cell volume).
The following compiled UDF named user_evap_condens_react
defines the reaction rate required to simulate evaporation or condensation
on the surface of droplets. Such a reaction can be formally described
by the following:
(2–14) |
Here, gas is a primary phase mixture of two species: and air. Droplets constitute the secondary phase and represent a mixture of one species - . Single-species mixtures are allowed in multiphase models.
The formulation for the reaction rate follows the model for particle evaporation that is defined in Droplet Vaporization (Law 2) in the Theory Guide.
#include "udf.h" /*Constants used in psat_h2o to calculate saturation pressure*/ #define PSAT_A 0.01 #define PSAT_TP 338.15 #define C_LOOP 8 #define H2O_PC 22.089E6 #define H2O_TC 647.286 /*user inputs*/ #define MAX_SPE_EQNS_PRIM 2 /*total number of species in primary phase*/ #define index_evap_primary 0 /*evaporating species index in primary phase*/ #define prim_index 0 /*index of primary phase*/ #define P_OPER 101325 /*operating pressure equal to GUI value*/ /*end of user inputs*/ /*************************************************************/ /* UDF for specifying an interfacial area density */ /*************************************************************/ double psat_h2o(double tsat) /* */ /* Computes saturation pressure of water vapor */ /* as function of temperature */ /* Equation is taken from THERMODYNAMIC PROPERTIES IN SI, */ /* by Reynolds, 1979 */ /* Returns pressure in PASCALS, given temperature in KELVIN */ { int i; double var1,sum1,ans1,psat; double constants[8]={-7.4192420, 2.97221E-1, -1.155286E-1, 8.68563E-3, 1.094098E-3, -4.39993E-3, 2.520658E-3, -5.218684E-4}; /* var1 is an expression that is used in the summation loop */ var1 = PSAT_A*(tsat-PSAT_TP); /* Compute summation loop */ i = 0; sum1 = 0.0; while (i < C_LOOP){ sum1+=constants[i]*pow(var1,i); ++i; } ans1 = sum1*(H2O_TC/tsat-1.0); /* compute exponential to determine result */ /* psat has units of Pascals */ psat = H2O_PC*exp(ans1); return psat; } DEFINE_HET_RXN_RATE(user_evap_condens_react, c, t, hr, mw, yi, rr, rr_t) { Thread **pt = THREAD_SUB_THREADS(t); Thread *tp = pt[0]; Thread *ts = pt[1]; int i; real concentration_evap_primary, accum = 0., mole_frac_evap_prim, concentration_sat ; real T_prim = C_T(c,tp); /*primary phase (gas) temperature*/ real T_sec = C_T(c,ts); /*secondary phase (droplet) temperature*/ real diam = C_PHASE_DIAMETER(c,ts); /*secondary phase diameter*/ real D_evap_prim = C_DIFF_EFF(c,tp,index_evap_primary) - 0.7*C_MU_T(c,tp)/C_R(c,tp); /*primary phase species turbulent diffusivity*/ real Re, Sc, Nu, urel, urelx,urely,urelz=0., mass_coeff, area_density, flux_evap ; if(Data_Valid_P()) { urelx = C_U(c,tp) - C_U(c,ts); urely = C_V(c,tp) - C_V(c,ts); #if RP_3D urelz = C_W(c,tp) - C_W(c,ts); #endif urel = sqrt(urelx*urelx + urely*urely + urelz*urelz); /*relative velocity*/ Re = urel * diam * C_R(c,tp) / C_MU_L(c,tp); Sc = C_MU_L(c,tp) / C_R(c,tp) / D_evap_prim ; Nu = 2. + 0.6 * pow(Re, 0.5)* pow(Sc, 0.333); mass_coeff = Nu * D_evap_prim / diam ; for (i=0; i < MAX_SPE_EQNS_PRIM ; i++) { accum = accum + C_YI(c,tp,i)/mw[i][prim_index]; } mole_frac_evap_prim = C_YI(c,tp,index_evap_primary) / mw[index_evap_primary][prim_index] / accum; concentration_evap_primary = mole_frac_evap_prim * P_OPER / UNIVERSAL_GAS_CONSTANT / T_prim ; concentration_sat = psat_h2o(T_sec)/UNIVERSAL_GAS_CONSTANT/T_sec ; area_density = 6. * C_VOF(c,ts) / diam ; flux_evap = mass_coeff * (concentration_sat - concentration_evap_primary) ; *rr = area_density * flux_evap ; } }
After the UDF that you have defined using DEFINE_HET_RXN_RATE
is
interpreted (Interpreting UDFs) or compiled (Compiling UDFs), the name of the argument that you supplied as
the first DEFINE
macro argument (for example,
user_evap_condens_react
) will become visible and selectable
under Reaction Rate Function in the Phase
Interaction > Heat, Mass, Reactions >
Reactions tab of the Multiphase Model dialog
box. (Note you will first need to specify the Total Number of
Reactions greater than 0.) See Hooking DEFINE_HET_RXN_RATE
UDFs for details.
You can use DEFINE_LINEARIZED_MASS_TRANSFER
when you want to model mass transfer in a multiphase problem. This
macro allows you to linearize the mass transfer source terms as well
as couple the interfacial mass transfer with flows. This is the recommend
UDF method for modeling mass transfer in multiphase flows.
You can linearize the mass transfer term for the calculation of the volume fraction equation in Ansys Fluent, such that
(2–15) |
where | |
= mass transfer rate | |
= volume fraction of the phase from which mass is transferred | |
= volume fraction of the phase to which mass is transferred | |
= mass transfer source which cannot be linearized to and | |
linearization coefficient related to | |
linearization coefficient related to |
To couple the mass transfer terms with flow transport equations, the derivative of the mass transfer rate to pressure is required to be computed and stored in a macro:
(2–16) |
Where and are the reference densities of the phases. Typically, they are the cell phase densities.
DEFINE_LINEARIZED_MASS_TRANSFER
(name
, c
, mixture_thread
, from_phase_index
, from_species_index
, to_phase_index
, to_species_index
, lin_from
, lin_to
)
Important: Note that all of the arguments to a DEFINE
macro need to be placed on the same line in your source code. Splitting
the DEFINE
statement onto several lines will
result in a compilation error.
Argument Type |
Description |
---|---|
|
UDF name. |
|
Index of cell on the thread pointed
to by |
|
Pointer to mixture-level thread. |
| Index of phase from which mass is transferred. |
|
ID of species from which mass is transferred (ID= -1 if phase does not have mixture material). |
|
Index of phase to which mass is transferred. |
|
ID of species to which mass is transferred (ID= -1 if phase does not have mixture material). |
|
Linearization term for the origin phase |
|
Linearization term for the destination phase |
Function returns
real
There are nine arguments to DEFINE_LINEARIZED_MASS_TRANSFER
: name
, c
, mixture_thread
, from_phase_index
, from_species_index
, to_phase_index
, to_species_index
, lin_from
, lin_to
. You supply name
, the name of the UDF. The variables c
, mixture_thread
, from_phase_index
, from_species_index
, to_phase_index
, to_species_index
, lin_from
, and lin_to
are passed by the Ansys Fluent solver
to your UDF. Your UDF will need to return the real value of the mass
transfer rate to the solver and the two linearized terms lin_from
and lin_to
.
Important: The linearization terms *lin_from
and *lin_to
are calculated based on the linearization coefficients and in Equation 2–15:
(2–17) |
(2–18) |
The derivative of mass transfer rate to pressure is stored in the above-mentioned macro in Equation 2–16.
The arguments from_species_index
and to_species_index
are relevant for multiphase species
transport problems only, and only if the respective phase has a mixture
material associated with it.
The following UDF, named cav_source
, specifies mass transfer source terms as a function of liquid vaporization
pressure and flow pressure.
Important: Note that in the example that follows, the DEFINE_LINEARIZED_MASS_TRANSFER
statement is broken up into three lines for the sake of readability.
In your source file, you must make sure that the DEFINE
statement is on one line only.
#include "udf.h" DEFINE_LINEARIZED_MASS_TRANSFER(cav_source,cell,thread,from_index,from_species_index, to_index, to_species_index, d_mdot_d_vof_from,d_mdot_d_vof_to) { real vof_nuc = RP_Get_Real("mp/cvt/cfx/vof-nuc"); real r_b = RP_Get_Real("mp/cvt/cfx/r-bubbles"); real F_evap = RP_Get_Real("mp/cvt/cfx/f-evap"); real F_cond = RP_Get_Real("mp/cvt/cfx/f-cond"); real c_evap = 3.0*F_evap*vof_nuc/r_b; real c_cond = 3.0*F_cond/r_b; real P_SAT = RP_Get_Real("mp/cvt/vapor-p"); Thread *liq = THREAD_SUB_THREAD(thread, from_index); Thread *vap = THREAD_SUB_THREAD(thread, to_index); real m_dot, dp, m_source; real p_op = RP_Get_Real ("operating-pressure"); real press = C_P(cell, thread) + p_op; real rho_l = C_R(cell,liq); real rho_v = C_R(cell,vap); real vof_l = C_VOF(cell,liq); real vof_v = C_VOF(cell,vap); real r_rho_lv = 1./rho_v - 1./rho_l; m_dot = 0.; m_source = 0.0; if (press <= P_SAT) { dp = P_SAT - press; dp = MAX(dp, 1e-4); m_dot = c_evap*rho_v*sqrt(2/3.0*dp/rho_l); m_source = m_dot*vof_l; *d_mdot_d_vof_from = m_dot; *d_mdot_d_vof_to = -m_dot; } else { dp = press - P_SAT; dp = MAX(dp, 1e-4); m_dot = -c_cond*rho_v*sqrt(2/3.0*dp/rho_l); m_source = m_dot*vof_v; *d_mdot_d_vof_from = m_dot; *d_mdot_d_vof_to = -m_dot; } /* ++++++++++ ds/dp term ++++++++++++++ */ if(NNULLP(THREAD_STORAGE(thread, SV_MT_DS_DP))) C_STORAGE_R(cell,thread,SV_MT_DS_DP) = ABS(r_rho_lv*m_source/(2*dp)); return m_source; }
After the UDF that you have defined using
DEFINE_LINEARIZED_MASS_TRANSFER
is interpreted (Interpreting UDFs) or compiled (Compiling UDFs), the name of the argument that you supplied as
the first DEFINE
macro argument will become visible and
selectable under Mass Transfer in the Phase
Interaction > Heat, Mass, Reactions >
Mass tab of the Multiphase Model dialog box
after you specify the Number of Mass Transfer Mechanisms. See Hooking DEFINE_LINEARIZED_MASS_TRANSFER
UDFs for details.
You can use DEFINE_MASS_TRANSFER
when
you want to model mass transfer in a multiphase problem. The mass
transfer rate specified using a DEFINE_MASS_TRANSFER
UDF is used to compute mass, momentum, energy, and species sources
for the phases involved in the mass transfer. For problems in which
species transport is enabled, the mass transfer will be from one species
in one phase, to another species in another phase. If one of the phases
does not have a mixture material associated with it, then the mass
transfer will be with the bulk fluid of that phase.
Note: By default, Fluent will attempt to hook mass transfer UDFs defined using
DEFINE_LINEARIZED_MASS_TRANSFER
. In order to hook a
DEFINE_MASS_TRANSFER
UDF to Fluent, you must first disable
the following text command:
solve/set/advanced/linearized-mass-transfer-udf
.
You may want to consider using the
DEFINE_LINEARIZED_MASS_TRANSFER
macro (
DEFINE_LINEARIZED_MASS_TRANSFER
) as it may provide a more
robust solution, even though the results may be the same when converged.
DEFINE_MASS_TRANSFER
(name
, c
, mixture_thread
, from_phase_index
, from_species_index
, to_phase_index
, to_species_index
)
Important: Note that all of the arguments to a DEFINE
macro need to be placed on the same line in your source code. Splitting
the DEFINE
statement onto several lines will
result in a compilation error.
Argument Type |
Description |
---|---|
|
UDF name. |
|
Index of cell on the thread pointed
to by |
|
Pointer to mixture-level thread. |
|
Index of phase from which mass is transferred. |
|
ID of species from which mass is transferred (ID= -1 if phase does not have mixture material). |
|
Index of phase to which mass is transferred. |
|
ID of species to which mass is transferred (ID= -1 if phase does not have mixture material). |
Function returns
real
There are seven arguments to DEFINE_MASS_TRANSFER
: name
, c
, mixture_thread
, from_phase_index
, from_species_index
, to_phase_index
, to_species_index
. You supply name
, the name of the UDF. The variables c
, mixture_thread
, from_phase_index
, from_species_index
, to_phase_index
, and to_species_index
are passed by the Ansys Fluent solver to your UDF. Your UDF will need
to return the real value of the mass transfer to the solver in the
units of kg//s.
Important: The arguments from_species_index
and to_species_index
are relevant for multiphase species
transport problems only, and only if the respective phase has a mixture
material associated with it.
The following UDF, named liq_gas_source
, specifies a simple
mass transfer coefficient based on saturation temperature:
Important: Note that in the example that follows, the
DEFINE_MASS_TRANSFER
statement is broken up into two lines
for the sake of readability. In your source file, you must make sure that the
DEFINE
statement is on one line only.
/* UDF to define a simple mass transfer based on Saturation Temperature. The "from" phase is the gas phase and the "to" phase is the liquid phase */ #include "udf.h" DEFINE_MASS_TRANSFER(liq_gas_source,cell,thread,from_index,from_species_index,to_index,to_species_index) { real m_lg; real T_SAT = 373.15; Thread *gas, *liq; gas = THREAD_SUB_THREAD(thread, from_index); liq = THREAD_SUB_THREAD(thread, to_index); m_lg = 0.0; if (C_T(cell, liq) > T_SAT) { /* Evaporating */ m_lg = -0.1*C_VOF(cell,liq)*C_R(cell,liq)* (C_T(cell,liq)-T_SAT)/T_SAT; } else if (C_T(cell, gas) < T_SAT) { /* Condensing */ m_lg = 0.1*C_VOF(cell,gas)*C_R(cell,gas)* (T_SAT-C_T(cell,gas))/T_SAT; } return (m_lg); }
In order to hook a DEFINE_MASS_TRANSFER
UDF to Fluent, you must
first disable the following text command:
solve/set/advanced/linearized-mass-transfer-udf
.
After the UDF that you have defined using DEFINE_MASS_TRANSFER
is
interpreted (Interpreting UDFs) or compiled (Compiling UDFs), the name of the argument that you supplied as
the first DEFINE
macro argument (for example,
liq_gas_source
) will become visible and selectable under
Mass Transfer in the Phase Interaction >
Heat, Mass, Reactions > Mass tab of the
Multiphase Model dialog box after you specify the Number
of Mass Transfer Mechanisms. See Hooking DEFINE_MASS_TRANSFER
UDFs for details.
You can use DEFINE_VECTOR_EXCHANGE_PROPERTY
to specify:
Custom slip velocities for the multiphase mixture model
Blending factors for flow regime modeling (supported with the algebraic slip model)
Turbulent dispersion forces for the multiphase Eulerian model.
DEFINE_VECTOR_EXCHANGE_PROPERTY
(name
,
c
, mixture_thread
,
phase_index_1
, phase_index_2
,
vector_result
)
Important: Note that all of the arguments to a DEFINE
macro need to be placed on the same line in your source code. Splitting
the DEFINE
statement onto several lines will
result in a compilation error.
Argument Type |
Description |
---|---|
|
UDF name. |
|
Cell index. |
|
Pointer to cell thread of mixture domain. |
|
Identifier of the first phase from the pair of phases for which you want to
customize a phase interaction property. An index of |
|
Identifier of the second phase from the pair of phases for which you want to
customize a phase interaction property. See |
|
Pointer to array containing vector properties such as slip velocity vector, list of blending factors, turbulent dispersion force |
Function returns
real
There are six arguments to DEFINE_VECTOR_EXCHANGE_PROPERTY
:
name
, c
,
mixture_thread
, phase_index_1
,
phase_index_2
, and vector_result
. You
supply name
, the name of the UDF. c
,
mixture_thread
, phase_index_1
,
phase_index_2
, and vector_result
are
variables that are passed by the Ansys Fluent solver to your UDF. Your UDF will need to set
the values referenced by the real
pointer to the
vector_result
array:
For the slip velocity vector, the components of the slip velocity vector are stored in
vector_result[0]
for 3D orvector_result[1]
for 2D.For the dispersion force, the components of the force vector are stored in
vector_result[0]
for 3D orvector_result[1]
for 2D.For the flow regime modeling, the blending factors are stored in
vector_result[0]
,vector_result[1]
,vector_result[2]
,vector_result[3]
.
Note: ANSYS Fluent passes values of phase identifiers in the same order as they appear in the Phase Pairs selection list in the Multiphase Model dialog box (Phase Interaction > Forces tab).
The following UDF, named custom_slip
, specifies a custom slip velocity in a two-phase mixture problem.
Important: Note that in the example that follows, the DEFINE_VECTOR_EXCHANGE_PROPERTY
statement is broken up into two lines for the sake of readability.
In your source file, you must make sure that the DEFINE
statement is on one line only.
/*************************************************************** UDF for a defining a custom slip velocity in a 2-phase mixture problem ****************************************************************/ #include "udf.h" DEFINE_VECTOR_EXCHANGE_PROPERTY(custom_slip,c,mixture_thread,phase_index_1, phase_index_2,vector_result) { real grav[2] = {0., -9.81}; real K = 5.e4; real pgrad_x, pgrad_y; Thread *pt, *st;/* thread pointers for primary and secondary phases*/ pt = THREAD_SUB_THREAD(mixture_thread, phase_index_1); st = THREAD_SUB_THREAD(mixture_thread, phase_index_2); /* at this point the phase threads are known for primary (0) and secondary(1) phases */ pgrad_x = C_DP(c,mixture_thread)[0]; pgrad_y = C_DP(c,mixture_thread)[1]; vector_result[0] = -(pgrad_x/K) +(((C_R(c, st)- C_R(c, pt))/K)* grav[0]); vector_result[1] = -(pgrad_y/K) +(((C_R(c, st)- C_R(c, pt))/K)* grav[1]); }
Important: Note that the pressure gradient macro C_DP
is now obsolete. A more current pressure gradient macro can be found
in Table 3.4: Macro for Cell Volume Defined in mem.h
.
The following UDF, named custom_td
, specifies a custom turbulent
dispersion force in a two-phase mixture problem (see Lopez de Bertodano Model in the Fluent Theory Guide for
details).
/* This example uses DEFINE_VECTOR_EXCHANGE_PROPERTY to implement the lopez-de-bertodano model */ /* and should be used for primary-secondary phase pair */ /* The Force to be returned is the force exerted on the dispersed phase */ /* Force = A * gradient of continuous phase - B * gradient of dispersed phase */ /* A and B are coefficients which vary depending on the turbulent dispersion models */ /* P_DRIFT_COEFF(c, tj) = coefficient of continuous phase volume fraction gradient = A */ /* S_DRIFT_COEFF(c, tj) = coefficient of dispersed phase volume fraction gradient = -B */ #include "udf.h" DEFINE_VECTOR_EXCHANGE_PROPERTY(custom_td, c, t, i, j, val) { /* i -- liquid-phase; j -- vapor-phase */ Thread **pt = THREAD_SUB_THREADS(t); Thread *ti = pt[i]; Thread *tj = pt[j]; real term; real rho_P = C_R(c,ti); real turb_k = (mp_ke_type==TURB_MP_KE_MIXTURE || mp_rsm_type==TURB_MP_RSM_MIXTURE)? C_K(c,t) : C_K(c,pt[P_PHASE]); term = -rho_P*turb_k; term *= C_VOLUME(c,t); NV_VS(val, =, C_VOF_G(c,tj),*,term); if(NNULLP(THREAD_STORAGE(tj,SV_MP_DRIFT_S_P_COEFF))) P_DRIFT_COEFF(c,tj) = 0.0; if(NNULLP(THREAD_STORAGE(tj,SV_MP_DRIFT_COEFF))) S_DRIFT_COEFF(c,tj) = term; }
The following UDF, named aiad
, defines custom blending
factors via exponential functions for the flow regime modeling (see Flow Regime Modeling in the Fluent Theory Guide for details).
/* This example uses exponential functions based blending factors for flow regime modeling */ #include "udf.h" #define COEFF 50 #define VF_LIMIT 0.3 DEFINE_VECTOR_EXCHANGE_PROPERTY(aiad,c,mixture_thread,phase_index_1, phase_index_2,vector_result) { Thread *pt, *st; int n; /* thread pointers for primary and secondary phases*/ pt = THREAD_SUB_THREAD(mixture_thread, phase_index_1); st = THREAD_SUB_THREAD(mixture_thread, phase_index_2); /* array initialization */ for (n = 0; n < 4; n++) vector_result[n] = 0.; real vof_1 = C_VOF(c, pt); real vof_2 = C_VOF(c, st); /* Calculating the blending factors */ /* if index_1 behaves as continuous and index_2 behaves as dispersed */ real Fcd = 1.0 / (1.0 + exp(COEFF * (vof_2 - VF_LIMIT))); /* if index_1 behaves as dispersed and index_2 behaves as continuous */ real Fdc = 1.0 / (1.0 + exp(COEFF * (vof_1 - VF_LIMIT))); /* if index_1 and index_2 both behave as continuous */ real Fcc = 1 - (Fcd + Fdc); /* if index_1 and index_2 both behave as dispersed */ real Fdd = 0.; vector_result[0] = Fcc; vector_result[1] = Fcd; vector_result[2] = Fdc; vector_result[3] = Fdd; }
After the UDF that you have defined using
DEFINE_VECTOR_EXCHANGE_PROPERTY
is interpreted (Interpreting UDFs) or compiled (Compiling UDFs), the name of the argument that you supplied as
the first DEFINE
macro argument (for example,
custom_slip
) will become visible and selectable in the
Phase Interaction tab of the Multiphase Model
dialog box in Ansys Fluent. See Hooking DEFINE_VECTOR_EXCHANGE_PROPERTY
UDFs for
details.