Changeset f3054ea in flexpart.git for src/concoutput_nest_mpi.f90
- Timestamp:
- Aug 27, 2020, 8:00:15 PM (4 years ago)
- Branches:
- GFS_025, dev
- Children:
- 4ab2fbf
- Parents:
- a756649
- File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
src/concoutput_nest_mpi.f90
r16b61a5 rf3054ea 58 58 !***************************************************************************** 59 59 60 61 60 use unc_mod 62 61 use point_mod … … 74 73 real :: sp_fact 75 74 real :: outnum,densityoutrecept(maxreceptor),xl,yl 75 ! RLT 76 real :: densitydryrecept(maxreceptor) 77 real :: factor_dryrecept(maxreceptor) 76 78 77 79 !real densityoutgrid(0:numxgrid-1,0:numygrid-1,numzgrid), … … 100 102 character :: adate*8,atime*6 101 103 character(len=3) :: anspec 104 logical :: lexist 102 105 integer :: mind 103 106 ! mind eso:added to ensure identical results between 2&3-fields versions … … 176 179 densityoutgrid(ix,jy,kz)=(rho(iix,jjy,kzz,mind)*dz1+ & 177 180 rho(iix,jjy,kzz-1,mind)*dz2)/dz 181 ! RLT 182 densitydrygrid(ix,jy,kz)=(rho_dry(iix,jjy,kzz,mind)*dz1+ & 183 rho_dry(iix,jjy,kzz-1,mind)*dz2)/dz 178 184 end do 179 185 end do … … 187 193 !densityoutrecept(i)=rho(iix,jjy,1,2) 188 194 densityoutrecept(i)=rho(iix,jjy,1,mind) 195 ! RLT 196 densitydryrecept(i)=rho_dry(iix,jjy,1,mind) 189 197 end do 190 198 199 ! RLT 200 ! conversion factor for output relative to dry air 201 factor_drygrid=densityoutgrid/densitydrygrid 202 factor_dryrecept=densityoutrecept/densitydryrecept 191 203 192 204 ! Output is different for forward and backward simulations … … 211 223 212 224 write(anspec,'(i3.3)') ks 225 226 if (DRYBKDEP.or.WETBKDEP) then !scavdep output 227 if (DRYBKDEP) & 228 open(unitoutgrid,file=path(2)(1:length(2))//'grid_drydep_nest_'//adate// & 229 atime//'_'//anspec,form='unformatted') 230 if (WETBKDEP) & 231 open(unitoutgrid,file=path(2)(1:length(2))//'grid_wetdep_nest_'//adate// & 232 atime//'_'//anspec,form='unformatted') 233 write(unitoutgrid) itime 234 else 213 235 if ((iout.eq.1).or.(iout.eq.3).or.(iout.eq.5)) then 214 236 if (ldirect.eq.1) then … … 223 245 write(unitoutgrid) itime 224 246 endif 247 endif 225 248 226 249 if ((iout.eq.2).or.(iout.eq.3)) then ! mixing ratio … … 553 576 end do 554 577 555 556 557 ! Reinitialization of grid 558 !************************* 559 560 do ks=1,nspec 561 do kp=1,maxpointspec_act 562 do i=1,numreceptor 563 creceptor(i,ks)=0. 564 end do 565 do jy=0,numygridn-1 566 do ix=0,numxgridn-1 567 do l=1,nclassunc 568 do nage=1,nageclass 569 do kz=1,numzgrid 570 griduncn(ix,jy,kz,ks,kp,l,nage)=0. 571 end do 572 end do 573 end do 574 end do 578 ! RLT Aug 2017 579 ! Write out conversion factor for dry air 580 inquire(file=path(2)(1:length(2))//'factor_drygrid_nest',exist=lexist) 581 if (lexist) then 582 ! open and append 583 open(unitoutfactor,file=path(2)(1:length(2))//'factor_drygrid_nest',form='unformatted',& 584 status='old',action='write',access='append') 585 else 586 ! create new 587 open(unitoutfactor,file=path(2)(1:length(2))//'factor_drygrid_nest',form='unformatted',& 588 status='new',action='write') 589 endif 590 sp_count_i=0 591 sp_count_r=0 592 sp_fact=-1. 593 sp_zer=.true. 594 do kz=1,numzgrid 595 do jy=0,numygridn-1 596 do ix=0,numxgridn-1 597 if (factor_drygrid(ix,jy,kz).gt.(1.+smallnum).or.factor_drygrid(ix,jy,kz).lt.(1.-smallnum)) then 598 if (sp_zer.eqv..true.) then ! first value not equal to one 599 sp_count_i=sp_count_i+1 600 sparse_dump_i(sp_count_i)= & 601 ix+jy*numxgridn+kz*numxgridn*numygridn 602 sp_zer=.false. 603 sp_fact=sp_fact*(-1.) 604 endif 605 sp_count_r=sp_count_r+1 606 sparse_dump_r(sp_count_r)= & 607 sp_fact*factor_drygrid(ix,jy,kz) 608 else ! factor is one 609 sp_zer=.true. 610 endif 575 611 end do 576 612 end do 577 613 end do 578 614 write(unitoutfactor) sp_count_i 615 write(unitoutfactor) (sparse_dump_i(i),i=1,sp_count_i) 616 write(unitoutfactor) sp_count_r 617 write(unitoutfactor) (sparse_dump_r(i),i=1,sp_count_r) 618 close(unitoutfactor) 619 620 ! Reinitialization of grid 621 !************************* 622 623 ! do ks=1,nspec 624 ! do kp=1,maxpointspec_act 625 ! do i=1,numreceptor 626 ! creceptor(i,ks)=0. 627 ! end do 628 ! do jy=0,numygridn-1 629 ! do ix=0,numxgridn-1 630 ! do l=1,nclassunc 631 ! do nage=1,nageclass 632 ! do kz=1,numzgrid 633 ! griduncn(ix,jy,kz,ks,kp,l,nage)=0. 634 ! end do 635 ! end do 636 ! end do 637 ! end do 638 ! end do 639 ! end do 640 ! end do 641 creceptor(:,:)=0. 642 griduncn(:,:,:,:,:,:,:)=0. 643 579 644 if (mp_measure_time) call mpif_mtime('iotime',1) 580 645 ! if (mp_measure_time) then
Note: See TracChangeset
for help on using the changeset viewer.