Skip to content
Snippets Groups Projects
Commit 40b96f8a authored by WAUTELET Philippe's avatar WAUTELET Philippe
Browse files

lfi2cdf: new option: '-c or --compress' to compress data

parent 2860ba37
No related branches found
No related tags found
No related merge requests found
subroutine LFI2CDFMAIN(hinfile,iiflen,ooutname,houtfile,ioflen,hvarlist,ivlen,olfi2cdf,olfilist,ohdf5,omerge,nb_levels,&
oreduceprecision)
oreduceprecision,ocompress,compress_level)
USE mode_util
IMPLICIT NONE
INTEGER :: iiflen, ioflen, ivlen
......@@ -7,8 +7,9 @@ subroutine LFI2CDFMAIN(hinfile,iiflen,ooutname,houtfile,ioflen,hvarlist,ivlen,o
CHARACTER(LEN=iiflen) :: hinfile
CHARACTER(LEN=ioflen) :: houtfile
CHARACTER(LEN=ivlen) :: hvarlist
LOGICAL :: ooutname, olfi2cdf, olfilist, ohdf5, omerge, oreduceprecision
LOGICAL :: ooutname, olfi2cdf, olfilist, ohdf5, omerge, oreduceprecision, ocompress
INTEGER :: compress_level
INTEGER :: ibuflen
INTEGER :: ilu
INTEGER :: inaf, ji
......@@ -44,7 +45,7 @@ subroutine LFI2CDFMAIN(hinfile,iiflen,ooutname,houtfile,ioflen,hvarlist,ivlen,o
!Standard treatment (one LFI file only)
IF (.not.omerge) THEN
CALL parse_lfi(ilu,hvarlist,inaf,tzreclist,ibuflen)
CALL def_ncdf(tzreclist,inaf,oreduceprecision,icdf_id,omerge)
CALL def_ncdf(tzreclist,inaf,oreduceprecision,icdf_id,omerge,ocompress,compress_level)
CALL fill_ncdf(ilu,icdf_id,tzreclist,inaf,ibuflen)
ELSE
!Treat several LFI files and merge into 1 NC file
......@@ -57,7 +58,7 @@ subroutine LFI2CDFMAIN(hinfile,iiflen,ooutname,houtfile,ioflen,hvarlist,ivlen,o
!Read 1st LFI file
CALL parse_lfi(ilu,hvarlist,inaf,tzreclist,ibuflen,current_level)
CALL def_ncdf(tzreclist,inaf,oreduceprecision,icdf_id,omerge)
CALL def_ncdf(tzreclist,inaf,oreduceprecision,icdf_id,omerge,ocompress,compress_level)
DO current_level = first_level,last_level
print *,'Treating level ',current_level
......
......@@ -274,12 +274,14 @@ CONTAINS
END IF
END SUBROUTINE HANDLE_ERR
SUBROUTINE def_ncdf(tpreclist,knaf,oreduceprecision,kcdf_id,omerge)
SUBROUTINE def_ncdf(tpreclist,knaf,oreduceprecision,kcdf_id,omerge,ocompress,compress_level)
TYPE(workfield),DIMENSION(:),INTENT(INOUT) :: tpreclist
INTEGER, INTENT(IN) :: knaf
LOGICAL, INTENT(IN) :: oreduceprecision
INTEGER, INTENT(OUT):: kcdf_id
LOGICAL, INTENT(IN) :: omerge
LOGICAL, INTENT(IN) :: ocompress
INTEGER, INTENT(IN) :: compress_level
INTEGER :: status
INTEGER :: ji
......@@ -383,6 +385,12 @@ CONTAINS
END SELECT
! Compress data (costly operation for the CPU)
IF (ocompress .AND. invdims>0) THEN
status = NF90_DEF_VAR_DEFLATE(kcdf_id,tpreclist(ji)%id,1,1,compress_level)
IF (status /= NF90_NOERR) CALL HANDLE_ERR(status,__LINE__)
END IF
! GRID attribute definition
status = NF90_PUT_ATT(kcdf_id,tpreclist(ji)%id,'GRID',tpreclist(ji)%grid)
IF (status /= NF90_NOERR) CALL HANDLE_ERR(status,__LINE__)
......
......@@ -5,7 +5,7 @@
#define BUFSIZE 4096
extern lfi2cdfmain_(char*, int*, int *, char*, int*, char*, int*, int*, int*, int*, int*, int*, int*);
extern lfi2cdfmain_(char*, int*, int *, char*, int*, char*, int*, int*, int*, int*, int*, int*, int*, int*, int*);
char *cleancomma(char *varlist)
{
......@@ -33,6 +33,7 @@ int main(int argc, char **argv)
int merge_flag, nb_levels;
int reduceprecision_flag;
int outname_flag;
int compress_flag, compress_level;
char *cmd, *infile;
int c;
char buff[BUFSIZE];
......@@ -53,6 +54,7 @@ int main(int argc, char **argv)
list_flag = 0;
hdf5_flag = 0;
reduceprecision_flag = 0;
compress_flag = 0;
p = buff;
*p = '\0';
......@@ -65,6 +67,7 @@ int main(int argc, char **argv)
static struct option long_options[] = {
{"cdf4", no_argument, 0, '4' },
{"compress", required_argument, 0, 'c' },
{"list", no_argument, 0, 'l' },
{"merge", required_argument, 0, 'm' },
{"output", required_argument, 0, 'o' },
......@@ -73,7 +76,7 @@ int main(int argc, char **argv)
{0, 0, 0, 0 }
};
c = getopt_long(argc, argv, "4lm:o:rv:",
c = getopt_long(argc, argv, "4c:lm:o:rv:",
long_options, &option_index);
if (c == -1)
break;
......@@ -85,6 +88,14 @@ int main(int argc, char **argv)
printf(" with arg %s", optarg);
printf("\n");
break;
case 'c':
compress_flag = 1;
compress_level = atoi(optarg);
if(compress_level<1 || compress_level>9) {
printf("Error: compression level should in the 1 to 9 interval\n");
exit(EXIT_FAILURE);
}
break;
case '4':
hdf5_flag = 1;
break;
......@@ -124,7 +135,7 @@ int main(int argc, char **argv)
}
if (optind == argc) {
printf("usage : lfi2cdf [--cdf4 -4] [-l] [-v --var var1[,...]] [-r --reduce-precision] [-m --merge number_of_z_levels] [-o --output output-file.nc] input-file.lfi\n");
printf("usage : lfi2cdf [--cdf4 -4] [-l] [-v --var var1[,...]] [-r --reduce-precision] [-m --merge number_of_z_levels] [-o --output output-file.nc] [-c --compress compression_level] input-file.lfi\n");
printf(" cdf2lfi [-o --output output-file.lfi] input-file.nc\n");
exit(EXIT_FAILURE);
}
......@@ -156,13 +167,19 @@ int main(int argc, char **argv)
olen = strlen(outfile);
}
/* Compression flag only supported if using netCDF4 */
if (hdf5_flag==0 && compress_flag==1) {
compress_flag = 0;
printf("Warning: compression is forced to disable (only supported from netCDF4).\n");
}
/*
printf("cmd=%s; inputfile=%s(%d); outputfile=%s(%d); varlistclean=%s with size : %d\n", cmd,
infile, ilen, outfile, olen, varlist, varlistlen);
*/
lfi2cdfmain_(infile, &ilen, &outname_flag, outfile, &olen, varlist, &varlistlen, &l2c_flag, &list_flag, &hdf5_flag, &merge_flag,
&nb_levels, &reduceprecision_flag);
&nb_levels, &reduceprecision_flag, &compress_flag, &compress_level);
exit(EXIT_SUCCESS);
}
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment