program sender_apple use mod_oasis use mpi implicit none integer :: i, kinfo, date integer :: comp_id, part_id, var_id integer :: n_points = 16 integer :: part_params(OASIS_Apple_Params), offset, local_size integer :: local_comm, local_comm_size, local_comm_rank integer :: icpl, coupl_comm, comm_size, comm_rank integer :: var_nodims(2) character(len=13) :: comp_name = "sender-apple" character(len=8) :: var_name = "FSENDOCN" real, allocatable :: field(:) print '(2A)', "Component name: ", comp_name call oasis_init_comp(comp_id, comp_name, kinfo) if(kinfo<0) call oasis_abort(comp_id, comp_name, & & "Error in oasis_init_comp: ", rcode=kinfo) print '(A,I0)', "Sender: Component ID: ", comp_id call oasis_get_localcomm(local_comm, kinfo) if(kinfo<0) call oasis_abort(comp_id, comp_name, & & "Error in oasis_get_localcomm: ", rcode=kinfo) call mpi_comm_size(local_comm, local_comm_size, kinfo) call mpi_comm_rank(local_comm, local_comm_rank, kinfo) print '(A,I0,A,I0)', "Sender: local_comm_rank = ",local_comm_rank, & & " of ",local_comm_size icpl = 1 if (local_comm_size .gt. 3) then if (local_comm_rank .ge. local_comm_size - 2) icpl = 0 end if call mpi_comm_split(local_comm, icpl, local_comm_rank, coupl_comm, kinfo) if(kinfo<0) call oasis_abort(comp_id, comp_name, & & "Error in mpi_comm_split: ", rcode=kinfo) if (icpl == 0) coupl_comm = MPI_COMM_NULL call oasis_set_couplcomm(coupl_comm, kinfo) if(kinfo<0) call oasis_abort(comp_id, comp_name, & & "Error in oasis_set_couplcomm: ", rcode=kinfo) if (icpl == 1) then call mpi_comm_size(coupl_comm, comm_size, kinfo) call mpi_comm_rank(coupl_comm, comm_rank, kinfo) print '(A,I0,A,I0)', "Sender: comm_rank = ",comm_rank, " of ",comm_size if ( mod(n_points,comm_size) /= 0) & & call oasis_abort(comp_id, comp_name, & & "Sender: comm_size has to divide n_points exaclty", rcode=kinfo) local_size=n_points/comm_size offset=comm_rank*local_size part_params(OASIS_Strategy) = OASIS_Apple part_params(OASIS_Offset) = offset part_params(OASIS_Length) = local_size call oasis_def_partition(part_id, part_params, kinfo) if(kinfo<0) call oasis_abort(comp_id, comp_name, & & "Error in oasis_def_partition: ", rcode=kinfo) print '(A,I0,A,I0)', "Sender rank(",comm_rank,"): part_id: ", part_id var_nodims=[1, 1] print '(A,I0,2A)', "Sender rank(",comm_rank,"): var_name: ", var_name call oasis_def_var(var_id, var_name, part_id, var_nodims, OASIS_OUT, & & [1], OASIS_REAL, kinfo) if(kinfo<0 .or. var_id<0) call oasis_abort(comp_id, comp_name, & & "Error in oasis_def_var: ", rcode=kinfo) print '(A,I0,A,I0)', "Sender rank(",comm_rank,"): var_id: ", var_id end if call oasis_enddef(kinfo) if(kinfo<0) call oasis_abort(comp_id, comp_name, & & "Error in oasis_enddef: ", rcode=kinfo) if (icpl == 1) then allocate(field(local_size)) field(:) = [(i+offset, i=1,local_size)] date=0 call oasis_put(var_id, date, field, kinfo) if(kinfo<0) call oasis_abort(comp_id, comp_name, & & "Error in oasis_put: ", rcode=kinfo) end if call oasis_terminate(kinfo) if(kinfo<0) call oasis_abort(comp_id, comp_name, & & "Error in oasis_terminate: ", rcode=kinfo) end program sender_apple