Reaktoro
A unified framework for modeling chemically reactive systems
ChemicalVector.hpp
1 // Reaktoro is a unified framework for modeling chemically reactive systems.
2 //
3 // Copyright (C) 2014-2018 Allan Leal
4 //
5 // This library is free software; you can redistribute it and/or
6 // modify it under the terms of the GNU Lesser General Public
7 // License as published by the Free Software Foundation; either
8 // version 2.1 of the License, or (at your option) any later version.
9 //
10 // This library is distributed in the hope that it will be useful,
11 // but WITHOUT ANY WARRANTY; without even the implied warranty of
12 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
13 // Lesser General Public License for more details.
14 //
15 // You should have received a copy of the GNU Lesser General Public License
16 // along with this library. If not, see <http://www.gnu.org/licenses/>.
17 
18 #pragma once
19 
20 // Reaktoro includes
21 #include <Reaktoro/Common/ChemicalScalar.hpp>
22 #include <Reaktoro/Common/ThermoScalar.hpp>
23 #include <Reaktoro/Common/ThermoVector.hpp>
24 #include <Reaktoro/Math/Matrix.hpp>
25 
26 namespace Reaktoro {
27 
28 // Forward declaration
29 template<typename V, typename T, typename P, typename N>
30 class ChemicalVectorBase;
31 
40 
43 
46 
49 template<typename V, typename T, typename P, typename N>
51 {
52 public:
54  V val;
55 
57  T ddT;
58 
60  P ddP;
61 
63  N ddn;
64 
67 
70  explicit ChemicalVectorBase(Index nspecies)
71  : ChemicalVectorBase(nspecies, nspecies) {}
72 
76  ChemicalVectorBase(Index nrows, Index nspecies)
77  : ChemicalVectorBase(zeros(nrows), zeros(nrows), zeros(nrows), zeros(nrows, nspecies)) {}
78 
83  ChemicalVectorBase(Index nrows, Index nspecies, double val)
84  : ChemicalVectorBase(constants(nrows, val), zeros(nrows), zeros(nrows), zeros(nrows, nspecies)) {}
85 
91  ChemicalVectorBase(const V& val, const T& ddT, const P& ddP, const N& ddn)
92  : val(val), ddT(ddT), ddP(ddP), ddn(ddn) {}
93 
95  template<typename VR, typename TR, typename PR, typename NR>
97  : val(other.val), ddT(other.ddT), ddP(other.ddP), ddn(other.ddn) {}
98 
100  template<typename VR, typename TR, typename PR, typename NR>
102  : val(other.val), ddT(other.ddT), ddP(other.ddP), ddn(other.ddn) {}
103 
105  auto size() const -> Index
106  {
107  return val.size();
108  }
109 
113  auto resize(Index nrows, Index nspecies) -> void
114  {
115  val = zeros(nrows);
116  ddT = zeros(nrows);
117  ddP = zeros(nrows);
118  ddn = zeros(nrows, nspecies);
119  }
120 
123  auto resize(Index nspecies) -> void
124  {
125  resize(nspecies, nspecies);
126  }
127 
129  template<typename VR, typename NR>
130  auto fill(const ChemicalScalarBase<VR, NR>& other) -> void
131  {
132  val.fill(other.val);
133  ddT.fill(other.ddT);
134  ddP.fill(other.ddP);
135  for(auto i = 0; i < ddn.rows(); ++i) ddn.row(i) = other.ddn;
136  }
137 
139  template<typename VR>
140  auto fill(const ThermoScalarBase<VR>& other) -> void
141  {
142  val.fill(other.val);
143  ddT.fill(other.ddT);
144  ddP.fill(other.ddP);
145  ddn.fill(0.0);
146  }
147 
149  auto fill(double value) -> void
150  {
151  val.fill(value);
152  ddT.fill(0.0);
153  ddP.fill(0.0);
154  ddn.fill(0.0);
155  }
156 
158  template<typename VR, typename TR, typename PR, typename NR>
160  {
161  val = other.val;
162  ddT = other.ddT;
163  ddP = other.ddP;
164  ddn = other.ddn;
165  return *this;
166  }
167 
169  template<typename VR, typename NR>
171  {
172  val.fill(other.val);
173  ddT.fill(other.ddT);
174  ddP.fill(other.ddP);
175  for(auto i = 0; i < ddn.rows(); ++i) ddn.row(i) = other.ddn;
176  return *this;
177  }
178 
180  template<typename VR>
182  {
183  val.fill(other.val);
184  ddT.fill(other.ddT);
185  ddP.fill(other.ddP);
186  ddn.fill(0.0);
187  return *this;
188  }
189 
191  auto operator=(double other) -> ChemicalVectorBase&
192  {
193  val.fill(other);
194  ddT.fill(0.0);
195  ddP.fill(0.0);
196  ddn.fill(0.0);
197  return *this;
198  }
199 
201  template<typename VR, typename TR, typename PR, typename NR>
203  {
204  val += other.val;
205  ddT += other.ddT;
206  ddP += other.ddP;
207  ddn += other.ddn;
208  return *this;
209  }
210 
212  template<typename VR, typename TR, typename PR>
214  {
215  val += other.val;
216  ddT += other.ddT;
217  ddP += other.ddP;
218  return *this;
219  }
220 
222  template<typename VR>
224  {
225  val.array() += other.val;
226  ddT.array() += other.ddT;
227  ddP.array() += other.ddP;
228  return *this;
229  }
230 
232  auto operator+=(double other) -> ChemicalVectorBase&
233  {
234  val.array() += other;
235  return *this;
236  }
237 
239  template<typename VR, typename TR, typename PR, typename NR>
241  {
242  val -= other.val;
243  ddT -= other.ddT;
244  ddP -= other.ddP;
245  ddn -= other.ddn;
246  return *this;
247  }
248 
250  template<typename VR, typename TR, typename PR>
252  {
253  val -= other.val;
254  ddT -= other.ddT;
255  ddP -= other.ddP;
256  return *this;
257  }
258 
260  template<typename VR>
262  {
263  val.array() -= other.val;
264  ddT.array() -= other.ddT;
265  ddP.array() -= other.ddP;
266  return *this;
267  }
268 
270  auto operator-=(double other) -> ChemicalVectorBase&
271  {
272  val.array() -= other;
273  return *this;
274  }
275 
277  template<typename VR, typename TR, typename PR, typename NR>
279  {
280  ddT = diag(val) * other.ddT + diag(other.val) * ddT;
281  ddP = diag(val) * other.ddP + diag(other.val) * ddP;
282  ddn = diag(val) * other.ddn + diag(other.val) * ddn;
283  val = diag(val) * other.val;
284  return *this;
285  }
286 
288  auto operator*=(double other) -> ChemicalVectorBase&
289  {
290  val *= other;
291  ddT *= other;
292  ddP *= other;
293  ddn *= other;
294  return *this;
295  }
296 
298  auto operator/=(double other) -> ChemicalVectorBase&
299  {
300  *this *= 1.0/other;
301  return *this;
302  }
303 
305  auto operator[](Index irow) -> ChemicalScalarBase<decltype(val[irow]), decltype(ddn.row(irow))>
306  {
307  return {val[irow], ddT[irow], ddP[irow], ddn.row(irow)};
308  }
309 
311  auto operator[](Index irow) const -> ChemicalScalarBase<decltype(val[irow]), decltype(ddn.row(irow))>
312  {
313  return {val[irow], ddT[irow], ddP[irow], ddn.row(irow)};
314  }
315 
319  auto view(Index irow, Index nrows) -> ChemicalVectorRef
320  {
321  return {rowsmap(val, irow, nrows), rowsmap(ddT, irow, nrows), rowsmap(ddP, irow, nrows), rowsmap(ddn, irow, nrows)};
322  }
323 
327  auto view(Index irow, Index nrows) const -> ChemicalVectorConstRef
328  {
329  return {rowsmap(val, irow, nrows), rowsmap(ddT, irow, nrows), rowsmap(ddP, irow, nrows), rowsmap(ddn, irow, nrows)};
330  }
331 
337  auto view(Index irow, Index icol, Index nrows, Index ncols) -> ChemicalVectorRef
338  {
339  return {rowsmap(val, irow, nrows), rowsmap(ddT, irow, nrows), rowsmap(ddP, irow, nrows), blockmap(ddn, irow, icol, nrows, ncols)};
340  }
341 
347  auto view(Index irow, Index icol, Index nrows, Index ncols) const -> ChemicalVectorConstRef
348  {
349  return {rowsmap(val, irow, nrows), rowsmap(ddT, irow, nrows), rowsmap(ddP, irow, nrows), blockmap(ddn, irow, icol, nrows, ncols)};
350  }
351 
353  explicit operator Vector() const
354  {
355  return val;
356  }
357 };
358 
360 template<typename V, typename T, typename P, typename N>
361 auto zeros(const ChemicalVectorBase<V,T,P,N>& v) -> ChemicalVectorBase<decltype(zeros(0)), decltype(zeros(0)), decltype(zeros(0)), decltype(zeros(0,0))>
362 {
363  const Index n = v.size();
364  return {zeros(n), zeros(n), zeros(n), zeros(n, n)};
365 }
366 
368 template<typename V, typename T, typename P, typename N>
369 auto ones(const ChemicalVectorBase<V,T,P,N>& v) -> ChemicalVectorBase<decltype(ones(0)), decltype(zeros(0)), decltype(zeros(0)), decltype(zeros(0,0))>
370 {
371  const Index n = v.size();
372  return { ones(n), zeros(n), zeros(n), zeros(n, n) };
373 }
374 
376 class Composition : public ChemicalVectorBase<VectorConstRef, decltype(zeros(0)), decltype(zeros(0)), decltype(identity(0,0))>
377 {
378 public:
380  using Base = ChemicalVectorBase<VectorConstRef, decltype(zeros(0)), decltype(zeros(0)), decltype(identity(0,0))>;
381 
383  Composition(VectorConstRef n) : Base(n, zeros(n.rows()), zeros(n.rows()), identity(n.rows(), n.rows())) {}
384 
386  auto operator[](Index irow) const -> ChemicalScalarBase<double, decltype(unitrow(val.size(), irow))>
387  {
388  return { val[irow], 0.0, 0.0, unitrow(val.size(), irow) };
389  }
390 };
391 
392 template<typename V, typename T, typename P, typename N>
393 auto operator<<(std::ostream& out, const ChemicalVectorBase<V,T,P,N>& r) -> std::ostream&
394 {
395  out << r.val;
396  return out;
397 }
398 
399 template<typename V, typename T, typename P, typename N>
400 auto operator+(const ChemicalVectorBase<V,T,P,N>& r) -> ChemicalVectorBase<V,T,P,N>
401 {
402  return r;
403 }
404 
405 template<typename V, typename T, typename P, typename N>
406 auto operator-(const ChemicalVectorBase<V,T,P,N>& r) -> ChemicalVectorBase<decltype(-r.val), decltype(-r.ddT), decltype(-r.ddP), decltype(-r.ddn)>
407 {
408  return {-r.val, -r.ddT, -r.ddP, -r.ddn};
409 }
410 
411 template<typename VL, typename TL, typename PL, typename NL, typename VR, typename TR, typename PR, typename NR>
412 auto operator+(const ChemicalVectorBase<VL,TL,PL,NL>& l, const ChemicalVectorBase<VR,TR,PR,NR>& r) -> ChemicalVectorBase<decltype(l.val + r.val), decltype(l.ddT + r.ddT), decltype(l.ddP + r.ddP), decltype(l.ddn + r.ddn)>
413 {
414  return {l.val + r.val, l.ddT + r.ddT, l.ddP + r.ddP, l.ddn + r.ddn};
415 }
416 
417 template<typename VL, typename TL, typename PL, typename NL, typename VR, typename TR, typename PR>
418 auto operator+(const ChemicalVectorBase<VL,TL,PL,NL>& l, const ThermoVectorBase<VR,TR,PR>& r) -> ChemicalVectorBase<decltype(l.val + r.val), decltype(l.ddT + r.ddT), decltype(l.ddP + r.ddP), decltype(l.ddn)>
419 {
420  return {l.val + r.val, l.ddT + r.ddT, l.ddP + r.ddP, l.ddn};
421 }
422 
423 template<typename VL, typename TL, typename PL, typename VR, typename TR, typename PR, typename NR>
424 auto operator+(const ThermoVectorBase<VL,TL,PL>& l, const ChemicalVectorBase<VR,TR,PR,NR>& r) -> decltype(r + l)
425 {
426  return r + l;
427 }
428 
429 template<typename VL, typename TL, typename PL, typename NL, typename VR>
430 auto operator+(const ChemicalVectorBase<VL,TL,PL,NL>& l, const ThermoScalarBase<VR>& r) -> ChemicalVectorBase<decltype(l.val + r.val*ones(l.size())), decltype(l.ddT + r.ddT*ones(l.size())), decltype(l.ddP + r.ddP*ones(l.size())), decltype(l.ddn)>
431 {
432  return {l.val + r.val*ones(l.size()), l.ddT + r.ddT*ones(l.size()), l.ddP + r.ddP*ones(l.size()), l.ddn};
433 }
434 
435 template<typename VL, typename VR, typename TR, typename PR, typename NR>
436 auto operator+(const ThermoScalarBase<VL>& l, const ChemicalVectorBase<VR,TR,PR,NR>& r) -> decltype(r + l)
437 {
438  return r + l;
439 }
440 
441 template<typename V, typename T, typename P, typename N, typename Derived>
442 auto operator+(const ChemicalVectorBase<V,T,P,N>& l, const Eigen::MatrixBase<Derived>& r) -> ChemicalVectorBase<decltype(l.val + r),T,P,N>
443 {
444  return {l.val + r, l.ddT, l.ddP, l.ddn};
445 }
446 
447 template<typename V, typename T, typename P, typename N, typename Derived>
448 auto operator+(const Eigen::MatrixBase<Derived>& l, const ChemicalVectorBase<V,T,P,N>& r) -> decltype(r + l)
449 {
450  return r + l;
451 }
452 
453 template<typename VL, typename TL, typename PL, typename NL, typename VR, typename TR, typename PR, typename NR>
454 auto operator-(const ChemicalVectorBase<VL,TL,PL,NL>& l, const ChemicalVectorBase<VR,TR,PR,NR>& r) -> ChemicalVectorBase<decltype(l.val - r.val), decltype(l.ddT - r.ddT), decltype(l.ddP - r.ddP), decltype(l.ddn - r.ddn)>
455 {
456  return {l.val - r.val, l.ddT - r.ddT, l.ddP - r.ddP, l.ddn - r.ddn};
457 }
458 
459 template<typename VL, typename TL, typename PL, typename NL, typename VR, typename TR, typename PR>
460 auto operator-(const ChemicalVectorBase<VL,TL,PL,NL>& l, const ThermoVectorBase<VR,TR,PR>& r) -> ChemicalVectorBase<decltype(l.val - r.val), decltype(l.ddT - r.ddT), decltype(l.ddP - r.ddP), decltype(l.ddn)>
461 {
462  return {l.val - r.val, l.ddT - r.ddT, l.ddP - r.ddP, l.ddn};
463 }
464 
465 template<typename VL, typename TL, typename PL, typename VR, typename TR, typename PR, typename NR>
466 auto operator-(const ThermoVectorBase<VL,TL,PL>& l, const ChemicalVectorBase<VR,TR,PR,NR>& r) -> decltype(-(r - l))
467 {
468  return -(r - l);
469 }
470 
471 template<typename VL, typename TL, typename PL, typename NL, typename VR>
472 auto operator-(const ChemicalVectorBase<VL,TL,PL,NL>& l, const ThermoScalarBase<VR>& r) -> ChemicalVectorBase<decltype(l.val - r.val*ones(l.size())), decltype(l.ddT - r.ddT*ones(l.size())), decltype(l.ddP - r.ddP*ones(l.size())), decltype(l.ddn)>
473 {
474  return {l.val - r.val*ones(l.size()), l.ddT - r.ddT*ones(l.size()), l.ddP - r.ddP*ones(l.size()), l.ddn};
475 }
476 
477 template<typename VL, typename VR, typename TR, typename PR, typename NR>
478 auto operator-(const ThermoScalarBase<VL>& l, const ChemicalVectorBase<VR,TR,PR,NR>& r) -> decltype(r + l)
479 {
480  return -(r - l);
481 }
482 
483 template<typename V, typename T, typename P, typename N, typename Derived>
484 auto operator-(const ChemicalVectorBase<V,T,P,N>& l, const Eigen::MatrixBase<Derived>& r) -> ChemicalVectorBase<decltype(l.val - r),T,P,N>
485 {
486  return {l.val - r, l.ddT, l.ddP, l.ddn};
487 }
488 
489 template<typename V, typename T, typename P, typename N, typename Derived>
490 auto operator-(const Eigen::MatrixBase<Derived>& l, const ChemicalVectorBase<V,T,P,N>& r) -> decltype(-(r - l))
491 {
492  return -(r - l);
493 }
494 
495 template<typename V, typename T, typename P, typename N>
496 auto operator*(double l, const ChemicalVectorBase<V,T,P,N>& r) -> ChemicalVectorBase<decltype(l * r.val), decltype(l * r.ddT), decltype(l * r.ddP), decltype(l * r.ddn)>
497 {
498  return {l * r.val, l * r.ddT, l * r.ddP, l * r.ddn};
499 }
500 
501 template<typename V, typename T, typename P, typename N>
502 auto operator*(const ChemicalVectorBase<V,T,P,N>& l, double r) -> decltype(r * l)
503 {
504  return r * l;
505 }
506 
507 template<typename VL, typename VR, typename T, typename P, typename N>
508 auto operator*(const ThermoScalarBase<VL>& l, const ChemicalVectorBase<VR,T,P,N>& r) -> ChemicalVectorBase<decltype(l.val * r.val), decltype(l.val * r.ddT + l.ddT * r.val), decltype(l.val * r.ddP + l.ddP * r.val), decltype(l.val * r.ddn)>
509 {
510  return {l.val * r.val, l.val * r.ddT + l.ddT * r.val, l.val * r.ddP + l.ddP * r.val, l.val * r.ddn};
511 }
512 
513 template<typename VL, typename VR, typename T, typename P, typename N>
514 auto operator*(const ChemicalVectorBase<VL,T,P,N>& l, const ThermoScalarBase<VR>& r) -> decltype(r * l)
515 {
516  return r * l;
517 }
518 
519 template<typename VL, typename TL, typename PL, typename NL, typename VR, typename NR>
520 auto operator*(const ChemicalVectorBase<VL,TL,PL,NL>& l, const ChemicalScalarBase<VR,NR>& r) -> ChemicalVectorBase<decltype(l.val * r.val), decltype(l.val * r.ddT + l.ddT * r.val), decltype(l.val * r.ddP + l.ddP * r.val), decltype(l.val * r.ddn + l.ddn * r.val)>
521 {
522  return {l.val * r.val, l.val * r.ddT + l.ddT * r.val, l.val * r.ddP + l.ddP * r.val, l.val * r.ddn + l.ddn * r.val};
523 }
524 
525 template<typename VL, typename NL, typename VR, typename TR, typename PR, typename NR>
526 auto operator*(const ChemicalScalarBase<VL,NL>& l, const ChemicalVectorBase<VR,TR,PR,NR>& r) -> decltype(r * l)
527 {
528  return r * l;
529 }
530 
531 template<typename VL, typename TL, typename PL, typename NL, typename VR, typename TR, typename PR, typename NR>
532 auto operator%(const ChemicalVectorBase<VL, TL, PL, NL>& l, const ChemicalVectorBase<VR, TR, PR, NR>& r) -> ChemicalVectorBase<decltype(diag(l.val) * r.val), decltype(diag(l.val) * r.ddT + diag(r.val) * l.ddT), decltype(diag(l.val) * r.ddP + diag(r.val) * l.ddP), decltype(diag(l.val) * r.ddn + diag(r.val) * l.ddn)>
533 {
534  return { diag(l.val) * r.val,
535  diag(l.val) * r.ddT + diag(r.val) * l.ddT,
536  diag(l.val) * r.ddP + diag(r.val) * l.ddP,
537  diag(l.val) * r.ddn + diag(r.val) * l.ddn };
538 }
539 
540 template<typename VL, typename TL, typename PL, typename VR, typename TR, typename PR, typename NR>
541 auto operator%(const ThermoVectorBase<VL, TL, PL>& l, const ChemicalVectorBase<VR, TR, PR, NR>& r) -> ChemicalVectorBase<decltype(diag(l.val) * r.val), decltype(diag(l.val) * r.ddT + diag(r.val) * l.ddT), decltype(diag(l.val) * r.ddP + diag(r.val) * l.ddP), decltype(diag(l.val) * r.ddn)>
542 {
543  return { diag(l.val) * r.val,
544  diag(l.val) * r.ddT + diag(r.val) * l.ddT,
545  diag(l.val) * r.ddP + diag(r.val) * l.ddP,
546  diag(l.val) * r.ddn };
547 }
548 
549 
550 template<typename VL, typename TL, typename PL, typename NL, typename VR, typename TR, typename PR>
551 auto operator%(const ChemicalVectorBase<VL,TL,PL,NL>& l, const ThermoVectorBase<VR,TR,PR>& r) -> decltype(r % l)
552 {
553  return r % l;
554 }
555 
556 template<typename V, typename T, typename P, typename N, typename Derived>
557 auto operator%(const Eigen::MatrixBase<Derived>& l, const ChemicalVectorBase<V,T,P,N>& r) -> ChemicalVectorBase<decltype(diag(l) * r.val), decltype(diag(l) * r.ddT), decltype(diag(l) * r.ddP), decltype(diag(l) * r.ddn)>
558 {
559  return {diag(l) * r.val, diag(l) * r.ddT, diag(l) * r.ddP, diag(l) * r.ddn};
560 }
561 
562 template<typename VL, typename TL, typename PL, typename NL, typename Derived>
563 auto operator%(const ChemicalVectorBase<VL,TL,PL,NL>& l, const Eigen::MatrixBase<Derived>& r) -> decltype(r % l)
564 {
565  return r % l;
566 }
567 
568 template<typename VL, typename TL, typename PL, typename NL, typename VR, typename TR, typename PR, typename NR>
569 auto operator/(const ChemicalVectorBase<VL,TL,PL,NL>& l, const ChemicalVectorBase<VR,TR,PR,NR>& r) -> ChemicalVector
570 {
571  const Vector tmp = 1.0/(r.val % r.val);
572  return {l.val/r.val,
573  diag(tmp) * (diag(r.val) * l.ddT - diag(l.val) * r.ddT),
574  diag(tmp) * (diag(r.val) * l.ddP - diag(l.val) * r.ddP),
575  diag(tmp) * (diag(r.val) * l.ddn - diag(l.val) * r.ddn)};
576 }
577 
578 template<typename VL, typename TL, typename PL, typename NL, typename VR, typename TR, typename PR, typename NR>
579 auto operator/(const ChemicalVectorBase<VL,TL,PL,NL>& l, const ThermoVectorBase<VR,TR,PR>& r) -> ChemicalVector
580 {
581  const Vector tmp = 1.0/(r.val % r.val);
582  return {l.val/r.val,
583  diag(tmp) * (diag(r.val) * l.ddT - diag(l.val) * r.ddT),
584  diag(tmp) * (diag(r.val) * l.ddP - diag(l.val) * r.ddP),
585  diag(tmp) * (diag(r.val) * l.ddn)};
586 }
587 
588 template<typename VL, typename TL, typename PL, typename NL, typename VR, typename NR>
589 auto operator/(const ChemicalVectorBase<VL,TL,PL,NL>& l, const ChemicalScalarBase<VR,NR>& r) -> ChemicalVector
590 {
591  const double tmp = 1.0/(r.val * r.val);
592  return {l.val/r.val,
593  tmp * (r.val * l.ddT - l.val * r.ddT),
594  tmp * (r.val * l.ddP - l.val * r.ddP),
595  tmp * (r.val * l.ddn - l.val * r.ddn)};
596 }
597 
598 template<typename VL, typename VR, typename T, typename P, typename N>
599 auto operator/(const ChemicalVectorBase<VL,T,P,N>& l, const ThermoScalarBase<VR>& r) -> ChemicalVectorBase<decltype(l.val/r.val), decltype(double() * (l.ddT * r.val - l.val * r.ddT)), decltype(double() * (l.ddP * r.val - l.val * r.ddP)), decltype(double() * (l.ddn * r.val))>
600 {
601  const double tmp = 1.0/(r.val * r.val);
602  return {l.val/r.val,
603  tmp * (l.ddT * r.val - l.val * r.ddT),
604  tmp * (l.ddP * r.val - l.val * r.ddP),
605  tmp * (l.ddn * r.val)};
606 }
607 
608 template<typename V, typename T, typename P, typename N>
609 auto operator/(const ChemicalVectorBase<V,T,P,N>& l, double r) -> decltype((1.0/r) * l)
610 {
611  return (1.0/r) * l;
612 }
613 
614 template<typename VL, typename TL, typename PL, typename NL, typename VR, typename TR, typename PR, typename NR>
615 auto operator<(const ChemicalVectorBase<VL,TL,PL,NL>& l, const ChemicalVectorBase<VR,TR,PR,NR>& r) -> bool
616 {
617  return l.val < r.val;
618 }
619 
620 template<typename VL, typename TL, typename PL, typename NL, typename VR, typename TR, typename PR, typename NR>
621 auto operator<=(const ChemicalVectorBase<VL,TL,PL,NL>& l, const ChemicalVectorBase<VR,TR,PR,NR>& r) -> bool
622 {
623  return l.val <= r.val;
624 }
625 
626 template<typename VL, typename TL, typename PL, typename NL, typename VR, typename TR, typename PR, typename NR>
627 auto operator>(const ChemicalVectorBase<VL,TL,PL,NL>& l, const ChemicalVectorBase<VR,TR,PR,NR>& r) -> bool
628 {
629  return l.val > r.val;
630 }
631 
632 template<typename VL, typename TL, typename PL, typename NL, typename VR, typename TR, typename PR, typename NR>
633 auto operator>=(const ChemicalVectorBase<VL,TL,PL,NL>& l, const ChemicalVectorBase<VR,TR,PR,NR>& r) -> bool
634 {
635  return l.val >= r.val;
636 }
637 
638 template<typename VL, typename TL, typename PL, typename NL, typename VR, typename TR, typename PR, typename NR>
639 auto operator==(const ChemicalVectorBase<VL,TL,PL,NL>& l, const ChemicalVectorBase<VR,TR,PR,NR>& r) -> bool
640 {
641  return l.val == r.val;
642 }
643 
644 template<typename VL, typename TL, typename PL, typename NL, typename VR, typename TR, typename PR, typename NR>
645 auto operator!=(const ChemicalVectorBase<VL,TL,PL,NL>& l, const ChemicalVectorBase<VR,TR,PR,NR>& r) -> bool
646 {
647  return l.val == r.val;
648 }
649 
651 template<typename V, typename T, typename P, typename N>
652 auto row(ChemicalVectorBase<V,T,P,N>& vec, Index irow) -> ChemicalScalarBase<decltype(vec.val[irow]), decltype(vec.ddn.row(irow))>
653 {
654  return {vec.val[irow], vec.ddT[irow], vec.ddP[irow], vec.ddn.row(irow)};
655 }
656 
658 template<typename V, typename T, typename P, typename N>
659 auto row(const ChemicalVectorBase<V,T,P,N>& vec, Index irow) -> ChemicalScalarBase<decltype(vec.val[irow]), decltype(vec.ddn.row(irow))>
660 {
661  return {vec.val[irow], vec.ddT[irow], vec.ddP[irow], vec.ddn.row(irow)};
662 }
663 
665 template<typename V, typename T, typename P, typename N>
666 auto row(ChemicalVectorBase<V,T,P,N>& vec, Index irow, Index icol, Index ncols) -> ChemicalScalarBase<decltype(vec.val[irow]), decltype(vec.ddn.row(irow).segment(icol, ncols))>
667 {
668  return {vec.val[irow], vec.ddT[irow], vec.ddP[irow], vec.ddn.row(irow).segment(icol, ncols)};
669 }
670 
672 template<typename V, typename T, typename P, typename N>
673 auto row(const ChemicalVectorBase<V,T,P,N>& vec, Index irow, Index icol, Index ncols) -> ChemicalScalarBase<decltype(vec.val[irow]), decltype(vec.ddn.row(irow).segment(icol, ncols))>
674 {
675  return {vec.val[irow], vec.ddT[irow], vec.ddP[irow], vec.ddn.row(irow).segment(icol, ncols)};
676 }
677 
679 template<typename V, typename T, typename P, typename N>
680 auto rows(ChemicalVectorBase<V,T,P,N>& vec, Index irow, Index nrows) -> ChemicalVectorBase<decltype(rows(vec.val, irow, nrows)), decltype(rows(vec.ddT, irow, nrows)), decltype(rows(vec.ddP, irow, nrows)), decltype(rows(vec.ddn, irow, nrows))>
681 {
682  return {rows(vec.val, irow, nrows), rows(vec.ddT, irow, nrows), rows(vec.ddP, irow, nrows), rows(vec.ddn, irow, nrows)};
683 }
684 
686 template<typename V, typename T, typename P, typename N>
687 auto rows(const ChemicalVectorBase<V,T,P,N>& vec, Index irow, Index nrows) -> ChemicalVectorBase<decltype(rows(vec.val, irow, nrows)), decltype(rows(vec.ddT, irow, nrows)), decltype(rows(vec.ddP, irow, nrows)), decltype(rows(vec.ddn, irow, nrows))>
688 {
689  return {rows(vec.val, irow, nrows), rows(vec.ddT, irow, nrows), rows(vec.ddP, irow, nrows), rows(vec.ddn, irow, nrows)};
690 }
691 
693 template<typename V, typename T, typename P, typename N>
694 auto rows(ChemicalVectorBase<V,T,P,N>& vec, Index irow, Index icol, Index nrows, Index ncols) -> ChemicalVectorBase<decltype(rows(vec.val, irow, nrows)), decltype(rows(vec.ddT, irow, nrows)), decltype(rows(vec.ddP, irow, nrows)), decltype(block(vec.ddn, irow, icol, nrows, ncols))>
695 {
696  return {rows(vec.val, irow, nrows), rows(vec.ddT, irow, nrows), rows(vec.ddP, irow, nrows), block(vec.ddn, irow, icol, nrows, ncols)};
697 }
698 
700 template<typename V, typename T, typename P, typename N>
701 auto rows(const ChemicalVectorBase<V,T,P,N>& vec, Index irow, Index icol, Index nrows, Index ncols) -> ChemicalVectorBase<decltype(rows(vec.val, irow, nrows)), decltype(rows(vec.ddT, irow, nrows)), decltype(rows(vec.ddP, irow, nrows)), decltype(block(vec.ddn, irow, icol, nrows, ncols))>
702 {
703  return {rows(vec.val, irow, nrows), rows(vec.ddT, irow, nrows), rows(vec.ddP, irow, nrows), block(vec.ddn, irow, icol, nrows, ncols)};
704 }
705 
707 template<typename V, typename T, typename P, typename N>
708 auto rows(ChemicalVectorBase<V,T,P,N>& vec, const Indices& irows) -> ChemicalVectorBase<decltype(rows(vec.val, irows)), decltype(rows(vec.ddT, irows)), decltype(rows(vec.ddP, irows)), decltype(rows(vec.ddn, irows))>
709 {
710  return {rows(vec.val, irows), rows(vec.ddT, irows), rows(vec.ddP, irows), rows(vec.ddn, irows)};
711 }
712 
714 template<typename V, typename T, typename P, typename N>
715 auto rows(const ChemicalVectorBase<V,T,P,N>& vec, const Indices& irows) -> ChemicalVectorBase<decltype(rows(vec.val, irows)), decltype(rows(vec.ddT, irows)), decltype(rows(vec.ddP, irows)), decltype(rows(vec.ddn, irows))>
716 {
717  return {rows(vec.val, irows), rows(vec.ddT, irows), rows(vec.ddP, irows), rows(vec.ddn, irows)};
718 }
719 
721 template<typename V, typename T, typename P, typename N>
722 auto rows(ChemicalVectorBase<V,T,P,N>& vec, const Indices& irows, const Indices& icols) -> ChemicalVectorBase<decltype(rows(vec.val, irows)), decltype(rows(vec.ddT, irows)), decltype(rows(vec.ddP, irows)), decltype(submatrix(vec.ddn, irows, icols))>
723 {
724  return {rows(vec.val, irows), rows(vec.ddT, irows), rows(vec.ddP, irows), submatrix(vec.ddn, irows, icols)};
725 }
726 
728 template<typename V, typename T, typename P, typename N>
729 auto rows(const ChemicalVectorBase<V,T,P,N>& vec, const Indices& irows, const Indices& icols) -> ChemicalVectorBase<decltype(rows(vec.val, irows)), decltype(rows(vec.ddT, irows)), decltype(rows(vec.ddP, irows)), decltype(submatrix(vec.ddn, irows, icols))>
730 {
731  return {rows(vec.val, irows), rows(vec.ddT, irows), rows(vec.ddP, irows), submatrix(vec.ddn, irows, icols)};
732 }
733 
734 template<typename V, typename T, typename P, typename N>
735 auto abs(const ChemicalVectorBase<V,T,P,N>& l) -> ChemicalVector
736 {
737  const Vector tmp1 = abs(l.val);
738  const Vector tmp2 = l.val/tmp1;
739  return {tmp1, diag(tmp2) * l.ddT, diag(tmp2) * l.ddP, diag(tmp2) * l.ddn};
740 }
741 
742 template<typename V, typename T, typename P, typename N>
743 auto sqrt(const ChemicalVectorBase<V,T,P,N>& l) -> ChemicalVector
744 {
745  const Vector tmp1 = sqrt(l.val);
746  const Vector tmp2 = 0.5 * tmp1/l.val;
747  return {tmp1, diag(tmp2) * l.ddT, diag(tmp2) * l.ddP, diag(tmp2) * l.ddn};
748 }
749 
750 template<typename V, typename T, typename P, typename N>
751 auto pow(const ChemicalVectorBase<V,T,P,N>& l, double power) -> ChemicalVector
752 {
753  const Vector tmp1 = pow(l.val, power);
754  const Vector tmp2 = power * tmp1/l.val;
755  return {tmp1, diag(tmp2) * l.ddT, diag(tmp2) * l.ddP, diag(tmp2) * l.ddn};
756 }
757 
758 template<typename V, typename T, typename P, typename N>
759 auto exp(const ChemicalVectorBase<V,T,P,N>& l) -> ChemicalVector
760 {
761  const Vector tmp = exp(l.val);
762  return {tmp, diag(tmp) * l.ddT, diag(tmp) * l.ddP, diag(tmp) * l.ddn};
763 }
764 
765 template<typename V, typename T, typename P, typename N>
766 auto log(const ChemicalVectorBase<V,T,P,N>& l) -> ChemicalVector
767 {
768  const Vector tmp1 = log(l.val);
769  const Vector tmp2 = 1.0/l.val;
770  return {tmp1, diag(tmp2) * l.ddT, diag(tmp2) * l.ddP, diag(tmp2) * l.ddn};
771 }
772 
773 template<typename V, typename T, typename P, typename N>
774 auto log10(const ChemicalVectorBase<V,T,P,N>& l) -> decltype(log(l)/double())
775 {
776  const double ln10 = 2.302585092994046;
777  return log(l)/ln10;
778 }
779 
780 template<typename V, typename T, typename P, typename N>
781 auto min(const ChemicalVectorBase<V,T,P,N>& r) -> double
782 {
783  return min(r.val);
784 }
785 
786 template<typename V, typename T, typename P, typename N>
787 auto max(const ChemicalVectorBase<V,T,P,N>& r) -> double
788 {
789  return max(r.val);
790 }
791 
792 template<typename V, typename T, typename P, typename N>
793 auto sum(const ChemicalVectorBase<V, T, P, N>& r) -> ChemicalScalarBase<double, decltype(r.ddn.colwise().sum())>
794 {
795  return { r.val.sum(), r.ddT.sum(), r.ddP.sum(), r.ddn.colwise().sum() };
796 }
797 
798 template<typename VL, typename TL, typename PL, typename NL, typename VR, typename TR, typename PR, typename NR>
799 auto dot(const ChemicalVectorBase<VL,TL,PL,NL>& l, const ChemicalVectorBase<VR,TR,PR,NR>& r) -> ChemicalScalarBase<double, decltype(tr(l.val) * r.ddn + tr(r.val) * l.ddn)>
800 {
801  return {tr(l.val) * r.val, tr(l.val) * r.ddT + tr(r.val) * l.ddT, tr(l.val) * r.ddP + tr(r.val) * l.ddP, tr(l.val) * r.ddn + tr(r.val) * l.ddn};
802 }
803 
804 template<typename VL, typename TL, typename PL, typename NL, typename VR, typename TR, typename PR, typename NR>
805 auto dot(const ThermoVectorBase<VL,TL,PL>& l, const ChemicalVectorBase<VR,TR,PR,NR>& r) -> ChemicalScalarBase<double, decltype(tr(l.val) * r.ddn)>
806 {
807  return {tr(l.val) * r.val, tr(l.val) * r.ddT + tr(r.val) * l.ddT, tr(l.val) * r.ddP + tr(r.val) * l.ddP, tr(l.val) * r.ddn};
808 }
809 
810 template<typename VL, typename TL, typename PL, typename NL, typename VR, typename TR, typename PR, typename NR>
811 auto dot(const ChemicalVectorBase<VR,TR,PR,NR>& l, const ThermoVectorBase<VL,TL,PL>& r) -> decltype(dot(r, l))
812 {
813  return dot(r, l);
814 }
815 
816 template<typename V, typename T, typename P, typename N, typename Derived>
817 auto dot(const Eigen::MatrixBase<Derived>& l, const ChemicalVectorBase<V,T,P,N>& r) -> ChemicalScalarBase<double, decltype(tr(l) * r.ddn)>
818 {
819  return {tr(l) * r.val, tr(l) * r.ddT, tr(l) * r.ddP, tr(l) * r.ddn};
820 }
821 
822 template<typename V, typename T, typename P, typename N, typename Derived>
823 auto dot(const ChemicalVectorBase<V,T,P,N>& l, const Eigen::MatrixBase<Derived>& r) -> decltype(dot(r, l))
824 {
825  return dot(r, l);
826 }
827 
828 } // namespace Reaktoro
auto view(Index irow, Index icol, Index nrows, Index ncols) -> ChemicalVectorRef
Return a view of an interval of the ChemicalVectorBase instance.
Definition: ChemicalVector.hpp:337
auto rowsmap(Eigen::Matrix< Scalar, Rows, Cols, Options, MaxRows, MaxCols > &mat, Index row, Index nrows) -> Eigen::Map< Eigen::Matrix< Scalar, Rows, Cols, Options, MaxRows, MaxCols >, Eigen::Unaligned, Eigen::Stride< Rows, Cols >>
Return a mapped view of a sequence of rows of a matrix.
Definition: Matrix.hpp:266
ChemicalVectorBase(Index nrows, Index nspecies, double val)
Construct a ChemicalVectorBase instance with given number of rows and species, and a constant value.
Definition: ChemicalVector.hpp:83
auto operator=(const ChemicalVectorBase< VR, TR, PR, NR > &other) -> ChemicalVectorBase &
Assign another ChemicalVectorBase instance to this.
Definition: ChemicalVector.hpp:159
auto operator-=(const ThermoScalarBase< VR > &other) -> ChemicalVectorBase &
Assign-subtraction of a ThermoScalarBase instance to this.
Definition: ChemicalVector.hpp:261
auto size() const -> Index
Return the number of rows in this ChemicalVectorBase instance.
Definition: ChemicalVector.hpp:105
auto fill(const ThermoScalarBase< VR > &other) -> void
Assign a ThermoScalarBase instance to this.
Definition: ChemicalVector.hpp:140
auto operator=(double other) -> ChemicalVectorBase &
Assign a scalar to this.
Definition: ChemicalVector.hpp:191
A template base class to represent a chemical scalar and its partial derivatives.
Definition: ChemicalScalar.hpp:49
V val
The value of the chemical scalar.
Definition: ChemicalScalar.hpp:52
auto submatrix(Eigen::MatrixBase< Derived > &mat, const Indices &irows, const Indices &icols) -> decltype(mat(irows, icols))
Return a view of some rows and columns of a matrix.
Definition: Matrix.hxx:155
ChemicalVectorBase(ChemicalVectorBase< VR, TR, PR, NR > &other)
Construct a ChemicalVectorBase instance from another.
Definition: ChemicalVector.hpp:96
auto row(ChemicalVectorBase< V, T, P, N > &vec, Index irow) -> ChemicalScalarBase< decltype(vec.val[irow]), decltype(vec.ddn.row(irow))>
Return a ChemicalScalarBase with reference to the chemical scalar in a given row.
Definition: ChemicalVector.hpp:652
ChemicalVectorBase(const ChemicalVectorBase< VR, TR, PR, NR > &other)
Construct a ChemicalVectorBase instance from another.
Definition: ChemicalVector.hpp:101
auto operator*=(double other) -> ChemicalVectorBase &
Assign-multiplication of a scalar to this.
Definition: ChemicalVector.hpp:288
auto operator/=(double other) -> ChemicalVectorBase &
Assign-division of a scalar to this.
Definition: ChemicalVector.hpp:298
auto operator+=(const ThermoVectorBase< VR, TR, PR > &other) -> ChemicalVectorBase &
Assign-addition of a ThermoVectorBase instance to this.
Definition: ChemicalVector.hpp:213
The namespace containing all components of the Reaktoro library.
Definition: ChemicalScalar.hpp:24
auto operator-=(const ThermoVectorBase< VR, TR, PR > &other) -> ChemicalVectorBase &
Assign-subtraction of a ThermoVectorBase instance to this.
Definition: ChemicalVector.hpp:251
auto identity(Index rows, Index cols) -> decltype(Matrix::Identity(rows, cols))
Return an expression of an identity matrix.
Definition: Matrix.hxx:77
A type that describes temperature in units of K.
Definition: ChemicalVector.hpp:377
auto zeros(const ChemicalVectorBase< V, T, P, N > &v) -> ChemicalVectorBase< decltype(zeros(0)), decltype(zeros(0)), decltype(zeros(0)), decltype(zeros(0, 0))>
Return a ChemicalVectorBase expression representing zeros with same dimension of given vector.
Definition: ChemicalVector.hpp:361
ChemicalVectorBase(const V &val, const T &ddT, const P &ddP, const N &ddn)
Construct a ChemicalVectorBase instance with given values and derivatives.
Definition: ChemicalVector.hpp:91
auto rows(ChemicalVectorBase< V, T, P, N > &vec, Index irow, Index nrows) -> ChemicalVectorBase< decltype(rows(vec.val, irow, nrows)), decltype(rows(vec.ddT, irow, nrows)), decltype(rows(vec.ddP, irow, nrows)), decltype(rows(vec.ddn, irow, nrows))>
Return a reference of a sequence of rows of this ChemicalVectorBase instance.
Definition: ChemicalVector.hpp:680
auto fill(const ChemicalScalarBase< VR, NR > &other) -> void
Assign a ChemicalScalarBase instance to this.
Definition: ChemicalVector.hpp:130
auto resize(Index nrows, Index nspecies) -> void
Resize this ChemicalVectorBase instance with new number of rows and number of species.
Definition: ChemicalVector.hpp:113
auto blockmap(Eigen::Matrix< Scalar, Rows, Cols, Options, MaxRows, MaxCols > &mat, Index row, Index col, Index nrows, Index ncols) -> Eigen::Map< Eigen::Matrix< Scalar, Rows, Cols, Options, MaxRows, MaxCols >, Eigen::Unaligned, Eigen::Stride< Rows, Cols >>
Return a block mapped view of a matrix.
Definition: Matrix.hpp:240
auto operator-=(const ChemicalVectorBase< VR, TR, PR, NR > &other) -> ChemicalVectorBase &
Assign-subtraction of a ChemicalVectorBase instance to this.
Definition: ChemicalVector.hpp:240
N ddn
The matrix of partial mole derivatives of the chemical scalars.
Definition: ChemicalVector.hpp:63
auto operator+=(const ChemicalVectorBase< VR, TR, PR, NR > &other) -> ChemicalVectorBase &
Assign-addition of a ChemicalVectorBase instance to this.
Definition: ChemicalVector.hpp:202
auto operator*=(const ChemicalVectorBase< VR, TR, PR, NR > &other) -> ChemicalVectorBase &
Assign-multiplication of a ChemicalVectorBase instance to this.
Definition: ChemicalVector.hpp:278
auto operator+=(double other) -> ChemicalVectorBase &
Assign-addition of a scalar to this.
Definition: ChemicalVector.hpp:232
T ddT
The vector of partial temperature derivatives of the chemical scalars.
Definition: ChemicalVector.hpp:57
auto operator=(const ChemicalScalarBase< VR, NR > &other) -> ChemicalVectorBase &
Assign a ChemicalScalarBase instance to this.
Definition: ChemicalVector.hpp:170
auto operator[](Index irow) -> ChemicalScalarBase< decltype(val[irow]), decltype(ddn.row(irow))>
Return a ChemicalScalarBase with reference to the chemical scalar in a given row.
Definition: ChemicalVector.hpp:305
ChemicalVectorBase(Index nspecies)
Construct a ChemicalVectorBase instance with number of rows equal to given number of species.
Definition: ChemicalVector.hpp:70
std::size_t Index
Define a type that represents an index.
Definition: Index.hpp:26
ChemicalVectorBase()
Construct a default ChemicalVectorBase instance.
Definition: ChemicalVector.hpp:66
Composition(VectorConstRef n)
Construct a Composition instance with given composition vector.
Definition: ChemicalVector.hpp:383
P ddP
The vector of partial pressure derivatives of the chemical scalars.
Definition: ChemicalVector.hpp:60
auto fill(double value) -> void
Assign a scalarsto this.
Definition: ChemicalVector.hpp:149
V val
The vector of chemical scalars.
Definition: ChemicalVector.hpp:54
auto resize(Index nspecies) -> void
Resize this ChemicalVectorBase instance with number of rows equal the number of species.
Definition: ChemicalVector.hpp:123
std::vector< Index > Indices
Define a type that represents a collection of indices.
Definition: Index.hpp:29
auto operator-=(double other) -> ChemicalVectorBase &
Assign-subtraction of a scalar to this.
Definition: ChemicalVector.hpp:270
auto view(Index irow, Index nrows) const -> ChemicalVectorConstRef
Return a view of an interval of the ChemicalVectorBase instance.
Definition: ChemicalVector.hpp:327
auto view(Index irow, Index nrows) -> ChemicalVectorRef
Return a view of an interval of the ChemicalVectorBase instance.
Definition: ChemicalVector.hpp:319
A template base class to represent a thermodynamic scalar and its partial derivatives.
Definition: ThermoScalar.hpp:45
auto operator[](Index irow) const -> ChemicalScalarBase< double, decltype(unitrow(val.size(), irow))>
Return a ChemicalScalarBase with const reference to the chemical scalar in a given row.
Definition: ChemicalVector.hpp:386
auto ones(const ChemicalVectorBase< V, T, P, N > &v) -> ChemicalVectorBase< decltype(ones(0)), decltype(zeros(0)), decltype(zeros(0)), decltype(zeros(0, 0))>
Return a ChemicalVectorBase expression representing ones with same dimension of given vector.
Definition: ChemicalVector.hpp:369
ChemicalVectorBase< Vector, Vector, Vector, Matrix > ChemicalVector
A type that represents a vector of chemical properties and their derivatives.
Definition: ChemicalVector.hpp:39
auto tr(Eigen::MatrixBase< Derived > &mat) -> decltype(mat.transpose())
Return the transpose of the matrix.
Definition: Matrix.hxx:167
auto operator=(const ThermoScalarBase< VR > &other) -> ChemicalVectorBase &
Assign a ThermoScalarBase instance to this.
Definition: ChemicalVector.hpp:181
auto operator+=(const ThermoScalarBase< VR > &other) -> ChemicalVectorBase &
Assign-addition of a ThermoScalarBase instance to this.
Definition: ChemicalVector.hpp:223
auto operator[](Index irow) const -> ChemicalScalarBase< decltype(val[irow]), decltype(ddn.row(irow))>
Return a ChemicalScalarBase with const reference to the chemical scalar in a given row.
Definition: ChemicalVector.hpp:311
auto block(Eigen::MatrixBase< Derived > &mat, Index irow, Index icol, Index nrows, Index ncols) -> decltype(mat.block(irow, icol, nrows, ncols))
Return a view of some rows and columns of a matrix.
Definition: Matrix.hxx:143
auto view(Index irow, Index icol, Index nrows, Index ncols) const -> ChemicalVectorConstRef
Return a view of an interval of the ChemicalVectorBase instance.
Definition: ChemicalVector.hpp:347
Eigen::Ref< const Eigen::VectorXd > VectorConstRef
< Alias to Eigen type Ref<VectorXd>.
Definition: Matrix.hpp:31
auto diag(const Eigen::MatrixBase< Derived > &vec) -> decltype(vec.asDiagonal())
Return a diagonal matrix representation of a vector.
Definition: Matrix.hxx:185
A template base class to represent a vector of thermodynamic scalars and their partial derivatives.
Definition: ThermoVector.hpp:48
ChemicalVectorBase(Index nrows, Index nspecies)
Construct a ChemicalVectorBase instance with given number of rows and species.
Definition: ChemicalVector.hpp:76