20 #ifndef EIGEN_SIMPLICIAL_CHOLESKY_IMPL_H
21 #define EIGEN_SIMPLICIAL_CHOLESKY_IMPL_H
30 template <
typename Scalar,
typename StorageIndex>
35 static constexpr StorageIndex
kEmpty = -1;
41 #ifndef EIGEN_NO_DEBUG
53 StorageIndex
back()
const {
58 #ifndef EIGEN_NO_DEBUG
75 StorageIndex
find(StorageIndex u)
const {
77 while (
m_set[u] != u) {
88 StorageIndex next =
m_set[u];
101 for (StorageIndex
j = 1;
j <
size;
j++) {
103 StorageIndex
i = it.index();
104 if (
i <
j) outerIndex[
i + 1]++;
107 std::partial_sum(outerIndex, outerIndex +
size + 1, outerIndex);
112 StorageIndex* innerIndex, StorageIndex*
tmp) {
115 for (StorageIndex
j = 1;
j <
size;
j++) {
117 StorageIndex
i = it.index();
119 StorageIndex
b = outerIndex[
i] +
tmp[
i];
140 for (StorageIndex
j = 1;
j <
size;
j++) {
142 StorageIndex
i = it.index();
144 StorageIndex
r = ancestor.
find(
i);
145 if (
r !=
j) parent[
r] =
j;
153 static void calc_lineage(
const StorageIndex
size,
const StorageIndex* parent, StorageIndex* firstChild,
154 StorageIndex* firstSibling) {
158 for (StorageIndex
j = 0;
j <
size;
j++) {
159 StorageIndex
p = parent[
j];
161 StorageIndex
c = firstChild[
p];
165 while (firstSibling[
c] !=
kEmpty)
c = firstSibling[
c];
172 static void calc_post(
const StorageIndex
size,
const StorageIndex* parent, StorageIndex* firstChild,
173 const StorageIndex* firstSibling, StorageIndex* post, StorageIndex*
dfs) {
175 for (StorageIndex
j = 0;
j <
size;
j++) {
176 if (parent[
j] !=
kEmpty)
continue;
180 while (!dfs_stack.
empty()) {
181 StorageIndex
i = dfs_stack.
back();
182 StorageIndex
c = firstChild[
i];
189 firstChild[
i] = firstSibling[
c];
203 static void calc_colcount(
const StorageIndex
size,
const StorageIndex* hadjOuter,
const StorageIndex* hadjInner,
204 const StorageIndex* parent, StorageIndex* prevLeaf, StorageIndex*
tmp,
205 const StorageIndex* post, StorageIndex* nonZerosPerCol,
bool doLDLT) {
207 std::fill_n(nonZerosPerCol,
size, 1);
208 for (StorageIndex
j = 0;
j <
size;
j++) {
209 StorageIndex
p = parent[
j];
211 if (
p !=
kEmpty) nonZerosPerCol[
p] = 0;
214 DisjointSet parentSet(
tmp,
size);
218 for (StorageIndex j_ = 0; j_ <
size; j_++) {
219 StorageIndex
j = post[j_];
220 nonZerosPerCol[
j] += hadjOuter[
j + 1] - hadjOuter[
j];
221 for (StorageIndex
k = hadjOuter[
j];
k < hadjOuter[
j + 1];
k++) {
222 StorageIndex
i = hadjInner[
k];
227 parentSet.compress(prev,
q);
232 StorageIndex
p = parent[
j];
233 if (
p !=
kEmpty) parentSet.compress(
j,
p);
236 for (StorageIndex
j = 0;
j <
size;
j++) {
237 StorageIndex
p = parent[
j];
238 if (
p !=
kEmpty) nonZerosPerCol[
p] += nonZerosPerCol[
j] - 1;
239 if (doLDLT) nonZerosPerCol[
j]--;
246 std::partial_sum(nonZerosPerCol, nonZerosPerCol +
size,
L.outerIndexPtr() + 1);
247 L.resizeNonZeros(
L.outerIndexPtr()[
size]);
252 VectorI& workSpace,
bool doLDLT) {
257 StorageIndex* tmp1 = workSpace.
data();
258 StorageIndex* tmp2 = workSpace.
data() +
size;
259 StorageIndex* tmp3 = workSpace.
data() + 2 *
size;
260 StorageIndex* tmp4 = workSpace.
data() + 3 *
size;
263 StorageIndex* hadj_outer =
L.outerIndexPtr();
279 template <
typename Derived>
288 m_isInitialized =
true;
290 m_analysisIsOk =
true;
291 m_factorizationIsOk =
false;
294 template <
typename Derived>
295 template <
bool DoLDLT,
bool NonHermitian>
300 eigen_assert(m_analysisIsOk &&
"You must first call analyzePattern()");
315 m_diag.resize(DoLDLT ?
size : 0);
322 nonZerosPerCol[
k] = 0;
326 y[
i] += getSymm(it.value());
328 for (len = 0; tags[
i] !=
k;
i = m_parent[
i]) {
332 while (len > 0) pattern[--top] = pattern[--len];
339 getDiag(
y[
k]) * m_shiftScale + m_shiftOffset;
341 for (; top <
size; ++top) {
349 l_ki = yi / getDiag(m_diag[
i]);
351 yi = l_ki = yi /
Lx[Lp[
i]];
353 Index p2 = Lp[
i] + nonZerosPerCol[
i];
355 for (
p = Lp[
i] + (DoLDLT ? 0 : 1);
p < p2; ++
p)
y[Li[
p]] -= getSymm(
Lx[
p]) * yi;
356 d -= getDiag(l_ki * getSymm(yi));
368 Index p = Lp[
k] + nonZerosPerCol[
k]++;
379 m_factorizationIsOk =
true;
AnnoyingScalar sqrt(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:134
Array< int, Dynamic, 1 > v
Definition: Array_initializer_list_vector_cxx11.cpp:1
int i
Definition: BiCGSTAB_step_by_step.cpp:9
MatrixXd L
Definition: LLT_example.cpp:6
#define eigen_assert(x)
Definition: Macros.h:910
int data[]
Definition: Map_placement_new.cpp:1
#define ei_declare_aligned_stack_constructed_variable(TYPE, NAME, SIZE, BUFFER)
Definition: Memory.h:806
float * p
Definition: Tutorial_Map_using.cpp:9
Scalar Scalar int size
Definition: benchVecAdd.cpp:17
Scalar * b
Definition: benchVecAdd.cpp:17
SCALAR Scalar
Definition: bench_gemm.cpp:45
NumTraits< Scalar >::Real RealScalar
Definition: bench_gemm.cpp:46
constexpr EIGEN_DEVICE_FUNC const Scalar * data() const
Definition: PlainObjectBase.h:273
EIGEN_DEVICE_FUNC constexpr EIGEN_STRONG_INLINE void resize(Index rows, Index cols)
Definition: PlainObjectBase.h:294
MatrixType::Scalar Scalar
Definition: SimplicialCholesky.h:59
MatrixType::StorageIndex StorageIndex
Definition: SimplicialCholesky.h:62
internal::traits< Derived >::DiagonalScalar DiagonalScalar
Definition: SimplicialCholesky.h:61
void analyzePattern_preordered(const CholMatrixType &a, bool doLDLT)
Definition: SimplicialCholesky_impl.h:280
void factorize_preordered(const CholMatrixType &a)
Definition: SimplicialCholesky_impl.h:296
Base::InnerIterator InnerIterator
Definition: SparseMatrix.h:138
float real
Definition: datatypes.h:10
@ NumericalIssue
Definition: Constants.h:442
@ Success
Definition: Constants.h:440
Scalar * y
Definition: level1_cplx_impl.h:128
const Scalar * a
Definition: level2_cplx_impl.h:32
Scalar * ap
Definition: level2_cplx_impl.h:161
char char char int int * k
Definition: level2_impl.h:374
Eigen::Matrix< Scalar, Dynamic, Dynamic, ColMajor > tmp
Definition: level3_impl.h:365
EIGEN_DEVICE_FUNC const Scalar & q
Definition: SpecialFunctionsImpl.h:2019
Namespace containing all symbols from the Eigen library.
Definition: bench_norm.cpp:70
squared absolute value
Definition: GlobalFunctions.h:87
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:83
double Lx
Length of domain in x direction.
Definition: periodic_load.cc:55
void dfs(Node *node_pt, Node *parent_node_pt, std::map< Node *, Vector< std::pair< NodeNodeConstraintElementType *, Node * >>> &tree, std::map< Node *, Colour > &colours, Vector< NodeNodeConstraintElementType * > &pruned_branches)
Definition: mortaring_helpers.h:320
int c
Definition: calibrate.py:100
Definition: Eigen_Colamd.h:49
void set(Container &c, Position position, const Value &value)
Definition: stdlist_overload.cpp:36
Definition: SimplicialCholesky_impl.h:71
DisjointSet(StorageIndex *set, StorageIndex size)
Definition: SimplicialCholesky_impl.h:73
void compress(StorageIndex u, StorageIndex v)
Definition: SimplicialCholesky_impl.h:84
StorageIndex find(StorageIndex u) const
Definition: SimplicialCholesky_impl.h:75
StorageIndex * m_set
Definition: SimplicialCholesky_impl.h:72
Definition: SimplicialCholesky_impl.h:38
Index size() const
Definition: SimplicialCholesky_impl.h:52
Index m_size
Definition: SimplicialCholesky_impl.h:40
Stack(StorageIndex *data, StorageIndex size, StorageIndex maxSize)
Definition: SimplicialCholesky_impl.h:43
StorageIndex * m_data
Definition: SimplicialCholesky_impl.h:39
bool empty() const
Definition: SimplicialCholesky_impl.h:51
void push(const StorageIndex &value)
Definition: SimplicialCholesky_impl.h:57
void pop()
Definition: SimplicialCholesky_impl.h:64
const Index m_maxSize
Definition: SimplicialCholesky_impl.h:42
StorageIndex back() const
Definition: SimplicialCholesky_impl.h:53
Definition: SimplicialCholesky_impl.h:31
static void calc_hadj_outer(const StorageIndex size, const CholMatrixType &ap, StorageIndex *outerIndex)
Definition: SimplicialCholesky_impl.h:100
static constexpr StorageIndex kEmpty
Definition: SimplicialCholesky_impl.h:35
static void calc_etree(const StorageIndex size, const CholMatrixType &ap, StorageIndex *parent, StorageIndex *tmp)
Definition: SimplicialCholesky_impl.h:135
static void run(const StorageIndex size, const CholMatrixType &ap, CholMatrixType &L, VectorI &parent, VectorI &workSpace, bool doLDLT)
Definition: SimplicialCholesky_impl.h:251
static void calc_colcount(const StorageIndex size, const StorageIndex *hadjOuter, const StorageIndex *hadjInner, const StorageIndex *parent, StorageIndex *prevLeaf, StorageIndex *tmp, const StorageIndex *post, StorageIndex *nonZerosPerCol, bool doLDLT)
Definition: SimplicialCholesky_impl.h:203
static void init_matrix(const StorageIndex size, const StorageIndex *nonZerosPerCol, CholMatrixType &L)
Definition: SimplicialCholesky_impl.h:244
static void calc_hadj_inner(const StorageIndex size, const CholMatrixType &ap, const StorageIndex *outerIndex, StorageIndex *innerIndex, StorageIndex *tmp)
Definition: SimplicialCholesky_impl.h:111
typename CholMatrixType::InnerIterator InnerIterator
Definition: SimplicialCholesky_impl.h:33
static void calc_post(const StorageIndex size, const StorageIndex *parent, StorageIndex *firstChild, const StorageIndex *firstSibling, StorageIndex *post, StorageIndex *dfs)
Definition: SimplicialCholesky_impl.h:172
static void calc_lineage(const StorageIndex size, const StorageIndex *parent, StorageIndex *firstChild, StorageIndex *firstSibling)
Definition: SimplicialCholesky_impl.h:153
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2
void run(const string &dir_name, LinearSolver *linear_solver_pt, const unsigned nel_1d, bool mess_up_order)
Definition: two_d_poisson_compare_solvers.cc:317