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