Changeset 442 for branches


Ignore:
Timestamp:
07/03/23 15:32:49 (11 months ago)
Author:
aquiquet
Message:

A few modifications on bmelt-beckmann-gcm, K_t and nbassins are in param_list now

File:
1 edited

Legend:

Unmodified
Added
Removed
  • branches/GRISLIv3/SOURCES/bmelt-beckmann-gcm_mod.f90

    r346 r442  
    2929  implicit none 
    3030 
    31   integer, parameter :: nbassins = 10 !< number of sectors in the ocean 
     31  real*8,parameter :: cpw    = 4182.d0 
     32  integer :: nbassins !< number of sectors in the ocean 
    3233  integer :: nzoc !< number of vertical levels in the ocean (read in netcdf T/S file) 
    3334  real*8 :: coef_OM !combined coefficient by DeConto et Pollard (m/yr/C2) 
     
    4546  real*8, dimension (nx,ny)      :: ice_draft       !< ice draft (S-H) 
    4647  character(len=200) :: TO_file, SO_file, Bassin_file   ! fichiers de forcage 
    47   logical, dimension(nbassins) :: mask_bassin_vide  ! mask true if bassin without any GCM data 
     48  logical, dimension(:), allocatable :: mask_bassin_vide  ! mask true if bassin without any GCM data 
    4849   
    4950contains 
     
    7071    integer, dimension(:,:,:), allocatable :: mask_gcm  ! mask were T S are defined 
    7172    integer, dimension(:,:,:), allocatable :: mask_gcm_init ! mask were T S are defined in GCM before interpolation 
    72     integer, dimension(nbassins) :: max_depth_bassin ! max level were T S are defined for every bassins 
     73    integer, dimension(:), allocatable :: max_depth_bassin ! max level were T S are defined for every bassins 
    7374    integer :: nmax ! to avoid inifite loop  
    7475    integer :: err 
    75     namelist/bmelt_beckmann_gcm_mod/TO_file,SO_file,Bassin_file,bmelt_empty 
     76    namelist/bmelt_beckmann_gcm_mod/K_t,TO_file,SO_file,Bassin_file,nbassins,bmelt_empty 
    7677     
    7778    rewind(num_param)        ! loop back at the beginning of the param_list.dat file 
     
    8182    write(num_rep_42,bmelt_beckmann_gcm_mod) 
    8283 
    83     K_t = 15.77 
     84    !K_t = 15.77 
    8485 
    8586    file_inputs=TRIM(DIRNAMEINP)//'Snapshots-GCM/'//TRIM(TO_file) 
     
    9798       end if 
    9899    end if 
     100    if (.not.allocated(mask_bassin_vide)) then 
     101      allocate(mask_bassin_vide(nbassins),stat=err) 
     102      if (err/=0) then 
     103          print *,"erreur a l'allocation du tableau mask_bassin_vide",err 
     104          stop 4 
     105       end if 
     106    end if 
    99107    if (.not.allocated(temp_ocean)) then 
    100108      allocate(temp_ocean(nx,ny,nzoc),stat=err) 
     
    136144      if (err/=0) then 
    137145          print *,"erreur a l'allocation du tableau mask_gcm_init",err 
     146          stop 4 
     147       end if 
     148    end if 
     149    if (.not.allocated(max_depth_bassin)) then 
     150      allocate(max_depth_bassin(nbassins),stat=err) 
     151      if (err/=0) then 
     152          print *,"erreur a l'allocation du tableau max_depth_bassin",err 
    138153          stop 4 
    139154       end if 
     
    178193         
    179194   depth_max_grid=maxval(depth_max(:,:)) ! max depth were T and S are defined in the grid 
    180     
     195 
    181196   !$OMP PARALLEL 
    182197   if (depth_max_grid.lt.nzoc) then 
     
    310325        !temp_ocean(:,:,:) = temp_ocean(:,:,:) - 273.15 
    311326         
    312     coef_OM = K_t * 0.01420418516 
     327    !coef_OM = K_t * 0.01420418516 
     328    coef_OM = K_t * row * cpw / (cl * ro) 
    313329 
    314330    if ( ubound(zoc,dim=1).ne.nzoc) then 
     
    394410       do i=1,nx 
    395411          if ( TO_draft(i,j) .gt. -9000.d0 ) then ! floating points 
    396             if (mask_bassin_vide(bassin(i,j)) .or. (bassin(i,j).eq.1)) then ! bassins without any GCM data or bassin 1 
     412            !if (mask_bassin_vide(bassin(i,j)) .or. (bassin(i,j).eq.1)) then ! bassins without any GCM data or bassin 1 
     413            !if (mask_bassin_vide(bassin(i,j)) .or. (bassin(i,j).eq.1)) then !bassins without any GCM data 
     414            if (mask_bassin_vide(bassin(i,j))) then !bassins without any GCM data 
    397415               bmelt(i,j) = bmelt_empty 
    398416            else  ! Bassins with GCM data 
     
    412430             ngr=0 
    413431             bmloc=0.d0 
    414              if (TO_draft(i+1,j) .le. -9000.d0 ) then 
     432             !if (TO_draft(i+1,j) .le. -9000.d0 ) then 
     433             if (flot(i+1,j)) then 
    415434                ngr=ngr+1 
    416435                bmloc=bmloc+bmelt(i+1,j) 
    417436             endif 
    418              if (TO_draft(i-1,j) .le. -9000.d0 ) then 
     437             !if (TO_draft(i-1,j) .le. -9000.d0 ) then 
     438             if (flot(i-1,j)) then 
    419439                ngr=ngr+1 
    420440                bmloc=bmloc+bmelt(i-1,j) 
    421441             endif 
    422              if (TO_draft(i,j+1) .le. -9000.d0 ) then 
     442             !if (TO_draft(i,j+1) .le. -9000.d0 ) then 
     443             if (flot(i,j+1)) then 
    423444                ngr=ngr+1 
    424445                bmloc=bmloc+bmelt(i,j+1) 
    425446             endif 
    426              if (TO_draft(i,j-1) .le. -9000.d0 ) then 
     447             !if (TO_draft(i,j-1) .le. -9000.d0 ) then 
     448             if (flot(i,j-1)) then 
    427449                ngr=ngr+1 
    428450                bmloc=bmloc+bmelt(i,j-1) 
    429451             endif 
     452             ngr=0   ! afq -- I do not allow subgrid ocean melt 
    430453             bmelt(i,j)=(ngr/4.)*bmloc+(1-ngr/4.)*bmelt(i,j)  
    431454          endif 
Note: See TracChangeset for help on using the changeset viewer.