Skip to content
Snippets Groups Projects
soil.f 38.4 KiB
Newer Older
Petra Lasch-Born's avatar
Petra Lasch-Born committed
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 675 676 677 678 679 680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 719 720 721 722 723 724 725 726 727 728 729 730 731 732 733 734 735 736 737 738 739 740 741 742 743 744 745 746 747 748 749 750 751 752 753 754 755 756 757 758 759 760 761 762 763 764 765 766 767 768 769 770 771 772 773 774 775 776 777 778 779 780 781 782 783 784 785 786 787 788 789 790 791 792 793 794 795 796 797 798 799 800 801 802 803 804 805 806 807 808 809 810 811 812 813 814 815 816 817 818 819 820 821 822 823 824 825 826 827 828 829 830 831 832 833 834 835 836 837 838 839 840 841 842 843 844 845 846 847 848 849 850 851 852 853 854 855 856 857 858 859 860 861 862 863 864 865 866 867 868 869 870 871 872 873 874 875 876 877 878 879 880 881 882 883 884 885 886 887 888 889 890 891 892 893 894 895 896 897 898 899 900 901 902 903 904 905 906 907 908 909 910 911 912 913 914 915 916 917 918 919 920 921 922 923 924 925 926 927 928 929 930 931 932 933 934 935 936 937 938 939 940 941 942 943 944 945 946 947 948 949 950 951 952 953 954 955 956 957 958 959 960 961 962 963 964 965 966 967 968 969 970 971 972 973 974 975 976 977 978 979 980 981 982 983 984 985 986 987 988 989 990 991 992 993 994 995 996 997 998 999 1000
!*****************************************************************!
!*                                                               *!
!*              4C (FORESEE) Simulation Model                    *!
!*                                                               *!
!*                                                               *!
!*                    Subroutines for:                           *!
!*                Soil and Water - Programs                      *!
!*                                                               *!
!*   contains:                                                   *!
!*   SOIL                                                        *!
!*   SOIL_INI                                                    *!    
!*   SOIL_WAT                                                    *!
!*   UPT_WAT                                                     *!
!*   FRED1 - ...11                                                *!
!*   TAKE_WAT                                                    *!
!*   BUCKET                                                      *!
!*   SNOWPACK                                                    *!
!*   HUM_ADD                                                     *!
!*   BC_APPL: application of biochar                             *!
!*                                                               *!
!*                  Copyright (C) 1996-2018                      *!
!*     Potsdam Institute for Climate Impact Reserach (PIK)       *!
!*          Authors and contributors see AUTHOR file             *!
!*  This file is part of 4C and is licensed under BSD-2-Clause   *!
!*                   See LICENSE file or under:                  *!
!*     http://www.https://opensource.org/licenses/BSD-2-Clause   *!
!*                           Contact:                            *!
!*       https://gitlab.pik-potsdam.de/foresee/4C                *!
!*                                                               *!
!*****************************************************************!
    
SUBROUTINE soil

!   Soil processes (frame)

use data_climate
use data_depo
use data_out 
use data_simul
use data_soil
use data_soil_cn

implicit none

call evapo
call intercep
call soil_wat
call soil_temp

if (flag_wurz .eq. 4 .or. flag_wurz .eq. 6) then
    call soil_stress	!calculate ground stress factors
    call cr_depth		!define root depth 
endif
call soil_cn
call root_eff

END subroutine soil

!**************************************************************

SUBROUTINE soil_ini

! Initialisation of soil data and parameters

use data_inter
use data_evapo
use data_out
use data_par
use data_simul
use data_soil
use data_soil_cn
use data_species
use data_stand

implicit none

integer i,j,k
real d_0, t_0
! Table of quarz and clay content (mass%) versus wlam
real, dimension(17) :: xwlam  = (/1.5, 1.15, 0.9, 0.67, 0.6,  0.5, 0.38, 0.37, 0.3, 0.29, 0.27, 0.26, 0.25, 0.24, 0.23, 0.22, 0.15/),  &
                       yquarz = (/93.0,85.0,80.0, 82.0,76.0, 64.0, 65.0, 51.0,60.0, 30.0, 14.0, 10.0, 12.0, 20.0, 30.0, 43.0, 23.0/),  &
                       yclay  = (/3.0, 3.0,  3.0, 12.0, 6.0,  6.0, 10.0,  4.0,21.0, 12.0, 10.0, 37.0, 15.0, 40.0, 30.0, 35.0, 55.0/) 
real value
real, dimension(nlay):: humush(nlay)
real, dimension(nlay):: xfcap, xwiltp, xpv   ! output of addition for water capacities

! estimation of soil layer values
d_0 = 0.
do j = 1, nlay
   t_0      = thick(j)
   mid(j)   = d_0 + 0.5*t_0
   d_0      = d_0 + t_0
   depth(j) = d_0
enddo

perc	= 0.
wupt_r	= 0.
wupt_ev	= 0.
thick_1 = thick(1)

select case (flag_soilin)
case (0,2)
    do i=1,nlay
        if (i .gt. 1) then
            call tab_int(xwlam, yquarz,  17, wlam(i), value)
            sandv(i)  = value / 100.   ! Mass% of mineral fraction
            call tab_int(xwlam, yclay,  17, wlam(i), value)
            clayv(i)  = value / 100.
            siltv(i)  = 1. - clayv(i) - sandv(i)
        else
            sandv(1)  = 0.0
            clayv(1)  = 0.0
            siltv(1)  = 0.0
        endif
    enddo

case (1,3,4)
    clayv  = clayv / 100.
    sandv  = sandv / 100.
    siltv  = 1. - clayv - sandv
    if ((sandv(1) .le. zero) .and. (clayv(1) .le. zero)) siltv(1) = 0.
    skelv  = skelv / 100.
    humusv = humusv / 100.
end select

! Settings for subroutine take_wat
skelfact  = 1.
pv        = skelfact * pv_v * thick * 0.1	       ! mm
wilt_p    = skelfact * wilt_p_v * thick * 0.1	   ! mm
field_cap = skelfact * f_cap_v * thick * 0.1	   ! mm
wats      = field_cap   ! mm
watvol    = f_cap_v                    ! vol%

n_ev_d = nlay
nlgrw  = nlay+1
do i=1,nlay
	if (w_ev_d .gt. depth(i)) n_ev_d = i
	if (grwlev .gt. depth(i)) nlgrw = i+1
    vol(i)    = thick(i) * 10000.  
enddo

! dry mass of first layer
dmass = vol * dens
rmass1 = dmass(1) - (C_hum(1) + C_opm(1)) / cpart    ! corection term of first layer

humush = humusv
if (2.*C_hum(1) .lt. humusv(1)*dmass(1)) then
    humusv(1) = C_hum(1) / (dmass(1) * cpart)
endif
do i=2, nlay
    humusv(i) = C_hum(i) / (dmass(i) * cpart)
enddo

if (flag_bc .gt. 0) y_bc_n = 1      ! actual number of biochar application

! calculation of additions for water capacities 
call hum_add(xfcap, xwiltp, xpv)

fcaph  = f_cap_v - xfcap
wiltph = wilt_p_v - xwiltp
pvh    = pv_v - xpv

! ground water
do i = nlgrw, nlay
    wats(i) = pv(i)
enddo

interc_can  = 0.
int_st_can  = 0.
interc_sveg = 0.
int_st_sveg = 0.
aev_i       = 0.




wat_tot = SUM(wats)
END subroutine soil_ini

!**************************************************************

SUBROUTINE soil_wat
use data_out
! soil water balance

use data_climate
use data_evapo
use data_inter
use data_out
use data_par 
use data_simul
use data_soil
use data_species
use data_stand

implicit none
 
real :: eva_dem		! evaporation demand of soil
real :: p_inf = 0.  ! infiltrated water
real :: pev		    ! local: soil evaporation
real :: watot, wetot    ! total water content at start and end
real :: wutot           ! total water uptake from the soil
real :: wutot_ev        ! total water uptake by soil evaporation
real :: wutot_r         ! total water uptake by roots
real, allocatable, dimension(:) :: upt   ! local array for uptake
real, external :: wat_new

real enr, wa, we, percolnl, snow_sm, buckdepth
integer j

allocate (upt(nlay))
wupt_ev = 0.
aev_s = 0.
prec_stand = MAX(0., prec - interc_can - interc_sveg)    ! stand precipitation	

if (flag_int .gt. 1000) then
    prec_stand = prec_stand * prec_stand_red / 100.
endif

call snowpack(snow_sm, p_inf, pev)

if (anz_coh .le. 0) pev = pet
eva_dem = MAX(0., p_inf) - pev		! evaporation demand of soil

!     if all stand precipitation is evaporated and there is still a demand
!     there is an uptake from soil layers (up to an certain depth)
if (eva_dem .lt. 0.) then
   if (snow .le. 0) call take_wat(eva_dem, cover) 
   aev_s = aev_s + p_inf + SUM(wupt_ev)
else
   aev_s = aev_s + pev
endif	! eva_dem

aet   = aev_s + aev_i
p_inf = MAX(eva_dem, 0.)
upt   = wupt_ev
watot = SUM(wats)			! total initial water content

do j = 1, nlgrw-1

  enr   = p_inf - upt(j)
  wa    = wats(j) - field_cap(j)
  we    = wat_new(wa, enr, j)
  p_inf = enr + wa - we

  perc(j) = MAX(p_inf, 0.)
  wats(j) = MAX(we+field_cap(j), wilt_p(j))

enddo

do j = nlgrw, nlay

  enr   = p_inf - upt(j)
  wa    = wats(j) - field_cap(j)
  we    = pv(j) - field_cap(j)   ! ground water level is constant!
  p_inf = enr + wa - we
  perc(j) = MAX(p_inf, 0.)
  wats(j) = MAX(we+field_cap(j), wilt_p(j))

enddo

if (flag_wred .le. 10) then
   call upt_wat
else
   call upt_wat1
endif

! root uptake balanced imediate after calculation
upt  = upt + wupt_r
wats = wats - wupt_r

watvol = 10.*wats/(thick * skelfact)   ! estimation for complete layer in Vol% without skeleton (only soil substrate)

! total water quantities
wetot    = SUM(wats)		! total final water content
wutot_ev = SUM(wupt_ev)		! total water uptake by soil evaporation
wutot_r  = SUM(wupt_r)		! total water uptake by roots
wutot    = wutot_ev + wutot_r	! total water uptake
aet      = aet + wutot_r        ! daily total aet

percolnl = perc(nlay)

trans_tree = 0.
trans_sveg = 0.
   zeig => pt%first
   do while (associated(zeig))
     if (zeig%coh%species .le. nspec_tree) then
        trans_tree = trans_tree + zeig%coh%supply 
     else
        trans_sveg = trans_sveg + zeig%coh%supply 
     endif
     zeig => zeig%next
   enddo  ! zeig (cohorts) 

! cumulative water quantities
perc_cum   = perc_cum + perc(nlay)
wupt_r_c   = wupt_r_c + wutot_r
wupt_e_c   = wupt_e_c + wutot_ev
wupt_cum   = wupt_cum + wutot
aet_cum    = aet_cum + aet
dew_cum    = dew_cum + dew_rime
tra_tr_cum = tra_tr_cum + trans_tree 
tra_sv_cum = tra_sv_cum + trans_sveg

call bucket(bucks_100, bucks_root, buckdepth)

! number of drought days per layer
where ((wats-0.2) .le. wilt_p) s_drought = s_drought+1

 if (flag_dayout .ge. 2) then
   write (unit_wat, '(2I5, 7F7.2, 24F8.2)') time_cur, iday, airtemp, prec, interc_can, int_st_can, &
                    interc_sveg, int_st_sveg, snow, snow_sm,  pet, trans_dem, &
                    pev, aev_s, aev_i, perc(nlay), watot, wetot, wutot, wutot_ev, wutot_r,&
                    trans_tree,trans_sveg, eva_dem, gp_can, aet, ceppot_can, ceppot_sveg
endif
deallocate (upt)
		    
END	subroutine soil_wat

!**************************************************************
SUBROUTINE upt_wat

! Water uptake by roots

use data_simul
use data_evapo
use data_soil
use data_stand 
use data_par
use data_species
use data_climate

implicit none

real, dimension(1:anz_coh) :: tr_dem   ! auxiliary arrays for cohorts
real wat_ava, hdem, frtrel, frtrel_1, hupt, hupt_c, totfrt_2, hv, demand_mistletoe_canopy
real wat_at    ! total available water per layer
real wat_ar    ! total available water per layer with uptake resistance
real hred      ! resistance coefficient
real, external :: fred1, fred2, fred3, fred4, fred5, fred6, fred7, fred11
integer i, ianz, j, nroot3


! Calculation of Water Demand
ianz = anz_coh

tr_dem = 0
hdem   = 0
hv     = pet-aev_i-pev_s
if (hv .lt. 0.) hv = 0.

select case (flag_eva)
case (0,1,3)
 if((pet .gt. 0.) .and. (hv .gt. 0.)) then
    trans_dem = hv * alfm * (1. - exp(-gp_tot/gpmax))  ! pet (potential evapotranspiration) is reduced by intereption evaporation and potential ground evaporation
 else
 trans_dem = 0.0 
 endif 
   if (gp_tot .gt. zero) then
      hdem   = trans_dem / gp_tot
   else
      hdem= 0.
   endif
case (8)  
 ! potential transpiration demand of each cohort
   if ((gp_tot .gt. zero) .and. (hv .gt. 0.)) then
      hdem   = (pet-aev_i-aev_s) / gp_tot
   else
      hdem= 0.
   endif

case (2,4)
   trans_dem = 0.

case (6,16,36)
 ! Eucalyptus
 hv = pet

 if((pet .gt. 0.) .and. (hv .gt. 0.)) then
    trans_dem = hv * alfm * (1. - exp(-gp_tot/gpmax)) 
 else
    trans_dem = 0.
 endif 

 ! preparation: potential transpiration demand of each cohort
   if (gp_tot .gt. zero) then
      hdem   = trans_dem / gp_tot
   else
      hdem= 0.
   endif

case (7,17,37)
 trans_dem = hv

 ! potential transpiration demand of each cohort
   if (gp_tot .gt. zero) then
      hdem   = trans_dem / gp_tot
   else
      hdem= 0.
   endif

end select

hdem = max(0., hdem)

! Distribution of total Demand into Demands of Cohorts

!extraction of demand of mistletoe cohort (for case flag eva = 1,3,6,7...)
zeig => pt%first
 do while (associated(zeig))
  if (zeig%coh%species.eq.nspec_tree+2) then
      demand_mistletoe_canopy=zeig%coh%gp * zeig%coh%ntreea * hdem
  end if
  zeig => zeig%next
 enddo
zeig => pt%first
i = 1
do while (associated(zeig))

  select case (flag_eva)
  case (0, 1, 3, 6, 7, 16, 17, 36, 37)

            !uppermost tree cohort (with flag mistletoe) gets additinal demand of mistletoe
            if (zeig%coh%mistletoe.eq.1) then
                zeig%coh%demand = zeig%coh%gp * zeig%coh%ntreea * hdem + demand_mistletoe_canopy
            elseif (zeig%coh%species.eq.nspec_tree+2) then                ! set demand of mistletoe to zero as it will be fullfilled by the tree
               zeig%coh%demand=0.                                         ! set to zero because demand has been added to the infested tree cohort
            else
                zeig%coh%demand = zeig%coh%gp * zeig%coh%ntreea * hdem    ! all other cohorts get their demand
            end if

  case (2,4)
 !uppermost tree cohort (with flag mistletoe) gets additinal demand, that of mistletoe
            if (zeig%coh%mistletoe.eq.1) then
                zeig%coh%demand = (max(0., zeig%coh%demand - zeig%coh%aev_i)  + demand_mistletoe_cohort)
            endif
            if (zeig%coh%species.eq.nspec_tree+2) then                ! set demand of mistletoe to zero as it will be fullfilled by the tree
                zeig%coh%demand=0.
            endif
            if (zeig%coh%mistletoe.ne.1 .AND. zeig%coh%species.ne.nspec_tree+2) then
               zeig%coh%demand = max(0., zeig%coh%demand - zeig%coh%aev_i)   
            end if
            trans_dem = trans_dem + zeig%coh%demand
    end select

  tr_dem(i) = zeig%coh%demand   ! demand of transpiration per cohort 
   i = i + 1
  zeig => zeig%next
enddo  ! zeig (cohorts)  

! Calculation of Water Supply
frtrel_1 = 1.

select case (flag_wurz) 
case (0)
  if (nroot_max .gt. 5) then 
     nroot3 = 5
  else
     nroot3 = nroot_max
  endif
case default
  nroot3=nroot_max
end select  




! layers with seedlings   
do j = 1,nroot3
   ! determination of resisctance coefficient
   
   select case (flag_wred)
   case(1)
      hred = fred1(j)
   case(2)
      hred = fred2(j)
   case(3)
      hred = fred3(j)
   case(4)
      hred = fred4(j)
   case(5)
      hred = 1.
   case(6)
      hred = 0.5
   case(7)
      hred = 0.25
   case(8)  
      hred = fred6(j)
   case(10)
      hred = fred7(j)
   case(11)
      hred = fred11(j)
   end select
   
   wat_res(j) = hred

   if (temps(j) .gt. -0.3)  then
      wat_at = max(wats(j) - wilt_p(j), 0.)     ! total available water per layer
      wat_ar = hred * wat_at                    ! total available water per layer with uptake resistance
	  hupt = 0.
   else
      wat_ar = 0.                              ! frost
      wat_at = 0.
      hupt   = 0.
   endif



! Distribution of Water Supply into the Cohorts
! Distribution of Fine Roots
  
   zeig => pt%first
   i = 1   ! cohort index
   do while (associated(zeig))
     if (zeig%coh%species .ne. nspec_tree+2) then     ! not for mistletoe
      frtrel = zeig%coh%frtrelc(j)
     	
     wat_ava = frtrel * wat_ar        ! available water per tree cohort and layer 

    	
	 if (wat_ava .ge. tr_dem(i)) then
         hupt_c    = tr_dem(i)
         tr_dem(i) = 0.
     else
         hupt_c    = wat_ava
         tr_dem(i) = tr_dem(i) - wat_ava
     endif
	     
     select case (flag_dis)    !xylem clogger disturbance 
     case (1,2)
     hupt_c = hupt_c * xylem_dis
     end select
	 
     xwatupt(i,j) = hupt_c   ! water uptake per cohorte and layer
     zeig%coh%supply = zeig%coh%supply + hupt_c
     if (zeig%coh%supply .lt.0.) then
     continue
     endif
     hupt = hupt + hupt_c
  

	 i = i + 1
    end if              ! exclusion of mistletoe
     zeig => zeig%next
   enddo  ! zeig (cohorts) 

   wupt_r(j) = hupt 

enddo    ! j

! layers without seedlings   

if (totfrt_p.gt.(seedlfrt+zero)) then
 totfrt_2 = 1./(totfrt_p-seedlfrt)

  do j = nroot3+1, nroot_max
   ! determination of resisctance coefficient
   select case (flag_wred)
   case(1)
      hred = fred1(j)
   case(2)
      hred = fred2(j)
   case(3)
      hred = fred3(j)
   case(4)
      hred = fred4(j)
   case(5)
      hred = 1.
   case(6)
      hred = 0.5
   case(7)
      hred = 0.25
   end select
   
   wat_res(j) = hred
   
   if (temps(j) .gt. -0.3)  then
      wat_at = max(wats(j) - wilt_p(j), 0.)     ! total available water per layer
      wat_ar = hred * wat_at                    ! total available water per layer with uptake resistance
      hupt   = 0.
   else
      wat_ar = 0.
   endif
   
   zeig => pt%first
   i = 1   ! cohort index
   do while (associated(zeig))
     frtrel = zeig%coh%frtrelc(j)		
     wat_ava = frtrel * wat_ar        ! available water per tree cohort and layer 
     if (wat_ava .ge. tr_dem(i)) then
         hupt_c    = tr_dem(i)
         tr_dem(i) = 0.
     else
         hupt_c    = wat_ava
         tr_dem(i) = tr_dem(i) - wat_ava
     endif
       
     select case (flag_dis)    !xylem clogger disturbance 
     case (1,2)
     hupt_c = hupt_c * xylem_dis
     end select
	 
	 xwatupt(i,j) = hupt_c
     zeig%coh%supply = zeig%coh%supply + hupt_c
     hupt = hupt + hupt_c
     i = i + 1
     zeig => zeig%next
   enddo  ! zeig (cohorts) 

   wupt_r(j) = hupt 
enddo     ! j
endif
END	subroutine upt_wat

!**************************************************************

SUBROUTINE upt_wat1

! Water uptake by roots
! 2. Version

use data_simul
use data_evapo
use data_soil
use data_stand 
use data_par

implicit none

real, dimension(1:anz_coh) :: tr_dem,frt_rel   ! help arrays for cohorts
real wat_ava, hdem, frtrel, frtrel_1, hupt, hupt_c, totfrt_2
real wat_at    ! total available water per layer
real wat_ar    ! total available water per layer with uptake resistance
real hred      ! resistance coefficient
real, external :: fred1, fred2, fred3, fred4, fred5, fred6, fred7, fred11

integer i, ianz, j, nroot3

ianz = anz_coh
tr_dem=0
trans_dem = (pet-aev_i) * alfm * (1. - exp(-gp_can/gpmax))  ! pet NOT reduced by ground evaporation
if (trans_dem .lt. 0.) trans_dem = 0.

! potential transpiration demand of each cohort
if (gp_can .gt. zero) then
   hdem   = trans_dem / gp_can
else
   hdem= 0.
endif

! Estimation of transpiration demand of tree cohorts and total fine root mass
! in layers with and without seedlings
zeig => pt%first
i = 1
do while (associated(zeig))
  select case (flag_eva)
  case (0, 1, 3)
    zeig%coh%demand = zeig%coh%gp * zeig%coh%ntreea * hdem  
  case (2)
    zeig%coh%demand = zeig%coh%demand - zeig%coh%aev_i
  end select
  tr_dem(i) = zeig%coh%demand
  i = i + 1
  zeig => zeig%next
enddo  ! zeig (cohorts)  

! uptake controlled by share of roots 
frtrel_1 = 1.

! layers with seedlings   
do j = 1,nroot_max
   ! determination of resisctance coefficient
   select case (flag_wred)
   case(1)
      hred = fred1(j)
   case(2)
      hred = fred2(j)
   case(3)
      hred = fred3(j)
   case(4)
      hred = fred4(j)
   case(5)
      hred = 1.
   case(6)
      hred = 0.5
   case(7)
      hred = 0.25
   case(8)  ! BKL, ArcEGMO
      hred = fred6(j)
   case(10)
      hred = fred7(j)
   end select

   wat_at = max(wats(j) - wilt_p(j), 0.)     ! total available water per layer
   wat_ar = hred * wat_at                    ! total available water per layer with uptake resistance
   hupt   = 0.
   
   zeig => pt%first
   i = 1   ! cohort index
   do while (associated(zeig))
     
     frtrel = zeig%coh%frtrel(j) * zeig%coh%x_frt * zeig%coh%ntreea * totfrt_1    
     wat_ava = frtrel * wat_ar        ! available water per tree cohort and layer 
     if (wat_ava .ge. tr_dem(i)) then
         hupt_c    = tr_dem(i)
         tr_dem(i) = 0.
     else
         hupt_c    = wat_ava
         tr_dem(i) = tr_dem(i) - wat_ava
     endif
	     
     select case (flag_dis)    !xylem clogger disturbance 
     case (1,2)
     hupt_c = hupt_c * xylem_dis
     end select
	 
     xwatupt(i,j) = hupt_c
     zeig%coh%supply = zeig%coh%supply + hupt_c
     hupt = hupt + hupt_c
     i = i + 1
     zeig => zeig%next
   enddo  ! zeig (cohorts) 

   wupt_r(j) = hupt 
enddo    ! j
END	subroutine upt_wat1

!**************************************************************

real FUNCTION fred1(j)

! Function for calculating uptake resistance
! from CHEN (1993)
! empirical relation between soil water content and resistance
! fred1=1 if (field_cap - 10%*field_cap) <= wats <= (field_cap + 10%*field_cap)

use data_par
use data_soil

implicit none
real hf, f09, f11, wc, diff
integer j

f09 = 0.9 * field_cap(j)
f11 = 1.1 * field_cap(j)
wc  = wats(j)

if (wc .lt. wilt_p(j)) then
   hf = 0.
else if (wc .lt. f09) then
   diff = f09-wilt_p(j)
   if (diff .lt. zero) diff = 0.001
   hf = 1. - (f09-wc) / diff
else if (wc .gt. f11) then
   diff = pv(j)-f11
   if (diff .lt. zero) diff = 0.001
   hf = 0.3 + 0.7 * (pv(j)-wc) / diff
   if (hf .lt. zero) hf = 0.001
else
   hf = 1.
endif
fred1 = hf
END	function fred1

!**************************************************************
 
real FUNCTION fred2(j)

! Function for calculating uptake resistance
! from Aber and Federer (f=0.04 fuer Wasser in cm)
! only 40% of total available water are plant available per day

implicit none

integer j

fred2 = 0.05 

END	function fred2

!**************************************************************
 
real FUNCTION fred3(j)

! Function for calculating uptake resistance
! from CHEN (1993)
! modified to a profile defined in fred 
! fred3 may be described by a function (old version): 
!     fred3 = 0.0004*j*j - 0.0107*j + 0.0735
! this case: set from a root profile, defined by input of root_fr

use data_par
use data_soil

implicit none
real hf, f09, f11, wc, diff
! hf is a reduction factor in dependence on water content
real fred(15)

integer j

! uptake reduction depending on water content
f09 = 0.9 * field_cap(j)
f11 = 1.1 * field_cap(j)
wc  = wats(j)

if (wc .lt. wilt_p(j)) then
   hf = 0.
else if (wc .lt. f09) then
   diff = f09-wilt_p(j)
   if (diff .lt. zero) diff = 0.001
   hf = 1. - (f09-wc) / diff
else if (wc .gt. f11) then
   diff = pv(j)-f11
   if (diff .lt. zero) diff = 0.001
   hf = 0.3 + 0.7 * (pv(j)-wc) / diff
   if (hf .lt. zero) hf = 0.001
else
   hf = 1.
endif

fred3 = root_fr(j) * hf 

END	function fred3

!**************************************************************
 
real FUNCTION fred4(j)

! Function for calculating uptake resistance
! modified to a profile defined in fred 
! profile at Beerenbusch

use data_soil

implicit none
real fred(15)

integer j

fred = (/ 0.0, 0.03, 0.03, 0.02, 0.02, 0.02, 0.02, 0.01, 0.01, 0.01,    &
              0.01, 0.01, 0.01, 0.01, 0.01 /)  ! fred fuer Beerenbusch
 
fred4 = fred(j) 

END	function fred4

!**************************************************************

real FUNCTION fred6(j)

! Function for calculating uptake resistance
! from Kloecking (2006) simular to fred1
! empirical relation between soil water content and resistance
! fred6=1 if field_cap  <= wats <= (field_cap + 10%*field_cap)

use data_soil

implicit none
real hf, f09, f11, wc
integer j

f09 = field_cap(j)
f11 = 1.1 * field_cap(j)
wc  = wats(j)

if (wc .le. wilt_p(j)) then
   hf = 0.
else if (wc .lt. f09) then
   hf = 0.1 + (0.9 *(wc-wilt_p(j)) / (f09-wilt_p(j)))
else if (wc .gt. f11) then
   hf = 0.3 + 0.7 * (pv(j)-wc) / (pv(j)-f11)
   if (hf .lt. 0.) hf = 0.001
else
   hf = 1.
endif

fred6 = hf

END	function fred6

!**************************************************************

real FUNCTION fred7(j)

! Function for calculating uptake resistance
! from CHEN (1993)
! empirical relation between soil water content and resistance
! fred1=1 if (field_cap - 10%*field_cap) <= wats <= (field_cap + 10%*field_cap)

use data_par
use data_soil

implicit none
real hf, f09, f11, wc, diff
integer j

f09 = 0.9 * field_cap(j)
f11 = 1.1 * field_cap(j)
wc  = wats(j)

if (wc .lt. wilt_p(j)) then
   hf = 0.
else if (wc .lt. f09) then
   diff = f09-wilt_p(j)
   if (diff .lt. zero) diff = 0.001
   hf = exp(-5.*(f09-wc) / diff)
else if (wc .gt. f11) then
   diff = pv(j)-f11
   if (diff .lt. zero) diff = 0.001
   hf = 0.3 + 0.7 * (pv(j)-wc) / diff
   if (hf .lt. zero) hf = 0.001
else
   hf = 1.
endif

fred7 = hf

END	function fred7

!**************************************************************

real FUNCTION fred11(j)

! Function for calculating uptake resistance, especially adapted for Mistletoe disturbance
! function after van Wijk, 2000

use data_par
use data_soil
implicit none
real hf, S, f11, wc, diff
integer j
f11 = 1.1 * field_cap(j)
wc  = wats(j)
if (wc .lt. wilt_p(j)) then
   hf = 0.
else if (wc .lt. field_cap(j)) then
   S=(field_cap(j)-wc)/(field_cap(j)-wilt_p(j))
   hf = exp(-30*S)                                !30 = strong reduction in water avail.
else if (wc .gt. f11) then
   diff = pv(j)-f11
   if (diff .lt. zero) diff = 0.001
   hf = 0.3 + 0.7 * (pv(j)-wc) / diff
   if (hf .lt. zero) hf = 0.001
else
   hf = 1.
endif
fred11 = hf
END function fred11

!**************************************************************

 





 
 SUBROUTINE take_wat(eva_dem, psi)   
      
! Estimation of water taking out for uncovered soil 
use data_soil
use data_simul

implicit none

!input:
real :: eva_dem	   ! evaporation demand
real :: psi        ! covering

integer i, ii, j, ntag  ! max. layer of taking out
real, allocatable, dimension(:) :: gj
real, external :: b_r, funcov
real diff, gj_j, depth_j, depth_n, rij, rmax, rr, rs, sr

allocate (gj(nlay))

do i=1,nlay
   wupt_ev(i)=0.0
   gj(i)=0.0 
enddo  
      
ntag    = 0
rmax    = 0.0
depth_n = depth(n_ev_d)
	
do i=1,n_ev_d
  rij = 0.0
  rr  = depth_n/depth(i)
  sr  = 0.0
  rs  = 0.0
	    
  do j=1,i
!         depth for uncovered take out
     depth_j = depth(j)
     gj(j)   = FUNCOV(w_ev_d, rs*rr, rr*depth_j)
     rs      = depth_j
     sr      = sr + gj(j)
  enddo		! i