00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011 subroutine psmile_store_send_info (var_id, id_trans_out, dir_index, &
00012 cpl_index, appl_index, ierror)
00013
00014
00015
00016 use PRISM_constants
00017
00018 use PSMILe, dummy_interface => PSMILe_store_send_info
00019
00020 implicit none
00021
00022
00023
00024 Integer, Intent (In) :: var_id
00025
00026
00027
00028 Integer, Intent (In) :: id_trans_out
00029
00030
00031
00032 Integer, Intent (In) :: cpl_index
00033
00034
00035
00036
00037 Integer, Intent (In) :: dir_index
00038
00039
00040
00041
00042 Integer, Intent (In) :: appl_index
00043
00044
00045
00046
00047
00048
00049
00050 integer, Intent (Out) :: ierror
00051
00052
00053
00054
00055
00056
00057
00058 Integer, Parameter :: send_coupler_index = 1
00059 Integer, Parameter :: send_direct_index = 2
00060 Integer, Parameter :: send_appl_index = 5
00061
00062
00063
00064 Type (GridFunction), Pointer :: field
00065 Type (Taskout_type), Pointer :: fieldout
00066
00067 Integer :: s_index
00068
00069 Integer :: task_id
00070
00071
00072
00073 Integer, parameter :: nerrp = 1
00074 Integer :: ierrp (nerrp)
00075
00076
00077
00078
00079
00080
00081
00082
00083
00084
00085
00086
00087
00088
00089
00090
00091
00092
00093 Character(len=len_cvs_string), save :: mycvs =
00094 '$Id: psmile_store_send_info.F90 2325 2010-04-21 15:00:07Z valcke $'
00095
00096
00097
00098
00099
00100 #ifdef VERBOSE
00101 print 9990, trim(ch_id), var_id
00102
00103 call psmile_flushstd
00104 #endif /* VERBOSE */
00105
00106 ierror = 0
00107
00108 field => Fields(var_id)
00109
00110 #ifdef PRISM_ASSERTION
00111
00112
00113
00114 if (var_id < 1 .or. var_id > Number_of_Fields_allocated) then
00115 print *, 'var_id', var_id
00116 call psmile_assert (__FILE__, __LINE__, "Incorrect var_id")
00117 endif
00118
00119 if ( field%status == PSMILe_Status_free .or. &
00120 field%status == PSMILe_Status_undefined ) then
00121 call psmile_assert (__FILE__, __LINE__, "Incorrect field")
00122 endif
00123 #endif
00124
00125
00126
00127 do task_id = 1, size(Fields(var_id)%Taskout)
00128 if ( Fields(var_id)%Taskout(task_id)%global_transi_id == id_trans_out ) exit
00129 enddo
00130
00131 if ( task_id > size(Fields(var_id)%Taskout) ) then
00132 ierror = PRISM_Error_Internal
00133 ierrp(1) = id_trans_out
00134 call psmile_error ( PRISM_Error_Arglist, &
00135 'No matching id was found id_trans_out ', &
00136 ierrp, 1, __FILE__, __LINE__ )
00137 return
00138 endif
00139
00140 fieldout => Fields(var_id)%Taskout(task_id)
00141
00142
00143
00144 if (cpl_index /= PRISM_UNDEFINED) then
00145 call psmile_get_exch_index (var_id, task_id, send_coupler_index, &
00146 s_index, ierror)
00147 if (ierror > 0) return
00148
00149 fieldout%send_coupler(s_index)%trans_out_id = id_trans_out
00150 fieldout%send_coupler(s_index)%send_info_index = cpl_index
00151 endif
00152
00153 if (dir_index /= PRISM_UNDEFINED) then
00154 call psmile_get_exch_index (var_id, task_id, send_direct_index, &
00155 s_index, ierror)
00156 if (ierror > 0) return
00157
00158 fieldout%send_direct(s_index)%trans_out_id = id_trans_out
00159 fieldout%send_direct(s_index)%send_info_index = dir_index
00160 endif
00161
00162 if (appl_index /= PRISM_UNDEFINED) then
00163 call psmile_get_exch_index (var_id, task_id, send_appl_index, &
00164 s_index, ierror)
00165 if (ierror > 0) return
00166
00167 fieldout%send_appl(s_index)%trans_out_id = id_trans_out
00168 fieldout%send_appl(s_index)%send_info_index = appl_index
00169 endif
00170
00171
00172
00173 #ifdef VERBOSE
00174 print 9980, trim(ch_id), ierror
00175
00176 call psmile_flushstd
00177 #endif /* VERBOSE */
00178
00179
00180
00181 #ifdef VERBOSE
00182
00183 9990 format (1x, a, ': psmile_store_send_info: var_id', i3)
00184 9980 format (1x, a, ': psmile_store_send_info: eof ierror =', i3)
00185
00186 #endif /* VERBOSE */
00187
00188 end subroutine PSMILe_store_send_info