source: trunk/LMDZ.COMMON/libf/evolution/io_netcdf.F90 @ 4090

Last change on this file since 4090 was 4090, checked in by jbclement, 5 days ago

PEM:
Update README and change "startpem.nc" into "startevol.nc".
JBC

File size: 41.2 KB
Line 
1MODULE io_netcdf
2!-----------------------------------------------------------------------
3! NAME
4!     io_netcdf
5!
6! DESCRIPTION
7!     Provides input/output procedures and parameters for the PEM.
8!
9! AUTHORS & DATE
10!     JB Clement, 12/2025
11!
12! NOTES
13!
14!-----------------------------------------------------------------------
15
16! DEPEDENCIES
17! -----------
18use numerics, only: dp, di, k4, minieps
19use netcdf,   only: nf90_double, nf90_noerr, nf90_strerror, nf90_write, nf90_nowrite,                &
20                    nf90_open, nf90_close, nf90_redef, nf90_enddef, nf90_inquire, nf90_max_var_dims, &
21                    nf90_inq_dimid, nf90_inquire_dimension, nf90_inq_varid, nf90_inquire_variable,   &
22                    nf90_def_var, nf90_get_var, nf90_put_var, nf90_get_att, nf90_put_att
23use stoppage, only: stop_clean
24
25! DECLARATION
26! -----------
27implicit none
28
29! PARAMETERS
30! ----------
31character(8),  parameter :: start_name     = 'start.nc'
32character(11), parameter :: start1D_name   = 'start1D.txt'
33character(10), parameter :: startfi_name   = 'startfi.nc'
34character(12), parameter :: startpem_name  = 'startevol.nc'
35character(19), parameter :: xios_day_name1 = 'Xoutdaily4pem_Y1.nc'  ! XIOS daily output file, second to last PCM year
36character(19), parameter :: xios_day_name2 = 'Xoutdaily4pem_Y2.nc'  ! XIOS daily output file, last PCM year
37character(20), parameter :: xios_yr_name1  = 'Xoutyearly4pem_Y1.nc' ! XIOS yearly output file, second to last PCM year
38character(20), parameter :: xios_yr_name2  = 'Xoutyearly4pem_Y2.nc' ! XIOS yearly output file, last PCM year
39character(11), parameter :: diagevol_name  = 'diagevol.nc'
40
41! VARIABLES
42! ---------
43logical(k4), private :: open_ncfile = .false.   ! Flag to true if a NetCDF file is already open
44logical(k4), private :: open_diagevol = .false. ! Flag to true if "diagevol.nc" is already open
45integer(di), private :: ncid_file               ! File ID
46integer(di), private :: varid                   ! Variable ID
47integer(di)          :: ncid_diagevol           ! File ID specific to "diagevol.nc"
48
49! INTERFACES
50! ----------
51interface get_var_nc
52    module procedure get_var0d_nc, get_var1d_nc, get_var2d_nc, get_var3d_nc, get_var4d_nc
53end interface get_var_nc
54
55interface put_var_nc
56    module procedure put_var0d_nc, put_var1d_nc, put_var2d_nc, put_var3d_nc, put_var4d_nc
57end interface put_var_nc
58
59interface check_valid_var_nc
60    module procedure check_valid_var0d_nc, check_valid_var1d_nc, check_valid_var2d_nc, check_valid_var3d_nc, check_valid_var4d_nc
61end interface check_valid_var_nc
62
63contains
64!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
65
66!=======================================================================
67SUBROUTINE check_nc(ierr,msg,found)
68!-----------------------------------------------------------------------
69! NAME
70!     check_nc
71!
72! DESCRIPTION
73!     NetCDF error handler.
74!
75! AUTHORS & DATE
76!     JB Clement, 12/2025
77!
78! NOTES
79!
80!-----------------------------------------------------------------------
81
82! DEPENDENCIES
83! ------------
84use display, only: print_err
85
86! ARGUMENTS
87! ---------
88integer(di),  intent(in)            :: ierr
89character(*), intent(in)            :: msg
90logical(k4),  intent(out), optional :: found
91
92! LOCAL VARIABLES
93! ---------------
94logical(k4) :: tmp_found
95
96! CODE
97! ----
98if (ierr == nf90_noerr) then
99    tmp_found = .true.
100else
101    tmp_found = .false.
102end if
103
104if (present(found)) then
105    found = tmp_found
106else
107    if (.not. tmp_found) then
108        call print_err(trim(nf90_strerror(ierr)))
109        call stop_clean(__FILE__,__LINE__,'NetCDF error when '//trim(msg),ierr)
110    end if
111end if
112
113END SUBROUTINE check_nc
114!=======================================================================
115
116!=======================================================================
117SUBROUTINE open_nc(filename,mode,itime)
118!-----------------------------------------------------------------------
119! NAME
120!     open_nc
121!
122! DESCRIPTION
123!     Open a netCDF file for reading.
124
125! AUTHORS & DATE
126!     JB Clement, 12/2025
127
128! NOTES
129!
130!-----------------------------------------------------------------------
131
132! DECLARATION
133! -----------
134implicit none
135
136! ARGUMENTS
137! ---------
138character(*), intent(in)            :: filename, mode
139integer(di),  intent(out), optional :: itime
140
141! CODE
142! ----
143! Diagevol logic
144if (adjustl(trim(filename)) == diagevol_name) then
145    if (adjustl(trim(mode)) == 'read') then
146        call stop_clean(__FILE__,__LINE__,'opening mode "read" cannot be used with '//trim(filename)//'"!',1)
147    else if (adjustl(trim(mode)) == 'write') then
148        if (.not. open_diagevol) then
149            call check_nc(nf90_open(trim(filename),nf90_write,ncid_diagevol),'opening file '//trim(filename)//' to write')
150            if (present(itime)) call get_next_itime_nc('Time',itime) ! Next time index
151            open_diagevol = .true.
152        end if
153        return
154    else
155        call stop_clean(__FILE__,__LINE__,'opening mode "'//mode//'" unknown!',1)
156    end if
157end if
158
159! Standard logic
160if (open_ncfile) then ! A file is already opened
161    call stop_clean(__FILE__,__LINE__,'a NetCDF file is already opened!',1)
162else
163    if (adjustl(trim(mode)) == 'read') then
164        call check_nc(nf90_open(trim(filename),nf90_nowrite,ncid_file),'opening file '//trim(filename)//' to read')
165    else if (adjustl(trim(mode)) == 'write') then
166        call check_nc(nf90_open(trim(filename),nf90_write,ncid_file),'opening file '//trim(filename)//' to write')
167        if (present(itime)) call get_next_itime_nc('Time',itime) ! Next time index
168    else
169        call stop_clean(__FILE__,__LINE__,'opening mode "'//mode//'" unknown!',1)
170    end if
171    open_ncfile = .true.
172end if
173
174END SUBROUTINE open_nc
175!=======================================================================
176
177!=======================================================================
178SUBROUTINE close_nc(filename)
179!-----------------------------------------------------------------------
180! NAME
181!     close_nc
182!
183! DESCRIPTION
184!     Open a netCDF file.
185
186! AUTHORS & DATE
187!     JB Clement, 12/2025
188
189! NOTES
190!
191!-----------------------------------------------------------------------
192
193! DECLARATION
194! -----------
195implicit none
196
197! ARGUMENTS
198! ---------
199character(*), intent(in) :: filename
200
201! CODE
202! ----
203if (adjustl(trim(filename)) == diagevol_name) then ! Diagevol logic
204    call check_nc(nf90_close(ncid_diagevol),'closing file '//trim(filename))
205    open_diagevol = .false.
206else ! Standard logic
207    call check_nc(nf90_close(ncid_file),'closing file '//trim(filename))
208    open_ncfile = .false.
209end if
210
211END SUBROUTINE close_nc
212!=======================================================================
213
214!=======================================================================
215SUBROUTINE get_dim_nc(dim_name,d,found)
216!-----------------------------------------------------------------------
217! NAME
218!     get_dim_nc
219!
220! DESCRIPTION
221!     Read a 0D variable from open NetCDF file.
222!
223! AUTHORS & DATE
224!     JB Clement, 12/2025
225!
226! NOTES
227!
228!-----------------------------------------------------------------------
229
230! DECLARATION
231! -----------
232implicit none
233
234! ARGUMENTS
235! ---------
236character(*), intent(in)            :: dim_name ! Variable name
237integer(di),  intent(out)           :: d        ! Output dimension
238logical(k4),  intent(out), optional :: found
239
240! LOCAL VARIABLES
241! ---------------
242integer(di) :: dimid ! Dimension ID
243
244! CODE
245! ----
246if (present(found)) then
247    call check_nc(nf90_inq_dimid(ncid_file,dim_name,dimid),'inquiring '//dim_name,found)
248    call check_nc(nf90_inquire_dimension(ncid_file,dimid,len = d),'getting '//dim_name,found)
249else
250    call check_nc(nf90_inq_dimid(ncid_file,dim_name,dimid),'inquiring '//dim_name)
251    call check_nc(nf90_inquire_dimension(ncid_file,dimid,len = d),'getting '//dim_name)
252end if
253
254END SUBROUTINE get_dim_nc
255!=======================================================================
256
257!=======================================================================
258SUBROUTINE def_var_nc(var_name,title,units,dimids)
259!-----------------------------------------------------------------------
260! NAME
261!     def_var_nc
262!
263! DESCRIPTION
264!     Define a variable into open NetCDF file.
265!
266! AUTHORS & DATE
267!     JB Clement, 01/2026
268!
269! NOTES
270!
271!-----------------------------------------------------------------------
272
273! DECLARATION
274! -----------
275implicit none
276
277! ARGUMENTS
278! ---------
279character(*),              intent(in) :: var_name, title, units ! Variable|title|units name
280integer(di), dimension(:), intent(in) :: dimids                 ! Dimensions IDs
281
282! LOCAL VARIABLES
283! ---------------
284logical(k4) :: found
285integer(di) :: ncid
286
287! CODE
288! ----
289! Diagevol logic priority over standard logic
290if (open_diagevol) then
291    ncid = ncid_diagevol
292else
293    ncid = ncid_file
294end if
295
296! Check if the variable exists
297call check_nc(nf90_inq_varid(ncid,var_name,varid),'inquiring '//var_name,found)
298if (found) return ! Variable is already defined
299
300! Enter define mode
301call check_nc(nf90_redef(ncid),'entering define mode')
302
303! Define variable
304call check_nc(nf90_def_var(ncid,var_name,NF90_DOUBLE,dimids,varid),'defining variable '//var_name)
305call check_nc(nf90_put_att(ncid,varid,'title',title),'putting title attribute for '//var_name)
306call check_nc(nf90_put_att(ncid,varid,'units',units),'putting units attribute for '//var_name)
307
308! Leave define mode
309call check_nc(nf90_enddef(ncid),'leaving define mode')
310
311END SUBROUTINE def_var_nc
312!=======================================================================
313
314!=======================================================================
315SUBROUTINE get_var0d_nc(var_name,var,found)
316!-----------------------------------------------------------------------
317! NAME
318!     get_var0d_nc
319!
320! DESCRIPTION
321!     Read a 0D variable from open NetCDF file.
322!
323! AUTHORS & DATE
324!     JB Clement, 12/2025
325!
326! NOTES
327!
328!-----------------------------------------------------------------------
329
330! DECLARATION
331! -----------
332implicit none
333
334! ARGUMENTS
335! ---------
336character(*), intent(in)            :: var_name ! Variable name
337real(dp),     intent(out)           :: var      ! Output variable
338logical(k4),  intent(out), optional :: found
339
340! CODE
341! ----
342if (present(found)) then
343    call check_nc(nf90_inq_varid(ncid_file,var_name,varid),'inquiring '//var_name,found)
344    call check_nc(nf90_get_var(ncid_file,varid,var),'getting '//var_name,found)
345else
346    call check_nc(nf90_inq_varid(ncid_file,var_name,varid),'inquiring '//var_name)
347    call check_nc(nf90_get_var(ncid_file,varid,var),'getting '//var_name)
348end if
349call check_valid_var_nc(var_name,var)
350
351END SUBROUTINE get_var0d_nc
352!=======================================================================
353
354!=======================================================================
355SUBROUTINE get_var1d_nc(var_name,var,found)
356!-----------------------------------------------------------------------
357! NAME
358!     get_var1d_nc
359!
360! DESCRIPTION
361!     Read a 1D variable from open NetCDF file.
362!
363! AUTHORS & DATE
364!     JB Clement, 12/2025
365!
366! NOTES
367!
368!-----------------------------------------------------------------------
369
370! DECLARATION
371! -----------
372implicit none
373
374! ARGUMENTS
375! ---------
376character(*),           intent(in)            :: var_name ! Variable name
377real(dp), dimension(:), intent(out)           :: var      ! Output variable
378logical(k4),            intent(out), optional :: found
379
380! CODE
381! ----
382if (present(found)) then
383    call check_nc(nf90_inq_varid(ncid_file,var_name,varid),'inquiring '//var_name,found)
384    call check_nc(nf90_get_var(ncid_file,varid,var),'getting '//var_name,found)
385else
386    call check_nc(nf90_inq_varid(ncid_file,var_name,varid),'inquiring '//var_name)
387    call check_nc(nf90_get_var(ncid_file,varid,var),'getting '//var_name)
388end if
389call check_valid_var_nc(var_name,var)
390
391END SUBROUTINE get_var1d_nc
392!=======================================================================
393
394!=======================================================================
395SUBROUTINE get_var2d_nc(var_name,var,found)
396!-----------------------------------------------------------------------
397! NAME
398!     get_var2d_nc
399!
400! DESCRIPTION
401!     Read a 2D variable from open NetCDF file.
402!
403! AUTHORS & DATE
404!     JB Clement, 12/2025
405!
406! NOTES
407!
408!-----------------------------------------------------------------------
409
410! DECLARATION
411! -----------
412implicit none
413
414! ARGUMENTS
415! ---------
416character(*),             intent(in)            :: var_name ! Variable name
417real(dp), dimension(:,:), intent(out)           :: var      ! Output variable
418logical(k4),              intent(out), optional :: found
419
420! CODE
421! ----
422if (present(found)) then
423    call check_nc(nf90_inq_varid(ncid_file,var_name,varid),'inquiring '//var_name,found)
424    call check_nc(nf90_get_var(ncid_file,varid,var),'getting '//var_name,found)
425else
426    call check_nc(nf90_inq_varid(ncid_file,var_name,varid),'inquiring '//var_name)
427    call check_nc(nf90_get_var(ncid_file,varid,var),'getting '//var_name)
428end if
429call check_valid_var_nc(var_name,var)
430
431END SUBROUTINE get_var2d_nc
432!=======================================================================
433
434!=======================================================================
435SUBROUTINE get_var3d_nc(var_name,var,found)
436!-----------------------------------------------------------------------
437! NAME
438!     get_var3d_nc
439!
440! DESCRIPTION
441!     Read a 3D variable from open NetCDF file.
442!
443! AUTHORS & DATE
444!     JB Clement, 12/2025
445!
446! NOTES
447!
448!-----------------------------------------------------------------------
449
450! DECLARATION
451! -----------
452implicit none
453
454! ARGUMENTS
455! ---------
456character(*),               intent(in)            :: var_name ! Variable name
457real(dp), dimension(:,:,:), intent(out)           :: var      ! Output variable
458logical(k4),                intent(out), optional :: found
459
460! CODE
461! ----
462if (present(found)) then
463    call check_nc(nf90_inq_varid(ncid_file,var_name,varid),'inquiring '//var_name,found)
464    call check_nc(nf90_get_var(ncid_file,varid,var),'getting '//var_name,found)
465else
466    call check_nc(nf90_inq_varid(ncid_file,var_name,varid),'inquiring '//var_name)
467    call check_nc(nf90_get_var(ncid_file,varid,var),'getting '//var_name)
468end if
469call check_valid_var_nc(var_name,var)
470
471END SUBROUTINE get_var3d_nc
472!=======================================================================
473
474!=======================================================================
475SUBROUTINE get_var4d_nc(var_name,var,found)
476!-----------------------------------------------------------------------
477! NAME
478!     get_var4d_nc
479!
480! DESCRIPTION
481!     Read a 4D variable from open NetCDF file.
482!
483! AUTHORS & DATE
484!     JB Clement, 12/2025
485!
486! NOTES
487!
488!-----------------------------------------------------------------------
489
490! DECLARATION
491! -----------
492implicit none
493
494! ARGUMENTS
495! ---------
496character(*),                 intent(in)            :: var_name ! Variable name
497real(dp), dimension(:,:,:,:), intent(out)           :: var      ! Output variable
498logical(k4),                  intent(out), optional :: found
499
500! CODE
501! ----
502if (present(found)) then
503    call check_nc(nf90_inq_varid(ncid_file,var_name,varid),'inquiring '//var_name,found)
504    call check_nc(nf90_get_var(ncid_file,varid,var),'getting '//var_name,found)
505else
506    call check_nc(nf90_inq_varid(ncid_file,var_name,varid),'inquiring '//var_name)
507    call check_nc(nf90_get_var(ncid_file,varid,var),'getting '//var_name)
508end if
509call check_valid_var_nc(var_name,var)
510
511END SUBROUTINE get_var4d_nc
512!=======================================================================
513
514!=======================================================================
515SUBROUTINE put_var0d_nc(var_name,var,itime)
516!-----------------------------------------------------------------------
517! NAME
518!     put_var0d_nc
519!
520! DESCRIPTION
521!     Write a 0D variable into open NetCDF file.
522!
523! AUTHORS & DATE
524!     JB Clement, 01/2026
525!
526! NOTES
527!
528!-----------------------------------------------------------------------
529
530! DECLARATION
531! -----------
532implicit none
533
534! ARGUMENTS
535! ---------
536character(*), intent(in)           :: var_name ! Variable name
537real(dp),     intent(in)           :: var      ! Input variable
538integer(di),  intent(in), optional :: itime    ! Current time index
539
540! LOCAL VARIABLES
541! ---------------
542integer(di)                               :: ndims, unlimdimid, ncid
543integer(di), dimension(nf90_max_var_dims) :: dimids
544logical(k4)                               :: has_time
545
546! CODE
547! ----
548! Diagevol logic priority over standard logic
549if (open_diagevol) then
550    ncid = ncid_diagevol
551else
552    ncid = ncid_file
553end if
554
555! Check if the variable holds a Time dimension (unlimited dim)
556call check_nc(nf90_inq_varid(ncid,var_name,varid),'inquiring '//var_name)
557call check_nc(nf90_inquire_variable(ncid,varid,ndims = ndims,dimids = dimids),'inquiring dims of '//var_name)
558if (ndims > 0) then
559    call check_nc(nf90_inquire(ncid,unlimitedDimId = unlimdimid),'inquiring unlimited dim')
560    has_time = (dimids(ndims) == unlimdimid)
561else
562    has_time = .false.
563end if
564
565if (present(itime)) then ! Time-dependent write
566    if (.not. has_time) call stop_clean(__FILE__,__LINE__,'itime present but variable has no Time dimension: '//var_name,1)
567    ! For 0D variable with time, just write at the time index
568    call check_nc(nf90_put_var(ncid,varid,var,start = (/itime/)),'putting '//var_name)
569else ! Static write
570    if (has_time) call stop_clean(__FILE__,__LINE__,'itime absent but variable is time-dependent: '//var_name,1)
571    call check_nc(nf90_put_var(ncid,varid,var),'putting '//var_name)
572end if
573
574END SUBROUTINE put_var0d_nc
575!=======================================================================
576
577!=======================================================================
578SUBROUTINE put_var1d_nc(var_name,var,itime)
579!-----------------------------------------------------------------------
580! NAME
581!     put_var1d_nc
582!
583! DESCRIPTION
584!     Write a 1D variable into open NetCDF file.
585!
586! AUTHORS & DATE
587!     JB Clement, 01/2026
588!
589! NOTES
590!
591!-----------------------------------------------------------------------
592
593! DECLARATION
594! -----------
595implicit none
596
597! ARGUMENTS
598! ---------
599character(*),           intent(in)           :: var_name ! Variable name
600real(dp), dimension(:), intent(in)           :: var      ! Input variable
601integer(di),            intent(in), optional :: itime    ! Current time index
602
603! LOCAL VARIABLES
604! ---------------
605integer(di)                               :: i, ndims, unlimdimid, ncid
606integer(di), dimension(nf90_max_var_dims) :: dimids
607integer(di), dimension(:), allocatable    :: strt, cnt
608logical(k4)                               :: has_time
609
610! CODE
611! ----
612! Diagevol logic priority over standard logic
613if (open_diagevol) then
614    ncid = ncid_diagevol
615else
616    ncid = ncid_file
617end if
618
619! Check if the variable holds a Time dimension (unlimited dim)
620call check_nc(nf90_inq_varid(ncid,var_name,varid),'inquiring '//var_name)
621call check_nc(nf90_inquire_variable(ncid,varid,ndims = ndims,dimids = dimids),'inquiring dims of '//var_name)
622if (ndims > 0) then
623    call check_nc(nf90_inquire(ncid,unlimitedDimId = unlimdimid),'inquiring unlimited dim')
624    has_time = (dimids(ndims) == unlimdimid)
625else
626    has_time = .false.
627end if
628
629if (present(itime)) then ! Time-dependent write
630    if (.not. has_time) call stop_clean(__FILE__,__LINE__,'itime present but variable has no Time dimension: '//var_name,1)
631    allocate(strt(ndims),cnt(ndims))
632    strt(:) = 1
633    cnt(:) = 1
634    strt(ndims) = itime
635    do i = 1,ndims - 1
636        cnt(i) = size(var,i)
637    end do
638    call check_nc(nf90_put_var(ncid,varid,var,start = strt,count = cnt),'putting '//var_name)
639    deallocate(strt,cnt)
640else ! Static write
641    if (has_time) call stop_clean(__FILE__,__LINE__,'itime absent but variable is time-dependent: '//var_name,1)
642    call check_nc(nf90_put_var(ncid,varid,var),'putting '//var_name)
643end if
644
645END SUBROUTINE put_var1d_nc
646!=======================================================================
647
648!=======================================================================
649SUBROUTINE put_var2d_nc(var_name,var,itime)
650!-----------------------------------------------------------------------
651! NAME
652!     put_var2d_nc
653!
654! DESCRIPTION
655!     Write a 2D variable into open NetCDF file.
656!
657! AUTHORS & DATE
658!     JB Clement, 01/2026
659!
660! NOTES
661!
662!-----------------------------------------------------------------------
663
664! DECLARATION
665! -----------
666implicit none
667
668! ARGUMENTS
669! ---------
670character(*),             intent(in)           :: var_name ! Variable name
671real(dp), dimension(:,:), intent(in)           :: var      ! Input variable
672integer(di),              intent(in), optional :: itime    ! Current time index
673
674! LOCAL VARIABLES
675! ---------------
676integer(di)                               :: i, ndims, unlimdimid, ncid
677integer(di), dimension(nf90_max_var_dims) :: dimids
678integer(di), dimension(:), allocatable    :: strt, cnt
679logical(k4)                               :: has_time
680
681! CODE
682! ----
683! Diagevol logic priority over standard logic
684if (open_diagevol) then
685    ncid = ncid_diagevol
686else
687    ncid = ncid_file
688end if
689
690! Check if the variable holds a Time dimension (unlimited dim)
691call check_nc(nf90_inq_varid(ncid,var_name,varid),'inquiring '//var_name)
692call check_nc(nf90_inquire_variable(ncid,varid,ndims = ndims,dimids = dimids),'inquiring dims of '//var_name)
693if (ndims > 0) then
694    call check_nc(nf90_inquire(ncid,unlimitedDimId = unlimdimid),'inquiring unlimited dim')
695    has_time = (dimids(ndims) == unlimdimid)
696else
697    has_time = .false.
698end if
699
700if (present(itime)) then ! Time-dependent write
701    if (.not. has_time) call stop_clean(__FILE__,__LINE__,'itime present but variable has no Time dimension: '//var_name,1)
702    allocate(strt(ndims),cnt(ndims))
703    strt(:) = 1
704    cnt(:) = 1
705    strt(ndims) = itime
706    do i = 1,ndims - 1
707        cnt(i) = size(var,i)
708    end do
709    call check_nc(nf90_put_var(ncid,varid,var,start = strt,count = cnt),'putting '//var_name)
710    deallocate(strt,cnt)
711else ! Static write
712    if (has_time) call stop_clean(__FILE__,__LINE__,'itime absent but variable is time-dependent: '//var_name,1)
713    call check_nc(nf90_put_var(ncid,varid,var),'putting '//var_name)
714end if
715
716END SUBROUTINE put_var2d_nc
717!=======================================================================
718
719!=======================================================================
720SUBROUTINE put_var3d_nc(var_name,var,itime)
721!-----------------------------------------------------------------------
722! NAME
723!     put_var3d_nc
724!
725! DESCRIPTION
726!     Write a 3D variable into open NetCDF file.
727!
728! AUTHORS & DATE
729!     JB Clement, 01/2026
730!
731! NOTES
732!
733!-----------------------------------------------------------------------
734
735! DECLARATION
736! -----------
737implicit none
738
739! ARGUMENTS
740! ---------
741character(*),               intent(in)           :: var_name ! Variable name
742real(dp), dimension(:,:,:), intent(in)           :: var      ! Input variable
743integer(di),                intent(in), optional :: itime    ! Current time index
744
745! LOCAL VARIABLES
746! ---------------
747integer(di)                               :: i, ndims, unlimdimid, ncid
748integer(di), dimension(nf90_max_var_dims) :: dimids
749integer(di), dimension(:), allocatable    :: strt, cnt
750logical(k4)                               :: has_time
751
752! CODE
753! ----
754! Diagevol logic priority over standard logic
755if (open_diagevol) then
756    ncid = ncid_diagevol
757else
758    ncid = ncid_file
759end if
760
761! Check if the variable holds a Time dimension (unlimited dim)
762call check_nc(nf90_inq_varid(ncid,var_name,varid),'inquiring '//var_name)
763call check_nc(nf90_inquire_variable(ncid,varid,ndims = ndims,dimids = dimids),'inquiring dims of '//var_name)
764if (ndims > 0) then
765    call check_nc(nf90_inquire(ncid,unlimitedDimId = unlimdimid),'inquiring unlimited dim')
766    has_time = (dimids(ndims) == unlimdimid)
767else
768    has_time = .false.
769end if
770
771if (present(itime)) then ! Time-dependent write
772    if (.not. has_time) call stop_clean(__FILE__,__LINE__,'itime present but variable has no Time dimension: '//var_name,1)
773    allocate(strt(ndims),cnt(ndims))
774    strt(:) = 1
775    cnt(:) = 1
776    strt(ndims) = itime
777    do i = 1,ndims - 1
778        cnt(i) = size(var,i)
779    end do
780    call check_nc(nf90_put_var(ncid,varid,var,start = strt,count = cnt),'putting '//var_name)
781    deallocate(strt,cnt)
782else ! Static write
783    if (has_time) call stop_clean(__FILE__,__LINE__,'itime absent but variable is time-dependent: '//var_name,1)
784    call check_nc(nf90_put_var(ncid,varid,var),'putting '//var_name)
785end if
786
787END SUBROUTINE put_var3d_nc
788!=======================================================================
789
790!=======================================================================
791SUBROUTINE put_var4d_nc(var_name,var,itime)
792!-----------------------------------------------------------------------
793! NAME
794!     put_var4d_nc
795!
796! DESCRIPTION
797!     Write a 4D variable into open NetCDF file.
798!
799! AUTHORS & DATE
800!     JB Clement, 01/2026
801!
802! NOTES
803!
804!-----------------------------------------------------------------------
805
806! DECLARATION
807! -----------
808implicit none
809
810! ARGUMENTS
811! ---------
812character(*),                 intent(in)           :: var_name ! Variable name
813real(dp), dimension(:,:,:,:), intent(in)           :: var      ! Input variable
814integer(di),                  intent(in), optional :: itime    ! Current time index
815
816! LOCAL VARIABLES
817! ---------------
818integer(di)                               :: i, ndims, unlimdimid, ncid
819integer(di), dimension(nf90_max_var_dims) :: dimids
820integer(di), dimension(:), allocatable    :: strt, cnt
821logical(k4)                               :: has_time
822
823! CODE
824! ----
825! Diagevol logic priority over standard logic
826if (open_diagevol) then
827    ncid = ncid_diagevol
828else
829    ncid = ncid_file
830end if
831
832! Check if the variable holds a Time dimension (unlimited dim)
833call check_nc(nf90_inq_varid(ncid,var_name,varid),'inquiring '//var_name)
834call check_nc(nf90_inquire_variable(ncid,varid,ndims = ndims,dimids = dimids),'inquiring dims of '//var_name)
835if (ndims > 0) then
836    call check_nc(nf90_inquire(ncid,unlimitedDimId = unlimdimid),'inquiring unlimited dim')
837    has_time = (dimids(ndims) == unlimdimid)
838else
839    has_time = .false.
840end if
841
842if (present(itime)) then ! Time-dependent write
843    if (.not. has_time) call stop_clean(__FILE__,__LINE__,'itime present but variable has no Time dimension: '//var_name,1)
844    allocate(strt(ndims),cnt(ndims))
845    strt(:) = 1
846    cnt(:) = 1
847    strt(ndims) = itime
848    do i = 1,ndims - 1
849        cnt(i) = size(var,i)
850    end do
851    call check_nc(nf90_put_var(ncid,varid,var,start = strt,count = cnt),'putting '//var_name)
852    deallocate(strt,cnt)
853else ! Static write
854    if (has_time) call stop_clean(__FILE__,__LINE__,'itime absent but variable is time-dependent: '//var_name,1)
855    call check_nc(nf90_put_var(ncid,varid,var),'putting '//var_name)
856end if
857
858END SUBROUTINE put_var4d_nc
859!=======================================================================
860
861!=======================================================================
862SUBROUTINE get_next_itime_nc(dim_name,itime)
863!-----------------------------------------------------------------------
864! NAME
865!     get_next_itime_nc
866!
867! DESCRIPTION
868!     Get the next time index in a NetCDF file to record variables.
869!
870! AUTHORS & DATE
871!     JB Clement, 01/2026
872!
873! NOTES
874!     In most cases, dim_name = 'Time'.
875!-----------------------------------------------------------------------
876
877! DECLARATION
878! -----------
879implicit none
880
881! ARGUMENTS
882! ---------
883character(*), intent(in)  :: dim_name ! Dimension name
884integer(di),  intent(out) :: itime    ! Time index
885
886! LOCAL VARIABLES
887! ---------------
888integer(di) :: length
889logical(k4) :: found
890
891! CODE
892! ----
893! Get dimension length
894call get_dim_nc(dim_name,length,found)
895if (.not. found) call stop_clean(__FILE__,__LINE__,'dimension '//dim_name//' not found',1)
896
897! Next time index
898itime = length + 1
899
900END SUBROUTINE get_next_itime_nc
901!=======================================================================
902
903!=======================================================================
904SUBROUTINE check_valid_var0d_nc(var_name,var)
905!-----------------------------------------------------------------------
906! NAME
907!     check_valid_var0d_nc
908!
909! DESCRIPTION
910!     Check the validity of a 0D variable read from a NetCDF file.
911!
912! AUTHORS & DATE
913!     JB Clement, 02/2026
914!
915! NOTES
916!
917!-----------------------------------------------------------------------
918
919! DEPENDENCIES
920! ------------
921use numerics, only: largest_nb
922use stoppage, only: stop_clean
923use utility,  only: real2str
924
925! DECLARATION
926! -----------
927implicit none
928
929! ARGUMENTS
930! ---------
931character(*), intent(in) :: var_name ! Variable name
932real(dp),     intent(in) :: var      ! Input variable
933
934! LOCAL VARIABLES
935! ---------------
936logical(k4)            :: has_fill, has_range, has_valid_min, has_valid_max
937real(dp)               :: fill_value, valid_min, valid_max
938real(dp), dimension(2) :: valid_range
939integer(di)            :: ncid
940
941! CODE
942! ----
943! Diagevol logic priority over standard logic
944if (open_diagevol) then
945    ncid = ncid_diagevol
946else
947    ncid = ncid_file
948end if
949
950! NaN
951if (var /= var) call stop_clean(__FILE__,__LINE__,'NaN detected in variable '//var_name//'!',1)
952
953! Infinite
954if (abs(var) > largest_nb) call stop_clean(__FILE__,__LINE__,'Infs detected in variable '//var_name//'!',1)
955
956! Fill value
957has_fill = .false.
958call check_nc(nf90_get_att(ncid,varid,"_FillValue",fill_value),'getting fill value',has_fill)
959if (.not. has_fill) call check_nc(nf90_get_att(ncid,varid,"missing_value",fill_value),'getting missing value',has_fill)
960if (has_fill) then
961    if (abs(var - fill_value) < minieps) call stop_clean(__FILE__,__LINE__,'Fill values ('//real2str(fill_value)//') detected in '//var_name//'!',1)
962end if
963
964! Valid range
965has_range = .false.
966call check_nc(nf90_get_att(ncid,varid,"valid_range",valid_range),'getting valid range',has_range)
967if (has_range) then
968    valid_min = valid_range(1)
969    valid_max = valid_range(2)
970else! valid_min / valid_max (fallback)
971    has_valid_min = .false.
972    has_valid_max = .false.
973    call check_nc(nf90_get_att(ncid,varid,"valid_min",valid_min),'getting valid min',has_valid_min)
974    call check_nc(nf90_get_att(ncid,varid,"valid_max",valid_max),'getting valid max',has_valid_max)
975    if (has_valid_min .or. has_valid_max) then
976        has_range = .true.
977        if (.not. has_valid_min) valid_min = -largest_nb
978        if (.not. has_valid_max) valid_max = largest_nb
979    end if
980end if
981if (has_range) then
982    if (var < valid_min .or. var > valid_max) call stop_clean(__FILE__,__LINE__,'Values outside valid range ('//real2str(valid_min)//','//real2str(valid_max)//') detected in '//var_name//'!',1)
983end if
984
985END SUBROUTINE check_valid_var0d_nc
986!=======================================================================
987
988!=======================================================================
989SUBROUTINE check_valid_var1d_nc(var_name,var)
990!-----------------------------------------------------------------------
991! NAME
992!     check_valid_var1d_nc
993!
994! DESCRIPTION
995!     Check the validity of a 1D variable read from a NetCDF file.
996!
997! AUTHORS & DATE
998!     JB Clement, 02/2026
999!
1000! NOTES
1001!
1002!-----------------------------------------------------------------------
1003
1004! DEPENDENCIES
1005! ------------
1006use numerics, only: largest_nb
1007use stoppage, only: stop_clean
1008use utility,  only: real2str
1009
1010! DECLARATION
1011! -----------
1012implicit none
1013
1014! ARGUMENTS
1015! ---------
1016character(*),           intent(in) :: var_name ! Variable name
1017real(dp), dimension(:), intent(in) :: var      ! Input variable
1018
1019! LOCAL VARIABLES
1020! ---------------
1021logical(k4)            :: has_fill, has_range, has_valid_min, has_valid_max
1022real(dp)               :: fill_value, valid_min, valid_max
1023real(dp), dimension(2) :: valid_range
1024integer(di)            :: ncid
1025
1026! CODE
1027! ----
1028! Diagevol logic priority over standard logic
1029if (open_diagevol) then
1030    ncid = ncid_diagevol
1031else
1032    ncid = ncid_file
1033end if
1034
1035! NaN
1036if (any(var /= var)) call stop_clean(__FILE__,__LINE__,'NaN detected in variable '//var_name//'!',1)
1037
1038! Infinite
1039if (any(abs(var) > largest_nb)) call stop_clean(__FILE__,__LINE__,'Infs detected in variable '//var_name//'!',1)
1040
1041! Fill value
1042has_fill = .false.
1043call check_nc(nf90_get_att(ncid,varid,"_FillValue",fill_value),'getting fill value',has_fill)
1044if (.not. has_fill) call check_nc(nf90_get_att(ncid,varid,"missing_value",fill_value),'getting missing value',has_fill)
1045if (has_fill) then
1046    if (any(abs(var - fill_value) < minieps)) call stop_clean(__FILE__,__LINE__,'Fill values ('//real2str(fill_value)//') detected in '//var_name//'!',1)
1047end if
1048
1049! Valid range
1050has_range = .false.
1051call check_nc(nf90_get_att(ncid,varid,"valid_range",valid_range),'getting valid range',has_range)
1052if (has_range) then
1053    valid_min = valid_range(1)
1054    valid_max = valid_range(2)
1055else! valid_min / valid_max (fallback)
1056    has_valid_min = .false.
1057    has_valid_max = .false.
1058    call check_nc(nf90_get_att(ncid,varid,"valid_min",valid_min),'getting valid min',has_valid_min)
1059    call check_nc(nf90_get_att(ncid,varid,"valid_max",valid_max),'getting valid max',has_valid_max)
1060    if (has_valid_min .or. has_valid_max) then
1061        has_range = .true.
1062        if (.not. has_valid_min) valid_min = -largest_nb
1063        if (.not. has_valid_max) valid_max = largest_nb
1064    end if
1065end if
1066if (has_range) then
1067    if (any(var < valid_min) .or. any(var > valid_max)) call stop_clean(__FILE__,__LINE__,'Values outside valid range ('//real2str(valid_min)//','//real2str(valid_max)//') detected in '//var_name//'!',1)
1068end if
1069
1070END SUBROUTINE check_valid_var1d_nc
1071!=======================================================================
1072
1073!=======================================================================
1074SUBROUTINE check_valid_var2d_nc(var_name,var)
1075!-----------------------------------------------------------------------
1076! NAME
1077!     check_valid_var2d_nc
1078!
1079! DESCRIPTION
1080!     Check the validity of a 1D variable read from a NetCDF file.
1081!
1082! AUTHORS & DATE
1083!     JB Clement, 02/2026
1084!
1085! NOTES
1086!
1087!-----------------------------------------------------------------------
1088
1089! DEPENDENCIES
1090! ------------
1091use numerics, only: largest_nb
1092use stoppage, only: stop_clean
1093use utility,  only: real2str
1094
1095! DECLARATION
1096! -----------
1097implicit none
1098
1099! ARGUMENTS
1100! ---------
1101character(*),             intent(in) :: var_name ! Variable name
1102real(dp), dimension(:,:), intent(in) :: var      ! Input variable
1103
1104! LOCAL VARIABLES
1105! ---------------
1106logical(k4)            :: has_fill, has_range, has_valid_min, has_valid_max
1107real(dp)               :: fill_value, valid_min, valid_max
1108real(dp), dimension(2) :: valid_range
1109integer(di)            :: ncid
1110
1111! CODE
1112! ----
1113! Diagevol logic priority over standard logic
1114if (open_diagevol) then
1115    ncid = ncid_diagevol
1116else
1117    ncid = ncid_file
1118end if
1119
1120! NaN
1121if (any(var /= var)) call stop_clean(__FILE__,__LINE__,'NaN detected in variable '//var_name//'!',1)
1122
1123! Infinite
1124if (any(abs(var) > largest_nb)) call stop_clean(__FILE__,__LINE__,'Infs detected in variable '//var_name//'!',1)
1125
1126! Fill value
1127has_fill = .false.
1128call check_nc(nf90_get_att(ncid,varid,"_FillValue",fill_value),'getting fill value',has_fill)
1129if (.not. has_fill) call check_nc(nf90_get_att(ncid,varid,"missing_value",fill_value),'getting missing value',has_fill)
1130if (has_fill) then
1131    if (any(abs(var - fill_value) < minieps)) call stop_clean(__FILE__,__LINE__,'Fill values ('//real2str(fill_value)//') detected in '//var_name//'!',1)
1132end if
1133
1134! Valid range
1135has_range = .false.
1136call check_nc(nf90_get_att(ncid,varid,"valid_range",valid_range),'getting valid range',has_range)
1137if (has_range) then
1138    valid_min = valid_range(1)
1139    valid_max = valid_range(2)
1140else! valid_min / valid_max (fallback)
1141    has_valid_min = .false.
1142    has_valid_max = .false.
1143    call check_nc(nf90_get_att(ncid,varid,"valid_min",valid_min),'getting valid min',has_valid_min)
1144    call check_nc(nf90_get_att(ncid,varid,"valid_max",valid_max),'getting valid max',has_valid_max)
1145    if (has_valid_min .or. has_valid_max) then
1146        has_range = .true.
1147        if (.not. has_valid_min) valid_min = -largest_nb
1148        if (.not. has_valid_max) valid_max = largest_nb
1149    end if
1150end if
1151if (has_range) then
1152    if (any(var < valid_min) .or. any(var > valid_max)) call stop_clean(__FILE__,__LINE__,'Values outside valid range ('//real2str(valid_min)//','//real2str(valid_max)//') detected in '//var_name//'!',1)
1153end if
1154
1155END SUBROUTINE check_valid_var2d_nc
1156!=======================================================================
1157
1158!=======================================================================
1159SUBROUTINE check_valid_var3d_nc(var_name,var)
1160!-----------------------------------------------------------------------
1161! NAME
1162!     check_valid_var3d_nc
1163!
1164! DESCRIPTION
1165!     Check the validity of a 1D variable read from a NetCDF file.
1166!
1167! AUTHORS & DATE
1168!     JB Clement, 02/2026
1169!
1170! NOTES
1171!
1172!-----------------------------------------------------------------------
1173
1174! DEPENDENCIES
1175! ------------
1176use numerics, only: largest_nb
1177use stoppage, only: stop_clean
1178use utility,  only: real2str
1179
1180! DECLARATION
1181! -----------
1182implicit none
1183
1184! ARGUMENTS
1185! ---------
1186character(*),               intent(in) :: var_name ! Variable name
1187real(dp), dimension(:,:,:), intent(in) :: var      ! Input variable
1188
1189! LOCAL VARIABLES
1190! ---------------
1191logical(k4)            :: has_fill, has_range, has_valid_min, has_valid_max
1192real(dp)               :: fill_value, valid_min, valid_max
1193real(dp), dimension(2) :: valid_range
1194integer(di)            :: ncid
1195
1196! CODE
1197! ----
1198! Diagevol logic priority over standard logic
1199if (open_diagevol) then
1200    ncid = ncid_diagevol
1201else
1202    ncid = ncid_file
1203end if
1204
1205! NaN
1206if (any(var /= var)) call stop_clean(__FILE__,__LINE__,'NaN detected in variable '//var_name//'!',1)
1207
1208! Infinite
1209if (any(abs(var) > largest_nb)) call stop_clean(__FILE__,__LINE__,'Infs detected in variable '//var_name//'!',1)
1210
1211! Fill value
1212has_fill = .false.
1213call check_nc(nf90_get_att(ncid,varid,"_FillValue",fill_value),'getting fill value',has_fill)
1214if (.not. has_fill) call check_nc(nf90_get_att(ncid,varid,"missing_value",fill_value),'getting missing value',has_fill)
1215if (has_fill) then
1216    if (any(abs(var - fill_value) < minieps)) call stop_clean(__FILE__,__LINE__,'Fill values ('//real2str(fill_value)//') detected in '//var_name//'!',1)
1217end if
1218
1219! Valid range
1220has_range = .false.
1221call check_nc(nf90_get_att(ncid,varid,"valid_range",valid_range),'getting valid range',has_range)
1222if (has_range) then
1223    valid_min = valid_range(1)
1224    valid_max = valid_range(2)
1225else! valid_min / valid_max (fallback)
1226    has_valid_min = .false.
1227    has_valid_max = .false.
1228    call check_nc(nf90_get_att(ncid,varid,"valid_min",valid_min),'getting valid min',has_valid_min)
1229    call check_nc(nf90_get_att(ncid,varid,"valid_max",valid_max),'getting valid max',has_valid_max)
1230    if (has_valid_min .or. has_valid_max) then
1231        has_range = .true.
1232        if (.not. has_valid_min) valid_min = -largest_nb
1233        if (.not. has_valid_max) valid_max = largest_nb
1234    end if
1235end if
1236if (has_range) then
1237    if (any(var < valid_min) .or. any(var > valid_max)) call stop_clean(__FILE__,__LINE__,'Values outside valid range ('//real2str(valid_min)//','//real2str(valid_max)//') detected in '//var_name//'!',1)
1238end if
1239
1240END SUBROUTINE check_valid_var3d_nc
1241!=======================================================================
1242
1243!=======================================================================
1244SUBROUTINE check_valid_var4d_nc(var_name,var)
1245!-----------------------------------------------------------------------
1246! NAME
1247!     check_valid_var4d_nc
1248!
1249! DESCRIPTION
1250!     Check the validity of a 1D variable read from a NetCDF file.
1251!
1252! AUTHORS & DATE
1253!     JB Clement, 02/2026
1254!
1255! NOTES
1256!
1257!-----------------------------------------------------------------------
1258
1259! DEPENDENCIES
1260! ------------
1261use numerics, only: largest_nb
1262use stoppage, only: stop_clean
1263use utility,  only: real2str
1264
1265! DECLARATION
1266! -----------
1267implicit none
1268
1269! ARGUMENTS
1270! ---------
1271character(*),                 intent(in) :: var_name ! Variable name
1272real(dp), dimension(:,:,:,:), intent(in) :: var      ! Input variable
1273
1274! LOCAL VARIABLES
1275! ---------------
1276logical(k4)            :: has_fill, has_range, has_valid_min, has_valid_max
1277real(dp)               :: fill_value, valid_min, valid_max
1278real(dp), dimension(2) :: valid_range
1279integer(di)            :: ncid
1280
1281! CODE
1282! ----
1283! Diagevol logic priority over standard logic
1284if (open_diagevol) then
1285    ncid = ncid_diagevol
1286else
1287    ncid = ncid_file
1288end if
1289
1290! NaN
1291if (any(var /= var)) call stop_clean(__FILE__,__LINE__,'NaN detected in variable '//var_name//'!',1)
1292
1293! Infinite
1294if (any(abs(var) > largest_nb)) call stop_clean(__FILE__,__LINE__,'Infs detected in variable '//var_name//'!',1)
1295
1296! Fill value
1297has_fill = .false.
1298call check_nc(nf90_get_att(ncid,varid,"_FillValue",fill_value),'getting fill value',has_fill)
1299if (.not. has_fill) call check_nc(nf90_get_att(ncid,varid,"missing_value",fill_value),'getting missing value',has_fill)
1300if (has_fill) then
1301    if (any(abs(var - fill_value) < minieps)) call stop_clean(__FILE__,__LINE__,'Fill values ('//real2str(fill_value)//') detected in '//var_name//'!',1)
1302end if
1303
1304! Valid range
1305has_range = .false.
1306call check_nc(nf90_get_att(ncid,varid,"valid_range",valid_range),'getting valid range',has_range)
1307if (has_range) then
1308    valid_min = valid_range(1)
1309    valid_max = valid_range(2)
1310else! valid_min / valid_max (fallback)
1311    has_valid_min = .false.
1312    has_valid_max = .false.
1313    call check_nc(nf90_get_att(ncid,varid,"valid_min",valid_min),'getting valid min',has_valid_min)
1314    call check_nc(nf90_get_att(ncid,varid,"valid_max",valid_max),'getting valid max',has_valid_max)
1315    if (has_valid_min .or. has_valid_max) then
1316        has_range = .true.
1317        if (.not. has_valid_min) valid_min = -largest_nb
1318        if (.not. has_valid_max) valid_max = largest_nb
1319    end if
1320end if
1321if (has_range) then
1322    if(any(var < valid_min) .or. any(var > valid_max)) call stop_clean(__FILE__,__LINE__,'Values outside valid range ('//real2str(valid_min)//','//real2str(valid_max)//') detected in '//var_name//'!',1)
1323end if
1324
1325END SUBROUTINE check_valid_var4d_nc
1326!=======================================================================
1327
1328END MODULE io_netcdf
Note: See TracBrowser for help on using the repository browser.