Skip to content
GitLab
Explore
Sign in
Primary navigation
Search or go to…
Project
M
Méso-NH code
Manage
Activity
Members
Labels
Plan
Issues
Issue boards
Milestones
Wiki
Code
Merge requests
Repository
Branches
Commits
Tags
Repository graph
Compare revisions
Snippets
Build
Pipelines
Jobs
Pipeline schedules
Artifacts
Deploy
Releases
Package Registry
Model registry
Operate
Environments
Terraform modules
Monitor
Incidents
Analyze
Value stream analytics
Contributor analytics
CI/CD analytics
Repository analytics
Model experiments
Help
Help
Support
GitLab documentation
Compare GitLab plans
Community forum
Contribute to GitLab
Provide feedback
Terms and privacy
Keyboard shortcuts
?
Snippets
Groups
Projects
Show more breadcrumbs
RODIER Quentin
Méso-NH code
Commits
87c049ab
Commit
87c049ab
authored
9 years ago
by
WAUTELET Philippe
Browse files
Options
Downloads
Patches
Plain Diff
lfi2cdf: -s/--split option (1 variable per output file) is now implemented
parent
a8ba3340
No related branches found
No related tags found
No related merge requests found
Changes
2
Hide whitespace changes
Inline
Side-by-side
Showing
2 changed files
tools/lfi2cdf/src/lfi2cdf.f90
+36
-10
36 additions, 10 deletions
tools/lfi2cdf/src/lfi2cdf.f90
tools/lfi2cdf/src/mode_util.f90
+130
-40
130 additions, 40 deletions
tools/lfi2cdf/src/mode_util.f90
with
166 additions
and
50 deletions
tools/lfi2cdf/src/lfi2cdf.f90
+
36
−
10
View file @
87c049ab
...
@@ -18,19 +18,41 @@ subroutine LFI2CDFMAIN(hinfile,iiflen,ooutname,houtfile,ioflen,hvarlist,ivlen,o
...
@@ -18,19 +18,41 @@ subroutine LFI2CDFMAIN(hinfile,iiflen,ooutname,houtfile,ioflen,hvarlist,ivlen,o
INTEGER
::
nbvar_calc
! number of variables to be computed from others
INTEGER
::
nbvar_calc
! number of variables to be computed from others
INTEGER
::
nbvar_tbw
! number of variables to be written
INTEGER
::
nbvar_tbw
! number of variables to be written
INTEGER
::
nbvar
! number of defined variables
INTEGER
::
nbvar
! number of defined variables
INTEGER
::
icdf_id
INTEGER
::
first_level
,
current_level
,
last_level
INTEGER
::
first_level
,
current_level
,
last_level
INTEGER
(
KIND
=
LFI_INT
)
::
iresp
,
iverb
,
inap
INTEGER
(
KIND
=
LFI_INT
)
::
iresp
,
iverb
,
inap
CHARACTER
(
LEN
=
3
)
::
suffix
CHARACTER
(
LEN
=
3
)
::
suffix
CHARACTER
(
LEN
=
iiflen
)
::
filename
CHARACTER
(
LEN
=
iiflen
)
::
filename
TYPE
(
cdf_files
)
::
cdffiles
TYPE
(
workfield
),
DIMENSION
(:),
POINTER
::
tzreclist
TYPE
(
workfield
),
DIMENSION
(:),
POINTER
::
tzreclist
cdffiles
%
nbfiles
=
0
cdffiles
%
opened
=
.FALSE.
!Remove level in the filename if merging LFI splitted files
!Remove level in the filename if merging LFI splitted files
if
(
omerge
.AND.
.NOT.
ooutname
)
then
if
(
.NOT.
ooutname
)
then
if
(
omerge
.AND.
.NOT.
osplit
)
then
houtfile
=
houtfile
(
1
:
len
(
houtfile
)
-9
)//
houtfile
(
len
(
houtfile
)
-3
:)
houtfile
=
houtfile
(
1
:
len
(
houtfile
)
-9
)//
houtfile
(
len
(
houtfile
)
-3
:)
end
if
if
(
.NOT.
omerge
.AND.
osplit
)
then
if
(
ohdf5
)
then
ji
=
4
else
ji
=
3
end
if
houtfile
=
houtfile
(
1
:
len
(
houtfile
)
-
ji
)
end
if
if
(
omerge
.AND.
osplit
)
then
if
(
ohdf5
)
then
ji
=
9
else
ji
=
8
end
if
houtfile
=
houtfile
(
1
:
len
(
houtfile
)
-
ji
)
end
if
end
if
end
if
CALL
OPEN_FILES
(
hinfile
,
houtfile
,
olfi2cdf
,
olfilist
,
ohdf5
,
i
cdf
_id
,
ilu
,
nbvar_lfi
)
CALL
OPEN_FILES
(
hinfile
,
houtfile
,
olfi2cdf
,
olfilist
,
ohdf5
,
cdf
files
,
ilu
,
nbvar_lfi
,
osplit
)
IF
(
olfilist
)
return
IF
(
olfilist
)
return
IF
(
olfi2cdf
)
THEN
IF
(
olfi2cdf
)
THEN
...
@@ -59,8 +81,10 @@ subroutine LFI2CDFMAIN(hinfile,iiflen,ooutname,houtfile,ioflen,hvarlist,ivlen,o
...
@@ -59,8 +81,10 @@ subroutine LFI2CDFMAIN(hinfile,iiflen,ooutname,houtfile,ioflen,hvarlist,ivlen,o
!Standard treatment (one LFI file only)
!Standard treatment (one LFI file only)
IF
(
.not.
omerge
)
THEN
IF
(
.not.
omerge
)
THEN
CALL
parse_lfi
(
ilu
,
hvarlist
,
nbvar_lfi
,
nbvar_tbr
,
nbvar_calc
,
nbvar_tbw
,
tzreclist
,
ibuflen
)
CALL
parse_lfi
(
ilu
,
hvarlist
,
nbvar_lfi
,
nbvar_tbr
,
nbvar_calc
,
nbvar_tbw
,
tzreclist
,
ibuflen
)
CALL
def_ncdf
(
tzreclist
,
nbvar
,
oreduceprecision
,
icdf_id
,
omerge
,
ocompress
,
compress_level
)
IF
(
osplit
)
call
open_split_ncfiles
(
houtfile
,
nbvar
,
tzreclist
,
cdffiles
,
ohdf5
)
CALL
fill_ncdf
(
ilu
,
icdf_id
,
tzreclist
,
nbvar
,
ibuflen
)
CALL
def_ncdf
(
tzreclist
,
nbvar
,
oreduceprecision
,
cdffiles
,
omerge
,
ocompress
,
compress_level
)
CALL
fill_ncdf
(
ilu
,
tzreclist
,
nbvar
,
ibuflen
,
cdffiles
)
ELSE
ELSE
!Treat several LFI files and merge into 1 NC file
!Treat several LFI files and merge into 1 NC file
iverb
=
0
!Verbosity level for LFI
iverb
=
0
!Verbosity level for LFI
...
@@ -72,7 +96,9 @@ subroutine LFI2CDFMAIN(hinfile,iiflen,ooutname,houtfile,ioflen,hvarlist,ivlen,o
...
@@ -72,7 +96,9 @@ subroutine LFI2CDFMAIN(hinfile,iiflen,ooutname,houtfile,ioflen,hvarlist,ivlen,o
!Read 1st LFI file
!Read 1st LFI file
CALL
parse_lfi
(
ilu
,
hvarlist
,
nbvar_lfi
,
nbvar_tbr
,
nbvar_calc
,
nbvar_tbw
,
tzreclist
,
ibuflen
,
current_level
)
CALL
parse_lfi
(
ilu
,
hvarlist
,
nbvar_lfi
,
nbvar_tbr
,
nbvar_calc
,
nbvar_tbw
,
tzreclist
,
ibuflen
,
current_level
)
CALL
def_ncdf
(
tzreclist
,
nbvar
,
oreduceprecision
,
icdf_id
,
omerge
,
ocompress
,
compress_level
)
IF
(
osplit
)
call
open_split_ncfiles
(
houtfile
,
nbvar
,
tzreclist
,
cdffiles
,
ohdf5
)
!Define NC variables
CALL
def_ncdf
(
tzreclist
,
nbvar
,
oreduceprecision
,
cdffiles
,
omerge
,
ocompress
,
compress_level
)
DO
current_level
=
first_level
,
last_level
DO
current_level
=
first_level
,
last_level
print
*
,
'Treating level '
,
current_level
print
*
,
'Treating level '
,
current_level
...
@@ -82,18 +108,18 @@ subroutine LFI2CDFMAIN(hinfile,iiflen,ooutname,houtfile,ioflen,hvarlist,ivlen,o
...
@@ -82,18 +108,18 @@ subroutine LFI2CDFMAIN(hinfile,iiflen,ooutname,houtfile,ioflen,hvarlist,ivlen,o
CALL
LFIOUV
(
iresp
,
ilu
,
ltrue
,
filename
,
'OLD'
,
lfalse
,
lfalse
,
iverb
,
inap
,
nbvar_lfi
)
CALL
LFIOUV
(
iresp
,
ilu
,
ltrue
,
filename
,
'OLD'
,
lfalse
,
lfalse
,
iverb
,
inap
,
nbvar_lfi
)
CALL
read_data_lfi
(
ilu
,
hvarlist
,
nbvar
,
tzreclist
,
ibuflen
,
current_level
)
CALL
read_data_lfi
(
ilu
,
hvarlist
,
nbvar
,
tzreclist
,
ibuflen
,
current_level
)
END
IF
END
IF
CALL
fill_ncdf
(
ilu
,
icdf_id
,
tzreclist
,
nbvar
,
ibuflen
,
current_level
)
CALL
fill_ncdf
(
ilu
,
tzreclist
,
nbvar
,
ibuflen
,
cdffiles
,
current_level
)
IF
(
current_level
/
=
last_level
)
CALL
LFIFER
(
iresp
,
ilu
,
'KEEP'
)
IF
(
current_level
/
=
last_level
)
CALL
LFIFER
(
iresp
,
ilu
,
'KEEP'
)
END
DO
END
DO
END
IF
END
IF
ELSE
ELSE
! Conversion NetCDF -> LFI
! Conversion NetCDF -> LFI
CALL
parse_cdf
(
i
cdf_id
,
tzreclist
,
ibuflen
)
CALL
parse_cdf
(
cdffiles
%
cdf_id
(
1
)
,
tzreclist
,
ibuflen
)
CALL
build_lfi
(
i
cdf_id
,
ilu
,
tzreclist
,
ibuflen
)
CALL
build_lfi
(
cdffiles
%
cdf_id
(
1
)
,
ilu
,
tzreclist
,
ibuflen
)
END
IF
END
IF
CALL
CLOSE_FILES
(
ilu
,
i
cdf
_id
)
CALL
CLOSE_FILES
(
ilu
,
cdf
files
,
osplit
)
end
subroutine
LFI2CDFMAIN
end
subroutine
LFI2CDFMAIN
This diff is collapsed.
Click to expand it.
tools/lfi2cdf/src/mode_util.f90
+
130
−
40
View file @
87c049ab
...
@@ -7,6 +7,15 @@ MODULE mode_util
...
@@ -7,6 +7,15 @@ MODULE mode_util
IMPLICIT
NONE
IMPLICIT
NONE
INTEGER
,
PARAMETER
::
MAXRAW
=
10
INTEGER
,
PARAMETER
::
MAXRAW
=
10
INTEGER
,
PARAMETER
::
MAXLEN
=
512
TYPE
cdf_files
INTEGER
::
nbfiles
LOGICAL
::
opened
INTEGER
,
DIMENSION
(:),
ALLOCATABLE
::
cdf_id
!ID of the netCDF file
INTEGER
,
DIMENSION
(:),
ALLOCATABLE
::
var_id
!position of the variable in the workfield structure
END
TYPE
cdf_files
TYPE
workfield
TYPE
workfield
CHARACTER
(
LEN
=
FM_FIELD_SIZE
)
::
name
! nom du champ
CHARACTER
(
LEN
=
FM_FIELD_SIZE
)
::
name
! nom du champ
...
@@ -388,17 +397,18 @@ END DO
...
@@ -388,17 +397,18 @@ END DO
END
IF
END
IF
END
SUBROUTINE
HANDLE_ERR
END
SUBROUTINE
HANDLE_ERR
SUBROUTINE
def_ncdf
(
tpreclist
,
nbvar
,
oreduceprecision
,
k
cdf
_id
,
omerge
,
ocompress
,
compress_level
)
SUBROUTINE
def_ncdf
(
tpreclist
,
nbvar
,
oreduceprecision
,
cdf
files
,
omerge
,
ocompress
,
compress_level
)
TYPE
(
workfield
),
DIMENSION
(:),
INTENT
(
INOUT
)
::
tpreclist
TYPE
(
workfield
),
DIMENSION
(:),
INTENT
(
INOUT
)
::
tpreclist
INTEGER
,
INTENT
(
IN
)
::
nbvar
INTEGER
,
INTENT
(
IN
)
::
nbvar
LOGICAL
,
INTENT
(
IN
)
::
oreduceprecision
LOGICAL
,
INTENT
(
IN
)
::
oreduceprecision
INTEGER
,
INTENT
(
OUT
)
::
k
cdf
_id
TYPE
(
cdf_files
),
INTENT
(
IN
)
::
cdf
files
LOGICAL
,
INTENT
(
IN
)
::
omerge
LOGICAL
,
INTENT
(
IN
)
::
omerge
LOGICAL
,
INTENT
(
IN
)
::
ocompress
LOGICAL
,
INTENT
(
IN
)
::
ocompress
INTEGER
,
INTENT
(
IN
)
::
compress_level
INTEGER
,
INTENT
(
IN
)
::
compress_level
INTEGER
::
status
INTEGER
::
status
INTEGER
::
ji
INTEGER
::
idx
,
ji
,
nbfiles
INTEGER
::
kcdf_id
TYPE
(
dimCDF
),
POINTER
::
tzdim
TYPE
(
dimCDF
),
POINTER
::
tzdim
INTEGER
::
invdims
INTEGER
::
invdims
INTEGER
::
type_float
INTEGER
::
type_float
...
@@ -406,29 +416,36 @@ END DO
...
@@ -406,29 +416,36 @@ END DO
CHARACTER
(
LEN
=
20
)
::
ycdfvar
CHARACTER
(
LEN
=
20
)
::
ycdfvar
nbfiles
=
cdffiles
%
nbfiles
IF
(
oreduceprecision
)
THEN
IF
(
oreduceprecision
)
THEN
type_float
=
NF90_REAL
type_float
=
NF90_REAL
ELSE
ELSE
type_float
=
NF90_DOUBLE
type_float
=
NF90_DOUBLE
END
IF
END
IF
DO
ji
=
1
,
nbfiles
kcdf_id
=
cdffiles
%
cdf_id
(
ji
)
! global attributes
! global attributes
status
=
NF90_PUT_ATT
(
kcdf_id
,
NF90_GLOBAL
,
'Title'
,
VERSION_ID
)
status
=
NF90_PUT_ATT
(
kcdf_id
,
NF90_GLOBAL
,
'Title'
,
VERSION_ID
)
IF
(
status
/
=
NF90_NOERR
)
CALL
HANDLE_ERR
(
status
,
__
LINE__
)
IF
(
status
/
=
NF90_NOERR
)
CALL
HANDLE_ERR
(
status
,
__
LINE__
)
! define DIMENSIONS
! define DIMENSIONS
tzdim
=>
first_DimCDF
()
tzdim
=>
first_DimCDF
()
DO
WHILE
(
ASSOCIATED
(
tzdim
))
DO
WHILE
(
ASSOCIATED
(
tzdim
))
IF
(
tzdim
%
create
)
THEN
IF
(
tzdim
%
create
)
THEN
status
=
NF90_DEF_DIM
(
kcdf_id
,
tzdim
%
name
,
tzdim
%
len
,
tzdim
%
id
)
status
=
NF90_DEF_DIM
(
kcdf_id
,
tzdim
%
name
,
tzdim
%
len
,
tzdim
%
id
)
IF
(
status
/
=
NF90_NOERR
)
CALL
HANDLE_ERR
(
status
,
__
LINE__
)
IF
(
status
/
=
NF90_NOERR
)
CALL
HANDLE_ERR
(
status
,
__
LINE__
)
END
IF
END
IF
tzdim
=>
tzdim
%
next
tzdim
=>
tzdim
%
next
END
DO
END
DO
END
DO
PRINT
*
,
'------------- NetCDF DEFINITION ---------------'
PRINT
*
,
'------------- NetCDF DEFINITION ---------------'
! define VARIABLES and ATTRIBUTES
! define VARIABLES and ATTRIBUTES
idx
=
1
DO
ji
=
1
,
nbvar
DO
ji
=
1
,
nbvar
IF
(
.NOT.
tpreclist
(
ji
)
%
tbw
)
CYCLE
IF
(
.NOT.
tpreclist
(
ji
)
%
tbw
)
CYCLE
...
@@ -470,6 +487,8 @@ END DO
...
@@ -470,6 +487,8 @@ END DO
!! ni les '.' remplaces par '--'
!! ni les '.' remplaces par '--'
ycdfvar
=
str_replace
(
ycdfvar
,
'.'
,
'--'
)
ycdfvar
=
str_replace
(
ycdfvar
,
'.'
,
'--'
)
if
(
nbfiles
>
1
)
kcdf_id
=
cdffiles
%
cdf_id
(
idx
)
SELECT
CASE
(
tpreclist
(
ji
)
%
TYPE
)
SELECT
CASE
(
tpreclist
(
ji
)
%
TYPE
)
CASE
(
TEXT
)
CASE
(
TEXT
)
! PRINT *,'TEXT : ',tpreclist(ji)%name
! PRINT *,'TEXT : ',tpreclist(ji)%name
...
@@ -514,29 +533,30 @@ END DO
...
@@ -514,29 +533,30 @@ END DO
status
=
NF90_PUT_ATT
(
kcdf_id
,
tpreclist
(
ji
)
%
id
,
'COMMENT'
,
trim
(
tpreclist
(
ji
)
%
comment
))
status
=
NF90_PUT_ATT
(
kcdf_id
,
tpreclist
(
ji
)
%
id
,
'COMMENT'
,
trim
(
tpreclist
(
ji
)
%
comment
))
IF
(
status
/
=
NF90_NOERR
)
CALL
HANDLE_ERR
(
status
,
__
LINE__
)
IF
(
status
/
=
NF90_NOERR
)
CALL
HANDLE_ERR
(
status
,
__
LINE__
)
idx
=
idx
+
1
END
DO
END
DO
status
=
NF90_ENDDEF
(
kcdf_id
)
DO
ji
=
1
,
nbfiles
IF
(
status
/
=
NF90_NOERR
)
CALL
HANDLE_ERR
(
status
,
__
LINE__
)
kcdf_id
=
cdffiles
%
cdf_id
(
ji
)
status
=
NF90_ENDDEF
(
kcdf_id
)
IF
(
status
/
=
NF90_NOERR
)
CALL
HANDLE_ERR
(
status
,
__
LINE__
)
END
DO
END
SUBROUTINE
def_ncdf
END
SUBROUTINE
def_ncdf
SUBROUTINE
fill_ncdf
(
klu
,
kcdf_id
,
tpreclist
,
knaf
,
kbuflen
,
current_level
)
SUBROUTINE
fill_ncdf
(
klu
,
tpreclist
,
knaf
,
kbuflen
,
cdffiles
,
current_level
)
INTEGER
,
INTENT
(
IN
)::
klu
INTEGER
,
INTENT
(
IN
)::
klu
INTEGER
,
INTENT
(
IN
)::
kcdf_id
TYPE
(
workfield
),
DIMENSION
(:),
INTENT
(
IN
)::
tpreclist
TYPE
(
workfield
),
DIMENSION
(:),
INTENT
(
IN
)::
tpreclist
INTEGER
,
INTENT
(
IN
)::
knaf
INTEGER
,
INTENT
(
IN
)::
knaf
INTEGER
,
INTENT
(
IN
)::
kbuflen
INTEGER
,
INTENT
(
IN
)::
kbuflen
TYPE
(
cdf_files
),
INTENT
(
IN
)::
cdffiles
INTEGER
,
INTENT
(
IN
),
OPTIONAL
::
current_level
INTEGER
,
INTENT
(
IN
),
OPTIONAL
::
current_level
#ifdef LOWMEM
#ifdef LOWMEM
INTEGER
(
KIND
=
8
),
DIMENSION
(:),
ALLOCATABLE
::
iwork
INTEGER
(
KIND
=
8
),
DIMENSION
(:),
ALLOCATABLE
::
iwork
#endif
#endif
INTEGER
::
ji
,
jj
INTEGER
::
idx
,
ji
,
jj
INTEGER
,
DIMENSION
(:),
ALLOCATABLE
::
itab
INTEGER
::
kcdf_id
REAL
(
KIND
=
8
),
DIMENSION
(:),
ALLOCATABLE
::
xtab
CHARACTER
,
DIMENSION
(:),
ALLOCATABLE
::
ytab
INTEGER
::
status
INTEGER
::
status
INTEGER
::
extent
,
ndims
INTEGER
::
extent
,
ndims
INTEGER
::
ich
INTEGER
::
ich
...
@@ -544,7 +564,12 @@ END DO
...
@@ -544,7 +564,12 @@ END DO
INTEGER
::
level
INTEGER
::
level
INTEGER
(
KIND
=
LFI_INT
)
::
iresp
,
ilu
,
ileng
,
ipos
INTEGER
(
KIND
=
LFI_INT
)
::
iresp
,
ilu
,
ileng
,
ipos
CHARACTER
(
LEN
=
4
)
::
suffix
CHARACTER
(
LEN
=
4
)
::
suffix
INTEGER
,
DIMENSION
(:),
ALLOCATABLE
::
itab
REAL
(
KIND
=
8
),
DIMENSION
(:),
ALLOCATABLE
::
xtab
CHARACTER
,
DIMENSION
(:),
ALLOCATABLE
::
ytab
kcdf_id
=
cdffiles
%
cdf_id
(
1
)
!
!
ilu
=
klu
ilu
=
klu
!
!
...
@@ -563,9 +588,12 @@ END DO
...
@@ -563,9 +588,12 @@ END DO
ALLOCATE
(
itab
(
kbuflen
))
ALLOCATE
(
itab
(
kbuflen
))
ALLOCATE
(
xtab
(
kbuflen
))
ALLOCATE
(
xtab
(
kbuflen
))
idx
=
1
DO
ji
=
1
,
knaf
DO
ji
=
1
,
knaf
IF
(
.NOT.
tpreclist
(
ji
)
%
tbw
)
CYCLE
IF
(
.NOT.
tpreclist
(
ji
)
%
tbw
)
CYCLE
IF
(
cdffiles
%
nbfiles
>
1
)
kcdf_id
=
cdffiles
%
cdf_id
(
idx
)
#if LOWMEM
#if LOWMEM
CALL
LFINFO
(
iresp
,
ilu
,
trim
(
tpreclist
(
ji
)
%
name
)//
trim
(
suffix
),
ileng
,
ipos
)
CALL
LFINFO
(
iresp
,
ilu
,
trim
(
tpreclist
(
ji
)
%
name
)//
trim
(
suffix
),
ileng
,
ipos
)
CALL
LFILEC
(
iresp
,
ilu
,
trim
(
tpreclist
(
ji
)
%
name
)//
trim
(
suffix
),
iwork
,
ileng
)
CALL
LFILEC
(
iresp
,
ilu
,
trim
(
tpreclist
(
ji
)
%
name
)//
trim
(
suffix
),
iwork
,
ileng
)
...
@@ -701,6 +729,7 @@ print *,'lowmem: not tested!!!!!' (to be compared with no low mem version)
...
@@ -701,6 +729,7 @@ print *,'lowmem: not tested!!!!!' (to be compared with no low mem version)
END
SELECT
END
SELECT
idx
=
idx
+
1
END
DO
END
DO
DEALLOCATE
(
itab
,
xtab
)
DEALLOCATE
(
itab
,
xtab
)
#if LOWMEM
#if LOWMEM
...
@@ -900,11 +929,12 @@ print *,'lowmem: not tested!!!!!' (to be compared with no low mem version)
...
@@ -900,11 +929,12 @@ print *,'lowmem: not tested!!!!!' (to be compared with no low mem version)
END
SUBROUTINE
build_lfi
END
SUBROUTINE
build_lfi
SUBROUTINE
OPEN_FILES
(
hinfile
,
houtfile
,
olfi2cdf
,
olfilist
,
ohdf5
,
k
cdf
_id
,
klu
,
knaf
)
SUBROUTINE
OPEN_FILES
(
hinfile
,
houtfile
,
olfi2cdf
,
olfilist
,
ohdf5
,
cdf
files
,
klu
,
knaf
,
osplit
)
LOGICAL
,
INTENT
(
IN
)
::
olfi2cdf
,
olfilist
,
ohdf5
LOGICAL
,
INTENT
(
IN
)
::
olfi2cdf
,
olfilist
,
ohdf5
,
osplit
CHARACTER
(
LEN
=*
),
INTENT
(
IN
)
::
hinfile
CHARACTER
(
LEN
=*
),
INTENT
(
IN
)
::
hinfile
CHARACTER
(
LEN
=*
),
INTENT
(
IN
)
::
houtfile
CHARACTER
(
LEN
=*
),
INTENT
(
IN
)
::
houtfile
INTEGER
,
INTENT
(
OUT
)
::
kcdf_id
,
klu
,
knaf
TYPE
(
cdf_files
)
,
INTENT
(
OUT
)
::
cdffiles
INTEGER
,
INTENT
(
OUT
)
::
klu
,
knaf
INTEGER
::
extindex
INTEGER
::
extindex
INTEGER
(
KIND
=
LFI_INT
)
::
ilu
,
iresp
,
iverb
,
inap
,
inaf
INTEGER
(
KIND
=
LFI_INT
)
::
ilu
,
iresp
,
iverb
,
inap
,
inaf
...
@@ -927,18 +957,23 @@ print *,'lowmem: not tested!!!!!' (to be compared with no low mem version)
...
@@ -927,18 +957,23 @@ print *,'lowmem: not tested!!!!!' (to be compared with no low mem version)
CALL
LFILAF
(
iresp
,
ilu
,
lfalse
)
CALL
LFILAF
(
iresp
,
ilu
,
lfalse
)
CALL
LFIFER
(
iresp
,
ilu
,
'KEEP'
)
CALL
LFIFER
(
iresp
,
ilu
,
'KEEP'
)
return
return
end
IF
END
IF
IF
(
ohdf5
)
THEN
IF
(
.NOT.
osplit
)
THEN
status
=
NF90_CREATE
(
houtfile
,
IOR
(
NF90_CLOBBER
,
NF90_NETCDF4
),
kcdf_id
)
cdffiles
%
nbfiles
=
1
ELSE
allocate
(
cdffiles
%
cdf_id
(
1
))
status
=
NF90_CREATE
(
houtfile
,
IOR
(
NF90_CLOBBER
,
NF90_64BIT_OFFSET
),
kcdf_id
)
end
IF
IF
(
ohdf5
)
THEN
status
=
NF90_CREATE
(
houtfile
,
IOR
(
NF90_CLOBBER
,
NF90_NETCDF4
),
cdffiles
%
cdf_id
(
1
))
ELSE
status
=
NF90_CREATE
(
houtfile
,
IOR
(
NF90_CLOBBER
,
NF90_64BIT_OFFSET
),
cdffiles
%
cdf_id
(
1
))
END
IF
IF
(
status
/
=
NF90_NOERR
)
CALL
HANDLE_ERR
(
status
,
__
LINE__
)
IF
(
status
/
=
NF90_NOERR
)
CALL
HANDLE_ERR
(
status
,
__
LINE__
)
cdffiles
%
opened
=
.TRUE.
status
=
NF90_SET_FILL
(
k
cdf_id
,
NF90_NOFILL
,
omode
)
status
=
NF90_SET_FILL
(
cdffiles
%
cdf_id
(
1
)
,
NF90_NOFILL
,
omode
)
IF
(
status
/
=
NF90_NOERR
)
CALL
HANDLE_ERR
(
status
,
__
LINE__
)
IF
(
status
/
=
NF90_NOERR
)
CALL
HANDLE_ERR
(
status
,
__
LINE__
)
!!$ SELECT CASE(omode)
!!$ SELECT CASE(omode)
!!$ CASE (NF90_FILL)
!!$ CASE (NF90_FILL)
!!$ PRINT *,'Ancien mode : NF90_FILL'
!!$ PRINT *,'Ancien mode : NF90_FILL'
...
@@ -947,11 +982,15 @@ print *,'lowmem: not tested!!!!!' (to be compared with no low mem version)
...
@@ -947,11 +982,15 @@ print *,'lowmem: not tested!!!!!' (to be compared with no low mem version)
!!$ CASE default
!!$ CASE default
!!$ PRINT *, 'Ancien mode : inconnu'
!!$ PRINT *, 'Ancien mode : inconnu'
!!$ END SELECT
!!$ END SELECT
END
IF
! .NOT.osplit
ELSE
ELSE
! Cas NetCDF -> LFI
! Cas NetCDF -> LFI
status
=
NF90_OPEN
(
hinfile
,
NF90_NOWRITE
,
kcdf_id
)
cdffiles
%
nbfiles
=
1
allocate
(
cdffiles
%
cdf_id
(
1
))
status
=
NF90_OPEN
(
hinfile
,
NF90_NOWRITE
,
cdffiles
%
cdf_id
(
1
))
IF
(
status
/
=
NF90_NOERR
)
CALL
HANDLE_ERR
(
status
,
__
LINE__
)
IF
(
status
/
=
NF90_NOERR
)
CALL
HANDLE_ERR
(
status
,
__
LINE__
)
cdffiles
%
opened
=
.TRUE.
inap
=
100
inap
=
100
CALL
LFIOUV
(
iresp
,
ilu
,
ltrue
,
houtfile
,
'NEW'
&
CALL
LFIOUV
(
iresp
,
ilu
,
ltrue
,
houtfile
,
'NEW'
&
...
@@ -964,21 +1003,72 @@ print *,'lowmem: not tested!!!!!' (to be compared with no low mem version)
...
@@ -964,21 +1003,72 @@ print *,'lowmem: not tested!!!!!' (to be compared with no low mem version)
PRINT
*
,
'--> Fichier converti : '
,
houtfile
PRINT
*
,
'--> Fichier converti : '
,
houtfile
END
SUBROUTINE
OPEN_FILES
END
SUBROUTINE
OPEN_FILES
SUBROUTINE
OPEN_SPLIT_NCFILES
(
houtfile
,
nbvar
,
tpreclist
,
cdffiles
,
ohdf5
)
CHARACTER
(
LEN
=*
),
INTENT
(
IN
)
::
houtfile
INTEGER
,
INTENT
(
IN
)
::
nbvar
TYPE
(
workfield
),
DIMENSION
(:),
INTENT
(
IN
)
::
tpreclist
TYPE
(
cdf_files
),
INTENT
(
INOUT
)
::
cdffiles
LOGICAL
,
INTENT
(
IN
)
::
ohdf5
INTEGER
::
ji
,
idx
INTEGER
::
status
INTEGER
::
omode
CHARACTER
(
LEN
=
MAXLEN
)
::
filename
cdffiles
%
nbfiles
=
0
DO
ji
=
1
,
nbvar
IF
(
tpreclist
(
ji
)
%
tbw
)
cdffiles
%
nbfiles
=
cdffiles
%
nbfiles
+
1
END
DO
allocate
(
cdffiles
%
cdf_id
(
cdffiles
%
nbfiles
))
allocate
(
cdffiles
%
var_id
(
cdffiles
%
nbfiles
))
idx
=
1
DO
ji
=
1
,
nbvar
IF
(
.NOT.
tpreclist
(
ji
)
%
tbw
)
CYCLE
cdffiles
%
var_id
(
idx
)
=
ji
IF
(
ohdf5
)
THEN
filename
=
trim
(
houtfile
)//
'.'
//
trim
(
tpreclist
(
ji
)
%
name
)//
'.nc4'
status
=
NF90_CREATE
(
trim
(
filename
),
IOR
(
NF90_CLOBBER
,
NF90_NETCDF4
),
cdffiles
%
cdf_id
(
idx
))
ELSE
filename
=
trim
(
houtfile
)//
'.'
//
trim
(
tpreclist
(
ji
)
%
name
)//
'.nc'
status
=
NF90_CREATE
(
trim
(
filename
),
IOR
(
NF90_CLOBBER
,
NF90_64BIT_OFFSET
),
cdffiles
%
cdf_id
(
idx
))
END
IF
IF
(
status
/
=
NF90_NOERR
)
CALL
HANDLE_ERR
(
status
,
__
LINE__
)
status
=
NF90_SET_FILL
(
cdffiles
%
cdf_id
(
idx
),
NF90_NOFILL
,
omode
)
IF
(
status
/
=
NF90_NOERR
)
CALL
HANDLE_ERR
(
status
,
__
LINE__
)
idx
=
idx
+
1
END
DO
cdffiles
%
opened
=
.TRUE.
END
SUBROUTINE
OPEN_SPLIT_NCFILES
SUBROUTINE
CLOSE_FILES
(
klu
,
kcdf_id
)
SUBROUTINE
CLOSE_FILES
(
klu
,
cdffiles
,
osplit
)
INTEGER
,
INTENT
(
IN
)
::
klu
,
kcdf_id
INTEGER
,
INTENT
(
IN
)
::
klu
TYPE
(
cdf_files
),
INTENT
(
INOUT
)
::
cdffiles
LOGICAl
,
INTENT
(
IN
)
::
osplit
INTEGER
(
KIND
=
LFI_INT
)
::
iresp
,
ilu
INTEGER
(
KIND
=
LFI_INT
)
::
iresp
,
ilu
INTEGER
::
status
INTEGER
::
ji
,
status
ilu
=
klu
ilu
=
klu
! close LFI file
! close LFI file
CALL
LFIFER
(
iresp
,
ilu
,
'KEEP'
)
CALL
LFIFER
(
iresp
,
ilu
,
'KEEP'
)
! close NetCDF file
! close NetCDF files
status
=
NF90_CLOSE
(
kcdf_id
)
DO
ji
=
1
,
cdffiles
%
nbfiles
IF
(
status
/
=
NF90_NOERR
)
CALL
HANDLE_ERR
(
status
,
__
LINE__
)
status
=
NF90_CLOSE
(
cdffiles
%
cdf_id
(
ji
))
IF
(
status
/
=
NF90_NOERR
)
CALL
HANDLE_ERR
(
status
,
__
LINE__
)
END
DO
cdffiles
%
opened
=
.false.
END
SUBROUTINE
CLOSE_files
END
SUBROUTINE
CLOSE_files
END
MODULE
mode_util
END
MODULE
mode_util
This diff is collapsed.
Click to expand it.
Preview
0%
Loading
Try again
or
attach a new file
.
Cancel
You are about to add
0
people
to the discussion. Proceed with caution.
Finish editing this message first!
Save comment
Cancel
Please
register
or
sign in
to comment