source: trunk/LMDZ.GENERIC/deftank/def_generic_giants/randomRestart/createRandomStart.f90 @ 2468

Last change on this file since 2468 was 1284, checked in by milmd, 10 years ago

LMDZ.GENERIC. added a tool to create start.nc files with perturbated horizontal wind and temperature.

File size: 9.7 KB
Line 
1program createRandomStart
2
3USE netcdf
4
5implicit none
6
7!=======================================================================
8! Copy a start.nc file and change ucov, vcov, teta variables:
9!       teta -> teta +/- dTeta (K)
10!       ucov -> (0 +/- duNat)*cu (m/s)
11!       vcov -> (0 +/- dvNat)*cv (m/s)
12! Perturbations are vertically the same.
13!=======================================================================
14
15!Input
16 character (len=100) :: fileNetcdfInput ! input start.nc netcdf file
17!Output
18 character (len=100) :: fileNetcdf ! output start.nc netcdf file
19!Usefull stuff
20integer :: nLat,nRlatv,nLong,nRlonu,nAlt,nTime ! dimensions  in netcdf file (redundant point in longitude)
21real, dimension (:,:,:,:), allocatable :: uCov,vCov,uNat,vNat,uNat2,vNat2,teta
22real, dimension (:,:), allocatable :: cu,cv
23real :: duNat,dvNat,dTeta
24!Netcdf declarations
25integer :: idFile,idStuff,idUcov,idVcov,idCu,idCv,idTeta
26integer :: ierror
27!Temporary variables
28 character (len=100) :: format,varTmp
29integer :: i,j,k,l,N=1,dimFilter
30logical :: isFile
31real :: alea
32 character (len=200) :: command
33
34
35print*,'input start.nc?'
36read*, fileNetcdfInput
37inquire(file=fileNetcdfInput,exist=isFile)
38if (isFile .eqv. .false.) then
39        stop 'input doesn t exist'
40end if
41
42print*,'output start.nc?'
43read*, fileNetcdf
44inquire(file=fileNetcdf,exist=isFile)
45if (fileNetcdfInput == fileNetcdf) then
46        !name of the input netcdf file back-up
47        varTmp=trim(fileNetcdfInput) // '.old'
48        inquire(file=varTmp,exist=isFile)
49        do while(isFile .eqv. .true.)
50                varTmp=trim(varTmp) // '.old'
51                print*,varTmp
52                inquire(file=varTmp,exist=isFile)
53        end do
54        !input netcdf file back-up
55        write(command,'(a6,a50,a1,a50)') 'cp -f ',trim(fileNetcdfInput),' ',trim(varTmp)
56        call system(command)
57else
58        inquire(file=fileNetcdf,exist=isFile)
59        if (isFile .eqv. .true.) then
60                print*,'output name existing. Want to overwrite? (y/n)'
61                read*,varTmp
62                if (varTmp == 'y') then
63                        write(command,'(a6,a50)') 'rm -f ',fileNetcdf
64                        print*,command
65                        call system(command)
66                else
67                        stop 'existing output not deleted'
68                end if
69        end if
70end if
71
72print*,'du (m/s)?'
73read*, duNat
74
75print*,'dv (m/s)?'
76read*, dvNat
77
78print*,'dT (K)?'
79read*, dTeta
80
81!print*,'dim average filter? (0 --> no | 1,2,3... --> average)'
82!read*, dimFilter
83dimFilter=0
84
85!Creation of the output netcdf file if different from input netcdf file
86if (fileNetcdfInput /= fileNetcdf) then
87        write(command,'(a6,a50,a1,a50)') 'cp -f ',fileNetcdfInput,' ',fileNetcdf
88        print*,command
89        call system(command)
90end if
91
92!**********
93!Netcdf file reading
94PRINT*,'**********'
95PRINT*,'Netcdf file reading'
96!Open
97ierror=nf90_open(fileNetcdf,NF90_WRITE,idFile)
98if(ierror /= nf90_noerr) then
99     print *, trim(nf90_strerror(ierror))
100     stop "Stopped open"
101end if
102!Read dimension latitude
103ierror=nf90_inq_dimid(idFile,'latitude',idStuff)
104if(ierror /= nf90_noerr) then
105     print *, trim(nf90_strerror(ierror))
106     stop "Stopped"
107end if
108ierror=nf90_inquire_dimension(idFile,idStuff,varTmp,nLat)
109if(ierror /= nf90_noerr) then
110     print *, trim(nf90_strerror(ierror))
111     stop "Stopped dim lat"
112end if
113!Read dimension latitude bis
114ierror=nf90_inq_dimid(idFile,'rlatv',idStuff)
115if(ierror /= nf90_noerr) then
116     print *, trim(nf90_strerror(ierror))
117     stop "Stopped"
118end if
119ierror=nf90_inquire_dimension(idFile,idStuff,varTmp,nRlatv)
120if(ierror /= nf90_noerr) then
121     print *, trim(nf90_strerror(ierror))
122     stop "Stopped dim rlatv"
123end if
124!Read dimension longitude
125ierror=nf90_inq_dimid(idFile,'longitude',idStuff)
126if(ierror /= nf90_noerr) then
127     print *, trim(nf90_strerror(ierror))
128     stop "Stopped"
129end if
130ierror=nf90_inquire_dimension(idFile,idStuff,varTmp,nLong)
131if(ierror /= nf90_noerr) then
132     print *, trim(nf90_strerror(ierror))
133     stop "Stopped dim long"
134end if
135!Read dimension longitude bis
136ierror=nf90_inq_dimid(idFile,'rlonu',idStuff)
137if(ierror /= nf90_noerr) then
138     print *, trim(nf90_strerror(ierror))
139     stop "Stopped"
140end if
141ierror=nf90_inquire_dimension(idFile,idStuff,varTmp,nRlonu)
142if(ierror /= nf90_noerr) then
143     print *, trim(nf90_strerror(ierror))
144     stop "Stopped dim rlonu"
145end if
146!Read dimension altitude
147ierror=nf90_inq_dimid(idFile,'altitude',idStuff)
148if(ierror /= nf90_noerr) then
149     print *, trim(nf90_strerror(ierror))
150     stop "Stopped"
151end if
152ierror=nf90_inquire_dimension(idFile,idStuff,varTmp,nAlt)
153if(ierror /= nf90_noerr) then
154     print *, trim(nf90_strerror(ierror))
155     stop "Stopped dim alt"
156end if
157!Read dimension time
158ierror=nf90_inq_dimid(idFile,'Time',idStuff)
159if(ierror /= nf90_noerr) then
160     print *, trim(nf90_strerror(ierror))
161     stop "Stopped"
162end if
163ierror=nf90_inquire_dimension(idFile,idStuff,varTmp,nTime)
164if(ierror /= nf90_noerr) then
165     print *, trim(nf90_strerror(ierror))
166     stop "Stopped dim time"
167end if
168!Read variables ucov and vcov
169allocate(uCov(nRlonu,nLat,nAlt,nTime))
170allocate(vCov(nLong,nRlatv,nAlt,nTime))
171ierror=nf90_inq_varid(idFile,"ucov",idUcov)
172if(ierror /= nf90_noerr) then
173     print *, trim(nf90_strerror(ierror))
174     stop "Stopped"
175end if
176ierror=nf90_get_var(idFile,idUcov,uCov)
177if(ierror /= nf90_noerr) then
178     print *, trim(nf90_strerror(ierror))
179     stop "Stopped get_var ucov"
180end if
181ierror=nf90_inq_varid(idFile,"vcov",idVcov)
182if(ierror /= nf90_noerr) then
183     print *, trim(nf90_strerror(ierror))
184     stop "Stopped"
185end if
186ierror=nf90_get_var(idFile,idVcov,vCov)
187if(ierror /= nf90_noerr) then
188     print *, trim(nf90_strerror(ierror))
189     stop "Stopped get_var vcov"
190end if
191!Read variables cu and cv
192allocate(cu(nRlonu,nLat))
193allocate(cv(nLong,nRlatv))
194ierror=nf90_inq_varid(idFile,"cu",idCu)
195if(ierror /= nf90_noerr) then
196     print *, trim(nf90_strerror(ierror))
197     stop "Stopped"
198end if
199ierror=nf90_get_var(idFile,idCu,cu)
200if(ierror /= nf90_noerr) then
201     print *, trim(nf90_strerror(ierror))
202     stop "Stopped get_var cu"
203end if
204ierror=nf90_inq_varid(idFile,"cv",idCv)
205if(ierror /= nf90_noerr) then
206     print *, trim(nf90_strerror(ierror))
207     stop "Stopped"
208end if
209ierror=nf90_get_var(idFile,idCv,cv)
210if(ierror /= nf90_noerr) then
211     print *, trim(nf90_strerror(ierror))
212     stop "Stopped get_var cv"
213end if
214!Read variables teta
215allocate(teta(nLong,nLat,nAlt,nTime))
216ierror=nf90_inq_varid(idFile,"teta",idTeta)
217if(ierror /= nf90_noerr) then
218     print *, trim(nf90_strerror(ierror))
219     stop "Stopped"
220end if
221ierror=nf90_get_var(idFile,idTeta,teta)
222if(ierror /= nf90_noerr) then
223     print *, trim(nf90_strerror(ierror))
224     stop "Stopped get_var teta"
225end if
226
227!**********
228!uNat
229allocate(uNat(nRlonu,nLat,nAlt,nTime))
230uNat(:,:,:,:)=0.0
231do i=1,nLong
232do j=1,nLat
233do k=1,nAlt
234do l=1,nTime
235        call random_number(alea)
236        uNat(i,j,k,l) = -duNat+2*duNat*alea
237end do
238end do
239end do
240end do
241!-/+180 longitude points must be the same
242uNat(1,:,:,:)=uNat(nLong,:,:,:)
243!**********
244!vNat
245allocate(vNat(nLong,nRlatv,nAlt,nTime))
246vNat(:,:,:,:)=0.0
247do i=1,nLong
248do j=1,nRlatv
249do k=1,nAlt
250do l=1,nTime
251        call random_number(alea)
252        vNat(i,j,k,l) = -dvNat+2*dvNat*alea
253end do
254end do
255end do
256end do
257!-/+180 longitude points must be the same
258vNat(1,:,:,:)=vNat(nLong,:,:,:)
259!**********
260!mean (if dimFilter=0 mean mean_filter has no effect)
261allocate(uNat2(nRlonu,nLat,nAlt,nTime))
262allocate(vNat2(nLong,nRlatv,nAlt,nTime))
263uNat2(:,:,:,:)=uNat(:,:,:,:)
264vNat2(:,:,:,:)=vNat(:,:,:,:)
265do k=1,nAlt
266do l=1,nTime
267        call mean_filter(uNat(:,:,k,l),nRlonu,nLat,dimFilter,uNat2(:,:,k,l))
268        call mean_filter(vNat(:,:,k,l),nLong,nRlatv,dimFilter,vNat2(:,:,k,l))
269end do
270end do
271!**********
272!new ucov and vcov
273do k=1,nAlt
274do l=1,nTime
275        uCov(:,:,k,l)=uNat2(:,:,k,l)*cu
276        vCov(:,:,k,l)=vNat2(:,:,k,l)*cv
277end do
278end do
279!**********
280!teta
281do i=1,nLong
282do j=1,nRlatv
283do k=1,nAlt
284do l=1,nTime
285        call random_number(alea)
286        teta(i,j,k,l) = teta(i,j,k,l)-dTeta+2*dTeta*alea
287end do
288end do
289end do
290end do
291teta(1,:,:,:)=teta(nLong,:,:,:)
292!**********
293!Write variable ucov and vcov
294ierror=nf90_put_var(idFile,idUcov,uCov) !,(1,1,1,1),(nLong,nLat,nAlt,nTime),
295if(ierror /= nf90_noerr) then
296     print *, trim(nf90_strerror(ierror))
297     stop "Stopped put ucov"
298end if
299ierror=nf90_put_var(idFile,idVcov,vCov) !,(1,1,1,1),(nLong,nLat,nAlt,nTime),
300if(ierror /= nf90_noerr) then
301     print *, trim(nf90_strerror(ierror))
302     stop "Stopped put vcov"
303end if
304!**********
305!Write variable teta
306ierror=nf90_put_var(idFile,idTeta,teta) !,(1,1,1,1),(nLong,nLat,nAlt,nTime),
307if(ierror /= nf90_noerr) then
308     print *, trim(nf90_strerror(ierror))
309     stop "Stopped put teta"
310end if
311!**********
312!Close
313ierror=nf90_close(idFile)
314if(ierror /= nf90_noerr) then
315     print *, trim(nf90_strerror(ierror))
316     stop "Stopped close"
317end if
318
319contains
320
321subroutine interface_user(fileNetcdf,duNat,dvNat,dTeta)
322
323implicit none
324
325 character (len=*), intent(out) :: fileNetcdf
326real, intent(out) :: duNat,dvNat,dTeta
327
328
329
330
331
332
333end subroutine interface_user
334
335subroutine mean_filter(matrix,m,n,dimFilter,matrixOut)
336
337implicit none
338
339integer, intent(in) :: m,n,dimFilter
340real, dimension(m,n), intent(in) :: matrix
341real, dimension(m,n), intent(out) :: matrixOut
342real, dimension(m+2*dimFilter,n+2*dimFilter) :: matrixTmp,matrixConv
343real, dimension(2*dimFilter+1,2*dimFilter+1) :: kernel
344integer :: mm,nn
345
346!matrix containing input matrix plus zero edges
347matrixTmp(:,:)=0.0
348mm=m+2*dimFilter
349nn=n+2*dimFilter
350do i=dimFilter+1,mm-dimFilter
351        do j=dimFilter+1,nn-dimFilter
352                matrixTmp(i,j)=matrix(i-dimFilter,j-dimFilter)
353        end do
354end do
355!print*,matrixTmp(5,:)
356
357!filter creation
358kernel(:,:) = 1.0/((2*dimFilter+1)*(2*dimFilter+1))
359!print*,kernel
360
361!matrix convolution with the filter
362do i=dimFilter+1,mm-dimFilter
363        do j=dimFilter+1,nn-dimFilter
364                matrixConv(i,j)=sum(matrixTmp(i-dimFilter:i+dimFilter,j-dimFilter:j+dimFilter)*kernel)
365        end do
366end do
367i=5
368j=5
369!print*,shape(matrixTmp(i-dimFilter:i+dimFilter,j-dimFilter:j+dimFilter))
370!print*,matrixConv(5,:)
371
372!extraction of submatrix without zeros edges
373matrixOut(:,:)=0.0
374do i=1,m
375        do j=1,n
376                matrixOut(i,j)=matrixConv(i+dimFilter,j+dimFilter)
377        end do
378end do
379
380end subroutine mean_filter
381
382end program createRandomStart
383
Note: See TracBrowser for help on using the repository browser.