!> \file mix-SIA-L1_mod.f90 !! Module qui defini le type d'association entre L1 et SIA !< !> \namespace resolmeca_SIA_L1 !! Module qui defini le type d'association entre L1 et SIA !! @note Test sur la variable globale i_resolmeca !! \author ... !! \date ... !! @note Used module !! @note - use module3D_phy !! !< module resolmeca_SIA_L1 ! module qui défini le type d'association entre L1 et SIA ! test sur la variable globale i_resolmeca use module3D_phy contains !> SUBROUTINE: init_resol_meca !! Routine pour l'initialisation de la resolution mecanique !> subroutine init_resol_meca namelist/meca_SIA_L1/i_resolmeca ! formats pour les ecritures dans 42 428 format(A) ! lecture des parametres du run block eaubasale1 !-------------------------------------------------------------------- rewind(num_param) ! pour revenir au debut du fichier param_list.dat read(num_param,meca_SIA_L1) write(num_rep_42,428)'!___________________________________________________________' write(num_rep_42,428) '&meca_SIA_L1 ! bloc resol_meca ' write(num_rep_42,*) write(num_rep_42,*) 'i_resolmeca = ', i_resolmeca write(num_rep_42,*)'/' write(num_rep_42,428) '! i_resolmeca type d association entre SIA et L1' write(num_rep_42,428) '! i_resolmeca=0 chacun dans sa zone' write(num_rep_42,428) '! i_resolmeca=1 dans les zones stream, addition si uxdef > uxL1 (MIS11 Cairns)' write(num_rep_42,428) '! i_resolmeca=2 addition systematique dans les zones stream ' write(num_rep_42,*) end subroutine init_resol_meca !> SUBROUTINE: mix_SIA_L1 !! Subroutine qui associe SIA et L1: la methode depend de i_resolmeca !> subroutine mix_SIA_L1 !$ USE OMP_LIB ! variables parallelisation openMP !!$ integer :: rang ,nb_taches ! subroutine qui associe SIA et L1 ! la methode depend de i_resolmeca ! use module_choix if (itracebug.eq.1) call tracebug(' routine mix_SIA_L1') debug_3D(:,:,9)=0. debug_3D(:,:,10)=0. debug_3D(:,:,11)=0. debug_3D(:,:,12)=0. debug_3D(:,:,13)=ubx(:,:) debug_3D(:,:,14)=uby(:,:) if (i_resolmeca.eq.1) then ! prend les vitesses de la SIA si elles s'averent plus fortes ! que celles du L1. On garde le L1 comme vitesse basale. ! version utilisee pour les runs MIS11 de Cairns !$OMP PARALLEL !$OMP DO do j=2,ny do i=2,nx if (flgzmx(i,j)) then debug_3D(i,j,11)=uxflgz(i,j) if (abs((uxdef(i,j)+ubx(i,j))).gt.abs(uxflgz(i,j))) then ! SIA plus grande que L1 uxbar(i,j)=uxflgz(i,j)+uxdef(i,j) ! on ajoute la deform. SIA ubx(i,j) = uxflgz(i,j) debug_3D(i,j,9)=uxdef(i,j) end if end if if (flgzmy(i,j)) then debug_3D(i,j,12)=uyflgz(i,j) if (abs((uydef(i,j)+uby(i,j))).gt.abs(uyflgz(i,j))) then ! SIA plus grande que L1 uybar(i,j)=uyflgz(i,j)+uydef(i,j) ! on ajoute la deform. SIA uby(i,j) = uyflgz(i,j) debug_3D(i,j,10)=uydef(i,j) end if end if end do end do !$OMP END DO !$OMP END PARALLEL else if (i_resolmeca.eq.2) then ! addition systematique ! !$OMP PARALLEL PRIVATE(rang,nb_taches) ! !$ rang=OMP_GET_THREAD_NUM() ! !$ nb_taches=OMP_GET_NUM_THREADS() ! !$OMP DO SCHEDULE(STATIC,NY/nb_taches) !$OMP PARALLEL !$OMP DO do j=2,ny do i=2,nx ! test sur tout le tableau : if (flgzmx(i,j)) then ! debug_3D(i,j,11)=uxflgz(i,j) uxbar(i,j)=uxflgz(i,j)+uxdef(i,j) ! on ajoute la deform. SIA ubx(i,j) = uxflgz(i,j) ! do k=1,nz ! ux(i,j,k)=ux(i,j,k)-ux(i,j,nz)+uxflgz(i,j) ! remplace le glissement par uxflgz ! end do ! debug_3D(i,j,9)=uxdef(i,j) endif if (flgzmy(i,j)) then ! debug_3D(i,j,12)=uyflgz(i,j) uybar(i,j)=uyflgz(i,j)+uydef(i,j) ! on ajoute la deform. SIA uby(i,j) = uyflgz(i,j) ! do k=1,nz ! uy(i,j,k)=uy(i,j,k)-uy(i,j,nz)+uyflgz(i,j) ! remplace le glissement par uyflgz ! end do ! debug_3D(i,j,10)=uydef(i,j) endif end do end do !$OMP END DO !$OMP END PARALLEL ! on ne recalcul ux que lorsque uxflgz est modifié soit tous les dtt if (isynchro.eq.1) then ! !$OMP PARALLEL PRIVATE(rang,nb_taches) ! !$ rang=OMP_GET_THREAD_NUM() ! !$ nb_taches=OMP_GET_NUM_THREADS() ! !$OMP DO SCHEDULE(STATIC,NY/nb_taches) !$OMP PARALLEL !$OMP DO do j=2,ny do i=2,nx if (flgzmx(i,j)) then do k=1,nz ux(i,j,k)=ux(i,j,k)-ux(i,j,nz)+uxflgz(i,j) ! remplace le glissement par uxflgz end do endif if (flgzmy(i,j)) then do k=1,nz uy(i,j,k)=uy(i,j,k)-uy(i,j,nz)+uyflgz(i,j) ! remplace le glissement par uyflgz end do endif enddo enddo !$OMP END DO !$OMP END PARALLEL endif else ! SIA et L1 en 2 zones separe : pas d'action (uxbar est deja remis a uxflgz dans diffusiv) endif end subroutine mix_SIA_L1 end module resolmeca_SIA_L1