<div dir="ltr"><div dir="ltr">On Fri, Dec 4, 2020 at 6:51 AM Roland Richter <<a href="mailto:roland.richter@ntnu.no">roland.richter@ntnu.no</a>> wrote:<br></div><div class="gmail_quote"><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex">
  

    
  
  <div>
    <p>Hei,</p>
    <p>is it possible to fill a dense distributed matrix from existing
      data in rank 0, distribute it to all involved nodes and retrieve
      it afterwards, such that it can be stored in a single matrix in
      rank 0 again? The background behind this question is the following
      thought process:</p>
    <ul>
      <li>I generate the matrix locally on node 0</li></ul></div></blockquote><div>I do not understand why you would ever want to do this. <br></div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div><ul>
      <li>I distribute it to all involved nodes</li></ul></div></blockquote><div>Just create your distributed matrix, but set the values from rank 0. <br></div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div><ul>
      <li>I use the matrix for data processing/process the matrix data
        on all nodes</li></ul></div></blockquote><div>I assume you change the matrix here. </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div><ul>
      <li>For storing the data I need it back on node 0 in a single
        matrix</li></ul></div></blockquote><div>You can use <a href="https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatCreateRedundantMatrix.html">https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatCreateRedundantMatrix.html</a></div><div>to easily create serial matrices from it. You can throw away any you do not want.</div><div><br></div><div>  Thanks,</div><div><br></div><div>     Matt </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div><ul>
    </ul>
    <p>For testing I wrote the following code (and using armadillo for
      generating the initial matrices within arma-namespace):</p>
    <p><i>    Mat C, F;</i><i><br>
      </i><i>    Vec x, y, z;</i><i><br>
      </i><i>    PetscViewer viewer;</i><i><br>
      </i><i>    PetscMPIInt rank, size;</i><i><br>
      </i><i>    PetscInitialize(&argc, &args, (char*) 0, help);</i><i><br>
      </i><i><br>
      </i><i>    MPI_Comm_size(PETSC_COMM_WORLD, &size);</i><i><br>
      </i><i>    MPI_Comm_rank(PETSC_COMM_WORLD, &rank);</i><i><br>
      </i><i><br>
      </i><i>    PetscPrintf(PETSC_COMM_WORLD,"Number of processors =
        %d, rank = %d\n", size, rank);</i><i><br>
      </i><i></i><i><br>
      </i><i>    PetscViewerCreate(PETSC_COMM_WORLD, &viewer);</i><i><br>
      </i><i>    PetscViewerSetType(viewer, PETSCVIEWERASCII);</i><i><br>
      </i><i>    arma::cx_mat local_mat, local_zero_mat;</i><i><br>
      </i><i>    const size_t matrix_size = 5;</i><i><br>
      </i><i>    if(rank == 0) {</i><i><br>
      </i><i>        local_mat =
        arma::randu<arma::cx_mat>(matrix_size, matrix_size);</i><i><br>
      </i><i>        local_zero_mat =
        arma::zeros<arma::cx_mat>(matrix_size, matrix_size);</i><i><br>
      </i><i>        MatCreateDense(PETSC_COMM_WORLD, matrix_size,
        matrix_size, PETSC_DECIDE, PETSC_DECIDE, local_mat.memptr(),
        &C);</i><i><br>
      </i><i>        MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY);</i><i><br>
      </i><i>        MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY);</i><i><br>
      </i><i>    }</i><i><br>
      </i><i>    MatCreateDense(PETSC_COMM_WORLD, matrix_size,
        matrix_size, PETSC_DECIDE, PETSC_DECIDE, NULL, &F);</i><i><br>
      </i><i>    if(rank == 0)</i><i><br>
      </i><i>        MatCopy(C, F, SAME_NONZERO_PATTERN);</i><i><br>
      </i><i>    MatAssemblyBegin(F, MAT_FINAL_ASSEMBLY);</i><i><br>
      </i><i>    MatAssemblyEnd(F, MAT_FINAL_ASSEMBLY);</i><i><br>
      </i><i>    if(rank == 0) {</i><i><br>
      </i><i>        MatView(C, viewer);</i><i><br>
      </i><i>        std::cout << local_mat << '\n';</i><i><br>
      </i><i>        std::cout << local_zero_mat << '\n';</i><i><br>
      </i><i>    }</i><i><br>
      </i><i>    //How can I retrieve the distributed data from F to C
        again when running on multiple nodes?</i><i><br>
      </i><i><br>
      </i><i>    if(rank == 0) {</i><i><br>
      </i><i>        MatDestroy(&C);</i><i><br>
      </i><i>    }</i><i><br>
      </i><i>    MatDestroy(&F);</i><i><br>
      </i><i>    PetscViewerDestroy(&viewer);</i></p>
    <p>Running the code with mpirun -n 1 works fine, but when running it
      with mpirun -n 2 it just stops after <i>PetscPrintf() </i>and
      hangs.</p>
    <p>Moreover, how can I retrieve the data from matrix F afterwards?
      For vectors I have the function/class VecScatter, but can I apply
      this for matrices (mpidense) too?</p>
    <p>Thanks!</p>
    <p>Cheers,</p>
    <p>Roland<br>
    </p>
  </div>

</blockquote></div><br clear="all"><div><br></div>-- <br><div dir="ltr" class="gmail_signature"><div dir="ltr"><div><div dir="ltr"><div><div dir="ltr"><div>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><br></div><div><a href="http://www.cse.buffalo.edu/~knepley/" target="_blank">https://www.cse.buffalo.edu/~knepley/</a><br></div></div></div></div></div></div></div></div>