source: trunk/NS3D_JMC/sous_routines.F90 @ 15

Last change on this file since 15 was 15, checked in by yculod, 17 years ago

Correction assignation dRho_dz dans read_dRho, variation remise suivant l'axe z

File size: 78.0 KB
Line 
1#include "config.h"
2
3module sous_routines
4
5    use parametres
6    use fft
7
8    implicit none
9
10    contains
11
12!**********************************************************************************
13!**********************************************************************************
14!**********    SOUS ROUTINES D'ENTREE SORTIES                            **********
15!**********************************************************************************
16!**********************************************************************************
17
18! ***************************************************************************
19        subroutine enter_param
20! ***************************************************************************
21!       lecture des tous les parametres qui sont dans le fichier 'data.in'
22! ***************************************************************************
23
24            character(len=29) :: bidon
25            character(len=100) :: lecture_reel, lecture_entier, lecture_logical
26            integer :: i, nx,ny,nz
27
28            lecture_reel = '(A29,E14.5)'
29            lecture_entier = '(A29,I10)'
30            lecture_logical = '(A29,L12)'
31
32            open (unit=1,file='data.in')
33
34
35
36            read(1,*)
37            read(1,*)
38            read(1,*)
39            read(1,*)
40
41            ! dimension du probleme (ces valeurs ne sont en fait pas reellement utilisees).
42            ! les dimensions du probleme sont definies dans le fichier parametres.F90 par le preprocesseur
43            ! en fonction des valeurs indiquees dans le fichier config.h
44            read(1,lecture_entier) bidon,nx
45            read(1,lecture_entier) bidon,ny
46            read(1,lecture_entier) bidon,nz
47
48            ! taille suivant x
49            read(1,lecture_reel) bidon,lx
50            read(1,lecture_reel) bidon,ly
51            read(1,lecture_reel) bidon,lz
52
53             ! Pas de temps et instant initial
54            read(1,lecture_reel) bidon,dt
55            read(1,lecture_reel) bidon,begin
56
57            ! nombre d'itérations
58            read(1,lecture_entier) bidon,itmax
59
60            ! Troncature
61            read(1,lecture_logical) bidon,truncate
62            read(1,lecture_entier) bidon,type_trunc
63            read(1,lecture_reel) bidon,rtrunc_x
64            read(1,lecture_reel) bidon,rtrunc_y
65            read(1,lecture_reel) bidon,rtrunc_z
66
67            read(1,*)
68            read(1,*)
69
70             ! Viscosité
71            read(1,lecture_reel) bidon,nu
72
73            ! stratification, fréquence de Brunt Vaisala, vorticité planétaire et nombre de Schmidt
74            read(1,lecture_logical) bidon,stratification
75            read(1,lecture_reel) bidon,xns
76            xnsqr=xns*xns
77            read(1,lecture_reel) bidon,schm
78            read(1,lecture_reel) bidon,omega2
79
80            read(1,*)
81            read(1,*)
82
83            ! type de run
84            read(1,lecture_logical) bidon,perturbe
85            read(1,lecture_logical) bidon,lin
86
87            read(1,*)
88            read(1,*)
89
90            ! reprise de run
91            read(1,lecture_logical) bidon,reprise_run
92
93            read(1,*)
94            read(1,*)
95
96            ! etat de base 2D
97            read(1,lecture_logical) bidon, base2D_sub
98            read(1,lecture_logical) bidon, base2D_read
99
100            read(1,*)
101            read(1,*)
102
103            ! vitesse initiale 2D
104            read(1,lecture_logical) bidon, velo2D_sub
105            read(1,lecture_logical) bidon, velo2D_read
106            read(1,lecture_logical) bidon,lnoise
107            read(1,lecture_reel) bidon,amplitude_noise
108
109            read(1,*)
110            read(1,*)
111
112             
113            !profil de densite
114            read(1,lecture_logical) bidon, dRho_read
115            read(1,lecture_reel) bidon, Rhob
116             
117            read(1,*)
118            read(1,*)
119
120            ! pba : deb
121            ! forcage
122            read(1,lecture_logical) bidon, forcage_sub
123            read(1,lecture_reel) bidon, A_ics
124            read(1,lecture_reel) bidon, freq
125
126            read(1,*)
127            read(1,*)
128
129            !pba fin
130
131            ! sorties
132            read(1,lecture_entier) bidon,periode_out_physique1
133            read(1,lecture_entier) bidon,periode_out_physique2
134            read(1,lecture_entier) bidon,periode_out_spectral1
135
136            ! parametres additionnels
137            read(1,*)
138            read(1,*)
139
140            ! Tourbillons gaussiens
141            read(1,*)
142            read(1,*)
143
144            read(1,lecture_entier) bidon,nb_tourbillons
145            do i=1,nb_tourbillons
146                read(1,*)
147                read(1,lecture_reel) bidon,xt(i)
148                read(1,lecture_reel) bidon,yt(i)
149                read(1,lecture_reel) bidon,gamma(i)
150                read(1,lecture_reel) bidon,at(i)
151            end do
152
153            close(1)
154
155        end subroutine
156
157! ***************************************************************************
158        subroutine out_param
159! ***************************************************************************
160!       affichage des parametres du run
161! ***************************************************************************
162
163            implicit none
164
165            character(len=100) :: aff_reel, aff_entier, aff_logical
166            integer :: i
167
168            ! descripteur des formats de sortie
169            aff_reel = '(A29,TR3,G12.5)'
170            aff_entier = '(A29,TR4,I7)'
171            aff_logical = '(A29,TR4,L1)'
172
173            write(*,*)
174            write(*,*) '---------------------------------------------------'
175            write(*,*) ' PROGRAMME NS3D'
176            write(*,*) '---------------------------------------------------'
177            write(*,*)
178            write(*,*)
179            ! Mode de calcul
180            if (lin) then
181                write(*,*) '~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~'
182                write(*,*) '@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@'
183                write(*,*) '@@@                                             @@@'
184                write(*,*) '        VERSION PERTURBATIVE  / LINEARISEE         '
185                write(*,*) '@@@                                             @@@'
186                write(*,*) '@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@'
187                write(*,*) '~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~'
188            else
189                 if (perturbe) then
190                    write(*,*) '~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~'
191                    write(*,*) '@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@'
192                    write(*,*) '@@@                                             @@@'
193                    write(*,*) '      VERSION PERTURBATIVE  / NON LINEARISEE       '
194                    write(*,*) '@@@                                             @@@'
195                    write(*,*) '@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@'
196                    write(*,*) '~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~'
197                 else
198                    write(*,*) '~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~'
199                    write(*,*) '@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@'
200                    write(*,*) '@@@                                             @@@'
201                    write(*,*) '            VERSION NON-PERTURBATIVE               '
202                    write(*,*) '@@@                                             @@@'
203                    write(*,*) '@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@'
204                    write(*,*) '~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~'
205                 end if
206            end if
207
208            write(*,*)
209            write(*,*)
210            write(*,*) '---------------------------------------------------'
211            write(*,*) 'PARAMETRES DE DISCRETISATION'
212            write(*,*) '---------------------------------------------------'
213            write(*,aff_entier) 'nx...........................',nxmax
214            write(*,aff_entier) 'ny...........................',nymax
215            write(*,aff_entier) 'nz...........................',nzmax
216            write(*,aff_reel) 'Lx...........................',lx
217            write(*,aff_reel) 'Ly...........................',ly
218            write(*,aff_reel) 'Lz...........................',lz
219            write(*,aff_reel) 'Dt...........................',dt
220            write(*,aff_reel) 'begin........................',begin
221            write(*,aff_entier) 'itmax........................',itmax
222            write(*,aff_logical) 'Troncature...................',truncate
223            write(*,aff_entier) '    carree:1 / elliptique:2..',type_trunc
224            write(*,aff_reel) '    rayon troncature x.......',rtrunc_x
225            write(*,aff_reel) '    rayon troncature y.......',rtrunc_y
226            write(*,aff_reel) '    rayon troncature z........',rtrunc_z
227
228            write(*,*)
229            write(*,*) '---------------------------------------------------'
230            write(*,*) 'PARAMETRES PHYSIQUES'
231            write(*,*) '---------------------------------------------------'
232            write(*,aff_reel) 'viscosite....................',nu
233            write(*,aff_logical) 'stratification...............',stratification
234            write(*,aff_reel) '    brunt_vaisala_frequency..',xns
235            write(*,aff_reel) '    Schmidt number...........',schm
236            write(*,aff_reel) 'Omega2.......................',omega2
237
238            write(*,*)
239            write(*,*) '---------------------------------------------------'
240            write(*,*) 'TYPE DE RUN'
241            write(*,*) '---------------------------------------------------'
242            write(*,aff_logical) 'en perturbation..............',perturbe
243            write(*,aff_logical) 'lineaire.....................',lin
244
245            write(*,*)
246            write(*,*) '---------------------------------------------------'
247            write(*,*) 'REPRISE DE RUN'
248            write(*,*) '---------------------------------------------------'
249            write(*,aff_logical) 'reprise de run...............',reprise_run
250
251            write(*,*)
252            write(*,*) '---------------------------------------------------'
253            write(*,*) 'ETAT DE BASE 2D'
254            write(*,*) '---------------------------------------------------'
255            write(*,aff_logical) 'call subroutine..............',base2D_sub
256            write(*,aff_logical) 'Lecture fichier..............',base2D_read
257
258
259            write(*,*)
260            write(*,*) '---------------------------------------------------'
261            write(*,*) 'VITESSE INITIALE 2D'
262            write(*,*) '---------------------------------------------------'
263            write(*,aff_logical) 'call subroutine..............',velo2D_sub
264            write(*,aff_logical) 'Lecture fichier..............',velo2D_read
265            write(*,aff_logical) 'Bruit blanc..................',lnoise
266            write(*,aff_reel) '    amplitude................',amplitude_noise
267           
268            write(*,*)
269            write(*,*) '---------------------------------------------------'
270            write(*,*) 'DENSITE INITIALE 1D'
271            write(*,*) '---------------------------------------------------'
272            write(*,aff_logical) 'Lecture fichier..............',dRho_read
273            write(*,aff_reel) 'densite de base..............',Rhob
274               
275       
276            write(*,*)
277            write(*,*)'---------------------------------------------------'
278            write(*,*)'FORCAGE EN VITESSE'
279            write(*,*)'---------------------------------------------------'
280            write(*,aff_logical)'call subroutine', forcage_sub
281            write(*,aff_reel)'Amplitude..................',A_ics
282            write(*,aff_reel)'Frequence..................',freq
283
284            write(*,*)
285            write(*,*) '---------------------------------------------------'
286            write(*,*) 'SORTIES'
287            write(*,*) '---------------------------------------------------'
288            write(*,aff_entier) 'periode out_physique1.........',periode_out_physique1
289            write(*,aff_entier) 'periode out_physique2........',periode_out_physique2
290            write(*,aff_entier) 'periode out_spectral1........',periode_out_spectral1
291
292            write(*,*)
293            write(*,*) '---------------------------------------------------'
294            write(*,*) 'PARAMETRES ADDITIONNELS'
295            write(*,*) '---------------------------------------------------'
296            write(*,*)
297            write(*,*) '** Tourbillons gaussiens **'
298            write(*,aff_entier) 'nb tourbillons...............',nb_tourbillons
299            do i=1,nb_tourbillons
300                write(*,aff_entier) 'Tourbillon...................',i
301                write(*,aff_reel) '    xt.......................',xt(i)
302                write(*,aff_reel) '    xt.......................',yt(i)
303                write(*,aff_reel) '    gamma....................',gamma(i)
304                write(*,aff_reel) '    at.......................',at(i)
305            end do
306
307            write(*,*)
308            write(*,*)
309
310        end subroutine
311
312
313
314! **********************************************************************
315        subroutine out_physique1(vox,voy,voz)
316! **********************************************************************
317!       Sortie effectuee de façon periodique et qui prend en argument le champ de vitesse
318!       dans l'espace physique.
319!       Ici on affiche le taux de croissance des perturbations ainsi que le temps cpu ecoule et
320!       une estimation du temps restant
321! **********************************************************************
322
323            implicit none
324
325            double precision, dimension(0:dimx_tab-1,0:dimy_tab-1,0:dimz_tab-1), intent(in) :: vox, voy, voz
326            double precision :: vquad, growthrate
327            double precision, save :: tps_ecoule = 0, tps_ecoule_prec_iteration = 0
328
329
330            ! on calcule le temps cpu ecoule depuis la rentree dans la boucle et une estimation du temps
331            ! restant (on se base sur le temps necessaire a la derniere iteration de la boucle)
332            call cpu_time(tps_temp)
333            tps_ecoule = tps_total + tps_temp
334
335            write(*,*)
336            write(*,'(A4,TR2,I6,TR3,A6,TR2,F10.3)') 'it =', it, 'time =', it*dt + begin
337                                                write(*,'(A37,TR2,I10,TR3,I10)') '   temps cpu en sec (ecoule/restant) =', int(tps_ecoule), &
338                                                           int((tps_ecoule-tps_ecoule_prec_iteration)/periode_out_physique1*(itmax-it))
339
340            tps_ecoule_prec_iteration = tps_ecoule
341
342
343            ! on calcule la vitesse quadratique moyenne des perturbations
344            vquad = sum(vox(0:nxmax-1, 0:nymax-1, 0:nzmax-1)**2 + voy(0:nxmax-1, 0:nymax-1, 0:nzmax-1)**2 + &
345                    voz(0:nxmax-1, 0:nymax-1, 0:nzmax-1)**2)
346
347            write(*,'(A32,TR2,G23.10)') '   vitesse quadratique moyenne =',vquad
348
349            ! on calcule le taux de croissance des perturbations entre deux appels à la fonction out_physique1
350            if (vquadanc>epsilo) then
351                growthrate=0.5*log(vquad/vquadanc)/(periode_out_physique1*dt)
352            else
353                growthrate=0.
354            endif
355
356            write(*,'(A34,TR2,G23.10)') '   taux de croissance instantane =', growthrate
357
358            vquadanc = vquad
359
360        end subroutine
361
362! *************************************************************************
363        subroutine out_physique2(vox,voy,voz,rho,fx,fy,fz)
364! *************************************************************************
365!       Sortie effectuee de façon periodique et qui prend en argument le champ de vitesse
366!       la densite et la vorticite dans l'espace physique.
367! **************************************************************************
368            double precision, dimension(0:dimx_tab-1,0:dimy_tab-1,0:dimz_tab-1), intent(in) :: vox,voy,voz,rho
369            double precision, dimension(0:dimx_tab-1,0:dimy_tab-1,0:dimz_tab-1), intent(in) :: fx,fy,fz
370
371            ! là on appelle simplement, une fonction qui écrit les champs dans un fichier.
372            call out_champs(vox,voy,voz,rho,fx,fy,fz)
373
374        end subroutine
375
376
377! *************************************************************************
378        subroutine out_spectral1(ax,ay,az)
379! *************************************************************************
380!       Fonctions qui sauvegarde l'energie quadratique moyenné suivant ky dans l'espace spectral
381!       On fait d'abord une moyenne suivant ky et on enregistre l'energie pour chaque mode kx et kz
382! *************************************************************************
383
384            implicit none
385
386            double precision, dimension(2,0:dimx_tab/2-1,0:dimy_tab-1,0:dimz_tab-1), intent(in) :: ax, ay, az
387            double precision, dimension(0:kxmax, 0:nzmax-1) :: energie_quadratique_mode
388            integer :: ix, iz
389
390
391            do iz=0,nzmax-1
392                do ix=0,kxmax
393                    energie_quadratique_mode(ix, iz) = sum(ax(1,ix,0:nymax-1,iz)**2+ax(2,ix,0:nymax-1,iz)**2+ &
394                        ay(1,ix,0:nymax-1,iz)**2+ ay(2,ix,0:nymax-1,iz)**2+ &
395                        az(1,ix,0:nymax-1,iz)**2+az(2,ix,0:nymax-1,iz)**2)/nymax
396                end do
397            end do
398
399
400            ! on l'écrit dans un fichier
401
402            ! si c'est le premier appel à la fonction out_spectral1, on cree le fichier de sortie et on ecrit toutes les
403            ! metadonnees.
404            ! On pourrait eventuellement tester si le fichier existe deja inquire(file="energie_quadratique_mode", EXIST=existe)
405            ! Cela permettrait de completer un fichier de sortie genere lors d'une precedente simulation...
406
407            ! S'il ne s'agit pas du premier appel, on ecrit simplement l'instant et les valeurs de l'energie_quadratique_mode
408            ! dans le fichier de sortie
409
410            ! premier appel a periode_out_spectral1: on cree le fichier et on ecrit les metadonnees
411            if (it == periode_out_spectral1) then
412                open(unit=76, file='energie_quadratique_mode',form='unformatted',status="replace",action="write")
413                write(76) 1, nxmax, nymax, nzmax, lx, ly, lz, dt, truncate, type_trunc, rtrunc_x, rtrunc_y, rtrunc_z, &
414                         nu, stratification, xns, schm, omega2, perturbe, lin
415                close(76)
416            end if
417
418            ! ce n'est pas le premier appel : on ajoute a la fin du fichier existant l'instant et les valeurs de l'energie
419            ! Si on ajoute l'argument action="write" dans la foncion write avec le compilateur Absoft Pro Fortran 8.2 sous Mac
420            ! l'enregistrement ne se fait pas correctement alors que ca marche sans probleme avec le compilateur Intel Fortran
421            ! ou le compilateur Fortran du NEC de l'IDRIS, bizarre, bizarre...
422            open(unit=76, file='energie_quadratique_mode',form='unformatted',status="old",position="append")
423            write(76) it*dt + begin
424            write(76) energie_quadratique_mode
425            close(76)
426
427        end subroutine
428
429
430! *************************************************************************
431        subroutine out_champs(vox,voy,voz,rho,fx,fy,fz)
432! *************************************************************************
433!       Assure la sortie des champs dans des fichiers
434! **************************************************************************
435
436            implicit none
437
438            double precision, dimension(0:dimx_tab-1,0:dimy_tab-1,0:dimz_tab-1), intent(in) :: vox,voy,voz,rho
439            double precision, dimension(0:dimx_tab-1,0:dimy_tab-1,0:dimz_tab-1), intent(in) :: fx,fy,fz
440            character(len=50) :: nomfic
441            character(len=13) :: sufx
442
443            call time2string(it*dt + begin,sufx)
444
445            ! on ecrit la vitesse, la densite et la vorticite dans un fichier
446
447            ! on écrit un certain nombre de metadonnees pour aider a identifier le fichier :
448            ! c'est utile pour eviter bien des erreurs lors des entrees-sorties et faciliter le post-traitement...
449            ! Cela permet notamment de s'assurer qu'on utilise le bon codage des données binaires (little-
450            ! endian ou big-endian) et que les tableaux ont bien les dimensions requises.
451
452            ! On écrit ensuite les sorties tableau par tableau dans le fichier final.
453
454            ! Attention, Fortran ajoute au début et à la fin de chaque enregistrement, à chaque
455            ! enregistrement binaire effectué avec l'instruction "write(,form='unformatted')", un entier
456            ! correspondant a la taille de l'enregitrement en octets. Il faut faire attention
457            ! et en tenir compte lors du traitement ulterieur éventuel des donnees (dans Matlab par ex)
458
459            ! Le format final du fichier de sortie est donc (i et j sont des entiers correspondant
460            ! à la taille de l'enregistrement qui respectivement suit et précède) :
461
462            ! i,1,nxmax,nymax,nzmax,lx,ly,lz,truncate,type_trunc,rtrunc_x,rtrunc_y,rtrunc_z,nu,stratification,
463            ! xns,schm,omega2,perturbe,lin,j,i,vox,j,i,voy,j,i,voz,j,i,rho,j,i,fx,j,i,fy,j,i,fz,j,
464            ! i,vbx,vby,vbz,j,i,wbx,wby,wbz,j
465
466            nomfic='velo_rho_vort.'//sufx
467           
468            open (76,file=nomfic,form='unformatted')
469            write(76) 1, nxmax, nymax, nzmax, lx, ly, lz, dt, truncate, type_trunc, rtrunc_x, rtrunc_y, rtrunc_z, &
470                      nu, stratification, xns, schm, omega2, perturbe, lin
471            write(76) it*dt + begin
472            write(76) vox(0:nxmax-1,0:nymax-1,0:nzmax-1)
473            write(76) voy(0:nxmax-1,0:nymax-1,0:nzmax-1)
474            write(76) voz(0:nxmax-1,0:nymax-1,0:nzmax-1)
475            write(76) rho(0:nxmax-1,0:nymax-1,0:nzmax-1)
476            write(76) fx(0:nxmax-1,0:nymax-1,0:nzmax-1)
477            write(76) fy(0:nxmax-1,0:nymax-1,0:nzmax-1)
478            write(76) fz(0:nxmax-1,0:nymax-1,0:nzmax-1)
479            write(76) vbx,vby,vbz
480            write(76) wbx,wby,wbz
481
482            close(76)
483
484        end subroutine
485
486! ***********************************************************
487        subroutine read_run(vox,voy,voz,rho,fx,fy,fz,vbx,vby,vbz,wbx,wby,wbz,begin)
488! ***********************************************************
489!       Lecture des champs d'un run precedent dans le fichier  nom_fichier
490! ***********************************************************
491
492            implicit none
493
494            double precision, dimension(0:dimx_tab-1,0:dimy_tab-1,0:dimz_tab-1),intent(out) :: vox,voy,voz,rho,fx,fy,fz
495            double precision, dimension(0:nxmax-1,0:nymax-1), intent(out) :: vbx,vby,vbz,wbx,wby,wbz
496            double precision, intent(out) :: begin
497            integer :: flag, nxread, nyread, nzread
498            double precision :: lxread, lyread, lzread
499            double precision :: dbl_poub
500            integer :: int_poub
501            logical :: log_poub
502
503            open (unit=8,file='velocity.init',form='unformatted')
504
505            ! on lit les metadonnees
506            read(8) flag,nxread,nyread,nzread,lxread,lyread,lzread,dbl_poub,log_poub,int_poub,dbl_poub,dbl_poub,dbl_poub, &
507                    dbl_poub,log_poub,dbl_poub,dbl_poub,dbl_poub,log_poub,log_poub
508            read(8) begin
509
510            ! on verifie que le format du fichier d'entree est correct. On ne verifie que les dimensions du probleme car ce sont
511            ! les seuls parametres qui ne peuvent vraiment pas etre changes entre deux simulations
512            if (flag /= 1 .or. nxread /= nxmax .or. nyread /= nymax .or. nzread /= nzmax .or. lxread /= lx .or. lyread /= ly) then
513                write(*,*) "ERREUR: le fichier de reprise de run n'a pas le bon format"
514                write(*,*) "ARRET DU PROGRAMME..."
515                stop
516            end if
517
518            read(8) vox(0:nxmax-1,0:nymax-1,0:nzmax-1)
519            read(8) voy(0:nxmax-1,0:nymax-1,0:nzmax-1)
520            read(8) voz(0:nxmax-1,0:nymax-1,0:nzmax-1)
521            read(8) rho(0:nxmax-1,0:nymax-1,0:nzmax-1)
522            read(8) fx(0:nxmax-1,0:nymax-1,0:nzmax-1)
523            read(8) fy(0:nxmax-1,0:nymax-1,0:nzmax-1)
524            read(8) fz(0:nxmax-1,0:nymax-1,0:nzmax-1)
525            read(8) vbx,vby,vbz
526            read(8) wbx,wby,wbz
527
528            close(8)
529       
530        end subroutine
531
532
533! on choisit quelle fonction de generation de l'etat de base 2D doit etre utilisee
534#if choix_base2D==1
535! ***********************************************************
536        subroutine gen_base2D(vbx,vby,vbz,wbx,wby,wbz,vox2,voy2,voz2)
537! ***********************************************************
538!       Generation de l'ecoulement de base
539!       Ecoulement de Stuart
540! ***********************************************************
541
542            implicit none
543
544            double precision, dimension(0:nxmax-1,0:nymax-1), intent(out) :: vbx,vby,vbz,wbx,wby,wbz
545            double precision, dimension(0:dimx_tab-1,0:dimy_tab-1,0:dimz_tab-1), intent(out) :: vox2, voy2, voz2
546            integer :: ix, iy
547            double precision ::rho_stuart
548
549            vbx = 0.
550            vby = 0.
551            vbz = 0.
552            wbx = 0.
553            wby = 0.
554            wbz = 0.
555            vox2 = 0.
556            voy2 = 0.
557            voz2 = 0.
558
559            ! on définit ici le parametre de concentration des tourbillons de stuart
560            rho_stuart = 0.25
561
562            do iy=0, nymax-1
563                do ix = 0,nxmax-1
564                    vbx(ix,iy) = sinh(2*(yy(iy)-ly/2))/(cosh(2* (yy(iy)-ly/2) ) - rho_stuart*cos(  2*xx(ix)  )  )
565                    vby(ix,iy) = - rho_stuart * sin( 2*xx(ix) )/( cosh( 2*( yy(iy)-ly/2 ) ) - rho_stuart*cos(2*xx(ix)))
566                    wbz(ix,iy) = 2*(rho_stuart**2-1)/(cosh(2*(yy(iy)-ly/2)) - rho_stuart*cos(2*xx(ix)))**2
567                end do
568            end do
569
570        end subroutine
571
572
573#elif choix_base2D==2
574! ***********************************************************
575        subroutine gen_base2D(vbx,vby,vbz,wbx,wby,wbz,vox2,voy2,voz2)
576! ***********************************************************
577!       Generation de l'ecoulement de base
578!       Couche de mélange avec un profil vbx=tanh(y)
579! ***********************************************************
580
581            implicit none
582
583            double precision, dimension(0:nxmax-1,0:nymax-1), intent(out) :: vbx,vby,vbz,wbx,wby,wbz
584            double precision, dimension(0:dimx_tab-1,0:dimy_tab-1,0:dimz_tab-1), intent(out) :: vox2, voy2, voz2
585            integer :: ix
586
587            vbx = 0.
588            vby = 0.
589            vbz = 0.
590            wbx = 0.
591            wby = 0.
592            wbz = 0.
593            vox2 = 0.
594            voy2 = 0.
595            voz2 = 0.
596
597            do ix = 0,nxmax-1
598                vbx(ix,:) = tanh(yy-ly/2)
599                wbz(ix,:) = -(1-tanh(yy-ly/2)**2)
600            end do
601
602        end subroutine
603
604
605#elif choix_base2D==3
606! ***********************************************************
607        subroutine gen_base2D(vbx,vby,vbz,wbx,wby,wbz,vox2,voy2,voz2)
608! ***********************************************************
609!       Generation de l'ecoulement de base
610!       Superposition de tourbillons dont les caracteristiques
611!       sont donnees dans le fichier de parametres
612! ***********************************************************
613
614            implicit none
615
616            double precision, dimension(0:nxmax-1,0:nymax-1), intent(out) :: vbx,vby,vbz,wbx,wby,wbz
617            double precision, dimension(0:dimx_tab-1,0:dimy_tab-1,0:dimz_tab-1), intent(out) :: vox2, voy2, voz2
618            double precision k2
619            integer num_tourbillon
620            integer ix, iy, iz
621
622            vbx = 0.
623            vby = 0.
624            vbz = 0.
625            wbx = 0.
626            wby = 0.
627            wbz = 0.
628            voz2 = 0.
629
630            ! on determine la vorticite verticale de l'ecoulement
631            do num_tourbillon=1,nb_tourbillons
632                do iz=0,nzmax-1
633                    do iy=0,nymax-1
634                        do ix=0,nxmax-1
635                            voz2(ix,iy,iz)=voz2(ix,iy,iz)+wz_gau((xx(ix)-xt(num_tourbillon)),(yy(iy)-yt(num_tourbillon)), &
636                            gamma(num_tourbillon),at(num_tourbillon))
637                        end do
638                    end do
639                end do
640            end do
641
642            ! initialisation de la vorticite verticale
643            wbz = voz2(0:nxmax-1,0:nymax-1,0)
644
645            !  inversion du laplacien 2D pour trouver la fct de courant 2D*1D
646            call fwd_fft(voz2,fft_work)
647
648            do iz=0,nzmax-1
649                do iy=0,nymax-1
650                    do ix=0,2*kxmax+1
651                        k2=max(ca2x(ix)+k2y(iy),epsilo)
652                        voz2(ix,iy,iz)= voz2(ix,iy,iz)/k2
653                    end do
654                end do
655            end do
656
657            !  vitesses en fct de la fct de courant 2D*1D
658            do iz=0,nzmax-1
659                do iy=0,nymax-1
660                    do ix=0,kxmax
661                        vox2(2*ix,iy,iz) = -ky(iy)*voz2(2*ix+1,iy,iz)
662                        vox2(2*ix+1,iy,iz) = ky(iy)*voz2(2*ix,iy,iz)
663                        voy2(2*ix,iy,iz) = kx(ix)*voz2(2*ix+1,iy,iz)
664                        voy2(2*ix+1,iy,iz) = -kx(ix)*voz2(2*ix,iy,iz)
665                    end do
666                end do
667            end do
668
669            !  passage dans l'espace physique
670            call bck_fft(vox2,fft_work)
671            call bck_fft(voy2,fft_work)
672
673            vbx = vox2(0:nxmax-1,0:nymax-1,0)
674            vby = voy2(0:nxmax-1,0:nymax-1,0)
675
676        end subroutine
677#endif
678
679
680! ********************************************************************************
681        subroutine read_base2D(vbx,vby,vbz,wbx,wby,wbz)
682! ********************************************************************************
683!       Subroutine qui permet de lire un champ de base en formatted ou unformatted
684! ********************************************************************************
685
686            implicit none
687
688            double precision, dimension(0:nxmax-1,0:nymax-1), intent(out) :: vbx, vby, vbz, wbx, wby, wbz
689            integer :: flag, nxread, nyread
690            double precision :: lxread, lyread
691
692            open (unit=88,file='base2D.init',form='unformatted')
693            read(88) flag, nxread, nyread, lxread, lyread
694            read(88) vbx, vby, vbz
695            read(88) wbx, wby, wbz
696            close(88)
697
698            ! on verifie que le format du fichier d'entree est correct
699            if (flag /= 1 .or. nxread /= nxmax .or. nyread /= nymax .or. lxread /= lx .or. lyread /= ly) then
700                write(*,*) "ERREUR: le fichier de reprise de run n'a pas le bon format"
701                write(*,*) "ARRET DU PROGRAMME..."
702                stop
703            end if
704
705    end subroutine
706
707
708
709! ********************************************************************************
710        subroutine read_velo2D(vox, voy, voz)
711! ********************************************************************************
712!       Subroutine qui permet de lire un champ 2D en unformatted
713! ********************************************************************************
714
715            implicit none
716
717            double precision, dimension(0:dimx_tab-1,0:dimy_tab-1,0:dimz_tab-1), intent(out) :: vox, voy, voz
718            double precision, dimension(0:nxmax-1,0:nymax-1) :: vi2Dx, vi2Dy, vi2Dz
719            integer :: flag, nxread, nyread
720            double precision :: lxread, lyread
721            integer :: ix,iy,iz
722
723            open (unit=88,file='velo2D.init',form='unformatted')
724            read(88) flag, nxread, nyread, lxread, lyread
725            read(88) vi2Dx,vi2Dy,vi2Dz
726            close(88)
727
728            ! on verifie que le format du fichier d'entree est correct
729            if (flag /= 1 .or. nxread /= nxmax .or. nyread /= nymax .or. lxread /= lx .or. lyread /= ly) then
730                write(*,*) "ERREUR: le fichier de reprise de run n'a pas le bon format"
731                write(*,*) "ARRET DU PROGRAMME..."
732                stop
733            end if
734
735            do iz=0,nzmax-1
736                do iy=0,nymax-1
737                    do ix=0,nxmax-1
738                        vox(ix,iy,iz)=vi2Dx(ix,iy)
739                        voy(ix,iy,iz)=vi2Dy(ix,iy)
740                        voz(ix,iy,iz)=vi2Dz(ix,iy)
741                    end do
742                end do
743            end do
744
745        end subroutine
746
747
748! ********************************************************************************
749        subroutine read_dRho(dRho_dz) 
750!ajout d'une sous routine de lecture d'un gradient de densite moyen non constant
751! ********************************************************************************
752!       Subroutine qui permet de lire un champ 2D en unformatted
753! ********************************************************************************
754
755            implicit none
756
757            double precision, dimension(0:dimx_tab-1,0:dimy_tab-1,0:dimz_tab-1), intent(out) :: dRho_dz
758            double precision, dimension(0:nymax-1) :: dRho_dzi
759            integer :: flag, nxread, nyread
760            double precision :: lxread, lyread
761            integer :: ix,iy,iz
762            write(*,*) 'boucle lecture rho'
763            open (unit=88,file='dRho_dz.init',status='old')
764            !read(88) flag, nxread, nyread, lxread, lyread
765            do iz=0,nzmax-1
766            read(88,*) dRho_dzi(iz)
767            enddo
768           
769            close(88)
770
771            ! on verifie que le format du fichier d'entree est correct
772            ! if (flag /= 1 .or. nyread /= nymax .or. lyread /= ly) then
773            !    write(*,*) "ERREUR: le fichier de reprise de run n'a pas le bon format"
774            !    write(*,*) "ARRET DU PROGRAMME..."
775            !    stop
776            !end if
777           
778            do iz=0,nzmax-1
779                do iy=0,nymax-1
780                    do ix=0,nxmax-1
781                        dRho_dz(ix,iy,iz)=dRho_dzi(iz)
782                    end do
783                end do
784            end do
785
786        end subroutine
787
788
789
790
791! ************************************************************************************
792       subroutine gen_velo2D(vox, voy, voz)
793! ************************************************************************************
794!       Generation d'un champ de vitesse correspondant a un dipole fait de
795!       tourbillons de Lamb-oseen.
796! ************************************************************************************
797
798            implicit none
799
800            double precision, dimension(0:dimx_tab-1,0:dimy_tab-1,0:dimz_tab-1), intent(out) :: vox,voy,voz
801            double precision k2
802            integer num_tourbillon
803            integer ix, iy, iz
804
805            vox = 0.
806            voy = 0.
807            voz = 0.
808
809            ! on determine la vorticite verticale de l'ecoulement
810            do num_tourbillon=1,nb_tourbillons
811                do iz=0,nzmax-1
812                    do iy=0,nymax-1
813                        do ix=0,nxmax-1
814                            voz(ix,iy,iz)=voz(ix,iy,iz)+wz_gau((xx(ix)-xt(num_tourbillon)),(yy(iy)-yt(num_tourbillon)), &
815                            gamma(num_tourbillon),at(num_tourbillon))
816                        end do
817                    end do
818                end do
819            end do
820
821            !  inversion du laplacien 2D pour trouver la fct de courant 2D*1D
822            call fwd_fft(voz,fft_work)
823
824            do iz=0,nzmax-1
825                do iy=0,nymax-1
826                    do ix=0,2*kxmax+1
827                        k2=max(ca2x(ix)+k2y(iy),epsilo)
828                        voz(ix,iy,iz)= voz(ix,iy,iz)/k2
829                    end do
830                end do
831            end do
832
833            !  vitesses en fct de la fct de courant 2D*1D
834            do iz=0,nzmax-1
835                do iy=0,nymax-1
836                    do ix=0,kxmax
837                        vox(2*ix,iy,iz) = -ky(iy)*voz(2*ix+1,iy,iz)
838                        vox(2*ix+1,iy,iz) = ky(iy)*voz(2*ix,iy,iz)
839                        voy(2*ix,iy,iz) = kx(ix)*voz(2*ix+1,iy,iz)
840                        voy(2*ix+1,iy,iz) = -kx(ix)*voz(2*ix,iy,iz)
841                    end do
842                end do
843            end do
844
845            !  passage dans l'espace physique
846            call bck_fft(vox,fft_work)
847            call bck_fft(voy,fft_work)
848            voz = 0
849
850        end subroutine
851
852! ***************************************************************
853        function wz_gau(x,y,gamma,a)
854! ***************************************************************
855
856            implicit none
857
858            double precision, intent(in) :: x, y, gamma,a
859            double precision :: wz_gau
860            double precision :: r2
861
862            r2 = x**2 + y**2
863            wz_gau =gamma/(pi*a**2)*exp(-r2/a**2)
864
865        end function
866
867
868! *************************************************************************
869        subroutine white_noise(vox,voy,voz)
870! *********************************************************
871
872            implicit none
873
874            double precision, dimension(0:dimx_tab-1,0:dimy_tab-1,0:dimz_tab-1), intent(inout) :: vox, voy, voz
875            double precision :: rnd_number, moy_x, moy_y, moy_z
876            integer :: ix, iy, iz
877
878            ! on initialise le generateur de nombre aleatoire avec une nouvelle "graine"
879            ! ne semble pas fontionner sur toutes les architectures (bizarre, bizarre...)
880
881! on desactive l'initialisation du generateur de nombre aleatoire en mode debugage de façon
882! a pouvoir comparer plus facilement les resultats de runs differents.
883#ifndef debug
884            call random_seed
885#endif
886
887            ! on determine la moyenne initiale des champs vox, voy, voz
888            moy_x = sum(vox(0:nxmax-1,0:nymax-1,0:nzmax-1))/(nxmax*nymax*nzmax)
889            moy_y = sum(voy(0:nxmax-1,0:nymax-1,0:nzmax-1))/(nxmax*nymax*nzmax)
890            moy_z = sum(voz(0:nxmax-1,0:nymax-1,0:nzmax-1))/(nxmax*nymax*nzmax)
891
892            ! on ajoute le bruit blanc dans l'espace physique
893            do iz=0,nzmax-1
894                do iy=0,nymax-1
895                    do ix=0,nxmax-1
896                        call random_number(rnd_number)
897                        vox(ix,iy,iz)=vox(ix,iy,iz)+amplitude_noise*rnd_number
898                        call random_number(rnd_number)
899                        voy(ix,iy,iz)=voy(ix,iy,iz)+amplitude_noise*rnd_number
900                        call random_number(rnd_number)
901                        voz(ix,iy,iz)=voz(ix,iy,iz)+amplitude_noise*rnd_number
902                    end do
903                end do
904            end do
905
906            ! on fait en sorte que la moyenne du bruit que l'on a ajouté soit nulle
907            ! (on utilise pour cela la difference des moyennes des champs vo avant et après ajout du bruit...)
908            moy_x = moy_x -  sum(vox(0:nxmax-1,0:nymax-1,0:nzmax-1))/(nxmax*nymax*nzmax)
909            moy_y = moy_y -  sum(voy(0:nxmax-1,0:nymax-1,0:nzmax-1))/(nxmax*nymax*nzmax)
910            moy_z = moy_z -  sum(voz(0:nxmax-1,0:nymax-1,0:nzmax-1))/(nxmax*nymax*nzmax)
911
912            vox = vox + moy_x
913            voy = voy + moy_y
914            voz = voz + moy_z
915
916        end subroutine
917
918! *************************************************************************
919        subroutine time2string(time,string_time)
920! *************************************************************************
921!       CONVERSION TIME TO STRING t=0000.000
922! *************************************************************************
923
924            implicit none
925
926            double precision, intent(in) :: time
927            character(len=13), intent(out) :: string_time
928
929            write(string_time,'(A2,I7.7,A1,I3.3)') 't=',int(time),'.',nint(1000*(time-int(time)))
930      end subroutine
931
932
933!**********************************************************************************
934!**********************************************************************************
935!**********    FIN DES ENTREE SORTIES                             *****************
936!**********************************************************************************
937!**********************************************************************************
938
939! *********************************************************************************
940        subroutine adams_bashforth(vox,voy,voz,rho,fx,fy,fz,frho,vox2,voy2,voz2,rho2)
941! *********************************************************************************
942
943        !  coeur du programme avance temporelle sur itmax
944        !  Time discretisation: second order Adams-Bashforth:
945        !  step t+dt: input fields are : V(k,t) and F(k,t-dt)
946        !             output fields are: V(k,t+dt) and F(k,t)
947
948            implicit none
949
950            double precision, dimension(0:dimx_tab-1,0:dimy_tab-1,0:dimz_tab-1), intent(inout) :: vox, voy, voz, rho
951            double precision, dimension(0:dimx_tab-1,0:dimy_tab-1,0:dimz_tab-1), intent(inout) :: fx, fy, fz, frho
952            double precision, dimension(0:dimx_tab-1,0:dimy_tab-1,0:dimz_tab-1), intent(inout) :: vox2, voy2, voz2, rho2
953
954            call cpu_time(tps_temp)
955            tps_copie = tps_copie - tps_temp
956            vox2 = vox
957            voy2 = voy
958            voz2 = voz
959            if (stratification) then
960                rho2 = rho
961            end if
962            call cpu_time(tps_temp)
963            tps_copie = tps_copie + tps_temp
964
965
966            ! les champs contiennent vo = V(t)  /  vo2 = V(t)  /   f = f(V(t-dt))
967
968            ! ************
969            ! Etape 1 / V(t*,k)=E*V(t,k)-(dt/2)*E**2*F(t-dt,k)
970            !           avec E=exp(-nu*k2*dt) des constantes precalculees
971            ! ************
972            call cpu_time(tps_temp)
973            tps_adams_step1 = tps_adams_step1 - tps_temp
974            call adams_step1(vox,voy,voz,rho,fx,fy,fz,frho)
975            call cpu_time(tps_temp)
976            tps_adams_step1 = tps_adams_step1 + tps_temp
977
978            ! les champs contiennent vo = V(t*)  /  vo2 = V(t)  /   f = f(V(t-dt))
979
980            ! ************
981            ! Etape 2 / V(t+dt) = V(t*)+(3.*dt/2.)*E*F(V(t))
982            !           avec E=exp(-nu*k2*dt) des constantes precalculees
983            ! ************
984
985            call nonli(vox2,voy2,voz2,rho2,fx,fy,fz,frho,dRho_dz)
986
987            call cpu_time(tps_temp)
988            tps_adams_step2 = tps_adams_step2 - tps_temp
989            call adams_step2(vox,voy,voz,rho,fx,fy,fz,frho)
990            call cpu_time(tps_temp)
991            tps_adams_step2 = tps_adams_step2 + tps_temp
992
993            ! les champs contiennent vo = V(t+dt)  /  vo2 = ##  /   f = f(V(t))
994
995        end subroutine
996
997
998! *************************************************************************
999        subroutine adams_step1(vox,voy,voz,rho,fx,fy,fz,frho)
1000! *************************************************************************
1001!       Première étape du schéma d'ordre deux d'Adams-Bashforth
1002! *************************************************************************
1003
1004
1005            implicit none
1006
1007            double precision, dimension(0:dimx_tab-1,0:dimy_tab-1,0:dimz_tab-1), intent(inout) :: vox, voy, voz, rho
1008            double precision, dimension(0:dimx_tab-1,0:dimy_tab-1,0:dimz_tab-1), intent(in) :: fx, fy, fz, frho
1009            double precision :: ee,ees,e2,e2s
1010            double precision :: temp_dt
1011            integer :: ix, iy, iz
1012
1013            temp_dt=dt/2.
1014
1015
1016            do iz=0,nzmax-1
1017                do iy=0,nymax-1
1018                    do ix=0,2*kxmax+1
1019#ifdef test_troncature
1020                        if (trunctab(ix,iy,iz)) then
1021#endif
1022                            ee=eex(ix)*eey(iy)*eez(iz)
1023                            ees=exs(ix)*eeys(iy)*eezs(iz)
1024                            e2=ee*ee
1025                            e2s=ees*ees
1026                            vox(ix,iy,iz)=ee*vox(ix,iy,iz)-temp_dt*e2*fx(ix,iy,iz)
1027                            voy(ix,iy,iz)=ee*voy(ix,iy,iz)-temp_dt*e2*fy(ix,iy,iz)
1028                            voz(ix,iy,iz)=ee*voz(ix,iy,iz)-temp_dt*e2*fz(ix,iy,iz)
1029                            rho(ix,iy,iz)=ees*rho(ix,iy,iz)-temp_dt*e2s*frho(ix,iy,iz)
1030#ifdef test_troncature
1031                        end if
1032#endif
1033                    end do
1034                end do
1035            end do
1036
1037
1038        end subroutine
1039
1040
1041! *************************************************************************
1042        subroutine adams_step2(vox,voy,voz,rho,fx,fy,fz,frho)
1043! *************************************************************************
1044!       Deuxième étape du schéma d'ordre deux d'Adams-Bashforth
1045! *************************************************************************
1046
1047            implicit none
1048
1049            double precision, dimension(0:dimx_tab-1,0:dimy_tab-1,0:dimz_tab-1), intent(inout) :: vox, voy, voz, rho
1050            double precision, dimension(0:dimx_tab-1,0:dimy_tab-1,0:dimz_tab-1), intent(in) :: fx, fy, fz, frho
1051            double precision :: ee,ees
1052            double precision :: temp_dt
1053            integer :: ix, iy, iz
1054
1055            temp_dt=3.*dt/2.
1056
1057
1058            do iz=0,nzmax-1
1059                do iy=0,nymax-1
1060                    do ix=0,2*kxmax+1
1061#ifdef test_troncature
1062                        if (trunctab(ix,iy,iz)) then
1063#endif
1064                            ee=eex(ix)*eey(iy)*eez(iz)
1065                            ees=exs(ix)*eeys(iy)*eezs(iz)
1066                            vox(ix,iy,iz)=vox(ix,iy,iz)+ee*temp_dt*fx(ix,iy,iz)
1067                            voy(ix,iy,iz)=voy(ix,iy,iz)+ee*temp_dt*fy(ix,iy,iz)
1068                            voz(ix,iy,iz)=voz(ix,iy,iz)+ee*temp_dt*fz(ix,iy,iz)
1069                            rho(ix,iy,iz)=rho(ix,iy,iz)+ees*temp_dt*frho(ix,iy,iz)
1070#ifdef test_troncature
1071                        end if
1072#endif
1073                    end do
1074                end do
1075            end do
1076
1077
1078        end subroutine
1079
1080
1081
1082! *************************************************************************
1083        subroutine runge_step2(vox,voy,voz,rho,fx,fy,fz,frho)
1084! *************************************************************************
1085!       schema de type leap-frog centre a t+dt/2.
1086!       Utilise comme etape du subroutine  raccro
1087!          qui effectue un pas Runge-Kutta a l'ordre deux pour redemarrer.
1088!          ( v(n+1/2)+1/2 - v(n+1/2)-1/2 )/(2.dt/2) = f( v(n+1/2) ) .
1089! **************************************************************************
1090
1091            implicit none
1092
1093            double precision, dimension(0:dimx_tab-1,0:dimy_tab-1,0:dimz_tab-1), intent(inout) :: vox, voy, voz, rho
1094            double precision, dimension(0:dimx_tab-1,0:dimy_tab-1,0:dimz_tab-1), intent(in) :: fx, fy, fz, frho
1095            double precision :: k2,ee,a,ees,as
1096            integer :: ix, iy, iz
1097
1098            do iz=0,nzmax-1
1099                do iy=0,nymax-1
1100                    do ix=0,2*kxmax+1
1101                        k2=max(ca2x(ix)+k2y(iy)+k2z(iz),epsilo)
1102                        ee=eex(ix)*eey(iy)*eez(iz)
1103                        a=(1.-ee)/(nu*k2)
1104                        ees=exs(ix)*eeys(iy)*eezs(iz)
1105                        as=(1.-ees)/(nu/schm*k2)
1106
1107                        vox(ix,iy,iz)=vox(ix,iy,iz)*ee+fx(ix,iy,iz)*a
1108                        voy(ix,iy,iz)=voy(ix,iy,iz)*ee+fy(ix,iy,iz)*a
1109                        voz(ix,iy,iz)=voz(ix,iy,iz)*ee+fz(ix,iy,iz)*a
1110                        rho(ix,iy,iz)=rho(ix,iy,iz)*ees+frho(ix,iy,iz)*as
1111                    end do
1112                end do
1113            end do
1114
1115        end subroutine
1116
1117
1118! *************************************************************************
1119        subroutine runge_step1(vox,voy,voz,rho,fx,fy,fz,frho)
1120! *************************************************************************
1121!       schema de type Euler.
1122!       Utilisé comme étape du subroutine  raccro
1123!       qui effectue un pas Runge-Kutta a l'ordre deux pour redemarrer.
1124!       ( v(n+1/2)-v(n))/(dt/2)  =  f( v(n) )
1125! *************************************************************************
1126
1127            implicit none
1128
1129            double precision, dimension(0:dimx_tab-1,0:dimy_tab-1,0:dimz_tab-1), intent(inout) :: vox, voy, voz, rho
1130            double precision, dimension(0:dimx_tab-1,0:dimy_tab-1,0:dimz_tab-1), intent(in) :: fx, fy, fz, frho
1131            double precision :: k2,ee,a,ees,as
1132            integer :: ix, iy, iz
1133
1134            do iz=0,nzmax-1
1135                do iy=0,nymax-1
1136                    do ix=0,2*kxmax+1
1137                        k2=max(ca2x(ix)+k2y(iy)+k2z(iz),epsilo)
1138                        ee=eex(ix)*eey(iy)*eez(iz)
1139                        a=(1.-ee)/(nu*k2)
1140                        ees=exs(ix)*eeys(iy)*eezs(iz)
1141                        as=(1.-ees)/(nu/schm*k2)
1142
1143                        vox(ix,iy,iz)=vox(ix,iy,iz)*ee+fx(ix,iy,iz)*a
1144                        voy(ix,iy,iz)=voy(ix,iy,iz)*ee+fy(ix,iy,iz)*a
1145                        voz(ix,iy,iz)=voz(ix,iy,iz)*ee+fz(ix,iy,iz)*a
1146                        rho(ix,iy,iz)=rho(ix,iy,iz)*ees+frho(ix,iy,iz)*as
1147                    end do
1148                end do
1149            end do
1150
1151        end subroutine
1152
1153! *************************************************************************
1154        subroutine rotv(ax,ay,az,bx,by,bz)
1155! *************************************************************************
1156!       calcul du rotationnel du champ (ax,ay,az).
1157! *************************************************************************
1158
1159            implicit none
1160
1161            double precision, dimension(2,0:dimx_tab/2-1,0:dimy_tab-1,0:dimz_tab-1), intent(in) :: ax, ay, az
1162            double precision, dimension(2,0:dimx_tab/2-1,0:dimy_tab-1,0:dimz_tab-1), intent(out) :: bx, by, bz
1163            integer :: ix, iy, iz
1164
1165
1166            do iz=0,nzmax-1
1167                do iy=0,nymax-1
1168                    do ix=0,kxmax
1169#ifndef mode_kz
1170                        bx(1,ix,iy,iz)=kz(iz)*ay(2,ix,iy,iz)-ky(iy)*az(2,ix,iy,iz)
1171                        bx(2,ix,iy,iz)=ky(iy)*az(1,ix,iy,iz)-kz(iz)*ay(1,ix,iy,iz)
1172                        by(1,ix,iy,iz)=kx(ix)*az(2,ix,iy,iz)-kz(iz)*ax(2,ix,iy,iz)
1173                        by(2,ix,iy,iz)=kz(iz)*ax(1,ix,iy,iz)-kx(ix)*az(1,ix,iy,iz)
1174#else
1175                        bx(1,ix,iy,iz)=-kz(iz)*ay(1,ix,iy,iz)-ky(iy)*az(2,ix,iy,iz)
1176                        bx(2,ix,iy,iz)=ky(iy)*az(1,ix,iy,iz)-kz(iz)*ay(2,ix,iy,iz)
1177                        by(1,ix,iy,iz)=kx(ix)*az(2,ix,iy,iz)+kz(iz)*ax(1,ix,iy,iz)
1178                        by(2,ix,iy,iz)=kz(iz)*ax(2,ix,iy,iz)-kx(ix)*az(1,ix,iy,iz)
1179#endif
1180
1181                        bz(1,ix,iy,iz)=ky(iy)*ax(2,ix,iy,iz)-kx(ix)*ay(2,ix,iy,iz)
1182                        bz(2,ix,iy,iz)=kx(ix)*ay(1,ix,iy,iz)-ky(iy)*ax(1,ix,iy,iz)
1183                     end do
1184                end do
1185            end do
1186
1187        end subroutine
1188
1189
1190! *************************************************************************
1191        subroutine trunc0(ax,ay,az,b)
1192! *************************************************************************
1193!       classical desaliasing in the spectral space.
1194! **************************************************************************
1195
1196            implicit none
1197
1198            double precision, dimension(0:dimx_tab-1,0:dimy_tab-1,0:dimz_tab-1), intent(inout) ::  ax, ay, az, b
1199            integer :: ix, iy, iz
1200
1201            do iz=0,nzmax-1
1202                do iy=0,nymax-1
1203                    do ix=0,2*kxmax+1
1204                        if(.not. trunctab(ix,iy,iz)) then
1205                            ax(ix,iy,iz) = 0
1206                            ay(ix,iy,iz) = 0
1207                            az(ix,iy,iz) = 0
1208                            b(ix,iy,iz) = 0
1209                        end if
1210                    end do
1211                end do
1212            end do
1213
1214        end subroutine
1215
1216! *************************************************************************
1217        subroutine init_trunctab(trunctab)
1218! *************************************************************************
1219!       Génération du tableau utilisé lors de la troncature
1220! **************************************************************************
1221
1222            implicit none
1223
1224            logical, dimension(0:dimx_tab-1,0:dimy_tab-1,0:dimz_tab-1), intent(out) :: trunctab
1225            double precision :: kx2max,ky2max,kz2max,dks
1226            integer :: ix, iy, iz
1227
1228            ! rayon du cercle: kmax*rtrunc
1229            kx2max=(rtrunc_x*nxmax/2*2*pi/lx)**2
1230            ky2max=(rtrunc_y*nymax/2*2*pi/ly)**2
1231
1232#ifndef mode_kz
1233            kz2max=(rtrunc_z*nzmax/2*2*pi/lz)**2
1234#else
1235            kz2max=rtrunc_z**2
1236#endif
1237
1238            trunctab = .false.
1239            do iz=0,nzmax-1
1240                do iy=0,nymax-1
1241                    do ix=0,2*kxmax+1
1242
1243                        if (type_trunc==1) then
1244                            dks = max(ca2x(ix)/kx2max, k2y(iy)/ky2max, k2z(iz)/kz2max)
1245                        else
1246                            dks = ca2x(ix)/kx2max + k2y(iy)/ky2max + k2z(iz)/kz2max
1247                        end if
1248
1249                        ! on rajoute un epsilo pour eviter que numériquement, un terme soit tronqué alors qu'il
1250                        ! ne devrait pas l'être.
1251                        if (dks <= 1. + epsilo) then
1252                            trunctab(ix,iy,iz)=.true.
1253                        end if
1254                    end do
1255                end do
1256            end do
1257
1258        end subroutine
1259
1260
1261! *************************************************************************
1262        subroutine vecpro2(ax,ay,az,bx,by,bz)
1263! *************************************************************************
1264!       produits vectoriels dans l'espace physique pour le code perturbatif.
1265!       inclut les termes avec l'ecoulement de base et coriolis
1266! **************************************************************************
1267
1268            implicit none
1269
1270            double precision, dimension(0:dimx_tab-1,0:dimy_tab-1,0:dimz_tab-1), intent(inout) :: bx, by, bz
1271            double precision, dimension(0:dimx_tab-1,0:dimy_tab-1,0:dimz_tab-1), intent(in) :: ax, ay, az
1272            double precision :: vbx0,vby0,vbz0,wbx0,wby0,wbz0
1273            double precision :: ax0,ay0,az0,bx0,by0,bz0
1274            integer :: ix, iy, iz
1275
1276            ! VERSION LINEARISEE
1277            if (lin) then
1278
1279                do iz=0,nzmax-1
1280                    do iy=0,nymax-1
1281                        do ix=0,nxmax-1
1282
1283                            vbx0=vbx(ix,iy)
1284                            vby0=vby(ix,iy)
1285                            vbz0=vbz(ix,iy)
1286                            wbx0=wbx(ix,iy)
1287                            wby0=wby(ix,iy)
1288                            wbz0=wbz(ix,iy)+omega2
1289
1290                            ax0=ax(ix,iy,iz)
1291                            ay0=ay(ix,iy,iz)
1292                            az0=az(ix,iy,iz)
1293                            bx0=bx(ix,iy,iz)
1294                            by0=by(ix,iy,iz)
1295                            bz0=bz(ix,iy,iz)
1296
1297                            bx(ix,iy,iz)=vby0*bz0-vbz0*by0 + ay0*wbz0-az0*wby0
1298                            by(ix,iy,iz)=vbz0*bx0-vbx0*bz0 + az0*wbx0-ax0*wbz0
1299                            bz(ix,iy,iz)=vbx0*by0-vby0*bx0 + ax0*wby0-ay0*wbx0
1300                        end do
1301                    end do
1302                end do
1303
1304            ! NON LINEAIRE
1305            else
1306                ! calcul en nonlineaire perturbatif :
1307                if (perturbe) then
1308
1309                    do iz=0,nzmax-1
1310                        do iy=0,nymax-1
1311                            do ix=0,nxmax-1
1312                                vbx0=vbx(ix,iy)
1313                                vby0=vby(ix,iy)
1314                                vbz0=vbz(ix,iy)
1315                                wbx0=wbx(ix,iy)
1316                                wby0=wby(ix,iy)
1317                                wbz0=wbz(ix,iy)+omega2
1318
1319                                ax0=ax(ix,iy,iz)
1320                                ay0=ay(ix,iy,iz)
1321                                az0=az(ix,iy,iz)
1322                                bx0=bx(ix,iy,iz)
1323                                by0=by(ix,iy,iz)
1324                                bz0=bz(ix,iy,iz)
1325
1326                                bx(ix,iy,iz)=vby0*bz0-vbz0*by0 + ay0*(wbz0+bz0) - az0*(wby0+by0)
1327                                by(ix,iy,iz)=vbz0*bx0-vbx0*bz0 + az0*(wbx0+bx0) - ax0*(wbz0+bz0)
1328                                bz(ix,iy,iz)=vbx0*by0-vby0*bx0 + ax0*(wby0+by0) - ay0*(wbx0+bx0)
1329                            end do
1330                        end do
1331                    end do
1332
1333                !  calcul en non perturbatif :
1334                else
1335
1336                    do iz=0,nzmax-1
1337                        do iy=0,nymax-1
1338                            do ix=0,nxmax-1
1339                                ax0=ax(ix,iy,iz)
1340                                ay0=ay(ix,iy,iz)
1341                                az0=az(ix,iy,iz)
1342                                bx0=bx(ix,iy,iz)
1343                                by0=by(ix,iy,iz)
1344                                bz0=bz(ix,iy,iz)
1345
1346                                bx(ix,iy,iz)=ay0*(omega2+bz0)-az0*by0
1347                                by(ix,iy,iz)=az0*bx0-ax0*(omega2+bz0)
1348                                bz(ix,iy,iz)=ax0*by0-ay0*bx0
1349                            end do
1350                        end do
1351                    end do
1352
1353                end if
1354            end if
1355
1356        end subroutine
1357
1358! pba: deb
1359! *************************************************************************
1360        subroutine forcage(bx,by,bz)
1361! *************************************************************************
1362!       ajout d'un terme de forcage
1363!       parametres du forcage entres ici provisoirement modifier
1364!       data.in etc pour entree plus propre
1365! **************************************************************************
1366
1367
1368            implicit none
1369
1370            double precision, dimension(0:dimx_tab-1,0:dimy_tab-1,0:dimz_tab-1), intent(inout) :: bx, by, bz
1371            double precision :: vbx0,vby0,vbz0,wbx0,wby0,wbz0
1372            double precision :: ax0,ay0,az0,bx0,by0,bz0
1373            double precision :: time, arg, sigma1, sigma2, FF, GG, acc
1374            integer :: ix, iy, iz, itime, i0, j0, k0    !itime a faire passer en parametre?
1375                ! itime=it0
1376
1377                ! itime=itime+1   ! penser a it0=it0+itmax(run precedent si restat
1378                ! print *,'temps:', itime
1379                !    A_ics=100.             !deplacement   a definir dans  data.in et parametres.F90
1380                !freq=2*pi/(3600*12.4)  !M2 a verifier idem
1381                arg=freq*it*Dt      !Dt est defini dans data.in ok
1382               ! print *,'arg:', arg
1383               ! print *,'A_ics:', A_ics  !acceleration
1384             
1385
1386            ! localisation centre de la source
1387
1388              i0=floor(0.5*nxmax)
1389              j0=floor(0.5*nymax)
1390              k0=floor(0.5*nzmax)
1391             
1392            ! enveloppe gaussienne
1393
1394!sigma1=4km, sigma2=100m     revoir les valeurs ci-dessous
1395              sigma1=20.   ! defini en unite de dx et dz
1396              sigma2=1.    !                    dy
1397
1398
1399                    do iz=0,nzmax-1
1400                        do iy=0,nymax-1
1401                            do ix=0,nxmax-1
1402
1403                                FF=exp(-((ix-i0)**2+(iy-j0)**2)/sigma1**2)
1404                                GG=exp(-(iz-k0)**2/sigma2**2)
1405
1406                                bx0=bx(ix,iy,iz)
1407                                by0=by(ix,iy,iz)
1408                                bz0=bz(ix,iy,iz)
1409
1410                                bx(ix,iy,iz)=bx0+A_ics*freq**2*FF*GG*sin(arg)
1411                                !  Acc=A_ics*FF*GG*sin(arg)
1412                                !print *,'acceleration:', Acc
1413                                by(ix,iy,iz)=by0
1414                                bz(ix,iy,iz)=bz0
1415
1416                            end do
1417                        end do
1418                    end do
1419
1420
1421
1422        end subroutine
1423
1424!pba : fin
1425
1426! *************************************************************************
1427        subroutine urho(ax,ay,az,b)
1428! *************************************************************************
1429!       term u.b pour l'equation de la densite (calcule dans l'espace physique)
1430! **************************************************************************
1431            implicit none
1432            double precision, dimension(0:dimx_tab-1,0:dimy_tab-1,0:dimz_tab-1), intent(inout) :: ax, ay, az, b
1433            integer :: ix, iy, iz 
1434
1435            ! VERSION LINEARISEE
1436            if (lin) then
1437                do iz=0,nzmax-1
1438                    do iy=0,nymax-1
1439                        do ix=0,nxmax-1
1440                            ax(ix,iy,iz)=b(ix,iy,iz)*vbx(ix,iy)
1441                            ay(ix,iy,iz)=b(ix,iy,iz)*vby(ix,iy)
1442                            az(ix,iy,iz)=b(ix,iy,iz)*vbz(ix,iy)
1443                        end do
1444                    end do
1445                end do
1446            ! NON LINEAIRE
1447            else
1448                ! calcul en nonlineaire perturbatif :
1449                if (perturbe) then
1450                    do iz=0,nzmax-1
1451                         do iy=0,nymax-1
1452                            do ix=0,nxmax-1
1453                                ax(ix,iy,iz)=b(ix,iy,iz)*(ax(ix,iy,iz)+vbx(ix,iy))
1454                                ay(ix,iy,iz)=b(ix,iy,iz)*(ay(ix,iy,iz)+vby(ix,iy))
1455                                az(ix,iy,iz)=b(ix,iy,iz)*(az(ix,iy,iz)+vbz(ix,iy))
1456                            end do
1457                        end do
1458                    end do
1459                !  calcul en non perturbatif :
1460                else
1461                    do iz=0,nzmax-1
1462                        do iy=0,nymax-1
1463                            do ix=0,nxmax-1
1464                                ax(ix,iy,iz)=b(ix,iy,iz)*ax(ix,iy,iz)
1465                                ay(ix,iy,iz)=b(ix,iy,iz)*ay(ix,iy,iz)
1466                                az(ix,iy,iz)=b(ix,iy,iz)*az(ix,iy,iz)
1467                            end do
1468                        end do
1469                    end do
1470                end if
1471            end if
1472
1473         end subroutine
1474
1475! *************************************************************************
1476        subroutine vozdRho_dz(a,b,c)
1477! *************************************************************************
1478!       term c=a.b pour l'equation de la densite (avec a=voz  b=dRho_dz) (calcule dans l'espace physique)
1479! **************************************************************************
1480            implicit none
1481            double precision, dimension(0:dimx_tab-1,0:dimy_tab-1,0:dimz_tab-1), intent(in) :: a, b
1482            double precision, dimension(0:dimx_tab-1,0:dimy_tab-1,0:dimz_tab-1), intent(out) :: c 
1483            integer :: ix, iy, iz 
1484
1485           
1486                    do iz=0,nzmax-1
1487                        do iy=0,nymax-1
1488                            do ix=0,nxmax-1
1489                                c(ix,iy,iz)=b(ix,iy,iz)*a(ix,iy,iz)
1490                               
1491                            end do
1492                        end do
1493                    end do
1494           
1495
1496         end subroutine
1497
1498
1499
1500
1501! *************************************************************************
1502        subroutine div_speciale(ax,ay,az,bz1,bz2)
1503! *************************************************************************
1504!       Cette sous-routine renvoie dans le champ bz2 le terme bz1 - div(ax,ay,az)
1505!       utile pour evaluer le terme non lineaire sur l'equation de la densite en economisant un champ de stockage...
1506! **************************************************************************
1507       
1508            implicit none
1509
1510            double precision ,dimension(2,0:dimx_tab/2-1,0:dimy_tab-1,0:dimz_tab-1), intent(in) :: ax, ay, az, bz1
1511            double precision ,dimension(2,0:dimx_tab/2-1,0:dimy_tab-1,0:dimz_tab-1), intent(inout) :: bz2
1512            integer :: ix, iy, iz
1513           
1514 
1515                       
1516
1517            do iz=0,nzmax-1
1518                do iy=0,nymax-1
1519                    do ix=0,kxmax
1520                       bz2(1,ix,iy,iz) = - bz1(1,ix,iy,iz) + (kx(ix)*ax(2,ix,iy,iz)+ky(iy)*ay(2,ix,iy,iz)+kz(iz)*az(2,ix,iy,iz))
1521                       bz2(2,ix,iy,iz) = - bz1(2,ix,iy,iz) - (kx(ix)*ax(1,ix,iy,iz)+ky(iy)*ay(1,ix,iy,iz)+kz(iz)*az(1,ix,iy,iz))
1522                    end do
1523                end do
1524            end do
1525   
1526        end subroutine
1527
1528! **************************************************************************
1529        subroutine init_k(kx,k2x,ky,k2y,kz,k2z,cax,ca2x)
1530! **************************************************************************
1531!       Définition des vecteurs d'ondes
1532! **************************************************************************
1533
1534            implicit none
1535
1536            double precision, dimension(0:kxmax), intent(out) :: kx, k2x
1537            double precision, dimension(0:nymax-1), intent(out) :: ky, k2y
1538            double precision, dimension(0:nzmax-1), intent(out) :: kz, k2z
1539            double precision, dimension(0:2*kxmax+1), intent(out) :: cax, ca2x
1540
1541#ifndef mode_kz
1542            integer :: ix, iy, iz
1543#else
1544            integer :: ix, iy
1545#endif
1546
1547            do ix=0,kxmax
1548                kx(ix)=ix*2.*pi/lx
1549            end do
1550
1551            do iy=0,nymax/2
1552                ky(iy)=iy*2.*pi/ly
1553            end do
1554
1555            do iy=nymax/2+1,nymax-1
1556                ky(iy)=(iy-nymax)*2.*pi/ly
1557            end do
1558
1559#ifndef mode_kz
1560            do iz=0,nzmax/2
1561                kz(iz)=iz*2.*pi/lz
1562            end do
1563
1564            do iz=nzmax/2+1,nzmax-1
1565                kz(iz)=(iz-nzmax)*2.*pi/lz
1566            end do
1567#else
1568            if (lz==0) then
1569                kz = 0
1570            else
1571                kz = 2.*pi/lz
1572            end if
1573#endif
1574
1575            k2x=kx*kx
1576            k2y=ky*ky
1577            k2z=kz*kz
1578
1579            ! définition alternative des vecteurs d'onde suivant x
1580            ! utile pour un certains nombre de sous-routines (améliore la vectorisation sur le NEC SX-5 de l'Idris...)
1581            do ix=0,kxmax
1582                cax(2*ix)=ix*2.*pi/lx
1583                cax(2*ix+1)=ix*2.*pi/lx
1584            end do
1585
1586            ca2x=cax*cax
1587
1588        end subroutine
1589
1590! *************************************************************************
1591        subroutine nonli(vox,voy,voz,rho,fx,fy,fz,frho,dRho_dz)
1592! *************************************************************************
1593!       calcul des termes non linéaires f(t) et frho(t) en spectral
1594!       en entrée v(t) et rho(t) en spectral
1595           
1596            implicit none
1597
1598            double precision, dimension(0:dimx_tab-1,0:dimy_tab-1,0:dimz_tab-1), intent(in) :: dRho_dz
1599            double precision, dimension(0:dimx_tab-1,0:dimy_tab-1,0:dimz_tab-1), intent(inout) :: vox, voy, voz, rho
1600            double precision, dimension(0:dimx_tab-1,0:dimy_tab-1,0:dimz_tab-1), intent(out) :: fx, fy, fz, frho
1601           
1602! necessaire au calcul d'une vitesse de translation
1603#ifdef calcul_ut
1604            double precision, dimension(0:dimx_tab-1,0:dimy_tab-1,0:dimz_tab-1) :: dwz_dx, dwz_dy
1605            integer :: ix, iy, iz
1606            double precision :: u0, erreur
1607#endif
1608
1609            ! sauvegarde de voz en sphérique dans le champ frho (ruse pour gagner de la place en memoire)
1610            !(cette ruse ne marche plus si on rentre un dRho_dz)
1611            ! nécessaire pour le cas stratifié
1612            !if (stratification) then
1613            !    call cpu_time(tps_temp)
1614            !    tps_copie = tps_copie - tps_temp
1615            !    frho=voz
1616            !    call cpu_time(tps_temp)
1617            !    tps_copie = tps_copie + tps_temp
1618            !end if
1619
1620            ! calcul de la vorticite de vo. Le resultat est placé dans f.
1621            call cpu_time(tps_temp)
1622            tps_rotv = tps_rotv - tps_temp
1623            call rotv(vox,voy,voz,fx,fy,fz)
1624            call cpu_time(tps_temp)
1625            tps_rotv = tps_rotv + tps_temp
1626
1627#ifdef calcul_ut
1628            ! on determine la derivee de la vorticite verticale par rapport a x et par rapport a y
1629            do iz=0, nzmax-1
1630                do iy=0, nymax-1
1631                    do ix=0, kxmax
1632                        dwz_dx(2*ix,iy,iz) = -kx(ix)*fz(2*ix+1,iy,iz)
1633                        dwz_dx(2*ix+1,iy,iz) = kx(ix)*fz(2*ix,iy,iz)
1634
1635                        dwz_dy(2*ix,iy,iz) = -ky(iy)*fz(2*ix+1,iy,iz)
1636                        dwz_dy(2*ix+1,iy,iz) = ky(iy)*fz(2*ix,iy,iz)
1637                    end do
1638                end do
1639            end do
1640#endif
1641
1642            ! on appelle ici éventuellement la fonction de sortie dans l'espace spectral out_spectral1
1643            call cpu_time(tps_temp)
1644            tps_out = tps_out - tps_temp
1645
1646            if (periode_out_spectral1 /= 0 .and. mod(it, periode_out_spectral1) == 0 .and. it /= 0) then
1647                call out_spectral1(vox, voy, voz)
1648            end if
1649
1650            call cpu_time(tps_temp)
1651            tps_out = tps_out + tps_temp
1652
1653            ! passage dans l'espace physique de la vitesse, de la vorticite et de la densite
1654            call cpu_time(tps_temp)
1655
1656            tps_bck_fft = tps_bck_fft - tps_temp
1657           
1658           
1659            call bck_fft(vox,fft_work)
1660            call bck_fft(voy,fft_work)
1661            call bck_fft(voz,fft_work)
1662            call bck_fft(fx,fft_work)
1663            call bck_fft(fy,fft_work)
1664            call bck_fft(fz,fft_work)
1665
1666            if (stratification) then
1667                call bck_fft(rho,fft_work)
1668            end if
1669
1670#ifdef calcul_ut
1671            ! on passe dans l'espace physique les termes dwz_dx et dwz_dy
1672            call bck_fft(dwz_dx,fft_work)
1673            call bck_fft(dwz_dy,fft_work)
1674#endif
1675
1676            call cpu_time(tps_temp)
1677            tps_bck_fft = tps_bck_fft + tps_temp
1678
1679            ! on appelle eventuellement les sorties periodiques dans l'espace physique
1680            ! il y a deux sorties periodiques dans l'espace physique qui peuvent avoir des periodes
1681            ! différentes
1682            call cpu_time(tps_temp)
1683            tps_out = tps_out - tps_temp
1684
1685            if (periode_out_physique1 /= 0 .and. mod(it, periode_out_physique1) == 0 .and. it /= 0) then
1686                call out_physique1(vox,voy,voz)
1687            end if
1688
1689            if (periode_out_physique2 /= 0 .and. mod(it, periode_out_physique2) == 0 .and. it /= 0) then
1690                call out_physique2(vox,voy,voz,rho,fx,fy,fz)
1691            end if
1692
1693            call cpu_time(tps_temp)
1694            tps_out = tps_out + tps_temp
1695
1696
1697            ! calcul de V x W dans l'espace physique
1698            call cpu_time(tps_temp)
1699            tps_vecpro2 = tps_vecpro2 - tps_temp
1700            call vecpro2(vox,voy,voz,fx,fy,fz)
1701!pba : deb
1702            call forcage(fx,fy,fz)
1703!pba : fin
1704            call cpu_time(tps_temp)
1705            tps_vecpro2 = tps_vecpro2 + tps_temp
1706
1707#ifdef calcul_ut
1708            ! on calcule une estimation de u0
1709            u0 = sum((vox(0:nxmax-1,0:nymax-1,0:nzmax-1)*dwz_dx(0:nxmax-1,0:nymax-1,0:nzmax-1)+voy(0:nxmax-1,0:nymax-1, &
1710                      0:nzmax-1)*dwz_dy(0:nxmax-1,0:nymax-1,0:nzmax-1))*dwz_dx(0:nxmax-1,0:nymax-1,0:nzmax-1))
1711
1712            u0 = u0/max(sum(dwz_dx(0:nxmax-1,0:nymax-1,0:nzmax-1)**2), epsilo)
1713
1714            ! on estime l'erreur que l'on commet
1715            erreur = sum(((vox(0:nxmax-1,0:nymax-1,0:nzmax-1)-u0)*dwz_dx(0:nxmax-1,0:nymax-1,0:nzmax-1) + &
1716                          voy(0:nxmax-1,0:nymax-1,0:nzmax-1)*dwz_dy(0:nxmax-1,0:nymax-1,0:nzmax-1))**2)/(nxmax*nymax*nzmax)
1717
1718            ! on imprime les resultats
1719            write(*,'(A21,TR2,G23.10)') '   Estimation de ut =', u0
1720            write(*,'(A34,TR2,G23.10)') '   Estimation de l erreur sur ut =', erreur
1721
1722#endif
1723
1724            ! calcul du terme non linéaire pour l'équation de la densité
1725            if (stratification) then
1726
1727                ! on calcul u*rho dans l'espace physique
1728                call cpu_time(tps_temp)
1729                tps_urho = tps_urho - tps_temp
1730               
1731                           
1732                call vozdRho_dz(voz,dRho_dz,wdRho)
1733                call urho(vox,voy,voz,rho)
1734                 
1735                call cpu_time(tps_temp)
1736                tps_urho = tps_urho + tps_temp
1737
1738                ! on passe u*rho et w*dRho_dz en spectral
1739                call cpu_time(tps_temp)
1740                tps_fwd_fft = tps_fwd_fft - tps_temp
1741
1742                call fwd_fft(vox,fft_work)
1743                call fwd_fft(voy,fft_work)
1744                call fwd_fft(voz,fft_work)
1745                call fwd_fft(wdRho,fft_work)
1746
1747                call cpu_time(tps_temp)
1748                tps_fwd_fft = tps_fwd_fft + tps_temp
1749
1750                ! la fonction div_speciale renvoie frho = voz*dRho_dz - u*grad(rho) (drho_dz= gradient de densité moyen rentré  au départ)
1751                !
1752                ! u*grad(rho) est calcule en prenant la divergence de u*rho parce que div(u) = 0
1753                ! de plus on a ruse en stockant voz dans le champ frho pour gagner un champ en memoire
1754                call cpu_time(tps_temp)
1755                tps_div_speciale = tps_div_speciale - tps_temp
1756                call div_speciale(vox,voy,voz,wdRho,frho)
1757                call cpu_time(tps_temp)
1758                tps_div_speciale = tps_div_speciale + tps_temp
1759
1760                ! ajout du terme de flottaison pour l'équation de la vitesse
1761                !if (Rho_read) !on a lu un profil de densite initial
1762                fz=fz-9.81*rho/Rhob 
1763                !else
1764                !fz=fz-xnsqr*rho  !si pas de profil de densité initial entré N est rescalé de manière à avoir dRho/dz=1         
1765                 
1766            end if
1767
1768            ! Retour dans l'espace spectral de f
1769            call cpu_time(tps_temp)
1770            tps_fwd_fft = tps_fwd_fft - tps_temp
1771
1772            call fwd_fft(fx,fft_work)
1773            call fwd_fft(fy,fft_work)
1774            call fwd_fft(fz,fft_work)
1775
1776            call cpu_time(tps_temp)
1777            tps_fwd_fft = tps_fwd_fft + tps_temp
1778
1779        end subroutine
1780
1781
1782! **************************************************************************
1783        subroutine projection(ax,ay,az)
1784! **************************************************************************
1785!  calcul de la projection : Pij( V x W )(k)  dans l'espace spectral.
1786! **************************************************************************
1787
1788            implicit none
1789
1790            double precision, dimension(2,0:dimx_tab/2-1,0:dimy_tab-1,0:dimz_tab-1), intent(inout) :: ax, ay, az
1791            double precision inv_k2, s1, s2
1792            integer :: ix, iy, iz
1793
1794
1795            do iz=0,nzmax-1
1796                do iy=0,nymax-1
1797                    do ix=0,kxmax
1798#ifdef test_troncature
1799                        if (trunctab(ix,iy,iz)) then
1800#endif
1801                            ! On utilise epsilo pour éviter que le module de k ne soit nul
1802                            inv_k2=1./max(k2x(ix)+k2y(iy)+k2z(iz),epsilo)
1803
1804#ifndef mode_kz
1805                            s1 =(kx(ix)*ax(1,ix,iy,iz)+ky(iy)*ay(1,ix,iy,iz)+kz(iz)*az(1,ix,iy,iz))*inv_k2
1806                            s2 =(kx(ix)*ax(2,ix,iy,iz)+ky(iy)*ay(2,ix,iy,iz)+kz(iz)*az(2,ix,iy,iz))*inv_k2
1807#else
1808                            s1 =(kx(ix)*ax(1,ix,iy,iz)+ky(iy)*ay(1,ix,iy,iz)-kz(iz)*az(2,ix,iy,iz))*inv_k2
1809                            s2 =(kx(ix)*ax(2,ix,iy,iz)+ky(iy)*ay(2,ix,iy,iz)+kz(iz)*az(1,ix,iy,iz))*inv_k2
1810#endif
1811
1812                            ax(1,ix,iy,iz)= ax(1,ix,iy,iz)-s1*kx(ix)
1813                            ax(2,ix,iy,iz)= ax(2,ix,iy,iz)-s2*kx(ix)
1814                            ay(1,ix,iy,iz)= ay(1,ix,iy,iz)-s1*ky(iy)
1815                            ay(2,ix,iy,iz)= ay(2,ix,iy,iz)-s2*ky(iy)
1816
1817#ifndef mode_kz
1818                            az(1,ix,iy,iz)= az(1,ix,iy,iz)-s1*kz(iz)
1819                            az(2,ix,iy,iz)= az(2,ix,iy,iz)-s2*kz(iz)
1820#else
1821                            az(1,ix,iy,iz)= az(1,ix,iy,iz)-s2*kz(iz)
1822                            az(2,ix,iy,iz)= az(2,ix,iy,iz)+s1*kz(iz)
1823#endif
1824
1825
1826#ifdef test_troncature
1827                        end if
1828#endif
1829                    end do
1830                end do
1831            end do
1832
1833
1834        end subroutine
1835
1836! ***************************************************************************
1837        subroutine precal2(dtl)
1838! ***************************************************************************
1839!       Précalcul des termes exponentiels pour l'intégration en temps
1840! ***************************************************************************
1841
1842            implicit none
1843
1844            double precision, intent(in) :: dtl
1845
1846            eex=exp(-nu*ca2x*dtl)
1847            exs=exp(-nu/schm*ca2x*dtl)
1848
1849            eey=exp(-nu*k2y*dtl)
1850            eeys=exp(-nu/schm*k2y*dtl)
1851
1852            eez=exp(-nu*k2z*dtl)
1853            eezs=exp(-nu/schm*k2z*dtl)
1854
1855        end subroutine
1856
1857! **********************************************************************
1858        subroutine raccro(vox,voy,voz,rho,fx,fy,fz,frho,vox2,voy2,voz2,rho2)
1859! **********************************************************************
1860!        Cette subroutine permet à partir du champs v(t) de déterminer
1861!        avec un schéma Runge-Kutta à l'ordre deux les champs v(t+dt) et f(t)
1862! **********************************************************************
1863
1864            implicit none
1865
1866            double precision, dimension(0:dimx_tab-1,0:dimy_tab-1,0:dimz_tab-1), intent(inout) :: vox, voy, voz, rho
1867            double precision, dimension(0:dimx_tab-1,0:dimy_tab-1,0:dimz_tab-1), intent(out) :: fx, fy, fz, frho
1868            double precision, dimension(0:dimx_tab-1,0:dimy_tab-1,0:dimz_tab-1), intent(inout) :: vox2, voy2, voz2, rho2
1869
1870            vox2=vox
1871            voy2=voy
1872            voz2=voz
1873            if (stratification) then
1874                rho2 = rho
1875            end if
1876
1877            ! les champs contiennent vo = V(t)  /  vo2 = V(t)  /   f = ##
1878
1879            ! ************
1880            ! Etape 1 / Schema Euler avant : (V(t+dt/2)-V(t))/(dt/2) = f(V(t))
1881            ! ************
1882
1883            ! calcul des termes non linéaires f(V(t=0))
1884            call nonli(vox,voy,voz,rho,fx,fy,fz,frho,dRho_dz)
1885
1886            vox=vox2
1887            voy=voy2
1888            voz=voz2
1889            if (stratification) then
1890                rho = rho2
1891            end if
1892
1893            ! les champs contiennent vo = V(t)  /  vo2 = V(t)  /   f = f(V(t))
1894
1895            ! precalcul des termes exponentiels associes au shéma Euler avant
1896            call precal2(dt/2.)
1897
1898            ! shema Euler avant proprement dit
1899            call runge_step1(vox,voy,voz,rho,fx,fy,fz,frho)
1900
1901            ! on projette et on tronque
1902            call projection(vox,voy,voz)
1903            if(truncate) then
1904                call trunc0(vox,voy,voz,rho)
1905            endif
1906
1907            ! les champs contiennent vo = V(t+dt/2)  /  vo2 = V(t)  /   f = f(V(t))
1908
1909            ! ************
1910            ! Etape 2 / Schema leap-frog centre en t+dt/2 : ( V(t+dt)-V(t) )/(dt) = f( V(t+dt/2) )
1911            ! ************
1912
1913            call nonli(vox,voy,voz,rho,fx,fy,fz,frho,dRho_dz)
1914
1915            vox=vox2
1916            voy=voy2
1917            voz=voz2
1918            rho=rho2
1919
1920            ! les champs contiennent vo = V(t)  /  vo2 = V(t)  /   f = f(V(t+dt/2))
1921
1922            ! precalcul des termes exponentiels associes au shéma Euler avant
1923            call precal2(dt)
1924
1925            ! shema leap-frog proprement dit
1926            call runge_step2(vox,voy,voz,rho,fx,fy,fz,frho)
1927
1928            ! les champs contiennent vo = V(t+dt)  /  vo2 = V(t)  /   f = f(V(t+dt/2))
1929
1930            ! ************
1931            ! Etape 3 / Recalcul du terme f(V(t)) qui sera utile par la suite pour le schéma d'Adams-Bashforth
1932            ! ************
1933            call nonli(vox2,voy2,voz2,rho2,fx,fy,fz,frho,dRho_dz)
1934
1935            ! les champs contiennent vo = V(t+dt)  /  vo2 = V(t)  /   f = f(V(t))
1936
1937        end subroutine
1938
1939! ***************************************************************************
1940        subroutine init_coord(xx,yy,zz)
1941! ***************************************************************************
1942!       Initialisations de tableaux de longueur
1943! ***************************************************************************
1944
1945            implicit none
1946
1947            double precision, dimension(0:nxmax-1), intent(out) :: xx
1948            double precision, dimension(0:nymax-1), intent(out) :: yy
1949            double precision, dimension(0:nzmax-1), intent(out) :: zz
1950            integer :: i
1951
1952            do i=0,nxmax-1
1953                xx(i)=i*lx/nxmax
1954            end do
1955
1956            do i=0,nymax-1
1957                yy(i)=i*ly/nymax
1958            end do
1959
1960            do i=0,nzmax-1
1961                zz(i)=i*lz/nzmax
1962            end do
1963
1964        end subroutine
1965
1966end module
Note: See TracBrowser for help on using the repository browser.