<html><body>
<p><font size="2" face="sans-serif">Keep in mind, BG/P is a 64-bit MPI_AINT on a 32-bit platform.</font><br>
<br>
<br>
<font size="2" face="sans-serif">_______________________________________________<br>
Douglas Miller                  BlueGene Messaging Development<br>
IBM Corp., Rochester, MN USA                     Bldg 030-2 A401<br>
dougmill@us.ibm.com               Douglas Miller/Rochester/IBM</font><br>
<br>
<img width="16" height="16" src="cid:1__=09BBF339DFF4A3CC8f9e8a93df938@us.ibm.com" border="0" alt="Inactive hide details for Rob Latham ---02/20/2012 09:33:42 AM---Rob Latham <robl@mcs.anl.gov>"><font size="2" color="#424282" face="sans-serif">Rob Latham ---02/20/2012 09:33:42 AM---Rob Latham <robl@mcs.anl.gov></font><br>
<br>

<table width="100%" border="0" cellspacing="0" cellpadding="0">
<tr valign="top"><td style="background-image:url(cid:2__=09BBF339DFF4A3CC8f9e8a93df938@us.ibm.com); background-repeat: no-repeat; " width="40%">
<ul style="padding-left: 72pt"><font size="1" face="sans-serif"><b>Rob Latham <robl@mcs.anl.gov></b></font><font size="1" face="sans-serif"> </font><br>
<font size="1" face="sans-serif">Sent by: mpi-forum-bounces@lists.mpi-forum.org</font>
<p><font size="1" face="sans-serif">02/20/2012 09:26 AM</font>
<table border="1">
<tr valign="top"><td width="168" bgcolor="#FFFFFF"><div align="center"><font size="1" face="sans-serif">Please respond to<br>
Main MPI Forum mailing list <mpi-forum@lists.mpi-forum.org></font></div></td></tr>
</table>
</ul>
</td><td width="60%">
<table width="100%" border="0" cellspacing="0" cellpadding="0">
<tr valign="top"><td width="1%"><img width="58" height="1" src="cid:3__=09BBF339DFF4A3CC8f9e8a93df938@us.ibm.com" border="0" alt=""><br>
<div align="right"><font size="1" face="sans-serif">To</font></div></td><td width="100%"><img width="1" height="1" src="cid:3__=09BBF339DFF4A3CC8f9e8a93df938@us.ibm.com" border="0" alt=""><br>

<ul style="padding-left: 7pt"><font size="1" face="sans-serif">Main MPI Forum mailing list <mpi-forum@lists.mpi-forum.org>, </font></ul>
</td></tr>

<tr valign="top"><td width="1%"><img width="58" height="1" src="cid:3__=09BBF339DFF4A3CC8f9e8a93df938@us.ibm.com" border="0" alt=""><br>
<div align="right"><font size="1" face="sans-serif">cc</font></div></td><td width="100%"><img width="1" height="1" src="cid:3__=09BBF339DFF4A3CC8f9e8a93df938@us.ibm.com" border="0" alt=""><br>

<ul style="padding-left: 7pt"><font size="1" face="sans-serif">Harald Klimach <h.klimach@grs-sim.de></font></ul>
</td></tr>

<tr valign="top"><td width="1%"><img width="58" height="1" src="cid:3__=09BBF339DFF4A3CC8f9e8a93df938@us.ibm.com" border="0" alt=""><br>
<div align="right"><font size="1" face="sans-serif">Subject</font></div></td><td width="100%"><img width="1" height="1" src="cid:3__=09BBF339DFF4A3CC8f9e8a93df938@us.ibm.com" border="0" alt=""><br>

<ul style="padding-left: 7pt"><font size="1" face="sans-serif">Re: [Mpi-forum] MPI-IO for mesh based simulations beyond 2 billion elements</font></ul>
</td></tr>
</table>

<table border="0" cellspacing="0" cellpadding="0">
<tr valign="top"><td width="58"><img width="1" height="1" src="cid:3__=09BBF339DFF4A3CC8f9e8a93df938@us.ibm.com" border="0" alt=""></td><td width="336"><img width="1" height="1" src="cid:3__=09BBF339DFF4A3CC8f9e8a93df938@us.ibm.com" border="0" alt=""></td></tr>
</table>
</td></tr>
</table>
<br>
<tt><font size="2">On Sun, Feb 19, 2012 at 10:51:42AM -0600, Rajeev Thakur wrote:<br>
> Harald,<br>
>             If all you have is a one-dimensional array, and each process reads or writes a large contiguous chunk of it, you don't need to use Type_create_subarray or even use collective I/O. You can just have each process call MPI_File_read_at or write_at with the right offset. The subarray type is more useful if you have a multidimensional array and need to read a subsection of it that is not contiguously located in the file.<br>
<br>
Hi Rajeev: that's not entirely true any longer.  On BlueGene,<br>
collective I/O will carry out aggregation and alignment to block<br>
boundaries.  On Lustre, collective I/O can conjure up a workload such<br>
that all operations from a particular client go to one server.  <br>
<br>
So, anything that allows collective I/O, even for large requests, we<br>
should try to encourage.<br>
<br>
32 bit MPI_Aint on 32 bit platforms isn't perfect (I think there's<br>
some problem in the Fortran cases?)  but it's allowed us to actually<br>
get work done on BlueGene.<br>
<br>
==rob<br>
<br>
> <br>
> Rajeev<br>
> <br>
> <br>
> On Feb 19, 2012, at 4:23 AM, Rolf Rabenseifner wrote:<br>
> <br>
> > Dear Harald, dear I/O and datatype committee members,<br>
> > <br>
> > Indeed, you need a workaround. Our LARGE COUNT #265 is not enough.<br>
> > The workaround (instead of MPI_TYPE_CREATE_SUBARRAY) is:<br>
> > <br>
> > 1. Call MPI_FILE_SET_VIEW(fh, 0, vectype, vectype, datarep, MPI_INFO_NULL, iError)<br>
> > 2. Call MPI_FILE_GET_TYPE_EXTENT(fh, vectype, vectype_byte_extent, iError)<br>
> > 3. Byte_elemOffset = elemOff (as in your source code) * vectype_byte_extent<br>
> >    Byte_GlobalSize = globElems (as in your source code) * vectype_byte_extent<br>
> > 4. Call MPI_Type_contiguous(locElems, vectype, chunktype, iError)<br>
> > 5. Call MPI_TYPE_CREATE_RESIZED(chunktype,-Byte_elemOffset,<br>
> >                                 Byte_GlobalSize, ftype)<br>
> >    (I hope that I understood thisroutine correctly <br>
> >     and here, a negativ lb is needed.)<br>
> > <br>
> > 6. call MPI_File_set_view( fh, 0_MPI_OFFSET_KIND, &<br>
> >> & vectype, ftype, "native", &<br>
> >> & MPI_INFO_NULL, iError )<br>
> > <br>
> > Any comments by specialists?<br>
> > <br>
> > Best regards<br>
> > Rolf<br>
> > <br>
> > ----- Original Message -----<br>
> >> From: "Harald Klimach" <h.klimach@grs-sim.de><br>
> >> To: "Rolf Rabenseifner" <rabenseifner@hlrs.de><br>
> >> Sent: Sunday, February 19, 2012 10:05:40 AM<br>
> >> Subject: MPI-IO for mesh based simulations beyond 2 billion elements<br>
> >> <br>
> >> Dear Rolf,<br>
> >> <br>
> >> Mesh based simulations are an important backbone of CFD. The mesh<br>
> >> sizes we are capable<br>
> >> to compute on new HPC systems steadily increases even for complex<br>
> >> geometrical layouts.<br>
> >> For simple geometries, much larger domains then 2 billion elements<br>
> >> where already done<br>
> >> 10 years ago, see<br>
> >> M. Yokokawa, K. Itakura, A. Uno, T. Ishihara and Y. Kaneda,<br>
> >> "16.4-TFlops Direct Numerical Simulation of Turbulence by a Fourier<br>
> >> Spectral Method on the Earth Simulator", Proceedings of the 2002<br>
> >> ACM/IEEE Conference on Supercomputing, Baltimore MD (2002).<br>
> >> To describe the spatial domain around complex geometries usually some<br>
> >> unstructured<br>
> >> approach has to be followed to allow a sparse representation. This<br>
> >> results in a one-dimensional<br>
> >> list of all elements which is distributed across all processes in a<br>
> >> parallel simulation.<br>
> >> Therefore the natural representation of all elements in each partition<br>
> >> is perfectly<br>
> >> described by an one-dimensional MPI subarray datatype.<br>
> >> This would also be the most natural choice to describe the MPI File<br>
> >> view for each<br>
> >> partition to enable dumping and reading of results and restart data.<br>
> >> However, due to the interface definition of MPI_Type_create_subarray,<br>
> >> the global<br>
> >> problem size would be limited to roughly 2 billion elements by MPI.<br>
> >> <br>
> >> Here is the actual implementation we are facing in our Code:<br>
> >> We use a Lattice-Boltzmann representation, where each element is<br>
> >> described by<br>
> >> 19 MPI_DOUBLE_PRECISION values:<br>
> >> <br>
> >> call MPI_Type_contiguous(19, MPI_DOUBLE_PRECISION, vectype, iError)<br>
> >> call MPI_Type_commit( vectype, iError )<br>
> >> <br>
> >> This type is then used to build up the subarray description of each<br>
> >> partition:<br>
> >> <br>
> >> call MPI_Type_create_subarray( 1, [ globElems ], [ locElems ], [<br>
> >> elemOff ], &<br>
> >> & MPI_ORDER_FORTRAN, vectype, ftype, iError)<br>
> >> call MPI_Type_commit( ftype, iError )<br>
> >> <br>
> >> Where globElems is the global number of Elements, locElems the number<br>
> >> of<br>
> >> Elements in the process local partition, and elemOff the Offset of the<br>
> >> process<br>
> >> local partition in terms of elements. This type is then used to set<br>
> >> the view<br>
> >> for the local process on a file opened for MPI-IO:<br>
> >> <br>
> >> call MPI_File_open( MPI_COMM_WORLD, &<br>
> >> & 'restart.dat', &<br>
> >> & MPI_MODE_WRONLY+MPI_MODE_CREATE, &<br>
> >> & MPI_INFO_NULL, fhandle, &<br>
> >> & iError )<br>
> >> <br>
> >> call MPI_File_set_view( fhandle, 0_MPI_OFFSET_KIND, &<br>
> >> & vectype, ftype, "native", &<br>
> >> & MPI_INFO_NULL, iError )<br>
> >> <br>
> >> Finally the data is written to disk with collective write:<br>
> >> <br>
> >> call MPI_File_write_all(fhandle, data, locElems, vectype, iostatus,<br>
> >> iError)<br>
> >> <br>
> >> The global mesh is sorted according to a space-filling curve and to<br>
> >> provide<br>
> >> mostly equal distribution on all processes, the local number of<br>
> >> elements is<br>
> >> computed by the following integer arithmetic:<br>
> >> <br>
> >> locElems = globElems / nProcs<br>
> >> if (myrank < mod(globElems, nProcs) locElems = locElems + 1<br>
> >> <br>
> >> Where nProcs is the number of processes in the communicator and myRank<br>
> >> the rank of the local partition in the communicator.<br>
> >> Note, that for differently sized elements, the workload varies and the<br>
> >> elements<br>
> >> shouldn't be distributed equally for work load balance. And with<br>
> >> different<br>
> >> partitioning strategies other distributions might arise. However in<br>
> >> any case it<br>
> >> is not generally possible to find common discriminators greater than<br>
> >> one<br>
> >> across all partitions.<br>
> >> Thus it is not easily possible to find a larger etype, which spans<br>
> >> more than a<br>
> >> single element.<br>
> >> <br>
> >> What is the recommended implementation for efficient parallel IO<br>
> >> within MPI<br>
> >> in such a scenario?<br>
> >> <br>
> >> Thanks a lot,<br>
> >> Harald<br>
> >> <br>
> >> --<br>
> >> Dipl.-Ing. Harald Klimach<br>
> >> <br>
> >> German Research School for<br>
> >> Simulation Sciences GmbH<br>
> >> Schinkelstr. 2a<br>
> >> 52062 Aachen | Germany<br>
> >> <br>
> >> Tel +49 241 / 80-99758<br>
> >> Fax +49 241 / 806-99758<br>
> >> Web </font></tt><tt><font size="2">www.grs-sim.de</font></tt><tt><font size="2"><br>
> >> JID: haraldkl@jabber.ccc.de<br>
> >> <br>
> >> Members: Forschungszentrum Jülich GmbH | RWTH Aachen University<br>
> >> Registered in the commercial register of the local court of<br>
> >> Düren (Amtsgericht Düren) under registration number HRB 5268<br>
> >> Registered office: Jülich<br>
> >> Executive board: Prof. Marek Behr Ph.D. | Dr. Norbert Drewes<br>
> > <br>
> > -- <br>
> > Dr. Rolf Rabenseifner . . . . . . . . . .. email rabenseifner@hlrs.de<br>
> > High Performance Computing Center (HLRS) . phone ++49(0)711/685-65530<br>
> > University of Stuttgart . . . . . . . . .. fax ++49(0)711 / 685-65832<br>
> > Head of Dpmt Parallel Computing . . . </font></tt><tt><font size="2">www.hlrs.de/people/rabenseifner</font></tt><tt><font size="2"><br>
> > Nobelstr. 19, D-70550 Stuttgart, Germany . (Office: Allmandring 30)<br>
> > <br>
> > _______________________________________________<br>
> > mpi-forum mailing list<br>
> > mpi-forum@lists.mpi-forum.org<br>
> > </font></tt><tt><font size="2"><a href="http://lists.mpi-forum.org/mailman/listinfo.cgi/mpi-forum">http://lists.mpi-forum.org/mailman/listinfo.cgi/mpi-forum</a></font></tt><tt><font size="2"><br>
> <br>
> <br>
> _______________________________________________<br>
> mpi-forum mailing list<br>
> mpi-forum@lists.mpi-forum.org<br>
> </font></tt><tt><font size="2"><a href="http://lists.mpi-forum.org/mailman/listinfo.cgi/mpi-forum">http://lists.mpi-forum.org/mailman/listinfo.cgi/mpi-forum</a></font></tt><tt><font size="2"><br>
<br>
-- <br>
Rob Latham<br>
Mathematics and Computer Science Division<br>
Argonne National Lab, IL USA<br>
_______________________________________________<br>
mpi-forum mailing list<br>
mpi-forum@lists.mpi-forum.org<br>
</font></tt><tt><font size="2"><a href="http://lists.mpi-forum.org/mailman/listinfo.cgi/mpi-forum">http://lists.mpi-forum.org/mailman/listinfo.cgi/mpi-forum</a></font></tt><tt><font size="2"><br>
<br>
</font></tt><br>
</body></html>