]> www.ginac.de Git - ginac.git/blobdiff - ginac/tensor.cpp
Synced to HEAD
[ginac.git] / ginac / tensor.cpp
index 48393665ad5fc0555add8a19bba463af14abf288..a8c83b65ab20edeb18258411bd1e42dfa434788d 100644 (file)
@@ -3,7 +3,7 @@
  *  Implementation of GiNaC's special tensors. */
 
 /*
- *  GiNaC Copyright (C) 1999-2003 Johannes Gutenberg University Mainz, Germany
+ *  GiNaC Copyright (C) 1999-2004 Johannes Gutenberg University Mainz, Germany
  *
  *  This program is free software; you can redistribute it and/or modify
  *  it under the terms of the GNU General Public License as published by
 #include "indexed.h"
 #include "symmetry.h"
 #include "relational.h"
+#include "operators.h"
 #include "lst.h"
 #include "numeric.h"
 #include "matrix.h"
-#include "print.h"
 #include "archive.h"
 #include "utils.h"
 
 namespace GiNaC {
 
 GINAC_IMPLEMENT_REGISTERED_CLASS(tensor, basic)
-GINAC_IMPLEMENT_REGISTERED_CLASS(tensdelta, tensor)
-GINAC_IMPLEMENT_REGISTERED_CLASS(tensmetric, tensor)
-GINAC_IMPLEMENT_REGISTERED_CLASS(minkmetric, tensmetric)
-GINAC_IMPLEMENT_REGISTERED_CLASS(spinmetric, tensmetric)
-GINAC_IMPLEMENT_REGISTERED_CLASS(tensepsilon, tensor)
+
+GINAC_IMPLEMENT_REGISTERED_CLASS_OPT(tensdelta, tensor,
+  print_func<print_dflt>(&tensdelta::do_print).
+  print_func<print_latex>(&tensdelta::do_print_latex))
+
+GINAC_IMPLEMENT_REGISTERED_CLASS_OPT(tensmetric, tensor,
+  print_func<print_dflt>(&tensmetric::do_print).
+  print_func<print_latex>(&tensmetric::do_print))
+
+GINAC_IMPLEMENT_REGISTERED_CLASS_OPT(minkmetric, tensmetric,
+  print_func<print_dflt>(&minkmetric::do_print).
+  print_func<print_latex>(&minkmetric::do_print_latex))
+
+GINAC_IMPLEMENT_REGISTERED_CLASS_OPT(spinmetric, tensmetric,
+  print_func<print_dflt>(&spinmetric::do_print).
+  print_func<print_latex>(&spinmetric::do_print_latex))
+
+GINAC_IMPLEMENT_REGISTERED_CLASS_OPT(tensepsilon, tensor,
+  print_func<print_dflt>(&tensepsilon::do_print).
+  print_func<print_latex>(&tensepsilon::do_print_latex))
 
 //////////
-// default ctor, dtor, copy ctor, assignment operator and helpers
+// constructors
 //////////
 
-DEFAULT_CTORS(tensor)
-DEFAULT_CTORS(tensdelta)
-DEFAULT_CTORS(tensmetric)
-DEFAULT_COPY(spinmetric)
-DEFAULT_DESTROY(spinmetric)
-DEFAULT_DESTROY(minkmetric)
-DEFAULT_DESTROY(tensepsilon)
+tensor::tensor() : inherited(TINFO_tensor)
+{
+       setflag(status_flags::evaluated | status_flags::expanded);
+}
+
+DEFAULT_CTOR(tensdelta)
+DEFAULT_CTOR(tensmetric)
 
 minkmetric::minkmetric() : pos_sig(false)
 {
@@ -72,12 +87,6 @@ minkmetric::minkmetric(bool ps) : pos_sig(ps)
        tinfo_key = TINFO_minkmetric;
 }
 
-void minkmetric::copy(const minkmetric & other)
-{
-       inherited::copy(other);
-       pos_sig = other.pos_sig;
-}
-
 tensepsilon::tensepsilon() : minkowski(false), pos_sig(false)
 {
        tinfo_key = TINFO_tensepsilon;
@@ -88,13 +97,6 @@ tensepsilon::tensepsilon(bool mink, bool ps) : minkowski(mink), pos_sig(ps)
        tinfo_key = TINFO_tensepsilon;
 }
 
-void tensepsilon::copy(const tensepsilon & other)
-{
-       inherited::copy(other);
-       minkowski = other.minkowski;
-       pos_sig = other.pos_sig;
-}
-
 //////////
 // archiving
 //////////
@@ -106,7 +108,7 @@ DEFAULT_ARCHIVING(spinmetric)
 DEFAULT_UNARCHIVE(minkmetric)
 DEFAULT_UNARCHIVE(tensepsilon)
 
-minkmetric::minkmetric(const archive_node &n, const lst &sym_lst) : inherited(n, sym_lst)
+minkmetric::minkmetric(const archive_node &n, lst &sym_lst) : inherited(n, sym_lst)
 {
        n.find_bool("pos_sig", pos_sig);
 }
@@ -117,7 +119,7 @@ void minkmetric::archive(archive_node &n) const
        n.add_bool("pos_sig", pos_sig);
 }
 
-tensepsilon::tensepsilon(const archive_node &n, const lst &sym_lst) : inherited(n, sym_lst)
+tensepsilon::tensepsilon(const archive_node &n, lst &sym_lst) : inherited(n, sym_lst)
 {
        n.find_bool("minkowski", minkowski);
        n.find_bool("pos_sig", pos_sig);
@@ -183,7 +185,10 @@ ex tensdelta::eval_indexed(const basic & i) const
        // dimension
        if (!i1.get_dim().is_equal(i2.get_dim())) {
                ex min_dim = i1.minimal_dim(i2);
-               return i.subs(lst(i1 == i1.replace_dim(min_dim), i2 == i2.replace_dim(min_dim)));
+               exmap m;
+               m[i1] = i1.replace_dim(min_dim);
+               m[i2] = i2.replace_dim(min_dim);
+               return i.subs(m, subs_options::no_pattern);
        }
 
        // Trace of delta tensor is the (effective) dimension of the space
@@ -224,7 +229,10 @@ ex tensmetric::eval_indexed(const basic & i) const
        // dimension
        if (!i1.get_dim().is_equal(i2.get_dim())) {
                ex min_dim = i1.minimal_dim(i2);
-               return i.subs(lst(i1 == i1.replace_dim(min_dim), i2 == i2.replace_dim(min_dim)));
+               exmap m;
+               m[i1] = i1.replace_dim(min_dim);
+               m[i2] = i2.replace_dim(min_dim);
+               return i.subs(m, subs_options::no_pattern);
        }
 
        // A metric tensor with one covariant and one contravariant index gets
@@ -312,15 +320,15 @@ ex tensepsilon::eval_indexed(const basic & i) const
                // a canonic order but we can't assume what exactly that order is)
                std::vector<int> v;
                v.reserve(i.nops() - 1);
-               for (unsigned j=1; j<i.nops(); j++)
+               for (size_t j=1; j<i.nops(); j++)
                        v.push_back(ex_to<numeric>(ex_to<idx>(i.op(j)).get_value()).to_int());
                int sign = permutation_sign(v.begin(), v.end());
 
                // In a Minkowski space, check for covariant indices
                if (minkowski) {
-                       for (unsigned j=1; j<i.nops(); j++) {
+                       for (size_t j=1; j<i.nops(); j++) {
                                const ex & x = i.op(j);
-                               if (!is_ex_of_type(x, varidx))
+                               if (!is_a<varidx>(x))
                                        throw(std::runtime_error("indices of epsilon tensor in Minkowski space must be of type varidx"));
                                if (ex_to<varidx>(x).is_covariant())
                                        if (ex_to<idx>(x).get_value().is_zero())
@@ -346,7 +354,7 @@ bool tensor::replace_contr_index(exvector::iterator self, exvector::iterator oth
 
 again:
        if (self_idx->is_symbolic()) {
-               for (unsigned i=1; i<other->nops(); i++) {
+               for (size_t i=1; i<other->nops(); i++) {
                        const idx &other_idx = ex_to<idx>(other->op(i));
                        if (is_dummy_pair(*self_idx, other_idx)) {
 
@@ -400,11 +408,11 @@ bool tensmetric::contract_with(exvector::iterator self, exvector::iterator other
 
        // If contracting with the delta tensor, let the delta do it
        // (don't raise/lower delta indices)
-       if (is_ex_of_type(other->op(0), tensdelta))
+       if (is_a<tensdelta>(other->op(0)))
                return false;
 
        // Replace the dummy index with this tensor's other index and remove
-       // the tensor (this is valid for contractions with all other tensors)
+       // the tensor
        return replace_contr_index(self, other);
 }
 
@@ -417,7 +425,7 @@ bool spinmetric::contract_with(exvector::iterator self, exvector::iterator other
        GINAC_ASSERT(is_a<spinmetric>(self->op(0)));
 
        // Contractions between spinor metrics
-       if (is_ex_of_type(other->op(0), spinmetric)) {
+       if (is_a<spinmetric>(other->op(0))) {
                const idx &self_i1 = ex_to<idx>(self->op(1));
                const idx &self_i2 = ex_to<idx>(self->op(2));
                const idx &other_i1 = ex_to<idx>(other->op(1));
@@ -450,7 +458,7 @@ bool spinmetric::contract_with(exvector::iterator self, exvector::iterator other
 
        // If contracting with the delta tensor, let the delta do it
        // (don't raise/lower delta indices)
-       if (is_ex_of_type(other->op(0), tensdelta))
+       if (is_a<tensdelta>(other->op(0)))
                return false;
 
        // Try to contract first index
@@ -461,7 +469,7 @@ bool spinmetric::contract_with(exvector::iterator self, exvector::iterator other
 
 again:
        if (self_idx->is_symbolic()) {
-               for (unsigned i=1; i<other->nops(); i++) {
+               for (size_t i=1; i<other->nops(); i++) {
                        const idx &other_idx = ex_to<idx>(other->op(i));
                        if (is_dummy_pair(*self_idx, other_idx)) {
 
@@ -494,15 +502,15 @@ bool tensepsilon::contract_with(exvector::iterator self, exvector::iterator othe
        GINAC_ASSERT(is_a<indexed>(*self));
        GINAC_ASSERT(is_a<indexed>(*other));
        GINAC_ASSERT(is_a<tensepsilon>(self->op(0)));
-       unsigned num = self->nops() - 1;
+       size_t num = self->nops() - 1;
 
-       if (is_ex_exactly_of_type(other->op(0), tensepsilon) && num+1 == other->nops()) {
+       if (is_exactly_a<tensepsilon>(other->op(0)) && num+1 == other->nops()) {
 
                // Contraction of two epsilon tensors is a determinant
                bool variance = is_a<varidx>(self->op(1));
                matrix M(num, num);
-               for (unsigned i=0; i<num; i++) {
-                       for (unsigned j=0; j<num; j++) {
+               for (size_t i=0; i<num; i++) {
+                       for (size_t j=0; j<num; j++) {
                                if (minkowski)
                                        M(i, j) = lorentz_g(self->op(i+1), other->op(j+1), pos_sig);
                                else if (variance)
@@ -526,41 +534,52 @@ bool tensepsilon::contract_with(exvector::iterator self, exvector::iterator othe
 
 ex delta_tensor(const ex & i1, const ex & i2)
 {
-       if (!is_ex_of_type(i1, idx) || !is_ex_of_type(i2, idx))
+       static ex delta = (new tensdelta)->setflag(status_flags::dynallocated);
+
+       if (!is_a<idx>(i1) || !is_a<idx>(i2))
                throw(std::invalid_argument("indices of delta tensor must be of type idx"));
 
-       return indexed(tensdelta(), sy_symm(), i1, i2);
+       return indexed(delta, symmetric2(), i1, i2);
 }
 
 ex metric_tensor(const ex & i1, const ex & i2)
 {
-       if (!is_ex_of_type(i1, varidx) || !is_ex_of_type(i2, varidx))
+       static ex metric = (new tensmetric)->setflag(status_flags::dynallocated);
+
+       if (!is_a<varidx>(i1) || !is_a<varidx>(i2))
                throw(std::invalid_argument("indices of metric tensor must be of type varidx"));
 
-       return indexed(tensmetric(), sy_symm(), i1, i2);
+       return indexed(metric, symmetric2(), i1, i2);
 }
 
 ex lorentz_g(const ex & i1, const ex & i2, bool pos_sig)
 {
-       if (!is_ex_of_type(i1, varidx) || !is_ex_of_type(i2, varidx))
+       static ex metric_neg = (new minkmetric(false))->setflag(status_flags::dynallocated);
+       static ex metric_pos = (new minkmetric(true))->setflag(status_flags::dynallocated);
+
+       if (!is_a<varidx>(i1) || !is_a<varidx>(i2))
                throw(std::invalid_argument("indices of metric tensor must be of type varidx"));
 
-       return indexed(minkmetric(pos_sig), sy_symm(), i1, i2);
+       return indexed(pos_sig ? metric_pos : metric_neg, symmetric2(), i1, i2);
 }
 
 ex spinor_metric(const ex & i1, const ex & i2)
 {
-       if (!is_ex_of_type(i1, spinidx) || !is_ex_of_type(i2, spinidx))
+       static ex metric = (new spinmetric)->setflag(status_flags::dynallocated);
+
+       if (!is_a<spinidx>(i1) || !is_a<spinidx>(i2))
                throw(std::invalid_argument("indices of spinor metric must be of type spinidx"));
        if (!ex_to<idx>(i1).get_dim().is_equal(2) || !ex_to<idx>(i2).get_dim().is_equal(2))
                throw(std::runtime_error("index dimension for spinor metric must be 2"));
 
-       return indexed(spinmetric(), sy_anti(), i1, i2);
+       return indexed(metric, antisymmetric2(), i1, i2);
 }
 
 ex epsilon_tensor(const ex & i1, const ex & i2)
 {
-       if (!is_ex_of_type(i1, idx) || !is_ex_of_type(i2, idx))
+       static ex epsilon = (new tensepsilon)->setflag(status_flags::dynallocated);
+
+       if (!is_a<idx>(i1) || !is_a<idx>(i2))
                throw(std::invalid_argument("indices of epsilon tensor must be of type idx"));
 
        ex dim = ex_to<idx>(i1).get_dim();
@@ -569,12 +588,14 @@ ex epsilon_tensor(const ex & i1, const ex & i2)
        if (!ex_to<idx>(i1).get_dim().is_equal(_ex2))
                throw(std::runtime_error("index dimension of epsilon tensor must match number of indices"));
 
-       return indexed(tensepsilon(), sy_anti(), i1, i2);
+       return indexed(epsilon, antisymmetric2(), i1, i2);
 }
 
 ex epsilon_tensor(const ex & i1, const ex & i2, const ex & i3)
 {
-       if (!is_ex_of_type(i1, idx) || !is_ex_of_type(i2, idx) || !is_ex_of_type(i3, idx))
+       static ex epsilon = (new tensepsilon)->setflag(status_flags::dynallocated);
+
+       if (!is_a<idx>(i1) || !is_a<idx>(i2) || !is_a<idx>(i3))
                throw(std::invalid_argument("indices of epsilon tensor must be of type idx"));
 
        ex dim = ex_to<idx>(i1).get_dim();
@@ -583,12 +604,15 @@ ex epsilon_tensor(const ex & i1, const ex & i2, const ex & i3)
        if (!ex_to<idx>(i1).get_dim().is_equal(_ex3))
                throw(std::runtime_error("index dimension of epsilon tensor must match number of indices"));
 
-       return indexed(tensepsilon(), sy_anti(), i1, i2, i3);
+       return indexed(epsilon, antisymmetric3(), i1, i2, i3);
 }
 
 ex lorentz_eps(const ex & i1, const ex & i2, const ex & i3, const ex & i4, bool pos_sig)
 {
-       if (!is_ex_of_type(i1, varidx) || !is_ex_of_type(i2, varidx) || !is_ex_of_type(i3, varidx) || !is_ex_of_type(i4, varidx))
+       static ex epsilon_neg = (new tensepsilon(true, false))->setflag(status_flags::dynallocated);
+       static ex epsilon_pos = (new tensepsilon(true, true))->setflag(status_flags::dynallocated);
+
+       if (!is_a<varidx>(i1) || !is_a<varidx>(i2) || !is_a<varidx>(i3) || !is_a<varidx>(i4))
                throw(std::invalid_argument("indices of Lorentz epsilon tensor must be of type varidx"));
 
        ex dim = ex_to<idx>(i1).get_dim();
@@ -597,7 +621,7 @@ ex lorentz_eps(const ex & i1, const ex & i2, const ex & i3, const ex & i4, bool
        if (!ex_to<idx>(i1).get_dim().is_equal(_ex4))
                throw(std::runtime_error("index dimension of epsilon tensor must match number of indices"));
 
-       return indexed(tensepsilon(true, pos_sig), sy_anti(), i1, i2, i3, i4);
+       return indexed(pos_sig ? epsilon_pos : epsilon_neg, antisymmetric4(), i1, i2, i3, i4);
 }
 
 } // namespace GiNaC