<div dir="ltr"><div class="gmail_quote"><div dir="ltr">On Mon, Oct 15, 2018 at 10:42 AM Yingjie Wu <<a href="mailto:yjwu16@gmail.com">yjwu16@gmail.com</a>> wrote:<br></div><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex"><div dir="ltr"><div dir="ltr"><div dir="ltr"><div dir="ltr"><div>Dear Petsc developer:           </div><div>Hi,</div><div>Thank you very much for your previous reply. </div><div>I recently wanted to modify my program to parallel version, and encountered some problems in modifying it.            </div><div>1. There are functions (read matrix) in the program that reads files,will they affect my parallelism? </div></div></div></div></div></blockquote><div><br></div><div>No, MatLoad works in parallel. This will work unchanged if you want the default layout. If you want a special partition,</div><div>you must call MatSetSizes() after MatCreate().</div><div> </div><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex"><div dir="ltr"><div dir="ltr"><div dir="ltr"><div dir="ltr"><div> The codes are as follows:           </div></div></div><blockquote style="margin:0px 0px 0px 40px;border:none;padding:0px"><div><div><div><font size="1">ierr = PetscViewerBinaryOpen (PETSC_COMM_WORLD, file, FILE_MODE_READ, &viewer); CHKERRQ (ierr);            </font></div></div></div><div><div><div><font size="1">ierr = MatCreate (PETSC_COMM_WORLD, &A1); CHKERRQ (ierr);            </font></div></div></div><div><div><div><font size="1">ierr = MatSetFromOptions (A1); CHKERRQ (ierr);            </font></div></div></div><div><div><div><font size="1">ierr = MatCreate (PETSC_COMM_WORLD, &A2); CHKERRQ (ierr);            </font></div></div></div><div><div><div><font size="1">ierr = MatSetFromOptions (A2); CHKERRQ (ierr);            </font></div></div></div><div><div><div><font size="1">ierr = MatLoad (A1, viewer); CHKERRQ (ierr);            </font></div></div></div><div><div><div><font size="1">ierr = MatLoad (A2, viewer); CHKERRQ (ierr);            </font></div></div></div><div><div><div><font size="1">ierr = PetscViewerDestroy (&viewer); CHKERRQ (ierr);  </font>          </div></div></div></blockquote><div dir="ltr"><div dir="ltr"><div>I read two matrix information from a binary file and wanted to use it on each processor (in fact, my goal was to use these two matrices to give initial values to the two field variables). The program can run in serial time. Now I want to change it to parallel program. What do I need to modify?            </div><div>2. Following the last question, I used the following code in giving initial value procedure FormInitialGuess():            </div></div></div></div><blockquote style="margin:0 0 0 40px;border:none;padding:0px"><div><div><div><div><font size="1">ierr = MatSeqAIJGetArray (A1, &initial_phi1); CHKERRQ (ierr);            </font></div></div></div></div><div><div><div><div><font size="1">ierr = MatSeqAIJGetArray (A2, &initial_phi2); CHKERRQ (ierr); </font>           </div></div></div></div></blockquote><div dir="ltr"><div dir="ltr"><div dir="ltr"><div>I found this function on manualpages, and I felt that it could fulfill my need to represent the elements of the matrix in arrays to give field variables an initial value in each grid. The matrix A1 and A2 above are actually the matrix information that was read from the file before. Now I want to modify it as a parallel program. How do I use matrix information to give initial values in parallel? (In p<span style="white-space:pre-wrap">rogram, field variables are implemented through DM because parallel computing and Ghost Value are supported</span>)</div></div></div></div></div></blockquote><div><br></div><div>I do not understand the use of matrices to initialize field values. Usually field values are stored on Vec objects, and this is</div><div>the philosophy of DM objects.</div><div><br></div><div>  Thanks,</div><div><br></div><div>     Matt</div><div> </div><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex"><div dir="ltr"><div dir="ltr"><div dir="ltr"><div dir="ltr"><div>Thanks,</div><div>Yingjie</div></div></div></div></div>
</blockquote></div><br clear="all"><div><br></div>-- <br><div dir="ltr" class="gmail_signature" data-smartmail="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>