1 | /*! |
---|
2 | \file domain_algorithm_interpolate_from_file.cpp |
---|
3 | \author Ha NGUYEN |
---|
4 | \since 09 Jul 2015 |
---|
5 | \date 15 Sep 2015 |
---|
6 | |
---|
7 | \brief Algorithm for interpolation on a domain. |
---|
8 | */ |
---|
9 | #include "domain_algorithm_interpolate.hpp" |
---|
10 | #include <unordered_map> |
---|
11 | #include "context.hpp" |
---|
12 | #include "context_client.hpp" |
---|
13 | #include "distribution_client.hpp" |
---|
14 | #include "client_server_mapping_distributed.hpp" |
---|
15 | #include "netcdf.hpp" |
---|
16 | #include "mapper.hpp" |
---|
17 | #include "mpi_tag.hpp" |
---|
18 | #include "domain.hpp" |
---|
19 | #include "grid_transformation_factory_impl.hpp" |
---|
20 | #include "interpolate_domain.hpp" |
---|
21 | #include "grid.hpp" |
---|
22 | using namespace ep_lib; |
---|
23 | |
---|
24 | namespace xios { |
---|
25 | CGenericAlgorithmTransformation* CDomainAlgorithmInterpolate::create(CGrid* gridDst, CGrid* gridSrc, |
---|
26 | CTransformation<CDomain>* transformation, |
---|
27 | int elementPositionInGrid, |
---|
28 | std::map<int, int>& elementPositionInGridSrc2ScalarPosition, |
---|
29 | std::map<int, int>& elementPositionInGridSrc2AxisPosition, |
---|
30 | std::map<int, int>& elementPositionInGridSrc2DomainPosition, |
---|
31 | std::map<int, int>& elementPositionInGridDst2ScalarPosition, |
---|
32 | std::map<int, int>& elementPositionInGridDst2AxisPosition, |
---|
33 | std::map<int, int>& elementPositionInGridDst2DomainPosition) |
---|
34 | { |
---|
35 | std::vector<CDomain*> domainListDestP = gridDst->getDomains(); |
---|
36 | std::vector<CDomain*> domainListSrcP = gridSrc->getDomains(); |
---|
37 | |
---|
38 | CInterpolateDomain* interpolateDomain = dynamic_cast<CInterpolateDomain*> (transformation); |
---|
39 | int domainDstIndex = elementPositionInGridDst2DomainPosition[elementPositionInGrid]; |
---|
40 | int domainSrcIndex = elementPositionInGridSrc2DomainPosition[elementPositionInGrid]; |
---|
41 | |
---|
42 | return (new CDomainAlgorithmInterpolate(domainListDestP[domainDstIndex], domainListSrcP[domainSrcIndex], interpolateDomain)); |
---|
43 | } |
---|
44 | |
---|
45 | bool CDomainAlgorithmInterpolate::registerTrans() |
---|
46 | { |
---|
47 | CGridTransformationFactory<CDomain>::registerTransformation(TRANS_INTERPOLATE_DOMAIN, create); |
---|
48 | } |
---|
49 | |
---|
50 | CDomainAlgorithmInterpolate::CDomainAlgorithmInterpolate(CDomain* domainDestination, CDomain* domainSource, CInterpolateDomain* interpDomain) |
---|
51 | : CDomainAlgorithmTransformation(domainDestination, domainSource), interpDomain_(interpDomain), writeToFile_(false), readFromFile_(false) |
---|
52 | { |
---|
53 | CContext* context = CContext::getCurrent(); |
---|
54 | interpDomain_->checkValid(domainSource); |
---|
55 | |
---|
56 | detectMissingValue = interpDomain_->detect_missing_value ; |
---|
57 | renormalize = interpDomain_->renormalize ; |
---|
58 | quantity = interpDomain_->quantity ; |
---|
59 | |
---|
60 | if (interpDomain_->read_write_convention == CInterpolateDomain::read_write_convention_attr::fortran) fortranConvention=true ; |
---|
61 | else fortranConvention=false ; |
---|
62 | |
---|
63 | fileToReadWrite_ = "xios_interpolation_weights_"; |
---|
64 | |
---|
65 | if (interpDomain_->weight_filename.isEmpty()) |
---|
66 | { |
---|
67 | fileToReadWrite_ += context->getId() + "_" + |
---|
68 | domainSource->getDomainOutputName() + "_" + |
---|
69 | domainDestination->getDomainOutputName() + ".nc"; |
---|
70 | } |
---|
71 | else |
---|
72 | fileToReadWrite_ = interpDomain_->weight_filename; |
---|
73 | |
---|
74 | ifstream f(fileToReadWrite_.c_str()); |
---|
75 | switch (interpDomain_->mode) |
---|
76 | { |
---|
77 | case CInterpolateDomain::mode_attr::read: |
---|
78 | readFromFile_ = true; |
---|
79 | break; |
---|
80 | case CInterpolateDomain::mode_attr::compute: |
---|
81 | readFromFile_ = false; |
---|
82 | break; |
---|
83 | case CInterpolateDomain::mode_attr::read_or_compute: |
---|
84 | if (!f.good()) |
---|
85 | readFromFile_ = false; |
---|
86 | else |
---|
87 | readFromFile_ = true; |
---|
88 | break; |
---|
89 | default: |
---|
90 | break; |
---|
91 | } |
---|
92 | |
---|
93 | writeToFile_ = interpDomain_->write_weight; |
---|
94 | |
---|
95 | } |
---|
96 | |
---|
97 | /*! |
---|
98 | Compute remap with integrated remap calculation module |
---|
99 | */ |
---|
100 | void CDomainAlgorithmInterpolate::computeRemap() |
---|
101 | { |
---|
102 | using namespace sphereRemap; |
---|
103 | |
---|
104 | CContext* context = CContext::getCurrent(); |
---|
105 | CContextClient* client=context->client; |
---|
106 | int clientRank = client->clientRank; |
---|
107 | int i, j, k, idx; |
---|
108 | std::vector<double> srcPole(3,0), dstPole(3,0); |
---|
109 | int orderInterp = interpDomain_->order.getValue(); |
---|
110 | |
---|
111 | |
---|
112 | const double poleValue = 90.0; |
---|
113 | const int constNVertex = 4; // Value by default number of vertex for rectangular domain |
---|
114 | int nVertexSrc, nVertexDest; |
---|
115 | nVertexSrc = nVertexDest = constNVertex; |
---|
116 | |
---|
117 | // First of all, try to retrieve the boundary values of domain source and domain destination |
---|
118 | int localDomainSrcSize = domainSrc_->i_index.numElements(); |
---|
119 | int niSrc = domainSrc_->ni.getValue(), njSrc = domainSrc_->nj.getValue(); |
---|
120 | bool hasBoundSrc = domainSrc_->hasBounds; |
---|
121 | if (hasBoundSrc) nVertexSrc = domainSrc_->nvertex.getValue(); |
---|
122 | CArray<double,2> boundsLonSrc(nVertexSrc,localDomainSrcSize); |
---|
123 | CArray<double,2> boundsLatSrc(nVertexSrc,localDomainSrcSize); |
---|
124 | |
---|
125 | if (domainSrc_->hasPole) srcPole[2] = 1; |
---|
126 | if (hasBoundSrc) // Suppose that domain source is curvilinear or unstructured |
---|
127 | { |
---|
128 | if (!domainSrc_->bounds_lon_2d.isEmpty()) |
---|
129 | { |
---|
130 | for (j = 0; j < njSrc; ++j) |
---|
131 | for (i = 0; i < niSrc; ++i) |
---|
132 | { |
---|
133 | k=j*niSrc+i; |
---|
134 | for(int n=0;n<nVertexSrc;++n) |
---|
135 | { |
---|
136 | boundsLonSrc(n,k) = domainSrc_->bounds_lon_2d(n,i,j); |
---|
137 | boundsLatSrc(n,k) = domainSrc_->bounds_lat_2d(n,i,j); |
---|
138 | } |
---|
139 | } |
---|
140 | } |
---|
141 | else |
---|
142 | { |
---|
143 | boundsLonSrc = domainSrc_->bounds_lon_1d; |
---|
144 | boundsLatSrc = domainSrc_->bounds_lat_1d; |
---|
145 | } |
---|
146 | } |
---|
147 | else // if domain source is rectilinear, not do anything now |
---|
148 | { |
---|
149 | CArray<double,1> lon_g ; |
---|
150 | CArray<double,1> lat_g ; |
---|
151 | |
---|
152 | if (!domainSrc_->lonvalue_1d.isEmpty() && !domainSrc_->latvalue_1d.isEmpty()) |
---|
153 | { |
---|
154 | domainSrc_->AllgatherRectilinearLonLat(domainSrc_->lonvalue_1d,domainSrc_->latvalue_1d, lon_g,lat_g) ; |
---|
155 | } |
---|
156 | else if (! domainSrc_->latvalue_rectilinear_read_from_file.isEmpty() && ! domainSrc_->lonvalue_rectilinear_read_from_file.isEmpty() ) |
---|
157 | { |
---|
158 | lat_g=domainSrc_->latvalue_rectilinear_read_from_file ; |
---|
159 | lon_g=domainSrc_->lonvalue_rectilinear_read_from_file ; |
---|
160 | } |
---|
161 | else if (!domainSrc_->lon_start.isEmpty() && !domainSrc_->lon_end.isEmpty() && |
---|
162 | !domainSrc_->lat_start.isEmpty() && !domainSrc_->lat_end.isEmpty()) |
---|
163 | { |
---|
164 | double step=(domainSrc_->lon_end-domainSrc_->lon_start)/domainSrc_->ni_glo ; |
---|
165 | for (int i=0; i<domainSrc_->ni_glo; ++i) lon_g(i)=domainSrc_->lon_start+i*step ; |
---|
166 | step=(domainSrc_->lat_end-domainSrc_->lat_start)/domainSrc_->nj_glo ; |
---|
167 | for (int i=0; i<domainSrc_->ni_glo; ++i) lat_g(i)=domainSrc_->lat_start+i*step ; |
---|
168 | } |
---|
169 | else ERROR("void CDomainAlgorithmInterpolate::computeRemap()",<<"Cannot compute bounds for rectilinear domain") ; |
---|
170 | |
---|
171 | nVertexSrc = constNVertex; |
---|
172 | domainSrc_->fillInRectilinearBoundLonLat(lon_g,lat_g, boundsLonSrc, boundsLatSrc); |
---|
173 | } |
---|
174 | |
---|
175 | std::map<int,std::vector<std::pair<int,double> > > interpMapValueNorthPole; |
---|
176 | std::map<int,std::vector<std::pair<int,double> > > interpMapValueSouthPole; |
---|
177 | |
---|
178 | int localDomainDestSize = domainDest_->i_index.numElements(); |
---|
179 | int niDest = domainDest_->ni.getValue(), njDest = domainDest_->nj.getValue(); |
---|
180 | bool hasBoundDest = domainDest_->hasBounds; |
---|
181 | if (hasBoundDest) nVertexDest = domainDest_->nvertex.getValue(); |
---|
182 | CArray<double,2> boundsLonDest(nVertexDest,localDomainDestSize); |
---|
183 | CArray<double,2> boundsLatDest(nVertexDest,localDomainDestSize); |
---|
184 | |
---|
185 | if (domainDest_->hasPole) dstPole[2] = 1; |
---|
186 | if (hasBoundDest) |
---|
187 | { |
---|
188 | if (!domainDest_->bounds_lon_2d.isEmpty()) |
---|
189 | { |
---|
190 | for (j = 0; j < njDest; ++j) |
---|
191 | for (i = 0; i < niDest; ++i) |
---|
192 | { |
---|
193 | k=j*niDest+i; |
---|
194 | for(int n=0;n<nVertexDest;++n) |
---|
195 | { |
---|
196 | boundsLonDest(n,k) = domainDest_->bounds_lon_2d(n,i,j); |
---|
197 | boundsLatDest(n,k) = domainDest_->bounds_lat_2d(n,i,j); |
---|
198 | } |
---|
199 | } |
---|
200 | } |
---|
201 | else |
---|
202 | { |
---|
203 | boundsLonDest = domainDest_->bounds_lon_1d; |
---|
204 | boundsLatDest = domainDest_->bounds_lat_1d; |
---|
205 | } |
---|
206 | } |
---|
207 | else |
---|
208 | { |
---|
209 | bool isNorthPole = false; |
---|
210 | bool isSouthPole = false; |
---|
211 | |
---|
212 | CArray<double,1> lon_g ; |
---|
213 | CArray<double,1> lat_g ; |
---|
214 | |
---|
215 | if (!domainDest_->lonvalue_1d.isEmpty() && !domainDest_->latvalue_1d.isEmpty()) |
---|
216 | { |
---|
217 | domainDest_->AllgatherRectilinearLonLat(domainDest_->lonvalue_1d,domainDest_->latvalue_1d, lon_g,lat_g) ; |
---|
218 | } |
---|
219 | else if (! domainDest_->latvalue_rectilinear_read_from_file.isEmpty() && ! domainDest_->lonvalue_rectilinear_read_from_file.isEmpty() ) |
---|
220 | { |
---|
221 | lat_g=domainDest_->latvalue_rectilinear_read_from_file ; |
---|
222 | lon_g=domainDest_->lonvalue_rectilinear_read_from_file ; |
---|
223 | } |
---|
224 | else if (!domainDest_->lon_start.isEmpty() && !domainDest_->lon_end.isEmpty() && |
---|
225 | !domainDest_->lat_start.isEmpty() && !domainDest_->lat_end.isEmpty()) |
---|
226 | { |
---|
227 | double step=(domainDest_->lon_end-domainDest_->lon_start)/domainDest_->ni_glo ; |
---|
228 | for(int i=0; i<domainDest_->ni_glo; ++i) lon_g(i)=domainDest_->lon_start+i*step ; |
---|
229 | step=(domainDest_->lat_end-domainDest_->lat_start)/domainDest_->nj_glo ; |
---|
230 | for(int i=0; i<domainDest_->ni_glo; ++i) lat_g(i)=domainDest_->lat_start+i*step ; |
---|
231 | } |
---|
232 | else ERROR("void CDomainAlgorithmInterpolate::computeRemap()",<<"Cannot compute bounds for rectilinear domain") ; |
---|
233 | |
---|
234 | if (std::abs(poleValue - std::abs(lat_g(0))) < NumTraits<double>::epsilon()) isNorthPole = true; |
---|
235 | if (std::abs(poleValue - std::abs(lat_g(domainDest_->nj_glo-1))) < NumTraits<double>::epsilon()) isSouthPole = true; |
---|
236 | |
---|
237 | |
---|
238 | |
---|
239 | |
---|
240 | if (isNorthPole && (0 == domainDest_->jbegin.getValue())) |
---|
241 | { |
---|
242 | int ibegin = domainDest_->ibegin.getValue(); |
---|
243 | for (i = 0; i < niDest; ++i) |
---|
244 | { |
---|
245 | interpMapValueNorthPole[i+ibegin]; |
---|
246 | } |
---|
247 | } |
---|
248 | |
---|
249 | if (isSouthPole && (domainDest_->nj_glo.getValue() == (domainDest_->jbegin.getValue() + njDest))) |
---|
250 | { |
---|
251 | int ibegin = domainDest_->ibegin.getValue(); |
---|
252 | int njGlo = domainDest_->nj_glo.getValue(); |
---|
253 | int niGlo = domainDest_->ni_glo.getValue(); |
---|
254 | for (i = 0; i < niDest; ++i) |
---|
255 | { |
---|
256 | k = (njGlo - 1)*niGlo + i + ibegin; |
---|
257 | interpMapValueSouthPole[k]; |
---|
258 | } |
---|
259 | } |
---|
260 | |
---|
261 | // Ok, fill in boundary values for rectangular domain |
---|
262 | nVertexDest = constNVertex; |
---|
263 | domainDest_->fillInRectilinearBoundLonLat(lon_g,lat_g, boundsLonDest, boundsLatDest); |
---|
264 | } |
---|
265 | |
---|
266 | |
---|
267 | |
---|
268 | // Ok, now use mapper to calculate |
---|
269 | int nSrcLocal = domainSrc_->i_index.numElements(); |
---|
270 | int nDstLocal = domainDest_->i_index.numElements(); |
---|
271 | long int * globalSrc = new long int [nSrcLocal]; |
---|
272 | long int * globalDst = new long int [nDstLocal]; |
---|
273 | |
---|
274 | long int globalIndex; |
---|
275 | int i_ind, j_ind; |
---|
276 | for (int idx = 0; idx < nSrcLocal; ++idx) |
---|
277 | { |
---|
278 | i_ind=domainSrc_->i_index(idx) ; |
---|
279 | j_ind=domainSrc_->j_index(idx) ; |
---|
280 | |
---|
281 | globalIndex = i_ind + j_ind * domainSrc_->ni_glo; |
---|
282 | globalSrc[idx] = globalIndex; |
---|
283 | } |
---|
284 | |
---|
285 | for (int idx = 0; idx < nDstLocal; ++idx) |
---|
286 | { |
---|
287 | i_ind=domainDest_->i_index(idx) ; |
---|
288 | j_ind=domainDest_->j_index(idx) ; |
---|
289 | |
---|
290 | globalIndex = i_ind + j_ind * domainDest_->ni_glo; |
---|
291 | globalDst[idx] = globalIndex; |
---|
292 | } |
---|
293 | |
---|
294 | |
---|
295 | // Calculate weight index |
---|
296 | Mapper mapper(client->intraComm); |
---|
297 | mapper.setVerbosity(PROGRESS) ; |
---|
298 | |
---|
299 | |
---|
300 | // supress masked data for the source |
---|
301 | int nSrcLocalUnmasked = 0 ; |
---|
302 | for (int idx=0 ; idx < nSrcLocal; idx++) if (domainSrc_->localMask(idx)) ++nSrcLocalUnmasked ; |
---|
303 | |
---|
304 | |
---|
305 | CArray<double,2> boundsLonSrcUnmasked(nVertexSrc,nSrcLocalUnmasked); |
---|
306 | CArray<double,2> boundsLatSrcUnmasked(nVertexSrc,nSrcLocalUnmasked); |
---|
307 | long int * globalSrcUnmasked = new long int [nSrcLocalUnmasked]; |
---|
308 | |
---|
309 | nSrcLocalUnmasked=0 ; |
---|
310 | for (int idx=0 ; idx < nSrcLocal; idx++) |
---|
311 | { |
---|
312 | if (domainSrc_->localMask(idx)) |
---|
313 | { |
---|
314 | for(int n=0;n<nVertexSrc;++n) |
---|
315 | { |
---|
316 | boundsLonSrcUnmasked(n,nSrcLocalUnmasked) = boundsLonSrc(n,idx) ; |
---|
317 | boundsLatSrcUnmasked(n,nSrcLocalUnmasked) = boundsLatSrc(n,idx) ; |
---|
318 | } |
---|
319 | globalSrcUnmasked[nSrcLocalUnmasked]=globalSrc[idx] ; |
---|
320 | ++nSrcLocalUnmasked ; |
---|
321 | } |
---|
322 | } |
---|
323 | |
---|
324 | |
---|
325 | int nDstLocalUnmasked = 0 ; |
---|
326 | for (int idx=0 ; idx < nDstLocal; idx++) if (domainDest_->localMask(idx)) ++nDstLocalUnmasked ; |
---|
327 | |
---|
328 | CArray<double,2> boundsLonDestUnmasked(nVertexDest,nDstLocalUnmasked); |
---|
329 | CArray<double,2> boundsLatDestUnmasked(nVertexDest,nDstLocalUnmasked); |
---|
330 | long int * globalDstUnmasked = new long int [nDstLocalUnmasked]; |
---|
331 | |
---|
332 | nDstLocalUnmasked=0 ; |
---|
333 | for (int idx=0 ; idx < nDstLocal; idx++) |
---|
334 | { |
---|
335 | if (domainDest_->localMask(idx)) |
---|
336 | { |
---|
337 | for(int n=0;n<nVertexDest;++n) |
---|
338 | { |
---|
339 | boundsLonDestUnmasked(n,nDstLocalUnmasked) = boundsLonDest(n,idx) ; |
---|
340 | boundsLatDestUnmasked(n,nDstLocalUnmasked) = boundsLatDest(n,idx) ; |
---|
341 | } |
---|
342 | globalDstUnmasked[nDstLocalUnmasked]=globalDst[idx] ; |
---|
343 | ++nDstLocalUnmasked ; |
---|
344 | } |
---|
345 | } |
---|
346 | |
---|
347 | mapper.setSourceMesh(boundsLonSrcUnmasked.dataFirst(), boundsLatSrcUnmasked.dataFirst(), nVertexSrc, nSrcLocalUnmasked, &srcPole[0], globalSrcUnmasked); |
---|
348 | mapper.setTargetMesh(boundsLonDestUnmasked.dataFirst(), boundsLatDestUnmasked.dataFirst(), nVertexDest, nDstLocalUnmasked, &dstPole[0], globalDstUnmasked); |
---|
349 | |
---|
350 | std::vector<double> timings = mapper.computeWeights(orderInterp,renormalize,quantity); |
---|
351 | |
---|
352 | std::map<int,std::vector<std::pair<int,double> > > interpMapValue; |
---|
353 | std::map<int,std::vector<std::pair<int,double> > >::const_iterator iteNorthPole = interpMapValueNorthPole.end(), |
---|
354 | iteSouthPole = interpMapValueSouthPole.end(); |
---|
355 | for (int idx = 0; idx < mapper.nWeights; ++idx) |
---|
356 | { |
---|
357 | interpMapValue[mapper.targetWeightId[idx]].push_back(make_pair(mapper.sourceWeightId[idx],mapper.remapMatrix[idx])); |
---|
358 | if (iteNorthPole != interpMapValueNorthPole.find(mapper.targetWeightId[idx])) |
---|
359 | { |
---|
360 | interpMapValueNorthPole[mapper.targetWeightId[idx]].push_back(make_pair(mapper.sourceWeightId[idx],mapper.remapMatrix[idx])); |
---|
361 | } |
---|
362 | |
---|
363 | if (iteSouthPole != interpMapValueSouthPole.find(mapper.targetWeightId[idx])) |
---|
364 | { |
---|
365 | interpMapValueSouthPole[mapper.targetWeightId[idx]].push_back(make_pair(mapper.sourceWeightId[idx],mapper.remapMatrix[idx])); |
---|
366 | } |
---|
367 | } |
---|
368 | int niGloDst = domainDest_->ni_glo.getValue(); |
---|
369 | processPole(interpMapValueNorthPole, niGloDst); |
---|
370 | processPole(interpMapValueSouthPole, niGloDst); |
---|
371 | |
---|
372 | if (!interpMapValueNorthPole.empty()) |
---|
373 | { |
---|
374 | std::map<int,std::vector<std::pair<int,double> > >::iterator itNorthPole = interpMapValueNorthPole.begin(); |
---|
375 | for (; itNorthPole != iteNorthPole; ++itNorthPole) |
---|
376 | { |
---|
377 | if (!(itNorthPole->second.empty())) |
---|
378 | itNorthPole->second.swap(interpMapValue[itNorthPole->first]); |
---|
379 | } |
---|
380 | } |
---|
381 | |
---|
382 | if (!interpMapValueSouthPole.empty()) |
---|
383 | { |
---|
384 | std::map<int,std::vector<std::pair<int,double> > >::iterator itSouthPole = interpMapValueSouthPole.begin(); |
---|
385 | for (; itSouthPole != iteSouthPole; ++itSouthPole) |
---|
386 | { |
---|
387 | if (!(itSouthPole->second.empty())) |
---|
388 | itSouthPole->second.swap(interpMapValue[itSouthPole->first]); |
---|
389 | } |
---|
390 | } |
---|
391 | |
---|
392 | if (writeToFile_ && !readFromFile_) writeRemapInfo(interpMapValue); |
---|
393 | // exchangeRemapInfo(interpMapValue); |
---|
394 | convertRemapInfo(interpMapValue) ; |
---|
395 | |
---|
396 | delete [] globalSrc; |
---|
397 | delete [] globalSrcUnmasked; |
---|
398 | delete [] globalDst; |
---|
399 | delete [] globalDstUnmasked; |
---|
400 | |
---|
401 | } |
---|
402 | |
---|
403 | void CDomainAlgorithmInterpolate::processPole(std::map<int,std::vector<std::pair<int,double> > >& interMapValuePole, |
---|
404 | int nbGlobalPointOnPole) |
---|
405 | { |
---|
406 | CContext* context = CContext::getCurrent(); |
---|
407 | CContextClient* client=context->client; |
---|
408 | |
---|
409 | ep_lib::MPI_Comm poleComme = MPI_COMM_NULL; |
---|
410 | ep_lib::MPI_Comm_split(client->intraComm, interMapValuePole.empty() ? 0 : 1, 0, &poleComme); |
---|
411 | if (poleComme!=MPI_COMM_NULL) |
---|
412 | { |
---|
413 | int nbClientPole; |
---|
414 | ep_lib::MPI_Comm_size(poleComme, &nbClientPole); |
---|
415 | |
---|
416 | std::map<int,std::vector<std::pair<int,double> > >::iterator itePole = interMapValuePole.end(), itPole, |
---|
417 | itbPole = interMapValuePole.begin(); |
---|
418 | |
---|
419 | int nbWeight = 0; |
---|
420 | for (itPole = itbPole; itPole != itePole; ++itPole) |
---|
421 | nbWeight += itPole->second.size(); |
---|
422 | |
---|
423 | std::vector<int> recvCount(nbClientPole,0); |
---|
424 | std::vector<int> displ(nbClientPole,0); |
---|
425 | ep_lib::MPI_Allgather(&nbWeight,1,MPI_INT,&recvCount[0],1,MPI_INT,poleComme) ; |
---|
426 | displ[0]=0; |
---|
427 | for(int n=1;n<nbClientPole;++n) displ[n]=displ[n-1]+recvCount[n-1] ; |
---|
428 | int recvSize=displ[nbClientPole-1]+recvCount[nbClientPole-1] ; |
---|
429 | |
---|
430 | std::vector<int> sendSourceIndexBuff(nbWeight); |
---|
431 | std::vector<double> sendSourceWeightBuff(nbWeight); |
---|
432 | int k = 0; |
---|
433 | for (itPole = itbPole; itPole != itePole; ++itPole) |
---|
434 | { |
---|
435 | for (int idx = 0; idx < itPole->second.size(); ++idx) |
---|
436 | { |
---|
437 | sendSourceIndexBuff[k] = (itPole->second)[idx].first; |
---|
438 | sendSourceWeightBuff[k] = (itPole->second)[idx].second; |
---|
439 | ++k; |
---|
440 | } |
---|
441 | } |
---|
442 | |
---|
443 | std::vector<int> recvSourceIndexBuff(recvSize); |
---|
444 | std::vector<double> recvSourceWeightBuff(recvSize); |
---|
445 | |
---|
446 | // Gather all index and weight for pole |
---|
447 | ep_lib::MPI_Allgatherv(&sendSourceIndexBuff[0],nbWeight,MPI_INT,&recvSourceIndexBuff[0],&recvCount[0],&displ[0],MPI_INT,poleComme); |
---|
448 | ep_lib::MPI_Allgatherv(&sendSourceWeightBuff[0],nbWeight,MPI_DOUBLE,&recvSourceWeightBuff[0],&recvCount[0],&displ[0],MPI_DOUBLE,poleComme); |
---|
449 | |
---|
450 | std::map<int,double> recvTemp; |
---|
451 | for (int idx = 0; idx < recvSize; ++idx) |
---|
452 | { |
---|
453 | if (recvTemp.end() != recvTemp.find(recvSourceIndexBuff[idx])) |
---|
454 | recvTemp[recvSourceIndexBuff[idx]] += recvSourceWeightBuff[idx]/nbGlobalPointOnPole; |
---|
455 | else |
---|
456 | recvTemp[recvSourceIndexBuff[idx]] = recvSourceWeightBuff[idx]/nbGlobalPointOnPole; |
---|
457 | } |
---|
458 | |
---|
459 | std::map<int,double>::const_iterator itRecvTemp, itbRecvTemp = recvTemp.begin(), iteRecvTemp = recvTemp.end(); |
---|
460 | |
---|
461 | for (itPole = itbPole; itPole != itePole; ++itPole) |
---|
462 | { |
---|
463 | itPole->second.clear(); |
---|
464 | for (itRecvTemp = itbRecvTemp; itRecvTemp != iteRecvTemp; ++itRecvTemp) |
---|
465 | itPole->second.push_back(make_pair(itRecvTemp->first, itRecvTemp->second)); |
---|
466 | } |
---|
467 | } |
---|
468 | |
---|
469 | } |
---|
470 | |
---|
471 | /*! |
---|
472 | Compute the index mapping between domain on grid source and one on grid destination |
---|
473 | */ |
---|
474 | void CDomainAlgorithmInterpolate::computeIndexSourceMapping_(const std::vector<CArray<double,1>* >& dataAuxInputs) |
---|
475 | { |
---|
476 | if (readFromFile_) |
---|
477 | readRemapInfo(); |
---|
478 | else |
---|
479 | { |
---|
480 | computeRemap(); |
---|
481 | } |
---|
482 | } |
---|
483 | |
---|
484 | void CDomainAlgorithmInterpolate::writeRemapInfo(std::map<int,std::vector<std::pair<int,double> > >& interpMapValue) |
---|
485 | { |
---|
486 | writeInterpolationInfo(fileToReadWrite_, interpMapValue); |
---|
487 | } |
---|
488 | |
---|
489 | void CDomainAlgorithmInterpolate::readRemapInfo() |
---|
490 | { |
---|
491 | std::map<int,std::vector<std::pair<int,double> > > interpMapValue; |
---|
492 | readInterpolationInfo(fileToReadWrite_, interpMapValue); |
---|
493 | |
---|
494 | exchangeRemapInfo(interpMapValue); |
---|
495 | } |
---|
496 | |
---|
497 | void CDomainAlgorithmInterpolate::convertRemapInfo(std::map<int,std::vector<std::pair<int,double> > >& interpMapValue) |
---|
498 | { |
---|
499 | CContext* context = CContext::getCurrent(); |
---|
500 | CContextClient* client=context->client; |
---|
501 | int clientRank = client->clientRank; |
---|
502 | |
---|
503 | this->transformationMapping_.resize(1); |
---|
504 | this->transformationWeight_.resize(1); |
---|
505 | |
---|
506 | TransformationIndexMap& transMap = this->transformationMapping_[0]; |
---|
507 | TransformationWeightMap& transWeight = this->transformationWeight_[0]; |
---|
508 | |
---|
509 | std::map<int,std::vector<std::pair<int,double> > >::const_iterator itb = interpMapValue.begin(), it, |
---|
510 | ite = interpMapValue.end(); |
---|
511 | |
---|
512 | for (it = itb; it != ite; ++it) |
---|
513 | { |
---|
514 | const std::vector<std::pair<int,double> >& tmp = it->second; |
---|
515 | for (int i = 0; i < tmp.size(); ++i) |
---|
516 | { |
---|
517 | transMap[it->first].push_back(tmp[i].first); |
---|
518 | transWeight[it->first].push_back(tmp[i].second); |
---|
519 | } |
---|
520 | } |
---|
521 | } |
---|
522 | |
---|
523 | /*! |
---|
524 | Read remap information from file then distribute it among clients |
---|
525 | */ |
---|
526 | void CDomainAlgorithmInterpolate::exchangeRemapInfo(std::map<int,std::vector<std::pair<int,double> > >& interpMapValue) |
---|
527 | { |
---|
528 | CContext* context = CContext::getCurrent(); |
---|
529 | CContextClient* client=context->client; |
---|
530 | int clientRank = client->clientRank; |
---|
531 | |
---|
532 | this->transformationMapping_.resize(1); |
---|
533 | this->transformationWeight_.resize(1); |
---|
534 | |
---|
535 | TransformationIndexMap& transMap = this->transformationMapping_[0]; |
---|
536 | TransformationWeightMap& transWeight = this->transformationWeight_[0]; |
---|
537 | |
---|
538 | std::unordered_map<size_t,int> globalIndexOfDomainDest; |
---|
539 | int ni = domainDest_->ni.getValue(); |
---|
540 | int nj = domainDest_->nj.getValue(); |
---|
541 | int ni_glo = domainDest_->ni_glo.getValue(); |
---|
542 | size_t globalIndex; |
---|
543 | int nIndexSize = domainDest_->i_index.numElements(), i_ind, j_ind; |
---|
544 | for (int idx = 0; idx < nIndexSize; ++idx) |
---|
545 | { |
---|
546 | i_ind=domainDest_->i_index(idx) ; |
---|
547 | j_ind=domainDest_->j_index(idx) ; |
---|
548 | |
---|
549 | globalIndex = i_ind + j_ind * ni_glo; |
---|
550 | globalIndexOfDomainDest[globalIndex] = clientRank; |
---|
551 | } |
---|
552 | |
---|
553 | CClientServerMappingDistributed domainIndexClientClientMapping(globalIndexOfDomainDest, |
---|
554 | client->intraComm, |
---|
555 | true); |
---|
556 | CArray<size_t,1> globalIndexInterp(interpMapValue.size()); |
---|
557 | std::map<int,std::vector<std::pair<int,double> > >::const_iterator itb = interpMapValue.begin(), it, |
---|
558 | ite = interpMapValue.end(); |
---|
559 | size_t globalIndexCount = 0; |
---|
560 | for (it = itb; it != ite; ++it) |
---|
561 | { |
---|
562 | globalIndexInterp(globalIndexCount) = it->first; |
---|
563 | ++globalIndexCount; |
---|
564 | } |
---|
565 | |
---|
566 | domainIndexClientClientMapping.computeServerIndexMapping(globalIndexInterp, client->serverSize); |
---|
567 | const CClientServerMapping::GlobalIndexMap& globalIndexInterpSendToClient = domainIndexClientClientMapping.getGlobalIndexOnServer(); |
---|
568 | |
---|
569 | //Inform each client number of index they will receive |
---|
570 | int nbClient = client->clientSize; |
---|
571 | int* sendBuff = new int[nbClient]; |
---|
572 | int* recvBuff = new int[nbClient]; |
---|
573 | for (int i = 0; i < nbClient; ++i) |
---|
574 | { |
---|
575 | sendBuff[i] = 0; |
---|
576 | recvBuff[i] = 0; |
---|
577 | } |
---|
578 | int sendBuffSize = 0; |
---|
579 | CClientServerMapping::GlobalIndexMap::const_iterator itbMap = globalIndexInterpSendToClient.begin(), itMap, |
---|
580 | iteMap = globalIndexInterpSendToClient.end(); |
---|
581 | for (itMap = itbMap; itMap != iteMap; ++itMap) |
---|
582 | { |
---|
583 | const std::vector<size_t>& tmp = itMap->second; |
---|
584 | int sizeIndex = 0, mapSize = (itMap->second).size(); |
---|
585 | for (int idx = 0; idx < mapSize; ++idx) |
---|
586 | { |
---|
587 | // sizeIndex += interpMapValue.at((itMap->second)[idx]).size(); |
---|
588 | sizeIndex += (interpMapValue[(int)(itMap->second)[idx]]).size(); |
---|
589 | } |
---|
590 | sendBuff[itMap->first] = sizeIndex; |
---|
591 | sendBuffSize += sizeIndex; |
---|
592 | } |
---|
593 | |
---|
594 | |
---|
595 | ep_lib::MPI_Allreduce(sendBuff, recvBuff, nbClient, MPI_INT, MPI_SUM, client->intraComm); |
---|
596 | |
---|
597 | int* sendIndexDestBuff = new int [sendBuffSize]; |
---|
598 | int* sendIndexSrcBuff = new int [sendBuffSize]; |
---|
599 | double* sendWeightBuff = new double [sendBuffSize]; |
---|
600 | |
---|
601 | std::vector<ep_lib::MPI_Request> sendRequest(3*globalIndexInterpSendToClient.size()); |
---|
602 | |
---|
603 | int sendOffSet = 0, l = 0; |
---|
604 | int position = 0; |
---|
605 | for (itMap = itbMap; itMap != iteMap; ++itMap) |
---|
606 | { |
---|
607 | const std::vector<size_t>& indexToSend = itMap->second; |
---|
608 | int mapSize = indexToSend.size(); |
---|
609 | int k = 0; |
---|
610 | for (int idx = 0; idx < mapSize; ++idx) |
---|
611 | { |
---|
612 | std::vector<std::pair<int,double> >& interpMap = interpMapValue[(int)indexToSend[idx]]; //interpMapValue.at(indexToSend[idx]); |
---|
613 | for (int i = 0; i < interpMap.size(); ++i) |
---|
614 | { |
---|
615 | sendIndexDestBuff[l] = indexToSend[idx]; |
---|
616 | sendIndexSrcBuff[l] = interpMap[i].first; |
---|
617 | sendWeightBuff[l] = interpMap[i].second; |
---|
618 | ++k; |
---|
619 | ++l; |
---|
620 | } |
---|
621 | } |
---|
622 | |
---|
623 | ep_lib::MPI_Isend(sendIndexDestBuff + sendOffSet, |
---|
624 | k, |
---|
625 | MPI_INT, |
---|
626 | itMap->first, |
---|
627 | MPI_DOMAIN_INTERPOLATION_DEST_INDEX, |
---|
628 | client->intraComm, |
---|
629 | &sendRequest[position++]); |
---|
630 | ep_lib::MPI_Isend(sendIndexSrcBuff + sendOffSet, |
---|
631 | k, |
---|
632 | MPI_INT, |
---|
633 | itMap->first, |
---|
634 | MPI_DOMAIN_INTERPOLATION_SRC_INDEX, |
---|
635 | client->intraComm, |
---|
636 | &sendRequest[position++]); |
---|
637 | ep_lib::MPI_Isend(sendWeightBuff + sendOffSet, |
---|
638 | k, |
---|
639 | MPI_DOUBLE, |
---|
640 | itMap->first, |
---|
641 | MPI_DOMAIN_INTERPOLATION_WEIGHT, |
---|
642 | client->intraComm, |
---|
643 | &sendRequest[position++]); |
---|
644 | sendOffSet += k; |
---|
645 | } |
---|
646 | |
---|
647 | int recvBuffSize = recvBuff[clientRank]; |
---|
648 | int* recvIndexDestBuff = new int [recvBuffSize]; |
---|
649 | int* recvIndexSrcBuff = new int [recvBuffSize]; |
---|
650 | double* recvWeightBuff = new double [recvBuffSize]; |
---|
651 | int receivedSize = 0; |
---|
652 | int clientSrcRank; |
---|
653 | while (receivedSize < recvBuffSize) |
---|
654 | { |
---|
655 | ep_lib::MPI_Status recvStatus; |
---|
656 | ep_lib::MPI_Recv((recvIndexDestBuff + receivedSize), |
---|
657 | recvBuffSize, |
---|
658 | MPI_INT, |
---|
659 | -2, |
---|
660 | MPI_DOMAIN_INTERPOLATION_DEST_INDEX, |
---|
661 | client->intraComm, |
---|
662 | &recvStatus); |
---|
663 | |
---|
664 | int countBuff = 0; |
---|
665 | ep_lib::MPI_Get_count(&recvStatus, MPI_INT, &countBuff); |
---|
666 | #ifdef _usingMPI |
---|
667 | clientSrcRank = recvStatus.MPI_SOURCE; |
---|
668 | #elif _usingEP |
---|
669 | clientSrcRank = recvStatus.ep_src; |
---|
670 | #endif |
---|
671 | |
---|
672 | ep_lib::MPI_Recv((recvIndexSrcBuff + receivedSize), |
---|
673 | recvBuffSize, |
---|
674 | MPI_INT, |
---|
675 | clientSrcRank, |
---|
676 | MPI_DOMAIN_INTERPOLATION_SRC_INDEX, |
---|
677 | client->intraComm, |
---|
678 | &recvStatus); |
---|
679 | |
---|
680 | ep_lib::MPI_Recv((recvWeightBuff + receivedSize), |
---|
681 | recvBuffSize, |
---|
682 | MPI_DOUBLE, |
---|
683 | clientSrcRank, |
---|
684 | MPI_DOMAIN_INTERPOLATION_WEIGHT, |
---|
685 | client->intraComm, |
---|
686 | &recvStatus); |
---|
687 | |
---|
688 | for (int idx = 0; idx < countBuff; ++idx) |
---|
689 | { |
---|
690 | transMap[*(recvIndexDestBuff + receivedSize + idx)].push_back(*(recvIndexSrcBuff + receivedSize + idx)); |
---|
691 | transWeight[*(recvIndexDestBuff + receivedSize + idx)].push_back(*(recvWeightBuff + receivedSize + idx)); |
---|
692 | } |
---|
693 | receivedSize += countBuff; |
---|
694 | } |
---|
695 | |
---|
696 | std::vector<ep_lib::MPI_Status> requestStatus(sendRequest.size()); |
---|
697 | ep_lib::MPI_Waitall(sendRequest.size(), &sendRequest[0], &requestStatus[0]); |
---|
698 | |
---|
699 | delete [] sendIndexDestBuff; |
---|
700 | delete [] sendIndexSrcBuff; |
---|
701 | delete [] sendWeightBuff; |
---|
702 | delete [] recvIndexDestBuff; |
---|
703 | delete [] recvIndexSrcBuff; |
---|
704 | delete [] recvWeightBuff; |
---|
705 | delete [] sendBuff; |
---|
706 | delete [] recvBuff; |
---|
707 | } |
---|
708 | |
---|
709 | /*! Redefined some functions of CONetCDF4 to make use of them */ |
---|
710 | CDomainAlgorithmInterpolate::WriteNetCdf::WriteNetCdf(const StdString& filename, const ep_lib::MPI_Comm comm) |
---|
711 | : CNc4DataOutput(NULL, filename, false, false, true, comm, false, true) {} |
---|
712 | |
---|
713 | CDomainAlgorithmInterpolate::WriteNetCdf::WriteNetCdf(const StdString& filename, bool exist, const ep_lib::MPI_Comm comm) |
---|
714 | : CNc4DataOutput(NULL, filename, exist, false, true, comm, false, true) {} |
---|
715 | |
---|
716 | int CDomainAlgorithmInterpolate::WriteNetCdf::addDimensionWrite(const StdString& name, |
---|
717 | const StdSize size) |
---|
718 | { |
---|
719 | return CONetCDF4::addDimension(name, size); |
---|
720 | } |
---|
721 | |
---|
722 | int CDomainAlgorithmInterpolate::WriteNetCdf::addVariableWrite(const StdString& name, nc_type type, |
---|
723 | const std::vector<StdString>& dim) |
---|
724 | { |
---|
725 | return CONetCDF4::addVariable(name, type, dim); |
---|
726 | } |
---|
727 | |
---|
728 | void CDomainAlgorithmInterpolate::WriteNetCdf::endDefinition() |
---|
729 | { |
---|
730 | CONetCDF4::definition_end(); |
---|
731 | } |
---|
732 | |
---|
733 | void CDomainAlgorithmInterpolate::WriteNetCdf::writeDataIndex(const CArray<int,1>& data, const StdString& name, |
---|
734 | bool collective, StdSize record, |
---|
735 | const std::vector<StdSize>* start, |
---|
736 | const std::vector<StdSize>* count) |
---|
737 | { |
---|
738 | CONetCDF4::writeData<int,1>(data, name, collective, record, start, count); |
---|
739 | } |
---|
740 | |
---|
741 | void CDomainAlgorithmInterpolate::WriteNetCdf::writeDataIndex(const CArray<double,1>& data, const StdString& name, |
---|
742 | bool collective, StdSize record, |
---|
743 | const std::vector<StdSize>* start, |
---|
744 | const std::vector<StdSize>* count) |
---|
745 | { |
---|
746 | CONetCDF4::writeData<double,1>(data, name, collective, record, start, count); |
---|
747 | } |
---|
748 | |
---|
749 | /* |
---|
750 | Write interpolation weights into a file |
---|
751 | \param [in] filename name of output file |
---|
752 | \param interpMapValue mapping of global index of domain destination and domain source as well as the corresponding weight |
---|
753 | */ |
---|
754 | void CDomainAlgorithmInterpolate::writeInterpolationInfo(std::string& filename, |
---|
755 | std::map<int,std::vector<std::pair<int,double> > >& interpMapValue) |
---|
756 | { |
---|
757 | CContext* context = CContext::getCurrent(); |
---|
758 | CContextClient* client=context->client; |
---|
759 | |
---|
760 | size_t n_src = domainSrc_->ni_glo * domainSrc_->nj_glo; |
---|
761 | size_t n_dst = domainDest_->ni_glo * domainDest_->nj_glo; |
---|
762 | |
---|
763 | long localNbWeight = 0; |
---|
764 | long globalNbWeight; |
---|
765 | long startIndex; |
---|
766 | typedef std::map<int,std::vector<std::pair<int,double> > > IndexRemap; |
---|
767 | IndexRemap::iterator itb = interpMapValue.begin(), it, |
---|
768 | ite = interpMapValue.end(); |
---|
769 | for (it = itb; it!=ite; ++it) |
---|
770 | { |
---|
771 | localNbWeight += (it->second).size(); |
---|
772 | } |
---|
773 | |
---|
774 | CArray<int,1> src_idx(localNbWeight); |
---|
775 | CArray<int,1> dst_idx(localNbWeight); |
---|
776 | CArray<double,1> weights(localNbWeight); |
---|
777 | |
---|
778 | int index = 0; |
---|
779 | int indexOffset=0 ; |
---|
780 | if (fortranConvention) indexOffset=1 ; |
---|
781 | for (it = itb; it !=ite; ++it) |
---|
782 | { |
---|
783 | std::vector<std::pair<int,double> >& tmp = it->second; |
---|
784 | for (int idx = 0; idx < tmp.size(); ++idx) |
---|
785 | { |
---|
786 | dst_idx(index) = it->first + indexOffset; |
---|
787 | src_idx(index) = tmp[idx].first + indexOffset; |
---|
788 | weights(index) = tmp[idx].second; |
---|
789 | ++index; |
---|
790 | } |
---|
791 | } |
---|
792 | |
---|
793 | ep_lib::MPI_Allreduce(&localNbWeight, &globalNbWeight, 1, MPI_LONG, MPI_SUM, client->intraComm); |
---|
794 | ep_lib::MPI_Scan(&localNbWeight, &startIndex, 1, MPI_LONG, MPI_SUM, client->intraComm); |
---|
795 | |
---|
796 | if (0 == globalNbWeight) |
---|
797 | { |
---|
798 | info << "There is no interpolation weights calculated between " |
---|
799 | << "domain source: " << domainSrc_->getDomainOutputName() |
---|
800 | << " and domain destination: " << domainDest_->getDomainOutputName() |
---|
801 | << std::endl; |
---|
802 | return; |
---|
803 | } |
---|
804 | |
---|
805 | std::vector<StdSize> start(1, startIndex - localNbWeight); |
---|
806 | std::vector<StdSize> count(1, localNbWeight); |
---|
807 | |
---|
808 | int my_rank_loc = client->intraComm->ep_comm_ptr->size_rank_info[1].first; |
---|
809 | int my_rank = client->intraComm->ep_comm_ptr->size_rank_info[0].first; |
---|
810 | |
---|
811 | |
---|
812 | |
---|
813 | WriteNetCdf *netCdfWriter; |
---|
814 | |
---|
815 | MPI_Barrier_local(client->intraComm); |
---|
816 | |
---|
817 | if(my_rank_loc==0) |
---|
818 | { |
---|
819 | info(100)<<"rank "<< my_rank <<" create weight info file"<< std::endl; |
---|
820 | |
---|
821 | WriteNetCdf my_writer(filename, client->intraComm); |
---|
822 | info(100)<<"rank "<< my_rank <<" file created"<< std::endl; |
---|
823 | netCdfWriter = &my_writer; |
---|
824 | |
---|
825 | // Define some dimensions |
---|
826 | netCdfWriter->addDimensionWrite("n_src", n_src); |
---|
827 | netCdfWriter->addDimensionWrite("n_dst", n_dst); |
---|
828 | netCdfWriter->addDimensionWrite("n_weight", globalNbWeight); |
---|
829 | info(100)<<"rank "<< my_rank <<" addDimensionWrite : n_src, n_dst, n_weight"<< std::endl; |
---|
830 | |
---|
831 | std::vector<StdString> dims(1,"n_weight"); |
---|
832 | |
---|
833 | // Add some variables |
---|
834 | netCdfWriter->addVariableWrite("src_idx", NC_INT, dims); |
---|
835 | netCdfWriter->addVariableWrite("dst_idx", NC_INT, dims); |
---|
836 | netCdfWriter->addVariableWrite("weight", NC_DOUBLE, dims); |
---|
837 | |
---|
838 | info(100)<<"rank "<< my_rank <<" addVariableWrite : src_idx, dst_idx, weight"<< std::endl; |
---|
839 | |
---|
840 | // End of definition |
---|
841 | netCdfWriter->endDefinition(); |
---|
842 | info(100)<<"rank "<< my_rank <<" endDefinition"<< std::endl; |
---|
843 | |
---|
844 | netCdfWriter->closeFile(); |
---|
845 | info(100)<<"rank "<< my_rank <<" file closed"<< std::endl; |
---|
846 | } |
---|
847 | |
---|
848 | MPI_Barrier_local(client->intraComm); |
---|
849 | |
---|
850 | #pragma omp critical (write_weight_data) |
---|
851 | { |
---|
852 | // open file |
---|
853 | info(100)<<"rank "<< my_rank <<" writing in weight info file"<< std::endl; |
---|
854 | |
---|
855 | WriteNetCdf my_writer(filename, true, client->intraComm); |
---|
856 | info(100)<<"rank "<< my_rank <<" file opened"<< std::endl; |
---|
857 | netCdfWriter = &my_writer; |
---|
858 | |
---|
859 | // // Write variables |
---|
860 | if (0 != localNbWeight) |
---|
861 | { |
---|
862 | netCdfWriter->writeDataIndex(src_idx, "src_idx", false, 0, &start, &count); |
---|
863 | netCdfWriter->writeDataIndex(dst_idx, "dst_idx", false, 0, &start, &count); |
---|
864 | netCdfWriter->writeDataIndex(weights, "weight", false, 0, &start, &count); |
---|
865 | |
---|
866 | info(100)<<"rank "<< my_rank <<" WriteDataIndex : src_idx, dst_idx, weight"<< std::endl; |
---|
867 | } |
---|
868 | |
---|
869 | netCdfWriter->closeFile(); |
---|
870 | info(100)<<"rank "<< my_rank <<" file closed"<< std::endl; |
---|
871 | |
---|
872 | } |
---|
873 | |
---|
874 | MPI_Barrier_local(client->intraComm); |
---|
875 | |
---|
876 | |
---|
877 | } |
---|
878 | |
---|
879 | /*! |
---|
880 | Read interpolation information from a file |
---|
881 | \param [in] filename interpolation file |
---|
882 | \param [in/out] interpMapValue Mapping between (global) index of domain on grid destination and |
---|
883 | corresponding global index of domain and associated weight value on grid source |
---|
884 | */ |
---|
885 | void CDomainAlgorithmInterpolate::readInterpolationInfo(std::string& filename, |
---|
886 | std::map<int,std::vector<std::pair<int,double> > >& interpMapValue) |
---|
887 | { |
---|
888 | int ncid ; |
---|
889 | int weightDimId ; |
---|
890 | size_t nbWeightGlo ; |
---|
891 | |
---|
892 | |
---|
893 | CContext* context = CContext::getCurrent(); |
---|
894 | CContextClient* client=context->client; |
---|
895 | int clientRank = client->clientRank; |
---|
896 | int clientSize = client->clientSize; |
---|
897 | |
---|
898 | |
---|
899 | { |
---|
900 | ifstream f(filename.c_str()); |
---|
901 | if (!f.good()) ERROR("void CDomainAlgorithmInterpolate::readInterpolationInfo", |
---|
902 | << "Attempt to read file weight :" << filename << " which doesn't seem to exist." << std::endl |
---|
903 | << "Please check this file "); |
---|
904 | } |
---|
905 | |
---|
906 | #pragma omp critical (read_weight_data) |
---|
907 | { |
---|
908 | nc_open(filename.c_str(),NC_NOWRITE, &ncid) ; |
---|
909 | nc_inq_dimid(ncid,"n_weight",&weightDimId) ; |
---|
910 | nc_inq_dimlen(ncid,weightDimId,&nbWeightGlo) ; |
---|
911 | |
---|
912 | size_t nbWeight ; |
---|
913 | size_t start ; |
---|
914 | size_t div = nbWeightGlo/clientSize ; |
---|
915 | size_t mod = nbWeightGlo%clientSize ; |
---|
916 | if (clientRank < mod) |
---|
917 | { |
---|
918 | nbWeight=div+1 ; |
---|
919 | start=clientRank*(div+1) ; |
---|
920 | } |
---|
921 | else |
---|
922 | { |
---|
923 | nbWeight=div ; |
---|
924 | start= mod * (div+1) + (clientRank-mod) * div ; |
---|
925 | } |
---|
926 | |
---|
927 | double* weight=new double[nbWeight] ; |
---|
928 | int weightId ; |
---|
929 | nc_inq_varid (ncid, "weight", &weightId) ; |
---|
930 | nc_get_vara_double(ncid, weightId, &start, &nbWeight, weight) ; |
---|
931 | |
---|
932 | long* srcIndex=new long[nbWeight] ; |
---|
933 | int srcIndexId ; |
---|
934 | nc_inq_varid (ncid, "src_idx", &srcIndexId) ; |
---|
935 | nc_get_vara_long(ncid, srcIndexId, &start, &nbWeight, srcIndex) ; |
---|
936 | |
---|
937 | long* dstIndex=new long[nbWeight] ; |
---|
938 | int dstIndexId ; |
---|
939 | nc_inq_varid (ncid, "dst_idx", &dstIndexId) ; |
---|
940 | nc_get_vara_long(ncid, dstIndexId, &start, &nbWeight, dstIndex) ; |
---|
941 | |
---|
942 | int indexOffset=0 ; |
---|
943 | if (fortranConvention) indexOffset=1 ; |
---|
944 | for(size_t ind=0; ind<nbWeight;++ind) |
---|
945 | interpMapValue[dstIndex[ind]-indexOffset].push_back(make_pair(srcIndex[ind]-indexOffset,weight[ind])); |
---|
946 | } |
---|
947 | } |
---|
948 | |
---|
949 | void CDomainAlgorithmInterpolate::apply(const std::vector<std::pair<int,double> >& localIndex, |
---|
950 | const double* dataInput, |
---|
951 | CArray<double,1>& dataOut, |
---|
952 | std::vector<bool>& flagInitial, |
---|
953 | bool ignoreMissingValue, bool firstPass ) |
---|
954 | { |
---|
955 | int nbLocalIndex = localIndex.size(); |
---|
956 | double defaultValue = std::numeric_limits<double>::quiet_NaN(); |
---|
957 | |
---|
958 | if (detectMissingValue) |
---|
959 | { |
---|
960 | if (firstPass && renormalize) |
---|
961 | { |
---|
962 | renormalizationFactor.resize(dataOut.numElements()) ; |
---|
963 | renormalizationFactor=1 ; |
---|
964 | } |
---|
965 | |
---|
966 | if (firstPass) |
---|
967 | { |
---|
968 | allMissing.resize(dataOut.numElements()) ; |
---|
969 | allMissing=true ; |
---|
970 | } |
---|
971 | |
---|
972 | for (int idx = 0; idx < nbLocalIndex; ++idx) |
---|
973 | { |
---|
974 | if (NumTraits<double>::isNan(*(dataInput + idx))) |
---|
975 | { |
---|
976 | allMissing(localIndex[idx].first) = allMissing(localIndex[idx].first) && true; |
---|
977 | if (renormalize) renormalizationFactor(localIndex[idx].first)-=localIndex[idx].second ; |
---|
978 | } |
---|
979 | else |
---|
980 | { |
---|
981 | dataOut(localIndex[idx].first) += *(dataInput + idx) * localIndex[idx].second; |
---|
982 | allMissing(localIndex[idx].first) = allMissing(localIndex[idx].first) && false; // Reset flag to indicate not all data source are nan |
---|
983 | } |
---|
984 | } |
---|
985 | |
---|
986 | } |
---|
987 | else |
---|
988 | { |
---|
989 | for (int idx = 0; idx < nbLocalIndex; ++idx) |
---|
990 | { |
---|
991 | dataOut(localIndex[idx].first) += *(dataInput + idx) * localIndex[idx].second; |
---|
992 | } |
---|
993 | } |
---|
994 | } |
---|
995 | |
---|
996 | void CDomainAlgorithmInterpolate::updateData(CArray<double,1>& dataOut) |
---|
997 | { |
---|
998 | if (detectMissingValue) |
---|
999 | { |
---|
1000 | double defaultValue = std::numeric_limits<double>::quiet_NaN(); |
---|
1001 | size_t nbIndex=dataOut.numElements() ; |
---|
1002 | |
---|
1003 | for (int idx = 0; idx < nbIndex; ++idx) |
---|
1004 | { |
---|
1005 | if (allMissing(idx)) dataOut(idx) = defaultValue; // If all data source are nan then data destination must be nan |
---|
1006 | } |
---|
1007 | |
---|
1008 | if (renormalize) |
---|
1009 | { |
---|
1010 | if (renormalizationFactor.numElements()>0) dataOut/=renormalizationFactor ; // In some case, process doesn't received any data for interpolation (mask) |
---|
1011 | // so renormalizationFactor is not initialized |
---|
1012 | } |
---|
1013 | } |
---|
1014 | } |
---|
1015 | |
---|
1016 | } |
---|