EMPIRE DA  v1.9.1
Data assimilation codes using EMPIRE communication
 All Classes Files Functions Variables Pages
Lorenz96_empire.f90
Go to the documentation of this file.
1 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2 ! Lorenz96_empire.f90 Implements Lorenz 1996 model with EMPIRE coupling
3 !
4 !The MIT License (MIT)
5 !
6 !Copyright (c) 2015 Philip A. Browne
7 !
8 !Permission is hereby granted, free of charge, to any person obtaining a copy
9 !of this software and associated documentation files (the "Software"), to deal
10 !in the Software without restriction, including without limitation the rights
11 !to use, copy, modify, merge, publish, distribute, sublicense, and/or sell
12 !copies of the Software, and to permit persons to whom the Software is
13 !furnished to do so, subject to the following conditions:
14 !
15 !The above copyright notice and this permission notice shall be included in all
16 !copies or substantial portions of the Software.
17 !
18 !THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR
19 !IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
20 !FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE
21 !AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
22 !LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM,
23 !OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN THE
24 !SOFTWARE.
25 !
26 !Email: p.browne@reading.ac.uk
27 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
28 
29 program lorenz96
30  implicit none
31  include 'mpif.h'
32  real(kind=kind(1.0D0)) :: dt=1.0d-2
33  real(kind=kind(1.0D0)), allocatable, dimension(:) :: x,k1,k2,k3,k4
34  real(kind=kind(1.0d0)) :: F=8.0d0
35  integer :: N=40
36  integer :: total_timesteps=100
37  integer :: t
38  integer :: mpi_err,mdl_id,cpl_root,cpl_mpi_comm
39  integer, dimension(MPI_STATUS_SIZE) :: mpi_status
40  integer :: ios
41  logical :: l96_exists
42 
43  namelist/l96/n,&
44  &total_timesteps,&
45  &f,&
46  &dt
47 
48  inquire(file='l96.nml', exist=l96_exists)
49  if(l96_exists) then
50  open(32,file='l96.nml',iostat=ios,action='read'&
51  &,status='old')
52  if(ios .ne. 0) stop 'Cannot open l96.nml'
53  read(32,nml=l96)
54  close(32)
55  end if
56 
57 
58 
59 
60 
61 
62 
63 
64  call initialise_mpi(mdl_id,cpl_root,cpl_mpi_comm)
65 
66 
67  allocate(x(n),k1(n),k2(n),k3(n),k4(n))
68  x = f
69  x(n/2) = f+0.05
70 
71 
72  if(mdl_id .eq. 0) then
73  call mpi_send(x,n,mpi_double_precision,cpl_root&
74  &,1,cpl_mpi_comm,mpi_err)
75  call mpi_recv(x,n,mpi_double_precision,cpl_root&
76  &,mpi_any_tag,cpl_mpi_comm,mpi_status,mpi_err)
77  end if
78 2 continue
79  print*,x
80 
81 
82  do t = 1,total_timesteps
83  k1 = g(x , n , f )
84  k2 = g(x +0.5d0 * dt * k1 , n , f )
85  k3 = g(x +0.5d0 * dt * k2 , n , f )
86  k4 = g(x + dt * k3 , n , f )
87  x = x + dt *( k1 + 2.0d0 *( k2 + k3 ) + k4 )/6.0d0
88 
89  if(mdl_id .eq. 0) then
90  call mpi_send(x,n,mpi_double_precision,cpl_root&
91  &,1,cpl_mpi_comm,mpi_err)
92  call mpi_recv(x,n,mpi_double_precision,cpl_root&
93  &,mpi_any_tag,cpl_mpi_comm,mpi_status,mpi_err)
94  end if
95  print*,x
96 
97  select case(mpi_status(mpi_tag))
98  case(2)
99  go to 2
100  case(3)
101  go to 3
102  case default
103  end select
104 
105  end do
106 3 continue
107  call mpi_finalize(mpi_err)
108 
109 contains
110  function g (x , N, F )
111  implicit none
112  real(kind=kind(1.0D0)),intent(in),dimension(0:N-1) :: x
113  integer, intent(in) :: N
114  real(kind=kind(1.0D0)),dimension(0:N-1) :: g
115  real(kind=kind(1.0D0)),intent(in) :: F
116  integer :: j
117  do j = 0,n-1
118  g(j) = ( x(modulo(j+1,n)) -x( modulo(j-2,n)) )*&
119  &x(modulo(j-1,n)) - x(j) + f
120  end do
121  end function g
122  subroutine initialise_mpi(mdl_id,cpl_root,cpl_mpi_comm)
123  implicit none
124  include 'mpif.h'
125  integer, intent(out) :: mdl_id,cpl_root,cpl_mpi_comm
126  integer :: mdl_num_proc=1
127  integer :: mpi_err,world_size,world_id
128  integer :: cpl_colour
129  integer :: particle_id,nens, da, nda
130  integer :: mdl_mpi_comm,mdlcolour
131  integer :: tmp_mdls_comm,models_id,models_size
132  call mpi_init(mpi_err)
133  da = 0
134  call mpi_comm_rank(mpi_comm_world,world_id,mpi_err)
135  call mpi_comm_size(mpi_comm_world,world_size,mpi_err)
136  call mpi_comm_split(mpi_comm_world,da,world_id,tmp_mdls_comm,mpi_err)
137  call mpi_comm_size(tmp_mdls_comm,models_size,mpi_err)
138  call mpi_comm_rank(tmp_mdls_comm,models_id, mpi_err)
139  mdlcolour = models_id/mdl_num_proc
140  call mpi_comm_split(tmp_mdls_comm,mdlcolour,models_id,mdl_mpi_comm,mpi_err)
141  call mpi_comm_rank(mdl_mpi_comm,mdl_id,mpi_err)
142  if(mdl_id .eq. 0) then
143  cpl_colour = 9999
144  else
145  cpl_colour = mpi_undefined
146  end if
147  call mpi_comm_split(mpi_comm_world,cpl_colour,mdlcolour,cpl_mpi_comm,mpi_err)
148  if(mdl_id .eq. 0) then
149  call mpi_comm_size(cpl_mpi_comm,nens,mpi_err)
150  call mpi_comm_rank(cpl_mpi_comm,particle_id,mpi_err)
151  nda = world_size-models_size;nens = nens - nda
152  cpl_root = ((nda*particle_id)/nens)+nens
153  else
154  cpl_root = -1
155  end if
156  end subroutine initialise_mpi
157 end program lorenz96
subroutine initialise_mpi(mdl_id, cpl_root, cpl_mpi_comm)
real(kind=kind(1.0d0)) function, dimension(n, 3) g(X, N, F, alpha, delta, epsilon, gamma)
program lorenz96