source: trunk/SOURCES/mix-SIA-L1_mod.f90 @ 25

Last change on this file since 25 was 4, checked in by dumas, 10 years ago

initial import GRISLI trunk

File size: 4.0 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
56! subroutine qui associe SIA et L1
57! la methode depend de i_resolmeca
58
59
60! use module_choix
61    if (itracebug.eq.1)  call tracebug(' routine mix_SIA_L1')
62
63debug_3D(:,:,9)=0.
64debug_3D(:,:,10)=0.
65debug_3D(:,:,11)=0.
66debug_3D(:,:,12)=0.
67debug_3D(:,:,13)=ubx(:,:)
68debug_3D(:,:,14)=uby(:,:)
69
70if (i_resolmeca.eq.1) then
71
72   ! prend les vitesses de la SIA si elles s'averent plus fortes
73   ! que celles du L1. On garde le L1 comme vitesse basale.
74   ! version utilisee pour les runs MIS11 de Cairns
75
76   do j=2,ny
77      do i=2,nx
78         if (flgzmx(i,j)) then
79            debug_3D(i,j,11)=uxflgz(i,j)
80            if (abs((uxdef(i,j)+ubx(i,j))).gt.abs(uxflgz(i,j))) then   ! SIA plus grande que L1
81               uxbar(i,j)=uxflgz(i,j)+uxdef(i,j)                       ! on ajoute la deform. SIA
82               ubx(i,j) = uxflgz(i,j)
83               debug_3D(i,j,9)=uxdef(i,j)             
84            end if
85         end if
86
87         if (flgzmy(i,j)) then
88            debug_3D(i,j,12)=uyflgz(i,j)
89            if (abs((uydef(i,j)+uby(i,j))).gt.abs(uyflgz(i,j))) then   ! SIA plus grande que L1
90               uybar(i,j)=uyflgz(i,j)+uydef(i,j)                       ! on ajoute la deform. SIA
91               uby(i,j) = uyflgz(i,j)
92               debug_3D(i,j,10)=uydef(i,j)
93            end if
94         end if
95      end do
96   end do
97
98else if (i_resolmeca.eq.2) then
99
100! addition systematique
101
102   do j=2,ny
103      do i=2,nx
104         if (flgzmx(i,j)) then
105            debug_3D(i,j,11)=uxflgz(i,j)
106            uxbar(i,j)=uxflgz(i,j)+uxdef(i,j)                       ! on ajoute la deform. SIA
107            ubx(i,j) = uxflgz(i,j)
108            do k=1,nz
109               ux(i,j,k)=ux(i,j,k)-ux(i,j,nz)+uxflgz(i,j)           ! remplace le glissement par uxflgz
110            end do
111
112            debug_3D(i,j,9)=uxdef(i,j)                   
113         end if
114
115         if (flgzmy(i,j)) then
116            debug_3D(i,j,12)=uyflgz(i,j)
117            uybar(i,j)=uyflgz(i,j)+uydef(i,j)                       ! on ajoute la deform. SIA
118            uby(i,j) = uyflgz(i,j)
119            do k=1,nz
120               uy(i,j,k)=uy(i,j,k)-uy(i,j,nz)+uyflgz(i,j)           ! remplace le glissement par uyflgz
121            end do
122
123            debug_3D(i,j,10)=uydef(i,j)             
124         end if
125      end do
126   end do
127
128
129
130else ! SIA et L1 en 2 zones separe : pas d'action (uxbar est deja remis a uxflgz dans diffusiv)
131
132
133endif
134
135end subroutine mix_SIA_L1
136
137end module resolmeca_SIA_L1
Note: See TracBrowser for help on using the repository browser.