Skip to content

Commit adbdeae

Browse files
committed
Merge branch 'atmosphere/physics_mmm_f2008_compliance' into release-v8.0.0 (PR #1086)
This merge corrects Fortran 2008 compliance issues in physics modules that were caught by the GNU Fortran compiler when using the -std=f2008 flag. Most of these issues are related to the use of non-standard kind/type specifications.
2 parents 644bb34 + 1a3eed0 commit adbdeae

File tree

3 files changed

+48
-35
lines changed

3 files changed

+48
-35
lines changed

src/core_atmosphere/physics/physics_mmm/cu_ntiedtke.F

Lines changed: 16 additions & 8 deletions
Original file line numberDiff line numberDiff line change
@@ -2527,9 +2527,10 @@ subroutine cudlfsn &
25272527
implicit none
25282528

25292529
!--- input arguments:
2530+
integer,intent(in):: klon
25302531
logical,intent(in),dimension(klon):: ldcum
25312532

2532-
integer,intent(in):: klev,klon
2533+
integer,intent(in):: klev
25332534
integer,intent(in),dimension(klon):: lndj
25342535
integer,intent(in),dimension(klon):: kcbot,kctop
25352536

@@ -2737,9 +2738,10 @@ subroutine cuddrafn &
27372738
implicit none
27382739

27392740
!--- input arguments:
2741+
integer,intent(in)::klon
27402742
logical,intent(in),dimension(klon):: lddraf
27412743

2742-
integer,intent(in)::klev,klon
2744+
integer,intent(in)::klev
27432745

27442746
real(kind=kind_phys),intent(in),dimension(klon,klev):: ptenh,pqenh,puen,pven
27452747
real(kind=kind_phys),intent(in),dimension(klon,klev):: pgeo,pmfu
@@ -2980,9 +2982,10 @@ subroutine cuflxn &
29802982
implicit none
29812983

29822984
!--- input arguments:
2985+
integer,intent(in):: klon
29832986
logical,intent(in),dimension(klon):: ldcum
29842987

2985-
integer,intent(in):: klev,klon
2988+
integer,intent(in):: klev
29862989
integer,intent(in),dimension(klon):: lndj
29872990
integer,intent(in),dimension(klon):: kcbot,kctop,kdtop
29882991

@@ -3237,9 +3240,10 @@ subroutine cudtdqn(klon,klev,ktopm2,kctop,kdtop,ldcum, &
32373240
implicit none
32383241

32393242
!--- input arguments:
3243+
integer,intent(in):: klon
32403244
logical,intent(in),dimension(klon):: ldcum,lddraf
32413245

3242-
integer,intent(in):: klon,klev,ktopm2
3246+
integer,intent(in):: klev,ktopm2
32433247
integer,intent(in),dimension(klon):: kctop,kdtop
32443248

32453249
real(kind=kind_phys),intent(in):: ztmst
@@ -3323,8 +3327,9 @@ subroutine cududvn(klon,klev,ktopm2,ktype,kcbot,kctop,ldcum, &
33233327
implicit none
33243328

33253329
!--- input arguments:
3330+
integer,intent(in):: klon
33263331
logical,intent(in),dimension(klon):: ldcum
3327-
integer,intent(in):: klon,klev,ktopm2
3332+
integer,intent(in):: klev,ktopm2
33283333
integer,intent(in),dimension(klon):: ktype,kcbot,kctop
33293334

33303335
real(kind=kind_phys),intent(in):: ztmst
@@ -3466,8 +3471,9 @@ subroutine cuadjtqn &
34663471
implicit none
34673472

34683473
!--- input arguments:
3474+
integer,intent(in):: klon
34693475
logical,intent(in),dimension(klon):: ldflag
3470-
integer,intent(in):: kcall,kk,klev,klon
3476+
integer,intent(in):: kcall,kk,klev
34713477

34723478
real(kind=kind_phys),intent(in),dimension(klon):: psp
34733479

@@ -3597,8 +3603,9 @@ subroutine cubasmcn &
35973603
! ----------------------------------------------------------------
35983604

35993605
!--- input arguments:
3606+
integer,intent(in):: klon
36003607
logical,intent(in),dimension(klon):: ldcum
3601-
integer,intent(in):: klon,kk,klev,klevm1
3608+
integer,intent(in):: kk,klev,klevm1
36023609

36033610
real(kind=kind_phys),intent(in),dimension(klon,klev):: pten,pqen,pqsen,pgeo,pverv
36043611
real(kind=kind_phys),intent(in),dimension(klon,klev):: puen,pven ! not used.
@@ -3657,9 +3664,10 @@ subroutine cuentrn(klon,klev,kk,kcbot,ldcum,ldwork, &
36573664

36583665
!--- input arguments:
36593666
logical,intent(in):: ldwork
3667+
integer,intent(in):: klon
36603668
logical,intent(in),dimension(klon):: ldcum
36613669

3662-
integer,intent(in):: klon,klev,kk
3670+
integer,intent(in):: klev,kk
36633671
integer,intent(in),dimension(klon):: kcbot
36643672

36653673
real(kind=kind_phys),intent(in),dimension(klon,klev):: pmfu

src/core_atmosphere/physics/physics_mmm/module_libmassv.F

Lines changed: 12 additions & 10 deletions
Original file line numberDiff line numberDiff line change
@@ -14,6 +14,8 @@ module module_libmassv
1414
module procedure vsqrt_s
1515
end interface
1616

17+
integer, parameter, private :: R4KIND = selected_real_kind(6)
18+
integer, parameter, private :: R8KIND = selected_real_kind(12)
1719

1820
contains
1921

@@ -22,14 +24,14 @@ module module_libmassv
2224
subroutine vrec_d(y,x,n)
2325
!=================================================================================================================
2426
integer,intent(in):: n
25-
real*8,dimension(*),intent(in):: x
26-
real*8,dimension(*),intent(out):: y
27+
real(kind=R8KIND),dimension(*),intent(in):: x
28+
real(kind=R8KIND),dimension(*),intent(out):: y
2729

2830
integer:: j
2931
!-----------------------------------------------------------------------------------------------------------------
3032

3133
do j=1,n
32-
y(j)=1.d0/x(j)
34+
y(j)=real(1.0,kind=R8KIND)/x(j)
3335
enddo
3436

3537
end subroutine vrec_d
@@ -38,14 +40,14 @@ end subroutine vrec_d
3840
subroutine vrec_s(y,x,n)
3941
!=================================================================================================================
4042
integer,intent(in):: n
41-
real*4,dimension(*),intent(in):: x
42-
real*4,dimension(*),intent(out):: y
43+
real(kind=R4KIND),dimension(*),intent(in):: x
44+
real(kind=R4KIND),dimension(*),intent(out):: y
4345

4446
integer:: j
4547
!-----------------------------------------------------------------------------------------------------------------
4648

4749
do j=1,n
48-
y(j)=1.e0/x(j)
50+
y(j)=real(1.0,kind=R4KIND)/x(j)
4951
enddo
5052

5153
end subroutine vrec_s
@@ -54,8 +56,8 @@ end subroutine vrec_s
5456
subroutine vsqrt_d(y,x,n)
5557
!=================================================================================================================
5658
integer,intent(in):: n
57-
real*8,dimension(*),intent(in):: x
58-
real*8,dimension(*),intent(out):: y
59+
real(kind=R8KIND),dimension(*),intent(in):: x
60+
real(kind=R8KIND),dimension(*),intent(out):: y
5961

6062
integer:: j
6163
!-----------------------------------------------------------------------------------------------------------------
@@ -71,8 +73,8 @@ subroutine vsqrt_s(y,x,n)
7173
!=================================================================================================================
7274

7375
integer,intent(in):: n
74-
real*4,dimension(*),intent(in):: x
75-
real*4,dimension(*),intent(out):: y
76+
real(kind=R4KIND),dimension(*),intent(in):: x
77+
real(kind=R4KIND),dimension(*),intent(out):: y
7678

7779
integer:: j
7880

src/core_atmosphere/physics/physics_mmm/mp_radar.F

Lines changed: 20 additions & 17 deletions
Original file line numberDiff line numberDiff line change
@@ -30,6 +30,9 @@ module mp_radar
3030
!.. calculations based on 50 individual size bins of the distributions.
3131
!+---+-----------------------------------------------------------------+
3232
33+
integer, parameter, private :: R4KIND = selected_real_kind(6)
34+
integer, parameter, private :: R8KIND = selected_real_kind(12)
35+
3336
integer,parameter,public:: nrbins = 50
3437
integer,parameter,public:: slen = 20
3538
character(len=slen), public:: &
@@ -38,7 +41,7 @@ module mp_radar
3841
mixingrulestring_g, matrixstring_g, inclusionstring_g, &
3942
hoststring_g, hostmatrixstring_g, hostinclusionstring_g
4043
41-
complex*16,public:: m_w_0, m_i_0
44+
complex(kind=R8KIND),public:: m_w_0, m_i_0
4245
4346
double precision,dimension(nrbins+1),public:: xxdx
4447
double precision,dimension(nrbins),public:: xxds,xdts,xxdg,xdtg
@@ -123,7 +126,7 @@ subroutine radar_init
123126
xxdx(1) = 100.d-6
124127
xxdx(nrbins+1) = 0.02d0
125128
do n = 2, nrbins
126-
xxdx(n) = dexp(dfloat(n-1)/dfloat(nrbins) &
129+
xxdx(n) = dexp(real(n-1,kind=R8KIND)/real(nrbins,kind=R8KIND) &
127130
* dlog(xxdx(nrbins+1)/xxdx(1)) +dlog(xxdx(1)))
128131
enddo
129132
do n = 1, nrbins
@@ -135,7 +138,7 @@ subroutine radar_init
135138
xxdx(1) = 100.d-6
136139
xxdx(nrbins+1) = 0.05d0
137140
do n = 2, nrbins
138-
xxdx(n) = dexp(dfloat(n-1)/dfloat(nrbins) &
141+
xxdx(n) = dexp(real(n-1,kind=R8KIND)/real(nrbins,kind=R8KIND) &
139142
* dlog(xxdx(nrbins+1)/xxdx(1)) +dlog(xxdx(1)))
140143
enddo
141144
do n = 1, nrbins
@@ -196,7 +199,7 @@ subroutine rayleigh_soak_wetgraupel(x_g,a_geo,b_geo,fmelt,meltratio_outside,m_w,
196199
character(len=*), intent(in):: mixingrule, matrix, inclusion, &
197200
host, hostmatrix, hostinclusion
198201
199-
complex*16,intent(in):: m_w, m_i
202+
complex(kind=R8KIND),intent(in):: m_w, m_i
200203
201204
double precision, intent(in):: x_g, a_geo, b_geo, fmelt, lambda, meltratio_outside
202205
@@ -206,7 +209,7 @@ subroutine rayleigh_soak_wetgraupel(x_g,a_geo,b_geo,fmelt,meltratio_outside,m_w,
206209
!--- local variables:
207210
integer:: error
208211
209-
complex*16:: m_core, m_air
212+
complex(kind=R8KIND):: m_core, m_air
210213
211214
double precision, parameter:: pix=3.1415926535897932384626434d0
212215
double precision:: d_large, d_g, rhog, x_w, xw_a, fm, fmgrenz, &
@@ -340,7 +343,7 @@ real(kind=kind_phys) function gammln(xx)
340343
end function gammln
341344
342345
!=================================================================================================================
343-
complex*16 function get_m_mix_nested (m_a, m_i, m_w, volair, &
346+
complex(kind=R8KIND) function get_m_mix_nested (m_a, m_i, m_w, volair, &
344347
volice, volwater, mixingrule, host, matrix, &
345348
inclusion, hostmatrix, hostinclusion, cumulerror)
346349
implicit none
@@ -350,7 +353,7 @@ complex*16 function get_m_mix_nested (m_a, m_i, m_w, volair, &
350353
character(len=*),intent(in):: mixingrule, host, matrix, &
351354
inclusion, hostmatrix, hostinclusion
352355
353-
complex*16,intent(in):: m_a, m_i, m_w
356+
complex(kind=R8KIND),intent(in):: m_a, m_i, m_w
354357
355358
double precision,intent(in):: volice, volair, volwater
356359
@@ -360,7 +363,7 @@ complex*16 function get_m_mix_nested (m_a, m_i, m_w, volair, &
360363
!--- local variables:
361364
integer:: error
362365
363-
complex*16:: mtmp
366+
complex(kind=R8KIND):: mtmp
364367
365368
double precision:: vol1, vol2
366369
@@ -481,7 +484,7 @@ complex*16 function get_m_mix_nested (m_a, m_i, m_w, volair, &
481484
end function get_m_mix_nested
482485
483486
!=================================================================================================================
484-
complex*16 function get_m_mix (m_a, m_i, m_w, volair, volice, &
487+
complex(kind=R8KIND) function get_m_mix (m_a, m_i, m_w, volair, volice, &
485488
volwater, mixingrule, matrix, inclusion, &
486489
error)
487490
implicit none
@@ -490,7 +493,7 @@ complex*16 function get_m_mix (m_a, m_i, m_w, volair, volice, &
490493
!--- input arguments:
491494
character(len=*),intent(in):: mixingrule, matrix, inclusion
492495
493-
complex*16, intent(in):: m_a, m_i, m_w
496+
complex(kind=R8KIND), intent(in):: m_a, m_i, m_w
494497
495498
double precision, intent(in):: volice, volair, volwater
496499
@@ -531,15 +534,15 @@ complex*16 function get_m_mix (m_a, m_i, m_w, volair, volice, &
531534
end function get_m_mix
532535
533536
!=================================================================================================================
534-
complex*16 function m_complex_maxwellgarnett(vol1, vol2, vol3, &
537+
complex(kind=R8KIND) function m_complex_maxwellgarnett(vol1, vol2, vol3, &
535538
m1, m2, m3, inclusion, error)
536539
implicit none
537540
!=================================================================================================================
538541
539542
!--- input arguments:
540543
character(len=*),intent(in):: inclusion
541544
542-
complex*16,intent(in):: m1,m2,m3
545+
complex(kind=R8KIND),intent(in):: m1,m2,m3
543546
544547
double precision,intent(in):: vol1,vol2,vol3
545548
@@ -548,7 +551,7 @@ complex*16 function m_complex_maxwellgarnett(vol1, vol2, vol3, &
548551
integer,intent(out):: error
549552
550553
!--- local variables:
551-
complex*16 :: beta2, beta3, m1t, m2t, m3t
554+
complex(kind=R8KIND) :: beta2, beta3, m1t, m2t, m3t
552555
553556
!-----------------------------------------------------------------------------------------------------------------
554557
@@ -576,7 +579,7 @@ complex*16 function m_complex_maxwellgarnett(vol1, vol2, vol3, &
576579
else
577580
write(radar_debug,*) 'M_COMPLEX_MAXWELLGARNETT: ', 'unknown inclusion: ', inclusion
578581
call physics_message(radar_debug)
579-
m_complex_maxwellgarnett=dcmplx(-999.99d0,-999.99d0)
582+
m_complex_maxwellgarnett=cmplx(-999.99d0,-999.99d0,kind=R8KIND)
580583
error = 1
581584
return
582585
endif
@@ -587,7 +590,7 @@ complex*16 function m_complex_maxwellgarnett(vol1, vol2, vol3, &
587590
end function m_complex_maxwellgarnett
588591
589592
!=================================================================================================================
590-
complex*16 function m_complex_water_ray(lambda,t)
593+
complex(kind=R8KIND) function m_complex_water_ray(lambda,t)
591594
implicit none
592595
!=================================================================================================================
593596
@@ -603,7 +606,7 @@ complex*16 function m_complex_water_ray(lambda,t)
603606
double precision,parameter:: pix=3.1415926535897932384626434d0
604607
double precision:: epsinf,epss,epsr,epsi
605608
double precision:: alpha,lambdas,sigma,nenner
606-
complex*16,parameter:: i = (0d0,1d0)
609+
complex(kind=R8KIND),parameter:: i = (0d0,1d0)
607610
608611
!-----------------------------------------------------------------------------------------------------------------
609612
@@ -627,7 +630,7 @@ complex*16 function m_complex_water_ray(lambda,t)
627630
end function m_complex_water_ray
628631
629632
!=================================================================================================================
630-
complex*16 function m_complex_ice_maetzler(lambda,t)
633+
complex(kind=R8KIND) function m_complex_ice_maetzler(lambda,t)
631634
implicit none
632635
!=================================================================================================================
633636

0 commit comments

Comments
 (0)