<html><head><meta http-equiv="Content-Type" content="text/html charset=iso-8859-1"></head><body style="word-wrap: break-word; -webkit-nbsp-mode: space; -webkit-line-break: after-white-space; ">I found a similar problem solving Laplace's equation: using MatCreate took forever, whereas using MatCreateAIJ instead the time for MatSetValues was essentially negligible. I set up MatCreate with a call to MatSetSizes.<div>...Peter</div><div><br><div><div>On May 17, 2013, at 2:50 PM, Matthew Knepley <<a href="mailto:knepley@gmail.com">knepley@gmail.com</a>> wrote:</div><br class="Apple-interchange-newline"><blockquote type="cite"><div dir="ltr">On Fri, May 17, 2013 at 3:40 PM, Joon hee Choi <span dir="ltr"><<a href="mailto:choi240@purdue.edu" target="_blank">choi240@purdue.edu</a>></span> wrote:<br><div class="gmail_extra"><div class="gmail_quote">
<blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left-width:1px;border-left-color:rgb(204,204,204);border-left-style:solid;padding-left:1ex">My petsc version is 3.3.5.<br>
And I made nnz using the number of each row. When I wrote the matrix using MatSetValues, I wrote from low row and column to high row and column. However, when I changed the order writing the matrix, it took much more time (up to 3hrs). So I am concerned that the slowness is because the big size of matrix (2.6*10^7, 1.248*10^15) is related to reading from or writing to the cache, ram, or hard.<br>
</blockquote><div><br></div><div style="">No, it is bad preallocation. This should be simple to fix. First, turn on errors</div><div style=""><br></div><div style="">  MatSetOption(<span style="font-family: Times; font-size: medium; ">MAT_NEW_NONZERO_ALLOCATION_ERR</span>, PETSC_TRUE)</div>
<div style=""><br></div><div style="">and second run with -info.</div><div style=""><br></div><div style="">  Thanks,</div><div style=""><br></div><div style="">     Matt</div><div> </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left-width:1px;border-left-color:rgb(204,204,204);border-left-style:solid;padding-left:1ex">

Anyway, the following code is the part that I read the data from a file and set up tuples and nnz:<br>
<br>
        FILE *fp = fopen("data.txt", "r");<br>
        while (fscanf(fp, "%d %d %d %d", &x, &y, &z, &v) == 4)<br>
        {<br>
                tups.push_back(std::tr1::make_tuple (x-1, z-1, y-1, v));<br>
                nzrow[i-1] += 1;<br>
                if (x > X) X = x;<br>
                if (y > Y) Y = y;<br>
                if (z > Z) Z = z;<br>
        }<br>
        fclose(fp);<br>
        PetscMalloc(X*sizeof(PetscInt), &nnz);<br>
        memset(nnz, 0, X);<br>
        for (itnz=nzrow.begin(); itnz!=nzrow.end(); ++itnz) {<br>
                nnz[itnz->first] = itnz->second;<br>
        }<br>
        sort(tups.begin(), tups.end());<br>
<br>
If my code is wrong, then please let me know.<br>
<br>
Thank you,<br>
Joon<br>
<br>
----- Original Message -----<br>
From: "Jed Brown" <<a href="mailto:jedbrown@mcs.anl.gov">jedbrown@mcs.anl.gov</a>><br>
To: "Joon hee Choi" <<a href="mailto:choi240@purdue.edu">choi240@purdue.edu</a>><br>
Cc: <a href="mailto:petsc-users@mcs.anl.gov">petsc-users@mcs.anl.gov</a><br>
Sent: Friday, May 17, 2013 3:24:09 PM<br>
Subject: Re: [petsc-users] 3-dimension to matrix<br>
<br>
Joon hee Choi <<a href="mailto:choi240@purdue.edu">choi240@purdue.edu</a>> writes:<br>
<br>
> Thank you for your fast reply. Your last comments looks like my first<br>
> code. The full size of the matrix is (X, Y*Z). Also, I used SEQAIJ as<br>
> matrix type and (X, Y) as block size. Also, I implemented<br>
> SeqAIJpreallocating with nnz. Nevertheless, it was very slow. The<br>
> Matrix-Set-Up part of the first code is as follows:<br>
><br>
>         sort(tups.begin(), tups.end());<br>
>         MatCreate(PETSC_COMM_SELF, &A);<br>
>         MatSetType(A, MATSEQAIJ);<br>
>         MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, X, Y*Z);<br>
>         MatSetBlockSizes(A, X, Y);<br>
>         MatSeqAIJSetPreallocation(A, PETSC_DEFAULT, nnz);<br>
><br>
>         sz = tups.size();<br>
>         for (i=0; i<sz; i++) {<br>
>                 x = std::tr1::get<0>(tups[i]);<br>
>                 y = std::tr1::get<2>(tups[i]) + std::tr1::get<1>(tups[i])*Y;<br>
>                 val = std::tr1::get<3>(tups[i]);<br>
>                 MatSetValues(A, 1, &x, 1, &y, &val, INSERT_VALUES);<br>
>         }<br>
>         MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);<br>
>         MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);<br>
><br>
><br>
> I used tuple (x, z, y, value) and vector of c++. I didn't get any<br>
> errors from this code. However, it took about 9 minutes in this<br>
> part.<br>
<br>
What version of PETSc?  Your preallocation was almost certainly not<br>
sufficient.<br>
</blockquote></div><br><br clear="all"><div><br></div>-- <br>What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.<br>
-- Norbert Wiener
</div></div>
</blockquote></div><br><div apple-content-edited="true">
<span class="Apple-style-span" style="border-collapse: separate; color: rgb(0, 0, 0); font-family: Helvetica; font-style: normal; font-variant: normal; font-weight: normal; letter-spacing: normal; line-height: normal; orphans: 2; text-align: -webkit-auto; text-indent: 0px; text-transform: none; white-space: normal; widows: 2; word-spacing: 0px; -webkit-border-horizontal-spacing: 0px; -webkit-border-vertical-spacing: 0px; -webkit-text-decorations-in-effect: none; -webkit-text-size-adjust: auto; -webkit-text-stroke-width: 0px; font-size: medium; "><span class="Apple-style-span" style="border-collapse: separate; color: rgb(0, 0, 0); font-family: Helvetica; font-style: normal; font-variant: normal; font-weight: normal; letter-spacing: normal; line-height: normal; orphans: 2; text-align: -webkit-auto; text-indent: 0px; text-transform: none; white-space: normal; widows: 2; word-spacing: 0px; -webkit-border-horizontal-spacing: 0px; -webkit-border-vertical-spacing: 0px; -webkit-text-decorations-in-effect: none; -webkit-text-size-adjust: auto; -webkit-text-stroke-width: 0px; font-size: medium; "><div style="word-wrap: break-word; -webkit-nbsp-mode: space; -webkit-line-break: after-white-space; "><span class="Apple-style-span" style="border-collapse: separate; color: rgb(0, 0, 0); font-family: Helvetica; font-style: normal; font-variant: normal; font-weight: normal; letter-spacing: normal; line-height: normal; orphans: 2; text-align: -webkit-auto; text-indent: 0px; text-transform: none; white-space: normal; widows: 2; word-spacing: 0px; -webkit-border-horizontal-spacing: 0px; -webkit-border-vertical-spacing: 0px; -webkit-text-decorations-in-effect: none; -webkit-text-size-adjust: auto; -webkit-text-stroke-width: 0px; font-size: medium; "><div style="word-wrap: break-word; -webkit-nbsp-mode: space; -webkit-line-break: after-white-space; "><div>________________</div><div>Peter Lichtner</div><div>Santa Fe, NM 87507</div><div>(505) 692-4029 (c)</div><div>OFM Research/LANL Guest Scientist</div></div></span></div></span></span>
</div>
<br></div></body></html>