void AqueousTransport::getMixDiffCoeffs(doublereal* const d) { update_T(); update_C(); // update the binary diffusion coefficients if necessary if (!m_bindiff_ok) { updateDiff_T(); } size_t k, j; doublereal mmw = m_thermo->meanMolecularWeight(); doublereal sumxw = 0.0, sum2; doublereal p = m_press; if (m_nsp == 1) { d[0] = m_bdiff(0,0) / p; } else { for (k = 0; k < m_nsp; k++) { sumxw += m_molefracs[k] * m_mw[k]; } for (k = 0; k < m_nsp; k++) { sum2 = 0.0; for (j = 0; j < m_nsp; j++) { if (j != k) { sum2 += m_molefracs[j] / m_bdiff(j,k); } } if (sum2 <= 0.0) { d[k] = m_bdiff(k,k) / p; } else { d[k] = (sumxw - m_molefracs[k] * m_mw[k])/(p * mmw * sum2); } } } }
void AqueousTransport::getSpeciesFluxesExt(size_t ldf, doublereal* const fluxes) { update_T(); update_C(); getMixDiffCoeffs(DATA_PTR(m_spwork)); const vector_fp& mw = m_thermo->molecularWeights(); const doublereal* y = m_thermo->massFractions(); doublereal rhon = m_thermo->molarDensity(); // Unroll wrt ndim vector_fp sum(m_nDim,0.0); for (size_t n = 0; n < m_nDim; n++) { for (size_t k = 0; k < m_nsp; k++) { fluxes[n*ldf + k] = -rhon * mw[k] * m_spwork[k] * m_Grad_X[n*m_nsp + k]; sum[n] += fluxes[n*ldf + k]; } } // add correction flux to enforce sum to zero for (size_t n = 0; n < m_nDim; n++) { for (size_t k = 0; k < m_nsp; k++) { fluxes[n*ldf + k] -= y[k]*sum[n]; } } }
void GasTransport::getMixDiffCoeffsMass(doublereal* const d) { update_T(); update_C(); // update the binary diffusion coefficients if necessary if (!m_bindiff_ok) { updateDiff_T(); } doublereal mmw = m_thermo->meanMolecularWeight(); doublereal p = m_thermo->pressure(); if (m_nsp == 1) { d[0] = m_bdiff(0,0) / p; } else { for (size_t k=0; k<m_nsp; k++) { double sum1 = 0.0; double sum2 = 0.0; for (size_t i=0; i<m_nsp; i++) { if (i==k) { continue; } sum1 += m_molefracs[i] / m_bdiff(k,i); sum2 += m_molefracs[i] * m_mw[i] / m_bdiff(k,i); } sum1 *= p; sum2 *= p * m_molefracs[k] / (mmw - m_mw[k]*m_molefracs[k]); d[k] = 1.0 / (sum1 + sum2); } } }
/* * Units for the returned fluxes are kg m-2 s-1. * * * The diffusive mass flux of species \e k is computed from * \f[ * \vec{j}_k = -n M_k D_k \nabla X_k. * \f] * * @param ndim Number of dimensions in the flux expressions * @param grad_T Gradient of the temperature * (length = ndim) * @param ldx Leading dimension of the grad_X array * (usually equal to m_nsp but not always) * @param grad_X Gradients of the mole fraction * Flat vector with the m_nsp in the inner loop. * length = ldx * ndim * @param ldf Leading dimension of the fluxes array * (usually equal to m_nsp but not always) * @param fluxes Output of the diffusive mass fluxes * Flat vector with the m_nsp in the inner loop. * length = ldx * ndim */ void MixTransport::getSpeciesFluxes(int ndim, const doublereal* grad_T, int ldx, const doublereal* grad_X, int ldf, doublereal* fluxes) { int n, k; update_T(); update_C(); getMixDiffCoeffs(DATA_PTR(m_spwork)); const array_fp& mw = m_thermo->molecularWeights(); const doublereal* y = m_thermo->massFractions(); doublereal rhon = m_thermo->molarDensity(); vector_fp sum(ndim,0.0); for (n = 0; n < ndim; n++) { for (k = 0; k < m_nsp; k++) { fluxes[n*ldf + k] = -rhon * mw[k] * m_spwork[k] * grad_X[n*ldx + k]; sum[n] += fluxes[n*ldf + k]; } } // add correction flux to enforce sum to zero for (n = 0; n < ndim; n++) { for (k = 0; k < m_nsp; k++) { fluxes[n*ldf + k] -= y[k]*sum[n]; } } }
void GasTransport::getMixDiffCoeffsMole(doublereal* const d) { update_T(); update_C(); // update the binary diffusion coefficients if necessary if (!m_bindiff_ok) { updateDiff_T(); } doublereal p = m_thermo->pressure(); if (m_nsp == 1) { d[0] = m_bdiff(0,0) / p; } else { for (size_t k = 0; k < m_nsp; k++) { double sum2 = 0.0; for (size_t j = 0; j < m_nsp; j++) { if (j != k) { sum2 += m_molefracs[j] / m_bdiff(j,k); } } if (sum2 <= 0.0) { d[k] = m_bdiff(k,k) / p; } else { d[k] = (1 - m_molefracs[k]) / (p * sum2); } } } }
double ApproxMixTransport::viscosity() { update_T(); update_C(); if (m_visc_ok) { return m_viscmix; } doublereal vismix = 0.0; // update m_visc and m_phi if necessary if (!m_viscwt_ok) { updateViscosity_T(); } for (size_t ik=0; ik<_kMajor.size(); ik++) { size_t k = _kMajor[ik]; double sum = 0; for (size_t ij=0; ij<_kMajor.size(); ij++) { size_t j = _kMajor[ij]; sum += m_molefracs[j]*m_phi(k,j); } vismix += m_molefracs[k]*m_visc[k] / sum; } m_visc_ok = true; m_viscmix = vismix; assert(m_viscmix > 0 && m_viscmix < 1e100); return vismix; }
/* * The viscosity is computed using the general mixture rules * specified in the variable compositionDepType_. * * Solvent-only: * \f[ * \mu = \mu_0 * \f] * Mixture-average: * \f[ * \mu = \sum_k {\mu_k X_k} * \f] * * Here \f$ \mu_k \f$ is the viscosity of pure species \e k. * * @see updateViscosity_T(); */ doublereal SimpleTransport::viscosity() { update_T(); update_C(); if (m_visc_mix_ok) { return m_viscmix; } // update m_viscSpecies[] if necessary if (!m_visc_temp_ok) { updateViscosity_T(); } if (compositionDepType_ == 0) { m_viscmix = m_viscSpecies[0]; } else if (compositionDepType_ == 1) { m_viscmix = 0.0; for (size_t k = 0; k < m_nsp; k++) { m_viscmix += m_viscSpecies[k] * m_molefracs[k]; } } m_visc_mix_ok = true; return m_viscmix; }
doublereal LiquidTransport::thermalConductivity() { update_T(); update_C(); if (!m_lambda_mix_ok) { m_lambda = m_lambdaMixModel->getMixTransProp(m_lambdaTempDep_Ns); m_cond_mix_ok = true; } return m_lambda; }
doublereal LiquidTransport::ionConductivity() { update_T(); update_C(); if (m_ionCond_mix_ok) { return m_ionCondmix; } ////// LiquidTranInteraction method m_ionCondmix = m_ionCondMixModel->getMixTransProp(m_ionCondTempDep_Ns); return m_ionCondmix; }
doublereal LiquidTransport::viscosity() { update_T(); update_C(); if (m_visc_mix_ok) { return m_viscmix; } ////// LiquidTranInteraction method m_viscmix = m_viscMixModel->getMixTransProp(m_viscTempDep_Ns); return m_viscmix; }
/* * Returns the simple diffusion coefficients input into the model. Nothing fancy here. */ void SimpleTransport::getMixDiffCoeffs(doublereal* const d) { update_T(); update_C(); // update the binary diffusion coefficients if necessary if (!m_diff_temp_ok) { updateDiff_T(); } for (size_t k = 0; k < m_nsp; k++) { d[k] = m_diffSpecies[k]; } }
void LiquidTransport::selfDiffusion(doublereal* const selfDiff) { update_T(); update_C(); if (!m_selfDiff_mix_ok) { for (size_t k = 0; k < m_nsp; k++) { m_selfDiffMix[k] = m_selfDiffMixModel[k]->getMixTransProp(m_selfDiffTempDep_Ns[k]); } } for (size_t k = 0; k < m_nsp; k++) { selfDiff[k] = m_selfDiffMix[k]; } }
/* * The thermal conductivity is computed from the following mixture rule: * \f[ * \lambda = 0.5 \left( \sum_k X_k \lambda_k + \frac{1}{\sum_k X_k/\lambda_k} \right) * \f] * * It's used to compute the flux of energy due to a thermal gradient * * \f[ * j_T = - \lambda \nabla T * \f] * * The flux of energy has units of energy (kg m2 /s2) per second per area. * * The units of lambda are W / m K which is equivalent to kg m / s^3 K. * * @return Returns the mixture thermal conductivity, with units of W/m/K */ doublereal MixTransport::thermalConductivity() { int k; update_T(); update_C(); if (!m_spcond_ok) updateCond_T(); if (!m_condmix_ok) { doublereal sum1 = 0.0, sum2 = 0.0; for (k = 0; k < m_nsp; k++) { sum1 += m_molefracs[k] * m_cond[k]; sum2 += m_molefracs[k] / m_cond[k]; } m_lambda = 0.5*(sum1 + 1.0/sum2); m_condmix_ok = true; } return m_lambda; }
void Plog::validate(const std::string& equation) { double T[] = {200.0, 500.0, 1000.0, 2000.0, 5000.0, 10000.0}; for (auto iter = pressures_.begin(); iter->first < 1000; iter++) { update_C(&iter->first); for (size_t i=0; i < 6; i++) { double k = updateRC(log(T[i]), 1.0/T[i]); if (!(k >= 0)) { // k is NaN. Increment the iterator so that the error // message will correctly indicate that the problematic rate // expression is at the higher of the adjacent pressures. throw CanteraError("Plog::validate", "Invalid rate coefficient for reaction '{}'\nat P = {}, T = {}", equation, std::exp((++iter)->first), T[i]); } } } }
doublereal AqueousTransport::thermalConductivity() { update_T(); update_C(); if (!m_spcond_ok) { updateCond_T(); } if (!m_condmix_ok) { doublereal sum1 = 0.0, sum2 = 0.0; for (size_t k = 0; k < m_nsp; k++) { sum1 += m_molefracs[k] * m_cond[k]; sum2 += m_molefracs[k] / m_cond[k]; } m_lambda = 0.5*(sum1 + 1.0/sum2); } return m_lambda; }
void ApproxMixTransport::getMixDiffCoeffs(double* const d) { update_T(); update_C(); // update the binary diffusion coefficients if necessary if (!m_bindiff_ok) { updateDiff_T(); } double mmw = m_thermo->meanMolecularWeight(); double sumxw = 0.0, sum2; double p = m_thermo->pressure(); if (m_nsp == 1) { d[0] = m_bdiff(0,0) / p; } else { for (size_t k = 0; k < m_nsp; k++) { sumxw += m_molefracs[k] * m_mw[k]; } for (size_t k = 0; k < m_nsp; k++) { sum2 = 0.0; if (m_molefracs[k] >= _threshold) { for (size_t j = 0; j < m_nsp; j++) { if (j != k) { sum2 += m_molefracs[j] / m_bdiff(j,k); } } } else { for (size_t j : _kMajor) { if (j != k) { sum2 += m_molefracs[j] / m_bdiff(j,k); } } } if (sum2 <= 0.0) { d[k] = m_bdiff(k,k) / p; } else { d[k] = (sumxw - m_molefracs[k] * m_mw[k])/(p * mmw * sum2); } } } }
void LiquidTransport::mobilityRatio(doublereal* mobRat) { update_T(); update_C(); // LiquidTranInteraction method if (!m_mobRat_mix_ok) { for (size_t k = 0; k < m_nsp2; k++) { if (m_mobRatMixModel[k]) { m_mobRatMix[k] = m_mobRatMixModel[k]->getMixTransProp(m_mobRatTempDep_Ns[k]); if (m_mobRatMix[k] > 0.0) { m_mobRatMix[k / m_nsp + m_nsp * (k % m_nsp)] = 1.0 / m_mobRatMix[k]; // Also must be off diagonal: k%(1+n)!=0, but then m_mobRatMixModel[k] shouldn't be initialized anyway } } } } for (size_t k = 0; k < m_nsp2; k++) { mobRat[k] = m_mobRatMix[k]; } }
/* * * The viscosity is computed using the Wilke mixture rule. * * \f[ * \mu = \sum_k \frac{\mu_k X_k}{\sum_j \Phi_{k,j} X_j}. * \f] * * Here \f$ \mu_k \f$ is the viscosity of pure species \e k, and * * \f[ * \Phi_{k,j} = \frac{\left[1 * + \sqrt{\left(\frac{\mu_k}{\mu_j}\sqrt{\frac{M_j}{M_k}}\right)}\right]^2} * {\sqrt{8}\sqrt{1 + M_k/M_j}} * \f] * * @see updateViscosity_T(); */ doublereal MixTransport::viscosity() { update_T(); update_C(); if (m_viscmix_ok) return m_viscmix; doublereal vismix = 0.0; int k; // update m_visc and m_phi if necessary if (!m_viscwt_ok) updateViscosity_T(); multiply(m_phi, DATA_PTR(m_molefracs), DATA_PTR(m_spwork)); for (k = 0; k < m_nsp; k++) { vismix += m_molefracs[k] * m_visc[k]/m_spwork[k]; //denom; } m_viscmix = vismix; return vismix; }
/* * The thermal is computed using the general mixture rules * specified in the variable compositionDepType_. * * Solvent-only: * \f[ * \lambda = \lambda_0 * \f] * Mixture-average: * \f[ * \lambda = \sum_k {\lambda_k X_k} * \f] * * Here \f$ \lambda_k \f$ is the thermal conductivity of pure species \e k. * * @see updateCond_T(); */ doublereal SimpleTransport::thermalConductivity() { update_T(); update_C(); if (!m_cond_temp_ok) { updateCond_T(); } if (!m_cond_mix_ok) { if (compositionDepType_ == 0) { m_lambda = m_condSpecies[0]; } else if (compositionDepType_ == 1) { m_lambda = 0.0; for (size_t k = 0; k < m_nsp; k++) { m_lambda += m_condSpecies[k] * m_molefracs[k]; } } m_cond_mix_ok = true; } return m_lambda; }
void Plog::validate(const ReactionData& rdata) { double T[] = {200.0, 500.0, 1000.0, 2000.0, 5000.0, 10000.0}; for (pressureIter iter = pressures_.begin(); iter->first < 1000; iter++) { update_C(&iter->first); for (size_t i=0; i < 6; i++) { double k = updateRC(log(T[i]), 1.0/T[i]); if (!(k >= 0)) { // k is NaN. Increment the iterator so that the error // message will correctly indicate that the problematic rate // expression is at the higher of the adjacent pressures. throw CanteraError("Plog::validate", "Invalid rate coefficient for reaction #" + int2str(rdata.number) + ":\n" + rdata.equation + "\n" + "at P = " + fp2str(std::exp((++iter)->first)) + ", T = " + fp2str(T[i])); } } } }
doublereal AqueousTransport::viscosity() { update_T(); update_C(); if (m_viscmix_ok) { return m_viscmix; } // update m_visc[] and m_phi[] if necessary if (!m_viscwt_ok) { updateViscosity_T(); } multiply(m_phi, DATA_PTR(m_molefracs), DATA_PTR(m_spwork)); m_viscmix = 0.0; for (size_t k = 0; k < m_nsp; k++) { m_viscmix += m_molefracs[k] * m_visc[k]/m_spwork[k]; //denom; } return m_viscmix; }
/* * * units = kg/m2/s * * Internally, gradients in the in mole fraction, temperature * and electrostatic potential contribute to the diffusive flux * * * The diffusive mass flux of species \e k is computed from the following * formula * * \f[ * j_k = - M_k z_k u^f_k F c_k \nabla \Psi - c M_k D_k \nabla X_k - Y_k V_c * \f] * * where V_c is the correction velocity * * \f[ * V_c = - \sum_j {M_k z_k u^f_k F c_k \nabla \Psi + c M_j D_j \nabla X_j} * \f] * * @param ldf stride of the fluxes array. Must be equal to * or greater than the number of species. * @param fluxes Vector of calculated fluxes */ void SimpleTransport::getSpeciesFluxesExt(size_t ldf, doublereal* fluxes) { AssertThrow(ldf >= m_nsp ,"SimpleTransport::getSpeciesFluxesExt: Stride must be greater than m_nsp"); update_T(); update_C(); getMixDiffCoeffs(DATA_PTR(m_spwork)); const vector_fp& mw = m_thermo->molecularWeights(); const doublereal* y = m_thermo->massFractions(); doublereal concTotal = m_thermo->molarDensity(); // Unroll wrt ndim if (doMigration_) { double FRT = ElectronCharge / (Boltzmann * m_temp); for (size_t n = 0; n < m_nDim; n++) { rhoVc[n] = 0.0; for (size_t k = 0; k < m_nsp; k++) { fluxes[n*ldf + k] = - concTotal * mw[k] * m_spwork[k] * (m_Grad_X[n*m_nsp + k] + FRT * m_molefracs[k] * m_chargeSpecies[k] * m_Grad_V[n]); rhoVc[n] += fluxes[n*ldf + k]; } } } else { for (size_t n = 0; n < m_nDim; n++) { rhoVc[n] = 0.0; for (size_t k = 0; k < m_nsp; k++) { fluxes[n*ldf + k] = - concTotal * mw[k] * m_spwork[k] * m_Grad_X[n*m_nsp + k]; rhoVc[n] += fluxes[n*ldf + k]; } } } if (m_velocityBasis == VB_MASSAVG) { for (size_t n = 0; n < m_nDim; n++) { rhoVc[n] = 0.0; for (size_t k = 0; k < m_nsp; k++) { rhoVc[n] += fluxes[n*ldf + k]; } } for (size_t n = 0; n < m_nDim; n++) { for (size_t k = 0; k < m_nsp; k++) { fluxes[n*ldf + k] -= y[k] * rhoVc[n]; } } } else if (m_velocityBasis == VB_MOLEAVG) { for (size_t n = 0; n < m_nDim; n++) { rhoVc[n] = 0.0; for (size_t k = 0; k < m_nsp; k++) { rhoVc[n] += fluxes[n*ldf + k] / mw[k]; } } for (size_t n = 0; n < m_nDim; n++) { for (size_t k = 0; k < m_nsp; k++) { fluxes[n*ldf + k] -= m_molefracs[k] * rhoVc[n] * mw[k]; } } } else if (m_velocityBasis >= 0) { for (size_t n = 0; n < m_nDim; n++) { rhoVc[n] = - fluxes[n*ldf + m_velocityBasis] / mw[m_velocityBasis]; for (size_t k = 0; k < m_nsp; k++) { rhoVc[n] += fluxes[n*ldf + k] / mw[k]; } } for (size_t n = 0; n < m_nDim; n++) { for (size_t k = 0; k < m_nsp; k++) { fluxes[n*ldf + k] -= m_molefracs[k] * rhoVc[n] * mw[k]; } fluxes[n*ldf + m_velocityBasis] = 0.0; } } else { throw CanteraError("SimpleTransport::getSpeciesFluxesExt()", "unknown velocity basis"); } }
void LiquidTransport::stefan_maxwell_solve() { doublereal tmp; m_B.resize(m_nsp, m_nDim, 0.0); m_A.resize(m_nsp, m_nsp, 0.0); //! grab a local copy of the molecular weights const vector_fp& M = m_thermo->molecularWeights(); //! grad a local copy of the ion molar volume (inverse total ion concentration) const doublereal vol = m_thermo->molarVolume(); /* * Update the temperature, concentrations and diffusion coefficients in the mixture. */ update_T(); update_C(); if (!m_diff_temp_ok) { updateDiff_T(); } double T = m_thermo->temperature(); update_Grad_lnAC(); m_thermo->getActivityCoefficients(DATA_PTR(m_actCoeff)); /* * Calculate the electrochemical potential gradient. This is the * driving force for relative diffusional transport. * * Here we calculate * * X_i * (grad (mu_i) + S_i grad T - M_i / dens * grad P * * This is Eqn. 13-1 p. 318 Newman. The original equation is from * Hershfeld, Curtis, and Bird. * * S_i is the partial molar entropy of species i. This term will cancel * out a lot of the grad T terms in grad (mu_i), therefore simplifying * the expression. * * Ok I think there may be many ways to do this. One way is to do it via basis * functions, at the nodes, as a function of the variables in the problem. * * For calculation of molality based thermo systems, we current get * the molar based values. This may change. * * Note, we have broken the symmetry of the matrix here, due to * considerations involving species concentrations going to zero. */ for (size_t a = 0; a < m_nDim; a++) { for (size_t i = 0; i < m_nsp; i++) { m_Grad_mu[a*m_nsp + i] = m_chargeSpecies[i] * Faraday * m_Grad_V[a] + GasConstant * T * m_Grad_lnAC[a*m_nsp+i]; } } if (m_thermo->activityConvention() == cAC_CONVENTION_MOLALITY) { int iSolvent = 0; double mwSolvent = m_thermo->molecularWeight(iSolvent); double mnaught = mwSolvent/ 1000.; double lnmnaught = log(mnaught); for (size_t a = 0; a < m_nDim; a++) { for (size_t i = 1; i < m_nsp; i++) { m_Grad_mu[a*m_nsp + i] -= m_molefracs[i] * GasConstant * m_Grad_T[a] * lnmnaught; } } } /* * Just for Note, m_A(i,j) refers to the ith row and jth column. * They are still fortran ordered, so that i varies fastest. */ double condSum1; const doublereal invRT = 1.0 / (GasConstant * T); switch (m_nDim) { case 1: /* 1-D approximation */ m_B(0,0) = 0.0; //equation for the reference velocity for (size_t j = 0; j < m_nsp; j++) { if (m_velocityBasis == VB_MOLEAVG) { m_A(0,j) = m_molefracs_tran[j]; } else if (m_velocityBasis == VB_MASSAVG) { m_A(0,j) = m_massfracs_tran[j]; } else if ((m_velocityBasis >= 0) && (m_velocityBasis < static_cast<int>(m_nsp))) { // use species number m_velocityBasis as reference velocity if (m_velocityBasis == static_cast<int>(j)) { m_A(0,j) = 1.0; } else { m_A(0,j) = 0.0; } } else { throw CanteraError("LiquidTransport::stefan_maxwell_solve", "Unknown reference velocity provided."); } } for (size_t i = 1; i < m_nsp; i++) { m_B(i,0) = m_Grad_mu[i] * invRT; m_A(i,i) = 0.0; for (size_t j = 0; j < m_nsp; j++) { if (j != i) { tmp = m_molefracs_tran[j] * m_bdiff(i,j); m_A(i,i) -= tmp; m_A(i,j) = tmp; } } } //! invert and solve the system Ax = b. Answer is in m_B solve(m_A, m_B); condSum1 = 0; for (size_t i = 0; i < m_nsp; i++) { condSum1 -= Faraday*m_chargeSpecies[i]*m_B(i,0)*m_molefracs_tran[i]/vol; } break; case 2: /* 2-D approximation */ m_B(0,0) = 0.0; m_B(0,1) = 0.0; //equation for the reference velocity for (size_t j = 0; j < m_nsp; j++) { if (m_velocityBasis == VB_MOLEAVG) { m_A(0,j) = m_molefracs_tran[j]; } else if (m_velocityBasis == VB_MASSAVG) { m_A(0,j) = m_massfracs_tran[j]; } else if ((m_velocityBasis >= 0) && (m_velocityBasis < static_cast<int>(m_nsp))) { // use species number m_velocityBasis as reference velocity if (m_velocityBasis == static_cast<int>(j)) { m_A(0,j) = 1.0; } else { m_A(0,j) = 0.0; } } else { throw CanteraError("LiquidTransport::stefan_maxwell_solve", "Unknown reference velocity provided."); } } for (size_t i = 1; i < m_nsp; i++) { m_B(i,0) = m_Grad_mu[i] * invRT; m_B(i,1) = m_Grad_mu[m_nsp + i] * invRT; m_A(i,i) = 0.0; for (size_t j = 0; j < m_nsp; j++) { if (j != i) { tmp = m_molefracs_tran[j] * m_bdiff(i,j); m_A(i,i) -= tmp; m_A(i,j) = tmp; } } } //! invert and solve the system Ax = b. Answer is in m_B solve(m_A, m_B); break; case 3: /* 3-D approximation */ m_B(0,0) = 0.0; m_B(0,1) = 0.0; m_B(0,2) = 0.0; //equation for the reference velocity for (size_t j = 0; j < m_nsp; j++) { if (m_velocityBasis == VB_MOLEAVG) { m_A(0,j) = m_molefracs_tran[j]; } else if (m_velocityBasis == VB_MASSAVG) { m_A(0,j) = m_massfracs_tran[j]; } else if ((m_velocityBasis >= 0) && (m_velocityBasis < static_cast<int>(m_nsp))) { // use species number m_velocityBasis as reference velocity if (m_velocityBasis == static_cast<int>(j)) { m_A(0,j) = 1.0; } else { m_A(0,j) = 0.0; } } else { throw CanteraError("LiquidTransport::stefan_maxwell_solve", "Unknown reference velocity provided."); } } for (size_t i = 1; i < m_nsp; i++) { m_B(i,0) = m_Grad_mu[i] * invRT; m_B(i,1) = m_Grad_mu[m_nsp + i] * invRT; m_B(i,2) = m_Grad_mu[2*m_nsp + i] * invRT; m_A(i,i) = 0.0; for (size_t j = 0; j < m_nsp; j++) { if (j != i) { tmp = m_molefracs_tran[j] * m_bdiff(i,j); m_A(i,i) -= tmp; m_A(i,j) = tmp; } } } //! invert and solve the system Ax = b. Answer is in m_B solve(m_A, m_B); break; default: printf("unimplemented\n"); throw CanteraError("routine", "not done"); break; } for (size_t a = 0; a < m_nDim; a++) { for (size_t j = 0; j < m_nsp; j++) { m_Vdiff(j,a) = m_B(j,a); m_flux(j,a) = concTot_ * M[j] * m_molefracs_tran[j] * m_B(j,a); } } }