cs205-lecture-examples

Example codes used during Harvard CS205 lectures
git clone https://git.0xfab.ch/cs205-lecture-examples.git
Log | Files | Refs | README | LICENSE

IO.h (4169B)


      1 #ifndef IO_H_C0WPQOQZ
      2 #define IO_H_C0WPQOQZ
      3 
      4 #include <fstream>
      5 #include <mpi.h>
      6 #include <sstream>
      7 #include <string>
      8 #include <vector>
      9 
     10 void initialize(std::vector<double> &v)
     11 {
     12     int rank, size;
     13     MPI_Comm_rank(MPI_COMM_WORLD, &rank);
     14     MPI_Comm_size(MPI_COMM_WORLD, &size);
     15     const int N = static_cast<int>(v.size()) - 2;
     16     const double dx = 1.0 / (size * N);
     17     const double x0 = rank * N * dx;
     18     for (int i = 0; i < N; ++i) {
     19         const double x = x0 + (i + 0.5) * dx;
     20         if (x > 0.4 && x < 0.6) {
     21             v[i + 1] = 1.0;
     22         } else {
     23             v[i + 1] = 0.0;
     24         }
     25     }
     26 }
     27 
     28 void dump_posix(const std::vector<double> &v)
     29 {
     30     int rank, size;
     31     MPI_Comm_rank(MPI_COMM_WORLD, &rank);
     32     MPI_Comm_size(MPI_COMM_WORLD, &size);
     33     const int N = static_cast<int>(v.size()) - 2; // ignore ghost cells
     34     std::vector<double> u(size * N);
     35     MPI_Gather(
     36         &v[1], N, MPI_DOUBLE, u.data(), N, MPI_DOUBLE, 0, MPI_COMM_WORLD);
     37     if (0 == rank) { // not parallel
     38         std::ofstream out("u.dat", std::ios::app);
     39         out << std::scientific;
     40         for (size_t i = 0; i < u.size(); ++i) {
     41             out << ' ' << u[i];
     42         }
     43         out << '\n';
     44     }
     45 }
     46 
     47 void dump_rank(const std::vector<double> &v)
     48 {
     49     int rank;
     50     MPI_Comm_rank(MPI_COMM_WORLD, &rank);
     51     std::ofstream out("u_rank" + std::to_string(rank) + ".dat", std::ios::app);
     52     out << std::scientific;
     53     for (int i = 1; i < (int)v.size() - 1; ++i) {
     54         out << ' ' << v[i]; // skip ghost cells
     55     }
     56     out << '\n';
     57 }
     58 
     59 void dump_mpi_binary(const std::vector<double> &v)
     60 {
     61     int rank;
     62     MPI_Comm_rank(MPI_COMM_WORLD, &rank);
     63     MPI_File fh; // MPI file handle
     64     MPI_Offset bytes, file_offset;
     65 
     66     const int N = static_cast<int>(v.size()) - 2; // ignore ghost cells
     67     bytes = N * sizeof(double);
     68     file_offset = rank * bytes;
     69 
     70     MPI_File_open(MPI_COMM_WORLD,
     71                   "u_mpi.bin",
     72                   MPI_MODE_CREATE | MPI_MODE_WRONLY,
     73                   MPI_INFO_NULL,
     74                   &fh);
     75     MPI_File_get_size(fh, &bytes); // append to file
     76     file_offset += bytes;
     77     MPI_File_write_at_all(
     78         fh, file_offset, &v[1], N, MPI_DOUBLE, MPI_STATUS_IGNORE);
     79     MPI_File_close(&fh);
     80 }
     81 
     82 void dump_mpi_ascii(const std::vector<double> &v)
     83 {
     84     int rank, size;
     85     MPI_Comm_rank(MPI_COMM_WORLD, &rank);
     86     MPI_Comm_size(MPI_COMM_WORLD, &size);
     87     std::ostringstream out;
     88     out << std::scientific;
     89     for (int i = 1; i < (int)v.size() - 1; ++i) {
     90         out << ' ' << v[i]; // skip ghost cells
     91     }
     92     if (rank == size - 1) {
     93         out << '\n';
     94     }
     95     const std::string ascii = out.str();
     96 
     97     MPI_File fh; // MPI file handle
     98     MPI_Offset bytes, file_offset = 0;
     99     bytes = ascii.size() * sizeof(char);
    100     MPI_Exscan(&bytes, &file_offset, 1, MPI_OFFSET, MPI_SUM, MPI_COMM_WORLD);
    101 
    102     MPI_File_open(MPI_COMM_WORLD,
    103                   "u_mpi.dat",
    104                   MPI_MODE_CREATE | MPI_MODE_WRONLY,
    105                   MPI_INFO_NULL,
    106                   &fh);
    107     MPI_File_get_size(fh, &bytes); // append to file
    108     file_offset += bytes;
    109     MPI_File_write_at_all(fh,
    110                           file_offset,
    111                           ascii.data(),
    112                           ascii.size(),
    113                           MPI_CHAR,
    114                           MPI_STATUS_IGNORE);
    115     MPI_File_close(&fh);
    116 }
    117 
    118 void dump_mpi_ascii_ordered(const std::vector<double> &v)
    119 {
    120     int rank, size;
    121     MPI_Comm_rank(MPI_COMM_WORLD, &rank);
    122     MPI_Comm_size(MPI_COMM_WORLD, &size);
    123     std::ostringstream out;
    124     out << std::scientific;
    125     for (int i = 1; i < (int)v.size() - 1; ++i) {
    126         out << ' ' << v[i]; // skip ghost cells
    127     }
    128     if (rank == size - 1) {
    129         out << '\n';
    130     }
    131     const std::string ascii = out.str();
    132 
    133     MPI_File fh; // MPI file handle
    134     MPI_File_open(MPI_COMM_WORLD,
    135                   "u_mpi_ordered.dat",
    136                   MPI_MODE_CREATE | MPI_MODE_WRONLY | MPI_MODE_APPEND,
    137                   MPI_INFO_NULL,
    138                   &fh);
    139     MPI_File_write_ordered(
    140         fh, ascii.data(), ascii.size(), MPI_CHAR, MPI_STATUS_IGNORE);
    141     MPI_File_close(&fh);
    142 }
    143 
    144 #endif /* IO_H_C0WPQOQZ */