33 real(kind=kind(1.0D0)) :: dt=2.5d-3
34 real(kind=kind(1.0D0)),
allocatable,
dimension(:,:) :: x,k1,k2,k3,k4
35 real(kind=kind(1.0d0)) :: F=8.0d0
36 real(kind=kind(1.0d0)) :: alpha=0.5d0
37 real(kind=kind(1.0d0)) :: delta=0.5d0
38 real(kind=kind(1.0d0)) :: epsilon=2.5d-3
39 real(kind=kind(1.0d0)) :: gamma=0.1d0
41 integer :: total_timesteps=1000
43 integer :: mpi_err,mdl_id,cpl_root,cpl_mpi_comm
44 integer,
dimension(MPI_STATUS_SIZE) :: mpi_status
57 inquire(file=
'l96.nml', exist=l96_exists)
59 open(32,file=
'l96.nml',iostat=ios,action=
'read'&
61 if(ios .ne. 0) stop
'Cannot open l96.nml'
72 allocate(x(n,3),k1(n,3),k2(n,3),k3(n,3),k4(n,3))
80 if(mdl_id .eq. 0)
then
81 call mpi_send(x(:,1:2),n*2,mpi_double_precision,cpl_root&
82 &,1,cpl_mpi_comm,mpi_err)
83 call mpi_recv(x(:,1:2),n*2,mpi_double_precision,cpl_root&
84 &,mpi_any_tag,cpl_mpi_comm,mpi_status,mpi_err)
89 do t = 1,total_timesteps
90 k1 =
g(x , n , f ,alpha,delta,epsilon,gamma)
91 k2 =
g(x +0.5d0 * dt * k1 , n , f ,alpha,delta,epsilon,gamma)
92 k3 =
g(x +0.5d0 * dt * k2 , n , f ,alpha,delta,epsilon,gamma)
93 k4 =
g(x + dt * k3 , n , f ,alpha,delta,epsilon,gamma)
94 x = x + dt *( k1 + 2.0d0 *( k2 + k3 ) + k4 )/6.0d0
96 if(mdl_id .eq. 0)
then
97 call mpi_send(x(:,1:2),n*2,mpi_double_precision,cpl_root&
98 &,1,cpl_mpi_comm,mpi_err)
99 call mpi_recv(x(:,1:2),n*2,mpi_double_precision,cpl_root&
100 &,mpi_any_tag,cpl_mpi_comm,mpi_status,mpi_err)
104 select case(mpi_status(mpi_tag))
114 call mpi_finalize(mpi_err)
117 function g (X , N, F ,alpha,delta,epsilon,gamma)
119 real(kind=kind(1.0D0)),
intent(in),
dimension(N,3) :: X
120 integer,
intent(in) :: N
121 real(kind=kind(1.0D0)),
dimension(N,3) :: g
122 real(kind=kind(1.0D0)),
intent(in) :: F
123 real(kind=kind(1.0D0)),
intent(in) :: alpha
124 real(kind=kind(1.0D0)),
intent(in) :: delta
125 real(kind=kind(1.0D0)),
intent(in) :: epsilon
126 real(kind=kind(1.0D0)),
intent(in) :: gamma
129 g(1,1) = (1.0d0-delta)*(x(2,1)-x(n-1,1) )*x(n,1) + delta*(x(n,1)*x(2,2)-x(n-1,1)*x(n,2)) - x(1,1) + f
131 g(1,3) = (-x(1,2) + alpha**2*(x(2,2)-2*x(1,2)+x(n,2)) + x(1,1) - gamma*epsilon**2*x(1,3))/(epsilon**2)
133 g(2,1) = (1.0d0-delta)*(x(3,1)-x(n,1) )*x(1,1) + delta*(x(1,1)*x(3,2)-x(n,1)*x(1,2)) - x(2,1) + f
135 g(2,3) = (-x(2,2) + alpha**2*(x(3,2)-2*x(2,2)+x(1,2)) + x(2,1) - gamma*epsilon**2*x(2,3))/(epsilon**2)
139 g(j,1) = (1.0d0-delta)*(x(j+1,1)-x(j-2,1))*x(j-1,1) + delta*(x(j-1,1)*x(j+1,2)-x(j-2,1)*x(j-1,2)) - x(j,1) + f
141 g(j,3) = (-x(j,2) + alpha**2*(x(j+1,2)-2*x(j,2)+x(j-1,2)) + x(j,1) - gamma*epsilon**2*x(j,3))/(epsilon**2)
144 g(n,1) = (1.0d0-delta)*(x(1,1)-x(n-2,1) )*x(n-1,1) + delta*(x(n-1,1)*x(1,2)-x(n-2,1)*x(n-1,2)) - x(n,1) + f
146 g(n,3) = (-x(n,2) + alpha**2*(x(1,2)-2*x(n,2)+x(n-1,2)) + x(n,1) - gamma*epsilon**2*x(n,3))/(epsilon**2)
154 integer,
intent(out) :: mdl_id,cpl_root,cpl_mpi_comm
155 integer :: mdl_num_proc=1
156 integer :: mpi_err,world_size,world_id
157 integer :: cpl_colour
158 integer :: particle_id,nens, da, nda
159 integer :: mdl_mpi_comm,mdlcolour
160 integer :: tmp_mdls_comm,models_id,models_size
161 call mpi_init(mpi_err)
163 call mpi_comm_rank(mpi_comm_world,world_id,mpi_err)
164 call mpi_comm_size(mpi_comm_world,world_size,mpi_err)
165 call mpi_comm_split(mpi_comm_world,da,world_id,tmp_mdls_comm,mpi_err)
166 call mpi_comm_size(tmp_mdls_comm,models_size,mpi_err)
167 call mpi_comm_rank(tmp_mdls_comm,models_id, mpi_err)
168 mdlcolour = models_id/mdl_num_proc
169 call mpi_comm_split(tmp_mdls_comm,mdlcolour,models_id,mdl_mpi_comm,mpi_err)
170 call mpi_comm_rank(mdl_mpi_comm,mdl_id,mpi_err)
171 if(mdl_id .eq. 0)
then
174 cpl_colour = mpi_undefined
176 call mpi_comm_split(mpi_comm_world,cpl_colour,mdlcolour,cpl_mpi_comm,mpi_err)
177 if(mdl_id .eq. 0)
then
178 call mpi_comm_size(cpl_mpi_comm,nens,mpi_err)
179 call mpi_comm_rank(cpl_mpi_comm,particle_id,mpi_err)
180 nda = world_size-models_size;nens = nens - nda
181 cpl_root = ((nda*particle_id)/nens)+nens
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_slow_fast