Changeset 2613 for XIOS3/trunk/src/io/onetcdf4.cpp
- Timestamp:
- 03/08/24 17:05:40 (3 months ago)
- File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
XIOS3/trunk/src/io/onetcdf4.cpp
r2600 r2613 478 478 normalizingWeight = userChunkingWeights[i]; 479 479 } 480 481 std::vector<double> chunkingRatioPerDims; // will store coefficients used to compute chunk size 482 double productRatios = 1; // last_coeff = pow( shrink_ratio / (product of all ratios), 1/countChunkingDims ) 483 for (int i=0;i<userChunkingWeights.size();i++) 480 if (normalizingWeight!=0) // no chunk for scalar 484 481 { 485 chunkingRatioPerDims.push_back( userChunkingWeights[i] / normalizingWeight ); 486 if (chunkingRatioPerDims[i]) productRatios *= chunkingRatioPerDims[i]; 487 } 488 for (int i=0;i<userChunkingWeights.size();i++) 489 { 490 chunkingRatioPerDims[i] *= pow( chunkingRatio / productRatios, 1./countChunkingDims ); 491 } 492 493 std::vector<double>::iterator itChunkingRatios = chunkingRatioPerDims.begin(); 494 //itId = dim.rbegin(); 495 double correctionFromPreviousDim = 1.; 496 for (vector<StdSize>::reverse_iterator itDim = dimsizes.rbegin(); itDim != dimsizes.rend(); ++itDim, ++itChunkingRatios, ++itId) 497 { 498 *itChunkingRatios *= correctionFromPreviousDim; 499 correctionFromPreviousDim = 1; 500 if (*itChunkingRatios > 1) // else target larger than size ! 482 std::vector<double> chunkingRatioPerDims; // will store coefficients used to compute chunk size 483 double productRatios = 1; // last_coeff = pow( shrink_ratio / (product of all ratios), 1/countChunkingDims ) 484 for (int i=0;i<userChunkingWeights.size();i++) 501 485 { 502 StdSize dimensionSize = *itDim; 503 //info(0) << *itId << " " << *itDim << " " << *itChunkingRatios << " " << (*itDim)/(*itChunkingRatios) << endl; 504 *itDim = ceil( *itDim / ceil(*itChunkingRatios) ); 505 correctionFromPreviousDim = *itChunkingRatios/ ((double)dimensionSize/(*itDim)); 486 chunkingRatioPerDims.push_back( userChunkingWeights[i] / normalizingWeight ); 487 if (chunkingRatioPerDims[i]) productRatios *= chunkingRatioPerDims[i]; 506 488 } 507 } 489 for (int i=0;i<userChunkingWeights.size();i++) 490 { 491 chunkingRatioPerDims[i] *= pow( chunkingRatio / productRatios, 1./countChunkingDims ); 492 } 493 494 std::vector<double>::iterator itChunkingRatios = chunkingRatioPerDims.begin(); 495 //itId = dim.rbegin(); 496 double correctionFromPreviousDim = 1.; 497 for (vector<StdSize>::reverse_iterator itDim = dimsizes.rbegin(); itDim != dimsizes.rend(); ++itDim, ++itChunkingRatios, ++itId) 498 { 499 *itChunkingRatios *= correctionFromPreviousDim; 500 correctionFromPreviousDim = 1; 501 if (*itChunkingRatios > 1) // else target larger than size ! 502 { 503 StdSize dimensionSize = *itDim; 504 //info(0) << *itId << " " << *itDim << " " << *itChunkingRatios << " " << (*itDim)/(*itChunkingRatios) << endl; 505 *itDim = ceil( *itDim / ceil(*itChunkingRatios) ); 506 correctionFromPreviousDim = *itChunkingRatios/ ((double)dimensionSize/(*itDim)); 507 } 508 } 509 } 508 510 int storageType = (0 == dimSize) ? NC_CONTIGUOUS : NC_CHUNKED; 509 511 CNetCdfInterface::defVarChunking(grpid, varid, storageType, &dimsizes[0]);
Note: See TracChangeset
for help on using the changeset viewer.