<html><head><meta http-equiv="Content-Type" content="text/html; charset=us-ascii"></head><body style="word-wrap: break-word; -webkit-nbsp-mode: space; line-break: after-white-space;" class=""><br class=""><div><br class=""><blockquote type="cite" class=""><div class="">Am 27.01.2021 um 17:30 schrieb Matthew Knepley <<a href="mailto:knepley@gmail.com" class="">knepley@gmail.com</a>>:</div><br class="Apple-interchange-newline"><div class=""><div dir="ltr" style="caret-color: rgb(0, 0, 0); font-family: Helvetica; font-size: 12px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none;" class=""><div dir="ltr" class="">On Wed, Jan 27, 2021 at 10:51 AM Viet H.Q.H. <<a href="mailto:hqhviet@tohoku.ac.jp" class="">hqhviet@tohoku.ac.jp</a>> wrote:<br class=""></div><div class="gmail_quote"><blockquote class="gmail_quote" style="margin: 0px 0px 0px 0.8ex; border-left-width: 1px; border-left-style: solid; border-left-color: rgb(204, 204, 204); padding-left: 1ex;"><div dir="ltr" class="">Dear Patrick and Matthew,<div class=""><br class=""><div class="">Thank you very much for your answers.</div><div class="">My test code is just the code that contains the inner product calculation and the calculation of the matrix-vector multiplication. It is easy to measure the time of each calculation. I'm going to revise the code to measure the time as it's done in Patrick's code. I will also switch to Mpich so as not to worry about the MPI environment settings.</div></div></div></blockquote></div></div></div></blockquote><blockquote type="cite" class=""><div class=""><div dir="ltr" style="caret-color: rgb(0, 0, 0); font-family: Helvetica; font-size: 12px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none;" class=""><div class="gmail_quote"><blockquote class="gmail_quote" style="margin: 0px 0px 0px 0.8ex; border-left-width: 1px; border-left-style: solid; border-left-color: rgb(204, 204, 204); padding-left: 1ex;"><div dir="ltr" class=""><div class=""><div class=""><br class="">With the current number of nodes (8), the reduction time for calculating the inner product is relatively small, but I think it would get bigger if I increased the number of nodes to 100 or more. In this case, the latency of the inner product calculation would increase, since the reduced values are transmitted via a complex network structure with cables and switches. Then, hopefully, the non-blocking communication can show its effectiveness.<br class=""></div></div></div></blockquote><div class=""><br class=""></div><div class="">This is very important to do _first_. It would probably only take you a day to measure the Allreduce time on your target, say the whole machine you run on.</div><div class=""><br class=""></div></div></div></div></blockquote>Seconded! A very simple performance model is to expect reductions to take C * log(P) time on P ranks, which is of course a slowly increasing function.</div><div><br class=""></div><div><br class=""><blockquote type="cite" class=""><div class=""><div dir="ltr" style="caret-color: rgb(0, 0, 0); font-family: Helvetica; font-size: 12px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; -webkit-text-stroke-width: 0px; text-decoration: none;" class=""><div class="gmail_quote"><div class="">  Thanks</div><div class=""><br class=""></div><div class="">     Matt</div><div class=""> </div><blockquote class="gmail_quote" style="margin: 0px 0px 0px 0.8ex; border-left-width: 1px; border-left-style: solid; border-left-color: rgb(204, 204, 204); padding-left: 1ex;"><div dir="ltr" class=""><div class=""><div class="">Best,</div><div class="">Viet</div></div><div class=""><br class=""></div></div><br class=""><div class="gmail_quote"><div dir="ltr" class="gmail_attr">On Wed, Jan 27, 2021 at 2:53 AM Patrick Sanan <<a href="mailto:patrick.sanan@gmail.com" target="_blank" class="">patrick.sanan@gmail.com</a>> wrote:<br class=""></div><blockquote class="gmail_quote" style="margin: 0px 0px 0px 0.8ex; border-left-width: 1px; border-left-style: solid; border-left-color: rgb(204, 204, 204); padding-left: 1ex;"><div class=""><br class=""><div class=""><br class=""><blockquote type="cite" class=""><div class="">Am 26.01.2021 um 12:01 schrieb Matthew Knepley <<a href="mailto:knepley@gmail.com" target="_blank" class="">knepley@gmail.com</a>>:</div><br class=""><div class=""><div dir="ltr" style="font-family: Helvetica; font-size: 12px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; text-decoration: none;" class=""><div dir="ltr" class="">On Mon, Jan 25, 2021 at 11:31 PM Viet H.Q.H. <<a href="mailto:hqhviet@tohoku.ac.jp" target="_blank" class="">hqhviet@tohoku.ac.jp</a>> wrote:<br class=""></div><div class="gmail_quote"><blockquote class="gmail_quote" style="margin: 0px 0px 0px 0.8ex; border-left-width: 1px; border-left-style: solid; border-left-color: rgb(204, 204, 204); padding-left: 1ex;"><div dir="ltr" class="">Dear Patrick Sanan,<div class=""><br class="">Thank you very much for your answer, especially for your code.<br class="">I was able to compile and run your code on 8 nodes with 20 processes per node. Below is the result<br class=""><br class=""><blockquote style="margin: 0px 0px 0px 40px; border: none; padding: 0px;" class=""><font color="#0000ff" class="">Testing with 160 MPI ranks<br class=""></font><font color="#0000ff" class="">reducing an array of size 32 (256 bytes)<br class=""></font><font color="#0000ff" class="">Running 5 burnin runs and 100 tests ...  Done.<br class=""></font><font color="#0000ff" class="">For 100 runs with 5 burnin runs, on 160 MPI processes, min/max times over all ranks:<br class=""></font><font color="#0000ff" class="">MPI timer resolution:         1.0000e-06 seconds<br class=""></font><font color="#0000ff" class="">MPI timer resolution/#trials: 1.0000e-08 seconds<br class=""></font><font color="#0000ff" class="">B.   Red. Only (min/max):    8.850098e-06 / 8.890629e-06 seconds<br class=""></font><font color="#0000ff" class="">N.B. Red. Only (min/max):    1.725912e-05 / 1.733065e-05 seconds<br class=""></font><font color="#0000ff" class="">Loc. Only (min/max):         2.364278e-04 / 2.374697e-04 seconds<br class=""></font><font color="#0000ff" class="">Blocking (min/max):          2.650309e-04 / 2.650595e-04 seconds<br class=""></font><font color="#0000ff" class="">Non-Blocking (min/max):      2.673984e-04 / 2.674508e-04 seconds<br class=""></font><font color="#0000ff" class="">Observe to see if the local time is enough to hide the reduction, and see if the reduction is indeed hidden</font></blockquote><br class="">It appears that the non-blocking computation with this test is no faster than the blocking computation.<br class="">I think I am missing some suitable Intel MPI environment settings.<br class="">I am now thinking about using MPICH, which does not require any environment settings for non-blocking computation.<br class="">Could you please let me know which MPI (MPICH or OpenMPI) you used in your tests?<br class=""></div></div></blockquote></div></div></div></blockquote>I used Cray MPI, which is based on MPICH. <br class=""><blockquote type="cite" class=""><div class=""><div dir="ltr" style="font-family: Helvetica; font-size: 12px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; text-decoration: none;" class=""><div class="gmail_quote"><div class=""><br class=""></div><div class="">Note that in the test, the cost of the reduction is only 3% of the total cost. Is that worth putting time into? Is the cost definitely</div><div class="">larger in your application?</div></div></div></div></blockquote>This is very important! Note that the reductions are quite fast compared to the local work. Even on a much larger number of ranks, they might still be quite fast. The local work here is just some useless computation (here tuned to be comparable to how long it took to do a reduction on a few thousand ranks at the time). As Matt says, having some estimate of how long it takes to apply your operator (matrix) and preconditioner is essential. Do you know how to estimate that for your application? <br class=""><blockquote type="cite" class=""><div class=""><div dir="ltr" style="font-family: Helvetica; font-size: 12px; font-style: normal; font-variant-caps: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px; text-decoration: none;" class=""><div class="gmail_quote"><div class=""><br class=""></div><div class="">   Thanks,</div><div class=""><br class=""></div><div class="">     <span class="Apple-converted-space"> </span>Matt</div><div class=""> </div><blockquote class="gmail_quote" style="margin: 0px 0px 0px 0.8ex; border-left-width: 1px; border-left-style: solid; border-left-color: rgb(204, 204, 204); padding-left: 1ex;"><div dir="ltr" class=""><div class="">Thanks again.<br class="">Viet<br class=""></div><div class=""><br class=""></div></div><br class=""><div class="gmail_quote"><div dir="ltr" class="gmail_attr">On Mon, Jan 25, 2021 at 7:47 PM Patrick Sanan <<a href="mailto:patrick.sanan@gmail.com" target="_blank" class="">patrick.sanan@gmail.com</a>> wrote:<br class=""></div><blockquote class="gmail_quote" style="margin: 0px 0px 0px 0.8ex; border-left-width: 1px; border-left-style: solid; border-left-color: rgb(204, 204, 204); padding-left: 1ex;"><div class="">Sorry about the delay in responding, but I'll add a couple of points here:<div class=""><br class=""></div><div class=""><br class=""></div><div class="">1) It's important to have some reason to believe that pipelining will actually help your problem. Pipelined Krylov methods work by overlapping reductions with operator and preconditioner applications. So, to see speedup, the time for a reduction needs to be comparable to the time for the operator/preconditioner application. This will only be true in some cases - typical cases are when you have a large number of ranks/nodes, a slow network, or very fast operator/preconditioner applications (assuming that these require the same time on each rank - it's an interesting case when they don't, but unless you say otherwise I'll assume this doesn't apply to your use case). <br class=""><div class=""><div class=""><br class=""></div><div class="">2) As you're discovering, simply ensuring that asynchronous progress works, at the pure MPI level, isn't as easy as it might be, as it's so dependent on the MPI implementation.</div><div class=""><br class=""></div><div class=""><br class=""></div><div class="">For both of these reasons, I suggest setting up a test that just directly uses MPI (which you can of course do from a PETSc-style code) and allows you to compare times for blocking and non-blocking reductions, overlapping some (useless) local work. You should make sure to run multiple iterations within the script, and also run the script multiple times on the cluster (bearing in mind that it's possible that the performance will be affected by other users of the system).</div><div class=""><br class=""></div><div class="">I attach an old script I found that I used to test some of these things, to give a more concrete idea of what I mean. Note that this was used early on in our own exploration of these topics so I'm only offering it to give an idea, not as a meaningful benchmark in its own right.</div><div class=""><br class=""></div><div class=""></div></div></div></div><div class=""><div class=""><div class=""><div class=""></div><div class=""><blockquote type="cite" class=""><div class="">Am 25.01.2021 um 09:17 schrieb Viet H.Q.H. <<a href="mailto:hqhviet@tohoku.ac.jp" target="_blank" class="">hqhviet@tohoku.ac.jp</a>>:</div><br class=""><div class=""><div dir="ltr" class=""><br class=""><div class=""><font class="">Dear Barry,</font></div><div class=""><font class=""><br class=""></font></div><div class=""><font class="">Thank you very much for your information.<br class=""><br class="">It seems complicated to set environment variables to allow asynchronous progress and pinning threads to cores when using Intel MPI.<br class=""><br class=""></font></div><blockquote style="margin: 0px 0px 0px 40px; border: none; padding: 0px;" class=""><div class=""><font color="#0000ff" class="">$ export I_MPI_ASYNC_PROGRESS = 1</font></div><div class=""><font color="#0000ff" class="">$ export I_MPI_ASYNC_PROGRESS_PIN = <CPU list></font></div></blockquote><div class=""><font class=""><br class=""><a href="https://techdecoded.intel.io/resources/hiding-communication-latency-using-mpi-3-non-blocking-collectives/" target="_blank" class="">https://techdecoded.intel.io/resources/hiding-communication-latency-using-mpi-3-non-blocking-collectives/</a><br class=""><br class="">I'm still not sure how to get an appropriate "CPU list" when running MPI with multiple nodes and multiple processes on one node.<br class=""></font></div><div class=""><span style="background-color: transparent;" class="">Best,</span><br class=""></div><div class=""><font class="">Viet.</font></div><div class=""><font class=""><br class=""></font></div><div class=""><br class=""></div><div class=""><br class=""></div></div><br class=""><div class="gmail_quote"><div dir="ltr" class="gmail_attr">On Sat, Jan 23, 2021 at 3:01 AM Barry Smith <<a href="mailto:bsmith@petsc.dev" target="_blank" class="">bsmith@petsc.dev</a>> wrote:<br class=""></div><blockquote class="gmail_quote" style="margin: 0px 0px 0px 0.8ex; border-left-width: 1px; border-left-style: solid; border-left-color: rgb(204, 204, 204); padding-left: 1ex;"><div class=""><div class=""><br class=""></div><a href="https://software.intel.com/content/www/us/en/develop/documentation/mpi-developer-guide-linux/top/additional-supported-features/asynchronous-progress-control.html" target="_blank" class="">https://software.intel.com/content/www/us/en/develop/documentation/mpi-developer-guide-linux/top/additional-supported-features/asynchronous-progress-control.html</a><div class=""><br class=""></div><div class="">It states "<span style="color: rgb(85, 85, 85); font-family: intel-clear, tahoma, Helvetica, helvetica, Arial, sans-serif; font-size: 16px; background-color: rgb(226, 231, 235);" class="">and a partial support for non-blocking collectives ( </span><span style="box-sizing: border-box; color: rgb(85, 85, 85); font-family: intel-clear, tahoma, Helvetica, helvetica, Arial, sans-serif; font-size: 16px;" class="">MPI_Ibcas</span><span style="color: rgb(85, 85, 85); font-family: intel-clear, tahoma, Helvetica, helvetica, Arial, sans-serif; font-size: 16px; background-color: rgb(226, 231, 235);" class=""> t, </span><span style="box-sizing: border-box; color: rgb(85, 85, 85); font-family: intel-clear, tahoma, Helvetica, helvetica, Arial, sans-serif; font-size: 16px;" class="">MPI_Ireduce</span><span style="color: rgb(85, 85, 85); font-family: intel-clear, tahoma, Helvetica, helvetica, Arial, sans-serif; font-size: 16px; background-color: rgb(226, 231, 235);" class=""> , and </span><span style="box-sizing: border-box; color: rgb(85, 85, 85); font-family: intel-clear, tahoma, Helvetica, helvetica, Arial, sans-serif; font-size: 16px;" class="">MPI_Iallreduce</span><span style="color: rgb(85, 85, 85); font-family: intel-clear, tahoma, Helvetica, helvetica, Arial, sans-serif; font-size: 16px; background-color: rgb(226, 231, 235);" class=""> )."  I do not know what partial support means but you can try setting the variables and see if that helps.</span></div><div class=""><font color="#555555" face="intel-clear, tahoma, Helvetica, helvetica, Arial, sans-serif" size="3" class=""><span style="background-color: rgb(226, 231, 235);" class=""><br class=""></span></font></div><div class=""><font color="#555555" face="intel-clear, tahoma, Helvetica, helvetica, Arial, sans-serif" size="3" class=""><span style="background-color: rgb(226, 231, 235);" class=""><br class=""></span></font><div class=""><br class=""><blockquote type="cite" class=""><div class="">On Jan 22, 2021, at 11:20 AM, Viet H.Q.H. <<a href="mailto:hqhviet@tohoku.ac.jp" target="_blank" class="">hqhviet@tohoku.ac.jp</a>> wrote:</div><br class=""><div class=""><div dir="ltr" class=""><div dir="ltr" class=""><br class=""><div class="">Dear Victor and Berry,<br class=""><br class="">Thank you so much for your answers.<br class=""><br class="">I fixed the code with the bug in the PetscCommSplitReductionBegin function as commented by Brave.<br class=""><br class="">   <font color="#0000ff" class=""><span class=""> </span> ierr = PetscCommSplitReductionBegin (PetscObjectComm ((PetscObject) u));</font><br class=""><br class="">It was also a mistake to set the vector size too small.<br class="">I just set a vector size of 100000000 and ran the code on 4 nodes with 2 processors per node. The result is as follows<br class=""><br class="">The time used for the asynchronous calculation:<span class=""> </span><font color="#ff0000" class="">0.022043</font><br class="">+ | u | = 10000.<br class="">The time used for the synchronous calculation:<span class=""> </span><font color="#ff0000" class="">0.016188</font><br class="">+ | b | = 10000.<br class=""></div><div class=""><br class=""></div><div class="">Asynchronous computation still takes a longer time.</div><div class=""><br class=""></div><div class="">I also confirmed that PETSC_HAVE_MPI_IALLREDUCE is defined in the file $PETSC_DIR/include/petscconf.h<br class=""><br class="">I built Petsc by using the following script</div><div class=""><br class=""></div><div class=""><font color="#0000ff" class="">#!/usr/bin/bash<br class="">set -e<br class="">DATE="21.01.18"<br class="">MPIIT_DIR="/work/A/intel/2018_update2/compilers_and_libraries_2018.2.199/linux/mpi/intel64"<br class="">MKL_DIR="/work/A/intel/2018_update2/compilers_and_libraries_2018.2.199/linux/mkl"<br class="">INSTL_DIR="${HOME}/local/petsc-3.14.3"<br class="">BUILD_DIR="${HOME}/tmp/petsc/build_${DATE}"<br class="">PETSC_DIR="${HOME}/tmp/petsc"<br class=""><br class="">cd ${PETSC_DIR}<br class="">./configure --force --prefix=${INSTL_DIR} --with-mpi-dir=${MPIIT_DIR}  --with-fortran-bindings=0 --with-mpiexe=${MPIIT_DIR}/bin/mpiexec --with-valgrind-dir=${HOME}/local/valgrind --with-blaslapack-dir=${MKL_DIR} --download-make --with-debugging=0 COPTFLAGS='-O3 -march=native -mtune=native' CXXOPTFLAGS='-O3 -march=native -mtune=native' FOPTFLAGS='-O3 -march=native -mtune=native'<br class=""><br class="">make PETSC_DIR=${HOME}/tmp/petsc PETSC_ARCH=arch-linux2-c-opt all<br class="">make PETSC_DIR=${HOME}/tmp/petsc PETSC_ARCH=arch-linux2-c-opt install </font><br class=""></div><div class=""><br class=""><br class="">Intel 2018 also complies with the MPI-3 standard.<br class=""><br class="">Are there specific settings for Intel MPI to obtain the performance of the MPI_IALLREDUCE function?<br class=""></div><div class=""><br class=""></div><div class="">Sincerely,</div><div class="">Viet.</div><div class=""><br class=""></div></div><br class=""><div class="gmail_quote"><div dir="ltr" class="gmail_attr">On Fri, Jan 22, 2021 at 11:20 AM Barry Smith <<a href="mailto:bsmith@petsc.dev" target="_blank" class="">bsmith@petsc.dev</a>> wrote:<br class=""></div><blockquote class="gmail_quote" style="margin: 0px 0px 0px 0.8ex; border-left-width: 1px; border-left-style: solid; border-left-color: rgb(204, 204, 204); padding-left: 1ex;"><div class=""><div class=""><br class=""></div><span style="color: rgb(0, 0, 255);" class=""> </span><span class=""><span class=""> </span>ierr = VecNormBegin(u,NORM_2,&norm1);</span><br class=""><span class="">   <span class=""> </span>ierr = PetscCommSplitReductionBegin(PetscObjectComm((PetscObject)Ax)); </span><div class=""><span class=""><br class=""></span></div><div class=""><span class="">How come you call this on Ax and not on u? For clarity, if nothing else, I think you should call it on u.</span></div><div class=""><span class=""><br class=""></span></div><div class=""><span class="">comb.c has </span></div><div class=""><span class=""><br class=""></span></div><div class=""><div class="">/*</div><div class="">     <span class=""> </span>Split phase global vector reductions with support for combining the</div><div class="">   communication portion of several operations. Using MPI-1.1 support only</div><div class=""><br class=""></div><div class="">     <span class=""> </span>The idea for this and much of the initial code is contributed by</div><div class="">   Victor Eijkhout.</div><div class=""><br class=""></div><div class="">       Usage:</div><div class="">             VecDotBegin(Vec,Vec,PetscScalar *);</div><div class="">             VecNormBegin(Vec,NormType,PetscReal *);</div><div class="">             ....</div><div class="">             VecDotEnd(Vec,Vec,PetscScalar *);</div><div class="">             VecNormEnd(Vec,NormType,PetscReal *);</div><div class=""><br class=""></div><div class="">       Limitations:</div><div class="">         - The order of the xxxEnd() functions MUST be in the same order</div><div class="">           as the xxxBegin(). There is extensive error checking to try to</div><div class="">           insure that the user calls the routines in the correct order</div><div class="">*/</div><div class=""><br class=""></div><div class="">#include <petsc/private/vecimpl.h>    /*I   "petscvec.h"    I*/</div><div class=""><br class=""></div><div class="">static PetscErrorCode MPIPetsc_Iallreduce(void *sendbuf,void *recvbuf,PetscMPIInt count,MPI_Datatype datatype,MPI_Op op,MPI_Comm comm,MPI_Request *request)</div><div class="">{</div><div class=""> <span class=""> </span>PETSC_UNUSED PetscErrorCode ierr;</div><div class=""><br class=""></div><div class=""> <span class=""> </span>PetscFunctionBegin;</div><div class="">#if defined(PETSC_HAVE_MPI_IALLREDUCE)</div><div class=""> <span class=""> </span>ierr = MPI_Iallreduce(sendbuf,recvbuf,count,datatype,op,comm,request);CHKERRMPI(ierr);</div><div class="">#elif defined(PETSC_HAVE_MPIX_IALLREDUCE)</div><div class=""> <span class=""> </span>ierr = MPIX_Iallreduce(sendbuf,recvbuf,count,datatype,op,comm,request);CHKERRQ(ierr);</div><div class="">#else</div><div class=""> <span class=""> </span>ierr = MPIU_Allreduce(sendbuf,recvbuf,count,datatype,op,comm);CHKERRQ(ierr);</div><div class=""> <span class=""> </span>*request = MPI_REQUEST_NULL;</div><div class="">#endif</div><div class=""> <span class=""> </span>PetscFunctionReturn(0);</div><div class="">}</div><div style="color: rgb(0, 0, 255);" class=""><br class=""></div></div><div class=""><font color="#0000ff" class=""><span class=""><br class=""></span></font><div class="">So first check if $PETSC_DIR/include/petscconf.h has </div><div class=""><br class=""></div><div class=""><span class="">PETSC_HAVE_MPI_IALLREDUCE</span></div><div class=""><span class=""><br class=""></span></div><div class=""><span class="">if it does not then the standard MPI reduce is called. </span></div><div class=""><span class=""><br class=""></span></div><div class=""><span class="">If this is set then any improvement depends on the implementation of iallreduce inside the MPI you are using. </span></div><div class=""><span class=""><br class=""></span></div><div class=""><span class="">Barry</span></div><div class=""><span class=""><br class=""></span></div><div class=""><font color="#0000ff" class=""><span class=""><br class=""></span></font><blockquote type="cite" class=""><div class="">On Jan 21, 2021, at 6:52 AM, Viet H.Q.H. <<a href="mailto:hqhviet@tohoku.ac.jp" target="_blank" class="">hqhviet@tohoku.ac.jp</a>> wrote:</div><br class=""><div class=""><div dir="ltr" class=""><div class=""><br class=""></div><div class="">Hello Petsc developers and supporters,<br class=""><br class="">I would like to confirm the performance of asynchronous computations of inner product computation overlapping with matrix-vector multiplication computation by the below code.<br class=""></div><div class=""><br class=""></div><div class=""><br class=""></div><font color="#0000ff" class=""> PetscLogDouble tt1,tt2;<br class="">   <span class=""> </span>KSP ksp;<br class="">   <span class=""> </span>//ierr = VecSet(c,one);<br class="">   <span class=""> </span>ierr = VecSet(c,one);<br class="">   <span class=""> </span>ierr = VecSet(u,one);<br class="">   <span class=""> </span>ierr = VecSet(b,one);<br class=""><br class="">   <span class=""> </span>ierr = KSPCreate(PETSC_COMM_WORLD,&ksp); CHKERRQ(ierr);<br class="">   <span class=""> </span>ierr = KSP_MatMult(ksp,A,x,Ax); CHKERRQ(ierr);<br class=""><br class=""><br class="">   <span class=""> </span>ierr = PetscTime(&tt1);CHKERRQ(ierr);<br class="">   <span class=""> </span>ierr = VecNormBegin(u,NORM_2,&norm1);<br class="">   <span class=""> </span>ierr = PetscCommSplitReductionBegin(PetscObjectComm((PetscObject)Ax));<span class=""> </span><br class="">   <span class=""> </span>ierr = KSP_MatMult(ksp,A,c,Ac);<span class=""> </span><br class="">   <span class=""> </span>ierr = VecNormEnd(u,NORM_2,&norm1);<br class="">   <span class=""> </span>ierr = PetscTime(&tt2);CHKERRQ(ierr);<br class=""><br class="">   <span class=""> </span>ierr = PetscPrintf(PETSC_COMM_WORLD, "The time used for the asynchronous calculation: %f\n",tt2-tt1); CHKERRQ(ierr);<br class="">   <span class=""> </span>ierr = PetscPrintf(PETSC_COMM_WORLD,"+ |u| =  %g\n",(double) norm1); CHKERRQ(ierr);<br class=""><br class=""><br class="">   <span class=""> </span>ierr = PetscTime(&tt1);CHKERRQ(ierr);<br class="">   <span class=""> </span>ierr = VecNorm(b,NORM_2,&norm2); CHKERRQ(ierr);<br class="">   <span class=""> </span>ierr = KSP_MatMult(ksp,A,c,Ac);<span class=""> </span><br class="">   <span class=""> </span>ierr = PetscTime(&tt2);CHKERRQ(ierr);<br class=""><br class="">   <span class=""> </span>ierr = PetscPrintf(PETSC_COMM_WORLD, "The time used for the synchronous calculation: %f\n",tt2-tt1); CHKERRQ(ierr);<br class="">   <span class=""> </span>ierr = PetscPrintf(PETSC_COMM_WORLD,"+ |b| =  %g\n",(double) norm2); CHKERRQ(ierr);<br class=""></font><div class=""><font color="#0000ff" class=""><br class=""></font></div><div class=""><br class=""></div><div class="">On a cluster with two or four nodes, the asynchronous computation is always much slower than synchronous computation.<br class=""><br class=""><font color="#ff0000" class="">The time used for the asynchronous calculation: 0.000203<br class="">+ |u| =  100.<br class="">The time used for the synchronous calculation: 0.000006<br class="">+ |b| =  100.</font><br class=""><br class="">Are there any necessary settings on MPI or Petsc to gain performance of asynchronous computation?<br class=""><br class="">Thank you very much for anything you can provide.<br class="">Sincerely,<br class="">Viet.<br class=""></div><div class=""><br class=""></div><div class=""><br class=""></div><div class=""><br class=""></div><div class=""><br class=""></div><div class=""><br class=""></div><div class=""><br class=""></div></div></div></blockquote></div><br class=""></div></div></blockquote></div></div></div></blockquote></div><br class=""></div></div></blockquote></div></div></blockquote></div><br class=""></div></div></div></blockquote></div></blockquote></div><br clear="all" class=""><div class=""><br class=""></div>--<span class=""> </span><br class=""><div dir="ltr" class=""><div dir="ltr" class=""><div class=""><div dir="ltr" class=""><div class=""><div dir="ltr" class=""><div class="">What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.<br class="">-- Norbert Wiener</div><div class=""><br class=""></div><div class=""><a href="http://www.cse.buffalo.edu/~knepley/" target="_blank" class="">https://www.cse.buffalo.edu/~knepley/</a></div></div></div></div></div></div></div></div></div></blockquote></div><br class=""></div></blockquote></div></blockquote></div><br clear="all" class=""><div class=""><br class=""></div>--<span class="Apple-converted-space"> </span><br class=""><div dir="ltr" class="gmail_signature"><div dir="ltr" class=""><div class=""><div dir="ltr" class=""><div class=""><div dir="ltr" class=""><div class="">What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.<br class="">-- Norbert Wiener</div><div class=""><br class=""></div><div class=""><a href="http://www.cse.buffalo.edu/~knepley/" target="_blank" class="">https://www.cse.buffalo.edu/~knepley/</a></div></div></div></div></div></div></div></div></div></blockquote></div><br class=""></body></html>