EMPIRE DA  v1.9.1
Data assimilation codes using EMPIRE communication
 All Classes Files Functions Variables Pages
fourdenvardata.f90
Go to the documentation of this file.
1 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2 !!! Time-stamp: <2015-10-13 13:51:31 pbrowne>
3 !!!
4 !!! module to store data for 4DEnVar
5 !!! Copyright (C) 2015 Philip A. Browne
6 !!!
7 !!! This program is free software: you can redistribute it and/or modify
8 !!! it under the terms of the GNU General Public License as published by
9 !!! the Free Software Foundation, either version 3 of the License, or
10 !!! (at your option) any later version.
11 !!!
12 !!! This program is distributed in the hope that it will be useful,
13 !!! but WITHOUT ANY WARRANTY; without even the implied warranty of
14 !!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15 !!! GNU General Public License for more details.
16 !!!
17 !!! You should have received a copy of the GNU General Public License
18 !!! along with this program. If not, see <http://www.gnu.org/licenses/>.
19 !!!
20 !!! Email: p.browne @ reading.ac.uk
21 !!! Mail: School of Mathematical and Physical Sciences,
22 !!! University of Reading,
23 !!! Reading, UK
24 !!! RG6 6BB
25 !!!
26 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
27 
28 
33  implicit none
34  integer :: m
35  real(kind=kind(1.0d0)), allocatable, dimension(:) :: xb
37  real(kind=kind(1.0d0)), allocatable, dimension(:,:) :: x0
41  real(kind=kind(1.0d0)), allocatable, dimension(:,:) :: xt
43 contains
45  use sizes, only : state_dim
46  use comms, only : cnt,nens,pfrank
47  m = nens-1
48  allocate(xb(state_dim))
49  if(pfrank .ne. 0) then
50  allocate(x0(state_dim,cnt))
51  allocate(xt(state_dim,cnt))
52  else
53  allocate(x0(state_dim,cnt-1)) !no perturbation associated
54  !with optimization solution
55  allocate(xt(state_dim,cnt))!here we include the optimization
56  !current solution as well as the perturbations
57  end if
58 
59  end subroutine allocate4denvardata
60 
62  subroutine read_background_term()
63  real(kind=kind(1.0d0)), dimension(3) :: tempxb
64  !read xb from file
65  print*,'READ xb from file not implemented'
66  !stop 10
67  !type in a normally distributed random vector
68  call normalrandomnumbers1d(0.0d0,1.0d0,3,tempxb)
69  ! tempxb = (/1.8368174447696750d-1,-1.3102740999288161d-2, &
70  ! &-5.9318375014562030d-1/)
71  call bhalf(1,tempxb,xb)
72  print*,'background guess pert = '
73  print*,xb
74 
75  xb = xb + (/-3.12346395, -3.12529803, 20.69823159/)
76  print*,'background guess xb = '
77  print*,xb
78  end subroutine read_background_term
79 
81  deallocate(xb)
82  end subroutine deallocate4denvardata
83 
90  use comms, only : cnt,pfrank
91  use sizes
92  real(kind=kind(1.0d0)),allocatable, dimension(:,:) :: tempx0
93 
94  if(pfrank .ne. 0) then
95  ! READ x0(particle) from file
96  print*,'READ FROM FILE NOT IMPLEMENTED'
97  !stop 7
98  allocate(tempx0(state_dim,cnt))
99  call normalrandomnumbers2d(0.0d0,1.0d0,state_dim,cnt,tempx0)
100  ! print*,'sup tempx0 = ',tempx0
101  call bhalf(cnt,tempx0,x0)
102  ! print*,'arse x0 = ',x0
103  deallocate(tempx0)
104  else!no perturbation
105  !associated with optimization solution
106 
107  ! READ x0(particle-1) from file
108  print*,'READ FROM FILE NOT IMPLEMENTED'
109  !stop 8
110  allocate(tempx0(state_dim,cnt-1))
111  call normalrandomnumbers2d(0.0d0,1.0d0,state_dim,cnt-1,tempx0)
112  ! print*,'sup tempx0 = ',tempx0
113  call bhalf(cnt-1,tempx0,x0)
114  ! print*,'arse x0 = ',x0
115  deallocate(tempx0)
116  end if
117  x0 = x0!/(real(m-1,kind(1.0d0))**0.5d0)
118 
119  print*,'read_ensemble_perturbation_matrix'
120  print*,x0
121  print*,'read_ensemble_perturbation_matrix'
122  print*,'x0(1,3) = ',x0(1,3)
123  print*,'x0(3,2) = ',x0(3,2)
124  end subroutine read_ensemble_perturbation_matrix
125 
126 
127 
128 end module fourdenvardata
subroutine normalrandomnumbers2d(mean, stdev, n, k, phi)
generate two dimensional Normal random numbers
Definition: gen_rand.f90:109
subroutine bhalf(nrhs, x, bx)
subroutine to take a full state vector x and return in state space.
Module containing EMPIRE coupling data.
Definition: comms.f90:57
module holding data specific for 4denvar, not var itself. this is necessary because of the difference...
subroutine read_ensemble_perturbation_matrix
subroutine to read in the ensemble perturbation matrix
Module that stores the dimension of observation and state spaces.
Definition: sizes.f90:29
subroutine normalrandomnumbers1d(mean, stdev, n, phi)
generate one dimension of Normal random numbers
Definition: gen_rand.f90:74
subroutine deallocate4denvardata
subroutine allocate4denvardata
subroutine read_background_term()
subroutine to read xb from file