TrioCFD 1.9.8
TrioCFD documentation
Loading...
Searching...
No Matches
Triple_Line_Model_FT_Disc.h
1/****************************************************************************
2* Copyright (c) 2019, CEA
3* All rights reserved.
4*
5* Redistribution and use in source and binary forms, with or without modification, are permitted provided that the following conditions are met:
6* 1. Redistributions of source code must retain the above copyright notice, this list of conditions and the following disclaimer.
7* 2. Redistributions in binary form must reproduce the above copyright notice, this list of conditions and the following disclaimer in the documentation and/or other materials provided with the distribution.
8* 3. Neither the name of the copyright holder nor the names of its contributors may be used to endorse or promote products derived from this software without specific prior written permission.
9*
10* THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE DISCLAIMED.
11* IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS;
12* OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
13*
14*****************************************************************************/
15
16#ifndef Triple_Line_Model_FT_Disc_included
17#define Triple_Line_Model_FT_Disc_included
18
19#include <Transport_Interfaces_FT_Disc.h>
20#include <Domaine_VDF.h>
21#include <FTd_tools.h>
22#include <TRUST_Ref.h>
23
27class Equation_base;
28class Param;
29
31{
32 Declare_instanciable_sans_constructeur( Triple_Line_Model_FT_Disc ) ;
33public :
35
36 void associer_pb(const Probleme_base&);
37 void initialize();
38 void completer();
39 void set_param(Param& p) const override;
40 double get_Qtcl(const int num_face);
41 const double& get_ls() const { return ls_; }
42 const double& get_ym() const { return ym_; }
43 const double& get_xm() const { return xm_; }
44 double get_sm() { return sm_; }
45 const double& get_initial_CL_xcoord() const { return initial_CL_xcoord_; }
46 const double& get_rhocpl() const { return rhocpl_; }
47 bool is_activated() const { return activated_; }
49 bool is_read_via_file() const { return read_via_file_; }
50 bool Rc_tcl_GridN() const { return Rc_tcl_GridN_; }
51 const double& Rc_inject() const { return Rc_inject_; }
52 const double& thetaC_tcl() const { return thetaC_tcl_; }
53 const double& t_injection() const { return t_injection_; }
54 double& t_injection() { return t_injection_; }
55 bool reinjection_tcl() const { return reinjection_tcl_; }
57 bool ready_inject_tcl() const { return ready_inject_tcl_; }
59 const double& tempC_tcl() const { return tempC_tcl_; }
60 const int& tag_tcl() const { return tag_tcl_; }
61 bool lissage_tcl() const { return lissage_tcl_; }
62
63 const bool& adjust_meso_ML() const
64 {
65 return adjust_meso_ML_;
66 }
67 const double& Ri() const
68 {
69 return Ri_;
70 }
71
72 double compute_capillary_number() const;
73 int get_any_tcl_face() const;
74
75 // Using the exact in/out intersections of the interface within the cell (works only in 2D, and may assume the wall in y-)
76 void get_in_out_coords(const Domaine_VDF& zvdf, const int elem, const double dt, DoubleTab& in_out, FTd_vecteur3& norm_elem, double& surface_tot);
77
78 // computes an approximate interface position (equivalent mean-plane)
79 // returns in_out table and the mean normal in elem, and the interface surface
80 void compute_approximate_interface_inout(const Domaine_VDF& zvdf, const int korient, const int elem, const int num_face, DoubleTab& in_out, FTd_vecteur3& norm_elem, double& surface_tot) const;
81
82 double compute_local_cos_theta(const Parcours_interface& parcours, const int num_face, const FTd_vecteur3& norm_elem) const;
83
84 // Computes the integral flux in a cell in the meso zone.
85 // in and out are the interface position on both sides of the cell.
86 double compute_Qint(const DoubleTab& in_out, const double theta_app_locs, const int num_face, double& Qmeso) const;
88 void compute_TCL_fluxes_in_all_boundary_cells(ArrOfInt& elems_with_CL_contrib, ArrOfInt& faces_with_CL_contrib, ArrOfDouble& mpoint_from_CL, ArrOfDouble& Q_from_CL);
89
90 void corriger_mpoint(DoubleTab& mpoint) const;
91 void corriger_secmem(const double coef, DoubleTab& secmem2) const;
92 void correct_wall_adjacent_temperature_according_to_TCL_fluxes(DoubleTab& temperature) const;
93 void set_wall_adjacent_temperature_according_to_TCL_model(DoubleTab& temperature) const;
94 void correct_TCL_energy_evolution(DoubleTab& temperature) const;
95 double get_theta_app(const int num_face);
96 void correct_theta_app_qtcl(double& theta_app_, double& qtcl, const int num_face_wall, const int num_som, const DoubleTab& vit) const;
97 const double& get_theta_app() const { return theta_app_; }
98 const double& get_Qtcl() const { return Qtcl_; }
99
100
102
103 inline ArrOfInt& elems();
104
105 inline const ArrOfInt& elems() const
106 {
107 assert(tag_tcl() == ref_eq_interf_->maillage_interface().get_mesh_tag());
108 return elems_;
109 }
110
111 inline ArrOfInt& boundary_faces() { return boundary_faces_; }
112
113 inline const ArrOfInt& boundary_faces() const
114 {
115 assert(tag_tcl() == ref_eq_interf_->maillage_interface().get_mesh_tag());
116 return boundary_faces_;
117 }
118
119 inline ArrOfDouble& mp() { return mp_; }
120 inline const ArrOfDouble& mp() const { return mp_; }
121 inline ArrOfDouble& Q() { return Q_; }
122 inline const ArrOfDouble& Q() const { return Q_; }
123
124protected:
125
126 bool activated_ = false;
127 bool deactivate_ = false;
130 int n_ext_meso_; // number of layers in the extension of the meso-zone.
132 double Qtcl_; //
133 double ls_; // The slip length...
135 double x_cl_; // position of the CL.
136 // The end of micro-region:
137 double xm_;
138 double ym_;
139 double sm_;
140 // End of meso region
141 double ymeso_;
142// double old_xcl_ ; // Former position of the contact line. Not nice, but needed to copute the CL velocity
143 double initial_CL_xcoord_; // Former position of the contact line. Not nice, but needed to copute the CL velocity
144 double kl_cond_; // We store the liquid conductivity for easy access.
145 double rhocpl_;
147 // Tabulation TCL model
150 DoubleTab tab_Mtcl_;
151 // total nb of columns
153 // num -column for temperature
155 // num column for theta app
157 // num -column for Qtcl
159
160 // parameters to force numerical breakup
161 // Radius of nucleate site; [in number of grids]
162 // Distance below which the theta_app will be replaced by a big contact angle thetaC_tcl_ [in degree]
163 // to force bubble pinching / necking
167 // parameters injection seed nucleate
168 // Size: Using the same Radius defined above (Rc_tcl_gridN_), [in number of grids]
169 // Use the average temperature in the same zone to tell if nucleate site is activted and inject bubbles
170 // Temperature beyond which the site is activeted [in K] T-Tsat
171 // in this case, we need the initial value for contact angle, thetaC_tcl_ will be used
172 // reinjection: if a re-injection of interface if nessaire
173 bool reinjection_tcl_ = false;
175 bool ready_inject_tcl_ = false;
176 // summaries:
177 // - when reinjection_tcl_ is 0; i.e., by default, the pinching breakup with be used,
178 // Rc_tcl_GridN_ thetaC_tcl_ control the breakup, thetaC_tcl_ should be large
179 // - when reinjection_tcl_ is 1; i.e., there is no the pinching breakup. The bubble goes away directly, we reinject seed
180 // Rc_tcl_GridN_ thetaC_tcl_ control the initial shape, thetaC_tcl_ should be small...
181
182 // if to distribute the Qtcl into the first facette
184
185 // if to correct of qmeso in Micro-layer
186 // par default, it should NOT be activated
187 // interface thermal resistence solid-fluid
188 bool adjust_meso_ML_ = false;
189 double Ri_;
190
191 // if to activate lissage of Qmicro and Theta_app
192 // over time period, [t_injection_, t]
193 bool lissage_tcl_ = false;
195
196 // Information on the TCL region :
197 // Note that the same elem may appear twice in the list, once for the micro contribution, once for the meso.
198 ArrOfInt elems_; // The list of elements containing either the TCL itself, or the meso domaine.
199 ArrOfInt boundary_faces_; // corresponding list of faces on the wall boundary.
200 ArrOfDouble mp_; // corresponding value to be assigned to each cell (as a mass flux m [unit?])
201 ArrOfDouble Q_; // corresponding value to be assigned to each cell (as a thermal flux Q [unit?])
202
203 // Variables for post-processing :
213
215
220};
221
223{
224 const Transport_Interfaces_FT_Disc& eq_transport = ref_eq_interf_.valeur();
225 const Maillage_FT_Disc& maillage = eq_transport.maillage_interface();
226 if (tag_tcl() == maillage.get_mesh_tag())
227 {
228 Cerr << "Why do you need to access elems_ in a non-const maner when tag_cl_= " << tag_tcl()
229 << " and maillage.get_mesh_tag()= " << maillage.get_mesh_tag()
230 << " are equal (which means by-the-way that the model should be up-to-date already!)" << finl;
232 }
233 return elems_;
234}
235
236#endif /* Triple_Line_Model_FT_Disc_included */
class Domaine_VDF
Definition Domaine_VDF.h:64
classe Equation_base Le role d'une equation est le calcul d'un ou plusieurs champs....
: class Maillage_FT_Disc Cette classe decrit un maillage:
int get_mesh_tag() const
return mesh_state_tag_
class Nom Une chaine de caractere pour nommer les objets de TRUST
Definition Nom.h:31
Objet_U()
Constructeur par defaut : attribue un numero d'identifiant unique a l'objet (object_id_),...
Definition Objet_U.cpp:55
Helper class to factorize the readOn method of Objet_U classes.
Definition Param.h:112
classe Probleme_base C'est un Probleme_U qui n'est pas un couplage.
static void exit(int exit_code=-1)
Routine de sortie de TRUST dans une region Kokkos.
Definition Process.cpp:455
const Maillage_FT_Disc & maillage_interface() const
void set_wall_adjacent_temperature_according_to_TCL_model(DoubleTab &temperature) const
const double & get_theta_app() const
OBS_PTR(Probleme_base) pb_base_
const double & get_initial_CL_xcoord() const
const ArrOfInt & boundary_faces() const
OBS_PTR(Convection_Diffusion_Temperature_FT_Disc) ref_eq_temp_
void corriger_secmem(const double coef, DoubleTab &secmem2) const
OBS_PTR(Transport_Interfaces_FT_Disc) ref_eq_interf_
void corriger_mpoint(DoubleTab &mpoint) const
void correct_TCL_energy_evolution(DoubleTab &temperature) const
void associer_pb(const Probleme_base &)
OBS_PTR(Navier_Stokes_FT_Disc) ref_ns_
const ArrOfDouble & mp() const
void correct_theta_app_qtcl(double &theta_app_, double &qtcl, const int num_face_wall, const int num_som, const DoubleTab &vit) const
double compute_local_cos_theta(const Parcours_interface &parcours, const int num_face, const FTd_vecteur3 &norm_elem) const
double compute_Qint(const DoubleTab &in_out, const double theta_app_locs, const int num_face, double &Qmeso) const
void set_param(Param &p) const override
void correct_wall_adjacent_temperature_according_to_TCL_fluxes(DoubleTab &temperature) const
void compute_approximate_interface_inout(const Domaine_VDF &zvdf, const int korient, const int elem, const int num_face, DoubleTab &in_out, FTd_vecteur3 &norm_elem, double &surface_tot) const
const ArrOfDouble & Q() const
void get_in_out_coords(const Domaine_VDF &zvdf, const int elem, const double dt, DoubleTab &in_out, FTd_vecteur3 &norm_elem, double &surface_tot)