- Timestamp:
- 01/10/17 14:36:29 (7 years ago)
- File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
XIOS/dev/dev_olga/src/transformation/domain_algorithm_interpolate.cpp
r978 r1021 48 48 49 49 CDomainAlgorithmInterpolate::CDomainAlgorithmInterpolate(CDomain* domainDestination, CDomain* domainSource, CInterpolateDomain* interpDomain) 50 : CDomainAlgorithmTransformation(domainDestination, domainSource), interpDomain_(interpDomain) 51 { 50 : CDomainAlgorithmTransformation(domainDestination, domainSource), interpDomain_(interpDomain), writeToFile_(false), readFromFile_(false) 51 { 52 CContext* context = CContext::getCurrent(); 52 53 interpDomain_->checkValid(domainSource); 54 fileToReadWrite_ = "xios_interpolation_weights_"; 55 56 if (interpDomain_->weight_filename.isEmpty()) 57 { 58 fileToReadWrite_ += context->getId() + "_" + 59 domainSource->getDomainOutputName() + "_" + 60 domainDestination->getDomainOutputName() + ".nc"; 61 } 62 else 63 fileToReadWrite_ = interpDomain_->weight_filename; 64 65 ifstream f(fileToReadWrite_.c_str()); 66 switch (interpDomain_->mode) 67 { 68 case CInterpolateDomain::mode_attr::read: 69 readFromFile_ = true; 70 break; 71 case CInterpolateDomain::mode_attr::compute: 72 readFromFile_ = false; 73 break; 74 case CInterpolateDomain::mode_attr::read_or_compute: 75 if (!f.good()) 76 readFromFile_ = false; 77 else 78 readFromFile_ = true; 79 break; 80 default: 81 break; 82 } 83 84 writeToFile_ = interpDomain_->write_weight; 85 53 86 } 54 87 … … 351 384 } 352 385 386 if (writeToFile_ && !readFromFile_) 387 writeRemapInfo(interpMapValue); 353 388 exchangeRemapInfo(interpMapValue); 354 389 … … 434 469 void CDomainAlgorithmInterpolate::computeIndexSourceMapping_(const std::vector<CArray<double,1>* >& dataAuxInputs) 435 470 { 436 if ( !interpDomain_->file.isEmpty())471 if (readFromFile_ && !writeToFile_) 437 472 readRemapInfo(); 438 473 else 439 computeRemap(); 474 { 475 computeRemap(); 476 } 477 } 478 479 void CDomainAlgorithmInterpolate::writeRemapInfo(std::map<int,std::vector<std::pair<int,double> > >& interpMapValue) 480 { 481 writeInterpolationInfo(fileToReadWrite_, interpMapValue); 440 482 } 441 483 442 484 void CDomainAlgorithmInterpolate::readRemapInfo() 443 { 444 CContext* context = CContext::getCurrent(); 445 CContextClient* client=context->client; 446 int clientRank = client->clientRank; 447 448 std::string filename = interpDomain_->file.getValue(); 485 { 449 486 std::map<int,std::vector<std::pair<int,double> > > interpMapValue; 450 readInterpolationInfo(file name, interpMapValue);487 readInterpolationInfo(fileToReadWrite_, interpMapValue); 451 488 452 489 exchangeRemapInfo(interpMapValue); … … 637 674 delete [] recvBuff; 638 675 } 676 677 /*! Redefined some functions of CONetCDF4 to make use of them */ 678 CDomainAlgorithmInterpolate::WriteNetCdf::WriteNetCdf(const StdString& filename, const MPI_Comm comm) 679 : CNc4DataOutput(filename, false, false, true, comm, false, true) {} 680 int CDomainAlgorithmInterpolate::WriteNetCdf::addDimensionWrite(const StdString& name, 681 const StdSize size) 682 { 683 CONetCDF4::addDimension(name, size); 684 } 685 686 int CDomainAlgorithmInterpolate::WriteNetCdf::addVariableWrite(const StdString& name, nc_type type, 687 const std::vector<StdString>& dim) 688 { 689 CONetCDF4::addVariable(name, type, dim); 690 } 691 692 void CDomainAlgorithmInterpolate::WriteNetCdf::writeDataIndex(const CArray<int,1>& data, const StdString& name, 693 bool collective, StdSize record, 694 const std::vector<StdSize>* start, 695 const std::vector<StdSize>* count) 696 { 697 CONetCDF4::writeData<int,1>(data, name, collective, record, start, count); 698 } 699 700 void CDomainAlgorithmInterpolate::WriteNetCdf::writeDataIndex(const CArray<double,1>& data, const StdString& name, 701 bool collective, StdSize record, 702 const std::vector<StdSize>* start, 703 const std::vector<StdSize>* count) 704 { 705 CONetCDF4::writeData<double,1>(data, name, collective, record, start, count); 706 } 707 708 /* 709 Write interpolation weights into a file 710 \param [in] filename name of output file 711 \param interpMapValue mapping of global index of domain destination and domain source as well as the corresponding weight 712 */ 713 void CDomainAlgorithmInterpolate::writeInterpolationInfo(std::string& filename, 714 std::map<int,std::vector<std::pair<int,double> > >& interpMapValue) 715 { 716 CContext* context = CContext::getCurrent(); 717 CContextClient* client=context->client; 718 719 size_t n_src = domainSrc_->ni_glo * domainSrc_->nj_glo; 720 size_t n_dst = domainDest_->ni_glo * domainDest_->nj_glo; 721 722 long localNbWeight = 0; 723 long globalNbWeight; 724 long startIndex; 725 typedef std::map<int,std::vector<std::pair<int,double> > > IndexRemap; 726 IndexRemap::iterator itb = interpMapValue.begin(), it, 727 ite = interpMapValue.end(); 728 for (it = itb; it!=ite; ++it) 729 { 730 localNbWeight += (it->second).size(); 731 } 732 733 CArray<int,1> src_idx(localNbWeight); 734 CArray<int,1> dst_idx(localNbWeight); 735 CArray<double,1> weights(localNbWeight); 736 737 int index = 0; 738 for (it = itb; it !=ite; ++it) 739 { 740 std::vector<std::pair<int,double> >& tmp = it->second; 741 for (int idx = 0; idx < tmp.size(); ++idx) 742 { 743 dst_idx(index) = it->first + 1; 744 src_idx(index) = tmp[idx].first + 1; 745 weights(index) = tmp[idx].second; 746 ++index; 747 } 748 } 749 750 MPI_Allreduce(&localNbWeight, &globalNbWeight, 1, MPI_LONG, MPI_SUM, client->intraComm); 751 MPI_Scan(&localNbWeight, &startIndex, 1, MPI_LONG, MPI_SUM, client->intraComm); 752 753 std::vector<StdSize> start(1, startIndex - localNbWeight); 754 std::vector<StdSize> count(1, localNbWeight); 755 756 WriteNetCdf netCdfWriter(filename, client->intraComm); 757 758 // netCdfWriter = CONetCDF4(filename, false, false, true, client->intraComm, false); 759 760 // Define some dimensions 761 netCdfWriter.addDimensionWrite("n_src", n_src); 762 netCdfWriter.addDimensionWrite("n_dst", n_dst); 763 netCdfWriter.addDimensionWrite("n_weight", globalNbWeight); 764 765 std::vector<StdString> dims(1,"n_weight"); 766 767 // Add some variables 768 netCdfWriter.addVariableWrite("src_idx", NC_INT, dims); 769 netCdfWriter.addVariableWrite("dst_idx", NC_INT, dims); 770 netCdfWriter.addVariableWrite("weight", NC_DOUBLE, dims); 771 772 // // Write variables 773 netCdfWriter.writeDataIndex(src_idx, "src_idx", true, 0, &start, &count); 774 netCdfWriter.writeDataIndex(dst_idx, "dst_idx", true, 0, &start, &count); 775 netCdfWriter.writeDataIndex(weights, "weight", true, 0, &start, &count); 776 777 netCdfWriter.closeFile(); 778 } 639 779 640 780 /*!
Note: See TracChangeset
for help on using the changeset viewer.