Mentions légales du service
Skip to content
GitLab
Explore
Sign in
Primary navigation
Search or go to…
Project
spm
Manage
Activity
Members
Labels
Plan
Issues
Issue boards
Milestones
Code
Merge requests
Repository
Branches
Commits
Tags
Repository graph
Compare revisions
Deploy
Releases
Container Registry
Model registry
Monitor
Incidents
Service Desk
Analyze
Value stream analytics
Contributor analytics
Repository analytics
Model experiments
Help
Help
Support
GitLab documentation
Compare GitLab plans
Community forum
Contribute to GitLab
Provide feedback
Keyboard shortcuts
?
Snippets
Groups
Projects
Show more breadcrumbs
solverstack
spm
Commits
c1be92c7
Commit
c1be92c7
authored
7 years ago
by
Mathieu Faverge
Browse files
Options
Downloads
Patches
Plain Diff
Exploit the new sort to perform inplace conversion from IJV to CS[rc] formats
parent
6395bcdf
No related branches found
No related tags found
2 merge requests
!3
Feature/inplace ijv to csx recup
,
!2
Exploit the new sort to perform inplace conversion from IJV to CS[rc] formats
Changes
3
Hide whitespace changes
Inline
Side-by-side
Showing
3 changed files
src/z_spm.c
+20
-5
20 additions, 5 deletions
src/z_spm.c
src/z_spm_convert_to_csc.c
+21
-56
21 additions, 56 deletions
src/z_spm_convert_to_csc.c
src/z_spm_convert_to_csr.c
+26
-55
26 additions, 55 deletions
src/z_spm_convert_to_csr.c
with
67 additions
and
116 deletions
src/z_spm.c
+
20
−
5
View file @
c1be92c7
...
@@ -23,10 +23,12 @@
...
@@ -23,10 +23,12 @@
*
*
* @ingroup spm_dev_check
* @ingroup spm_dev_check
*
*
* @brief This routine sorts the subarray of edges of each vertex in a
* @brief This routine sorts the spm matrix.
* centralized spm stored in CSC or CSR format.
*
*
* Nothing is performed if IJV format is used.
* For the CSC and CSR formats, the subarray of edges for each vertex are sorted.
* For the IJV format, the edges are storted first by column indexes, and then
* by row indexes. To perform a sort first by row, second by column, please swap
* the colptr and rowptr of the structure before calling the subroutine.
*
*
* @warning This function should NOT be called if dof is greater than 1.
* @warning This function should NOT be called if dof is greater than 1.
*
*
...
@@ -44,13 +46,13 @@ z_spmSort( spmatrix_t *spm )
...
@@ -44,13 +46,13 @@ z_spmSort( spmatrix_t *spm )
spm_int_t
*
colptr
=
spm
->
colptr
;
spm_int_t
*
colptr
=
spm
->
colptr
;
spm_int_t
*
rowptr
=
spm
->
rowptr
;
spm_int_t
*
rowptr
=
spm
->
rowptr
;
spm_complex64_t
*
values
=
spm
->
values
;
spm_complex64_t
*
values
=
spm
->
values
;
void
*
sortptr
[
2
];
void
*
sortptr
[
3
];
spm_int_t
n
=
spm
->
n
;
spm_int_t
n
=
spm
->
n
;
spm_int_t
i
,
size
;
spm_int_t
i
,
size
;
(
void
)
sortptr
;
(
void
)
sortptr
;
if
(
spm
->
dof
>
1
){
if
(
spm
->
dof
>
1
){
fprintf
(
stderr
,
"z_spmSort: Number of dof (%d)
greater than
one not supported
\n
"
,
(
int
)
spm
->
dof
);
fprintf
(
stderr
,
"z_spmSort: Number of dof (%d)
different from
one not supported
\n
"
,
(
int
)
spm
->
dof
);
exit
(
1
);
exit
(
1
);
}
}
...
@@ -87,6 +89,19 @@ z_spmSort( spmatrix_t *spm )
...
@@ -87,6 +89,19 @@ z_spmSort( spmatrix_t *spm )
values
+=
size
;
values
+=
size
;
}
}
}
}
else
if
(
spm
->
fmttype
==
SpmIJV
)
{
size
=
spm
->
nnz
;
sortptr
[
0
]
=
colptr
;
sortptr
[
1
]
=
rowptr
;
#if defined(PRECISION_p)
spmIntSort2Asc2
(
sortptr
,
size
);
#else
sortptr
[
2
]
=
values
;
z_spmIntIntFltSortAsc
(
sortptr
,
size
);
#endif
}
}
}
/**
/**
...
...
This diff is collapsed.
Click to expand it.
src/z_spm_convert_to_csc.c
+
21
−
56
View file @
c1be92c7
...
@@ -40,12 +40,8 @@
...
@@ -40,12 +40,8 @@
int
int
z_spmConvertIJV2CSC
(
spmatrix_t
*
spm
)
z_spmConvertIJV2CSC
(
spmatrix_t
*
spm
)
{
{
#if !defined(PRECISION_p)
spm_int_t
*
newcol
,
*
oldcol
;
spm_complex64_t
*
navals
=
NULL
;
spm_int_t
i
,
tmp
,
baseval
,
total
;
spm_complex64_t
*
oavals
=
NULL
;
#endif
spm_int_t
*
spmptx
,
*
otmp
;
spm_int_t
i
,
j
,
tmp
,
baseval
,
total
;
spmatrix_t
oldspm
;
spmatrix_t
oldspm
;
/* Backup the input */
/* Backup the input */
...
@@ -56,66 +52,35 @@ z_spmConvertIJV2CSC( spmatrix_t *spm )
...
@@ -56,66 +52,35 @@ z_spmConvertIJV2CSC( spmatrix_t *spm )
*/
*/
baseval
=
spmFindBase
(
spm
);
baseval
=
spmFindBase
(
spm
);
/* Compute the new colptr */
/*
* Sort the IJV structure by column/row indexes
*/
z_spmSort
(
spm
);
/* Allocate and compute the new colptr */
spm
->
colptr
=
(
spm_int_t
*
)
calloc
(
spm
->
n
+
1
,
sizeof
(
spm_int_t
));
spm
->
colptr
=
(
spm_int_t
*
)
calloc
(
spm
->
n
+
1
,
sizeof
(
spm_int_t
));
/* Compute the number of edges per row */
/* Compute the number of edges per row */
spmptx
=
spm
->
colptr
-
baseval
;
newcol
=
spm
->
colptr
-
baseval
;
o
tmp
=
oldspm
.
colptr
;
o
ldcol
=
oldspm
.
colptr
;
for
(
i
=
0
;
i
<
spm
->
nnz
;
i
++
,
o
tmp
++
)
for
(
i
=
0
;
i
<
spm
->
nnz
;
i
++
,
o
ldcol
++
)
{
{
spmptx
[
*
otmp
]
++
;
newcol
[
*
oldcol
]
++
;
}
}
/*
Compute the indexes in C numbering for
the
f
ol
lowing sort
*/
/*
Update
the
c
ol
ptr
*/
total
=
0
;
total
=
baseval
;
spmptx
=
spm
->
colptr
;
newcol
=
spm
->
colptr
;
for
(
i
=
0
;
i
<
(
spm
->
n
+
1
);
i
++
,
spmptx
++
)
for
(
i
=
0
;
i
<
(
spm
->
n
+
1
);
i
++
,
newcol
++
)
{
{
tmp
=
*
spmptx
;
tmp
=
*
newcol
;
*
spmptx
=
total
;
*
newcol
=
total
;
total
+=
tmp
;
total
+=
tmp
;
}
}
assert
(
total
==
spm
->
nnz
);
assert
(
(
total
-
baseval
)
==
spm
->
nnz
);
/* Sort the rows and avals arrays by column */
spm
->
rowptr
=
malloc
(
spm
->
nnz
*
sizeof
(
spm_int_t
));
#if defined(PRECISION_p)
spm
->
values
=
NULL
;
#else
spm
->
values
=
malloc
(
spm
->
nnz
*
sizeof
(
spm_complex64_t
));
navals
=
(
spm_complex64_t
*
)(
spm
->
values
);
oavals
=
(
spm_complex64_t
*
)(
oldspm
.
values
);
#endif
for
(
j
=
0
;
j
<
spm
->
nnz
;
j
++
)
{
i
=
oldspm
.
colptr
[
j
]
-
baseval
;
spm
->
rowptr
[
spm
->
colptr
[
i
]
]
=
oldspm
.
rowptr
[
j
];
#if !defined(PRECISION_p)
navals
[
spm
->
colptr
[
i
]
]
=
oavals
[
j
];
#endif
(
spm
->
colptr
[
i
])
++
;
assert
(
spm
->
colptr
[
i
]
<=
spm
->
colptr
[
i
+
1
]
);
}
/* Rebuild the colptr with the correct baseval */
tmp
=
spm
->
colptr
[
0
];
spm
->
colptr
[
0
]
=
baseval
;
spmptx
=
spm
->
colptr
+
1
;
for
(
i
=
1
;
i
<
(
spm
->
n
+
1
);
i
++
,
spmptx
++
)
{
total
=
*
spmptx
;
*
spmptx
=
tmp
+
baseval
;
tmp
=
total
;
}
assert
(
spm
->
colptr
[
spm
->
n
]
==
(
spm
->
nnz
+
baseval
)
);
oldspm
.
rowptr
=
NULL
;
oldspm
.
values
=
NULL
;
spmExit
(
&
oldspm
);
spmExit
(
&
oldspm
);
spm
->
fmttype
=
SpmCSC
;
spm
->
fmttype
=
SpmCSC
;
...
...
This diff is collapsed.
Click to expand it.
src/z_spm_convert_to_csr.c
+
26
−
55
View file @
c1be92c7
...
@@ -13,6 +13,7 @@
...
@@ -13,6 +13,7 @@
* @date 2015-01-01
* @date 2015-01-01
*
*
* @precisions normal z -> c d s p
* @precisions normal z -> c d s p
*
**/
**/
#include
"common.h"
#include
"common.h"
#include
"z_spm.h"
#include
"z_spm.h"
...
@@ -123,12 +124,8 @@ z_spmConvertCSC2CSR( spmatrix_t *spm )
...
@@ -123,12 +124,8 @@ z_spmConvertCSC2CSR( spmatrix_t *spm )
int
int
z_spmConvertIJV2CSR
(
spmatrix_t
*
spm
)
z_spmConvertIJV2CSR
(
spmatrix_t
*
spm
)
{
{
#if !defined(PRECISION_p)
spm_int_t
*
newrow
,
*
oldrow
;
spm_complex64_t
*
navals
=
NULL
;
spm_int_t
i
,
tmp
,
baseval
,
total
;
spm_complex64_t
*
oavals
=
NULL
;
#endif
spm_int_t
*
spmptx
,
*
otmp
;
spm_int_t
i
,
j
,
tmp
,
baseval
,
total
;
spmatrix_t
oldspm
;
spmatrix_t
oldspm
;
/* Backup the input */
/* Backup the input */
...
@@ -139,66 +136,40 @@ z_spmConvertIJV2CSR( spmatrix_t *spm )
...
@@ -139,66 +136,40 @@ z_spmConvertIJV2CSR( spmatrix_t *spm )
*/
*/
baseval
=
spmFindBase
(
spm
);
baseval
=
spmFindBase
(
spm
);
/*
* Sort the IJV structure by column/row indexes
*/
newrow
=
spm
->
colptr
;
spm
->
colptr
=
spm
->
rowptr
;
spm
->
rowptr
=
newrow
;
z_spmSort
(
spm
);
spm
->
rowptr
=
spm
->
colptr
;
spm
->
colptr
=
newrow
;
/* Compute the new rowptr */
/* Compute the new rowptr */
spm
->
rowptr
=
(
spm_int_t
*
)
calloc
(
spm
->
n
+
1
,
sizeof
(
spm_int_t
));
spm
->
rowptr
=
(
spm_int_t
*
)
calloc
(
spm
->
n
+
1
,
sizeof
(
spm_int_t
));
/* Compute the number of edges per row */
/* Compute the number of edges per row */
spmptx
=
spm
->
rowptr
-
baseval
;
newrow
=
spm
->
rowptr
-
baseval
;
o
tmp
=
oldspm
.
rowptr
;
o
ldrow
=
oldspm
.
rowptr
;
for
(
i
=
0
;
i
<
spm
->
nnz
;
i
++
,
o
tmp
++
)
for
(
i
=
0
;
i
<
spm
->
nnz
;
i
++
,
o
ldrow
++
)
{
{
spmptx
[
*
otmp
]
++
;
newrow
[
*
oldrow
]
++
;
}
}
/*
Compu
te the
indexes in C numbering for the following sort
*/
/*
Upda
te the
rowptr
*/
total
=
0
;
total
=
baseval
;
spmptx
=
spm
->
rowptr
;
newrow
=
spm
->
rowptr
;
for
(
i
=
0
;
i
<
(
spm
->
n
+
1
);
i
++
,
spmptx
++
)
for
(
i
=
0
;
i
<
(
spm
->
n
+
1
);
i
++
,
newrow
++
)
{
{
tmp
=
*
spmptx
;
tmp
=
*
newrow
;
*
spmptx
=
total
;
*
newrow
=
total
;
total
+=
tmp
;
total
+=
tmp
;
}
}
assert
(
total
==
spm
->
nnz
);
assert
(
(
total
-
baseval
)
==
spm
->
nnz
);
/* Sort the colptr and avals arrays by rows */
spm
->
colptr
=
malloc
(
spm
->
nnz
*
sizeof
(
spm_int_t
));
#if defined(PRECISION_p)
spm
->
values
=
NULL
;
#else
spm
->
values
=
malloc
(
spm
->
nnz
*
sizeof
(
spm_complex64_t
));
navals
=
(
spm_complex64_t
*
)(
spm
->
values
);
oavals
=
(
spm_complex64_t
*
)(
oldspm
.
values
);
#endif
for
(
j
=
0
;
j
<
spm
->
nnz
;
j
++
)
{
i
=
oldspm
.
rowptr
[
j
]
-
baseval
;
spm
->
colptr
[
spm
->
rowptr
[
i
]
]
=
oldspm
.
colptr
[
j
];
#if !defined(PRECISION_p)
navals
[
spm
->
rowptr
[
i
]
]
=
oavals
[
j
];
#endif
(
spm
->
rowptr
[
i
])
++
;
assert
(
spm
->
rowptr
[
i
]
<=
spm
->
rowptr
[
i
+
1
]
);
}
/* Rebuild the rows (rowptr) with the correct baseval */
tmp
=
spm
->
rowptr
[
0
];
spm
->
rowptr
[
0
]
=
baseval
;
spmptx
=
spm
->
rowptr
+
1
;
for
(
i
=
1
;
i
<
(
spm
->
n
+
1
);
i
++
,
spmptx
++
)
{
total
=
*
spmptx
;
*
spmptx
=
tmp
+
baseval
;
tmp
=
total
;
}
assert
(
spm
->
rowptr
[
spm
->
n
]
==
(
spm
->
nnz
+
baseval
)
);
oldspm
.
colptr
=
NULL
;
oldspm
.
values
=
NULL
;
spmExit
(
&
oldspm
);
spmExit
(
&
oldspm
);
spm
->
fmttype
=
SpmCSR
;
spm
->
fmttype
=
SpmCSR
;
...
...
This diff is collapsed.
Click to expand it.
Preview
0%
Loading
Try again
or
attach a new file
.
Cancel
You are about to add
0
people
to the discussion. Proceed with caution.
Finish editing this message first!
Save comment
Cancel
Please
register
or
sign in
to comment