EMPIRE DA  v1.9.1
Data assimilation codes using EMPIRE communication
 All Classes Files Functions Variables Pages
linear_empire_vader.f90
Go to the documentation of this file.
1 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2 ! linear_empire_vader.f90 Implements a linear model with EMPIRE
3 ! coupling extended to include reverse communication (VADER)
4 !
5 !
6 !The MIT License (MIT)
7 !
8 !Copyright (c) 2015 Philip A. Browne
9 !
10 !Permission is hereby granted, free of charge, to any person obtaining a copy
11 !of this software and associated documentation files (the "Software"), to deal
12 !in the Software without restriction, including without limitation the rights
13 !to use, copy, modify, merge, publish, distribute, sublicense, and/or sell
14 !copies of the Software, and to permit persons to whom the Software is
15 !furnished to do so, subject to the following conditions:
16 !
17 !The above copyright notice and this permission notice shall be included in all
18 !copies or substantial portions of the Software.
19 !
20 !THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR
21 !IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
22 !FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE
23 !AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
24 !LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM,
25 !OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN THE
26 !SOFTWARE.
27 !
28 !Email: p.browne@reading.ac.uk
29 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
30 
44 program linear
45  implicit none
46  include 'mpif.h'
47  real(kind=kind(1.0D0)), dimension(:), allocatable :: x
48  integer :: i,n,maxt
49  integer, dimension(2) :: data
50  integer :: mpi_err,mdl_id,cpl_root,cpl_mpi_comm
51  integer, dimension(MPI_STATUS_SIZE) :: mpi_status
52 
53  ! SET UP EMPIRE COMMUNICATORS
54  call initialise_mpi(mdl_id,cpl_root,cpl_mpi_comm)
55 
56 
57  ! GET THE COUPLET FROM THE DATA ASSIMILATION CODE
58  call mpi_recv(data,2,mpi_integer,cpl_root,mpi_any_tag,cpl_mpi_comm&
59  &,mpi_status,mpi_err)
60 
61  ! PUT THE DATA INTO THE RIGHT VARIABLES AND ALLOCATE
62  n = data(1)
63  maxt = data(2)
64  allocate(x(n))
65  x = 1.0d0
66 
67 
68  ! DO THE INITIAL SEND AND RECIEVES FROM EMPIRE
69  call mpi_send(x,n,mpi_double_precision,cpl_root&
70  &,1,cpl_mpi_comm,mpi_err)
71  call mpi_recv(x,n,mpi_double_precision,cpl_root&
72  &,mpi_any_tag,cpl_mpi_comm,mpi_status,mpi_err)
73 2 continue
74 
75  write(6,*) 0,x
76  ! START THE TIMESTEP LOOP
77  do i = 1,maxt
78 
79  ! CALL THE SIMPLE LINEAR MODEL AND UPDATE MODEL STATE
80  x = f(n,x)
81 
82  ! DO THE TIMESTEP LOOP SEND AND RECIEVES WITH EMPIRE
83  call mpi_send(x,n,mpi_double_precision,cpl_root&
84  &,1,cpl_mpi_comm,mpi_err)
85  call mpi_recv(x,n,mpi_double_precision,cpl_root&
86  &,mpi_any_tag,cpl_mpi_comm,mpi_status,mpi_err)
87 
88  write(6,*) i,x
89  if(mpi_status(mpi_tag) .eq. 1) then
90  go to 1 !simply continue code
91  elseif(mpi_status(mpi_tag) .eq. 2) then
92  go to 2 !restart code
93  elseif(mpi_status(mpi_tag) .eq. 3) then
94  go to 3 !end code
95  else
96  print*,'Linear model error: unknown MPI_TAG: ',mpi_status(mpi_tag)
97  stop '-1'
98  end if
99 1 continue
100 
101 
102  !END TIMESTEP LOOP
103  end do
104 3 continue
105 
106  ! END CODE
107  deallocate(x)
108  call mpi_finalize(mpi_err)
109 contains
110 
111  ! DEFINE A SIMPLE LINEAR MODEL
112  function f (n,x)
113  implicit none
114  integer, intent(in) :: n
115  real(kind=kind(1.0D0)),intent(in),dimension (n) :: x
116  real(kind=kind(1.0D0)),dimension(n) :: f
117  f = x
118  end function f
119 
120  ! STANDARD EMPIRE INITIALISATION SUBROUTINE
121  subroutine initialise_mpi(mdl_id,cpl_root,cpl_mpi_comm)
122  implicit none
123  include 'mpif.h'
124  integer, intent(out) :: mdl_id,cpl_root,cpl_mpi_comm
125  integer :: mdl_num_proc=1
126  integer :: mpi_err,world_size,world_id
127  integer :: cpl_colour
128  integer :: particle_id,nens, da, nda
129  integer :: mdl_mpi_comm,mdlcolour
130  integer :: tmp_mdls_comm,models_id,models_size
131  call mpi_init(mpi_err)
132  da = 0
133  call mpi_comm_rank(mpi_comm_world,world_id,mpi_err)
134  call mpi_comm_size(mpi_comm_world,world_size,mpi_err)
135  call mpi_comm_split(mpi_comm_world,da,world_id,tmp_mdls_comm,mpi_err)
136  call mpi_comm_size(tmp_mdls_comm,models_size,mpi_err)
137  call mpi_comm_rank(tmp_mdls_comm,models_id, mpi_err)
138  mdlcolour = models_id/mdl_num_proc
139  call mpi_comm_split(tmp_mdls_comm,mdlcolour,models_id,mdl_mpi_comm,mpi_err)
140  call mpi_comm_rank(mdl_mpi_comm,mdl_id,mpi_err)
141  if(mdl_id .eq. 0) then
142  cpl_colour = 9999
143  else
144  cpl_colour = mpi_undefined
145  end if
146  call mpi_comm_split(mpi_comm_world,cpl_colour,mdlcolour,cpl_mpi_comm,mpi_err)
147  if(mdl_id .eq. 0) then
148  call mpi_comm_size(cpl_mpi_comm,nens,mpi_err)
149  call mpi_comm_rank(cpl_mpi_comm,particle_id,mpi_err)
150  nda = world_size-models_size;nens = nens - nda
151  cpl_root = ((nda*particle_id)/nens)+nens
152  if(nda ==0) cpl_root = particle_id
153  else
154  cpl_root = -1
155  end if
156  end subroutine initialise_mpi
157 
158 end program linear
subroutine initialise_mpi(mdl_id, cpl_root, cpl_mpi_comm)
real(kind=kind(1.0d0)) function, dimension(n) f(n, x)
program linear
program to implement a simple linear model of no use to anyone but for testing and debugging purposes...