Changeset 3485 for trunk/LMDZ.GENERIC


Ignore:
Timestamp:
Oct 24, 2024, 1:10:52 PM (2 months ago)
Author:
afalco
Message:

Generic PCM: abort_physic() in soil.F instead of stop().
AF

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.GENERIC/libf/phystd/soil.F

    r1543 r3485  
    1616!
    1717!  Purpose: Compute soil temperature using an implict 1st order scheme
    18 ! 
    19 !  Note: depths of layers and mid-layers, soil thermal inertia and 
     18!
     19!  Note: depths of layers and mid-layers, soil thermal inertia and
    2020!        heat capacity are commons in comsoil.h
    2121!-----------------------------------------------------------------------
     
    2525!  ---------
    2626!  inputs:
    27       integer,intent(in) :: ngrid       ! number of (horizontal) grid-points 
    28       integer,intent(in) :: nsoil       ! number of soil layers 
    29       logical,intent(in) :: firstcall ! identifier for initialization call 
     27      integer,intent(in) :: ngrid       ! number of (horizontal) grid-points
     28      integer,intent(in) :: nsoil       ! number of soil layers
     29      logical,intent(in) :: firstcall ! identifier for initialization call
    3030      logical,intent(in) :: lastcall
    3131      real,intent(in) :: therm_i(ngrid,nsoil) ! thermal inertia
     
    4646      real,save :: mu
    4747!$OMP THREADPRIVATE(mthermdiff,thermdiff,coefq,coefd,alph,beta,mu)
    48            
     48
    4949! local variables:
    5050      integer ig,ik
     
    5656      if (firstcall) then
    5757      ! note: firstcall is set to .true. or .false. by the caller
    58       !       and not changed by soil.F 
     58      !       and not changed by soil.F
    5959
    6060      ALLOCATE(mthermdiff(ngrid,0:nsoil-1)) ! mid-layer thermal diffusivity
     
    100100          coefd(ig,ik)=thermdiff(ig,ik)/(mlayer(ik)-mlayer(ik-1))
    101101        enddo
    102        
     102
    103103        ! alph_{N-1}
    104104        alph(ig,nsoil-1)=coefd(ig,nsoil-1)/
     
    120120      !write(*,*)'soil: ig=',ig,' capcal(ig)=',capcal(ig)
    121121      enddo ! of do ig=1,ngrid
    122      
     122
    123123      ! Additional checks: is the vertical discretization sufficient
    124124      ! to resolve diurnal and annual waves?
     
    139139     &                diurnal_skin
    140140          write(*,*) " change soil layer distribution (comsoil_h.F90)"
    141           stop
     141          call abort_physic("soil",
     142     &          "change soil layer distribution (comsoil_h.F90)",1)
    142143        endif
    143144        if (2.*annual_skin>layer(nsoil)) then
     
    150151     &                annual_skin
    151152          write(*,*) " change soil layer distribution (comsoil_h.F90)"
    152           stop
     153          call abort_physic("soil",
     154     &          "change soil layer distribution (comsoil_h.F90)",1)
    153155        endif
    154156      enddo ! of do ig=1,ngrid
    155      
     157
    156158      else ! of if (firstcall)
    157159
     
    171173        enddo
    172174      enddo
    173      
     175
    174176      endif! of if (firstcall)
    175177
Note: See TracChangeset for help on using the changeset viewer.