Cantera 2.6.0
BasisOptimize.cpp
Go to the documentation of this file.
1/**
2 * @file BasisOptimize.cpp Functions which calculation optimized basis of the
3 * stoichiometric coefficient matrix (see /ref equil functions)
4 */
5
6// This file is part of Cantera. See License.txt in the top-level directory or
7// at https://cantera.org/license.txt for license and copyright information.
8
11#include "cantera/base/global.h"
12
13using namespace std;
14
15namespace Cantera
16{
18static const double USEDBEFORE = -1;
19
20size_t BasisOptimize(int* usedZeroedSpecies, bool doFormRxn, MultiPhase* mphase,
21 std::vector<size_t>& orderVectorSpecies,
22 std::vector<size_t>& orderVectorElements,
23 vector_fp& formRxnMatrix)
24{
25 // Get the total number of elements defined in the multiphase object
26 size_t ne = mphase->nElements();
27
28 // Get the total number of species in the multiphase object
29 size_t nspecies = mphase->nSpecies();
30
31 // Perhaps, initialize the element ordering
32 if (orderVectorElements.size() < ne) {
33 orderVectorElements.resize(ne);
34 iota(orderVectorElements.begin(), orderVectorElements.end(), 0);
35 }
36
37 // Perhaps, initialize the species ordering
38 if (orderVectorSpecies.size() != nspecies) {
39 orderVectorSpecies.resize(nspecies);
40 iota(orderVectorSpecies.begin(), orderVectorSpecies.end(), 0);
41 }
42
43 if (BasisOptimize_print_lvl >= 1) {
44 writelog(" ");
45 writeline('-', 77);
46 writelog(" --- Subroutine BASOPT called to ");
47 writelog("calculate the number of components and ");
48 writelog("evaluate the formation matrix\n");
50 writelog(" ---\n");
51
52 writelog(" --- Formula Matrix used in BASOPT calculation\n");
53 writelog(" --- Species | Order | ");
54 for (size_t j = 0; j < ne; j++) {
55 size_t jj = orderVectorElements[j];
56 writelog(" {:>4.4s}({:1d})", mphase->elementName(jj), j);
57 }
58 writelog("\n");
59 for (size_t k = 0; k < nspecies; k++) {
60 size_t kk = orderVectorSpecies[k];
61 writelog(" --- {:>11.11s} | {:4d} |",
62 mphase->speciesName(kk), k);
63 for (size_t j = 0; j < ne; j++) {
64 size_t jj = orderVectorElements[j];
65 double num = mphase->nAtoms(kk,jj);
66 writelogf("%6.1g ", num);
67 }
68 writelog("\n");
69 }
70 writelog(" --- \n");
71 }
72 }
73
74 // Calculate the maximum value of the number of components possible. It's
75 // equal to the minimum of the number of elements and the number of total
76 // species.
77 size_t nComponents = std::min(ne, nspecies);
78 size_t nNonComponents = nspecies - nComponents;
79
80 // Set this return variable to false
81 *usedZeroedSpecies = false;
82
83 // Create an array of mole numbers
84 vector_fp molNum(nspecies,0.0);
85 mphase->getMoles(molNum.data());
86
87 // Other workspace
88 DenseMatrix sm(ne, ne);
89 vector_fp ss(ne, 0.0);
90 vector_fp sa(ne, 0.0);
91 if (formRxnMatrix.size() < nspecies*ne) {
92 formRxnMatrix.resize(nspecies*ne, 0.0);
93 }
94
95 // For debugging purposes keep an unmodified copy of the array.
96 vector_fp molNumBase = molNum;
97 double molSave = 0.0;
98 size_t jr = 0;
99
100 // Top of a loop of some sort based on the index JR. JR is the current
101 // number of component species found.
102 while (jr < nComponents) {
103 // Top of another loop point based on finding a linearly independent
104 // species
105 size_t k = npos;
106 while (true) {
107 // Search the remaining part of the mole number vector, molNum for
108 // the largest remaining species. Return its identity. kk is the raw
109 // number. k is the orderVectorSpecies index.
110 size_t kk = max_element(molNum.begin(), molNum.end()) - molNum.begin();
111 size_t j;
112 for (j = 0; j < nspecies; j++) {
113 if (orderVectorSpecies[j] == kk) {
114 k = j;
115 break;
116 }
117 }
118 if (j == nspecies) {
119 throw CanteraError("BasisOptimize", "orderVectorSpecies contains an error");
120 }
121
122 if (molNum[kk] == 0.0) {
123 *usedZeroedSpecies = true;
124 }
125 // If the largest molNum is negative, then we are done.
126 if (molNum[kk] == USEDBEFORE) {
127 nComponents = jr;
128 nNonComponents = nspecies - nComponents;
129 break;
130 }
131
132 // Assign a small negative number to the component that we have
133 // just found, in order to take it out of further consideration.
134 molSave = molNum[kk];
135 molNum[kk] = USEDBEFORE;
136
137 // CHECK LINEAR INDEPENDENCE WITH PREVIOUS SPECIES
138
139 // Modified Gram-Schmidt Method, p. 202 Dalquist
140 // QR factorization of a matrix without row pivoting.
141 size_t jl = jr;
142 for (j = 0; j < ne; ++j) {
143 size_t jj = orderVectorElements[j];
144 sm(j, jr) = mphase->nAtoms(kk,jj);
145 }
146 if (jl > 0) {
147 // Compute the coefficients of JA column of the the upper
148 // triangular R matrix, SS(J) = R_J_JR (this is slightly
149 // different than Dalquist) R_JA_JA = 1
150 for (j = 0; j < jl; ++j) {
151 ss[j] = 0.0;
152 for (size_t i = 0; i < ne; ++i) {
153 ss[j] += sm(i, jr) * sm(i, j);
154 }
155 ss[j] /= sa[j];
156 }
157
158 // Now make the new column, (*,JR), orthogonal to the previous
159 // columns
160 for (j = 0; j < jl; ++j) {
161 for (size_t i = 0; i < ne; ++i) {
162 sm(i, jr) -= ss[j] * sm(i, j);
163 }
164 }
165 }
166
167 // Find the new length of the new column in Q.
168 // It will be used in the denominator in future row calcs.
169 sa[jr] = 0.0;
170 for (size_t ml = 0; ml < ne; ++ml) {
171 sa[jr] += pow(sm(ml, jr), 2);
172 }
173
174 // IF NORM OF NEW ROW .LT. 1E-3 REJECT
175 if (sa[jr] > 1.0e-6) {
176 break;
177 }
178 }
179
180 // REARRANGE THE DATA
181 if (jr != k) {
182 if (BasisOptimize_print_lvl >= 1) {
183 size_t kk = orderVectorSpecies[k];
184 writelogf(" --- %-12.12s", mphase->speciesName(kk));
185 size_t jj = orderVectorSpecies[jr];
186 writelogf("(%9.2g) replaces %-12.12s",
187 molSave, mphase->speciesName(jj));
188 writelogf("(%9.2g) as component %3d\n", molNum[jj], jr);
189 }
190 std::swap(orderVectorSpecies[jr], orderVectorSpecies[k]);
191 }
192
193 // If we haven't found enough components, go back and find some more
194 jr++;
195 }
196
197 if (! doFormRxn) {
198 return nComponents;
199 }
200
201 // EVALUATE THE STOICHIOMETRY
202 //
203 // Formulate the matrix problem for the stoichiometric
204 // coefficients. CX + B = 0
205 //
206 // C will be an nc x nc matrix made up of the formula vectors for the
207 // components. Each component's formula vector is a column. The rows are the
208 // elements.
209 //
210 // n RHS's will be solved for. Thus, B is an nc x n matrix.
211 //
212 // BIG PROBLEM 1/21/99:
213 //
214 // This algorithm makes the assumption that the first nc rows of the formula
215 // matrix aren't rank deficient. However, this might not be the case. For
216 // example, assume that the first element in FormulaMatrix[] is argon.
217 // Assume that no species in the matrix problem actually includes argon.
218 // Then, the first row in sm[], below will be identically zero. bleh.
219 //
220 // What needs to be done is to perform a rearrangement of the ELEMENTS ->
221 // that is, rearrange, FormulaMatrix, sp, and gai, such that the first nc
222 // elements form in combination with the nc components create an invertible
223 // sm[]. not a small project, but very doable.
224 //
225 // An alternative would be to turn the matrix problem below into an ne x nc
226 // problem, and do QR elimination instead of Gauss-Jordan elimination.
227 //
228 // Note the rearrangement of elements need only be done once in the problem.
229 // It's actually very similar to the top of this program with ne being the
230 // species and nc being the elements!!
231
232 sm.resize(nComponents, nComponents);
233 for (size_t k = 0; k < nComponents; ++k) {
234 size_t kk = orderVectorSpecies[k];
235 for (size_t j = 0; j < nComponents; ++j) {
236 size_t jj = orderVectorElements[j];
237 sm(j, k) = mphase->nAtoms(kk, jj);
238 }
239 }
240
241 for (size_t i = 0; i < nNonComponents; ++i) {
242 size_t k = nComponents + i;
243 size_t kk = orderVectorSpecies[k];
244 for (size_t j = 0; j < nComponents; ++j) {
245 size_t jj = orderVectorElements[j];
246 formRxnMatrix[j + i * ne] = - mphase->nAtoms(kk, jj);
247 }
248 }
249
250 // // Use LU factorization to calculate the reaction matrix
251 solve(sm, formRxnMatrix.data(), nNonComponents, ne);
252
253 if (BasisOptimize_print_lvl >= 1) {
254 writelog(" ---\n");
255 writelogf(" --- Number of Components = %d\n", nComponents);
256 writelog(" --- Formula Matrix:\n");
257 writelog(" --- Components: ");
258 for (size_t k = 0; k < nComponents; k++) {
259 size_t kk = orderVectorSpecies[k];
260 writelogf(" %3d (%3d) ", k, kk);
261 }
262 writelog("\n --- Components Moles: ");
263 for (size_t k = 0; k < nComponents; k++) {
264 size_t kk = orderVectorSpecies[k];
265 writelogf("%-11.3g", molNumBase[kk]);
266 }
267 writelog("\n --- NonComponent | Moles | ");
268 for (size_t i = 0; i < nComponents; i++) {
269 size_t kk = orderVectorSpecies[i];
270 writelogf("%-11.10s", mphase->speciesName(kk));
271 }
272 writelog("\n");
273
274 for (size_t i = 0; i < nNonComponents; i++) {
275 size_t k = i + nComponents;
276 size_t kk = orderVectorSpecies[k];
277 writelogf(" --- %3d (%3d) ", k, kk);
278 writelogf("%-10.10s", mphase->speciesName(kk));
279 writelogf("|%10.3g|", molNumBase[kk]);
280 // Print the negative of formRxnMatrix[]; it's easier to interpret.
281 for (size_t j = 0; j < nComponents; j++) {
282 writelogf(" %6.2f", - formRxnMatrix[j + i * ne]);
283 }
284 writelog("\n");
285 }
286 writelog(" ");
287 writeline('-', 77);
288 }
289
290 return nComponents;
291} // basopt()
292
293
294void ElemRearrange(size_t nComponents, const vector_fp& elementAbundances,
295 MultiPhase* mphase,
296 std::vector<size_t>& orderVectorSpecies,
297 std::vector<size_t>& orderVectorElements)
298{
299 size_t nelements = mphase->nElements();
300 // Get the total number of species in the multiphase object
301 size_t nspecies = mphase->nSpecies();
302
303 if (BasisOptimize_print_lvl > 0) {
304 writelog(" ");
305 writeline('-', 77);
306 writelog(" --- Subroutine ElemRearrange() called to ");
307 writelog("check stoich. coefficient matrix\n");
308 writelog(" --- and to rearrange the element ordering once\n");
309 }
310
311 // Perhaps, initialize the element ordering
312 if (orderVectorElements.size() < nelements) {
313 orderVectorElements.resize(nelements);
314 for (size_t j = 0; j < nelements; j++) {
315 orderVectorElements[j] = j;
316 }
317 }
318
319 // Perhaps, initialize the species ordering. However, this is dangerous, as
320 // this ordering is assumed to yield the component species for the problem
321 if (orderVectorSpecies.size() != nspecies) {
322 orderVectorSpecies.resize(nspecies);
323 for (size_t k = 0; k < nspecies; k++) {
324 orderVectorSpecies[k] = k;
325 }
326 }
327
328 // If the elementAbundances aren't input, just create a fake one based on
329 // summing the column of the stoich matrix. This will force elements with
330 // zero species to the end of the element ordering.
331 vector_fp eAbund(nelements,0.0);
332 if (elementAbundances.size() != nelements) {
333 for (size_t j = 0; j < nelements; j++) {
334 eAbund[j] = 0.0;
335 for (size_t k = 0; k < nspecies; k++) {
336 eAbund[j] += fabs(mphase->nAtoms(k, j));
337 }
338 }
339 } else {
340 copy(elementAbundances.begin(), elementAbundances.end(),
341 eAbund.begin());
342 }
343
344 vector_fp sa(nelements,0.0);
345 vector_fp ss(nelements,0.0);
346 vector_fp sm(nelements*nelements,0.0);
347
348 // Top of a loop of some sort based on the index JR. JR is the current
349 // number independent elements found.
350 size_t jr = 0;
351 while (jr < nComponents) {
352 // Top of another loop point based on finding a linearly independent
353 // element
354 size_t k = nelements;
355 while (true) {
356 // Search the element vector. We first locate elements that are
357 // present in any amount. Then, we locate elements that are not
358 // present in any amount. Return its identity in K.
359 k = nelements;
360 size_t kk;
361 for (size_t ielem = jr; ielem < nelements; ielem++) {
362 kk = orderVectorElements[ielem];
363 if (eAbund[kk] != USEDBEFORE && eAbund[kk] > 0.0) {
364 k = ielem;
365 break;
366 }
367 }
368 for (size_t ielem = jr; ielem < nelements; ielem++) {
369 kk = orderVectorElements[ielem];
370 if (eAbund[kk] != USEDBEFORE) {
371 k = ielem;
372 break;
373 }
374 }
375
376 if (k == nelements) {
377 // When we are here, there is an error usually.
378 // We haven't found the number of elements necessary.
379 if (BasisOptimize_print_lvl > 0) {
380 writelogf("Error exit: returning with nComponents = %d\n", jr);
381 }
382 throw CanteraError("ElemRearrange", "Required number of elements not found.");
383 }
384
385 // Assign a large negative number to the element that we have
386 // just found, in order to take it out of further consideration.
387 eAbund[kk] = USEDBEFORE;
388
389 // CHECK LINEAR INDEPENDENCE OF CURRENT FORMULA MATRIX
390 // LINE WITH PREVIOUS LINES OF THE FORMULA MATRIX
391
392 // Modified Gram-Schmidt Method, p. 202 Dalquist
393 // QR factorization of a matrix without row pivoting.
394 size_t jl = jr;
395
396 // Fill in the row for the current element, k, under consideration
397 // The row will contain the Formula matrix value for that element
398 // with respect to the vector of component species. (note j and k
399 // indices are flipped compared to the previous routine)
400 for (size_t j = 0; j < nComponents; ++j) {
401 size_t jj = orderVectorSpecies[j];
402 kk = orderVectorElements[k];
403 sm[j + jr*nComponents] = mphase->nAtoms(jj,kk);
404 }
405 if (jl > 0) {
406 // Compute the coefficients of JA column of the the upper
407 // triangular R matrix, SS(J) = R_J_JR (this is slightly
408 // different than Dalquist) R_JA_JA = 1
409 for (size_t j = 0; j < jl; ++j) {
410 ss[j] = 0.0;
411 for (size_t i = 0; i < nComponents; ++i) {
412 ss[j] += sm[i + jr*nComponents] * sm[i + j*nComponents];
413 }
414 ss[j] /= sa[j];
415 }
416
417 // Now make the new column, (*,JR), orthogonal to the
418 // previous columns
419 for (size_t j = 0; j < jl; ++j) {
420 for (size_t i = 0; i < nComponents; ++i) {
421 sm[i + jr*nComponents] -= ss[j] * sm[i + j*nComponents];
422 }
423 }
424 }
425
426 // Find the new length of the new column in Q.
427 // It will be used in the denominator in future row calcs.
428 sa[jr] = 0.0;
429 for (size_t ml = 0; ml < nComponents; ++ml) {
430 double tmp = sm[ml + jr*nComponents];
431 sa[jr] += tmp * tmp;
432 }
433 // IF NORM OF NEW ROW .LT. 1E-6 REJECT
434 if (sa[jr] > 1.0e-6) {
435 break;
436 }
437 }
438 // REARRANGE THE DATA
439 if (jr != k) {
440 if (BasisOptimize_print_lvl > 0) {
441 size_t kk = orderVectorElements[k];
442 writelog(" --- ");
443 writelogf("%-2.2s", mphase->elementName(kk));
444 writelog("replaces ");
445 kk = orderVectorElements[jr];
446 writelogf("%-2.2s", mphase->elementName(kk));
447 writelogf(" as element %3d\n", jr);
448 }
449 std::swap(orderVectorElements[jr], orderVectorElements[k]);
450 }
451
452 // If we haven't found enough components, go back and find some more
453 jr++;
454 };
455}
456
457}
Headers for the MultiPhase object that is used to set up multiphase equilibrium problems (see Classes...
vector_fp & data()
Return a reference to the data vector.
Definition: Array.h:218
Base class for exceptions thrown by Cantera classes.
Definition: ctexceptions.h:61
A class for full (non-sparse) matrices with Fortran-compatible data storage, which adds matrix operat...
Definition: DenseMatrix.h:51
void resize(size_t n, size_t m, doublereal v=0.0)
Resize the matrix.
Definition: DenseMatrix.cpp:70
A class for multiphase mixtures.
Definition: MultiPhase.h:60
void getMoles(doublereal *molNum) const
Get the mole numbers of all species in the multiphase object.
Definition: MultiPhase.cpp:372
size_t nSpecies() const
Number of species, summed over all phases.
Definition: MultiPhase.h:126
std::string speciesName(const size_t kGlob) const
Name of species with global index kGlob.
Definition: MultiPhase.cpp:761
std::string elementName(size_t m) const
Returns the name of the global element m.
Definition: MultiPhase.cpp:732
doublereal nAtoms(const size_t kGlob, const size_t mGlob) const
Returns the Number of atoms of global element mGlob in global species kGlob.
Definition: MultiPhase.cpp:766
size_t nElements() const
Number of elements.
Definition: MultiPhase.h:100
This file contains definitions for utility functions and text for modules, inputfiles,...
size_t BasisOptimize(int *usedZeroedSpecies, bool doFormRxn, MultiPhase *mphase, std::vector< size_t > &orderVectorSpecies, std::vector< size_t > &orderVectorElements, vector_fp &formRxnMatrix)
Choose the optimum basis of species for the equilibrium calculations.
void ElemRearrange(size_t nComponents, const vector_fp &elementAbundances, MultiPhase *mphase, std::vector< size_t > &orderVectorSpecies, std::vector< size_t > &orderVectorElements)
Handles the potential rearrangement of the constraint equations represented by the Formula Matrix.
void writelog(const std::string &fmt, const Args &... args)
Write a formatted message to the screen.
Definition: global.h:164
void writelogf(const char *fmt, const Args &... args)
Write a formatted message to the screen.
Definition: global.h:185
Namespace for the Cantera kernel.
Definition: AnyMap.h:29
const size_t npos
index returned by functions to indicate "no position"
Definition: ct_defs.h:192
std::vector< double > vector_fp
Turn on the use of stl vectors for the basic array type within cantera Vector of doubles.
Definition: ct_defs.h:184
int BasisOptimize_print_lvl
External int that is used to turn on debug printing for the BasisOptimize program.
int solve(DenseMatrix &A, double *b, size_t nrhs=1, size_t ldb=0)
Solve Ax = b. Array b is overwritten on exit with x.
Various templated functions that carry out common vector operations (see Templated Utility Functions)...