Skip to content
Snippets Groups Projects
aust_manag.f 14.2 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
!*****************************************************************!
!*                                                               *!
!*                     4C (FORESEE)                              *!
!*                                                               *!
!*                                                               *!
!*                Subroutines for:                               *!
!*                Austrian  management                           *!
!*                 contains:                                     *!
!*                SR aust_ini                                    *!
!*                SR aust_manag                                  *!
!*                SR plant_aust                                  *!
!*                SR calc_rel_class                              *!
!*                                                               *!
!*                  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 aust_ini
 use data_manag
 use data_species
 use data_simul
 use data_stand
 implicit none

 integer :: manag_unit,i, ih1,ih2,ios,ih4, flp , flag_help
 character(len=150) :: filename
 logical :: ex
 character ::text
 real    :: hp, ih3

 manag_unit=getunit()
 filename = manfile(ip)

 allocate(thin_flag1(nspec_tree))
 flag_help = 0
 thin_flag1=-1
 thin_dead = 1
 allocate(yman(1000))
 allocate(dbh_clm(1000))
 allocate(rem_clm(1000))
 allocate(spec_man(1000))
 allocate(act(1000))
 allocate(rel_part(1000))
yman = 0
 dbh_clm = 0
 rem_clm = 0.
 spec_man = 0
 act = 0
 rel_part = 0
 flp = 0
  call testfile(filename,ex)
  open(manag_unit,file=trim(filename))
! read head of data-file
 do
    read(manag_unit,*) text
    if(text .ne. 's')then

       backspace(manag_unit);exit
    endif
 enddo

 i=1
 do
    read(manag_unit,*,iostat=ios) ih1,ih2, ih3, hp,ih4
    if(ios<0) exit
      yman(i) = ih1            ! year of treatment
      if(ih2.eq.1) then
! Fichte/ Spruce
          spec_man(i) = 2
      else if(ih2.eq.2) then
! Kiefer/ Pine
          spec_man(i) = 3
      else if(ih2.eq.3) then
! Eiche/ oak
          spec_man(i) = 4
      else if(ih2.eq.4) then
          spec_man(i) = 1
      end if

        ! species  number
      act(i) = ih4
      if(ih1.ne.-999 ) then
         if(flp.eq.0) then
            dbh_clm(i) = int(ih3)    ! dbh-cluss number for treatment
            rem_clm(i) = hp          ! removal of biomass
            i = i+1
         else
            act(i) = ih4
            rel_part(i) = ih3
            rem_clm(i) = 0
            i = i+1
         end if
      else
       if(i.eq.1) thin_dead = 0
       flp = 1
       backspace(manag_unit)
      end if
 end do
 num_man = i-1

 close(manag_unit)
END SUBROUTINE aust_ini
    
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!

SUBROUTINE aust_manag

 use data_manag
 use data_stand
 use data_simul
 use data_species
 use data_par

 implicit none
 integer    :: i,j,hcl, ha, taxnr,k,l, help_fl,helpz, helps
 real, dimension(5)   :: rel_biom, harv_biom
 real, dimension(5)   :: num_ccl, contr
 real,dimension(5)    :: help_rem_clm
 integer,dimension(5) :: help_rel_dbh, hrd
 real                 :: stump_dw, stump_v, hrb
 ha =0
 rel_biom = 0.
 num_ccl =0
 harv_biom = 0.
 help_rel_dbh = 0
 help_rem_clm = 0.
 hrd = 0
 helpz = 0
 helps = 0
  call calc_rel_class

! calculation of stem biomass of relative dbh-class
  zeig=>pt%first
  do
    if(.not.associated(zeig)) exit
    if(zeig%coh%diam.ne.0) then
      hcl = zeig%coh%rel_dbh_cl
      if(hcl.ne.0) then
          num_ccl(hcl)= num_ccl(hcl) +1
          rel_biom(hcl)= rel_biom(hcl) + (zeig%coh%x_sap + zeig%coh%x_hrt)*zeig%coh%ntreea
      end if
     end if
       zeig=>zeig%next
  end do

 do l=1,nspecies
  help_rel_dbh = 0
  help_rem_clm = 0.
  helpz = 0
  helps = 0
! calculation of stem biomass of relative dbh-class
  zeig=>pt%first
  do
    if(.not.associated(zeig)) exit
    if(zeig%coh%diam.ne.0.and.zeig%coh%species.eq.l) then
      hcl = zeig%coh%rel_dbh_cl
      if(hcl.ne.0) then
          num_ccl(hcl)= num_ccl(hcl) +1
          rel_biom(hcl)= rel_biom(hcl) + (zeig%coh%x_sap + zeig%coh%x_hrt)*zeig%coh%ntreea
      end if
     end if
       zeig=>zeig%next
  end do

  hrd=0
   do i=1,num_man

    if(yman(i).eq.time) then

      if(act(i) .eq.1.and.spec_man(i).eq.l) then
         zeig=>pt%first
         do
           if(.not.associated(zeig)) exit
           if(zeig%coh%diam.ne.0) then
            if(zeig%coh%species.eq.l) then
               hrd(zeig%coh%rel_dbh_cl)= 1
            end if
           end if
           zeig=>zeig%next
         end do
         help_rel_dbh(dbh_clm(i)) = 1
         help_rem_clm(dbh_clm(i)) = rem_clm(i)
       end if ! act(i)
    end if !yman(i)
   end do ! num_man
      do j=1,5
         if(help_rel_dbh(j).eq.1.and.hrd(j).eq.0) then
                if(j.eq.1.) then
                     do k=2,5
                       if(hrd(k).ne.0) then
                         help_rem_clm(k) = help_rem_clm(k) + help_rem_clm(j)
                         help_rel_dbh(k)=1
                         exit
                       end if
                     end do
                else if (j.eq.5.) then
                       do k= 4,1,-1
                           if(hrd(k).eq.1) then
                             help_rem_clm(k) = help_rem_clm(k) + help_rem_clm(j)
                             help_rel_dbh(k) = 1
                             exit
                           endif
                        end do
                 else
                    do k=j,5
                        if(hrd(k).eq.1) then
                             help_rem_clm(k) = help_rem_clm(k) + help_rem_clm(j)*0.5
                              help_rel_dbh(k) = 1
                              exit
                        end if
                    end do
                    do k=j,1,-1
                        if(hrd(k).eq.1) then
                             help_rem_clm(k) = help_rem_clm(k) + help_rem_clm(j)*0.5
                              help_rel_dbh(k) = 1
                              exit
                        end if
                    end do
                 end if
          help_rel_dbh(j) = 0
          help_rem_clm(j) = 0.
          end if

      end do
! thinning
  help_fl = 0
   do i=1,num_man

    if(yman(i).eq.time.and.help_fl.eq.0) then
                do k=1,5
                 helps = helps +  help_rel_dbh(k)
               end do

       help_fl=1
       zeig=>pt%first
       do
           if(.not.associated(zeig)) exit
           if(zeig%coh%diam.ne.0.and.zeig%coh%species.eq.l) then


               do k=1,5
                  if(zeig%coh%rel_dbh_cl.eq.k.and.help_rel_dbh(k).eq.1) then
                     if(help_rem_clm(k).gt.1.) help_rem_clm(k) = 1.
                     if( help_rem_clm(k) .eq. 1.)then

                         if(zeig%coh%underst.eq.0.and.zeig%coh%x_age.gt. 20) ha=int(help_rem_clm(k)* zeig%coh%ntreea+0.5)
                         helpz = helpz +1
                     else
                         ha=int(help_rem_clm(k)* zeig%coh%ntreea+0.5)

                     end if
                     if(ha.lt.1) ha = 1
                     if(help_rem_clm(k) .ne.1) then
                          harv_biom(k) = harv_biom(k) + ha* (zeig%coh%x_sap + zeig%coh%x_hrt)
                          hrb =  help_rem_clm(k)* rel_biom(k)
                          if(harv_biom(k).eq.rel_biom(k)) then
                               ha = ha -1
                          end if
                      end if

                     zeig%coh%ntreea = zeig%coh%ntreea - ha
                     zeig%coh%nta = zeig%coh%ntreea
                     zeig%coh%ntreem = zeig%coh%ntreem + ha

                  end if
               end do  ! k loop
            end if
            zeig=>zeig%next
      end do   ! zeig loop

    end if
  end do ! num_man

 if(helps.gt.0.and.helpz.ge.helps) then
    zeig=>pt%first
    do
           if(.not.associated(zeig)) exit
           if(zeig%coh%species.eq.l.and.zeig%coh%underst.eq.1) then
                   zeig%coh%underst = 0
           end if
           zeig => zeig%next
    end do
 end if

 write(9898,*) time, 'totbio', rel_biom
 write(9898,*) time, 'harvbio', harv_biom
  do i=1,5
  if(rel_biom(i).ne.0.) then
      contr(i) = harv_biom(i)/rel_biom(i)
  else
      contr(i) = 0.
  end if
  end do
 write(9898,*) time,l, contr

  rel_biom = 0.
  harv_biom = 0.
 end do ! nspecies

 ! planting

 do i=1,num_man
   if(yman(i).eq.time.and.act(i).ne.1) then
      call plant_aust(i)

    end if  ! act
 end do

stump_sum = 0
 zeig=>pt%first

 do
   if(.not.associated(zeig)) exit
   taxnr=zeig%coh%species

   if(zeig%coh%ntreem>0)then
! all parts without stems of trees are input for litter
         zeig%coh%litC_fol = zeig%coh%litC_fol + zeig%coh%ntreem*(1.-spar(taxnr)%psf)*zeig%coh%x_fol*cpart
         zeig%coh%litN_fol = zeig%coh%litN_fol + zeig%coh%ntreem*((1.-spar(taxnr)%psf)*zeig%coh%x_fol*cpart)/spar(taxnr)%cnr_fol
         zeig%coh%litC_frt = zeig%coh%litC_frt + zeig%coh%ntreem*zeig%coh%x_frt*cpart
         zeig%coh%litN_frt = zeig%coh%litN_frt + zeig%coh%ntreem*zeig%coh%x_frt*cpart/spar(taxnr)%cnr_frt
         zeig%coh%litC_tb = zeig%coh%litC_tb + zeig%coh%ntreem*zeig%coh%x_tb*cpart
         zeig%coh%litN_tb = zeig%coh%litN_tb + zeig%coh%ntreem*zeig%coh%x_tb*cpart/spar(taxnr)%cnr_tbc
         zeig%coh%litC_crt = zeig%coh%litC_crt + zeig%coh%ntreem*zeig%coh%x_crt*cpart
         zeig%coh%litN_crt = zeig%coh%litN_crt + zeig%coh%ntreem*zeig%coh%x_crt*cpart/spar(taxnr)%cnr_crt

! stumps into stem litter
          call stump( zeig%coh%x_ahb, zeig%coh%asapw,zeig%coh%dcrb,zeig%coh%x_hbole,  &
                      zeig%coh%height, taxnr,stump_v, stump_dw)
          zeig%coh%litC_stem = zeig%coh%litC_stem +  zeig%coh%ntreem*stump_dw*cpart
          zeig%coh%litN_stem = zeig%coh%litC_stem/spar(taxnr)%cnr_stem
          stump_sum = stump_sum + zeig%coh%ntreem*stump_dw
   endif
 zeig=>zeig%next
 enddo

 sumvsab = 0.
 sumvsab_m3 = 0.
 svar%sumvsab = 0.

 zeig=>pt%first
   do while (associated(zeig)) 

     ns = zeig%coh%species
     sumvsab          = sumvsab + zeig%coh%ntreem*(zeig%coh%x_sap + zeig%coh%x_hrt)
     sumvsab_m3       = sumvsab_m3 +  zeig%coh%ntreem*(zeig%coh%x_sap + zeig%coh%x_hrt)/(spar(ns)%prhos*1000000)
     svar(ns)%sumvsab = svar(ns)%sumvsab +  zeig%coh%ntreem*(zeig%coh%x_sap + zeig%coh%x_hrt)

     zeig=>zeig%next

   end do
  sumvsab = sumvsab *  10000./kpatchsize                 ! kg/ha
  sumvsab_m3 = sumvsab_m3 *  10000./kpatchsize           ! kg/ha


  do k = 1, nspec_tree
    svar(k)%sumvsab = svar(k)%sumvsab * 10000./kpatchsize           ! kg/ha
  end do
! cumulative harvested stem mass
  cumsumvsab = cumsumvsab + sumvsab

 if(thin_dead.ne.0) then
     call class_man
 end if

END SUBROUTINE aust_manag

!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
    
SUBROUTINE plant_aust(mp)
use data_manag
use data_plant
use data_species
use data_stand

implicit none
integer     :: fl_plant, i, nplant,taxid,mp
 real       :: age,          &
               pl_height,    &
               sdev,         &
               plhmin

  infspec = 0
  npl_mix = 0
 fl_plant = act(mp)
 select case(fl_plant)
    case(2)
      infspec(2) = 1
      npl_mix(2) = 2500
      
    case(3)
      infspec(2) = 1
      npl_mix(2) = 10000
      
    case(4)
      infspec(3) = 1
      npl_mix(3) = 5000

    case(5)
      infspec(3) = 1
      npl_mix(3) = 2000

    case(6)
      infspec(1) = 1
      npl_mix(1) = 500

    case(7)
      infspec(1) = 1
      npl_mix(1) = 5000

    case(8)
      infspec(4) = 1
      npl_mix(4) = 5000

    case(9)
      infspec(1) = 1
      npl_mix(1) = 1000

      infspec(4) = 1
      npl_mix(4) = 3500

    case(10)
      infspec(3) = 1
      npl_mix(3) = 2500

      infspec(4) = 1
      npl_mix(4) = 2500

    case(11)
      infspec(3) = 1
      npl_mix(3) = 2500

      infspec(1) = 1
      npl_mix(1) = 2500

    case(12)
      infspec(3) = 1
      npl_mix(3) = 7000

    case(13)
      infspec(4) = 1
      npl_mix(4) = 2500

 end select
       do i = 1,nspec_tree
         if (infspec(i).eq.1) then

               taxid = i
! data for Austria
               age = pl_age(taxid)
               pl_height = plant_height(taxid)
               plhmin = plant_hmin(taxid)
               nplant = rel_part(mp)*nint(npl_mix(taxid)*kpatchsize/10000)
               sdev = hsdev(taxid)
               call gener_coh(taxid, age, pl_height, plhmin, nplant,sdev)

          end if
      end do  ! i

END SUBROUTINE plant_aust

!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!    
    
SUBROUTINE calc_rel_class
 use data_manag
 use data_stand
 use data_species

 implicit none
 integer    :: nrmax, i, j, k, adm
 real,dimension(10)       :: maxdbh, mindbh,class_wd

 integer  :: nrmin
 real     :: help, help_h1
 class_wd =0.
 maxdbh = 0.
 mindbh = 0.

do  j= 1,nspecies
  call max_dbh(nrmax,help,adm, j)
  call min_dbh(nrmin,help_h1,adm, j)
  zeig=>pt%first

  do
    if(.not.associated(zeig)) exit

       if(zeig%coh%ident.eq.nrmax.and.zeig%coh%species.eq.j) then
          maxdbh(j) = help

       else if(zeig%coh%ident.eq.nrmin.and.zeig%coh%species.eq.j) then
          mindbh(j) = help_h1

       end if

     zeig=>zeig%next
  end do
end do

 do j=1,nspecies

     class_wd(j) =  (maxdbh(j)-mindbh(j))/5
     k = 5
    zeig=>pt%first

    do
       if(.not.associated(zeig)) exit
          if(zeig%coh%species.eq.j.and. zeig%coh%diam.gt.0) then

               do i=1,k
                 if(zeig%coh%diam.ge.(mindbh(j)+class_wd(j)*(i-1)).and.zeig%coh%diam.lt.(mindbh(j)+class_wd(j)*i)) then
                      zeig%coh%rel_dbh_cl = i
                      exit
                 else if (zeig%coh%diam.eq.maxdbh(j)) then
                       zeig%coh%rel_dbh_cl = 5
                 end if
               end do
           end if

         zeig=>zeig%next
    end do
 end do
END SUBROUTINE calc_rel_class