Changeset 1014
- Timestamp:
- 01/04/17 17:09:50 (6 years ago)
- Location:
- XIOS/trunk
- Files:
-
- 6 edited
Legend:
- Unmodified
- Added
- Removed
-
XIOS/trunk/inputs/REMAP/iodef.xml
r982 r1014 109 109 <domain_group id="domain_dst"> 110 110 <domain id="dst_domain"> 111 <interpolate_domain file="weight.nc" />111 <interpolate_domain write_weight="true" /> 112 112 </domain> 113 113 <domain id="dst_domain_regular_pole" ni_glo="90" nj_glo="45" type="rectilinear"> 114 114 <generate_rectilinear_domain id="domain_regular_pole"/> 115 <interpolate_domain />115 <interpolate_domain write_weight="false"/> 116 116 <zoom_domain ibegin="0" ni="45" jbegin="0" nj="45" /> 117 117 </domain> -
XIOS/trunk/src/config/interpolate_domain_attribute.conf
r1004 r1014 1 1 /* GLOBAL */ 2 DECLARE_ATTRIBUTE(StdString, file)3 2 DECLARE_ATTRIBUTE(int, order) 4 3 DECLARE_ATTRIBUTE(bool, renormalize) -
XIOS/trunk/src/interface/c_attr/icinterpolate_domain_attr.cpp
r1005 r1014 17 17 { 18 18 typedef xios::CInterpolateDomain* interpolate_domain_Ptr; 19 20 void cxios_set_interpolate_domain_file(interpolate_domain_Ptr interpolate_domain_hdl, const char * file, int file_size)21 {22 std::string file_str;23 if (!cstr2string(file, file_size, file_str)) return;24 CTimer::get("XIOS").resume();25 interpolate_domain_hdl->file.setValue(file_str);26 CTimer::get("XIOS").suspend();27 }28 29 void cxios_get_interpolate_domain_file(interpolate_domain_Ptr interpolate_domain_hdl, char * file, int file_size)30 {31 CTimer::get("XIOS").resume();32 if (!string_copy(interpolate_domain_hdl->file.getInheritedValue(), file, file_size))33 ERROR("void cxios_get_interpolate_domain_file(interpolate_domain_Ptr interpolate_domain_hdl, char * file, int file_size)", << "Input string is too short");34 CTimer::get("XIOS").suspend();35 }36 37 bool cxios_is_defined_interpolate_domain_file(interpolate_domain_Ptr interpolate_domain_hdl)38 {39 CTimer::get("XIOS").resume();40 bool isDefined = interpolate_domain_hdl->file.hasInheritedValue();41 CTimer::get("XIOS").suspend();42 return isDefined;43 }44 45 19 46 20 void cxios_set_interpolate_domain_mode(interpolate_domain_Ptr interpolate_domain_hdl, const char * mode, int mode_size) -
XIOS/trunk/src/node/interpolate_domain.cpp
r1004 r1014 55 55 { 56 56 case mode_attr::read: 57 if (this-> file.isEmpty())57 if (this->weight_filename.isEmpty()) 58 58 { 59 59 if (!this->write_weight) … … 64 64 else 65 65 { 66 weightFile = this-> file;66 weightFile = this->weight_filename; 67 67 ifstream f(weightFile.c_str()); 68 68 if (!f.good()) … … 75 75 break; 76 76 case mode_attr::read_or_compute: 77 if (!this-> file.isEmpty() && !this->write_weight)77 if (!this->weight_filename.isEmpty() && !this->write_weight) 78 78 { 79 weightFile = this-> file;79 weightFile = this->weight_filename; 80 80 ifstream f(weightFile.c_str()); 81 81 if (!f.good()) -
XIOS/trunk/src/transformation/domain_algorithm_interpolate.cpp
r1004 r1014 681 681 const StdSize size) 682 682 { 683 CONetCDF4::addDimension(name, size);683 return CONetCDF4::addDimension(name, size); 684 684 } 685 685 … … 687 687 const std::vector<StdString>& dim) 688 688 { 689 CONetCDF4::addVariable(name, type, dim); 689 return CONetCDF4::addVariable(name, type, dim); 690 } 691 692 void CDomainAlgorithmInterpolate::WriteNetCdf::endDefinition() 693 { 694 CONetCDF4::definition_end(); 690 695 } 691 696 … … 751 756 MPI_Scan(&localNbWeight, &startIndex, 1, MPI_LONG, MPI_SUM, client->intraComm); 752 757 758 if (0 == globalNbWeight) 759 { 760 info << "There is no interpolation weights calculated between " 761 << "domain source: " << domainSrc_->getDomainOutputName() 762 << " and domain destination: " << domainDest_->getDomainOutputName() 763 << std::endl; 764 return; 765 } 766 753 767 std::vector<StdSize> start(1, startIndex - localNbWeight); 754 768 std::vector<StdSize> count(1, localNbWeight); 755 756 WriteNetCdf netCdfWriter(filename, client->intraComm); 757 758 // netCdfWriter = CONetCDF4(filename, false, false, true, client->intraComm, false); 769 770 WriteNetCdf netCdfWriter(filename, client->intraComm); 759 771 760 772 // Define some dimensions … … 770 782 netCdfWriter.addVariableWrite("weight", NC_DOUBLE, dims); 771 783 784 // End of definition 785 netCdfWriter.endDefinition(); 786 772 787 // // 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); 788 if (0 != localNbWeight) 789 { 790 netCdfWriter.writeDataIndex(src_idx, "src_idx", false, 0, &start, &count); 791 netCdfWriter.writeDataIndex(dst_idx, "dst_idx", false, 0, &start, &count); 792 netCdfWriter.writeDataIndex(weights, "weight", false, 0, &start, &count); 793 } 776 794 777 795 netCdfWriter.closeFile(); -
XIOS/trunk/src/transformation/domain_algorithm_interpolate.hpp
r1004 r1014 59 59 int addVariableWrite(const StdString& name, nc_type type, 60 60 const std::vector<StdString>& dim); 61 void endDefinition(); 61 62 void writeDataIndex(const CArray<int,1>& data, const StdString& name, 62 63 bool collective, StdSize record,
Note: See TracChangeset
for help on using the changeset viewer.