source: trunk/SOURCES/Temperature-routines/icetemp_mod.f90 @ 334

Last change on this file since 334 was 271, checked in by aquiquet, 5 years ago

Further computations of mass conservation in double precision for small timesteps

File size: 9.2 KB
Line 
1!> \file icetemp_mod.f90
2!! Interfaces of functions and subroutines for temperature calculation
3!<
4
5!> \namespace icetempmod
6!! Interfaces of functions and subroutines for the temperature calculation
7!<
8
9module icetempmod
10
11contains
12
13  !> SUBROUTINE: init_icetemp
14  !! Initialise the temperature calculation
15  !! Used modules:
16  !! - Icetemp_declar
17  !>
18
19  Subroutine init_icetemp(Num_rep_42_m)
20    Use Icetemp_declar, only : Nfracq,Iq,Ecart_phid
21
22    Implicit None
23
24    !< Arguments
25    Integer,                       Intent(In)::Num_rep_42_m !< Id Of Reponse
26    ! Character(Len=80) :: Filin
27    !________________________________________________________________
28    ! La Temperature Du Point De Congelation De L'Eau De Mer
29    ! Est Tiree De Jenkins (1991), Tb Est En Degres C
30    ! _______________________________________________________________
31
32    !  Nouvelle Formulation De Fracq
33    !  Nfracq=15  ! Ancienne Formulation
34
35    Nfracq=1
36
37!    Write(Num_rep_42_m,*) 'Calcul Des Temperatures'
38!    Write(Num_rep_42_m,*) 'Nfracq=',Nfracq
39
40    Iq=1
41!    If (Iq.Eq.1) Write(Num_rep_42_m,*) 'Chaleur Demi Maille, Iq=',Iq
42!    If (Iq.Eq.2) Write(Num_rep_42_m,*) 'Chaleur Centree, Iq=',Iq
43!    If (Iq.Eq.3) Write(Num_rep_42_m,*) 'Chaleur Pente, Iq=',Iq
44!    If (Iq.Eq.4) Write(Num_rep_42_m,*) 'Chaleur U_taub, Iq=',Iq
45!    If (Iq.Eq.5) Write(Num_rep_42_m,*) 'Chaleur U_taub_stag, Iq=',Iq
46!    If (Iq.Eq.6) Write(Num_rep_42_m,*) 'Chaleur Q_sia_stag,Iq=',Iq
47!    If (Iq.Eq.7) Write(Num_rep_42_m,*) 'Chaleur Q_all_stag,Iq=',Iq
48!    Write(Num_rep_42_m,*)
49
50    Ecart_phid=0.5
51!    Write(Num_rep_42_m,*) ' Si Base Froide, La Chaleur De Glissement Est En Exp(Delta T * Ecart_phid)'
52!    Write(Num_rep_42_m,*) 'Ecart_phid=',Ecart_phid
53!    Write(Num_rep_42_m,*) '---------------------------------------------------------------------------'
54
55  End Subroutine init_icetemp
56 
57
58  !-----------------------------------------------------------------------------------------------------------------------
59
60  !> SUBROUTINE: icetemp
61  !! Calculate the temperature
62  !!
63  !! Used modules:
64  !!  -  use Icetemp_declar
65  !>
66
67  subroutine icetemp
68    !$ use OMP_LIB
69    use Icetemp_declar
70
71
72    Implicit None 
73   
74!~   integer :: t1,t2,ir
75!~   real                           :: temps, t_cpu_0, t_cpu_1, t_cpu, norme 
76!~   
77!~    ! Temps CPU de calcul initial.
78!~    call cpu_time(t_cpu_0)
79!~    ! Temps elapsed de reference.
80!~    call system_clock(count=t1, count_rate=ir)
81   
82   
83         !$OMP PARALLEL
84         !$OMP WORKSHARE
85    ! init
86    Aa=0.   
87    Bb=0.             
88    Cc=1.             
89    Rr=0.             
90    Hh=0.             
91    Tdot=0
92    Chal2_x=0.
93    Chal2_y=0.
94    Chal2_z=0.
95    Chal2_xy=0.
96    Chaldef_maj=0.
97    Chalglissx=0.
98    Chalglissy=0.
99   
100    !Instruction
101    Dee=1./(Nz-1) 
102    Fracq=(1.-(1.-Dee/2.)**Nfracq)/Nfracq
103    Fracq=(1.-(1.-Dee/2.)**5.)/5.
104    Ee(1)=0.
105    Ee(Nz)=1.
106    !$OMP END WORKSHARE
107   
108    !$OMP DO
109    Do K=1,Nz
110       If ((K.Ne.1).And.(K.Ne.Nz)) Ee(K)=(K-1.)/(Nz-1.)
111    End Do
112    !$OMP END DO
113    !$OMP END PARALLEL
114
115    !   Variables Dependant Du Pas De Temps Dtt
116    Ctm=Dtt*Cm/Rom/Cpm/Dzm/Dzm
117    Dttdx=Dtt/Dx
118    Dx11=1./Dx
119
120    If (Itracebug.Eq.1) Write(num_tracebug,*)' Entree Dans Routine Ice_temp Time=',Time
121
122    !________________________________________________________________
123    ! La Temperature Du Point De Congelation De L'Eau De Mer
124    ! Est Tiree De Jenkins (1991), Tb Est En Degres C
125    ! _______________________________________________________________
126    If (Itracebug.Eq.1) Write(num_tracebug,*)' avant appel  Routine Temp_mer'
127    Call Temp_mer
128
129    ! Appel Aux Proprietes Thermiques De La Glace : Conductivite, Capacite Claorifique -Tpmp
130
131    If (Itracebug.Eq.1) Write(num_tracebug,*)' avant appel  Routine Thermal_property'
132    Call Thermal_prop_icetemp
133
134    If (Itracebug.Eq.1) Write(num_tracebug,*)' avant appel  Routine Qprod'
135    Call Qprod_ice(Iq)
136
137
138    ! Dans Le Socle : Les elements de La Matrice Tridiag Sont Invariants Dans L'Espace
139    !$OMP PARALLEL
140    !$OMP DO
141    Do K=Nz+1,nz+nzm-1
142       Aa(K)=-Ctm
143       Bb(K)=1.+2.*Ctm
144       Cc(K)=-Ctm
145    End Do
146         !$OMP END DO
147
148    ! Conditions Aux Limites Pour Le Socle
149    Aa(nz+nzm)=-1.
150    Bb(nz+nzm)=1.
151    Cc(nz+nzm)=0.
152
153    ! Calcul D'Advection : On Traite Les Tableaux En Une Seule Ligne De Commandes
154    ! les tableaux iadvec valent 1 quand la vitesse consideree est upwind et 0 sinon
155
156    !         Advection Selon X
157    !         -----------------
158    !$OMP WORKSHARE
159    Iadvec_w(:,:)=Nint(Max(Sign(1.0,Uxbar(:,:)),0.0))
160    Iadvec_e(:,:)=Nint(Abs(Min(Sign(1.0,Uxbar(:,:)),0.0)))
161
162    !         Advection Selon Y
163    !         -----------------
164
165    Iadvec_s(:,:)=Nint(Max(Sign(1.0,Uybar(:,:)),0.0)) 
166    Iadvec_n(:,:)=Nint(Abs(Min(Sign(1.0,Uybar(:,:)),0.0)))
167    !$OMP END WORKSHARE
168    ! Boucle De Resolution Des Temperatures Colonne Par Colonne
169         !$OMP SINGLE
170    If (Itracebug.Eq.1) Write(num_tracebug,*)' avant appel  boucle Advec_icetemp'
171    !$OMP END SINGLE
172
173         !$OMP DO PRIVATE(Deh22,Duu,Aa,Bb,Cc,Rr,Hh,Dou)
174    Do J=2,Ny-1
175       Do I=2,Nx-1
176
177 !  If (Itracebug.Eq.1) Write(num_tracebug,*)' avant appel  Advect_icetemp  i,j', i,j
178          call Advec_icetemp(Nz,Nzm,Nn,Chaldef_maj(I,J,:),Cp(I,J,:),& 
179               Advecx(I,J,:),Advecy(I,J,:),Advec(I,J,:),Ct(I,J,:),Iadvec_w(I,J),Iadvec_e(I+1,J),Iadvec_s(I,J),Iadvec_n(I,J+1), &
180                                        Ux(i,j,:),Ux(i+1,j,:),Uy(i,j,:),Uy(i,j+1,:),T(i,j,:),T(i-1,j,:),T(i+1,j,:),T(i,j-1,:),T(i,j+1,:), &
181                                        real(H(i,j)),Ts(i,j),Aa,Bb,Cc,Rr,Dx11,Dou,DTT,Dee,Uzr(i,j,:),Dttdx)
182
183          ! variables de calcul dans la glace
184          Deh22=2.*Dee*H(I,J)**2
185          Duu=Dtt/2./Dee
186          Rr(nz+nzm)=-Dzm*Ghf(I,J)/Cm  ! depend de I,J Pour Un Flux Geothermique Variable
187
188!   If (Itracebug.Eq.1) Write(num_tracebug,*)' avant appel  temp_col  i,j', i,j
189          Call temp_col(Nz,Nzm,Nn,T3d_new(i,j,:),Ct(i,j,:),&
190               flot(i,j),real(H(i,j)),Tbmer(i,j),Ibase(i,j),T(i,j,:),Tpmp(i,j,:),Aa,Bb,Cc,Rr,Hh,Ncond,Dee,Cm,&
191               Dzm,Phid(i,j),Ghf(i,j),Ifail1,Tbdot(i,j),DTT,Dou,Ts(i,j))
192
193          ! Test Permettant D'Ecrire Les Variables Quand Il Y A Une Erreur
194          ! Ifail1 = 1 : Erreur Dans Tridiag
195          If (Ifail1.Eq.1) Then
196             Write(6,*)'Erreur Dans Tridiag'
197             Write(6,*)'I = ',I,' J = ',J
198             Write(6,*) 'A, B, C, R, U'
199             Do K=1,nz+nzm
200                Write(6,*)'K = ',K
201                !Write(6,*) Aa(K),Bb(K),Cc(K),Rr(K),Hh(K)
202             End Do
203             Stop
204          End If
205       End Do
206    End Do
207         !$OMP END DO
208
209
210         !$OMP DO
211    Do J=1,Ny
212       Do I=1,Nx
213          H1(I,J)=H(I,J)
214          B1(I,J)=B(I,J)
215          Tpmp(I,J,1)=0.
216       End Do
217    End Do
218    !$OMP END DO
219
220    ! Attribution Finale De La Temperature
221    ! Limitation a 3 deg De Variation Par Dtt Et Pas Plus Froid Que -70
222         !$OMP DO COLLAPSE(2) PRIVATE(Tdot)
223    Do K=1,nz+nzm
224       Do J=1,Ny
225          Do I=1,Nx
226             If (H(I,J).Gt.10.) Then
227
228                Tdot(K)=T3d_new(I,J,K)-T(I,J,K)
229
230                If ((K.Gt.1).And.(Tdot(K).Lt.-3.)) Tdot(K)=-3.
231                If ((K.Gt.1).And.(Tdot(K).Gt.3.)) Tdot(K)=3.
232
233                T(I,J,K)=T(I,J,K)+Tdot(K)
234
235                If (T(I,J,K).Lt.-70.) Then
236                   T(I,J,K)=-70.
237                Endif
238             Endif
239
240          End Do
241       End Do
242    End Do
243    !$OMP END DO
244
245         !$OMP DO COLLAPSE(2)
246    Do K=2,Nz
247       Do J=1,Ny
248          Do I=1,Nx
249             If (T(I,J,K).Gt.Tpmp(I,J,K)) Then
250                T(I,J,K)=Tpmp(I,J,K)
251                Ibase(I,J)=2
252             Endif
253          End Do
254       End Do
255    End Do
256    !$OMP END DO
257
258
259    !     Temperature Sur Les Bords :
260    !     -------------------------
261 
262    !$OMP DO COLLAPSE(2)
263    Do K=1,Nz
264       Do J=1,Ny
265          !        T(Nx,J,K)= T(Nx-1,J,K)
266          !        T(1,J,K)= T(2,J,K)
267          T(Nx,J,K)= Ts(Nx,J)
268          T(1,J,K)= Ts(1,J)
269       End Do
270       !
271       Do I=1,Nx
272          !        T(I,1,K)=T(I,2,K)
273          !        T(I,Ny,K)=T(I,Ny-1,K)
274          T(I,1,K)=Ts(I,1)
275          T(I,Ny,K)=Ts(I,Ny)
276       End Do
277    End Do
278    !$OMP END DO
279
280    !$OMP DO COLLAPSE(2)
281    Do K=1,Nz
282       Do J=1,Ny
283          Do I=1,Nx
284             If (H(I,J).le.1.) then
285                T(I,J,K)=Ts(I,J)
286                if ((k.ge.nz).and.(flot(i,j))) T(I,J,K) =Tbmer(I,J)
287             end If
288               
289          End Do
290       End Do
291    End Do
292    !$OMP END DO
293    !$OMP END PARALLEL
294 
295!~   ! Temps elapsed final
296!~   call system_clock(count=t2, count_rate=ir)
297!~   temps=real(t2 - t1,kind=4)/real(ir,kind=4)
298!~   ! Temps CPU de calcul final
299!~   call cpu_time(t_cpu_1)
300!~   t_cpu = t_cpu_1 - t_cpu_0
301 
302!~   ! Impression du resultat.
303!~   print '(//,3X,"Valeurs de nx et ny : ",I5,I5/,              &
304!~            & 3X,"Temps elapsed       : ",1PE10.3," sec.",/, &
305!~            & 3X,"Temps CPU           : ",1PE10.3," sec.",/, &
306!~            & 3X,"Norme (PB si /= 0)  : ",1PE10.3,//)', &
307!~            nx,ny,temps,t_cpu,norme
308
309    where (flot(:,:))
310      debug_3D(:,:,122) = T(:,:,nz)+273.15
311    elsewhere
312      debug_3D(:,:,120) = T(:,:,nz)+273.15
313    endwhere
314   
315    If (Itracebug.Eq.1) Write(num_tracebug,*)' fin routine icetemp'
316    return
317  End Subroutine icetemp
318
319  !-----------------------------------------------------------------------------------------------------------------------
320
321end module icetempmod
Note: See TracBrowser for help on using the repository browser.