EMPIRE DA  v1.9.1
Data assimilation codes using EMPIRE communication
 All Classes Files Functions Variables Pages
perturb_particle.f90
Go to the documentation of this file.
1 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2 !!! Time-stamp: <2016-10-18 15:12:04 pbrowne>
3 !!!
4 !!! Collection of routines to perturb states
5 !!! Copyright (C) 2014 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 
30 subroutine perturb_particle(x)
31  use output_empire, only : emp_e
32  use sizes
33  use comms
34  use pf_control
35  integer, parameter :: rk=kind(1.0d0)
36  real(kind=rk), dimension(state_dim), intent(inout) :: x
37  real(kind=rk), dimension(state_dim) :: rdom,y,kgain
38  character(256) :: filename
39 
40  select case(pf%init)
41  case('P')
42  call normalrandomnumbers1d(0.0d0,1.0d0,state_dim,rdom)
43  call qhalf(1,rdom,y)
44  rdom = 0.0_rk
45  kgain = 0.0_rk
46  call update_state(rdom,x,kgain,y)
47  x = rdom
48  case('N')
49  call normalrandomnumbers1d(0.0d0,1.0d0,state_dim,rdom)
50  kgain = 0.0_rk
51  call update_state(y,x,kgain,rdom)
52  x = y
53  case('R')
54  !get ensemble member from the restart folder
55  write(filename,'(A,i2.2,A)') 'rstrt/',pfrank,'.state'
56  call get_state(x,filename)
57  case('S')
58  !get ensemble member from the start folder
59  if(pf%gen_data) then
60  write(filename,'(A,i2.2,A)') 'start/',32,'.state'
61  else
62  write(filename,'(A,i2.2,A)') 'start/',pfrank,'.state'
63  print*,'pf #',pfrank,' starting from ',filename
64  end if
65  call get_state(x,filename)
66  case('B')
67  call normalrandomnumbers1d(0.0d0,1.0d0,state_dim,rdom)
68  call bhalf(1,rdom,y)
69  rdom = 0.0_rk
70  kgain = 0.0_rk
71  call update_state(rdom,x,kgain,y)
72  x = rdom
73  case('U')
74  call user_perturb_particle(state_dim,x)
75  case('Z')
76  !no perturbation. x remains as is
77  case default
78  write(emp_e,*) 'ERROR: incorrect pf%init selected in perturb_particle: ',pf%init
79  stop
80  end select
81 
82 
83 end subroutine perturb_particle
84 
subroutine user_perturb_particle(n, x)
Subroutine to perturb state vector as defined by the user governed by the init option.
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 that stores the information about the outputting from empire.
subroutine update_state(state, fpsi, kgain, betan)
Subroutine to update the state.
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 get_state(state, filename)
subroutine to read the state vector from a named file as an unformatted fortran file ...
Definition: data_io.f90:258
subroutine perturb_particle(x)
Subroutine to perturb state vector governed by the init option.
module pf_control holds all the information to control the the main program
Definition: pf_control.f90:29
subroutine qhalf(nrhs, x, Qx)
subroutine to take a full state vector x and return in state space.