[Commits] [svn:einsteintoolkit] EOS_Omni/trunk/ (Rev. 65)

rhaas at tapir.caltech.edu rhaas at tapir.caltech.edu
Thu Aug 9 16:52:00 CDT 2012


User: rhaas
Date: 2012/08/09 04:52 PM

Modified:
 /trunk/
  configuration.ccl, param.ccl
 /trunk/src/nuc_eos/
  readtable.c
 /trunk/test/
  eos_omni.par

Log:
 add option read_table_on_single_process to have a named processor
 read the full table and MPI_Bcast it to the other processors
 
 From: Roland Haas <roland.haas at physics.gatech.edu>

File Changes:

Directory: /trunk/src/nuc_eos/
==============================

File [modified]: readtable.c
Delta lines: +54 -41
===================================================================
--- trunk/src/nuc_eos/readtable.c	2012-08-01 20:35:47 UTC (rev 64)
+++ trunk/src/nuc_eos/readtable.c	2012-08-09 21:52:00 UTC (rev 65)
@@ -8,20 +8,28 @@
 #include "cctk_Arguments.h"
 #include "cctk_Functions.h"
 
+// mini NoMPI
+#ifdef HAVE_CAPABILITY_MPI
+#include <mpi.h>
+#define BCAST(buffer, size) MPI_Bcast(buffer, size, MPI_BYTE, my_reader_process, MPI_COMM_WORLD)
+#else
+#define BCAST(buffer, size) do { /* do nothing */ } while(0)
+#endif
+
 // Interface of function for fortran module eostable
 void CCTK_FNAME(setup_eosmodule)
                 (CCTK_INT* , CCTK_INT* , CCTK_INT*, CCTK_REAL*,
                  CCTK_REAL*, CCTK_REAL*, CCTK_REAL*, CCTK_REAL*);
 
-// Catch HDF5 errors
+// call HDF5 if we are an IO processor, handle errors from HDF5
 #define HDF5_ERROR(fn_call)                                              \
-  do {                                                                   \
+  if (doIO) {                                                            \
     int _error_code = fn_call;                                           \
     if (_error_code < 0)                                                 \
       CCTK_VWarn (CCTK_WARN_ABORT, __LINE__, __FILE__, CCTK_THORNSTRING, \
                   "HDF5 call '%s' returned error code %d",               \
                   #fn_call, _error_code);                                \
-  } while (0)
+  }
 
 static int file_is_readable(const char* filename);
 static int file_is_readable(const char* filename)
@@ -48,9 +56,17 @@
   CCTK_Info(CCTK_THORNSTRING,nuceos_table_name);
   CCTK_Info(CCTK_THORNSTRING,"*******************************");
 
+  CCTK_INT my_reader_process = reader_process;
+  if (my_reader_process < 0 || my_reader_process >= CCTK_nProcs(cctkGH))
+  {
+    CCTK_VWarn(CCTK_WARN_COMPLAIN, __LINE__, __FILE__, CCTK_THORNSTRING,
+               "Requested IO process %d out of range. Reverting to process 0.", my_reader_process);
+    my_reader_process = 0;
+  }
+  const int doIO = !read_table_on_single_process || CCTK_MyProc(cctkGH) == my_reader_process;
 
-  hid_t file;
-  if (!file_is_readable(nuceos_table_name))
+  hid_t file = -1;
+  if (doIO && !file_is_readable(nuceos_table_name))
     CCTK_VWarn(CCTK_WARN_ABORT, __LINE__, __FILE__, CCTK_THORNSTRING,
                "Could not read nuceos_table_name'%s'",
                nuceos_table_name);
@@ -58,27 +74,27 @@
 
 // Use these two defines to easily read in a lot of variables in the same way
 // The first reads in one variable of a given type completely
-#define READ_EOS_HDF5(NAME,VAR,TYPE,MEM)                                      \
+#define READ_BCAST_EOS_HDF5(NAME,VAR,TYPE,NELEMS)                             \
   do {                                                                        \
     hid_t dataset;                                                            \
     HDF5_ERROR(dataset = H5Dopen(file, NAME));                                \
-    HDF5_ERROR(H5Dread(dataset, TYPE, MEM, H5S_ALL, H5P_DEFAULT, VAR));       \
+    HDF5_ERROR(H5Dread(dataset, TYPE, H5S_ALL, H5S_ALL, H5P_DEFAULT, VAR));   \
+    if (read_table_on_single_process)                                         \
+      BCAST (VAR, sizeof(*(VAR))*(NELEMS));                                   \
     HDF5_ERROR(H5Dclose(dataset));                                            \
   } while (0)
 // The second reads a given variable into a hyperslab of the alltables array
-#define READ_EOSTABLE_HDF5(NAME,OFF)                                     \
-  do {                                                                   \
-    hsize_t offset[2]     = {OFF,0};                                     \
-    H5Sselect_hyperslab(mem3, H5S_SELECT_SET, offset, NULL, var3, NULL); \
-    READ_EOS_HDF5(NAME,alltables,H5T_NATIVE_DOUBLE,mem3);                \
+#define READ_BCAST_EOSTABLE_HDF5(NAME,OFF,DIMS)                                        \
+  do {                                                                                 \
+    READ_BCAST_EOS_HDF5(NAME,&alltables[(OFF)*(DIMS)[1]],H5T_NATIVE_DOUBLE,(DIMS)[1]); \
   } while (0)
 
   // Read size of tables
   CCTK_INT nrho, ntemp, nye;
 
-  READ_EOS_HDF5("pointsrho",  &nrho,  H5T_NATIVE_INT, H5S_ALL);
-  READ_EOS_HDF5("pointstemp", &ntemp, H5T_NATIVE_INT, H5S_ALL);
-  READ_EOS_HDF5("pointsye",   &nye,   H5T_NATIVE_INT, H5S_ALL);
+  READ_BCAST_EOS_HDF5("pointsrho",  &nrho,  H5T_NATIVE_INT, 1);
+  READ_BCAST_EOS_HDF5("pointstemp", &ntemp, H5T_NATIVE_INT, 1);
+  READ_BCAST_EOS_HDF5("pointsye",   &nye,   H5T_NATIVE_INT, 1);
 
 
   // Allocate memory for tables
@@ -95,42 +111,39 @@
 
   // Prepare HDF5 to read hyperslabs into alltables
   hsize_t table_dims[2] = {19, nrho * ntemp * nye};
-  hsize_t var3[2]       = { 1, nrho * ntemp * nye};
-  hid_t mem3 =  H5Screate_simple(2, table_dims, NULL);
 
   // Read alltables
-  READ_EOSTABLE_HDF5("logpress",  0);
-  READ_EOSTABLE_HDF5("logenergy", 1);
-  READ_EOSTABLE_HDF5("entropy",   2);
-  READ_EOSTABLE_HDF5("munu",      3);
-  READ_EOSTABLE_HDF5("cs2",       4);
-  READ_EOSTABLE_HDF5("dedt",      5);
-  READ_EOSTABLE_HDF5("dpdrhoe",   6);
-  READ_EOSTABLE_HDF5("dpderho",   7);
+  READ_BCAST_EOSTABLE_HDF5("logpress",  0, table_dims);
+  READ_BCAST_EOSTABLE_HDF5("logenergy", 1, table_dims);
+  READ_BCAST_EOSTABLE_HDF5("entropy",   2, table_dims);
+  READ_BCAST_EOSTABLE_HDF5("munu",      3, table_dims);
+  READ_BCAST_EOSTABLE_HDF5("cs2",       4, table_dims);
+  READ_BCAST_EOSTABLE_HDF5("dedt",      5, table_dims);
+  READ_BCAST_EOSTABLE_HDF5("dpdrhoe",   6, table_dims);
+  READ_BCAST_EOSTABLE_HDF5("dpderho",   7, table_dims);
   // chemical potentials
-  READ_EOSTABLE_HDF5("muhat",     8);
-  READ_EOSTABLE_HDF5("mu_e",      9);
-  READ_EOSTABLE_HDF5("mu_p",     10);
-  READ_EOSTABLE_HDF5("mu_n",     11);
+  READ_BCAST_EOSTABLE_HDF5("muhat",     8, table_dims);
+  READ_BCAST_EOSTABLE_HDF5("mu_e",      9, table_dims);
+  READ_BCAST_EOSTABLE_HDF5("mu_p",     10, table_dims);
+  READ_BCAST_EOSTABLE_HDF5("mu_n",     11, table_dims);
   // compositions
-  READ_EOSTABLE_HDF5("Xa",       12);
-  READ_EOSTABLE_HDF5("Xh",       13);
-  READ_EOSTABLE_HDF5("Xn",       14);
-  READ_EOSTABLE_HDF5("Xp",       15);
+  READ_BCAST_EOSTABLE_HDF5("Xa",       12, table_dims);
+  READ_BCAST_EOSTABLE_HDF5("Xh",       13, table_dims);
+  READ_BCAST_EOSTABLE_HDF5("Xn",       14, table_dims);
+  READ_BCAST_EOSTABLE_HDF5("Xp",       15, table_dims);
   // average nucleus
-  READ_EOSTABLE_HDF5("Abar",     16);
-  READ_EOSTABLE_HDF5("Zbar",     17);
+  READ_BCAST_EOSTABLE_HDF5("Abar",     16, table_dims);
+  READ_BCAST_EOSTABLE_HDF5("Zbar",     17, table_dims);
   // Gamma
-  READ_EOSTABLE_HDF5("gamma",    18);
+  READ_BCAST_EOSTABLE_HDF5("gamma",    18, table_dims);
 
 
   // Read additional tables and variables
-  READ_EOS_HDF5("logrho",       logrho,        H5T_NATIVE_DOUBLE, H5S_ALL);
-  READ_EOS_HDF5("logtemp",      logtemp,       H5T_NATIVE_DOUBLE, H5S_ALL);
-  READ_EOS_HDF5("ye",           ye,            H5T_NATIVE_DOUBLE, H5S_ALL);
-  READ_EOS_HDF5("energy_shift", &energy_shift, H5T_NATIVE_DOUBLE, H5S_ALL);
+  READ_BCAST_EOS_HDF5("logrho",       logrho,        H5T_NATIVE_DOUBLE, nrho);
+  READ_BCAST_EOS_HDF5("logtemp",      logtemp,       H5T_NATIVE_DOUBLE, ntemp);
+  READ_BCAST_EOS_HDF5("ye",           ye,            H5T_NATIVE_DOUBLE, nye);
+  READ_BCAST_EOS_HDF5("energy_shift", &energy_shift, H5T_NATIVE_DOUBLE, 1);
 
-  HDF5_ERROR(H5Sclose(mem3));
   HDF5_ERROR(H5Fclose(file));
 
   // Give all values to fortran - which will store pointers to them, so don't

Directory: /trunk/test/
=======================

File [modified]: eos_omni.par
Delta lines: +3 -0
===================================================================
--- trunk/test/eos_omni.par	2012-08-01 20:35:47 UTC (rev 64)
+++ trunk/test/eos_omni.par	2012-08-09 21:52:00 UTC (rev 65)
@@ -5,6 +5,9 @@
 EOS_Omni::dump_nuceos_table = yes
 EOS_Omni::dump_nuceos_table_name = "TEST_TABLEEOS_rho4_temp4_ye4_version1.0_20120514.asc"
 
+EOS_Omni::read_table_on_single_process = yes
+EOS_Omni::reader_process = 1
+
 IO::out_dir			= $parfile
 
 #IOASCII::out1D_every = 1

Directory: /trunk/
==================

File [modified]: configuration.ccl
Delta lines: +3 -0
===================================================================
--- trunk/configuration.ccl	2012-08-01 20:35:47 UTC (rev 64)
+++ trunk/configuration.ccl	2012-08-09 21:52:00 UTC (rev 65)
@@ -5,3 +5,6 @@
   LANG
 }
 
+OPTIONAL MPI
+{
+}

File [modified]: param.ccl
Delta lines: +9 -0
===================================================================
--- trunk/param.ccl	2012-08-01 20:35:47 UTC (rev 64)
+++ trunk/param.ccl	2012-08-09 21:52:00 UTC (rev 65)
@@ -80,6 +80,15 @@
  .* :: "Can be anything"
 } "blah.h5"
 
+BOOLEAN read_table_on_single_process "read table on one process only and bcast data" STEERABLE=ALWAYS
+{
+} "no"
+
+INT reader_process "read table on this process and bcast data" STEERABLE=ALWAYS
+{
+  0:* :: "read on process N"
+} 0
+
 shares: IO
 
 USES STRING out_dir



More information about the Commits mailing list