局部与非局部突触耦合链接下神经元网络中的行波


Traveling chimera pattern in a neuronal network under local gap and nonlocal coupling

阅读

复现

figure 1

figure 1

  • 黑线是局部耦合,通过电紧张电位间隙链接到它的最近邻居。红虚线是非局部耦合,通过化学突触耦合到远的邻居节点。
module HR
    implicit none
    real,parameter :: h=0.01
    integer,parameter :: N=100,MaxT=398700
    integer,allocatable :: neighbour_matrix_2(:,:),neighbour_matrix_2P(:,:)
    real :: x(N),y(N),z(N)
contains

    subroutine x0_y0_z0()
        implicit none
        integer :: k
        real :: x1,x2,x3
        call random_seed()
        do k=1,N,1
            call random_number(x1)
            call random_number(x2)
            call random_number(x3)
            x(k)=0.001*(k-N/2.0)+x1
            y(k)=0.002*(k-N/2.0)+x2
            z(k)=0.003*(k-N/2.0)+x3
            write(10,*) x(k),y(k),z(k)
        end do
    end subroutine x0_y0_z0

    subroutine fnf(p)
        implicit none
        integer :: i,j,p,t
        real :: func_x(N),func_y(N),func_z(N),a,b,c,d,k1,k2,II,JJ,CC,r,s,x_R,v_s,lambda,thet_s
        k1=0.0
        k2=0.0
        a=1.0
        b=3.0
        c=1.0
        d=5.0
        r=0.01
        s=5.0
        v_s=2.0
        x_R=-1.6
        lambda=10.0
        thet_s=-0.25
        II=0.0
        do t=1,MaxT,1
            if(mod(t,5000)==0) then
                write(*,*) t
            end if
            do i=1,N,1
                JJ=0.0
                C=0.0
                do j=1,N,1
                    JJ=JJ+neighbour_matrix_2(i,j)*(x(j)-x(i))
                    C=C+neighbour_matrix_2P(i,j)*(1.0/(1.0+exp(-lambda*(x(j)-thet_s))))
                end do
                func_x(i)=x(i)+h*(y(i)-a*x(i)**3+b*x(i)**2-z(i)+II+k1*JJ+((-k2/(2.0*p-2.0))*(v_s-x(i))*CC))
                func_y(i)=y(i)+h*(c-d*x(i)**2-y(i))
                func_z(i)=z(i)+h*(r*(s*(x(i)-x_R)-z(i)))
                if(mod(t,50)==0.and.t>353700) then
                    write(20,*) i,t*h,x(i)
                end if
                if(i==1.and.t>353700) then
                    write(30,*) t*h,x(i)
                    write(70,*) x(i),y(i),z(i)
                    write(*,*) t*h,x(i)
                end if
                if(i==50.and.t>353700) then
                    write(40,*) t*h,x(i)
                end if
            end do
            x=func_x
            y=func_y
            z=func_z
        end do
        return
    end subroutine fnf


    !建立网状结构

    !k=2
    subroutine neighbour_2(K)
        implicit none
        integer :: K !单节点K个邻居
        integer :: i,j,number
        allocate(neighbour_matrix_2(N,N))
        !初始化网络矩阵
        neighbour_matrix_2=0
        do i=1,N-1,1
            do j=i+1,N,1
                if(abs(i-j)<=K/2.or.abs(i-j)>=N-K/2) then
                    neighbour_matrix_2(i,j)=1
                    neighbour_matrix_2(j,i)=1
                end if
            end do
        end do
        !输出矩阵
        do i=1,N,1
            do j=1,N,1
                write(50,"(I2)",advance='no') neighbour_matrix_2(i,j)
            end do
            write(50,*)
        end do
        return
    end subroutine neighbour_2
    !k=2*p+2
    subroutine neighbour_2P(K)
        implicit none
        integer :: K !单节点K个邻居
        integer :: i,j,number
        allocate(neighbour_matrix_2P(N,N))
        !初始化网络矩阵
        neighbour_matrix_2P=0
        do i=1,N-1,1
            do j=i+1,N,1
                if((abs(i-j)<=K/2.or.abs(i-j)>=N-K/2).and.abs(i-j)/=1) then
                    neighbour_matrix_2P(i,j)=1
                    neighbour_matrix_2P(j,i)=1
                end if
            end do
        end do
        !输出矩阵
        do i=1,N,1
            do j=1,N,1
                write(60,"(I2)",advance='no') neighbour_matrix_2P(i,j)
            end do
            write(60,*)
        end do
        return
    end subroutine neighbour_2P

end module HR



program main
    use HR
    implicit none
    integer :: p
    p=40
    open(10,file="data_x0_y0_z0.txt")
    open(20,file="data_i_t_x.txt")
    open(30,file="data_t_x5.txt")
    open(40,file="data_t_x50.txt")
    open(50,file="neighbour_matrix_2.txt")
    open(60,file="neighbour_matrix_2P.txt")
    open(70,file="x.txt")
    call x0_y0_z0()
    call neighbour_2(2)
    call neighbour_2P(2*p+2)
    call fnf(p)
    deallocate(neighbour_matrix_2)
    deallocate(neighbour_matrix_2P)
    close(10)
    close(20)
    close(30)
    close(40)
    close(50)
    close(60)
    close(70)
end program main

figure 2

figure 2

module HR
    implicit none
    real,parameter :: h=0.005
    integer,parameter :: MaxT=800000
    real :: x,y,z
contains

    subroutine x0_y0_z0()
        implicit none
        call random_seed()
        call random_number(x)
        call random_number(y)
        call random_number(z)
    end subroutine x0_y0_z0

    subroutine fnf()
        implicit none
        integer :: t
        real :: func_x,func_y,func_z,a,b,c,d,r,s,x_R,II
        a=1.0
        b=3.0
        c=1.0
        d=5.0
        r=0.01
        s=5.0
        x_R=-1.6
        II=3.8
        do t=1,MaxT,1
            if(mod(t,5000)==0) then
                write(*,*) t
            end if
            func_x=x+h*(y-a*x**3+b*x**2-z+II)
            func_y=y+h*(c-d*x**2-y)
            func_z=z+h*(r*(s*(x-x_R)-z))
            x=func_x
            y=func_y
            z=func_z
            if(t>353700) then
                write(30,*) t*h,x
                write(*,*) t*h,x
            end if
        end do
        return
    end subroutine fnf

end module HR



program main
    use HR
    implicit none
    open(10,file="data_x0_y0_z0.txt")
    open(20,file="data_i_t_x.txt")
    open(30,file="data_t_x.txt")
    call x0_y0_z0()
    call fnf()
    close(10)
    close(20)
    close(30)
    close(40)
end program main

figure 2

figure 2

figure 3

figure 3

module HR
    implicit none
    real,parameter :: h=0.01
    integer,parameter :: N=100,MaxT=400000
    integer,allocatable :: neighbour_matrix_2(:,:),neighbour_matrix_2P(:,:)
    real :: x(N),y(N),z(N)
contains

    subroutine x0_y0_z0()
        implicit none
        integer :: k
        real :: x1,x2,x3
        call random_seed()
        do k=1,N,1
            call random_number(x1)
            call random_number(x2)
            call random_number(x3)
            x(k)=0.001*(k-N/2.0)+x1
            y(k)=0.002*(k-N/2.0)+x2
            z(k)=0.003*(k-N/2.0)+x3
            write(10,*) x(k),y(k),z(k)
        end do
    end subroutine x0_y0_z0

    subroutine fnf(p)
        implicit none
        integer :: i,j,p,t
        real :: func_x(N),func_y(N),func_z(N),a,b,c,d,k1,k2,II,JJ,CC,r,s,x_R,v_s,lambda,thet_s
        k1=3.0
        k2=2.0
        a=1.0
        b=3.0
        c=1.0
        d=5.0
        r=0.01
        s=5.0
        v_s=2.0
        x_R=-1.6
        lambda=10.0
        thet_s=-0.25
        II=3.8
        do t=1,MaxT,1
            if(mod(t,5000)==0) then
                write(*,*) t
            end if
            do i=1,N,1
                JJ=0.0
                CC=0.0
                do j=1,N,1
                    JJ=JJ+neighbour_matrix_2(i,j)*(x(j)-x(i))
                    CC=CC+neighbour_matrix_2P(i,j)*(1.0/(1.0+exp(-lambda*(x(j)-thet_s))))
                end do
                func_x(i)=x(i)+h*(y(i)-a*x(i)**3+b*x(i)**2-z(i)+II+k1*JJ+((-k2/(2.0*p-2.0))*(v_s-x(i))*CC))
                func_y(i)=y(i)+h*(c-d*x(i)**2-y(i))
                func_z(i)=z(i)+h*(r*(s*(x(i)-x_R)-z(i)))
                if(mod(t,50)==0.and.t>300000) then
                    write(20,*) i,t*h,x(i)
                end if
                if(i==1.and.t>350000) then
                    write(30,*) t*h,x(i)
                    write(*,*) t*h,x(i)
                end if
                if(t==350000) then
                    write(40,*) i,x(i)
                end if
                if(t==400000) then
                    write(50,*) i,x(i)
                end if
            end do
            x=func_x
            y=func_y
            z=func_z
        end do
        return
    end subroutine fnf


    !建立网状结构

    !k=2
    subroutine neighbour_2(K)
        implicit none
        integer :: K !单节点K个邻居
        integer :: i,j,number
        allocate(neighbour_matrix_2(N,N))
        !初始化网络矩阵
        neighbour_matrix_2=0
        do i=1,N-1,1
            do j=i+1,N,1
                if(abs(i-j)<=K/2.or.abs(i-j)>=N-K/2) then
                    neighbour_matrix_2(i,j)=1
                    neighbour_matrix_2(j,i)=1
                end if
            end do
        end do
        return
    end subroutine neighbour_2
    !k=2*p+2
    subroutine neighbour_2P(K)
        implicit none
        integer :: K !单节点K个邻居
        integer :: i,j,number
        allocate(neighbour_matrix_2P(N,N))
        !初始化网络矩阵
        neighbour_matrix_2P=0
        do i=1,N-1,1
            do j=i+1,N,1
                if((abs(i-j)<=K/2.or.abs(i-j)>=N-K/2).and.abs(i-j)/=1) then
                    neighbour_matrix_2P(i,j)=1
                    neighbour_matrix_2P(j,i)=1
                end if
            end do
        end do
        return
    end subroutine neighbour_2P

end module HR



program main
    use HR
    implicit none
    integer :: p
    p=40
    open(10,file="data_x0_y0_z0.txt")
    open(20,file="data_i_t_x.txt")
    open(30,file="data_t_x5.txt")
    open(40,file="i_x1.txt")
    open(50,file="i_x2.txt")
    call x0_y0_z0()
    call neighbour_2(2)
    call neighbour_2P(2*p+2)
    call fnf(p)
    deallocate(neighbour_matrix_2)
    deallocate(neighbour_matrix_2P)
    close(10)
    close(20)
    close(30)
    close(40)
    close(50)
end program main

figure 3_A

figure 3_A

figure 3_B

figure 3_B

figure 4

figure 4

module HR
    implicit none
    real,parameter :: h=0.01
    integer,parameter :: N=100,MaxT=400000
    integer,allocatable :: neighbour_matrix_2(:,:),neighbour_matrix_2P(:,:)
    real :: x(N),y(N),z(N)
contains

    subroutine x0_y0_z0()
        implicit none
        integer :: k
        real :: x1,x2,x3
        call random_seed()
        do k=1,N,1
            call random_number(x1)
            call random_number(x2)
            call random_number(x3)
            x(k)=0.001*(k-N/2.0)+x1
            y(k)=0.002*(k-N/2.0)+x2
            z(k)=0.003*(k-N/2.0)+x3
            write(10,*) x(k),y(k),z(k)
        end do
    end subroutine x0_y0_z0

    subroutine fnf(p)
        implicit none
        integer :: i,j,p,t
        real :: func_x(N),func_y(N),func_z(N),a,b,c,d,k1,k2,II,JJ,CC,r,s,x_R,v_s,lambda,thet_s
        k1=3.0
        k2=4.0
        a=1.0
        b=3.0
        c=1.0
        d=5.0
        r=0.01
        s=5.0
        v_s=2.0
        x_R=-1.6
        lambda=10.0
        thet_s=-0.25
        II=4.0
        do t=1,MaxT,1
            if(mod(t,5000)==0) then
                write(*,*) t
            end if
            do i=1,N,1
                JJ=0.0
                CC=0.0
                do j=1,N,1
                    JJ=JJ+neighbour_matrix_2(i,j)*(x(j)-x(i))
                    CC=CC+neighbour_matrix_2P(i,j)*(1.0/(1.0+exp(-lambda*(x(j)-thet_s))))
                end do
                func_x(i)=x(i)+h*(y(i)-a*x(i)**3+b*x(i)**2-z(i)+II+k1*JJ+((-k2/(2.0*p-2.0))*(v_s-x(i))*CC))
                func_y(i)=y(i)+h*(c-d*x(i)**2-y(i))
                func_z(i)=z(i)+h*(r*(s*(x(i)-x_R)-z(i)))
                if(mod(t,50)==0.and.t>300000) then
                    write(20,*) i,t*h,x(i)
                end if
                if(i==1.and.t>350000) then
                    write(30,*) t*h,x(i)
                    write(*,*) t*h,x(i)
                end if
                if(t==350000) then
                    write(40,*) i,x(i)
                end if
                if(t==400000) then
                    write(50,*) i,x(i)
                end if
            end do
            x=func_x
            y=func_y
            z=func_z
        end do
        return
    end subroutine fnf


    !建立网状结构

    !k=2
    subroutine neighbour_2(K)
        implicit none
        integer :: K !单节点K个邻居
        integer :: i,j,number
        allocate(neighbour_matrix_2(N,N))
        !初始化网络矩阵
        neighbour_matrix_2=0
        do i=1,N-1,1
            do j=i+1,N,1
                if(abs(i-j)<=K/2.or.abs(i-j)>=N-K/2) then
                    neighbour_matrix_2(i,j)=1
                    neighbour_matrix_2(j,i)=1
                end if
            end do
        end do
        return
    end subroutine neighbour_2
    !k=2*p+2
    subroutine neighbour_2P(K)
        implicit none
        integer :: K !单节点K个邻居
        integer :: i,j,number
        allocate(neighbour_matrix_2P(N,N))
        !初始化网络矩阵
        neighbour_matrix_2P=0
        do i=1,N-1,1
            do j=i+1,N,1
                if((abs(i-j)<=K/2.or.abs(i-j)>=N-K/2).and.abs(i-j)/=1) then
                    neighbour_matrix_2P(i,j)=1
                    neighbour_matrix_2P(j,i)=1
                end if
            end do
        end do
        return
    end subroutine neighbour_2P

end module HR



program main
    use HR
    implicit none
    integer :: p
    p=40
    open(10,file="data_x0_y0_z0.txt")
    open(20,file="data_i_t_x.txt")
    open(30,file="data_t_x5.txt")
    open(40,file="i_x1.txt")
    open(50,file="i_x2.txt")
    call x0_y0_z0()
    call neighbour_2(2)
    call neighbour_2P(2*p+2)
    call fnf(p)
    deallocate(neighbour_matrix_2)
    deallocate(neighbour_matrix_2P)
    close(10)
    close(20)
    close(30)
    close(40)
    close(50)
end program main

figure 4_A

figure 4_A

figure 4_B

figure 4_B

figure 4_C

figure 4_C

figure 5

figure 5

figure 5_A

figure 5_A

figure 5_B

figure 5_B


文章作者: rep-rebirth
版权声明: 本博客所有文章除特別声明外,均采用 CC BY 4.0 许可协议。转载请注明来源 rep-rebirth !
评论
评论
  目录