Changeset 1791 for LMDZ5/trunk/libf/phylmd/iophy.F90
- Timestamp:
- Jul 17, 2013, 12:20:19 PM (11 years ago)
- File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
LMDZ5/trunk/libf/phylmd/iophy.F90
r1775 r1791 6 6 ! abd REAL,private,allocatable,dimension(:),save :: io_lat 7 7 ! abd REAL,private,allocatable,dimension(:),save :: io_lon 8 REAL,allocatable,dimension(:),save :: io_lat 9 REAL,allocatable,dimension(:),save :: io_lon 10 INTEGER, save :: phys_domain_id 11 INTEGER, save :: npstn 12 INTEGER, allocatable, dimension(:), save :: nptabij 8 REAL,ALLOCATABLE,DIMENSION(:),SAVE :: io_lat 9 REAL,ALLOCATABLE,DIMENSION(:),SAVE :: io_lon 10 INTEGER, SAVE :: phys_domain_id 11 INTEGER, SAVE :: npstn 12 INTEGER, ALLOCATABLE, DIMENSION(:), SAVE :: nptabij 13 INTEGER, SAVE :: itau_iophy 14 15 !$OMP THREADPRIVATE(itau_iophy) 13 16 14 17 INTERFACE histwrite_phy 15 MODULE PROCEDURE histwrite2d_phy,histwrite3d_phy 18 MODULE PROCEDURE histwrite2d_phy,histwrite3d_phy,histwrite2d_phy_old,histwrite3d_phy_old 16 19 END INTERFACE 17 20 … … 322 325 end subroutine histbeg_phy_points 323 326 324 subroutine histwrite2d_phy(nid,lpoint,name,itau,field) 325 USE dimphy 326 USE mod_phys_lmdz_para 327 SUBROUTINE histwrite2d_phy_old(nid,lpoint,name,itau,field) 328 USE dimphy 329 USE mod_phys_lmdz_para 330 USE phys_output_var_mod 327 331 USE ioipsl 328 implicit none332 IMPLICIT NONE 329 333 include 'dimensions.h' 330 334 include 'iniprint.h' … … 341 345 integer :: ip 342 346 real,allocatable,dimension(:) :: fieldok 347 343 348 344 349 IF (size(field)/=klon) CALL abort_gcm('iophy::histwrite2d','Field first dimension not equal to klon',1) … … 379 384 deallocate(index2d) 380 385 deallocate(fieldok) 381 !$OMP END MASTER 382 end subroutine histwrite2d_phy 383 384 subroutine histwrite3d_phy(nid,lpoint,name,itau,field) 385 USE dimphy 386 USE mod_phys_lmdz_para 386 !$OMP END MASTER 387 388 389 end subroutine histwrite2d_phy_old 390 391 subroutine histwrite3d_phy_old(nid,lpoint,name,itau,field) 392 USE dimphy 393 USE mod_phys_lmdz_para 394 USE phys_output_var_mod 387 395 388 396 use ioipsl … … 401 409 INTEGER, ALLOCATABLE, dimension(:) :: index3d 402 410 real,allocatable, dimension(:,:) :: fieldok 411 403 412 404 413 IF (size(field,1)/=klon) CALL abort_gcm('iophy::histwrite3d','Field first dimension not equal to klon',1) … … 450 459 deallocate(fieldok) 451 460 !$OMP END MASTER 452 end subroutine histwrite3d_phy 461 462 end subroutine histwrite3d_phy_old 463 464 465 ! ug NOUVELLE VERSION DES WRITE AVEC LA BOUCLE DO RENTREE 466 SUBROUTINE histwrite2d_phy(var,field, STD_iff) 467 USE dimphy 468 USE mod_phys_lmdz_para 469 USE ioipsl 470 !Pour avoir nfiles, nidfiles tout ça tout ça... 471 USE phys_output_var_mod 472 473 474 475 #ifdef CPP_XIOS 476 ! USE WXIOS 477 #endif 478 479 IMPLICIT NONE 480 include 'dimensions.h' 481 482 ! integer,intent(in) :: nid 483 ! logical,intent(in) :: lpoint 484 ! character*(*), intent(IN) :: name 485 ! integer, intent(in) :: itau 486 ! real,dimension(:),intent(in) :: field 487 488 TYPE(ctrl_out), INTENT(IN) :: var 489 REAL, DIMENSION(:), INTENT(IN) :: field 490 INTEGER, INTENT(IN), OPTIONAL :: STD_iff ! ug RUSTINE POUR LES STD LEVS..... 491 492 INTEGER :: iff, iff_beg, iff_end 493 494 REAL,dimension(klon_mpi) :: buffer_omp 495 INTEGER, allocatable, dimension(:) :: index2d 496 REAL :: Field2d(iim,jj_nb) 497 498 INTEGER :: ip 499 REAL, ALLOCATABLE, DIMENSION(:) :: fieldok 500 501 ! ug RUSTINE POUR LES STD LEVS..... 502 IF (PRESENT(STD_iff)) THEN 503 iff_beg = STD_iff 504 iff_end = STD_iff 505 ELSE 506 iff_beg = 1 507 iff_end = nfiles 508 END IF 509 510 IF (size(field)/=klon) CALL abort_gcm('iophy::histwrite2d','Field first dimension not equal to klon',1) 511 512 CALL Gather_omp(field,buffer_omp) 513 !$OMP MASTER 514 CALL grid1Dto2D_mpi(buffer_omp,Field2d) 515 516 ! La boucle sur les fichiers: 517 DO iff=iff_beg, iff_end 518 IF (var%flag(iff) <= lev_files(iff) .AND. clef_files(iff)) THEN 519 520 IF(.NOT.clef_stations(iff)) THEN 521 ALLOCATE(index2d(iim*jj_nb)) 522 ALLOCATE(fieldok(iim*jj_nb)) 523 524 CALL histwrite(nid_files(iff),var%name,itau_iophy,Field2d,iim*jj_nb,index2d) 525 #ifdef CPP_XIOS 526 ! IF (iff .EQ. 1) THEN 527 ! CALL wxios_write_2D(var%name, Field2d) 528 ! ENDIF 529 #endif 530 ELSE 531 ALLOCATE(fieldok(npstn)) 532 ALLOCATE(index2d(npstn)) 533 534 IF (is_sequential) THEN 535 ! klon_mpi_begin=1 536 ! klon_mpi_end=klon 537 DO ip=1, npstn 538 fieldok(ip)=buffer_omp(nptabij(ip)) 539 ENDDO 540 ELSE 541 DO ip=1, npstn 542 ! print*,'histwrite2d is_sequential npstn ip name nptabij',npstn,ip,name,nptabij(ip) 543 IF(nptabij(ip).GE.klon_mpi_begin.AND. & 544 nptabij(ip).LE.klon_mpi_end) THEN 545 fieldok(ip)=buffer_omp(nptabij(ip)-klon_mpi_begin+1) 546 ENDIF 547 ENDDO 548 ENDIF 549 550 CALL histwrite(nid_files(iff),var%name,itau_iophy,fieldok,npstn,index2d) 551 ENDIF 552 553 deallocate(index2d) 554 deallocate(fieldok) 555 ENDIF !levfiles 556 ENDDO 557 !$OMP END MASTER 558 559 END SUBROUTINE histwrite2d_phy 560 561 562 ! ug NOUVELLE VERSION DES WRITE AVEC LA BOUCLE DO RENTREE 563 SUBROUTINE histwrite3d_phy(var, field) 564 USE dimphy 565 USE mod_phys_lmdz_para 566 567 use ioipsl 568 !Pour avoir nfiles, nidfiles tout ça tout ça... 569 USE phys_output_var_mod 570 571 572 #ifdef CPP_XIOS 573 ! USE WXIOS 574 #endif 575 576 577 IMPLICIT NONE 578 include 'dimensions.h' 579 580 ! integer,intent(in) :: nid 581 ! logical,intent(in) :: lpoint 582 ! character*(*), intent(IN) :: name 583 ! integer, intent(in) :: itau 584 ! real,dimension(:,:),intent(in) :: field ! --> field(klon,:) 585 586 TYPE(ctrl_out), INTENT(IN) :: var 587 REAL, DIMENSION(:,:), INTENT(IN) :: field ! --> field(klon,:) 588 589 590 REAL,DIMENSION(klon_mpi,SIZE(field,2)) :: buffer_omp 591 REAL :: Field3d(iim,jj_nb,SIZE(field,2)) 592 INTEGER :: ip, n, nlev, iff 593 INTEGER, ALLOCATABLE, DIMENSION(:) :: index3d 594 REAL,ALLOCATABLE, DIMENSION(:,:) :: fieldok 595 596 IF (size(field,1)/=klon) CALL abort_gcm('iophy::histwrite3d','Field first dimension not equal to klon',1) 597 nlev=size(field,2) 598 599 ! print*,'hist3d_phy mpi_rank npstn=',mpi_rank,npstn 600 601 ! DO ip=1, npstn 602 ! print*,'hist3d_phy mpi_rank nptabij',mpi_rank,nptabij(ip) 603 ! ENDDO 604 605 CALL Gather_omp(field,buffer_omp) 606 !$OMP MASTER 607 CALL grid1Dto2D_mpi(buffer_omp,field3d) 608 609 610 ! BOUCLE SUR LES FICHIERS 611 DO iff=1, nfiles 612 IF (var%flag(iff) <= lev_files(iff) .AND. clef_files(iff)) THEN 613 IF (.NOT.clef_stations(iff)) THEN 614 ALLOCATE(index3d(iim*jj_nb*nlev)) 615 ALLOCATE(fieldok(iim*jj_nb,nlev)) 616 CALL histwrite(nid_files(iff),var%name,itau_iophy,Field3d,iim*jj_nb*nlev,index3d) 617 #ifdef CPP_XIOS 618 ! IF (iff .EQ. 1) THEN 619 ! CALL wxios_write_3D(var%name, Field3d(:,:,1:klev)) 620 ! ENDIF 621 #endif 622 623 ELSE 624 nlev=size(field,2) 625 ALLOCATE(index3d(npstn*nlev)) 626 ALLOCATE(fieldok(npstn,nlev)) 627 628 IF (is_sequential) THEN 629 ! klon_mpi_begin=1 630 ! klon_mpi_end=klon 631 DO n=1, nlev 632 DO ip=1, npstn 633 fieldok(ip,n)=buffer_omp(nptabij(ip),n) 634 ENDDO 635 ENDDO 636 ELSE 637 DO n=1, nlev 638 DO ip=1, npstn 639 IF(nptabij(ip).GE.klon_mpi_begin.AND. & 640 nptabij(ip).LE.klon_mpi_end) THEN 641 fieldok(ip,n)=buffer_omp(nptabij(ip)-klon_mpi_begin+1,n) 642 ENDIF 643 ENDDO 644 ENDDO 645 ENDIF 646 CALL histwrite(nid_files(iff),var%name,itau_iophy,fieldok,npstn*nlev,index3d) 647 ENDIF 648 deallocate(index3d) 649 deallocate(fieldok) 650 ENDIF 651 ENDDO 652 !$OMP END MASTER 653 END SUBROUTINE histwrite3d_phy 453 654 454 655 end module iophy
Note: See TracChangeset
for help on using the changeset viewer.