Changeset 3990 for trunk/LMDZ.GENERIC/libf
- Timestamp:
- Dec 12, 2025, 5:36:11 PM (3 months ago)
- Location:
- trunk/LMDZ.GENERIC/libf/phystd
- Files:
-
- 3 edited
-
ave_stelspec.F90 (modified) (8 diffs)
-
physiq_mod.F90 (modified) (1 diff)
-
setspv.F90 (modified) (3 diffs)
Legend:
- Unmodified
- Added
- Removed
-
trunk/LMDZ.GENERIC/libf/phystd/ave_stelspec.F90
r3942 r3990 1 module ave_stelspec_mod 2 3 implicit none 4 5 contains 6 1 7 subroutine ave_stelspec(STELLAR) 2 8 … … 38 44 integer ifine,band 39 45 40 real,allocatable,save :: lam(:),stel_f(:) !read by master 46 real,allocatable,save :: lam(:),stel_f(:) ! read by master thread 47 ! but used by all threads 41 48 real lamm,lamp 42 49 real dl … … 50 57 51 58 integer :: ios ! file opening/reading status 59 logical :: file_exists 52 60 53 61 STELLAR(:)=0.0 … … 93 101 write(*,*) 'Please use ',1,' and only ',1,' header line in ',trim(stelspec_file) 94 102 95 ! Opening file103 ! Check the target file is there 96 104 file_path = trim(datadir)//'/stellar_spectra/'//stelspec_file 97 105 print*, 'stellar flux : ', file_path 98 OPEN(UNIT=110,FILE=file_path,STATUS='old',iostat=ios)106 inquire(FILE=file_path,EXIST=file_exists) 99 107 100 if ( ios /= 0) THEN108 if (.not.file_exists) THEN 101 109 write(*,*)'Beware that startype is now deprecated, you should use ' 102 110 write(*,*)'stelspec_file and tstellar to define the input stellar spectrum.' … … 116 124 end if 117 125 126 !$OMP MASTER 127 ! Open the file 128 OPEN(UNIT=110,FILE=file_path,STATUS='old',iostat=ios) 118 129 ! Get number of line in the file 119 130 READ(110,*) ! skip first line header just in case … … 127 138 READ(110,*) ! skip first line header just in case 128 139 129 !$OMP MASTER130 140 allocate(lam(Nfine),stel_f(Nfine)) 131 141 … … 158 168 !$OMP BARRIER 159 169 !$OMP MASTER 160 if (allocated(lam))deallocate(lam)161 if (allocated(stel_f))deallocate(stel_f)170 deallocate(lam) 171 deallocate(stel_f) 162 172 !$OMP END MASTER 163 173 !$OMP BARRIER … … 165 175 166 176 end subroutine ave_stelspec 177 178 end module ave_stelspec_mod -
trunk/LMDZ.GENERIC/libf/phystd/physiq_mod.F90
r3928 r3990 24 24 use radcommon_h, only: sigma, glat, grav, BWNV, WNOI, DWNI, DWNV, WNOV 25 25 use suaer_corrk_mod, only: suaer_corrk 26 use setspv_mod, only: setspv 26 27 use radii_mod, only: h2o_reffrad, co2_reffrad 27 28 use aerosol_mod, only: iniaerosol, iaero_co2, iaero_h2o -
trunk/LMDZ.GENERIC/libf/phystd/setspv.F90
r3893 r3990 1 module setspv_mod 2 3 implicit none 4 5 contains 6 1 7 subroutine setspv 2 8 … … 27 33 use datafile_mod, only: datadir 28 34 use callkeys_mod, only: Fat1AU 35 use ave_stelspec_mod, only: ave_stelspec 29 36 30 37 implicit none … … 132 139 print*,' ' 133 140 134 RETURN135 141 END subroutine setspv 142 143 end module setspv_mod 144
Note: See TracChangeset
for help on using the changeset viewer.
