source: trunk/LMDZ.COMMON/libf/evolution/iostart_pem.F90 @ 3991

Last change on this file since 3991 was 3991, checked in by jbclement, 6 days ago

PEM:
Apply documentation template everywhere: standardized headers format with short description, separators between functions/subroutines, normalized code sections, aligned dependencies/arguments/variables declaration.
JBC

File size: 52.8 KB
Line 
1MODULE iostart_pem
2!-----------------------------------------------------------------------
3! NAME
4!     iostart_pem
5!
6! DESCRIPTION
7!     Read and write specific netcdf files for the PEM (start and restart)
8!     Provides interfaces for field and variable I/O operations
9!
10! AUTHORS & DATE
11!     L. Lange
12!     JB Clement, 2023-2025
13!
14! NOTES
15!     Inspired by iostart from the PCM.
16!
17!-----------------------------------------------------------------------
18
19    ! DECLARATION
20! -----------
21implicit none
22
23! MODULE VARIABLES
24! ----------------
25PRIVATE
26    INTEGER,SAVE :: nid_start ! NetCDF file identifier for startfi.nc file
27    INTEGER,SAVE :: nid_restart ! NetCDF file identifier for restartfi.nc file
28
29    ! restartfi.nc file dimension identifiers: (see open_restartphy())
30    INTEGER,SAVE :: idim1  ! "index" dimension
31    INTEGER,SAVE :: idim2  ! "physical_points" dimension
32    INTEGER,SAVE :: idim3  ! "subsurface_layers" dimension
33    INTEGER,SAVE :: idim4  ! "nlayer_plus_1" dimension
34    INTEGER,SAVE :: idim5  ! "number_of_advected_fields" dimension
35    INTEGER,SAVE :: idim6  ! "nlayer" dimension
36    INTEGER,SAVE :: idim7  ! "Time" dimension
37    INTEGER,SAVE :: idim8  ! "nslope" dimension
38    INTEGER,SAVE :: idim9  ! "inter slope" dimension
39    INTEGER,SAVE :: idim10 ! "nb_str_max" dimension
40    INTEGER,SAVE :: timeindex ! current time index (for time-dependent fields)
41    INTEGER,PARAMETER :: length=100 ! size of tab_cntrl array
42
43    INTERFACE get_field
44      MODULE PROCEDURE Get_field_r1,Get_field_r2,Get_field_r3
45    END INTERFACE get_field
46
47    INTERFACE get_var
48      MODULE PROCEDURE get_var_r0,Get_var_r1,Get_var_r2,Get_var_r3
49    END INTERFACE get_var
50
51    INTERFACE put_field
52      MODULE PROCEDURE put_field_r1,put_field_r2,put_field_r3
53    END INTERFACE put_field
54
55    INTERFACE put_var
56      MODULE PROCEDURE put_var_r0,put_var_r1,put_var_r2,put_var_r3
57    END INTERFACE put_var
58
59    PUBLIC nid_start, length
60    PUBLIC get_field,get_var,put_field,put_var
61    PUBLIC inquire_dimension, inquire_dimension_length
62    PUBLIC inquire_field, inquire_field_ndims
63    PUBLIC open_startphy,close_startphy,open_restartphy,close_restartphy
64
65CONTAINS
66!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
67
68!=======================================================================
69  SUBROUTINE open_startphy(filename)
70!-----------------------------------------------------------------------
71! NAME
72!     open_startphy
73!
74! DESCRIPTION
75!     Open the startphy netCDF file for reading.
76
77! AUTHORS & DATE
78!     L. Lange
79!     JB Clement, 2023-2025
80
81! NOTES
82!
83!-----------------------------------------------------------------------
84
85! DEPENDENCIES
86! ------------
87  USE netcdf,             only: NF90_OPEN, NF90_NOERR, NF90_NOWRITE, nf90_strerror
88  USE mod_phys_lmdz_para, only: is_master, bcast
89
90! DECLARATION
91! -----------
92implicit none
93
94! ARGUMENTS
95! ---------
96    CHARACTER(LEN=*) :: filename
97
98! LOCAL VARIABLES
99! ---------------
100    INTEGER :: ierr
101
102! CODE
103! ----
104    IF (is_master) THEN
105      ierr = NF90_OPEN (filename, NF90_NOWRITE, nid_start)
106      IF (ierr.NE.NF90_NOERR) THEN
107        write(*,*)'open_startphy: problem opening file '//trim(filename)
108        write(*,*)trim(nf90_strerror(ierr))
109        CALL abort_physic("open_startphy","Cannot open file",1)
110      ENDIF
111    ENDIF
112
113    CALL bcast(nid_start) ! tell all procs about nid_start
114
115  END SUBROUTINE open_startphy
116!=======================================================================
117
118!=======================================================================
119  SUBROUTINE close_startphy
120!-----------------------------------------------------------------------
121! NAME
122!     close_startphy
123!
124! DESCRIPTION
125!     Close the startphy netCDF file.
126
127! AUTHORS & DATE
128!     L. Lange
129!     JB Clement, 2023-2025
130
131! NOTES
132!
133!-----------------------------------------------------------------------
134
135! DEPENDENCIES
136! ------------
137  USE netcdf,             only: NF90_CLOSE
138  USE mod_phys_lmdz_para, only: is_master
139
140! DECLARATION
141! -----------
142implicit none
143
144! LOCAL VARIABLES
145! ---------------
146    INTEGER :: ierr
147
148! CODE
149! ----
150    IF (is_master) THEN
151        ierr = NF90_CLOSE (nid_start)
152    ENDIF
153
154  END SUBROUTINE close_startphy
155!=======================================================================
156
157!=======================================================================
158  FUNCTION inquire_field(Field_name)
159!-----------------------------------------------------------------------
160! NAME
161!     inquire_field
162!
163! DESCRIPTION
164!     Check if a given field is present in the input file.
165
166! AUTHORS & DATE
167!     L. Lange
168!     JB Clement, 2023-2025
169
170! NOTES
171!
172!-----------------------------------------------------------------------
173
174! DEPENDENCIES
175! ------------
176  USE netcdf,             only: NF90_INQ_VARID, NF90_NOERR
177  USE mod_phys_lmdz_para, only: is_master, bcast
178
179! DECLARATION
180! -----------
181implicit none
182
183! ARGUMENTS
184! ---------
185    CHARACTER(LEN=*), INTENT(IN) :: Field_name
186
187! LOCAL VARIABLES
188! ---------------
189    LOGICAL :: inquire_field
190    INTEGER :: varid
191    INTEGER :: ierr
192
193! CODE
194! ----
195    IF (is_master) THEN
196      ierr=NF90_INQ_VARID(nid_start,Field_name,varid)
197      IF (ierr==NF90_NOERR) THEN
198        Inquire_field=.TRUE.
199      ELSE
200        Inquire_field=.FALSE.
201      ENDIF
202    ENDIF
203
204    CALL bcast(inquire_field)
205
206  END FUNCTION inquire_field
207!=======================================================================
208
209!=======================================================================
210  FUNCTION inquire_field_ndims(Field_name)
211!-----------------------------------------------------------------------
212! NAME
213!     inquire_field_ndims
214!
215! DESCRIPTION
216!     Give the number of dimensions of "Field_name" stored in the input file.
217
218! AUTHORS & DATE
219!     L. Lange
220!     JB Clement, 2023-2025
221
222! NOTES
223!
224!-----------------------------------------------------------------------
225
226! DEPENDENCIES
227! ------------
228  USE netcdf,             only: nf90_inq_varid, nf90_inquire_variable, &
229                                NF90_NOERR, nf90_strerror
230  USE mod_phys_lmdz_para, only: is_master, bcast
231
232! DECLARATION
233! -----------
234implicit none
235
236! ARGUMENTS
237! ---------
238    CHARACTER(LEN=*), INTENT(IN) :: Field_name
239
240! LOCAL VARIABLES
241! ---------------
242    INTEGER :: inquire_field_ndims
243    INTEGER :: varid
244    INTEGER :: ierr
245
246! CODE
247! ----
248    IF (is_master) THEN
249      ierr=nf90_inq_varid(nid_start,Field_name,varid)
250      ierr=nf90_inquire_variable(nid_start,varid,&
251                                  ndims=inquire_field_ndims)
252      IF (ierr.NE.NF90_NOERR) THEN
253        write(*,*)'inquire_field_ndims: problem obtaining ndims of '&
254                  //trim(field_name)
255        write(*,*)trim(nf90_strerror(ierr))
256        CALL abort_physic("inquire_field_ndims","Failed to get ndims",1)
257      ENDIF
258    ENDIF
259
260    CALL bcast(inquire_field_ndims)
261
262  END FUNCTION inquire_field_ndims
263!=======================================================================
264
265!=======================================================================
266  FUNCTION inquire_dimension(Field_name)
267!-----------------------------------------------------------------------
268! NAME
269!     inquire_dimension
270!
271! DESCRIPTION
272!     Check if a given dimension is present in the input file.
273
274! AUTHORS & DATE
275!     L. Lange
276!     JB Clement, 2023-2025
277
278! NOTES
279!
280!-----------------------------------------------------------------------
281
282! DEPENDENCIES
283! ------------
284  USE netcdf, only: nf90_inq_dimid, NF90_NOERR
285  USE mod_phys_lmdz_para, only: is_master, bcast
286
287! DECLARATION
288! -----------
289implicit none
290
291! ARGUMENTS
292! ---------
293    CHARACTER(LEN=*), INTENT(IN) :: Field_name
294
295! LOCAL VARIABLES
296! ---------------
297    LOGICAL :: inquire_dimension
298    INTEGER :: varid
299    INTEGER :: ierr
300
301! CODE
302! ----
303    IF (is_master) THEN
304      ierr=NF90_INQ_DIMID(nid_start,Field_name,varid)
305      IF (ierr==NF90_NOERR) THEN
306        Inquire_dimension=.TRUE.
307      ELSE
308        Inquire_dimension=.FALSE.
309      ENDIF
310    ENDIF
311
312    CALL bcast(inquire_dimension)
313
314  END FUNCTION inquire_dimension
315!=======================================================================
316
317!=======================================================================
318  FUNCTION inquire_dimension_length(Field_name)
319!-----------------------------------------------------------------------
320! NAME
321!     inquire_dimension_length
322!
323! DESCRIPTION
324!     Give the length of the "Field_name" dimension stored in the input file.
325
326! AUTHORS & DATE
327!     L. Lange
328!     JB Clement, 2023-2025
329
330! NOTES
331!
332!-----------------------------------------------------------------------
333
334! DEPENDENCIES
335! ------------
336  USE netcdf, only: nf90_inquire_dimension, nf90_inq_dimid, &
337                    NF90_NOERR, nf90_strerror
338  USE mod_phys_lmdz_para, only: is_master, bcast
339
340! DECLARATION
341! -----------
342implicit none
343
344! ARGUMENTS
345! ---------
346    CHARACTER(LEN=*), INTENT(IN) :: Field_name
347
348! LOCAL VARIABLES
349! ---------------
350    INTEGER :: inquire_dimension_length
351    INTEGER :: varid
352    INTEGER :: ierr
353
354! CODE
355! ----
356    IF (is_master) THEN
357      ierr=nf90_inq_dimid(nid_start,Field_name,varid)
358      ierr=nf90_inquire_dimension(nid_start,varid,&
359                                  len=inquire_dimension_length)
360      IF (ierr.NE.NF90_NOERR) THEN
361        write(*,*)'inquire_field_length: problem obtaining length of '&
362                  //trim(field_name)
363        write(*,*)trim(nf90_strerror(ierr))
364        CALL abort_physic("inquire_field_ndims","Failed to get length",1)
365      ENDIF
366    ENDIF
367
368    CALL bcast(inquire_dimension_length)
369
370  END FUNCTION inquire_dimension_length
371!=======================================================================
372
373!=======================================================================
374  SUBROUTINE Get_Field_r1(field_name,field,found,timeindex)
375!-----------------------------------------------------------------------
376! NAME
377!     Get_Field_r1
378!
379! DESCRIPTION
380!     Read a surface field (1D) from the start file.
381
382! AUTHORS & DATE
383!     L. Lange
384!     JB Clement, 2023-2025
385
386! NOTES
387!
388!-----------------------------------------------------------------------
389
390! DEPENDENCIES
391! ------------
392  use mod_grid_phy_lmdz, only: klon_glo
393
394! DECLARATION
395! -----------
396implicit none
397
398! ARGUMENTS
399! ---------
400  CHARACTER(LEN=*), INTENT(IN)            :: field_name
401  REAL,             INTENT(INOUT)         :: field(:)
402  LOGICAL,          INTENT(OUT), OPTIONAL :: found
403  INTEGER,          INTENT(IN),  OPTIONAL :: timeindex
404
405! LOCAL VARIABLES
406! ---------------
407    integer :: edges(4), corners(4)
408
409! CODE
410! ----
411    edges(1)=klon_glo
412    edges(2:4)=1
413    corners(1:4)=1
414    if (PRESENT(timeindex)) then
415      corners(2)=timeindex
416    endif
417
418    IF (PRESENT(found)) THEN
419      CALL Get_field_rgen(field_name,field,1,corners,edges,found)
420    ELSE
421      CALL Get_field_rgen(field_name,field,1,corners,edges)
422    ENDIF
423
424  END SUBROUTINE Get_Field_r1
425!=======================================================================
426
427!=======================================================================
428  SUBROUTINE Get_Field_r2(field_name,field,found,timeindex)
429!-----------------------------------------------------------------------
430! NAME
431!     Get_Field_r2
432!
433! DESCRIPTION
434!     Read a "3D" horizontal-vertical field (2D) from the start file.
435
436! AUTHORS & DATE
437!     L. Lange
438!     JB Clement, 2023-2025
439
440! NOTES
441!
442!-----------------------------------------------------------------------
443
444! DEPENDENCIES
445! ------------
446  use mod_grid_phy_lmdz, only: klon_glo
447
448! DECLARATION
449! -----------
450implicit none
451
452! ARGUMENTS
453! ---------
454  CHARACTER(LEN=*), INTENT(IN)            :: field_name
455  REAL,             INTENT(INOUT)         :: field(:,:)
456  LOGICAL,          INTENT(OUT), OPTIONAL :: found
457  INTEGER,          INTENT(IN),  OPTIONAL :: timeindex
458
459! LOCAL VARIABLES
460! ---------------
461    integer :: edges(4), corners(4)
462
463! CODE
464! ----
465    edges(1)=klon_glo
466    edges(2)=size(field,2)
467    edges(3:4)=1
468    corners(1:4)=1
469    if (PRESENT(timeindex)) then
470      corners(3)=timeindex
471    endif
472
473    IF (PRESENT(found)) THEN
474      CALL Get_field_rgen(field_name,field,size(field,2),&
475                          corners,edges,found)
476    ELSE
477      CALL Get_field_rgen(field_name,field,size(field,2),&
478                          corners,edges)
479    ENDIF
480
481
482  END SUBROUTINE Get_Field_r2
483!=======================================================================
484
485!=======================================================================
486  SUBROUTINE Get_Field_r3(field_name,field,found,timeindex)
487!-----------------------------------------------------------------------
488! NAME
489!     Get_Field_r3
490!
491! DESCRIPTION
492!     Read a "4D" field surf/alt/?? (3D) from the start file.
493
494! AUTHORS & DATE
495!     L. Lange
496!     JB Clement, 2023-2025
497
498! NOTES
499!
500!-----------------------------------------------------------------------
501
502! DEPENDENCIES
503! ------------
504  use mod_grid_phy_lmdz, only: klon_glo
505
506! DECLARATION
507! -----------
508implicit none
509
510! ARGUMENTS
511! ---------
512  CHARACTER(LEN=*), INTENT(IN)            :: field_name
513  REAL,             INTENT(INOUT)         :: field(:,:,:)
514  LOGICAL,          INTENT(OUT), OPTIONAL :: found
515  INTEGER,          INTENT(IN),  OPTIONAL :: timeindex
516
517! LOCAL VARIABLES
518! ---------------
519    integer :: edges(4), corners(4)
520
521! CODE
522! ----
523    edges(1)=klon_glo
524    edges(2)=size(field,2)
525    edges(3)=size(field,3)
526    edges(4)=1
527    corners(1:4)=1
528    if (PRESENT(timeindex)) then
529      corners(4)=timeindex
530    endif
531
532    IF (PRESENT(found)) THEN
533      CALL Get_field_rgen(field_name,field,size(field,2)*size(field,3),&
534                          corners,edges,found)
535    ELSE
536      CALL Get_field_rgen(field_name,field,size(field,2)*size(field,3),&
537                          corners,edges)
538    ENDIF
539
540  END SUBROUTINE Get_Field_r3
541!=======================================================================
542
543!=======================================================================
544  SUBROUTINE Get_field_rgen(field_name,field,field_size, &
545                            corners,edges,found)
546!-----------------------------------------------------------------------
547! NAME
548!     Get_field_rgen
549!
550! DESCRIPTION
551!     Generic subroutine to read fields from start file with scatter to processes.
552
553! AUTHORS & DATE
554!     L. Lange
555!     JB Clement, 2023-2025
556
557! NOTES
558!
559!-----------------------------------------------------------------------
560
561! DEPENDENCIES
562! ------------
563  USE netcdf
564  USE dimphy
565  USE mod_grid_phy_lmdz
566  USE mod_phys_lmdz_para
567
568! DECLARATION
569! -----------
570implicit none
571
572! ARGUMENTS
573! ---------
574    CHARACTER(LEN=*), INTENT(IN)            :: field_name
575    INTEGER,          INTENT(IN)            :: field_size
576    REAL,             INTENT(OUT)           :: field(klon,field_size)
577    INTEGER,          INTENT(IN)            :: corners(4)
578    INTEGER,          INTENT(IN)            :: edges(4)
579    LOGICAL,          INTENT(OUT), OPTIONAL :: found
580
581! LOCAL VARIABLES
582! ---------------
583    REAL    :: field_glo(klon_glo,field_size)
584    LOGICAL :: tmp_found
585    INTEGER :: varid
586    INTEGER :: ierr
587
588! CODE
589! ----
590    IF (is_master) THEN
591
592      ierr=NF90_INQ_VARID(nid_start,field_name,varid)
593
594      IF (ierr==NF90_NOERR) THEN
595        CALL body(field_glo)
596        tmp_found=.TRUE.
597      ELSE
598        tmp_found=.FALSE.
599      ENDIF
600
601    ENDIF
602
603    CALL bcast(tmp_found)
604
605    IF (tmp_found) THEN
606      CALL scatter(field_glo,field)
607    ENDIF
608
609    IF (PRESENT(found)) THEN
610      found=tmp_found
611    ELSE
612      IF (.NOT. tmp_found) THEN
613        write(*,*) 'get_field_rgen: Field <'//field_name//'> not found'
614        CALL abort_physic("get_field_rgen","Failed to get field",1)
615      ENDIF
616    ENDIF
617
618
619        CONTAINS
620    !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
621
622    !=======================================================================
623      SUBROUTINE body(field_glo)
624    !-----------------------------------------------------------------------
625    ! NAME
626    !     body
627    !
628    ! DESCRIPTION
629    !     Read the requested field from the start file into the global array.
630    !
631    ! AUTHORS & DATE
632    !     L. Lange
633    !     JB Clement, 2023-2025
634    !
635    ! NOTES
636    !     Helper contained in Get_field_rgen.
637    !-----------------------------------------------------------------------
638
639    ! DECLARATION
640    ! -----------
641      ! DECLARATION
642! -----------
643implicit none
644
645    ! ARGUMENTS
646    ! ---------
647      REAL, INTENT(OUT) :: field_glo(klon_glo*field_size) ! Buffer for global field
648
649    ! CODE
650    ! ----
651          ierr=NF90_GET_VAR(nid_start,varid,field_glo,corners,edges)
652          IF (ierr/=NF90_NOERR) THEN
653         ! La variable exist dans le fichier mais la lecture a echouee.
654         write(*,*) 'get_field_rgen: Failed reading <'//field_name//'>'
655
656!           IF (field_name=='CLWCON' .OR. field_name=='RNEBCON' .OR. field_name=='RATQS') THEN
657!              ! Essaye de lire le variable sur surface uniqument, comme fait avant
658!              field_glo(:)=0.
659!              ierr=NF90_GET_VAR(nid_start,varid,field_glo(1:klon_glo))
660!              IF (ierr/=NF90_NOERR) THEN
661!                 write(*,*) 'phyetat0: Lecture echouee aussi en 2D pour <'//field_name//'>'
662!                 CALL stop_clean
663!              ELSE
664!                 write(*,*) 'phyetat0: La variable <'//field_name//'> lu sur surface seulement'!, selon ancien format, le reste mis a zero'
665!              END IF
666!           ELSE
667              CALL abort_physic("get_field_rgen","Failed to read field",1)
668!           ENDIF
669         ENDIF
670
671     END SUBROUTINE body
672
673  END SUBROUTINE Get_field_rgen
674!=======================================================================
675
676!=======================================================================
677  SUBROUTINE get_var_r0(var_name,var,found)
678!-----------------------------------------------------------------------
679! NAME
680!     get_var_r0
681!
682! DESCRIPTION
683!     Get a scalar variable from input file.
684!
685! AUTHORS & DATE
686!     L. Lange
687!     JB Clement, 2023-2025
688!
689! NOTES
690!
691!-----------------------------------------------------------------------
692
693! DECLARATION
694! -----------
695implicit none
696
697! ARGUMENTS
698! ---------
699  CHARACTER(LEN=*), INTENT(IN)            :: var_name
700  REAL,             INTENT(INOUT)         :: var
701  LOGICAL,          INTENT(OUT), OPTIONAL :: found
702
703! LOCAL VARIABLES
704! ---------------
705    REAL :: varout(1)
706
707! CODE
708! ----
709    IF (PRESENT(found)) THEN
710      CALL Get_var_rgen(var_name,varout,size(varout),found)
711    ELSE
712      CALL Get_var_rgen(var_name,varout,size(varout))
713    ENDIF
714    var=varout(1)
715
716  END SUBROUTINE get_var_r0
717!=======================================================================
718
719!=======================================================================
720  SUBROUTINE get_var_r1(var_name,var,found)
721!-----------------------------------------------------------------------
722! NAME
723!     get_var_r1
724!
725! DESCRIPTION
726!     Get a vector (1D) variable from input file.
727!
728! AUTHORS & DATE
729!     L. Lange
730!     JB Clement, 2023-2025
731!
732! NOTES
733!
734!-----------------------------------------------------------------------
735
736! DECLARATION
737! -----------
738implicit none
739
740! ARGUMENTS
741! ---------
742  CHARACTER(LEN=*), INTENT(IN)            :: var_name
743  REAL,             INTENT(INOUT)         :: var(:)
744  LOGICAL,          INTENT(OUT), OPTIONAL :: found
745
746! CODE
747! ----
748    IF (PRESENT(found)) THEN
749      CALL Get_var_rgen(var_name,var,size(var),found)
750    ELSE
751      CALL Get_var_rgen(var_name,var,size(var))
752    ENDIF
753
754  END SUBROUTINE get_var_r1
755!=======================================================================
756
757!=======================================================================
758  SUBROUTINE get_var_r2(var_name,var,found)
759!-----------------------------------------------------------------------
760! NAME
761!     get_var_r2
762!
763! DESCRIPTION
764!     Get a 2D field variable from input file.
765!
766! AUTHORS & DATE
767!     L. Lange
768!     JB Clement, 2023-2025
769!
770! NOTES
771!
772!-----------------------------------------------------------------------
773
774! DECLARATION
775! -----------
776implicit none
777
778! ARGUMENTS
779! ---------
780  CHARACTER(LEN=*), INTENT(IN)            :: var_name
781  REAL,             INTENT(OUT)           :: var(:,:)
782  LOGICAL,          INTENT(OUT), OPTIONAL :: found
783
784! CODE
785! ----
786    IF (PRESENT(found)) THEN
787      CALL Get_var_rgen(var_name,var,size(var),found)
788    ELSE
789      CALL Get_var_rgen(var_name,var,size(var))
790    ENDIF
791
792  END SUBROUTINE get_var_r2
793!=======================================================================
794
795!=======================================================================
796  SUBROUTINE get_var_r3(var_name,var,found)
797!-----------------------------------------------------------------------
798! NAME
799!     get_var_r3
800!
801! DESCRIPTION
802!     Get a 3D field variable from input file.
803!
804! AUTHORS & DATE
805!     L. Lange
806!     JB Clement, 2023-2025
807!
808! NOTES
809!
810!-----------------------------------------------------------------------
811
812! DECLARATION
813! -----------
814implicit none
815
816! ARGUMENTS
817! ---------
818  CHARACTER(LEN=*), INTENT(IN)            :: var_name
819  REAL,             INTENT(INOUT)         :: var(:,:,:)
820  LOGICAL,          INTENT(OUT), OPTIONAL :: found
821
822! CODE
823! ----
824    IF (PRESENT(found)) THEN
825      CALL Get_var_rgen(var_name,var,size(var),found)
826    ELSE
827      CALL Get_var_rgen(var_name,var,size(var))
828    ENDIF
829
830  END SUBROUTINE get_var_r3
831!=======================================================================
832
833!=======================================================================
834  SUBROUTINE Get_var_rgen(var_name,var,var_size,found)
835!-----------------------------------------------------------------------
836! NAME
837!     Get_var_rgen
838!
839! DESCRIPTION
840!     Generic subroutine to read variables from input file.
841!
842! AUTHORS & DATE
843!     L. Lange
844!     JB Clement, 2023-2025
845!
846! NOTES
847!
848!-----------------------------------------------------------------------
849
850! DEPENDENCIES
851! ------------
852  USE netcdf
853  USE dimphy
854  USE mod_grid_phy_lmdz
855  USE mod_phys_lmdz_para
856
857! DECLARATION
858! -----------
859implicit none
860
861! ARGUMENTS
862! ---------
863    CHARACTER(LEN=*), INTENT(IN)            :: var_name
864    INTEGER,          INTENT(IN)            :: var_size
865    REAL,             INTENT(OUT)           :: var(var_size)
866    LOGICAL,          INTENT(OUT), OPTIONAL :: found
867
868! LOCAL VARIABLES
869! ---------------
870    LOGICAL :: tmp_found
871    INTEGER :: varid
872    INTEGER :: ierr
873
874! CODE
875! ----
876    IF (is_mpi_root .AND. is_omp_root) THEN
877
878      ierr=NF90_INQ_VARID(nid_start,var_name,varid)
879
880      IF (ierr==NF90_NOERR) THEN
881        ierr=NF90_GET_VAR(nid_start,varid,var)
882        IF (ierr/=NF90_NOERR) THEN
883          write(*,*) 'phyetat0: Failed loading <'//trim(var_name)//'>'
884          CALL abort_physic("get_var_rgen","Failed to read variable",1)
885        ENDIF
886        tmp_found=.TRUE.
887      ELSE
888        tmp_found=.FALSE.
889      ENDIF
890
891    ENDIF
892
893    CALL bcast(tmp_found)
894
895    IF (tmp_found) THEN
896      CALL bcast(var)
897    ENDIF
898
899    IF (PRESENT(found)) THEN
900      found=tmp_found
901    ELSE
902      IF (.NOT. tmp_found) THEN
903        write(*,*) 'phyetat0: Variable <'//trim(var_name)//'> not found'
904        CALL abort_physic("get_var_rgen","Failed to read variable",1)
905      ENDIF
906    ENDIF
907
908  END SUBROUTINE Get_var_rgen
909!=======================================================================
910
911!=======================================================================
912  SUBROUTINE open_restartphy(filename)
913!-----------------------------------------------------------------------
914! NAME
915!     open_restartphy
916!
917! DESCRIPTION
918!     Open (or create) the restart netCDF file for writing.
919!     Define dimensions and global attributes on first call.
920!
921! AUTHORS & DATE
922!     L. Lange
923!     JB Clement, 2023-2025
924!
925! NOTES
926!
927!-----------------------------------------------------------------------
928
929! DEPENDENCIES
930! ------------
931  USE netcdf, only: NF90_CREATE, NF90_CLOBBER, NF90_64BIT_OFFSET, &
932                    NF90_NOERR, nf90_strerror, &
933                    NF90_PUT_ATT, NF90_GLOBAL, NF90_DEF_DIM, &
934                    NF90_UNLIMITED, NF90_ENDDEF, &
935                    NF90_WRITE, NF90_OPEN
936  USE mod_phys_lmdz_para, only: is_master
937  USE mod_grid_phy_lmdz,  only: klon_glo
938  USE dimphy,             only: klev, klevp1
939  USE tracer_mod,         only: nqmx
940  USE soil,               only:  nsoilmx_PEM
941  USE comslope_mod,       only: nslope
942  use layered_deposits,   only: nb_str_max
943
944! DECLARATION
945! -----------
946implicit none
947
948! ARGUMENTS
949! ---------
950  CHARACTER(LEN=*), INTENT(IN) :: filename
951
952! LOCAL VARIABLES
953! ---------------
954    INTEGER           :: ierr
955    LOGICAL,SAVE      :: already_created=.false.
956
957! CODE
958! ----
959    IF (is_master) THEN
960      if (.not.already_created) then
961        ! At the very first call, create the file
962        ierr=NF90_CREATE(filename,IOR(NF90_CLOBBER,NF90_64BIT_OFFSET), &
963                          nid_restart)
964        IF (ierr/=NF90_NOERR) THEN
965          write(*,*)'open_restartphy: problem creating file '//trim(filename)
966          write(*,*)trim(nf90_strerror(ierr))
967          CALL abort_physic("open_restartphy","Failed creating file",1)
968        ENDIF
969        already_created=.true.
970      else
971        ! Just open the file
972        ierr=NF90_OPEN(filename,NF90_WRITE,nid_restart)
973        IF (ierr/=NF90_NOERR) THEN
974          write(*,*)'open_restartphy: problem opening file '//trim(filename)
975          write(*,*)trim(nf90_strerror(ierr))
976          CALL abort_physic("open_restartphy","Failed opening file",1)
977        ENDIF
978        return
979      endif ! of if (.not.already_created)
980
981      ierr=NF90_PUT_ATT(nid_restart,NF90_GLOBAL,"title",&
982                        "Physics start file")
983      IF (ierr/=NF90_NOERR) THEN
984        write(*,*)'open_restartphy: problem writing title '
985        write(*,*)trim(nf90_strerror(ierr))
986      ENDIF
987
988      ierr=NF90_DEF_DIM(nid_restart,"index",length,idim1)
989      IF (ierr/=NF90_NOERR) THEN
990        write(*,*)'open_restartphy: problem defining index dimension '
991        write(*,*)trim(nf90_strerror(ierr))
992        CALL abort_physic("open_restartphy","Failed defining index",1)
993      ENDIF
994
995      ierr=NF90_DEF_DIM(nid_restart,"physical_points",klon_glo,idim2)
996      IF (ierr/=NF90_NOERR) THEN
997        write(*,*)'open_restartphy: problem defining physical_points dimension '
998        write(*,*)trim(nf90_strerror(ierr))
999        CALL abort_physic("open_restartphy","Failed defining physical_points",1)
1000      ENDIF
1001
1002      ierr=NF90_DEF_DIM(nid_restart,"subsurface_layers",nsoilmx_PEM,idim3)
1003      IF (ierr/=NF90_NOERR) THEN
1004        write(*,*)'open_restartphy: problem defining subsurface_layers dimension '
1005        write(*,*)trim(nf90_strerror(ierr))
1006        CALL abort_physic("open_restartphy","Failed defining subsurface_layers",1)
1007      ENDIF
1008
1009      ierr=NF90_DEF_DIM(nid_restart,"nlayer_plus_1",klevp1,idim4)
1010      IF (ierr/=NF90_NOERR) THEN
1011        write(*,*)'open_restartphy: problem defining nlayer_plus_1 dimension '
1012        write(*,*)trim(nf90_strerror(ierr))
1013        CALL abort_physic("open_restartphy","Failed defining nlayer_plus_1",1)
1014      ENDIF
1015
1016      if (nqmx.ne.0) then
1017        ierr=NF90_DEF_DIM(nid_restart,"number_of_advected_fields",nqmx,idim5)
1018      else
1019        ! pretend nqmx=1 because 0 size implies unlimited dimension for NetCDF
1020        ierr=NF90_DEF_DIM(nid_restart,"number_of_advected_fields",1,idim5)
1021      endif
1022      IF (ierr/=NF90_NOERR) THEN
1023        write(*,*)'open_restartphy: problem defining number_of_advected_fields dimension '
1024        write(*,*)trim(nf90_strerror(ierr))
1025        CALL abort_physic("open_restartphy","Failed defining number_of_advected_fields",1)
1026      ENDIF
1027
1028      ierr=NF90_DEF_DIM(nid_restart,"nlayer",klev,idim6)
1029      IF (ierr/=NF90_NOERR) THEN
1030        write(*,*)'open_restartphy: problem defining nlayer dimension '
1031        write(*,*)trim(nf90_strerror(ierr))
1032        CALL abort_physic("open_restartphy","Failed defining nlayer",1)
1033      ENDIF
1034
1035      ierr=NF90_DEF_DIM(nid_restart,"Time",NF90_UNLIMITED,idim7)
1036      IF (ierr/=NF90_NOERR) THEN
1037        write(*,*)'open_restartphy: problem defining Time dimension '
1038        write(*,*)trim(nf90_strerror(ierr))
1039        CALL abort_physic("open_restartphy","Failed defining Time",1)
1040      ENDIF
1041
1042      ierr=NF90_DEF_DIM(nid_restart,"nslope",nslope,idim8)
1043      IF (ierr/=NF90_NOERR) THEN
1044        write(*,*)'phyredem: problem defining nslope dimension'
1045        write(*,*)trim(nf90_strerror(ierr))
1046        CALL ABORT
1047      ENDIF
1048
1049      ierr=NF90_DEF_DIM(nid_restart,"inter slope",nslope+1,idim9)
1050      IF (ierr/=NF90_NOERR) THEN
1051        write(*,*)'phyredem: problem defining inter slope dimension'
1052        write(*,*)trim(nf90_strerror(ierr))
1053        CALL ABORT
1054      ENDIF
1055
1056      ierr=NF90_DEF_DIM(nid_restart,"nb_str_max",nb_str_max,idim10)
1057      IF (ierr/=NF90_NOERR) THEN
1058        write(*,*)'phyredem: problem defining nb_str_max dimension'
1059        write(*,*)trim(nf90_strerror(ierr))
1060        CALL ABORT
1061      ENDIF
1062
1063      ierr=NF90_ENDDEF(nid_restart)
1064      IF (ierr/=NF90_NOERR) THEN
1065        write(*,*)'open_restartphy: problem ending definition mode '
1066        write(*,*)trim(nf90_strerror(ierr))
1067        CALL abort_physic("open_restartphy","Failed ending definition mode",1)
1068      ENDIF
1069    ENDIF
1070
1071  END SUBROUTINE open_restartphy
1072!=======================================================================
1073
1074!=======================================================================
1075  SUBROUTINE close_restartphy
1076!-----------------------------------------------------------------------
1077! NAME
1078!     close_restartphy
1079!
1080! DESCRIPTION
1081!     Close the restart netCDF file.
1082!
1083! AUTHORS & DATE
1084!     L. Lange
1085!     JB Clement, 2023-2025
1086!
1087! NOTES
1088!
1089!-----------------------------------------------------------------------
1090
1091! DEPENDENCIES
1092! ------------
1093  USE netcdf,             only: NF90_CLOSE
1094  USE mod_phys_lmdz_para, only: is_master
1095
1096! DECLARATION
1097! -----------
1098implicit none
1099
1100! LOCAL VARIABLES
1101! ---------------
1102    INTEGER :: ierr
1103
1104! CODE
1105! ----
1106    IF (is_master) THEN
1107      ierr = NF90_CLOSE (nid_restart)
1108    ENDIF
1109
1110  END SUBROUTINE close_restartphy
1111!=======================================================================
1112
1113!=======================================================================
1114  SUBROUTINE put_field_r1(field_name,title,field,time)
1115!-----------------------------------------------------------------------
1116! NAME
1117!     put_field_r1
1118!
1119! DESCRIPTION
1120!     Write a surface field to the restart file.
1121!
1122! AUTHORS & DATE
1123!     L. Lange
1124!     JB Clement, 2023-2025
1125!
1126! NOTES
1127!
1128!-----------------------------------------------------------------------
1129
1130! DECLARATION
1131! -----------
1132implicit none
1133
1134
1135! ARGUMENTS
1136! ---------
1137  CHARACTER(LEN=*), INTENT(IN)           :: field_name
1138  CHARACTER(LEN=*), INTENT(IN)           :: title
1139  REAL,             INTENT(IN)           :: field(:)
1140  REAL,             INTENT(IN), OPTIONAL :: time
1141
1142! CODE
1143! ----
1144  IF (present(time)) THEN
1145    ! if timeindex is present, it is a time-dependent variable
1146    CALL put_field_rgen(field_name,title,field,1,time)
1147  ELSE
1148    CALL put_field_rgen(field_name,title,field,1)
1149  ENDIF
1150
1151  END SUBROUTINE put_field_r1
1152!=======================================================================
1153
1154!=======================================================================
1155  SUBROUTINE put_field_r2(field_name,title,field,time)
1156!-----------------------------------------------------------------------
1157! NAME
1158!     put_field_r2
1159!
1160! DESCRIPTION
1161!     Write a "3D" horizontal-vertical field (2D) to the restart file.
1162
1163! AUTHORS & DATE
1164!     L. Lange
1165!     JB Clement, 2023-2025
1166
1167! NOTES
1168!
1169!-----------------------------------------------------------------------
1170
1171! DECLARATION
1172! -----------
1173implicit none
1174
1175! ARGUMENTS
1176! ---------
1177  CHARACTER(LEN=*), INTENT(IN)           :: field_name
1178  CHARACTER(LEN=*), INTENT(IN)           :: title
1179  REAL,             INTENT(IN)           :: field(:,:)
1180  REAL,             INTENT(IN), OPTIONAL :: time
1181
1182! CODE
1183! ----
1184  IF (present(time)) THEN
1185    ! if timeindex is present, it is a time-dependent variable
1186    CALL put_field_rgen(field_name,title,field,size(field,2),time)
1187  ELSE
1188    CALL put_field_rgen(field_name,title,field,size(field,2))
1189  ENDIF
1190
1191  END SUBROUTINE put_field_r2
1192!=======================================================================
1193
1194!=======================================================================
1195  SUBROUTINE put_field_r3(field_name,title,field,time)
1196!-----------------------------------------------------------------------
1197! NAME
1198!     put_field_r3
1199!
1200! DESCRIPTION
1201!     Write a "4D" field surf/alt/?? (3D) to the restart file.
1202
1203! AUTHORS & DATE
1204!     L. Lange
1205!     JB Clement, 2023-2025
1206
1207! NOTES
1208!
1209!-----------------------------------------------------------------------
1210
1211! DECLARATION
1212! -----------
1213implicit none
1214
1215! ARGUMENTS
1216! ---------
1217  CHARACTER(LEN=*), INTENT(IN)           :: field_name
1218  CHARACTER(LEN=*), INTENT(IN)           :: title
1219  REAL,             INTENT(IN)           :: field(:,:,:)
1220  REAL,             INTENT(IN), OPTIONAL :: time
1221
1222! CODE
1223! ----
1224  IF (present(time)) THEN
1225    ! if timeindex is present, it is a time-dependent variable
1226    CALL put_field_rgen(field_name,title,field,size(field,2)*size(field,3),&
1227                        time)
1228  ELSE
1229    CALL put_field_rgen(field_name,title,field,size(field,2)*size(field,3))
1230  ENDIF
1231
1232  END SUBROUTINE put_field_r3
1233!=======================================================================
1234
1235!=======================================================================
1236  SUBROUTINE put_field_rgen(field_name,title,field,field_size,time)
1237!-----------------------------------------------------------------------
1238! NAME
1239!     put_field_rgen
1240!
1241! DESCRIPTION
1242!     Generic subroutine to write fields to restart file with gather from processes
1243
1244! AUTHORS & DATE
1245!     L. Lange
1246!     JB Clement, 2023-2025
1247
1248! NOTES
1249!
1250!-----------------------------------------------------------------------
1251
1252! DEPENDENCIES
1253! ------------
1254  USE netcdf
1255  USE dimphy
1256  USE soil,         only: nsoilmx_PEM
1257  USE comslope_mod, ONLY: nslope
1258  USE mod_grid_phy_lmdz
1259  USE mod_phys_lmdz_para
1260
1261! DECLARATION
1262! -----------
1263implicit none
1264
1265! ARGUMENTS
1266! ---------
1267  CHARACTER(LEN=*),INTENT(IN)    :: field_name
1268  CHARACTER(LEN=*),INTENT(IN)    :: title
1269  INTEGER,INTENT(IN)             :: field_size
1270  REAL,INTENT(IN)                :: field(klon,field_size)
1271  REAL,OPTIONAL,INTENT(IN)       :: time
1272
1273! LOCAL VARIABLES
1274! ---------------
1275  REAL                           :: field_glo(klon_glo,field_size)
1276  INTEGER                        :: ierr
1277  INTEGER                        :: nvarid
1278
1279! CODE
1280! ----
1281    CALL gather(field,field_glo)
1282
1283    IF (is_master) THEN
1284
1285     IF (index(field_name,"stratif") == 0) then
1286      IF (field_size==1) THEN
1287        ! input is a 1D "surface field" array
1288        if (.not.present(time)) then ! for a time-independent field
1289          ierr=NF90_REDEF(nid_restart)
1290#ifdef NC_DOUBLE
1291          ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_DOUBLE,&
1292                            (/idim2/),nvarid)
1293#else
1294          ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_FLOAT,&
1295                            (/idim2/),nvarid)
1296#endif
1297          if (ierr.ne.NF90_NOERR) then
1298            write(*,*)"put_field_rgen error: failed to define "//trim(field_name)
1299            write(*,*)trim(nf90_strerror(ierr))
1300          endif
1301          IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title)
1302          ierr=NF90_ENDDEF(nid_restart)
1303          ierr=NF90_PUT_VAR(nid_restart,nvarid,field_glo)
1304        else
1305          ! check if the variable has already been defined:
1306          ierr=NF90_INQ_VARID(nid_restart,field_name,nvarid)
1307          if (ierr/=NF90_NOERR) then ! variable not found, define it
1308            ierr=NF90_REDEF(nid_restart)
1309#ifdef NC_DOUBLE
1310            ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_DOUBLE,&
1311                              (/idim2,idim7/),nvarid)
1312#else
1313            ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_FLOAT,&
1314                              (/idim2,idim7/),nvarid)
1315#endif
1316            if (ierr.ne.NF90_NOERR) then
1317              write(*,*)"put_field_rgen error: failed to define "//trim(field_name)
1318              write(*,*)trim(nf90_strerror(ierr))
1319            endif
1320            IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title)
1321            ierr=NF90_ENDDEF(nid_restart)
1322          endif
1323          ! Write the variable
1324          ierr=NF90_PUT_VAR(nid_restart,nvarid,field_glo,&
1325                            start=(/1,timeindex/))
1326        endif ! of if (.not.present(timeindex))
1327
1328      ELSE IF (field_size==klev) THEN
1329        ! input is a 2D "atmospheric field" array
1330        if (.not.present(time)) then ! for a time-independent field
1331          ierr=NF90_REDEF(nid_restart)
1332#ifdef NC_DOUBLE
1333          ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_DOUBLE,&
1334                            (/idim2,idim6/),nvarid)
1335#else
1336          ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_FLOAT,&
1337                            (/idim2,idim6/),nvarid)
1338#endif
1339          if (ierr.ne.NF90_NOERR) then
1340            write(*,*)"put_field_rgen error: failed to define "//trim(field_name)
1341            write(*,*)trim(nf90_strerror(ierr))
1342          endif
1343          IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title)
1344          ierr=NF90_ENDDEF(nid_restart)
1345          ierr=NF90_PUT_VAR(nid_restart,nvarid,field_glo)
1346        else
1347          ! check if the variable has already been defined:
1348          ierr=NF90_INQ_VARID(nid_restart,field_name,nvarid)
1349          if (ierr/=NF90_NOERR) then ! variable not found, define it
1350            ierr=NF90_REDEF(nid_restart)
1351#ifdef NC_DOUBLE
1352            ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_DOUBLE,&
1353                              (/idim2,idim6,idim7/),nvarid)
1354#else
1355            ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_FLOAT,&
1356                              (/idim2,idim6,idim7/),nvarid)
1357#endif
1358            if (ierr.ne.NF90_NOERR) then
1359              write(*,*)"put_field_rgen error: failed to define "//trim(field_name)
1360              write(*,*)trim(nf90_strerror(ierr))
1361            endif
1362            IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title)
1363            ierr=NF90_ENDDEF(nid_restart)
1364          endif
1365          ! Write the variable
1366          ierr=NF90_PUT_VAR(nid_restart,nvarid,field_glo,&
1367                            start=(/1,1,timeindex/))
1368        endif ! of if (.not.present(time))
1369
1370      ELSE IF (field_size==klevp1) THEN
1371        ! input is a 2D "interlayer atmospheric field" array
1372        if (.not.present(time)) then ! for a time-independent field
1373          ierr=NF90_REDEF(nid_restart)
1374#ifdef NC_DOUBLE
1375          ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_DOUBLE,&
1376                            (/idim2,idim4/),nvarid)
1377#else
1378          ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_FLOAT,&
1379                            (/idim2,idim4/),nvarid)
1380#endif
1381          if (ierr.ne.NF90_NOERR) then
1382            write(*,*)"put_field_rgen error: failed to define "//trim(field_name)
1383            write(*,*)trim(nf90_strerror(ierr))
1384          endif
1385          IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title)
1386          ierr = NF90_ENDDEF(nid_restart)
1387          ierr = NF90_PUT_VAR(nid_restart,nvarid,field_glo)
1388        else
1389          ! check if the variable has already been defined:
1390          ierr=NF90_INQ_VARID(nid_restart,field_name,nvarid)
1391          if (ierr/=NF90_NOERR) then ! variable not found, define it
1392            ierr=NF90_REDEF(nid_restart)
1393#ifdef NC_DOUBLE
1394            ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_DOUBLE,&
1395                              (/idim2,idim4,idim7/),nvarid)
1396#else
1397            ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_FLOAT,&
1398                              (/idim2,idim4,idim7/),nvarid)
1399#endif
1400            if (ierr.ne.NF90_NOERR) then
1401              write(*,*)"put_field_rgen error: failed to define "//trim(field_name)
1402              write(*,*)trim(nf90_strerror(ierr))
1403            endif
1404            IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title)
1405            ierr=NF90_ENDDEF(nid_restart)
1406          endif
1407          ! Write the variable
1408          ierr=NF90_PUT_VAR(nid_restart,nvarid,field_glo,&
1409                            start=(/1,1,timeindex/))
1410        endif ! of if (.not.present(timeindex))
1411
1412      ELSE IF (field_size==nsoilmx_PEM) THEN
1413        ! input is a 2D "subsurface field" array
1414        if (.not.present(time)) then ! for a time-independent field
1415          ierr = NF90_REDEF(nid_restart)
1416#ifdef NC_DOUBLE
1417          ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_DOUBLE,&
1418                            (/idim2,idim3/),nvarid)
1419#else
1420          ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_FLOAT,&
1421                            (/idim2,idim3/),nvarid)
1422#endif
1423          if (ierr.ne.NF90_NOERR) then
1424            write(*,*)"put_field_rgen error: failed to define "//trim(field_name)
1425            write(*,*)trim(nf90_strerror(ierr))
1426          endif
1427          IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title)
1428          ierr = NF90_ENDDEF(nid_restart)
1429          ierr = NF90_PUT_VAR(nid_restart,nvarid,field_glo)
1430        else
1431          ! check if the variable has already been defined:
1432          ierr=NF90_INQ_VARID(nid_restart,field_name,nvarid)
1433          if (ierr/=NF90_NOERR) then ! variable not found, define it
1434            ierr=NF90_REDEF(nid_restart)
1435#ifdef NC_DOUBLE
1436            ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_DOUBLE,&
1437                              (/idim2,idim3,idim7/),nvarid)
1438#else
1439            ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_FLOAT,&
1440                              (/idim2,idim3,idim7/),nvarid)
1441#endif
1442           if (ierr.ne.NF90_NOERR) then
1443              write(*,*)"put_field_rgen error: failed to define "//trim(field_name)
1444              write(*,*)trim(nf90_strerror(ierr))
1445            endif
1446            IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title)
1447            ierr=NF90_ENDDEF(nid_restart)
1448          endif
1449          ! Write the variable
1450          ierr=NF90_PUT_VAR(nid_restart,nvarid,field_glo,&
1451                            start=(/1,1,timeindex/))
1452
1453        endif ! of if (.not.present(time))
1454
1455
1456      ELSE IF (field_size==nslope) THEN
1457        ! input is a 2D "subsurface field" array
1458        if (.not.present(time)) then ! for a time-independent field
1459          ierr = NF90_REDEF(nid_restart)
1460#ifdef NC_DOUBLE
1461          ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_DOUBLE,&
1462                            (/idim2,idim8/),nvarid)
1463#else
1464          ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_FLOAT,&
1465                            (/idim2,idim8/),nvarid)
1466#endif
1467          if (ierr.ne.NF90_NOERR) then
1468            write(*,*)"put_field_rgen error: failed to define"//trim(field_name)
1469            write(*,*)trim(nf90_strerror(ierr))
1470          endif
1471          IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title)
1472          ierr = NF90_ENDDEF(nid_restart)
1473          ierr = NF90_PUT_VAR(nid_restart,nvarid,field_glo)
1474        else
1475          ! check if the variable has already been defined:
1476          ierr=NF90_INQ_VARID(nid_restart,field_name,nvarid)
1477          if (ierr/=NF90_NOERR) then ! variable not found, define it
1478            ierr=NF90_REDEF(nid_restart)
1479#ifdef NC_DOUBLE
1480            ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_DOUBLE,&
1481                              (/idim2,idim8,idim7/),nvarid)
1482#else
1483            ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_FLOAT,&
1484                              (/idim2,idim8,idim7/),nvarid)
1485#endif
1486           if (ierr.ne.NF90_NOERR) then
1487              write(*,*)"put_field_rgen error: failed to define"//trim(field_name)
1488              write(*,*)trim(nf90_strerror(ierr))
1489            endif
1490            IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title)
1491            ierr=NF90_ENDDEF(nid_restart)
1492          endif
1493          ! Write the variable
1494          ierr=NF90_PUT_VAR(nid_restart,nvarid,field_glo,&
1495                            start=(/1,1,timeindex/))
1496
1497        endif ! of if (.not.present(time))
1498
1499      ELSE
1500        write(*,*) "Error phyredem(put_field_rgen) : wrong dimension for ",trim(field_name)
1501        write(*,*) "  field_size =",field_size
1502        CALL abort_physic("put_field_rgen","wrong field dimension",1)
1503      ENDIF
1504     ELSE
1505        ! input is a 2D "stratification" array
1506        if (.not.present(time)) then ! for a time-independent field
1507          ierr = NF90_REDEF(nid_restart)
1508#ifdef NC_DOUBLE
1509          ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_DOUBLE,&
1510                            (/idim2,idim8/),nvarid)
1511#else
1512          ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_FLOAT,&
1513                            (/idim2,idim8/),nvarid)
1514#endif
1515          if (ierr.ne.NF90_NOERR) then
1516            write(*,*)"put_field_rgen error: failed to define"//trim(field_name)
1517            write(*,*)trim(nf90_strerror(ierr))
1518          endif
1519          IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title)
1520          ierr = NF90_ENDDEF(nid_restart)
1521          ierr = NF90_PUT_VAR(nid_restart,nvarid,field_glo)
1522        else
1523          ! check if the variable has already been defined:
1524          ierr=NF90_INQ_VARID(nid_restart,field_name,nvarid)
1525          if (ierr/=NF90_NOERR) then ! variable not found, define it
1526            ierr=NF90_REDEF(nid_restart)
1527#ifdef NC_DOUBLE
1528            ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_DOUBLE,&
1529                              (/idim2,idim10,idim7/),nvarid)
1530#else
1531            ierr=NF90_DEF_VAR(nid_restart,field_name,NF90_FLOAT,&
1532                              (/idim2,idim10,idim7/),nvarid)
1533#endif
1534           if (ierr.ne.NF90_NOERR) then
1535              write(*,*)"put_field_rgen error: failed to define"//trim(field_name)
1536              write(*,*)trim(nf90_strerror(ierr))
1537            endif
1538            IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title)
1539            ierr=NF90_ENDDEF(nid_restart)
1540          endif
1541          ! Write the variable
1542          ierr=NF90_PUT_VAR(nid_restart,nvarid,field_glo,&
1543                            start=(/1,1,timeindex/))
1544
1545        endif ! of if (.not.present(time))
1546     ENDIF
1547
1548      ! Check the writting of field to file went OK
1549      if (ierr.ne.NF90_NOERR) then
1550        write(*,*) " Error phyredem(put_field_rgen) : failed writing ",trim(field_name)
1551        write(*,*)trim(nf90_strerror(ierr))
1552        CALL abort_physic("put_field_rgen","Failed writing field",1)
1553      endif
1554
1555    ENDIF ! of IF (is_master)
1556
1557  END SUBROUTINE put_field_rgen
1558!=======================================================================
1559
1560!=======================================================================
1561  SUBROUTINE put_var_r0(var_name,title,var)
1562!-----------------------------------------------------------------------
1563! NAME
1564!     put_var_r0
1565!
1566! DESCRIPTION
1567!     Write a scalar variable to the restart file
1568
1569! AUTHORS & DATE
1570!     L. Lange
1571!     JB Clement, 2023-2025
1572
1573! NOTES
1574!
1575!-----------------------------------------------------------------------
1576
1577! DECLARATION
1578! -----------
1579implicit none
1580
1581! ARGUMENTS
1582! ---------
1583  CHARACTER(LEN=*), INTENT(IN) :: var_name
1584  CHARACTER(LEN=*), INTENT(IN) :: title
1585  REAL,             INTENT(IN) :: var
1586
1587! LOCAL VARIABLES
1588! ---------------
1589     REAL :: varin(1)
1590
1591! CODE
1592! ----
1593     varin(1)=var
1594
1595     CALL put_var_rgen(var_name,title,varin,size(varin))
1596
1597  END SUBROUTINE put_var_r0
1598!=======================================================================
1599
1600!=======================================================================
1601  SUBROUTINE put_var_r1(var_name,title,var)
1602!-----------------------------------------------------------------------
1603! NAME
1604!    put_var_r1
1605!
1606! DESCRIPTION
1607!     Write a vector (1D) variable to the restart file
1608
1609! AUTHORS & DATE
1610!     L. Lange
1611!     JB Clement, 2023-2025
1612
1613! NOTES
1614!
1615!-----------------------------------------------------------------------
1616
1617! DECLARATION
1618! -----------
1619implicit none
1620
1621! ARGUMENTS
1622! ---------
1623  CHARACTER(LEN=*), INTENT(IN) :: var_name
1624  CHARACTER(LEN=*), INTENT(IN) :: title
1625  REAL,             INTENT(IN) :: var(:)
1626
1627! CODE
1628! ----
1629  CALL put_var_rgen(var_name,title,var,size(var))
1630
1631  END SUBROUTINE put_var_r1
1632!=======================================================================
1633
1634!=======================================================================
1635  SUBROUTINE put_var_r2(var_name,title,var)
1636!-----------------------------------------------------------------------
1637! NAME
1638!     put_var_r2
1639!
1640! DESCRIPTION
1641!     Write a 3D field variable to the restart file
1642
1643! AUTHORS & DATE
1644!     L. Lange
1645!     JB Clement, 2023-2025
1646
1647! NOTES
1648!
1649!-----------------------------------------------------------------------
1650
1651! DECLARATION
1652! -----------
1653implicit none
1654
1655! ARGUMENTS
1656! ---------
1657  CHARACTER(LEN=*), INTENT(IN) :: var_name
1658  CHARACTER(LEN=*), INTENT(IN) :: title
1659  REAL,             INTENT(IN) :: var(:,:)
1660
1661! CODE
1662! ----
1663  CALL put_var_rgen(var_name,title,var,size(var))
1664
1665  END SUBROUTINE put_var_r2
1666!=======================================================================
1667
1668!=======================================================================
1669  SUBROUTINE put_var_r3(var_name,title,var)
1670!-----------------------------------------------------------------------
1671! NAME
1672!     put_var_r3
1673!
1674! DESCRIPTION
1675!     Write a 3D field variable to the restart file
1676!
1677! AUTHORS & DATE
1678!     L. Lange
1679!     JB Clement, 2023-2025
1680!
1681! NOTES
1682!
1683!-----------------------------------------------------------------------
1684
1685! DECLARATION
1686! -----------
1687implicit none
1688
1689! ARGUMENTS
1690! ---------
1691  CHARACTER(LEN=*), INTENT(IN) :: var_name
1692  CHARACTER(LEN=*), INTENT(IN) :: title
1693  REAL,             INTENT(IN) :: var(:,:,:)
1694
1695! CODE
1696! ----
1697  CALL put_var_rgen(var_name,title,var,size(var))
1698
1699  END SUBROUTINE put_var_r3
1700!=======================================================================
1701
1702!=======================================================================
1703  SUBROUTINE put_var_rgen(var_name,title,var,var_size)
1704!-----------------------------------------------------------------------
1705! NAME
1706!     put_var_rgen
1707!
1708! DESCRIPTION
1709!     Generic subroutine to write variables to the restart file
1710
1711! AUTHORS & DATE
1712!     L. Lange
1713!     JB Clement, 2023-2025
1714
1715! NOTES
1716!
1717!-----------------------------------------------------------------------
1718
1719! DEPENDENCIES
1720! ------------
1721  USE netcdf, only: NF90_REDEF, NF90_DEF_VAR, NF90_ENDDEF, NF90_PUT_VAR, &
1722                    NF90_FLOAT, NF90_DOUBLE, &
1723                    NF90_PUT_ATT, NF90_NOERR, nf90_strerror, &
1724                    nf90_inq_dimid, nf90_inquire_dimension, NF90_INQ_VARID
1725  USE soil,               only:  nsoilmx_PEM
1726  USE comslope_mod,       only: nslope
1727  USE mod_phys_lmdz_para, only: is_master
1728
1729! DECLARATION
1730! -----------
1731implicit none
1732
1733! ARGUMENTS
1734! ---------
1735  CHARACTER(LEN=*), INTENT(IN) :: var_name
1736  CHARACTER(LEN=*), INTENT(IN) :: title
1737  INTEGER,          INTENT(IN) :: var_size
1738  REAL,             INTENT(IN) :: var(var_size)
1739
1740! LOCAL VARIABLES
1741! ---------------
1742  INTEGER       :: ierr
1743  INTEGER       :: nvarid
1744  INTEGER       :: idim1d
1745  LOGICAL, SAVE :: firsttime = .true.
1746
1747! CODE
1748! ----
1749    IF (is_master) THEN
1750
1751      IF (var_name=="Time") THEN
1752        ! Very specific case of "Time" variable
1753        if (firsttime) then
1754          ! Create the "Time variable"
1755          ierr=NF90_REDEF(nid_restart)
1756#ifdef NC_DOUBLE
1757          ierr=NF90_DEF_VAR(nid_restart,var_name,NF90_DOUBLE,&
1758                            (/idim7/),nvarid)
1759#else
1760          ierr=NF90_DEF_VAR(nid_restart,var_name,NF90_FLOAT,&
1761                            (/idim7/),nvarid)
1762#endif
1763          IF (LEN_TRIM(title) > 0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title)
1764          ierr=NF90_ENDDEF(nid_restart)
1765
1766          firsttime=.false.
1767        endif
1768        ! Append "time" value
1769        ! get current length of "Time" dimension
1770        ierr=nf90_inq_dimid(nid_restart,var_name,idim1d)
1771        ierr=nf90_inquire_dimension(nid_restart,idim1d,len=timeindex)
1772        timeindex=timeindex+1
1773        ierr=NF90_INQ_VARID(nid_restart,var_name,nvarid)
1774        ierr=NF90_PUT_VAR(nid_restart,nvarid,var,&
1775                           start=(/timeindex/))
1776        IF (ierr/=NF90_NOERR) THEN
1777          write(*,*)'put_var_rgen: problem writing Time'
1778          write(*,*)trim(nf90_strerror(ierr))
1779          CALL abort_physic("get_var_rgen","Failed to write Time",1)
1780        ENDIF
1781        return ! nothing left to do
1782      ELSEIF (var_size==length) THEN
1783        ! We know it is a "controle" kind of 1D array
1784        idim1d=idim1
1785      ELSEIF (var_size==nsoilmx_PEM) THEN
1786        ! We know it is an "mlayer" kind of 1D array
1787        idim1d=idim3
1788      ELSEIF (var_size==nslope+1) THEN
1789        ! We know it is an "inter slope" kind of 1D array
1790        idim1d=idim9
1791      ELSEIF (var_name == "nb_str_max") THEN
1792        ! We know it is a kind of stratification array
1793        idim1d = idim10
1794      ELSE
1795        write(*,*) "put_var_rgen error : wrong dimension"
1796        write(*,*) "  var_size =",var_size
1797        CALL abort_physic("get_var_rgen","Wrong variable dimension",1)
1798
1799      ENDIF ! of IF (var_size==length) THEN
1800
1801      ! Swich to NetCDF define mode
1802      ierr=NF90_REDEF (nid_restart)
1803      ! Define the variable
1804#ifdef NC_DOUBLE
1805      ierr=NF90_DEF_VAR(nid_restart,var_name,NF90_DOUBLE,(/idim1d/),nvarid)
1806#else
1807      ierr=NF90_DEF_VAR(nid_restart,var_name,NF90_FLOAT,(/idim1d/),nvarid)
1808#endif
1809      ! Add a "title" attribute
1810      IF (LEN_TRIM(title)>0) ierr=NF90_PUT_ATT(nid_restart,nvarid,"title",title)
1811      ! Swich out of define mode
1812      ierr=NF90_ENDDEF(nid_restart)
1813      ! Write variable to file
1814      ierr=NF90_PUT_VAR(nid_restart,nvarid,var)
1815      IF (ierr/=NF90_NOERR) THEN
1816        write(*,*)'put_var_rgen: problem writing '//trim(var_name)
1817        write(*,*)trim(nf90_strerror(ierr))
1818        CALL abort_physic("get_var_rgen","Failed writing variable",1)
1819      ENDIF
1820    ENDIF ! of IF (is_master)
1821
1822  END SUBROUTINE put_var_rgen
1823!=======================================================================
1824
1825END MODULE iostart_pem
Note: See TracBrowser for help on using the repository browser.