On 8/17/22 03:10, chong liu wrote:

I modified Step-74 based on the error_estimation part of Step-50. I found it can work for the attached step-74-mpi, while it cannot work for the attached step-74-mpi-error. The only difference is the location of the output command as the attached figure 1 shows. This is weird. The error information states that one vector is out of bounds as shown in figure 2.

The location of the print message isn't the problem -- it's just that because you don't end the
  std::cout << "...";
call with std::endl, the text is put into a buffer but never printed to screen. You should run the program in a debugger and check where the erroneous access comes from. I bet that if the two programs are otherwise the same, they should both fail.


In addition, there are three points I would like to ask

 1. The direct solver cannot work for the modified MPI program. I changed it
    to the iterative solver (solver_cg same as Step-40) since I am not
    familiar with the MPI direct solver. Could you give me some suggestions on
    the direct solver for MPI?

There isn't a good option. There are some parallel direct solvers in both the PETSc and Trilinos (for which you can probably find information by searching the mailing list archives), but at the end of the day, if the problem becomes sufficiently large, even parallel direct solvers cannot compete.


 2. Does not ConvergenceTable support the parallel data output? I found that
    the first parameter for convergencetable.write_text() is std::out. How can
    I modify it to pcout for MPI?

I'm not sure the class supports this, but you can always put a
  if (Utilities::MPI::this_mpi_process(...) == 0)
in front of the place where you generate output.


 3. I guess the l1_norm() calculation for a vector should be modified. For
    example, the code std::sqrt(energy_norm_square_per_cell.l1_norm()) should
    be modified to std::sqrt
    
<https://nam10.safelinks.protection.outlook.com/?url=https%3A%2F%2Fwww.dealii.org%2Fcurrent%2Fdoxygen%2Fdeal.II%2Fnamespacestd.html%23aa2def4dc0bd5b6d07af7aff83feba9b2&data=05%7C01%7CWolfgang.Bangerth%40colostate.edu%7C8676b13756db4734d87408da803061ae%7Cafb58802ff7a4bb1ab21367ff2ecfc8b%7C0%7C0%7C637963243499566483%7CUnknown%7CTWFpbGZsb3d8eyJWIjoiMC4wLjAwMDAiLCJQIjoiV2luMzIiLCJBTiI6Ik1haWwiLCJXVCI6Mn0%3D%7C2000%7C%7C%7C&sdata=G7i0Xg9AQIQkzIUAodl745glNO9sRfgnok1piSbkatM%3D&reserved=0>(Utilities::MPI::sum
    
<https://nam10.safelinks.protection.outlook.com/?url=https%3A%2F%2Fwww.dealii.org%2Fcurrent%2Fdoxygen%2Fdeal.II%2FnamespaceUtilities_1_1MPI.html%23ab544a3bf3301a6dd3e705ee352c5551b&data=05%7C01%7CWolfgang.Bangerth%40colostate.edu%7C8676b13756db4734d87408da803061ae%7Cafb58802ff7a4bb1ab21367ff2ecfc8b%7C0%7C0%7C637963243499566483%7CUnknown%7CTWFpbGZsb3d8eyJWIjoiMC4wLjAwMDAiLCJQIjoiV2luMzIiLCJBTiI6Ik1haWwiLCJXVCI6Mn0%3D%7C2000%7C%7C%7C&sdata=dZzc08aCivKuCOOuxnaeEr0ZBSnMB3qRzuVUNDTShvY%3D&reserved=0>(estimated_error_square_per_cell.l1_norm(),
    mpi_communicator)).

Yes, something like this.

Best
 W.


--
------------------------------------------------------------------------
Wolfgang Bangerth          email:                 [email protected]
                           www: http://www.math.colostate.edu/~bangerth/

--
The deal.II project is located at http://www.dealii.org/
For mailing list/forum options, see 
https://groups.google.com/d/forum/dealii?hl=en
--- You received this message because you are subscribed to the Google Groups "deal.II User Group" group.
To unsubscribe from this group and stop receiving emails from it, send an email 
to [email protected].
To view this discussion on the web visit 
https://groups.google.com/d/msgid/dealii/0dc890b4-4592-23da-095d-bd5d88ae613a%40colostate.edu.

Reply via email to