psmile_get_method_handle.F90
Go to the documentation of this file.00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011 subroutine psmile_get_method_handle ( grid_id, method_id, ierror)
00012
00013
00014
00015 use PRISM_constants
00016
00017 use PSMILe, dummy_interface => PSMILe_Get_method_handle
00018
00019 implicit none
00020
00021
00022
00023 integer, Intent (In) :: grid_id
00024
00025
00026
00027
00028
00029
00030 integer, Intent (Out) :: method_id
00031
00032
00033
00034 integer, Intent (Out) :: ierror
00035
00036
00037
00038
00039
00040
00041
00042
00043 integer :: i, new_dim
00044
00045 integer, parameter :: nerrp = 2
00046 integer :: ierrp (nerrp)
00047
00048 type (Method), Dimension(:), Pointer :: New_Methods
00049
00050
00051
00052
00053
00054
00055
00056
00057
00058
00059
00060
00061
00062
00063
00064
00065
00066
00067
00068
00069 Character(len=len_cvs_string), save :: mycvs =
00070 '$Id: psmile_get_method_handle.F90 2325 2010-04-21 15:00:07Z valcke $'
00071
00072
00073
00074 #ifdef VERBOSE
00075 print *, trim(ch_id), ': PSMILe_Get_method_handle: start'
00076
00077 call psmile_flushstd
00078 #endif /* VERBOSE */
00079
00080
00081
00082 ierror = 0
00083
00084
00085
00086 do i = 1, Number_of_Methods_allocated
00087 if (Methods(i)%status == PSMILe_status_free) exit
00088 end do
00089
00090 if (i <= Number_of_Methods_allocated) then
00091
00092 method_id = i
00093
00094 else
00095
00096
00097
00098 new_dim = Number_of_Methods_allocated + 8
00099
00100 Allocate (New_Methods (new_dim), STAT = ierror)
00101 if (ierror > 0) then
00102 ierrp (1) = ierror
00103 ierrp (2) = new_dim
00104 ierror = PRISM_Error_Alloc
00105
00106 call psmile_error ( ierror, 'New_Methods', &
00107 ierrp, 2, __FILE__, __LINE__ )
00108 return
00109 endif
00110
00111 New_Methods (1:Number_of_Methods_allocated) = &
00112 Methods (1:Number_of_Methods_allocated)
00113
00114 New_Methods(Number_of_Methods_allocated+1:new_dim)%status = &
00115 PSMILe_status_free
00116
00117 Nullify ( New_Methods(method_id)%halo_pointer )
00118
00119 do i = Number_of_Methods_allocated+1, new_dim
00120 Nullify ( New_Methods(i)%send_infos_direct )
00121 Nullify ( New_Methods(i)%send_infos_coupler )
00122 Nullify ( New_Methods(i)%recv_infos_direct )
00123 Nullify ( New_Methods(i)%recv_infos_coupler )
00124 Nullify ( New_Methods(i)%coords_pointer )
00125 Nullify ( New_Methods(i)%subgrid_pointer )
00126 Nullify ( New_Methods(i)%vector_pointer )
00127 Nullify ( New_Methods(i)%gauss2_real(1)%vector)
00128 Nullify ( New_Methods(i)%gauss2_real(2)%vector)
00129 Nullify ( New_Methods(i)%gauss2_dble(1)%vector)
00130 Nullify ( New_Methods(i)%gauss2_dble(2)%vector)
00131 #if defined ( PRISM_QUAD_TYPE )
00132 Nullify ( New_Methods(i)%gauss2_quad(1)%vector)
00133 Nullify ( New_Methods(i)%gauss2_quad(2)%vector)
00134 #endif
00135 enddo
00136
00137
00138
00139 Deallocate (Methods, STAT = ierror)
00140 if (ierror > 0) then
00141 ierrp (1) = ierror
00142 ierror = PRISM_Error_Dealloc
00143
00144 call psmile_error ( ierror, 'Methods', &
00145 ierrp, 1, __FILE__, __LINE__ )
00146 return
00147 endif
00148
00149
00150
00151 Methods => New_Methods
00152
00153 method_id = Number_of_Methods_allocated + 1
00154
00155 Number_of_Methods_allocated = new_dim
00156 endif
00157
00158
00159
00160 Methods(method_id)%next_method_in_grid = PRISM_UNDEFINED
00161
00162 do i = 1, Number_of_Methods_allocated
00163 if (Methods(i)%grid_id == grid_id .and. &
00164 Methods(i)%next_method_in_grid == PRISM_UNDEFINED .and. &
00165 Methods(i)%status /= PSMILe_status_free) exit
00166 end do
00167
00168 if (i > Number_of_Methods_allocated) then
00169
00170
00171
00172 Methods(method_id)%previous_method_in_grid = PRISM_UNDEFINED
00173 else
00174
00175
00176
00177 Methods(i)%next_method_in_grid = method_id
00178 endif
00179
00180
00181
00182 Methods(method_id)%status = PSMILe_status_defined
00183 Methods(method_id)%grid_id = grid_id
00184
00185
00186
00187 Methods(method_id)%used_for_coupling = .false.
00188
00189
00190
00191 Methods(method_id)%n_send_info_direct = 0
00192 Methods(method_id)%n_send_info_coupler = 0
00193 Methods(method_id)%n_recv_info_direct = 0
00194 Methods(method_id)%n_recv_info_coupler = 0
00195 Methods(method_id)%n_send_info_appl = 0
00196
00197 Methods(method_id)%n_alloc_send_direct = 0
00198 Methods(method_id)%n_alloc_send_coupler = 0
00199 Methods(method_id)%n_alloc_recv_direct = 0
00200 Methods(method_id)%n_alloc_recv_coupler = 0
00201 Methods(method_id)%n_alloc_send_appl = 0
00202
00203 #ifdef VERBOSE
00204 print *, trim(ch_id), ': PSMILe_Get_method_handle: eof handle, ierror =', &
00205 method_id, ierror
00206
00207 call psmile_flushstd
00208 #endif /* VERBOSE */
00209
00210 end subroutine PSMILe_Get_method_handle