MoDeNa  1.0
Software framework facilitating sequential multi-scale modelling
gasMixtureConductivity.py
1 
31 
32 
42 
43 import os
44 from modena import *
45 from fireworks.utilities.fw_utilities import explicit_serialize
46 import gasConductivity
47 
48 
49 weightedAverageCode=r'''
50 #include "modena.h"
51 #include "math.h"
52 #include "stdio.h"
53 
54 void gasMixtureConductivity
55 (
56  const modena_model_t* model,
57  const double* inputs,
58  double *outputs
59 )
60 {
61  {% block variables %}{% endblock %}
62 
63  double kgasmix=0; // gas mixture conductivity
64  double sumx=0;
65  int i;
66  for (i=0;i<x_size;i++) {
67  sumx = sumx + x[i];
68  }
69  for (i=0;i<x_size;i++) {
70  kgasmix = kgasmix + gas_thermal_conductivity[i]*x[i]/sumx;
71  }
72  outputs[0] = kgasmix;
73 }
74 '''
75 
78 dohrnCode=r'''
79 #include "modena.h"
80 #include "math.h"
81 #include "stdio.h"
82 
83 void gasMixtureConductivity
84 (
85  const modena_model_t* model,
86  const double* inputs,
87  double *outputs
88 )
89 {
90  {% block variables %}{% endblock %}
91 
92  double kgasmix=0; // gas mixture conductivity
93  double Rg=8.314; // gas constant
94  double y=0;
95  double xi[4]; // normalized molar fractions
96  double k[4]; // thermal conductivities
97  double Tc[4]={304.17,511.7,154.58,126.19}; // critical temperatures
98  double pc[4]={7.386e6,45.1e5,5.043e6,3.396e6}; // critical pressures
99  double M[4]={44e-3,70e-3,32e-3,28e-3}; // molar masses
100  double A[4][4];
101  double gam[4];
102  int i,j;
103  for (i=0;i<x_size;i++) {
104  y = y + x[i];
105  k[i]=gas_thermal_conductivity[i];
106  }
107  for (i=0;i<x_size;i++) {
108  xi[i] = x[i]/y; // normalize molar fractions
109  gam[i]=210*pow(Tc[i]*pow(M[i],3)/pow(pc[i],4),1.0/6.0);
110  }
111  for (i=0;i<x_size;i++) {
112  for (j=0;j<x_size;j++) {
113  y=gam[j]*(exp(0.0464*T/Tc[i])-exp(-0.2412*T/Tc[i]))/\
114  gam[i]/(exp(0.0464*T/Tc[j])-exp(-0.2412*T/Tc[j]));
115  A[i][j]=pow(1+pow(y,0.5)*pow(M[i]/M[j],0.25),2)/\
116  pow(8*(1+M[i]/M[j]),0.5);
117  }
118  }
119  for (i=0;i<x_size;i++) {
120  y=0;
121  for (j=0;j<x_size;j++) {
122  y=y+xi[j]*A[i][j];
123  }
124  kgasmix = kgasmix + xi[i]*k[i]/y;
125  }
126  outputs[0] = kgasmix;
127 }
128 '''
129 
132 lindsayBromleyCode=r'''
133 #include "modena.h"
134 #include "math.h"
135 #include "stdio.h"
136 
137 void gasMixtureConductivity
138 (
139  const modena_model_t* model,
140  const double* inputs,
141  double *outputs
142 )
143 {
144  {% block variables %}{% endblock %}
145 
146  double kgasmix=0; // gas mixture conductivity
147  double Rg=8.314; // gas constant
148  double y=0;
149  double xi[4]; // normalized molar fractions
150  double k[4]; // thermal conductivities
151  double Tb[4]={194.75,322.4,90.19,77.36}; // boiling point temperatures
152  double cp[4]; // thermal capacities at constant pressure
153  double M[4]={44e-3,70e-3,32e-3,28e-3}; // molar masses
154  double S[4]; // Sutherland constants
155  double gam[4]; // heat capacity ratio
156  double cv[4]; // thermal capacities at constant volume
157  double A[4][4];
158  double t=T/1e3;
159  int i,j;
160  cp[0]=24.99735+55.18696*t-33.69137*t*t+7.948387*t*t*t-0.136638/t/t; // CO2
161  cp[1]=-25.6132057+226.4176882*t+574.2688767*t*t-670.5517907*t*t*t+\
162  0.6765321/t/t; // CyP
163  cp[2]=31.32234-20.23531*t+57.86644*t*t-36.50624*t*t*t-0.007374/t/t; // O2
164  cp[3]=28.98641+1.853978*t-9.647459*t*t+16.63537*t*t*t+0.000117/t/t; // N2
165  for (i=0;i<x_size;i++) {
166  y = y + x[i];
167  k[i]=gas_thermal_conductivity[i];
168  }
169  for (i=0;i<x_size;i++) {
170  // normalize molar fractions
171  xi[i] = x[i]/y;
172  S[i] = 1.5*Tb[i];
173  cv[i] = cp[i] - Rg;
174  gam[i] = cp[i]/cv[i];
175  }
176  for (i=0;i<x_size;i++) {
177  for (j=0;j<x_size;j++) {
178  y=k[i]/k[j]*cp[j]/cp[i]*(9-5/gam[j])/(9-5/gam[i]);
179  A[i][j]=0.25*pow(1+pow(y*pow(M[j]/M[i],0.75)*(T+S[i])/\
180  (T+S[j]),0.5),2)*(T+sqrt(S[i]*S[j]))/(T+S[j]);
181  }
182  }
183  for (i=0;i<x_size;i++) {
184  y=0;
185  for (j=0;j<x_size;j++) {
186  y=y+xi[j]*A[i][j];
187  }
188  kgasmix = kgasmix + xi[i]*k[i]/y;
189  }
190  outputs[0] = kgasmix;
191 }
192 '''
193 
194 
197 f_gasMixtureConductivity = CFunction(
198  Ccode=dohrnCode,
199  # These are global bounds for the function
200  inputs={
201  'T': {'min': 273, 'max': 550},
202  'x': {'index': gasConductivity.species, 'min': 0, 'max': 1},
203  'gas_thermal_conductivity': {'index': gasConductivity.species,
204  'min': 0, 'max': 1},
205  },
206  outputs={
207  'gasMixtureConductivity': {'min': 0, 'max': +9e99, 'argPos': 0},
208  },
209  parameters={
210  'param0': {'min': -9e99, 'max': +9e99, 'argPos': 0},
211  'param1': {'min': -9e99, 'max': +9e99, 'argPos': 1},
212  'param2': {'min': -9e99, 'max': +9e99, 'argPos': 2},
213  },
214 )
215 
216 
219 m_gasMixtureConductivity = ForwardMappingModel(
220  _id='gasMixtureConductivity',
221  surrogateFunction=f_gasMixtureConductivity,
222  substituteModels=[gasConductivity.m_CO2_thermal_conductivity,\
223  gasConductivity.m_CyP_thermal_conductivity,\
224  gasConductivity.m_O2_thermal_conductivity,\
225  gasConductivity.m_N2_thermal_conductivity],
226  parameters=[1, 1, 1],
227 )