!*****************************************************************! !* *! !* 4C (FORESEE) *! !* *! !* *! !* Subroutines for: *! !* Aspen management *! !* contains: *! !* SR aspman_ini *! !* SR asp_manag *! !* SR asp_sprout *! !* SR asp_pruning *! !* *! !* 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 aspman_ini use data_manag use data_species use data_simul use data_stand use data_par implicit none integer :: manag_unit,i, ios character(len=150) :: filename logical :: ex character ::text manag_unit=getunit() filename = manfile(ip) allocate(thin_flag1(nspec_tree)) thin_flag1 = -1 allocate(yman(100)) allocate(rel_part(100)) yman = 0 rel_part = 0 call testfile(filename,ex) open(manag_unit,file=trim(filename)) ! read head of data-file do read(manag_unit,*) text if(text .ne. '#')then backspace(manag_unit);exit endif enddo i = 1 do read(manag_unit,*,iostat=ios) yman(i), rel_part(i) if(ios < 0) exit i = i+1 end do num_man = i-1 close(manag_unit) end subroutine aspman_ini !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! subroutine asp_manag use data_manag use data_simul implicit none integer :: i do i=1,num_man if(yman(i).eq.time) then call asp_pruning if(i.ne.num_man) then call asp_sprout flag_sprout = 1 end if end if end do end subroutine asp_manag !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! subroutine asp_sprout use data_manag use data_species use data_simul use data_stand use data_par use data_help use data_soil use data_tsort implicit none integer :: taxnr, i, j, nsp, acoh REAL :: shoot real :: faktor REAL :: x1,x2,xacc,h_root, root REAL :: rtflsp, stump_dw, stump_v, rtbis TYPE(cohort) ::tree_ini real, dimension(:), save, allocatable :: treea, crt, frt, stumpw integer, dimension(:), save, allocatable :: spectyp, cohid ! distribution of coarse root matter of coppice shoots real, dimension(6) :: fac_rob=(/0.0666, 0.1332, 0.1998, 0.2664,0.334, 0./) external weight1 external rtflsp external rtbis allocate ( treea(anz_coh), crt(anz_coh), frt(anz_coh), spectyp(anz_coh), cohid(anz_coh), stumpw(anz_coh)) if(flag_reg.eq.18) then nsprout = 5 end if i = 1 zeig=>pt%first do if(.not.associated(zeig)) exit if(zeig%coh%ntreem.ne.0.and. zeig%coh%ntreea.eq.0.and. zeig%coh%x_crt.ne.0) then treea(i) = zeig%coh%ntreem taxnr = zeig%coh%species crt(i) = zeig%coh%x_crt frt(i) = zeig%coh%x_frt spectyp(i) = zeig%coh%species cohid(i) = zeig%coh%ident call stump( zeig%coh%x_ahb, zeig%coh%asapw,zeig%coh%dcrb,zeig%coh%x_hbole, & zeig%coh%height, taxnr,stump_v, stump_dw) stumpw(i) = stump_dw i = i+1 end if zeig=>zeig%next end do acoh = i-1 do i =1, acoh if(flag_reg.eq.15) then faktor = 0.25 else faktor = fac_rob(1) end if do j = 1, nsprout tree_ini%species = spectyp(i) nsp = spectyp(i) hnspec = nsp h_root = faktor * (crt(i)*0.3 + stumpw(i)* 0.5) max_coh= max_coh +1 call coh_initial (tree_ini) tree_ini%ident = max_coh tree_ini%x_age = 1 tree_ini%ntreea = treea(i) tree_ini%nta = treea(i) mschelp = h_root x1 = 0. x2 = 0.1 xacc = (1.0e-10) * (x1+x2)/2 root = rtbis(weight1,x1,x2,xacc) tree_ini%x_sap = root shoot = root*1000. ! [g] tree_ini%x_fol= (spar(nsp)%seeda*(tree_ini%x_sap** spar(nsp)%seedb)) ![kg] ! [kg] tree_ini%x_frt = faktor * frt(i) ! [kg] tree_ini%med_sla = spar(nsp)%psla_min + spar(nsp)%psla_a*0.5 tree_ini%t_leaf = tree_ini%med_sla* tree_ini%x_fol ! [m-2] tree_ini%ca_ini = tree_ini%t_leaf IF(spar(tree_ini%species)%Phmodel==1) THEN tree_ini%P=0 tree_ini%I=1 ELSE tree_ini%P=0 tree_ini%I=0 tree_ini%Tcrit=0 END IF zeig=>pt%first do if(.not.associated(zeig)) exit if(zeig%coh%ident.eq. cohid(i)) then tree_ini%rooteff = zeig%coh%rooteff exit end if zeig=>zeig%next end do ! tranformation of shoot biomass kg --> mg if(nsp.ne.2)tree_ini%height = spar(nsp)%pheight1*(shoot*1000.)**spar(nsp)%pheight2 ! [cm] calculated from shoot biomass (mg) if(tree_ini%height.eq.0.) then nsp = nsp end if ! bole height from stump tree_ini%x_hbole = stoh(nsp) IF(tree_ini%ntreea.ne.0.) then IF (.not. associated(pt%first)) THEN ALLOCATE (pt%first) pt%first%coh = tree_ini NULLIFY(pt%first%next) ELSE ALLOCATE(zeig) zeig%coh = tree_ini zeig%next => pt%first pt%first => zeig END IF anz_coh=anz_coh+1 END IF if(flag_reg.eq.15) then faktor = faktor + 0.0833333 else faktor = fac_rob(j+1) end if end do ! j, nsprouts end do ! i deallocate ( treea, crt, frt, spectyp,cohid, stumpw) end subroutine asp_sprout subroutine asp_pruning use data_manag use data_species use data_simul use data_stand use data_par implicit none integer :: taxnr, j zeig=>pt%first do if(.not.associated(zeig)) exit zeig%coh%ntreem = zeig%coh%ntreea zeig%coh%ntreea = 0 zeig%coh%nta = 0. zeig=>zeig%next end do ! calculation of total dry mass of all harvested trees (stem + twigs and branches) sumNPP = 0 sumvsab = 0. sumvsab_m3 = 0. svar%sumvsab = 0. zeig=>pt%first do if(.not.associated(zeig)) exit ns = zeig%coh%species sumvsab = sumvsab + zeig%coh%ntreem*(zeig%coh%x_sap + zeig%coh%x_hrt + zeig%coh%x_tb) sumvsab_m3 = sumvsab_m3 + zeig%coh%ntreem*(zeig%coh%x_sap + zeig%coh%x_hrt+zeig%coh%x_tb)/(spar(ns)%prhos*1000000) svar(ns)%sumvsab = svar(ns)%sumvsab + zeig%coh%ntreem*(zeig%coh%x_sap + zeig%coh%x_hrt + zeig%coh%x_tb) sumnpp = sumnpp + zeig%coh%ntreem*zeig%coh%npp zeig=>zeig%next end do sumvsab_m3 = sumvsab_m3 * 10000./kpatchsize ! kg/ha sumvsab = sumvsab * 10000./kpatchsize ! kg/ha do j = 1, nspec_tree svar(j)%sumvsab = svar(j)%sumvsab * 10000./kpatchsize end do ! cumulative harvested stem mass cumsumvsab = cumsumvsab + sumvsab ! litter pools ! adding biomasses to litter pools depending on stage of stand 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 endif zeig=>zeig%next enddo end subroutine asp_pruning