source: branches/iLoveclim/SOURCES/mix-SIA-L1_mod.f90 @ 146

Last change on this file since 146 was 146, checked in by aquiquet, 7 years ago

Grisli-iLoveclim branch: merged to trunk at revision 145

File size: 5.8 KB
Line 
1!> \file mix-SIA-L1_mod.f90
2!! Module qui defini le type d'association entre L1 et SIA
3!<
4
5!> \namespace resolmeca_SIA_L1
6!! Module qui defini le type d'association entre L1 et SIA
7!! @note Test sur la variable globale i_resolmeca
8!! \author ...
9!! \date ...
10!! @note Used module
11!! @note   - use module3D_phy
12!!
13!<
14module resolmeca_SIA_L1
15
16! module qui défini le type d'association entre L1 et SIA
17! test sur la variable globale i_resolmeca
18
19use module3D_phy
20
21contains
22!> SUBROUTINE: init_resol_meca
23!! Routine pour l'initialisation de la resolution mecanique
24!>
25subroutine init_resol_meca
26
27namelist/meca_SIA_L1/i_resolmeca
28! formats pour les ecritures dans 42
29
30428 format(A)
31
32! lecture des parametres du run                      block eaubasale1
33!--------------------------------------------------------------------
34rewind(num_param)        ! pour revenir au debut du fichier param_list.dat
35read(num_param,meca_SIA_L1)
36
37
38write(num_rep_42,428)'!___________________________________________________________' 
39write(num_rep_42,428) '&meca_SIA_L1             ! bloc resol_meca '
40write(num_rep_42,*)
41write(num_rep_42,*) 'i_resolmeca = ', i_resolmeca
42write(num_rep_42,*)'/'                           
43write(num_rep_42,428) '! i_resolmeca type d association entre SIA et L1'
44write(num_rep_42,428) '! i_resolmeca=0  chacun dans sa zone'
45write(num_rep_42,428) '! i_resolmeca=1  dans les zones stream, addition si uxdef > uxL1 (MIS11 Cairns)'
46write(num_rep_42,428) '! i_resolmeca=2  addition systematique dans les zones stream '
47write(num_rep_42,*)
48
49end subroutine init_resol_meca
50!> SUBROUTINE: mix_SIA_L1
51!! Subroutine qui associe SIA et L1: la methode depend de i_resolmeca
52!>
53subroutine mix_SIA_L1
54
55!$ USE OMP_LIB
56! variables parallelisation openMP
57!!$  integer :: rang ,nb_taches
58
59! subroutine qui associe SIA et L1
60! la methode depend de i_resolmeca
61
62
63! use module_choix
64if (itracebug.eq.1)  call tracebug(' routine mix_SIA_L1')
65
66debug_3D(:,:,9)=0.
67debug_3D(:,:,10)=0.
68debug_3D(:,:,11)=0.
69debug_3D(:,:,12)=0.
70debug_3D(:,:,13)=ubx(:,:)
71debug_3D(:,:,14)=uby(:,:)
72
73if (i_resolmeca.eq.1) then
74
75   ! prend les vitesses de la SIA si elles s'averent plus fortes
76   ! que celles du L1. On garde le L1 comme vitesse basale.
77   ! version utilisee pour les runs MIS11 de Cairns
78        !$OMP PARALLEL
79   !$OMP DO
80   do j=2,ny
81      do i=2,nx
82         if (flgzmx(i,j)) then
83            debug_3D(i,j,11)=uxflgz(i,j)
84            if (abs((uxdef(i,j)+ubx(i,j))).gt.abs(uxflgz(i,j))) then   ! SIA plus grande que L1
85               uxbar(i,j)=uxflgz(i,j)+uxdef(i,j)                       ! on ajoute la deform. SIA
86               ubx(i,j) = uxflgz(i,j)
87               debug_3D(i,j,9)=uxdef(i,j)             
88            end if
89         end if
90
91         if (flgzmy(i,j)) then
92            debug_3D(i,j,12)=uyflgz(i,j)
93            if (abs((uydef(i,j)+uby(i,j))).gt.abs(uyflgz(i,j))) then   ! SIA plus grande que L1
94               uybar(i,j)=uyflgz(i,j)+uydef(i,j)                       ! on ajoute la deform. SIA
95               uby(i,j) = uyflgz(i,j)
96               debug_3D(i,j,10)=uydef(i,j)
97            end if
98         end if
99      end do
100   end do
101   !$OMP END DO
102        !$OMP END PARALLEL
103
104else if (i_resolmeca.eq.2) then
105
106! addition systematique
107!       !$OMP PARALLEL PRIVATE(rang,nb_taches)
108!       !$ rang=OMP_GET_THREAD_NUM()
109!       !$ nb_taches=OMP_GET_NUM_THREADS()
110!       !$OMP DO SCHEDULE(STATIC,NY/nb_taches)
111        !$OMP PARALLEL
112   !$OMP DO
113   do j=2,ny
114      do i=2,nx
115! test sur tout le tableau :
116         if (flgzmx(i,j)) then 
117!            debug_3D(i,j,11)=uxflgz(i,j)
118            uxbar(i,j)=uxflgz(i,j)+uxdef(i,j)                       ! on ajoute la deform. SIA
119            ubx(i,j) = uxflgz(i,j)
120!            do k=1,nz
121!               ux(i,j,k)=ux(i,j,k)-ux(i,j,nz)+uxflgz(i,j)           ! remplace le glissement par uxflgz
122!            end do
123
124!            debug_3D(i,j,9)=uxdef(i,j)                   
125         endif
126
127         if (flgzmy(i,j)) then 
128!            debug_3D(i,j,12)=uyflgz(i,j)
129            uybar(i,j)=uyflgz(i,j)+uydef(i,j)                       ! on ajoute la deform. SIA
130            uby(i,j) = uyflgz(i,j)
131!            do k=1,nz
132!               uy(i,j,k)=uy(i,j,k)-uy(i,j,nz)+uyflgz(i,j)           ! remplace le glissement par uyflgz
133!            end do
134
135!            debug_3D(i,j,10)=uydef(i,j)             
136         endif
137      end do
138   end do
139        !$OMP END DO
140        !$OMP END PARALLEL
141
142! on ne recalcul ux que lorsque uxflgz est modifié soit tous les dtt
143        if (isynchro.eq.1) then
144!       !$OMP PARALLEL PRIVATE(rang,nb_taches)
145!       !$ rang=OMP_GET_THREAD_NUM()
146!       !$ nb_taches=OMP_GET_NUM_THREADS()
147!       !$OMP DO SCHEDULE(STATIC,NY/nb_taches)
148        !$OMP PARALLEL
149   !$OMP DO
150        do j=2,ny
151                do i=2,nx
152        if (flgzmx(i,j)) then
153                        do k=1,nz
154               ux(i,j,k)=ux(i,j,k)-ux(i,j,nz)+uxflgz(i,j)           ! remplace le glissement par uxflgz
155            end do
156         endif
157         if (flgzmy(i,j)) then
158                        do k=1,nz
159               uy(i,j,k)=uy(i,j,k)-uy(i,j,nz)+uyflgz(i,j)           ! remplace le glissement par uyflgz
160            end do
161         endif
162       enddo
163   enddo
164        !$OMP END DO
165        !$OMP END PARALLEL
166endif
167
168else ! SIA et L1 en 2 zones separe : pas d'action (uxbar est deja remis a uxflgz dans diffusiv)
169
170
171endif
172
173
174! afq -- initMIP outputs:
175
176  !$OMP WORKSHARE
177  debug_3d(:,:,111) = ( (ux  (:,:,1)  + eoshift( ux(:,:,1) ,shift=1,boundary=0.,dim=1))/2.  )/secyear
178  debug_3d(:,:,112) = ( (uy  (:,:,1)  + eoshift( uy(:,:,1) ,shift=1,boundary=0.,dim=2))/2.  )/secyear
179  debug_3d(:,:,113) = ( uzr (:,:,1) * ice(:,:)                                              )/secyear
180  debug_3d(:,:,114) = ( (ux  (:,:,nz) + eoshift( ux(:,:,nz),shift=1,boundary=0.,dim=1))/2.  )/secyear
181  debug_3d(:,:,115) = ( (uy  (:,:,nz) + eoshift( uy(:,:,nz),shift=1,boundary=0.,dim=2))/2.  )/secyear
182  debug_3d(:,:,116) = ( uzr (:,:,nz) * ice(:,:)                                             )/secyear
183  !$OMP END WORKSHARE
184
185
186end subroutine mix_SIA_L1
187
188end module resolmeca_SIA_L1
Note: See TracBrowser for help on using the repository browser.