FV3 Bundle
generate_locations.f90
Go to the documentation of this file.
1 ! (C) Copyright 2009-2016 ECMWF.
2 !
3 ! This software is licensed under the terms of the Apache Licence Version 2.0
4 ! which can be obtained at http://www.apache.org/licenses/LICENSE-2.0.
5 ! In applying this licence, ECMWF does not waive the privileges and immunities
6 ! granted to it by virtue of its status as an intergovernmental organisation nor
7 ! does it submit to any jurisdiction.
8 
9 subroutine generate_locations(c_conf,nlocs,ntimes,bgn,step,times,obsloc)
10 
11 use iso_c_binding
12 use config_mod
13 use datetime_mod
14 use duration_mod
16 use kinds
17 
18 implicit none
19 type(c_ptr), intent(in) :: c_conf
20 integer, intent(in) :: nlocs
21 integer, intent(in) :: ntimes
22 type(datetime), intent(in) :: bgn
23 type(duration), intent(in) :: step
24 type(datetime), intent(inout) :: times(nlocs*ntimes)
25 type(obs_vect), intent(inout) :: obsloc
26 
27 integer :: ijk, jobs, intinv, iobs, jstep, ii, jj, kk, ij, nx, ny
28 real(kind=kind_real) :: goldinv, dx, dy
29 real(kind=kind_real), allocatable :: xxyyzz(:,:)
30 type(datetime) :: now
31 
32 nx = config_get_int(c_conf, "nx");
33 ny = config_get_int(c_conf, "ny");
34 dx=1.0_kind_real/real(nx,kind_real)
35 dy=1.0_kind_real/real(ny,kind_real)
36 
37 goldinv = 0.5_kind_real*(sqrt(5.0_kind_real)-1.0_kind_real) ! 1/(golden ratio)
38 intinv = nint(goldinv*real(nx*ny*2,kind_real))
39 
40 allocate(xxyyzz(3,nlocs))
41 ijk=0
42 do jobs=1,nlocs
43  ijk=ijk+intinv
44  ijk=modulo(ijk,2*nx*(ny-2))
45  kk=ijk/(nx*(ny-2))
46  ij=ijk-kk*nx*(ny-2)
47  jj=ij/nx
48  ii=ij-jj*nx
49  jj=jj+1
50  if (ii<0 .or. ii>nx-1) call abor1_ftn ("generate_locations: error ii")
51  if (jj<1 .or. jj>ny-2) call abor1_ftn ("generate_locations: error jj")
52  if (kk<0 .or. kk>1) call abor1_ftn ("generate_locations: error kk")
53  xxyyzz(1,jobs)=real(ii,kind_real)*dx
54  xxyyzz(2,jobs)=real(jj,kind_real)*dy
55  xxyyzz(3,jobs)=real(kk+1,kind_real)
56 enddo
57 
58 call obsvec_setup(obsloc,3,nlocs*ntimes)
59 
60 now = bgn
61 iobs=0
62 do jstep=1,ntimes
63  do jobs=1,nlocs
64  iobs=iobs+1
65  times(iobs)=now
66  obsloc%values(:,iobs)=xxyyzz(:,jobs)
67  enddo
68  call datetime_update(now,step)
69 enddo
70 
71 call datetime_delete(now)
72 deallocate(xxyyzz)
73 
74 end subroutine generate_locations
subroutine generate_locations(c_conf, nlocs, ntimes, bgn, step, times, obsloc)
subroutine, public obsvec_setup(self, nc, no)
Fortran module handling observation vectors.