Changeset 4065 for trunk/LMDZ.COMMON/libf/evolution/frost.F90
- Timestamp:
- Feb 12, 2026, 9:09:12 AM (2 weeks ago)
- File:
-
- 1 moved
-
trunk/LMDZ.COMMON/libf/evolution/frost.F90 (moved) (moved from trunk/LMDZ.COMMON/libf/evolution/metamorphism.F90) (6 diffs)
Legend:
- Unmodified
- Added
- Removed
-
trunk/LMDZ.COMMON/libf/evolution/frost.F90
r4064 r4065 1 MODULE metamorphism2 !----------------------------------------------------------------------- 3 ! NAME 4 ! metamorphism1 MODULE frost 2 !----------------------------------------------------------------------- 3 ! NAME 4 ! frost 5 5 ! 6 6 ! DESCRIPTION … … 14 14 !----------------------------------------------------------------------- 15 15 16 ! DECLARATION 17 ! ----------- 18 implicit none 19 20 ! MODULE VARIABLES 21 ! ---------------- 16 ! DEPENDENCIES 17 ! ------------ 18 use numerics, only: dp, di 19 20 ! DECLARATION 21 ! ----------- 22 implicit none 23 24 ! VARIABLES 25 ! --------- 22 26 ! Different types of frost retained by the PEM to give back to the PCM at the end 23 real, dimension(:,:), allocatable :: h2o_frost4PCM 24 real, dimension(:,:), allocatable :: co2_frost4PCM 25 26 ! Indices for frost taken from the PCM 27 integer :: iPCM_h2ofrost, iPCM_co2frost 27 real(dp), dimension(:,:), allocatable :: h2o_frost4PCM 28 real(dp), dimension(:,:), allocatable :: co2_frost4PCM 29 30 ! PARAMETERS 31 ! ---------- 32 ! Different types of frost in the PCM at the beginning [Pa] 33 real(dp), dimension(:,:), allocatable, protected :: h2ofrost_PCM 34 real(dp), dimension(:,:), allocatable, protected :: co2frost_PCM 28 35 29 36 contains … … 31 38 32 39 !======================================================================= 33 SUBROUTINE ini_frost_id(nqtot,noms) 34 !----------------------------------------------------------------------- 35 ! NAME 36 ! ini_frost_id 37 ! 38 ! DESCRIPTION 39 ! Initialize frost indices from PCM variable names. 40 SUBROUTINE ini_frost() 41 !----------------------------------------------------------------------- 42 ! NAME 43 ! ini_frost 44 ! 45 ! DESCRIPTION 46 ! Initialize the parameters of module 'frost'. 47 ! 48 ! AUTHORS & DATE 49 ! JB Clement, 12/2025 50 ! 51 ! NOTES 52 ! 53 !----------------------------------------------------------------------- 54 55 ! DEPENDENCIES 56 ! ------------ 57 use geometry, only: ngrid, nslope 58 59 ! DECLARATION 60 ! ----------- 61 implicit none 62 63 ! CODE 64 ! ---- 65 if (.not. allocated(h2ofrost_PCM)) allocate(h2ofrost_PCM(ngrid,nslope)) 66 if (.not. allocated(co2frost_PCM)) allocate(co2frost_PCM(ngrid,nslope)) 67 if (.not. allocated(h2o_frost4PCM)) allocate(h2o_frost4PCM(ngrid,nslope)) 68 if (.not. allocated(co2_frost4PCM)) allocate(co2_frost4PCM(ngrid,nslope)) 69 70 END SUBROUTINE ini_frost 71 !======================================================================= 72 73 !======================================================================= 74 SUBROUTINE end_frost() 75 !----------------------------------------------------------------------- 76 ! NAME 77 ! end_frost 78 ! 79 ! DESCRIPTION 80 ! Deallocate frost arrays. 81 ! 82 ! AUTHORS & DATE 83 ! JB Clement, 12/2025 84 ! 85 ! NOTES 86 ! 87 !----------------------------------------------------------------------- 88 89 ! DECLARATION 90 ! ----------- 91 implicit none 92 93 ! CODE 94 ! ---- 95 if (allocated(h2ofrost_PCM)) deallocate(h2ofrost_PCM) 96 if (allocated(co2frost_PCM)) deallocate(co2frost_PCM) 97 if (allocated(h2o_frost4PCM)) deallocate(h2o_frost4PCM) 98 if (allocated(co2_frost4PCM)) deallocate(co2_frost4PCM) 99 100 END SUBROUTINE end_frost 101 !======================================================================= 102 103 !======================================================================= 104 SUBROUTINE set_h2ofrost_PCM(h2ofrost_PCM_in) 105 !----------------------------------------------------------------------- 106 ! NAME 107 ! set_h2ofrost_PCM 108 ! 109 ! DESCRIPTION 110 ! Setter for 'h2ofrost_PCM'. 40 111 ! 41 112 ! AUTHORS & DATE … … 52 123 ! ARGUMENTS 53 124 ! --------- 54 integer, intent(in) :: nqtot 55 character(*), dimension(nqtot), intent(in) :: noms 56 57 ! LOCAL VARIABLES 58 ! --------------- 59 integer :: i 60 61 ! CODE 62 ! ---- 63 ! Initialization 64 iPCM_h2ofrost = -1 65 iPCM_co2frost = -1 66 67 ! Getting the index 68 do i = 1,nqtot 69 if (trim(noms(i)) == "h2o_ice") iPCM_h2ofrost = i 70 if (trim(noms(i)) == "co2") iPCM_co2frost = i 71 enddo 72 73 ! Checking if everything has been found 74 if (iPCM_h2ofrost < 0) error stop 'ini_frost_id: H2O frost index not found!' 75 if (iPCM_co2frost < 0) error stop 'ini_frost_id: CO2 frost index not found!' 76 77 END SUBROUTINE ini_frost_id 78 !======================================================================= 79 80 !======================================================================= 81 SUBROUTINE compute_frost(ngrid,nslope,h2ofrost_PCM,min_h2ofrost,co2frost_PCM,min_co2frost) 82 !----------------------------------------------------------------------- 83 ! NAME 84 ! compute_frost 85 ! 86 ! DESCRIPTION 87 ! Compute the frost to give back to the PCM. 125 real(dp), dimension(:,:), intent(in) :: h2ofrost_PCM_in 126 127 ! CODE 128 ! ---- 129 h2ofrost_PCM(:,:) = h2ofrost_PCM_in(:,:) 130 131 END SUBROUTINE set_h2ofrost_PCM 132 !======================================================================= 133 134 !======================================================================= 135 SUBROUTINE set_co2frost_PCM(co2frost_PCM_in) 136 !----------------------------------------------------------------------- 137 ! NAME 138 ! set_co2frost_PCM 139 ! 140 ! DESCRIPTION 141 ! Setter for 'co2frost_PCM'. 142 ! 143 ! AUTHORS & DATE 144 ! JB Clement, 12/2025 145 ! 146 ! NOTES 147 ! 148 !----------------------------------------------------------------------- 149 150 ! DECLARATION 151 ! ----------- 152 implicit none 153 154 ! ARGUMENTS 155 ! --------- 156 real(dp), dimension(:,:), intent(in) :: co2frost_PCM_in 157 158 ! CODE 159 ! ---- 160 co2frost_PCM(:,:) = co2frost_PCM_in(:,:) 161 162 END SUBROUTINE set_co2frost_PCM 163 !======================================================================= 164 165 !======================================================================= 166 SUBROUTINE compute_frost4PCM(min_h2ofrost,min_co2frost) 167 !----------------------------------------------------------------------- 168 ! NAME 169 ! compute_frost4PCM 170 ! 171 ! DESCRIPTION 172 ! Compute the frost to give back to the PCM (metamorphism). 88 173 ! 89 174 ! AUTHORS & DATE … … 95 180 !----------------------------------------------------------------------- 96 181 182 ! DEPENDENCIES 183 ! ------------ 184 use display, only: print_msg 185 97 186 ! DECLARATION 98 187 ! ----------- … … 101 190 ! ARGUMENTS 102 191 ! --------- 103 integer, intent(in) :: ngrid, nslope 104 real, dimension(ngrid,nslope), intent(in) :: h2ofrost_PCM, min_h2ofrost, co2frost_PCM, min_co2frost 105 106 ! CODE 107 ! ---- 108 write(*,*) '> Computing frost to give back to the PCM' 109 110 ! Allocation 111 call ini_frost(ngrid,nslope) 112 113 ! Initialization 114 h2o_frost4PCM(:,:) = 0. 115 co2_frost4PCM(:,:) = 0. 116 117 ! Computation: frost for the PEM is the extra part of the PCM frost above the yearly minimum 118 where (h2ofrost_PCM(:,:) > 0.) h2o_frost4PCM(:,:) = h2ofrost_PCM(:,:) - min_h2ofrost(:,:) 119 where (co2frost_PCM(:,:) > 0.) co2_frost4PCM(:,:) = co2frost_PCM(:,:) - min_co2frost(:,:) 120 121 END SUBROUTINE compute_frost 122 !======================================================================= 123 124 !======================================================================= 125 SUBROUTINE set_frost4PCM(PCMfrost) 126 !----------------------------------------------------------------------- 127 ! NAME 128 ! set_frost4PCM 129 ! 130 ! DESCRIPTION 131 ! Reconstruct frost for the PCM from PEM computations. 132 ! 133 ! AUTHORS & DATE 134 ! JB Clement, 12/2025 135 ! 136 ! NOTES 137 ! 138 !----------------------------------------------------------------------- 139 140 ! DECLARATION 141 ! ----------- 142 implicit none 143 144 ! ARGUMENTS 145 ! --------- 146 real, dimension(:,:,:), intent(inout) :: PCMfrost 147 148 ! CODE 149 ! ---- 150 write(*,*) '> Reconstructing frost for the PCM' 151 PCMfrost(:,iPCM_h2ofrost,:) = h2o_frost4PCM(:,:) 152 PCMfrost(:,iPCM_co2frost,:) = co2_frost4PCM(:,:) 153 154 ! Deallocation 155 call end_frost() 156 157 END SUBROUTINE set_frost4PCM 158 !======================================================================= 159 160 !======================================================================= 161 SUBROUTINE ini_frost(ngrid,nslope) 162 !----------------------------------------------------------------------- 163 ! NAME 164 ! ini_frost 165 ! 166 ! DESCRIPTION 167 ! Initialize frost arrays. 168 ! 169 ! AUTHORS & DATE 170 ! JB Clement, 12/2025 171 ! 172 ! NOTES 173 ! 174 !----------------------------------------------------------------------- 175 176 ! DECLARATION 177 ! ----------- 178 implicit none 179 180 ! ARGUMENTS 181 ! --------- 182 integer, intent(in) :: ngrid, nslope 183 184 ! CODE 185 ! ---- 186 if (.not. allocated(h2o_frost4PCM)) allocate(h2o_frost4PCM(ngrid,nslope)) 187 if (.not. allocated(co2_frost4PCM)) allocate(co2_frost4PCM(ngrid,nslope)) 188 189 END SUBROUTINE ini_frost 190 !======================================================================= 191 192 !======================================================================= 193 SUBROUTINE end_frost() 194 !----------------------------------------------------------------------- 195 ! NAME 196 ! end_frost 197 ! 198 ! DESCRIPTION 199 ! Deallocate frost arrays. 200 ! 201 ! AUTHORS & DATE 202 ! JB Clement, 12/2025 203 ! 204 ! NOTES 205 ! 206 !----------------------------------------------------------------------- 207 208 ! DECLARATION 209 ! ----------- 210 implicit none 211 212 ! CODE 213 ! ---- 214 if (allocated(h2o_frost4PCM)) deallocate(h2o_frost4PCM) 215 if (allocated(co2_frost4PCM)) deallocate(co2_frost4PCM) 216 217 END SUBROUTINE end_frost 218 !======================================================================= 219 220 END MODULE metamorphism 192 real(dp), dimension(:,:), intent(in) :: min_h2ofrost, min_co2frost 193 194 ! CODE 195 ! ---- 196 call print_msg('> Computing frost to give back to the PCM (metamorphism)') 197 198 h2o_frost4PCM(:,:) = 0._dp 199 co2_frost4PCM(:,:) = 0._dp 200 where (h2ofrost_PCM(:,:) > 0._dp) h2o_frost4PCM(:,:) = h2ofrost_PCM(:,:) - min_h2ofrost(:,:) 201 where (co2frost_PCM(:,:) > 0._dp) co2_frost4PCM(:,:) = co2frost_PCM(:,:) - min_co2frost(:,:) 202 203 END SUBROUTINE compute_frost4PCM 204 !======================================================================= 205 206 END MODULE frost
Note: See TracChangeset
for help on using the changeset viewer.
