forked from LRZ-BADW/DPEcho
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathDomain.cpp
176 lines (155 loc) · 9.64 KB
/
Domain.cpp
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
// Copyright(C) 2020 Fabio Baruffa, Intel Corp.
// Copyright(C) 2021 Salvatore Cielo, LRZ
//
// Licensed under the Apache License, Version 2.0 (the "License"); you may not use this file except in compliance with the
// License. You may obtain a copy of the License at http://www.apache.org/licenses/LICENSE-2.0
//
// Unless required by applicable law or agreed to in writing, software distributed under the License is distributed on an
// "AS IS" BASIS, WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied. See the License for the specific
// language governing permissions and limitations under the License.
#include "Domain.hpp"
Domain::~Domain( ){ free(bufL, qq); free(bufR, qq);
Log->setPar(false); *Log+4<<TAG<<"Domain removed and BCex buffers deallocated."; Log->fl();
}
Domain::Domain(mysycl::queue q, size_t bufMax) {
Log = Logger::getInstance(); Log->setPar(false);
int nprocs = Log->getNumProcs(), myrank = Log->getMyRank(), reorder = 0;
qq = q; // Use the queue to allocate the buffers
bufL= malloc_shared<field>(FLD_TOT*bufMax, qq); if(!bufL){Log->Error("%s Cannot allocate bufL.", TAG );}
bufR= malloc_shared<field>(FLD_TOT*bufMax, qq); if(!bufR){Log->Error("%s Cannot allocate bufR.", TAG );}
// Assign default values before reading parameters
boxMin_[0] = boxMin_[1] = boxMin_[2] =-0.5;
boxMax_[0] = boxMax_[1] = boxMax_[2] = 0.5;
bcType_[0] = bcType_[1] = bcType_[2] = BCPER;
#ifdef MPICODE
cartDims_[0] = cartDims_[1] = cartDims_[2] = 0;
#else
cartDims_[0] = cartDims_[1] = cartDims_[2] = 1;
#endif
// Now read the parameters
std::string key, value; std::ifstream inFile("echo.par"); *Log+0<<TAG<<" Reading input: ";
while (std::getline(inFile, key, ' ') && std::getline(inFile, value) ){
if(!key.compare("xMin" )){ boxMin_[0] = std::stof(value); *Log<<"\n\txMin "<<boxMin_[0] <<" "; continue;}
if(!key.compare("yMin" )){ boxMin_[1] = std::stof(value); *Log<<"\n\tyMin "<<boxMin_[1] <<" "; continue;}
if(!key.compare("zMin" )){ boxMin_[2] = std::stof(value); *Log<<"\n\tzMin "<<boxMin_[2] <<" "; continue;}
if(!key.compare("xMax" )){ boxMax_[0] = std::stof(value); *Log<<"\n\txMax "<<boxMax_[0] <<" "; continue;}
if(!key.compare("yMax" )){ boxMax_[1] = std::stof(value); *Log<<"\n\tyMax "<<boxMax_[1] <<" "; continue;}
if(!key.compare("zMax" )){ boxMax_[2] = std::stof(value); *Log<<"\n\tzMax "<<boxMax_[2] <<" "; continue;}
if(!key.compare("bcTypex")){ bcType_[0] = std::stoi(value); *Log<<"\n\tbcTypex "<<bcType_[0] <<" "; continue;}
if(!key.compare("bcTypey")){ bcType_[1] = std::stoi(value); *Log<<"\n\tbcTypey "<<bcType_[1] <<" "; continue;}
if(!key.compare("bcTypez")){ bcType_[2] = std::stoi(value); *Log<<"\n\tbcTypez "<<bcType_[2] <<" "; continue;}
#ifdef MPICODE
if(!key.compare("Rx")){ cartDims_[0] = std::stoi(value); *Log<<"\n\tcartDimx "<<cartDims_[0]<<" "; continue;}
if(!key.compare("Ry")){ cartDims_[1] = std::stoi(value); *Log<<"\n\tcartDimy "<<cartDims_[1]<<" "; continue;}
if(!key.compare("Rz")){ cartDims_[2] = std::stoi(value); *Log<<"\n\tcartDimz "<<cartDims_[2]<<" "; continue;}
#else
cartCoords_[0] = cartCoords_[1] = cartCoords_[2] = 0;
#endif
}; Log->fl();
boxSize_[0] = boxMax_[0]-boxMin_[0]; boxSize_[1] = boxMax_[1]-boxMin_[1]; boxSize_[2] = boxMax_[2]-boxMin_[2];
cartPeriodic_[0] = cartPeriodic_[1] = cartPeriodic_[2] = 1; // TODO Always periodic for now
#ifdef MPICODE
// Create cartesian domain
MPI_Dims_create(nprocs, 3, cartDims_);
MPI_Cart_create(MPI_COMM_WORLD, 3, cartDims_, cartPeriodic_, reorder, &cartComm_);
// Check that domain is consistent
unsigned const prod = cartDims_[0] * cartDims_[1] * cartDims_[2];
if ( (prod) && (prod != nprocs) ){
*Log+0<<TAG<<"Aborting. cartDims: "<<prod<<" don't match MPI ranks :"<<nprocs;
Log->fl(); Log->Error("");
}
MPI_Cart_coords(cartComm_, myrank, 3, cartCoords_);
int neighCoords[3]; // One of those MPI messy interfaces. Bad code. At least it scales! - SC
neighCoords[1] = cartCoords_[1]; neighCoords[2] = cartCoords_[2]; // X-neighbors
neighCoords[0] = cartCoords_[0]-1; MPI_Cart_rank(cartComm_, neighCoords, neighRankPrev_+0);
neighCoords[0] = cartCoords_[0]+1; MPI_Cart_rank(cartComm_, neighCoords, neighRankNext_+0);
neighCoords[0] = cartCoords_[0]; neighCoords[2] = cartCoords_[2]; // Y-neighbors
neighCoords[1] = cartCoords_[1]-1; MPI_Cart_rank(cartComm_, neighCoords, neighRankPrev_+1);
neighCoords[1] = cartCoords_[1]+1; MPI_Cart_rank(cartComm_, neighCoords, neighRankNext_+1);
neighCoords[0] = cartCoords_[0]; neighCoords[1] = cartCoords_[1]; // Z-neighbors
neighCoords[2] = cartCoords_[2]-1; MPI_Cart_rank(cartComm_, neighCoords, neighRankPrev_+2);
neighCoords[2] = cartCoords_[2]+1; MPI_Cart_rank(cartComm_, neighCoords, neighRankNext_+2);
#endif
// Initialize phisical (local!) dimensions
locSize_[0] = boxSize_[0]/cartDims_[0]; locMin_[0] = boxMin_[0]+locSize_[0]*cartCoords_[0]; locMax_[0] = locMin_[0]+locSize_[0];
locSize_[1] = boxSize_[1]/cartDims_[1]; locMin_[1] = boxMin_[1]+locSize_[1]*cartCoords_[1]; locMax_[1] = locMin_[1]+locSize_[1];
locSize_[2] = boxSize_[2]/cartDims_[2]; locMin_[2] = boxMin_[2]+locSize_[2]*cartCoords_[2]; locMax_[2] = locMin_[2]+locSize_[2];
*Log+0<<TAG<<"Domain created!"; Log->fl(); cartInfo();
}
void Domain::cartInfo() {
Log->setPar(false);
*Log+1<<TAG<<"3D Domain: " <<cartDims_[0]<<" "<<cartDims_[1]<<" "<<cartDims_[2]; Log->fl();
Log->setPar(true);
*Log+1<<TAG <<"CartCoords: (" <<cartCoords_ [0] <<" "<< cartCoords_ [1] <<" "<< cartCoords_ [2] <<") "
<<"CartPeriodic: (" <<cartPeriodic (0) <<" "<< cartPeriodic (1) <<" "<< cartPeriodic (2) <<")"; Log->fl();
#ifdef MPICODE
*Log+10<<TAG<<"CartPrev: (" <<neighRankPrev_[0] <<" "<< neighRankPrev_[1] <<" "<< neighRankPrev_[2] <<") "
<<"CartNext: (" <<neighRankNext_[0] <<" "<< neighRankNext_[1] <<" "<< neighRankNext_[2] <<") "; Log->fl();
int neighCoords_[3];
for (int i = 0; i < 3; ++i){
*Log+10<<TAG <<"Prev / This / Next ["<<i<<"]" ;
MPI_Cart_coords(cartComm_, neighRankPrev_[i], 3, neighCoords_);
*Log<<"("<< neighCoords_ [0] <<" "<< neighCoords_ [1] <<" "<< neighCoords_ [2] <<") ";
*Log<<"("<< cartCoords_ [0] <<" "<< cartCoords_ [1] <<" "<< cartCoords_ [2] <<") ";
MPI_Cart_coords(cartComm_, neighRankNext_[i], 3, neighCoords_);
*Log<<"("<< neighCoords_ [0] <<" "<< neighCoords_ [1] <<" "<< neighCoords_ [2] <<") ";
Log->fl();
}
#endif
}
void Domain::boxInfo() {
Log->setPar(false);
*Log+1<<TAG<<"Box size ("<<boxSize_[0]<<" "<<boxSize_[1]<<" "<<boxSize_[2]
<<") From (" <<boxMin_ [0]<<" "<<boxMin_ [1]<<" "<<boxMin_ [2]
<<") To (" <<boxMax_ [0]<<" "<<boxMax_ [1]<<" "<<boxMax_ [2]<<")"; Log->fl();
}
void Domain::locInfo() {
Log->setPar(true);
*Log+1<<TAG<<"Local size ("<<locSize_[0]<<" "<<locSize_[1]<<" "<<locSize_[2]
<<") From (" <<locMin_ [0]<<" "<<locMin_ [1]<<" "<<locMin_ [2]
<<") To (" <<locMax_ [0]<<" "<<locMax_ [1]<<" "<<locMax_ [2]<<")"; Log->fl();
}
// Variables are passed, so any array can be used.
// ACHTUNG:
// - It does fields one by one; must call it FLD_TOT times, or use the wrapper "fillAllHalos"
// - It assumes periodic, w or w/o MPI. For clarity;
// - Periodic? Call me and you're done.
// - Non-periodic MPI? Call me + appropriate BC fill.
// - Non-periodic serial? Don't call me, only BC fill.
void Domain::BCex(int myDir, Grid gr, field_array &v){ // gr is the usual local grid.
#ifdef MPICODE
int myRank; MPI_Comm_rank(cartComm_, &myRank);
MPI_Status status;
#endif
// Each direction may have a different buffer size. We may also store them, whatever.
int nBuf[]= {gr.nh[0], gr.nh[1], gr.nh[2]}; nBuf[myDir] = gr.h[myDir];
int nOff[]= { 0, 0, 0}; nOff[myDir] = gr.h[myDir]; // nOff = gr.h left corners out
range rBuf = range( nBuf[0], nBuf[1], nBuf[2]); auto sBuf = rBuf.size();
Log->setPar(true); *Log+8<<TAG<<" Filling buffers ..."; Log->fl();
// On the parallel_for: capture the buffers separately.
// Also: cumbersome to calculate rLoc (quite small anyway); with item<3> runtime decides.
qq.parallel_for( rBuf, [=, bufL=this->bufL, bufR=this->bufR ](item<3> it){
int iBufL = it.get_linear_id() , iBufR = sBuf -1 -iBufL; // The same, if we start from the end
int iVL = globLinId(it, gr.nh, nOff), iVR = gr.nht -1 -iVL ;
for(int iVar=0; iVar<FLD_TOT; ++iVar){
bufL[iVar*sBuf+iBufL] = v[iVar][iVL];
bufR[iVar*sBuf+iBufR] = v[iVar][iVR];
}
}).wait_and_throw();
#ifdef MPICODE // In serial, it's wasteful to allocate the buffers. As if someone cared about serial cases. -SC
*Log+8<<TAG<<"L+R MPI_Sendrrecv_replace ... "; Log->fl();
MPI_Sendrecv_replace(bufL,sBuf*FLD_TOT,MPI_FIELD,neighRankPrev_[myDir],myRank,neighRankNext_[myDir],MPI_ANY_TAG,cartComm_,&status);
MPI_Sendrecv_replace(bufR,sBuf*FLD_TOT,MPI_FIELD,neighRankNext_[myDir],myRank,neighRankPrev_[myDir],MPI_ANY_TAG,cartComm_,&status);
#endif
nOff[myDir] = 0; // To write, we start from 0
Log->setPar(true); *Log+8<<TAG<<" Recopying from buffers..."; Log->fl();
qq.parallel_for( rBuf, [=, bufL=this->bufL, bufR=this->bufR](item<3> it){ // v -> WHindex
int iBufL = it.get_linear_id( ), iBufR = sBuf -1 -iBufL; // The same, if we start from the end
int iVL = globLinId(it, gr.nh, nOff), iVR = gr.nht -1 -iVL ;
for(int iVar=0; iVar<FLD_TOT; ++iVar){
v[iVar][iVR] = bufL[iVar*sBuf+iBufR]; // ...besides the flipped assignments
v[iVar][iVL] = bufR[iVar*sBuf+iBufL]; // ACHTUNG: Must reverse both L<-->R and the indexes in them!
}
}).wait_and_throw();
return;
}