FMS 2025.01-dev
Flexible Modeling System
Loading...
Searching...
No Matches
mpp_util_nocomm.inc
1! -*-f90-*-
2
3!***********************************************************************
4!* GNU Lesser General Public License
5!*
6!* This file is part of the GFDL Flexible Modeling System (FMS).
7!*
8!* FMS is free software: you can redistribute it and/or modify it under
9!* the terms of the GNU Lesser General Public License as published by
10!* the Free Software Foundation, either version 3 of the License, or (at
11!* your option) any later version.
12!*
13!* FMS is distributed in the hope that it will be useful, but WITHOUT
14!* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
15!* FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
16!* for more details.
17!*
18!* You should have received a copy of the GNU Lesser General Public
19!* License along with FMS. If not, see <http://www.gnu.org/licenses/>.
20!***********************************************************************
21!> @file
22!> @brief Utility routines for parallelization, non-mpi version
23
24!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
25! !
26! MISCELLANEOUS UTILITIES: mpp_error !
27! !
28!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
29
30subroutine mpp_error_basic( errortype, errormsg )
31 !a very basic error handler
32 !uses ABORT and FLUSH calls, may need to use cpp to rename
33 integer, intent(in) :: errortype
34 character(len=*), intent(in), optional :: errormsg
35 character(len=512) :: text
36 logical :: opened
37 integer :: istat, errunit, outunit
38
39 if( .NOT.module_is_initialized )call abort()
40
41 select case( errortype )
42 case(note)
43 text = 'NOTE' !just FYI
44 case(warning)
45 text = 'WARNING' !probable error
46 case(fatal)
47 text = 'FATAL' !fatal error
48 case default
49 text = 'WARNING: non-existent errortype (must be NOTE|WARNING|FATAL)'
50 end select
51
52 if( npes.GT.1 )write( text,'(a,i5)' )trim(text)//' from PE', pe !this is the mpp part
53 if( PRESENT(errormsg) )text = trim(text)//': '//trim(errormsg)
54
55 errunit = stderr()
56 outunit = stdout()
57
58 select case( errortype )
59 case(note)
60 write( outunit,'(a)' )trim(text)
61 case default
62 write( errunit,'(/a/)' )trim(text)
63 write( outunit,'(/a/)' )trim(text)
64 if( errortype.EQ.fatal .OR. warnings_are_fatal )then
65 FLUSH(outunit)
66 call abort() !automatically calls traceback on Cray systems
67 end if
68 end select
69
70 error_state = errortype
71 return
72end subroutine mpp_error_basic
73
74!#####################################################################
75!> Makes a PE set out of a PE list. A PE list is an ordered list of PEs
76!! a PE set is a triad (start,log2stride,size) for SHMEM, an a communicator for MPI
77!! if stride is non-uniform or not a power of 2,
78!! will return error (not required for MPI but enforced for uniformity)
79function get_peset(pelist)
80 integer :: get_peset
81 integer, intent(in), optional :: pelist(:)
82
83 if( .NOT.PRESENT(pelist) )then !set it to current_peset_num
84 get_peset = current_peset_num; return
85 end if
86
87 get_peset = 0
88
89 return
90
91end function get_peset
92
93!#######################################################################
94!> synchronize PEs in list
95subroutine mpp_sync( pelist, check )
96 integer, intent(in), optional :: pelist(:)
97 integer, intent(in), optional :: check
98
99 return
100end subroutine mpp_sync
101
102!#######################################################################
103!> This is to check if current PE's outstanding puts are complete
104!! but we can't use shmem_fence because we are actually waiting for
105!! a remote PE to complete its get
106subroutine mpp_sync_self( pelist, check, request, msg_size, msg_type )
107 integer, intent(in), optional :: pelist(:)
108 integer, intent(in), optional :: check
109 integer, intent(inout), optional :: request(:)
110 integer, intent(in ), optional :: msg_size(:)
111 integer, intent(in ), optional :: msg_type(:)
112
113
114 return
115end subroutine mpp_sync_self
subroutine mpp_sync_self(pelist, check, request, msg_size, msg_type)
This is to check if current PE's outstanding puts are complete but we can't use shmem_fence because w...
subroutine mpp_error_basic(errortype, errormsg)
A very basic error handler uses ABORT and FLUSH calls, may need to use cpp to rename.
integer function get_peset(pelist)
Makes a PE set out of a PE list. A PE list is an ordered list of PEs a PE set is a triad (start,...
subroutine mpp_sync(pelist, do_self)
Synchronize PEs in list.