Changeset 4081 for trunk/LMDZ.GENERIC
- Timestamp:
- Feb 24, 2026, 9:59:11 AM (4 days ago)
- Location:
- trunk/LMDZ.GENERIC
- Files:
-
- 10 edited
-
changelog.txt (modified) (1 diff)
-
libf/phygeneric/aerosol_optical_properties_averaging.F (modified) (1 diff)
-
libf/phygeneric/physiq_mod.F90 (modified) (1 diff)
-
libf/phygeneric/rad_blackbody.F (modified) (1 diff)
-
libf/phygeneric/rad_correlatedk_init_stellar.F90 (modified) (2 diffs)
-
libf/phygeneric/rad_correlatedk_init_thermal.F90 (modified) (5 diffs)
-
libf/phygeneric/rad_correlatedk_online_recombination.F90 (modified) (2 diffs)
-
libf/phygeneric/rad_correlatedk_rayleigh_scattering_opacity.F90 (modified) (1 diff)
-
libf/phygeneric/rad_correlatedk_stellar_spectrum.F90 (modified) (7 diffs)
-
libf/phygeneric/radcommon_h.F90 (modified) (1 diff)
Legend:
- Unmodified
- Added
- Removed
-
trunk/LMDZ.GENERIC/changelog.txt
r4077 r4081 2179 2179 == 18/02/2026 == MT 2180 2180 Changes (renaming, merging) of radiative routine in phygeneric directory, following the hackathon (04/02/2026) 2181 2182 == 24/02/2026 == EM 2183 OpenMP bug fix in "rad_correlatedk_stellar_spectrum", all intermediate computations 2184 should be done by all OpenMp threads. 2185 While at it did some cleanup: 2186 - added some OMP threadprivate statements in radcommon_h.F90 (not alwyas 2187 necessary, but best practice is that all saved variables be threadprivate) 2188 - turned rad_blackbody.F into a module and modernized routines. 2189 - use clearer stategy (wrt OpenMP) in "rad_correlatedk_init_thermal.F90" and 2190 "rad_correlatedk_init_stellar.F90": have the master read in file and data 2191 and then broadcast to all cores. 2192 -
trunk/LMDZ.GENERIC/libf/phygeneric/aerosol_optical_properties_averaging.F
r4077 r4081 4 4 & epir,omegir,gir,qref,omegaref) 5 5 6 USE rad_blackbody_mod, ONLY: rad_blackbody_planck_law_wavelength 6 7 7 8 IMPLICIT NONE -
trunk/LMDZ.GENERIC/libf/phygeneric/physiq_mod.F90
r4079 r4081 25 25 use rad_correlatedk_ini_aerosol_mod, only: rad_correlatedk_ini_aerosol 26 26 use rad_correlatedk_init_stellar_mod, only: rad_correlatedk_init_stellar 27 use rad_correlatedk_init_thermal_mod, only: rad_correlatedk_init_thermal 27 28 use aerosol_radius, only: aerosol_radius_h2o_liquid_ice_mixture, aerosol_radius_co2 28 29 use aerosol_global_variables , only: aerosol_init, iaero_co2, iaero_h2o -
trunk/LMDZ.GENERIC/libf/phygeneric/rad_blackbody.F
r4077 r4081 1 module rad_blackbody_mod 2 3 implicit none 4 5 contains 6 1 7 subroutine rad_blackbody_planck_law_wavelength(blalong,blat,blae) 2 8 3 implicit double precision (a-h,o-z) 9 implicit none 10 11 double precision, intent(in) :: blalong 12 double precision, intent(in) :: blat 13 double precision, intent(out) :: blae 14 15 ! double precision :: sigma,pi,c0,h,cbol,rind,c,c1,c2 4 16 5 17 ! physical constants 6 sigma=5.670374D-87 pi=datan(1.d0)*4.d08 c0=2.997925d+089 h=6.62607d-3410 cbol=1.380649d-2311 rind=1.d012 c=c0/rind13 c1=h*(c**2)14 c2=h*c/cbol18 double precision, parameter :: sigma=5.670374D-8 19 double precision, parameter :: pi=atan(1.d0)*4.d0 20 double precision, parameter :: c0=2.997925d+08 21 double precision, parameter :: h=6.62607d-34 22 double precision, parameter :: cbol=1.380649d-23 23 double precision, parameter :: rind=1.d0 24 double precision, parameter :: c=c0/rind 25 double precision, parameter :: c1=h*(c**2) 26 double precision, parameter :: c2=h*c/cbol 15 27 16 28 17 blae=2.d0*pi*c1/blalong**5/( dexp(c2/blalong/blat)-1.d0)29 blae=2.d0*pi*c1/blalong**5/(exp(c2/blalong/blat)-1.d0) 18 30 19 20 return 21 end 31 end subroutine rad_blackbody_planck_law_wavelength 22 32 23 33 subroutine rad_blackbody_planck_law_wavenumber(blalong,blat,blae) 24 34 25 implicit double precision (a-h,o-z) 35 implicit none 36 37 double precision, intent(in) :: blalong 38 double precision, intent(in) :: blat 39 double precision, intent(out) :: blae 40 41 ! double precision :: sigma,pi,c0,h,cbol,rind,c,c1,c2 26 42 27 43 ! physical constants 28 sigma=5.670374D-829 pi=datan(1.d0)*4.d030 c0=2.997925d+0831 h=6.62607d-3432 cbol=1.380649d-2333 rind=1.d034 c=c0/rind35 c1=h*(c**2)36 c2=h*c/cbol44 double precision, parameter :: sigma=5.670374D-8 45 double precision, parameter :: pi=atan(1.d0)*4.d0 46 double precision, parameter :: c0=2.997925d+08 47 double precision, parameter :: h=6.62607d-34 48 double precision, parameter :: cbol=1.380649d-23 49 double precision, parameter :: rind=1.d0 50 double precision, parameter :: c=c0/rind 51 double precision, parameter :: c1=h*(c**2) 52 double precision, parameter :: c2=h*c/cbol 37 53 38 54 39 blae=2.d0*pi*c1*blalong**3/( dexp(c2*blalong/blat)-1.d0)55 blae=2.d0*pi*c1*blalong**3/(exp(c2*blalong/blat)-1.d0) 40 56 41 42 return43 end 57 end subroutine rad_blackbody_planck_law_wavenumber 58 59 end module rad_blackbody_mod -
trunk/LMDZ.GENERIC/libf/phygeneric/rad_correlatedk_init_stellar.F90
r4077 r4081 34 34 use callkeys_mod, only: Fat1AU 35 35 use rad_correlatedk_stellar_spectrum_mod, only: rad_correlatedk_stellar_spectrum 36 use mod_phys_lmdz_para, only : is_master, bcast 36 37 37 38 implicit none … … 74 75 endif 75 76 76 !$OMP MASTER77 nb=078 ierr=079 ! check that the file contains the right number of bands80 open(131,file=file_path,form='formatted')81 read(131,*,iostat=ierr) file_entries82 do while (ierr==0)77 if (is_master) then ! only the master needs to read in BWNV(:) from file 78 nb=0 79 ierr=0 80 ! check that the file contains the right number of bands 81 open(131,file=file_path,form='formatted') 82 read(131,*,iostat=ierr) file_entries 83 do while (ierr==0) 83 84 read(131,*,iostat=ierr) dummy 84 85 if (ierr==0) nb=nb+1 85 enddo86 close(131)86 enddo 87 close(131) 87 88 88 write(*,*) 'rad_correlatedk_init_stellar: L_NSPECTV = ',L_NSPECTV, 'in the model '89 write(*,*) ' there are ',nb, 'entries in ',TRIM(file_path)90 if(nb.ne.L_NSPECTV) then89 write(*,*) 'rad_correlatedk_init_stellar: L_NSPECTV = ',L_NSPECTV, 'in the model ' 90 write(*,*) ' there are ',nb, 'entries in ',TRIM(file_path) 91 if(nb.ne.L_NSPECTV) then 91 92 write(*,*) 'MISMATCH !! I stop here' 92 93 call abort_physic("rad_correlatedk_init_stellar","The number of entries in narrowbands_VI.in does not match with L_NSPECTV",1) 93 endif94 endif 94 95 95 ! load and display the data96 open(111,file=file_path,form='formatted')97 read(111,*)98 do M=1,L_NSPECTV-196 ! load and display the data 97 open(111,file=file_path,form='formatted') 98 read(111,*) 99 do M=1,L_NSPECTV-1 99 100 read(111,*) BWNV(M) 100 end do 101 read(111,*) lastband 102 close(111) 103 BWNV(L_NSPECTV) =lastband(1) 104 BWNV(L_NSPECTV+1)=lastband(2) 105 !$OMP END MASTER 106 !$OMP BARRIER 101 end do 102 read(111,*) lastband 103 close(111) 104 BWNV(L_NSPECTV) =lastband(1) 105 BWNV(L_NSPECTV+1)=lastband(2) 107 106 108 print*,'rad_correlatedk_init_stellar: VI band limits:'109 do M=1,L_NSPECTV+1107 print*,'rad_correlatedk_init_stellar: VI band limits:' 108 do M=1,L_NSPECTV+1 110 109 print*,m,'-->',BWNV(M),' cm^-1' 111 end do 112 print*,' ' 110 end do 111 print*,' ' 112 endif ! of if (is_master) 113 114 ! Broadcast BWNV to all cores 115 call bcast(BWNV) 113 116 114 117 ! Set up mean wavenumbers and wavenumber deltas. Units of -
trunk/LMDZ.GENERIC/libf/phygeneric/rad_correlatedk_init_thermal.F90
r4077 r4081 1 module rad_correlatedk_init_thermal_mod 2 3 implicit none 4 5 contains 6 1 7 subroutine rad_correlatedk_init_thermal 2 8 … … 26 32 use datafile_mod, only: datadir 27 33 use comcstfi_mod, only: pi 34 use mod_phys_lmdz_para, only : is_master, bcast 28 35 29 36 implicit none … … 90 97 endif 91 98 92 !$OMP MASTER93 nb=094 ierr=095 ! check that the file contains the right number of bands96 open(131,file=file_path,form='formatted')97 read(131,*,iostat=ierr) file_entries98 do while (ierr==0)99 if (is_master) then ! only the master needs to read in BWNI(:) from file 100 nb=0 101 ierr=0 102 ! check that the file contains the right number of bands 103 open(131,file=file_path,form='formatted') 104 read(131,*,iostat=ierr) file_entries 105 do while (ierr==0) 99 106 read(131,*,iostat=ierr) dummy 100 107 ! write(*,*) 'rad_correlatedk_init_thermal: file_entries:',dummy,'ierr=',ierr 101 108 if (ierr==0) nb=nb+1 102 enddo103 close(131)104 105 write(*,*) 'rad_correlatedk_init_thermal: L_NSPECTI = ',L_NSPECTI, 'in the model '106 write(*,*) ' there are ',nb, 'entries in ',TRIM(file_path)107 if(nb.ne.L_NSPECTI) then109 enddo 110 close(131) 111 112 write(*,*) 'rad_correlatedk_init_thermal: L_NSPECTI = ',L_NSPECTI, 'in the model ' 113 write(*,*) ' there are ',nb, 'entries in ',TRIM(file_path) 114 if(nb.ne.L_NSPECTI) then 108 115 write(*,*) 'MISMATCH !! I stop here' 109 116 call abort_physic("rad_correlatedk_init_thermal","The number of entries in narrowbands_IR.in does not match with L_NSPECTI",1) 110 endif111 112 ! load and display the data113 open(111,file=file_path,form='formatted')114 read(111,*)115 do M=1,L_NSPECTI-1117 endif 118 119 ! load and display the data 120 open(111,file=file_path,form='formatted') 121 read(111,*) 122 do M=1,L_NSPECTI-1 116 123 read(111,*) BWNI(M) 117 end do 118 read(111,*) lastband 119 close(111) 120 BWNI(L_NSPECTI) =lastband(1) 121 BWNI(L_NSPECTI+1)=lastband(2) 122 !$OMP END MASTER 123 !$OMP BARRIER 124 125 print*,'' 126 print*,'rad_correlatedk_init_thermal: IR band limits:' 127 do M=1,L_NSPECTI+1 124 end do 125 read(111,*) lastband 126 close(111) 127 BWNI(L_NSPECTI) =lastband(1) 128 BWNI(L_NSPECTI+1)=lastband(2) 129 130 print*,'' 131 print*,'rad_correlatedk_init_thermal: IR band limits:' 132 do M=1,L_NSPECTI+1 128 133 print*,m,'-->',BWNI(M),' cm^-1' 129 end do 134 end do 135 endif ! of if (is_master) 136 137 ! Broadcast BWNI to all cores 138 call bcast(BWNI) 130 139 131 140 ! Set up mean wavenumbers and wavenumber deltas. Units of … … 151 160 print*,'T = ',dble(NTstart)/NTfac, ' to ',dble(NTstop)/NTfac,' K.' 152 161 153 IF(.NOT.ALLOCATED(planckir))ALLOCATE(planckir(L_NSPECTI,NTstop-NTstart+1))162 ALLOCATE(planckir(L_NSPECTI,NTstop-NTstart+1)) 154 163 155 164 do NW=1,L_NSPECTI … … 223 232 endif 224 233 225 return226 234 end subroutine rad_correlatedk_init_thermal 235 236 end module rad_correlatedk_init_thermal_mod -
trunk/LMDZ.GENERIC/libf/phygeneric/rad_correlatedk_online_recombination.F90
r4077 r4081 38 38 ! Following arrays are allocated in rad_correlatedk_recombination_setup (excepted otf2tot_idx, in rad_correlatedk_recombination_init) and deallocated in callcork lastcall 39 39 REAL*8, save, DIMENSION(:,:,:,:,:), ALLOCATABLE :: gasi_recomb, gasv_recomb 40 !$OMP THREADPRIVATE(gasi_recomb,gasv_recomb) 40 41 REAL*8, save, DIMENSION(:,:,:,:,:), ALLOCATABLE :: gasi_otf, gasv_otf 42 !$OMP THREADPRIVATE(gasi_otf,gasv_otf) 41 43 REAL*8, save, DIMENSION(:), ALLOCATABLE :: w_cum 42 44 REAL*8, save, DIMENSION(:), ALLOCATABLE :: wtwospec 45 !$OMP THREADPRIVATE(w_cum,wtwospec) 43 46 44 47 INTEGER, save, DIMENSION(:), ALLOCATABLE :: otf2tot_idx … … 47 50 48 51 INTEGER, save, DIMENSION(:), ALLOCATABLE :: permut_idx 49 !$OMP THREADPRIVATE( gasi_recomb,gasv_recomb,w_cum,otf2tot_idx,rcb2tot_idx,otf2rcb_idx,permut_idx)52 !$OMP THREADPRIVATE(otf2tot_idx,rcb2tot_idx,otf2rcb_idx,permut_idx) 50 53 51 54 CONTAINS -
trunk/LMDZ.GENERIC/libf/phygeneric/rad_correlatedk_rayleigh_scattering_opacity.F90
r4077 r4081 38 38 use comcstfi_mod, only: g, pi 39 39 use callkeys_mod, only: strictboundrayleigh 40 use rad_blackbody_mod, only: rad_blackbody_planck_law_wavelength 40 41 41 42 implicit none -
trunk/LMDZ.GENERIC/libf/phygeneric/rad_correlatedk_stellar_spectrum.F90
r4077 r4081 35 35 use callkeys_mod, only: stelbbody,stelTbb 36 36 use ioipsl_getin_p_mod, only: getin_p 37 use rad_blackbody_mod, only: rad_blackbody_planck_law_wavelength 38 use mod_phys_lmdz_para, only : is_master, bcast 37 39 38 40 implicit none 39 41 40 real*8 STELLAR(L_NSPECTV)42 real*8,intent(out) :: STELLAR(L_NSPECTV) 41 43 42 44 integer Nfine … … 44 46 integer ifine,band 45 47 46 real,allocatable ,save :: lam(:),stel_f(:) ! read by master thread47 ! but used by all threads48 real,allocatable :: lam(:),stel_f(:) ! read by master 49 48 50 real lamm,lamp 49 51 real dl … … 61 63 STELLAR(:)=0.0 62 64 63 print*,'enter ave_stellspec' 65 if (is_master) print*,'entering rad_correlatedk_stellar_spectrum' 66 64 67 if(stelbbody)then 65 68 tstellar=stelTbb … … 89 92 end if 90 93 91 write(*,*) "Input stellar temperature is:" 92 write(*,*) "tstellar = ",tstellar 94 if (is_master) then 95 write(*,*) "Input stellar temperature is:" 96 write(*,*) "tstellar = ",tstellar 97 endif 93 98 94 99 ! load high resolution stellar data … … 97 102 call getin_p("stelspec_file",stelspec_file) ! default path 98 103 99 write(*,*) "Input stellar spectrum file is:" 100 write(*,*) "stelspec_file = ",trim(stelspec_file) 101 write(*,*) 'Please use ',1,' and only ',1,' header line in ',trim(stelspec_file) 104 if (is_master) then ! only the master needs to do this 105 write(*,*) "Input stellar spectrum file is:" 106 write(*,*) "stelspec_file = ",trim(stelspec_file) 107 write(*,*) 'Please use ',1,' and only ',1,' header line in ',trim(stelspec_file) 102 108 103 ! Check the target file is there104 file_path = trim(datadir)//'/stellar_spectra/'//stelspec_file105 print*, 'stellar flux : ', file_path106 inquire(FILE=file_path,EXIST=file_exists)109 ! Check the target file is there 110 file_path = trim(datadir)//'/stellar_spectra/'//stelspec_file 111 print*, 'stellar flux : ', file_path 112 inquire(FILE=file_path,EXIST=file_exists) 107 113 108 if (.not.file_exists) THEN114 if (.not.file_exists) THEN 109 115 write(*,*)'Beware that startype is now deprecated, you should use ' 110 116 write(*,*)'stelspec_file and tstellar to define the input stellar spectrum.' … … 122 128 write(*,*)'https://web.lmd.jussieu.fr/~lmdz/planets/generic/datagcm/stellar_spectra/' 123 129 call abort_physic("rad_correlatedk_stellar_spectrum", "Unable to read stellar flux file", 1) 124 end if130 end if 125 131 126 !$OMP MASTER 127 ! Open the file 128 OPEN(UNIT=110,FILE=file_path,STATUS='old',iostat=ios) 129 ! Get number of line in the file 130 READ(110,*) ! skip first line header just in case 131 Nfine = 0 132 do 132 ! Open the file 133 OPEN(UNIT=110,FILE=file_path,STATUS='old',iostat=ios) 134 ! Get number of line in the file 135 READ(110,*) ! skip first line header just in case 136 Nfine = 0 137 do 133 138 read(110,*,iostat=ios) 134 139 if (ios<0) exit 135 140 Nfine = Nfine + 1 136 end do137 rewind(110) ! Rewind file after counting lines138 READ(110,*) ! skip first line header just in case141 end do 142 rewind(110) ! Rewind file after counting lines 143 READ(110,*) ! skip first line header just in case 139 144 140 allocate(lam(Nfine),stel_f(Nfine)) 145 ! allocate arrays 146 allocate(lam(Nfine),stel_f(Nfine)) 141 147 142 do ifine=1,Nfine 148 ! read arrays 149 do ifine=1,Nfine 143 150 read(110,*) lam(ifine), stel_f(ifine) ! lam [um] stel_f [per unit of wavelength] (integrated and normalized by Fat1AU) 144 enddo151 enddo 145 152 146 !$OMP END MASTER 147 !$OMP BARRIER 153 endif ! of if(is_master) 154 155 ! brodcast Nfine, allocate arrays for all except master 156 call bcast(Nfine) 157 if (.not.is_master) allocate(lam(Nfine),stel_f(Nfine)) 158 ! broadcast arrays 159 call bcast(lam) 160 call bcast(stel_f) 148 161 162 149 163 ! sum data by band 150 164 band=1 … … 166 180 167 181 STELLAR(1:L_NSPECTV)=STELLAR(1:L_NSPECTV)/sum(STELLAR(1:L_NSPECTV)) 168 !$OMP BARRIER 169 !$OMP MASTER 170 deallocate(lam)182 183 ! cleanup: deallocate arrays 184 deallocate(lam) 171 185 deallocate(stel_f) 172 !$OMP END MASTER 173 !$OMP BARRIER 186 174 187 endif !stelbbody 175 188 -
trunk/LMDZ.GENERIC/libf/phygeneric/radcommon_h.F90
r4077 r4081 63 63 ! 64 64 65 REAL*8 BWNI(L_NSPECTI+1), WNOI(L_NSPECTI), DWNI(L_NSPECTI), WAVEI(L_NSPECTI) !BWNI read by master in rad_correlatedk_init_thermal 66 REAL*8 BWNV(L_NSPECTV+1), WNOV(L_NSPECTV), DWNV(L_NSPECTV), WAVEV(L_NSPECTV) !BWNV read by master in rad_correlatedk_init_stellar 65 REAL*8 BWNI(L_NSPECTI+1) !BWNI read by master in rad_correlatedk_init_thermal 66 !$OMP THREADPRIVATE(BWNI) 67 REAL*8 WNOI(L_NSPECTI), DWNI(L_NSPECTI), WAVEI(L_NSPECTI) 68 !$OMP THREADPRIVATE(WNOI,DWNI,WAVEI) 69 REAL*8 BWNV(L_NSPECTV+1) !BWNV read by master in rad_correlatedk_init_stellar 70 !$OMP THREADPRIVATE(BWNV) 71 REAL*8 WNOV(L_NSPECTV), DWNV(L_NSPECTV), WAVEV(L_NSPECTV) 67 72 REAL*8 STELLARF(L_NSPECTV) 68 !$OMP THREADPRIVATE(WNOI,DWNI,WAVEI,& 69 !$OMP WNOV,DWNV,WAVEV,& 70 !$OMP STELLARF) 73 !$OMP THREADPRIVATE(WNOV,DWNV,WAVEV,STELLARF) 71 74 72 75 REAL*8 blami(L_NSPECTI+1)
Note: See TracChangeset
for help on using the changeset viewer.
