MoDeNa  1.0
Software framework facilitating sequential multi-scale modelling
TSSource.H
Go to the documentation of this file.
1 
14 double dXOHdt,XOH_C,XOH_Old,dWdt,XW_C,XW_Old,dLdt,alphaFO,
16  rhoFO,C_TOT, L_l_C, L_l_Old, T_C,T_O, wBA_g_val, xBL_val, rhoPolyMixture;
17 
18 // Calling the model for density reaction mixture
19 size_t Tm_pos =
20  modena_model_inputs_argPos(density_reaction_mixturemodel, "T");
21 size_t XOHv_pos =
22  modena_model_inputs_argPos(density_reaction_mixturemodel, "XOH");
23 modena_model_argPos_check(density_reaction_mixturemodel);
24 
25 forAll(mesh.C(), celli)
26 {
27  alphaFO = alpha2[celli];
28  rhoFO = rho_foam[celli];
29 
30  XOH_C = XOH[celli];
31  XOH_Old = XOH.oldTime()[celli];
32  dXOHdt = (mag(XOH_C-XOH_Old))/(max(runTime.deltaTValue(),ROOTVSMALL));
33 
34  XW_C = XW[celli];
35  XW_Old = XW.oldTime()[celli];
36  dWdt = (mag(XW_C-XW_Old))/(max(runTime.deltaTValue(),ROOTVSMALL));
37 
38  T_C = TS[celli];
39  T_O = TS.oldTime()[celli];
40 
41  wBA_g_val = wBA_g[celli];
42 
43  if (PBESwitch)
44  {
45  L_l_C = wBA_l[celli];
46  L_l_Old = wBA_l.oldTime().oldTime()[celli];
47  dLdt = (L_l_C - L_l_Old)/(max(runTime.deltaTValue(),ROOTVSMALL));
48  }
49  else
50  {
51  double lmax;
52  lmax = LliqMax(T_C);
53  if (lmax < L0)
54  {
55  dLdt =
56  (
57  -1.0*ddT_LliqMax(T_C)*((T_C-T_O)
58  /(max(runTime.deltaTValue(),ROOTVSMALL)))
59  );
60  }
61  else
62  {
63  dLdt = 0.0;
64  }
65  }
66 
67  if (liquidMixtureDensitySurrogate)
68  {
69  // set input vector
70  modena_inputs_set(inputs_den, Tm_pos, T_C);
71  modena_inputs_set(inputs_den, XOHv_pos, XOH_C);
72  // call the model
73  int ret_den =
75  (
76  density_reaction_mixturemodel,
77  inputs_den,
78  outputs_den
79  );
80  if (ret_den != 0)
81  {
82  modena_inputs_destroy (inputs_den);
83  modena_outputs_destroy (outputs_den);
84  modena_model_destroy (density_reaction_mixturemodel);
85  exit(ret_den);
86  }
87  rhoPolyMixture = modena_outputs_get(outputs_den, 0);
88  }
89  else
90  {
91  rhoPolyMixture = rhoPoly;
92  }
93 
94  if (blowingAgent == "n-pentane")
95  {
96  C_TOT = C_Poly;
97  if (alphaFO > 0.5 && XOH[celli] < XOH_Gel)
98  {
99  TSSource[celli] =
100  (
101  ((-DH_OH*COH_0*dXOHdt)/(rhoPolyMixture*C_TOT)
102  + (-DH_W*CW_0*dWdt)/(rhoPolyMixture*C_TOT)
103  + (latenth*dLdt)/(C_TOT))
104  );
105  }
106  else if (alphaFO > 0.5 && XOH[celli] >= XOH_Gel)
107  {
108  TSSource[celli] = (latenth*dLdt)/(C_TOT);
109  }
110  else
111  {
112  TSSource[celli] = ROOTVSMALL;
113  }
114  }
115 
116  if (blowingAgent == "R-11")
117  {
118  if (wBA_g_val > 0.0)
119  {
120  xBL_val = xBL(T_C, dxdT);
121 
122  C_TOT =
123  (
124  C_Poly + (-(M_B/M_NCO)*
125  (1.0/(Foam::pow((1.0 - xBL_val),2)))*(dxdT))*latenth
126  );
127 
128  if (alphaFO > 0.5 && XOH[celli] < XOH_Gel)
129  {
130  TSSource[celli] =
131  (
132  ((-DH_OH*COH_0*dXOHdt)/(rhoPolyMixture*C_TOT)
133  + (-DH_W*CW_0*dWdt)/(rhoPolyMixture*C_TOT))
134  );
135  }
136  else if (alphaFO > 0.5 && XOH[celli] >= XOH_Gel)
137  {
138  TSSource[celli] =
139  (
140  ((-(M_B/M_NCO)*
141  (1.0/(Foam::pow((1.0 - xBL_val),2)))*(dxdT))*latenth)
142  );
143  }
144  else
145  {
146  TSSource[celli] = ROOTVSMALL;
147  }
148  }
149  else
150  {
151  C_TOT = C_Poly;
152  if (alphaFO > 0.5 && XOH[celli] < XOH_Gel)
153  {
154  TSSource[celli] =
155  (
156  ((-DH_OH*COH_0*dXOHdt)/(rhoPolyMixture*C_TOT)
157  + (-DH_W*CW_0*dWdt)/(rhoPolyMixture*C_TOT))
158  );
159  }
160  else
161  {
162  TSSource[celli] = ROOTVSMALL;
163  }
164  }
165  }
166  if (blowingAgent == "no")
167  {
168  C_TOT = C_Poly;
169  if (alphaFO > 0.5 && XOH[celli] < XOH_Gel)
170  {
171  TSSource[celli] =
172  (
173  ((-DH_OH*COH_0*dXOHdt)/(rhoPolyMixture*C_TOT)
174  + (-DH_W*CW_0*dWdt)/(rhoPolyMixture*C_TOT))
175  );
176  }
177  else
178  {
179  TSSource[celli] = ROOTVSMALL;
180  }
181  }
182  if (thermalConductivitySurrogateSwitch)
183  {
184  if (alphaFO > 0.5)
185  {
186  // surrogate model for thermal conductivity for CO2, cyclo-pentane
187  modena_inputs_set(inputs_thermalConductivity, porosity_Pos, max(1.0 - rhoFO/rhoPolyMixture,0.0));
188  scalar radiusBubble;
189  radiusBubble = bubbleRadius(mZero[celli], mOne[celli]);
190  modena_inputs_set(inputs_thermalConductivity, cell_size_Pos, (2.0*radiusBubble));
191  modena_inputs_set(inputs_thermalConductivity, temp_Pos, TS[celli]);
192  scalar pp_CO2, pp_BA;
193  pp_CO2 =
195  (
196  M_CO2, M_B, surfaceTension, wCO2_g[celli],
197  wBA_g[celli], p[celli], radiusBubble
198  );
199  pp_BA =
201  (
202  M_B, M_CO2, surfaceTension, wBA_g[celli], wCO2_g[celli],
203  p[celli], radiusBubble
204  );
205  modena_inputs_set(inputs_thermalConductivity, X_CO2_Pos, (pp_CO2/(pp_BA+pp_CO2)));
206  modena_inputs_set(inputs_thermalConductivity, X_O2_Pos, 0.0);
207  modena_inputs_set(inputs_thermalConductivity, X_N2_Pos, 0.0);
208  modena_inputs_set(inputs_thermalConductivity, X_Cyp_Pos, (pp_BA/(pp_BA+pp_CO2)));
209  // input for strut content
210  modena_inputs_set(inputs_strutContent, rho_foam_Pos, rhoFO);
211  int ret_strutContent = modena_model_call (strutContentmodel, inputs_strutContent, outputs_strutContent);
212  scalar st_c;
213  st_c = modena_outputs_get(outputs_strutContent, 0);
214  modena_inputs_set(inputs_thermalConductivity, strut_c_Pos, st_c);
215  int ret_thermalConductivitymodel = modena_model_call (thermalConductivitymodel, inputs_thermalConductivity, outputs_thermalConductivity);
216  if(modena_error_occurred())
217  {
218  exit(modena_error());
219  }
220  thermalConductivity[celli] = modena_outputs_get(outputs_thermalConductivity, 0);
221  thermalDiff_foam[celli] = (thermalConductivity[celli]/(rhoFO*C_Poly));
222  }
223  else
224  {
225  thermalDiff_foam[celli] = thermalDiffusivityGas(T_C);
226  }
227  thermalDiff_gas[celli] = thermalDiffusivityGas(T_C);
228  }
229  else
230  {
231  if (blowingAgent == "n-pentane")
232  {
233  if (alphaFO > 0.5)
234  {
235  if (rhoFO > 48)
236  {
237  thermalConductivity[celli] =
238  (
239  (scalar(8.7006e-8)*Foam::pow(rhoFO,2)
240  + scalar(8.4674e-5)*rhoFO + scalar(1.16e-2))
241  );
242  }
243  else
244  {
245  thermalConductivity[celli] =
246  (
247  (scalar(9.3738e-6)*Foam::pow(rhoFO,2)
248  - scalar(7.3511e-4)*rhoFO+scalar(2.956e-2))
249  );
250  }
251  thermalDiff_foam[celli] = (thermalConductivity[celli]/(rhoFO*C_Poly));
252  }
253  else
254  {
255  thermalDiff_foam[celli] = thermalDiffusivityGas(T_C);
256  }
257 
258  thermalDiff_gas[celli] = thermalDiffusivityGas(T_C);
259  }
260  else
261  {
262  if (alphaFO > 0.5)
263  {
264  thermalConductivity[celli] = 0.03;
265  thermalDiff_foam[celli] =
266  (thermalConductivity[celli]/(rhoFO*C_Poly));
267  }
268  else
269  {
270  thermalDiff_foam[celli] = thermalDiffusivityGas(T_C);
271  }
272  thermalDiff_gas[celli] = thermalDiffusivityGas(T_C);
273  }
274  }
275  thermalDiffusivity[celli] =
276  alpha1[celli]*thermalDiff_gas[celli]
277  + alpha2[celli]*thermalDiff_foam[celli];
278 }
279 thermalDiffusivity.correctBoundaryConditions();
280 ///@endcond
size_t modena_model_inputs_argPos(const modena_model_t *self, const char *name)
Function determining position of an argument in the input vector.
Definition: model.c:366
size_t X_CO2_Pos
memory allocation for the partial pressure of CO2 as the input of thermal conductivity surrogate mode...
Definition: modenaData.H:194
size_t temp_Pos
memory allocation for the temperature as the input of thermal conductivity surrogate model ...
Definition: modenaData.H:191
size_t porosity_Pos
memory allocation for the gas volume fraction as the input of thermal conductivity surrogate model ...
Definition: modenaData.H:182
size_t X_N2_Pos
memory allocation for the partial pressure of N2 as the input of thermal conductivity surrogate model...
Definition: modenaData.H:203
size_t cell_size_Pos
memory allocation for the bubble size as the input of thermal conductivity surrogate model ...
Definition: modenaData.H:185
modena_model_t * thermalConductivitymodel
pointer to the surrogate model for thermal conductivity
Definition: modenaData.H:22
size_t X_Cyp_Pos
memory allocation for the partial pressure of cyclo-pentane as the input of thermal conductivity surr...
Definition: modenaData.H:197
double partialPressureCO2(const state_type &y)
partial pressure of CO2
double surfaceTension
required for the computation of partial pressure
double M_CO2
Molecular mass of carbon dioxide, kg/kmol.
double M_NCO
Molecular weight of NCO, kg/kmol.
int modena_model_call(modena_model_t *self, modena_inputs_t *inputs, modena_outputs_t *outputs)
Function calling the surrogate model and checking for errors.
Definition: model.c:553
double DH_OH
Reaction heat for the gelling reaction, J/mol.
double ddT_LliqMax(double &)
derivative of LliqMax with respect to temperature
double C_TOT
Total specifc heat of the mixture.
double bubbleRadius(const double m0, const double m1)
radius of bubbles based on the moments
Definition: bubbleRadius.h:28
size_t rho_foam_Pos
memory allocation for the foam density as the input of thermal conductivity surrogate model ...
Definition: modenaData.H:179
double LliqMax(double &)
maximum allowable amount of liquid blowing agent (n-pentane) in liquid
double rhoPoly
density of the liquid polymer, kg/m3
double thermalDiffusivityGas(double &T)
thermal diffusivity of gas as a function of temperature
void modena_model_destroy(modena_model_t *self)
Function deallocating the memory allocated for the surrogate model.
Definition: model.c:665
double M_B
Molecular mass of blowing agent, kg/kmol.
double xBL(double &T, double &dxdT)
mole fraction of blowing agent (R-11) in liquid polymer
double partialPressureBA(const state_type &y)
partial pressure of the physical blowing agent
double L0
Initial weight fraction of blowing agent in the liquid, -.
size_t strut_c_Pos
memory allocation for the strut content as the input of thermal conductivity surrogate model ...
Definition: modenaData.H:188
size_t X_O2_Pos
memory allocation for the partial pressure of O2 as the input of thermal conductivity surrogate model...
Definition: modenaData.H:200
void modena_model_argPos_check(const modena_model_t *self)
Function checking that the user has queried all input positions.
Definition: model.c:408
double DH_W
Reaction heat for the blowing reaction, J/mol.