[mpich-discuss] MPICH does not work for quad precision

Jeff Hammond jeff.science at gmail.com
Sun Aug 12 16:11:18 CDT 2018


Is there a reason you do not want to use mpi_reduce with a user-defined
reduction?  Perhaps there is something that can be improved about that
feature set.

Jeff

On Thu, Aug 9, 2018 at 5:14 PM Jiancang Zhuang <zhuangjc at ism.ac.jp> wrote:

> Thanks a lot.
>
> In fact, I have tried MPI_REAL16. The same thing happened.  mpi_gather
> does not work for it either.
> I have decided to rewrite a new subroutine to replace mpi_reduce.
>
> Best regards,
>
>
>
> On Fri, Aug 10, 2018 at 12:38 AM Jeff Hammond <jeff.science at gmail.com>
> wrote:
>
>> The specific issue here is that MPI_LONG_DOUBLE does not map to REAL*16.
>> You should try MPI_REAL16.  I suspect it generates an error telling you
>> that it is not supported.
>>
>> MPI_LONG_DOUBLE corresponds to C "long double", which is an awful type
>> that means a bunch of different things depending on your hardware and
>> compiler.  If you are using x86 Linux, it should be the x87 80-bit type.
>> See https://en.wikipedia.org/wiki/Long_double for details.
>>
>> Jeff
>>
>> On Thu, Aug 9, 2018 at 8:33 AM, Jeff Hammond <jeff.science at gmail.com>
>> wrote:
>>
>>> There is a challenge here because there is no standard type in C
>>> corresponding to REAL*16.  Either the reduction operation needs to be
>>> written in Fortran or MPICH needs to figure out the compiler-dependent
>>> equivalent of REAL*16 that works in C.  While GCC __float128 and Intel
>>> _Quad might be equivalent, this is not a rigorous assumption.
>>>
>>> I recommend that you write your own user-defined reduction for REAL*16
>>> with the reduction operation callback in Fortran.
>>>
>>> Jeff
>>>
>>> On Thu, Aug 9, 2018 at 2:31 AM, Jiancang Zhuang <zhuangjc at ism.ac.jp>
>>> wrote:
>>>
>>>> I have found that the fortran version mpi_reduce does not work for
>>>> real*18. This can be shown by the following program. I have not test the C
>>>> version of mpi_reduce.
>>>>
>>>>
>>>> c-----------------------Fortran code begins
>>>> -----------------------------
>>>>       implicit real*8 (a-h, o-z)
>>>>       include 'mpif.h'
>>>>       real*16 h1, h
>>>>
>>>>
>>>>       call mpi_init(ierr)
>>>>       call mpi_comm_size(mpi_comm_world, nprocs, ierr)
>>>>       call mpi_comm_rank(mpi_comm_world, myrank, ierr)
>>>>
>>>>
>>>>       h1 = (myrank+4) *2.00000000000000
>>>>       write(*,'(''before reduce --'', i4,2f12.8)')myrank, h1,h
>>>>
>>>>
>>>>       call mpi_reduce(h1,h,1,mpi_long_double,mpi_sum,0,
>>>>      &    mpi_comm_world,ierr)
>>>>       write(*,'(''after reduce --'',i4,2f12.8)')myrank, h1,h
>>>>
>>>>       call mpi_bcast(h,1,mpi_long_double,0,
>>>>      &    mpi_comm_world,ierr)
>>>>
>>>>       write(*,'(''bcastvalue -- '',i4,2f12.8)')myrank, h1,h
>>>>
>>>>       call mpi_finalize(ierr)
>>>>        end
>>>>
>>>> c-----------------------Fortran code begins
>>>> -----------------------------
>>>>
>>>>
>>>> $ mpif77 a.f -o a.out
>>>> $ mpirun -np 3 ./a.out
>>>> before reduce --   1 10.00000000  0.00000000
>>>> after reduce --   1 10.00000000  0.00000000
>>>> before reduce --   2 12.00000000  0.00000000
>>>> before reduce --   0  8.00000000  0.00000000
>>>> after reduce --   2 12.00000000  0.00000000
>>>> after reduce --   0  8.00000000  8.00000000
>>>> bcastvalue --    0  8.00000000  8.00000000
>>>> bcastvalue --    1 10.00000000  8.00000000
>>>> bcastvalue --    2 12.00000000  8.00000000
>>>>
>>>>
>>>>
>>>> _______________________________________________
>>>> discuss mailing list     discuss at mpich.org
>>>> To manage subscription options or unsubscribe:
>>>> https://lists.mpich.org/mailman/listinfo/discuss
>>>>
>>>>
>>>
>>>
>>> --
>>> Jeff Hammond
>>> jeff.science at gmail.com
>>> http://jeffhammond.github.io/
>>>
>>
>>
>>
>> --
>> Jeff Hammond
>> jeff.science at gmail.com
>> http://jeffhammond.github.io/
>> _______________________________________________
>> discuss mailing list     discuss at mpich.org
>> To manage subscription options or unsubscribe:
>> https://lists.mpich.org/mailman/listinfo/discuss
>>
>
>
> Address: 10-3 Midori-cho, Tachikawa, Tokyo 190-8562, Japan
> Phone: +81-50-5533-8532
> Fax:     +81-42-526-4335
> email: zhuangjc at ism.ac.jp
> homepage:  http://bemlar.ism.ac.jp/zhuang
>
>  ----------------------------------------------------------------------------------------
> _______________________________________________
> discuss mailing list     discuss at mpich.org
> To manage subscription options or unsubscribe:
> https://lists.mpich.org/mailman/listinfo/discuss
>
-- 
Jeff Hammond
jeff.science at gmail.com
http://jeffhammond.github.io/
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mpich.org/pipermail/discuss/attachments/20180812/a215a3f8/attachment.html>
-------------- next part --------------
_______________________________________________
discuss mailing list     discuss at mpich.org
To manage subscription options or unsubscribe:
https://lists.mpich.org/mailman/listinfo/discuss


More information about the discuss mailing list