source: trunk/LMDZ.MARS/libf/phymars/paleoclimate_mod.F90 @ 3625

Last change on this file since 3625 was 3625, checked in by emillour, 2 days ago

Mars PCM:
OpenMP fix, zdqsdif_ssi_tot is a saved variable and has to be THREADPRIVATE.
EM

File size: 2.1 KB
Line 
1MODULE paleoclimate_mod
2!=======================================================================
3!   subject: Module dedicated to paleoclimates studies
4!   --------
5!
6!   author: LL, 06/2023
7!   ------
8!
9!=======================================================================
10
11implicit none
12
13logical :: paleoclimate ! False by default, is activate  for paleoclimates specific processes (e.g., lag layer) is initialized in conf_phys
14!$OMP THREADPRIVATE(paleoclimate)
15
16    real,    allocatable, dimension(:,:) :: h2o_ice_depth         ! Thickness of the lag before H2O ice [m]
17    real,    allocatable, dimension(:,:) :: lag_co2_ice           ! Thickness of the lag before CO2 ice [m]
18    real,    allocatable, dimension(:,:) :: zdqsdif_ssi_tot       ! Total flux of the interactions with SSI (kg/m^2/s^-1)
19!$OMP THREADPRIVATE(h2o_ice_depth,lag_co2_ice,zdqsdif_ssi_tot)
20    real,    allocatable, dimension(:,:) :: d_coef                ! Diffusion coefficent
21    logical                              :: lag_layer             ! Does lag layer is present?
22    logical                              :: include_waterbuoyancy ! Include the effect of water buoyancy when computing the sublimation of water ice ?
23!$OMP THREADPRIVATE(d_coef,lag_layer,include_waterbuoyancy)
24
25!=======================================================================
26contains
27!=======================================================================
28
29SUBROUTINE ini_paleoclimate_h(ngrid,nslope)
30
31implicit none
32
33integer, intent(in) :: ngrid  ! number of atmospheric columns
34integer, intent(in) :: nslope ! number of slope within a mesh
35
36allocate(h2o_ice_depth(ngrid,nslope))
37allocate(lag_co2_ice(ngrid,nslope))
38allocate(zdqsdif_ssi_tot(ngrid,nslope))
39allocate(d_coef(ngrid,nslope))
40
41END SUBROUTINE ini_paleoclimate_h
42
43!=======================================================================
44SUBROUTINE end_paleoclimate_h
45
46implicit none
47
48if (allocated(d_coef)) deallocate(d_coef)
49if (allocated(h2o_ice_depth)) deallocate(h2o_ice_depth)
50if (allocated(zdqsdif_ssi_tot)) deallocate(zdqsdif_ssi_tot)
51if (allocated(lag_co2_ice)) deallocate(lag_co2_ice)
52
53END SUBROUTINE end_paleoclimate_h
54
55END MODULE paleoclimate_mod
Note: See TracBrowser for help on using the repository browser.