blob: b83fb82ba6027a73e7441d12c00ffd1b3d03d169 [file] [log] [blame]
Narayan Kamathc981c482012-11-02 10:59:05 +00001// This file is part of Eigen, a lightweight C++ template library
2// for linear algebra.
3//
4// Copyright (C) 2011 Gael Guennebaud <gael.guennebaud@inria.fr>
5//
6// This Source Code Form is subject to the terms of the Mozilla
7// Public License v. 2.0. If a copy of the MPL was not distributed
8// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
9
10
11// Various sanity tests with exceptions:
12// - no memory leak when a custom scalar type trow an exceptions
13// - todo: complete the list of tests!
14
15#define EIGEN_STACK_ALLOCATION_LIMIT 100000000
16
17#include "main.h"
18
19struct my_exception
20{
21 my_exception() {}
22 ~my_exception() {}
23};
24
25class ScalarWithExceptions
26{
27 public:
28 ScalarWithExceptions() { init(); }
29 ScalarWithExceptions(const float& _v) { init(); *v = _v; }
30 ScalarWithExceptions(const ScalarWithExceptions& other) { init(); *v = *(other.v); }
31 ~ScalarWithExceptions() {
32 delete v;
33 instances--;
34 }
35
36 void init() {
37 v = new float;
38 instances++;
39 }
40
41 ScalarWithExceptions operator+(const ScalarWithExceptions& other) const
42 {
43 countdown--;
44 if(countdown<=0)
45 throw my_exception();
46 return ScalarWithExceptions(*v+*other.v);
47 }
48
49 ScalarWithExceptions operator-(const ScalarWithExceptions& other) const
50 { return ScalarWithExceptions(*v-*other.v); }
51
52 ScalarWithExceptions operator*(const ScalarWithExceptions& other) const
53 { return ScalarWithExceptions((*v)*(*other.v)); }
54
55 ScalarWithExceptions& operator+=(const ScalarWithExceptions& other)
56 { *v+=*other.v; return *this; }
57 ScalarWithExceptions& operator-=(const ScalarWithExceptions& other)
58 { *v-=*other.v; return *this; }
59 ScalarWithExceptions& operator=(const ScalarWithExceptions& other)
60 { *v = *(other.v); return *this; }
61
62 bool operator==(const ScalarWithExceptions& other) const
63 { return *v==*other.v; }
64 bool operator!=(const ScalarWithExceptions& other) const
65 { return *v!=*other.v; }
66
67 float* v;
68 static int instances;
69 static int countdown;
70};
71
Carlos Hernandez7faaa9f2014-08-05 17:53:32 -070072ScalarWithExceptions real(const ScalarWithExceptions &x) { return x; }
73ScalarWithExceptions imag(const ScalarWithExceptions & ) { return 0; }
74ScalarWithExceptions conj(const ScalarWithExceptions &x) { return x; }
75
Narayan Kamathc981c482012-11-02 10:59:05 +000076int ScalarWithExceptions::instances = 0;
77int ScalarWithExceptions::countdown = 0;
78
79
80#define CHECK_MEMLEAK(OP) { \
81 ScalarWithExceptions::countdown = 100; \
82 int before = ScalarWithExceptions::instances; \
83 bool exception_thrown = false; \
84 try { OP; } \
85 catch (my_exception) { \
86 exception_thrown = true; \
87 VERIFY(ScalarWithExceptions::instances==before && "memory leak detected in " && EIGEN_MAKESTRING(OP)); \
88 } \
89 VERIFY(exception_thrown && " no exception thrown in " && EIGEN_MAKESTRING(OP)); \
90 }
91
92void memoryleak()
93{
94 typedef Eigen::Matrix<ScalarWithExceptions,Dynamic,1> VectorType;
95 typedef Eigen::Matrix<ScalarWithExceptions,Dynamic,Dynamic> MatrixType;
96
97 {
98 int n = 50;
99 VectorType v0(n), v1(n);
100 MatrixType m0(n,n), m1(n,n), m2(n,n);
101 v0.setOnes(); v1.setOnes();
102 m0.setOnes(); m1.setOnes(); m2.setOnes();
103 CHECK_MEMLEAK(v0 = m0 * m1 * v1);
104 CHECK_MEMLEAK(m2 = m0 * m1 * m2);
105 CHECK_MEMLEAK((v0+v1).dot(v0+v1));
106 }
107 VERIFY(ScalarWithExceptions::instances==0 && "global memory leak detected in " && EIGEN_MAKESTRING(OP)); \
108}
109
110void test_exceptions()
111{
112 CALL_SUBTEST( memoryleak() );
113}