source: branches/iLoveclim/SOURCES/out_cptr_mod.f90

Last change on this file was 287, checked in by aquiquet, 5 years ago

Trunk merged to iLoveclim branch at revision 286

File size: 19.5 KB
Line 
1!> \file out_cptr_mod.f90
2!! Module avec les routines d'ecriture et de lecture de fichiers cptr
3!<
4
5!> \namespace out_cptr
6!! This module gathers routines to read and write cptr files
7!! \author Vincent Peyaud
8!! \date janvier 2006
9!!  @note  attention :  en ecriture fichier runname.cptr
10!!  @note               en lecture  fichier runname.CPTR
11!!  @note               il faut renommer le fichier !
12!<
13!----------------------------------------------------------
14
15
16!> subroutine: read_recovery
17
18!!Reprise d'un fichier compteur
19!! Used module:  - use module3D_phy
20!!               - use netcdf
21!!               - use io_netcdf_grisli
22!!               - use tracer_vars
23!!
24!
25!>
26
27subroutine read_recovery(icpt)           !---- reprise d'un fichier compteur -----
28
29  use module3D_phy
30  use netcdf
31  use io_netcdf_grisli
32  use tracer_vars                  !aurel neem
33  implicit none
34
35  character(len=20)                  :: titre
36  integer                            :: nzzm
37  integer                            :: icpt  ! icompteur 1->  tout (topo, T, U)
38                                              !           2 -> tout sauf topo
39                                              !           3 -> tout sauf topo et eau basale
40 
41
42 real,dimension(nx,ny)              :: bidon         !< sert à sauter les lectures de Ss,H,B,...
43                             
44  ! Variables specifiques pour les fichier .nc
45 
46  real*8, dimension(:,:),   pointer  :: tab => null()         !< tableau 2d real ecrit dans le fichier
47  real*8, dimension(:,:,:), pointer  :: tab1  => null()       !< tableau 3d real
48  real*8, dimension(:,:,:), pointer  :: tab1T => null()       !< tableau 3d real pour la temperature
49  real*8, dimension(:),     pointer  :: liste_time => null()
50
51! rappel de f90
52! reprcptr est un string (ici le nom du fichier recovery)
53! index (reprcptr,'.cptr')    renvoie :
54!                 0, si 'cptr' n'est pas dans reprcptr
55!                 la position de la premiere occurence de cptr sinon
56
57  if (ncdf_type.eq.0) call lect_netcdf_type         !< pour lire la valeur de netcdf_type (machine dependant)
58
59
60ascii_nc:  if ((index(reprcptr,'.cptr').ne.0) .or.(index(reprcptr,'.CPTR').ne.0))  then   !reprise fichier cptr
61 
62     if (itracebug.eq.1)  call tracebug(' Entree dans routine read_recovery : ascii')
63
64     open(num_forc,file=trim(reprcptr))
65     read(num_forc,*) time
66
67
68! temps de la reprise
69!-----------------------
70     if (tbegin.gt.1.d9) then       !    si tbegin est tres grand , on reprend le temps du cptr
71        tbegin=time
72     else
73        time=tbegin
74     endif
75
76
77     read(num_forc,*) titre
78     read(num_forc,*) ni,nzz,nzzm,nxx,nyy
79     read(num_forc,*)
80
81
82     if ((nxx.ne.nx).or.(nyy.ne.ny).or.(nzz.ne.nz)   &
83          .or.(nzzm.ne.nzm)) write(6,*) 'attention', &
84          'les tailles de tableaux ne sont pas compatibles'
85
86     !     lecture des tableaux 2D
87     read(num_forc,*)
88     read(num_forc,*)
89
90lect_topo:     if (icpt.eq.1) then                     ! lecture topo seulement si icompteur = 1 
91        read(num_forc,*) S
92        read(num_forc,*)
93
94        read(num_forc,*) H
95        read(num_forc,*)
96
97        read(num_forc,*) B
98        read(num_forc,*)
99
100        read(num_forc,*) Bsoc
101        read(num_forc,*)
102
103     else
104        read(num_forc,*) bidon
105        read(num_forc,*)
106
107        read(num_forc,*) bidon
108        read(num_forc,*)
109
110        read(num_forc,*) bidon
111        read(num_forc,*)
112
113        read(num_forc,*) bidon
114        read(num_forc,*)
115
116     end if lect_topo
117
118                                                       ! lecture temperature dans tous les cas     
119     read(num_forc,*) Ibase
120     read(num_forc,*)
121
122     read(num_forc,*) bmelt
123     read(num_forc,*)
124
125lect_Hwater: if  (icpt.ne.3) then                      ! pas de lecture Hwater si icompteur = 3 
126        read(num_forc,*) Hwater
127        read(num_forc,*)
128     else
129        read(num_forc,*) bidon
130        read(num_forc,*)
131        Hwater(:,:)=0.
132
133     end if lect_Hwater
134
135                                                       ! lecture vitesses dans tous les cas.
136        read(num_forc,*) Uxbar
137        read(num_forc,*)
138
139        read(num_forc,*) Uybar
140
141
142     !     tableaux 3d
143
144     read(num_forc,*)
145     read(num_forc,*)
146
147     read(num_forc,*,end=22) T
14822   continue
149
150     read(num_forc,*,end=23)
151     read(num_forc,*,end=23)
152
153     read(num_forc,*,end=23) Ux
15423   continue
155
156     read(num_forc,*,end=24)
157     read(num_forc,*,end=24)
158
159     read(num_forc,*,end=24) Uy
16024   continue
161
162! aurel pour faire des cptr avec les traceurs :
163      read(num_forc,*,end=25)
164      read(num_forc,*,end=25)
165      read(num_forc,*,end=25) Xdep_out  !xdepk ou xdep ? a voir
16625    continue
167
168      read(num_forc,*,end=26)
169      read(num_forc,*,end=26)
170      read(num_forc,*,end=26) Ydep_out
17126    continue
172
173      read(num_forc,*,end=27)
174      read(num_forc,*,end=27)
175      read(num_forc,*,end=27) tdep_out
17627    continue
177      tdep_out(:,:,:)=tdep_out(:,:,:)+time
178
179     close(num_forc)
180
181  else if ((index(reprcptr,'.nc').ne.0)) then      !     Modif hassine pour la reprise d'un fichier .nc
182
183     if (itracebug.eq.1)  call tracebug(' Entree dans routine read_recovery : netcdf')
184
185
186        if (.not.associated(liste_time)) then
187           allocate(liste_time(1)) 
188        end if
189        if (.not.associated(tab)) then
190           allocate(tab(nx,ny),tab1(nx,ny,nz),tab1T(nx,ny,nz+nzm))
191        end if
192
193        call Read_ncdf_var('temps',TRIM(reprcptr),liste_time)
194
195! temps de la reprise
196!-----------------------
197     if (tbegin.gt.1.d9) then       !    si tbegin est tres grand , on reprend le temps du .nc
198        tbegin=liste_time(1)
199        time = tbegin
200     else
201        time = tbegin
202     endif
203
204     if (itracebug.eq.1)   write(num_tracebug,*) tbegin, liste_time(1)
205
206
207!        TBEGIN=liste_time(1)-liste_time(1)
208       
209        call Read_ncdf_dim('x' ,TRIM(reprcptr),nxx)
210        call Read_ncdf_dim('y' ,TRIM(reprcptr),nyy)
211        call Read_ncdf_dim('z' ,TRIM(reprcptr),nzz)
212        call Read_ncdf_dim('zm',TRIM(reprcptr),nzzm)
213
214        if ((nxx.ne.nx).or.(nyy.ne.ny).or.(nzz.ne.nz)       &
215             .or.(nzzm-nzz.ne.nzm)) write(6,*) 'attention', &
216             'les tailles de tableaux ne sont pas compatibles'
217
218        !     lecture des tableaux 2D
219        ! -----------------------------
220
221        if (icpt.eq.1) then                             ! lecture topo seulement si icompteur = 1
222           call Read_ncdf_var('S',TRIM(reprcptr),tab)
223           S(:,:)=tab(:,:)
224           call Read_ncdf_var('H',TRIM(reprcptr),tab)
225           H(:,:)=tab(:,:)
226           call Read_ncdf_var('B',TRIM(reprcptr),tab)
227           B(:,:)=tab(:,:)
228           call Read_ncdf_var('BSOC',TRIM(reprcptr),tab)
229           Bsoc(:,:)=tab(:,:)
230        end if
231
232                                                       ! lecture temperature dans tous les cas     
233        call Read_ncdf_var('IBASE',TRIM(reprcptr),tab)
234        Ibase(:,:)=tab(:,:)
235        call Read_ncdf_var('BMELT',TRIM(reprcptr),tab)
236        Bmelt(:,:)=tab(:,:)
237
238        if (icpt.ne.3) then                            ! pas de lecture Hwater si icompteur = 3
239           call Read_ncdf_var('HWATER',TRIM(reprcptr),tab)
240           Hwater(:,:)=tab(:,:)
241        else
242           Hwater(:,:)=0.
243        end if
244
245                                                       ! lecture vitesses dans tous les cas.
246        call Read_ncdf_var('UXBAR',TRIM(reprcptr),tab)
247        UXBAR(:,:)=tab(:,:)
248        call Read_ncdf_var('UYBAR',TRIM(reprcptr),tab)
249        UYBAR(:,:)=tab(:,:)
250
251        do i=1,nx-1
252           do j=1,ny-1
253              oldu(i,j) = sqrt(((uxbar(i,j)+uxbar(i+1,j))/2)**2 &
254                   + ((uybar(i,j)+uybar(i,j+1))/2.)**2)
255           end do
256        end do
257
258        !     tableaux 3D
259        ! -----------------------------
260
261        call Read_ncdf_var('T',TRIM(reprcptr),tab1T)
262        T(:,:,:)=tab1T(:,:,:)
263
264        call Read_ncdf_var('UX',TRIM(reprcptr),tab1)
265        Ux(:,:,:)=tab1(:,:,:)
266
267        call Read_ncdf_var('UY',TRIM(reprcptr),tab1)
268        Uy(:,:,:)=tab1(:,:,:)
269     
270        if (itracebug.eq.1)  call tracebug('avant lectures variables traceur')
271
272
273        call Read_ncdf_var('XDEP',TRIM(reprcptr),tab1)  !aurel neem
274        Xdep_out(:,:,:)=tab1(:,:,:)
275
276        call Read_ncdf_var('YDEP',TRIM(reprcptr),tab1)
277        Ydep_out(:,:,:)=tab1(:,:,:)
278
279        call Read_ncdf_var('TDEP',TRIM(reprcptr),tab1)
280        Tdep_out(:,:,:)=tab1(:,:,:)+time
281 
282     end if ascii_nc
283
284
285     do i=1,nx-1
286        do j=1,ny-1
287           oldu(i,j) = sqrt(((Uxbar(i,j)+Uxbar(i+1,j))/2)**2 &
288                + ((Uybar(i,j)+Uybar(i,j+1))/2.)**2)
289        end do
290     end do
291
292
293     !    moyenne de l'epaisseur et pentes
294     do j=2,ny
295        do i=2,nx
296           Hmy(i,j)=(H(i,j)+H(i,j-1))/2.
297           Hmx(i,j)=(H(i,j)+H(i-1,j))/2.
298           Sdx(i,j)=(S(i,j)-S(i-1,j))/dx
299           Sdy(i,j)=(S(i,j)-S(i,j-1))/dy
300        end do
301     end do
302
303     do i=2,nx
304        Hmx(i,1)=(H(i,1)+H(i-1,1))/2.
305     end do
306     do j=2,ny
307        Hmy(1,j)=(H(1,j)+H(1,j-1))/2.
308     end do
309
310     ! vitesse a la base
311     Uxflgz(:,:) = Ux(:,:,nz)
312     Uyflgz(:,:) = Uy(:,:,nz)
313     
314     where (flot(:,:))
315       debug_3D(:,:,122) = T(:,:,nz)+273.15
316     elsewhere
317       debug_3D(:,:,120) = T(:,:,nz)+273.15
318     endwhere
319
320end subroutine read_recovery
321
322
323!> SUBROUTINE: read_no_recovery
324!!Pas de reprise d'un fichier compteur
325!>
326subroutine  read_no_recovery
327  use module3D_phy
328  implicit none
329  Real,Parameter :: Dzm=600        !< Grid Step In Mantle
330  !Prop Thermique
331  Real,Parameter ::  Cm=1.04e8 
332
333!!!pas reprise : il faut initier les temperatures ds le socle
334
335  !     temperature dans le socle lineaire avec le gradient geothermique
336  do k=nz+1,nz+nzm
337     do j=1,ny
338        do i=1,nx
339           T(i,j,k)=T(i,j,nz)-dzm*(k-nz)*ghf(i,j)/cm
340        end do
341     end do
342  end do
343
344
345end subroutine  read_no_recovery
346
347!----------------------------------------------------------
348
349!----------------------------------------------------------
350
351!> SUBROUTINE: out_recovery
352!! Sortie d'un fichier de reprise
353!! Attention : fichier en sortie .cptr ou fichier .nc
354!! used modules:   - use module3D_phy
355!!                 - use netcdf
356!!                 - use io_netcdf_grisli
357!!                 - use util_recovery
358!!                 - use tracer_vars
359!>
360
361subroutine out_recovery(isortie)
362
363  !   ---------------------------------------------------------------
364  !    SORTIE COMPTEUR
365  !   
366  !   
367  !     sortie qui stocke dans le fichier runname.cptr les variables
368  !     S, H, B, HDOT, BDOT, et T  plusieurs fois dans le programme
369  !     (chaque nouvelle sortie ecrase la precedente) =>
370  !     on peut faire repartir le programme a partir de ce fichier,
371  !     et donc du dernier pas de temps stocke
372  !
373  !     Maintenant, il y a deux formats de fichiers possibles :
374  !     ascii (isortie = 1) et netcdf (isortie = 2).
375  !   ---------------------------------------------------------------
376
377  use module3D_phy
378  use netcdf
379  use io_netcdf_grisli
380  use util_recovery
381  use tracer_vars      ! aurel neem
382
383  implicit none
384
385  character(len=80) :: filin
386  integer ncid,status
387  real*8, dimension(:,:), pointer   :: tab => null()   !< tableau 2d real ecrit dans le fichier
388  real*8, dimension(:,:,:), pointer :: tab1 => null()  !< tableau 3d real
389  real*8, dimension(:,:,:), pointer :: tab1T => null() !< tableau 3d real pour la temperature
390  character(len=20),dimension(3)    :: dimnames2d      !< dimensions pour netcdf pour tableau 2d pour les noeud majeur
391  character(len=20),dimension(4)    :: dimnames3d      !< pour 3d troisieme dim est nz
392  character(len=20),dimension(4)    :: dimnames3dT     !< pour 3d troisieme dim est nz+nzm
393!  integer                           :: nrecs=1         !< compteur pour les enregistrements temps des variables 
394!  integer                           :: idef=0          !< pour savoir si la variable a ete definie ou non
395  integer                           :: isortie         !< pour le choix de type de fichier en sortie
396
397  if (itracebug.eq.1)  call tracebug(' Entree dans routine out_recovery')
398
399  call testout_recovery(filin)
400  if (itracebug.eq.1)  call tracebug(' Entree dans routine out_recovery 2')
401
402  if (logic_out) then
403     Select Case (Isortie)
404     Case (1)                                       ! sortie de type ascii 
405
406        filin=trim(filin)//'.cptr'
407        open(num_file1,file=trim(filin))
408
409        write(num_file1,*) time_out(1), '    TIME '
410        write(num_file1,*) geoplace
411        write(num_file1,*) nx*ny,nz,nzm,nx,ny 
412        write(num_file1,*)
413
414        !     ecriture des tableaux 2D
415        write(num_file1,*)' Tableaux 2D'
416        write(num_file1,*)' Surface '
417        write(num_file1,*) S
418
419        write(num_file1,*)' Epaisseur '
420        write(num_file1,*) H
421
422        write(num_file1,*)' Base glace '
423        write(num_file1,*) B
424
425        write(num_file1,*)' Socle '
426        write(num_file1,*) Bsoc
427
428        write(num_file1,*)' Ibase '
429        write(num_file1,*) Ibase
430
431        write(num_file1,*)' Bmelt '
432        write(num_file1,*) Bmelt
433
434        write(num_file1,*)' Hwater '
435        write(num_file1,*) Hwater
436
437        write(num_file1,*)' Uxbar '
438        write(num_file1,*) Uxbar
439
440        write(num_file1,*)' Uybar '
441        write(num_file1,*) Uybar
442
443        !     tableau 3D
444
445        write(num_file1,*)
446        write(num_file1,*) 'Temperature (3D y compris le socle)'
447        write(num_file1,*) T(:,:,:)
448        write(num_file1,*)
449
450        write(num_file1,*) 'vitesse selon x'
451        write(num_file1,*) Ux(:,:,:)
452        write(num_file1,*)
453
454        write(num_file1,*) 'vitesse selon y'
455        write(num_file1,*) Uy(:,:,:)
456
457!aurel neem : traceurs dans le cptr
458
459        write(num_file1,*)
460        write(num_file1,*) 'origine de la glace en x'
461        write(num_file1,*) xdep(:,:,:)
462
463        write(num_file1,*)
464        write(num_file1,*) 'origine de la glace en y'
465        write(num_file1,*) ydep(:,:,:)
466
467        write(num_file1,*)
468        write(num_file1,*) 'age de la glace'
469        !tab1(:,:,:) = tdep(:,:,:) - time
470        write(num_file1,*) tdep(:,:,:) - time           ! note: on decide d'ecrire le temps 'relatif'
471
472        close(num_file1)
473
474
475     Case (2)                                     ! sortie de type netcdf
476
477        if (.not.associated(tab)) then
478           allocate(tab(nx,ny),tab1(nx,ny,nz),tab1T(nx,ny,nz+nzm))
479        end if
480        filin=trim(filin)//'.nc'
481
482! creation du fichier
483!------------------------
484
485        if (ncdf_type.eq.32) then
486           status  = nf90_create(TRIM(filin),NF90_WRITE,ncid)    ! en 32 bits
487        else if (ncdf_type.eq.64) then
488           status  = nf90_create(trim(filin),and(nf90_write,nf90_64bit_offset),ncid)   ! sur r2d2
489        else
490           write(6,*)'pb de lecture de netcdf_type dans out_cptr:', ncdf_type
491        endif
492
493        status  = nf90_close(ncid)                            ! fermeture
494
495
496        call initfile(nx,ny,nz,nzm,dimnames2d,dimnames3d,dimnames3dT,filin)   ! initialise le fichier nc
497
498        call write_ncdf_var('temps','time',trim(filin),time_out,'double')     ! ecrit le temps du snapshot
499
500!       variables 2D
501
502        tab(:,:)= S(:,:)
503        call write_ncdf_var('S',dimnames2d,trim(filin),tab,'double')          ! ecrit S
504
505        tab(:,:)= H(:,:)
506        call write_ncdf_var('H',dimnames2d,trim(filin),tab,'double')         
507
508        tab(:,:)= B(:,:)
509        call write_ncdf_var('B',dimnames2d,trim(filin),tab,'double')
510
511        tab(:,:)= Bsoc(:,:)
512        call write_ncdf_var('BSOC',dimnames2d,trim(filin),tab,'double')
513
514        tab(:,:)= IBASE(:,:)
515        call write_ncdf_var('IBASE',dimnames2d,trim(filin),tab,'double')
516
517        tab(:,:)= Bmelt(:,:)
518        call write_ncdf_var('BMELT',dimnames2d,trim(filin),tab,'double')
519
520        tab(:,:)= Hwater(:,:)
521        call write_ncdf_var('HWATER',dimnames2d,trim(filin),tab,'double')
522
523        tab(:,:)= UXBAR(:,:)
524        call write_ncdf_var('UXBAR',dimnames2d,trim(filin),tab,'double')
525
526        tab(:,:)= UYBAR(:,:)
527        call write_ncdf_var('UYBAR',dimnames2d,trim(filin),tab,'double')
528
529!       variables 3D
530
531        tab1T(:,:,:)= T(:,:,:)
532        call write_ncdf_var('T',dimnames3dT,trim(filin),tab1T,'double')       ! ecrit T
533
534        tab1(:,:,:)= UX(:,:,:)
535        call write_ncdf_var('UX',dimnames3d,trim(filin),tab1,'double')
536
537        tab1(:,:,:)= UY(:,:,:)
538        call write_ncdf_var('UY',dimnames3d,trim(filin),tab1,'double')
539
540!       variables 3D  pour les traceurs (Aurel NEEM)
541!       attention ce sont les tabelaux _out qui ont les bonnes dimensions.
542
543        tab1(:,:,:)=XDEP_out(:,:,:)
544        call write_ncdf_var('XDEP',dimnames3d,trim(filin),tab1,'double')  !aurel neem
545
546        tab1(:,:,:)=YDEP_out(:,:,:)
547        call write_ncdf_var('YDEP',dimnames3d,trim(filin),tab1,'double')
548
549        tab1(:,:,:)=TDEP_out(:,:,:) - time
550        call write_ncdf_var('TDEP',dimnames3d,trim(filin),tab1,'double')
551
552
553     End Select
554  End if
555end subroutine out_recovery
556
557!> subroutine initfile
558!! Initialise the netcdf file
559!! Used module:   - use netcdf
560!!                - use io_netcdf_grisli
561!! @param nxx                dimension along x
562!! @param nyy                dimension along y
563!! @param nzz                dimension along z
564!! @param nzmm               dimension along z for T
565!! @param fil_sortie         name of the file to initialise
566!! @param dimnames2d         dimensions for netcdf
567!>
568subroutine initfile(nxx,nyy,nzz,nzmm,&
569     dimnames2d,dimnames3d,dimnames3dT,file)
570
571  use netcdf
572  use io_netcdf_grisli
573  implicit none
574  integer,intent(in)             :: nxx               !< dimension along x
575  integer,intent(in)             :: nyy               !< dimension along y
576  integer,intent(in)             :: nzz               !< dimension along z
577  integer,intent(in)             :: nzmm              !< dimension along z in socle
578  character(len=20),dimension(3) :: dimnames2d        !< dimensions pour netcdf pour tableau 2d pour les noeud majeur
579  character(len=20),dimension(4) :: dimnames3d        !< pour 3d troisieme dim est nz
580  character(len=20),dimension(4) :: dimnames3dT       !< pour 3d troisieme dim est nz+nzm
581  character(len=*)               ::file               !< name of the file to init
582  ! initialisation
583  call write_ncdf_dim('x',trim(file),nxx)             !< dimensions des tableaux
584  call write_ncdf_dim('y',trim(file),nyy)
585  call write_ncdf_dim('z',trim(file),nzz)
586  call write_ncdf_dim('zm',trim(file),nzz+nzmm)
587  call write_ncdf_dim('time',trim(file),0)
588  dimnames2d(1)='x'
589  dimnames2d(2)='y'
590  dimnames2d(3)='time' 
591
592  dimnames3d(1)='x'
593  dimnames3d(2)='y'
594  dimnames3d(3)='z'
595  dimnames3d(4)='time' 
596
597  dimnames3dT(1)='x'
598  dimnames3dT(2)='y'
599  dimnames3dT(3)='zm'
600  dimnames3dT(4)='time' 
601
602end subroutine initfile
603
604
605!> SUBROUTINE: symetry_cptr
606!! Subroutine utilisee dans les experiences axysymetriques pour repartir d'un cptr dans lequel
607!! tous les champs sont axysmetriques ou symetriques par rapport a un axe
608!! Used module:      - use module3D_phy
609!! @param    iaxe    symetrie par rapport a iaxe
610!! @param    jaxe    symetrie par rapport a jaxe
611!>
612
613subroutine symetry_cptr(jaxe)
614  use module3D_phy
615  use tracer_vars
616  implicit none
617  integer :: jaxe
618  integer :: jsym
619
620  !symetrique par rapport à jaxe la référence est le bas. Pour l'instant seulement lui
621  jsym=0
622  do j=jaxe+1,ny
623     jsym=jsym+1
624
625     do i=1,nx
626        S(i,j)=S(i,jaxe-jsym)
627        H(i,j)=H(i,jaxe-jsym)
628        B(i,j)=B(i,jaxe-jsym)
629        Bsoc(i,j)=Bsoc(i,jaxe-jsym)
630        Ibase(i,j)=Ibase(i,jaxe-jsym)
631        Bmelt(i,j)=bmelt(i,jaxe-jsym)
632        hwater(i,j)=hwater(i,jaxe-jsym)
633        uxbar(i,j)=uxbar(i,jaxe-jsym)
634        uybar(i,j)=-uybar(i,jaxe-jsym+1)        ! attention different des autres
635        ! a cause des staggered grids
636
637        T(i,j,:)=T(i,jaxe-jsym,:)
638        ux(i,j,:)=ux(i,jaxe-jsym,:)
639        uy(i,j,:)=-uy(i,jaxe-jsym+1,:)          ! different des autres
640
641        XDEP(i,j,:)=XDEP(i,jaxe-jsym,:)
642        YDEP(i,j,:)=YDEP(i,jaxe-jsym,:)
643        TDEP(i,j,:)=TDEP(i,jaxe-jsym,:)
644
645     end do
646  end do
647  return
648end subroutine symetry_cptr
Note: See TracBrowser for help on using the repository browser.