Skip to content
Snippets Groups Projects
Commit cffd1657 authored by kar's avatar kar Committed by Christopher Subich
Browse files

initial work to refactoring QSP

parent 6e22524e
No related branches found
No related tags found
No related merge requests found
...@@ -58,175 +58,261 @@ public: ...@@ -58,175 +58,261 @@ public:
~QSPVector() { free(data); } ~QSPVector() { free(data); }
}; };
int index(int row, int col, int rows, int cols) { enum QSPType {
assert(row >= 0 && row < rows); QSP_u,
assert(col >= 0 && col < cols); QSP_v,
return (row * cols) + col; QSP_w,
} QSP_ke,
QSP_temp,
QSP_rho,
QSP_salinity,
};
void QSPCount(const TArrayn::DTArray &t, const TArrayn::DTArray &u, struct QSPOptions {
const TArrayn::DTArray &v, const TArrayn::DTArray &w, int NS;
const char T1_name, const char S1_name, const int NS, int NT;
const int NT, double T1_max, double S1_max, double T1_min, string filename;
double S1_min, const int Nx, const int Ny, const int Nz, double S1_max;
string filename, const int plotnum, bool mapped, double S1_min;
TArrayn::DTArray *xgrid, TArrayn::DTArray *ygrid, double T1_max;
TArrayn::DTArray *zgrid) { double T1_min;
string T1_name;
string S1_name;
};
int local_rank; struct QSPData {
MPI_Comm_rank(MPI_COMM_WORLD, &local_rank); TArrayn::DTArray *u;
const TArrayn::DTArray *T1_ptr = NULL, *S1_ptr = NULL; TArrayn::DTArray *v;
TArrayn::DTArray *w;
TArrayn::DTArray *temp;
TArrayn::DTArray *rho;
TArrayn::DTArray *salinity;
TArrayn::DTArray *xgrid;
TArrayn::DTArray *ygrid;
TArrayn::DTArray *zgrid;
int Nx;
int Ny;
int Nz;
int plotnum;
bool mapped;
};
switch (T1_name) { void QSP_write(int local_rank, const QSPVector &local_hist,
case 'u': const QSPOptions &qsp_options, const QSPData &qsp_data) {
T1_ptr = &u; if (local_rank == 0) {
break; QSPVector glob_hist(qsp_options.NS, qsp_options.NT);
case 'v': MPI_Reduce(local_hist.raw(), glob_hist.raw(), // send and receive buffers
T1_ptr = &v; qsp_options.NS * qsp_options.NT, // Count
break; MPI_DOUBLE, // datatype
case 'w': MPI_SUM, 0, // Reduction operator and root process #
T1_ptr = &w; MPI_COMM_WORLD); // Communicator
break; string filename = qsp_options.filename + "." +
case 'k': boost::lexical_cast<string>(qsp_data.plotnum) + ".csv";
break; std::fstream outfile;
case 't': outfile.open(filename.c_str(), std::ios_base::out);
T1_ptr = &t; if (outfile.is_open()) {
break; outfile << T1_max << ',' << T1_min << ',' << S1_max << ',' << S1_min;
case 'T': for (int i = 4; i < NT; i++) {
T1_ptr = &t; outfile << ',' << 0;
break; }
default: outfile << std::endl;
return; for (int ii = 0; ii < NS; ii++) {
outfile << glob_hist(ii, 0);
for (int jj = 1; jj < NT; jj++) {
outfile << ',' << glob_hist(ii, jj);
}
outfile << std::endl;
}
}
} else {
MPI_Reduce(local_hist.raw(), NULL, // send and receive buffers
NS * NT, MPI_DOUBLE, // count and datatype
MPI_SUM, 0, // Reduction operator and root process
MPI_COMM_WORLD); // Communicator
}
}
QSPType QSPConvert(const std::string &name) {
QSPType converted_type;
if (name.compare("u") == 0) {
converted_type = QSP_u;
} else if (name.compare("v") == 0) {
converted_type = QSP_v;
} else if (name.compare("w") == 0) {
converted_type = QSP_w;
} else if (name.compare("ke") == 0) {
converted_type = QSP_ke;
} else if (name.compare("temp") == 0) {
converted_type = QSP_temp;
} else if (name.compare("rho") == 0) {
converted_type = QSP_rho;
} else if (name.compare("salinity") == 0) {
converted_type = QSP_salinity;
} }
return converted_type;
}
TArrayn::DTArray *QSPPtr(const QSPData &qsp_data, const QSPType &type) {
TArrayn::DTArray *ptr = NULL;
switch (S1_name) { switch (type) {
case 'u': case QSP_u:
S1_ptr = &u; ptr = qsp_data.u;
break; break;
case 'v': case QSP_v:
S1_ptr = &v; ptr = qsp_data.v;
break; break;
case 'w': case QSP_w:
S1_ptr = &w; ptr = qsp_data.w;
break; break;
case 'k': case QSP_ke: // This is an odd case.
break; break;
case 't': case QSP_rho:
S1_ptr = &t; ptr = qsp_data.rho;
break; break;
case 'T': case QSP_temp:
S1_ptr = &t; ptr = qsp_data.temp;
break;
case QSP_salinity:
ptr = qsp_data.salinity;
break; break;
default:
return;
} }
int i_low = u.lbound(blitz::firstDim); return ptr;
int j_low = u.lbound(blitz::secondDim); }
int k_low = u.lbound(blitz::thirdDim);
int i_high = u.ubound(blitz::firstDim);
int j_high = u.ubound(blitz::secondDim);
int k_high = u.ubound(blitz::thirdDim);
// This block calculates the global min/max values, in case the user void QSPMaxMin(const QSPType &T1_type, const QSPType &S1_type,
// didn't want to specify them. TArrayn::DTArray *T1_ptr, TArrayn::DTArray *S1_ptr,
QSPOptions &qsp_options, const QSPData &qsp_data) {
double double_max = std::numeric_limits<double>::max(); double double_max = std::numeric_limits<double>::max();
double double_min = std::numeric_limits<double>::min();
if (T1_max == double_max || S1_max == double_max || S1_max == double_min ||
S1_min == double_min) {
// If One of the variables is K.E or rho, we need to hand-roll the max/min
// since, in general, the index of max(u) is the same as the index of max(v)
if (T1_name == 'k' || T1_name == 't' || S1_name == 'k' || S1_name == 't') {
double ke_max = -double_max, ke_min = double_max;
double rho_max = -double_max, rho_min = double_max;
// Main hand-rolled loop // If One of the variables is K.E or rho, we need to hand-roll the max / min
for (int i = i_low; i <= i_high; i++) { // since, in general, the index of max(u) is the same as the index of max(v)
for (int j = j_low; j <= j_high; j++) { if (T1_type == QSP_ke || T1_type == QSP_rho || S1_type == QSP_ke ||
for (int k = k_low; k <= k_high; k++) { S1_type == QSP_rho) {
double ke_current = 0, tmp = 0; double ke_max = -double_max, ke_min = double_max;
if (Nx > 1) { double rho_max = -double_max, rho_min = double_max;
tmp = u(i, j, k); // Main hand-rolled loop
ke_current += tmp * tmp; for (int i = i_low; i <= i_high; i++) {
} for (int j = j_low; j <= j_high; j++) {
for (int k = k_low; k <= k_high; k++) {
double tmp;
if (T1_type == QSP_ke || S1_type == QSP_ke) {
double ke_current = 0;
tmp = (*qsp_data.u)(i, j, k);
ke_current += tmp * tmp;
if (Ny > 1) { if (Ny > 1) {
tmp = v(i, j, k); tmp = (*qsp_data.v)(i, j, k);
ke_current += tmp * tmp;
}
if (Nz > 1) {
tmp = w(i, j, k);
ke_current += tmp * tmp; ke_current += tmp * tmp;
} }
tmp = (*qsp_data.w)(i, j, k);
ke_current += tmp * tmp;
ke_current = 0.5 * ke_current; ke_current = 0.5 * ke_current;
double rho_current = eqn_of_state_t(t(i, j, k));
ke_max = ke_current > ke_max ? ke_current : ke_max; ke_max = ke_current > ke_max ? ke_current : ke_max;
ke_min = ke_current < ke_min ? ke_current : ke_min; ke_min = ke_current < ke_min ? ke_current : ke_min;
}
if (T1_type == QSP_rho || S1_type == QSP_rho) {
double rho_current = eqn_of_state_t(t(i, j, k));
rho_max = rho_current > rho_max ? rho_current : rho_max; rho_max = rho_current > rho_max ? rho_current : rho_max;
rho_min = rho_current < rho_min ? rho_current : rho_min; rho_min = rho_current < rho_min ? rho_current : rho_min;
} }
} }
} }
}
double glob_ke_max, glob_ke_min, glob_rho_max, glob_rho_min;
MPI_Allreduce(&ke_max, &glob_ke_max, 1, MPI_DOUBLE, MPI_MAX,
MPI_COMM_WORLD);
MPI_Allreduce(&ke_min, &glob_ke_min, 1, MPI_DOUBLE, MPI_MIN,
MPI_COMM_WORLD);
MPI_Allreduce(&rho_max, &glob_rho_max, 1, MPI_DOUBLE, MPI_MAX,
MPI_COMM_WORLD);
MPI_Allreduce(&rho_min, &glob_rho_min, 1, MPI_DOUBLE, MPI_MIN,
MPI_COMM_WORLD);
switch (T1_type) {
case QSP_ke:
qsp_options.T1_max = glob_ke_max;
qsp_options.T1_min = glob_ke_min;
break;
case QSP_rho:
qsp_options.T1_max = glob_rho_max;
qsp_options.T1_min = glob_rho_min;
break;
default:
qsp_options.T1_max = psmax(max(*T1_ptr));
qsp_options.T1_min = psmin(min(*T1_ptr));
break;
}
switch (S1_type) {
case QSP_ke:
qsp_options.S1_max = glob_ke_max;
qsp_options.S1_min = glob_ke_min;
break;
case QSP_rho:
qsp_options.S1_max = glob_rho_max;
qsp_options.S1_min = glob_rho_min;
break;
default:
qsp_options.S1_max = psmax(max(*S1_ptr));
qsp_options.S1_min = psmin(min(*S1_ptr));
break;
}
} else { // !(cond1 || cond2) == !cond1 && !cond2
qsp_options.S1_max = psmax(max(*S1_ptr));
qsp_options.S1_min = psmin(min(*S1_ptr));
qsp_options.T1_max = psmax(max(*T1_ptr));
qsp_options.T1_min = psmin(min(*T1_ptr));
}
}
double glob_ke_max, glob_ke_min, glob_rho_max, glob_rho_min; void QSPCount(QSPOptions qsp_options, QSPData qsp_data) {
MPI_Allreduce(&ke_max, &glob_ke_max, 1, MPI_DOUBLE, MPI_MAX,
MPI_COMM_WORLD);
MPI_Allreduce(&ke_min, &glob_ke_min, 1, MPI_DOUBLE, MPI_MIN,
MPI_COMM_WORLD);
MPI_Allreduce(&rho_max, &glob_rho_max, 1, MPI_DOUBLE, MPI_MAX,
MPI_COMM_WORLD);
MPI_Allreduce(&rho_min, &glob_rho_min, 1, MPI_DOUBLE, MPI_MIN,
MPI_COMM_WORLD);
switch (T1_name) { int local_rank;
case 'k': MPI_Comm_rank(MPI_COMM_WORLD, &local_rank);
T1_max = glob_ke_max;
T1_min = glob_ke_min;
break;
case 't':
T1_max = glob_rho_max;
T1_min = glob_rho_min;
break;
default:
T1_max = psmax(max(*T1_ptr));
T1_min = psmin(min(*T1_ptr));
break;
}
switch (S1_name) { // Find out what
case 'k': QSPType S1_type = QSPConvert(qsp_options.S1_name);
S1_max = glob_ke_max; QSPType T1_type = QSPConvert(qsp_options.T1_name);
S1_min = glob_ke_min; TArrayn::DTArray *S1_ptr = QSPPtr(qsp_data, S1_type);
break; TArrayn::DTArray *T1_ptr = QSPPtr(qsp_data, T1_type);
case 't':
S1_max = glob_rho_max; if ((!S1_ptr && S1_type != QSP_ke) || (!T1_ptr && T1_type != QSP_ke)) {
S1_min = glob_rho_min; std::cout << "Not enough data was provided for the requested tracer. "
break; "Aborting...\n";
default: return;
S1_max = psmax(max(*S1_ptr));
S1_min = psmin(min(*S1_ptr));
break;
}
} else { // !(cond1 || cond2) == !cond1 && !cond2
S1_max = psmax(max(*S1_ptr));
S1_min = psmin(min(*S1_ptr));
T1_max = psmax(max(*T1_ptr));
T1_min = psmin(min(*T1_ptr));
}
} }
double hS = (S1_max - S1_min) / (double)NS; int i_low, j_low, k_low, i_high, j_high, k_high;
double hT = (T1_max - T1_min) / (double)NT; TArrayn::DTArray *temp_ptr;
if (S1_ptr) { // If S1 is not ke
temp_ptr = S1_ptr;
} else { // If S1 is ke we know u must exist
temp_ptr = qsp_data.u;
}
int i_low = S1_ptr->lbound(blitz::firstDim);
int j_low = S1_ptr->lbound(blitz::secondDim);
int k_low = S1_ptr->lbound(blitz::thirdDim);
int i_high = S1_ptr->ubound(blitz::firstDim);
int j_high = S1_ptr->ubound(blitz::secondDim);
int k_high = S1_ptr->ubound(blitz::thirdDim);
double double_max = std::numeric_limits<double>::max();
if (qsp_options.T1_max == double_max || qsp_options.S1_max == double_max ||
qsp_options.T1_min == -double_max || qsp_options.S1_min == -double_max) {
QSPMaxMin(T1_type, S1_type, T1_ptr, S1_ptr, qsp_options, qsp_data);
}
double hS = (qsp_options.S1_max - qsp_options.S1_min) / (double)NS;
double hT = (qsp_options.T1_max - qsp_options.T1_min) / (double)NT;
double hS_inv = 1 / hS; double hS_inv = 1 / hS;
double hT_inv = 1 / hT; double hT_inv = 1 / hT;
QSPVector local_hist(NS, NT); QSPVector local_hist(qsp_options.NS, qsp_options.NT);
QSPVector global_z_max(Nx, Ny); QSPVector global_z_max(qsp_data.Nx, qsp_data.Ny);
QSPVector global_z_min(Nx, Ny); QSPVector global_z_min(qsp_data.Nx, qsp_data.Ny);
// Find the range of Lz values per 2D-slice // Find the range of Lz values per 2D-slice
if (mapped) { if (qsp_data.mapped) {
QSPVector local_z_max(Nx, Ny); QSPVector local_z_max(qsp_data.Nx, qsp_data.Ny);
QSPVector local_z_min(Nx, Ny); QSPVector local_z_min(qsp_data.Nx, qsp_data.Ny);
// We are slicing as if we are doing zgrid[i, j, :] // We are slicing as if we are doing zgrid[i, j, :]
for (int ii = i_low; ii <= i_high; ii++) { for (int ii = i_low; ii <= i_high; ii++) {
for (int jj = j_low; jj <= j_high; jj++) { for (int jj = j_low; jj <= j_high; jj++) {
...@@ -236,7 +322,7 @@ void QSPCount(const TArrayn::DTArray &t, const TArrayn::DTArray &u, ...@@ -236,7 +322,7 @@ void QSPCount(const TArrayn::DTArray &t, const TArrayn::DTArray &u,
double tmp_z_max = -tmp_z_min; double tmp_z_max = -tmp_z_min;
double tmp; double tmp;
for (int kk = k_low; kk <= k_high; kk++) { for (int kk = k_low; kk <= k_high; kk++) {
tmp = (*zgrid)(ii, jj, kk); tmp = (*qsp_data.zgrid)(ii, jj, kk);
if (tmp > tmp_z_max) { if (tmp > tmp_z_max) {
tmp_z_max = tmp; tmp_z_max = tmp;
} else if (tmp < tmp_z_min) { } else if (tmp < tmp_z_min) {
...@@ -259,64 +345,57 @@ void QSPCount(const TArrayn::DTArray &t, const TArrayn::DTArray &u, ...@@ -259,64 +345,57 @@ void QSPCount(const TArrayn::DTArray &t, const TArrayn::DTArray &u,
for (int j = j_low; j <= j_high; j++) { for (int j = j_low; j <= j_high; j++) {
for (int k = k_low; k <= k_high; k++) { for (int k = k_low; k <= k_high; k++) {
if (T1_name == 'k') { switch (T1_type) {
case QSP_ke:
Tval = 0; Tval = 0;
if (Nx > 1) { tmp = (*qsp_data.u)(i, j, k);
tmp = u(i, j, k); Tval += tmp * tmp;
Tval += tmp * tmp; tmp = (*qsp_data.w)(i, j, k);
} Tval += tmp * tmp;
if (Ny > 1) { if (Ny > 1) {
tmp = v(i, j, k); tmp = (*qsp_data.v)(i, j, k);
Tval += tmp * tmp;
}
if (Nz > 1) {
tmp = w(i, j, k);
Tval += tmp * tmp; Tval += tmp * tmp;
} }
Tval = 0.5 * Tval; Tval = 0.5 * Tval;
} else if (T1_name == 't') { break;
case QSP_rho:
tmp = (*T1_ptr)(i, j, k); tmp = (*T1_ptr)(i, j, k);
Tval = eqn_of_state_t(tmp); Tval = eqn_of_state_t(tmp);
} else { break;
default:
Tval = (*T1_ptr)(i, j, k); Tval = (*T1_ptr)(i, j, k);
} break;
int idxT = floor((Tval - T1_min) * hT_inv);
if (idxT < 0) {
idxT = 0;
} else if (idxT >= NT) {
idxT = 0;
} }
if (S1_name == 'k') { switch (S1_type) {
case QSP_ke:
Sval = 0; Sval = 0;
if (Nx > 1) { tmp = (*qsp_data.u)(i, j, k);
tmp = u(i, j, k); Sval += tmp * tmp;
Sval += tmp * tmp; tmp = (*qsp_data.w)(i, j, k);
} Sval += tmp * tmp;
if (Ny > 1) { if (Ny > 1) {
tmp = v(i, j, k); tmp = (*qsp_data.v)(i, j, k);
Sval += tmp * tmp;
}
if (Nz > 1) {
tmp = w(i, j, k);
Sval += tmp * tmp; Sval += tmp * tmp;
} }
Sval = 0.5 * Sval; Sval = 0.5 * Sval;
} else if (S1_name == 't') { break;
case QSP_rho:
tmp = (*S1_ptr)(i, j, k); tmp = (*S1_ptr)(i, j, k);
Sval = eqn_of_state_t(tmp); Sval = eqn_of_state_t(tmp);
} else { break;
default:
Sval = (*S1_ptr)(i, j, k); Sval = (*S1_ptr)(i, j, k);
break;
} }
int idxS = floor((Sval - S1_min) * hS_inv); int idxS = floor((Sval - S1_min) * hS_inv);
if (idxS < 0) { int idxT = floor((Tval - T1_min) * hT_inv);
idxS = 0; idxS = std::max(std::min(idxS, qsp_options.NS), 0);
} else if (idxS >= NS) { idxT = std::max(std::min(idxT, qsp_options.NT), 0);
idxS = 0;
}
double volume_weight; double volume_weight;
if (mapped) { if (qsp_data.mapped) {
// Calculate the Lz range // Calculate the Lz range
double Lzmax_now = global_z_max(i, j); double Lzmax_now = global_z_max(i, j);
double Lzmin_now = global_z_min(i, j); double Lzmin_now = global_z_min(i, j);
...@@ -346,40 +425,11 @@ void QSPCount(const TArrayn::DTArray &t, const TArrayn::DTArray &u, ...@@ -346,40 +425,11 @@ void QSPCount(const TArrayn::DTArray &t, const TArrayn::DTArray &u,
volume_weight = 1.0; volume_weight = 1.0;
} }
// local_hist[index(idxS, idxT, NS, NT)] += volume_weight;
local_hist(idxS, idxT) += volume_weight; local_hist(idxS, idxT) += volume_weight;
} }
} }
} }
MPI_Barrier(MPI_COMM_WORLD); // Wait for everyone to finish MPI_Barrier(MPI_COMM_WORLD); // Wait for everyone to finish
if (local_rank == 0) { QSP_write(local_rank, local_hist, qsp_options, qsp_data);
QSPVector glob_hist(NS, NT);
MPI_Reduce(local_hist.raw(), glob_hist.raw(), // send and receive buffers
NS * NT, MPI_DOUBLE, // count and datatype
MPI_SUM, 0, // Reduction operator and root process #
MPI_COMM_WORLD); // Communicator
filename = filename + "." + boost::lexical_cast<string>(plotnum) + ".csv";
std::fstream outfile;
outfile.open(filename.c_str(), std::ios_base::out);
if (outfile.is_open()) {
outfile << T1_max << ',' << T1_min << ',' << S1_max << ',' << S1_min;
for (int i = 4; i < NT; i++) {
outfile << ',' << 0;
}
outfile << std::endl;
for (int ii = 0; ii < NS; ii++) {
outfile << glob_hist[index(ii, 0, NS, NT)];
for (int jj = 1; jj < NT; jj++) {
outfile << ',' << glob_hist[index(ii, jj, NS, NT)];
}
outfile << std::endl;
}
}
} else {
MPI_Reduce(local_hist.raw(), NULL, // send and receive buffers
NS * NT, MPI_DOUBLE, // count and datatype
MPI_SUM, 0, // Reduction operator and root process #
MPI_COMM_WORLD); // Communicator
}
} }
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment