source: branches/GRISLIv3/SOURCES/out_cptr_mod.f90 @ 447

Last change on this file since 447 was 447, checked in by aquiquet, 9 months ago

Cleaning branch: use only statements in out_cptr and strain_rate.

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