EMPIRE DA  v1.9.1
Data assimilation codes using EMPIRE communication
 All Classes Files Functions Variables Pages
three_d_var.f90
Go to the documentation of this file.
1 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2 !!! Time-stamp: <2016-10-18 15:39:57 pbrowne>
3 !!!
4 !!! Program to implement 3DVar
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 
29 subroutine three_d_var(x)
30  use output_empire, only : emp_e
31  use threedvar_data
32  use sizes
33  use var_data
34  implicit none
35  integer, parameter :: rk = kind(1.0d0)
36  real(kind=rk), dimension(state_dim), intent(inout) :: x
37 
38  vardata%n = state_dim
39 
40  if(allocated(xb)) deallocate(xb)
41  allocate(xb(state_dim))
42  xb = x
43  vardata%total_timesteps = 1
44  call allocate_vardata
45 
46  vardata%x0 = x
47 
48 
49  select case (vardata%opt_method)
50  case('cg')
51 
52  call subroutine_cg(vardata%cg_method,vardata%n,vardata%cg_eps,vardata%x0)
53 
54  case('lbfgs')
55  call lbfgs_sub(vardata%n,vardata%lbfgs_factr,vardata%lbfgs_pgtol,vardata%x0)
56 
57  case('lbfgsb')
59  call lbfgsb_sub(vardata%n,vardata%lbfgs_factr,vardata%lbfgs_pgtol,&
60  vardata%x0,vardata%nbd,vardata%l,vardata%u)
61 
62  case default
63  write(emp_e,*) 'three_d_var ERROR: vardata%opt_method incorrect. Stopping'
64  stop '-78'
65  end select
66 
67  x = vardata%x0
69 
70 end subroutine three_d_var
subroutine allocate_vardata
subroutine to allocate space for 4denvar
Definition: var_data.f90:328
subroutine deallocate_vardata
subroutine to deallocate space for 4denvar
Definition: var_data.f90:343
Module that stores the information about the outputting from empire.
subroutine lbfgsb_sub(n, factr_in, pgtol_in, x, nbd, l, u)
Limited memory BFGS bound constrained optimization code as callable subroutine.
Definition: lbfgsb_sub.f90:210
Module that stores the dimension of observation and state spaces.
Definition: sizes.f90:29
module to store stuff for 3DVar
subroutine three_d_var(x)
Definition: three_d_var.f90:29
subroutine subroutine_cg(method, n, epsin, x)
Nonlinear Conjugate gradient method as callable subroutine.
Definition: cgsub.f90:27
subroutine read_lbfgsb_bounds
subroutine to somehow read in bounds data
Definition: var_data.f90:353
subroutine lbfgs_sub(n, factr_in, pgtol_in, x)
Limited memory BFGS unconstrained optimization code as callable subroutine.
Definition: lbfgs_sub.f90:198
module holding data for variational problems
Definition: var_data.f90:29