00001
00002
00003
00004
00005
00006
00007
00008 subroutine prismtrs_store_link_cnsrv(add1, add2, lcl_weights, id_action, id_epio_id, num_wts)
00009
00010
00011
00012
00013
00014
00015 USE PRISMDrv, dummy_INTERFACE => PRISMTrs_store_link_cnsrv
00016 IMPLICIT NONE
00017 integer, intent(in) ::
00018 add1,
00019 add2,
00020 id_action,
00021 id_epio_id,
00022 num_wts
00023
00024
00025 DOUBLE PRECISION, dimension(2*num_wts), intent(in) ::
00026 lcl_weights
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047
00048
00049
00050 CHARACTER(LEN=len_cvs_string), SAVE :: mycvs =
00051 '$Id: prismtrs_store_link_cnsrv.F90 2082 2009-10-21 13:31:19Z hanke $'
00052
00053 integer :: nlink, min_link, max_link
00054 integer :: il_num_links_map1
00055 integer, dimension(:,:), allocatable, save ::
00056 link_add1,
00057 link_add2
00058
00059 #ifdef DEBUGX
00060 PRINT *, '| | | | | | | Enter PRISMTrs_store_link_cnsrv'
00061 call psmile_flushstd
00062 #endif
00063
00064
00065
00066
00067
00068
00069
00070 #ifdef DEBUGX
00071 WRITE(*,*) 'id_action, lcl_weights =', id_action, lcl_weights(:)
00072 call psmile_flushstd
00073 #endif
00074 if (all(lcl_weights == zero)) return
00075
00076 IF (id_action == 2) THEN
00077 DEALLOCATE(link_add1, link_add2)
00078 RETURN
00079 ENDIF
00080
00081
00082
00083
00084
00085
00086 if (id_action == 0) then
00087
00088
00089
00090 allocate(link_add1(2,Drv_Epios(id_epio_id)%src_size), &
00091 link_add2(2,Drv_Epios(id_epio_id)%tgt_size))
00092 link_add1 = 0
00093 link_add2 = 0
00094
00095 min_link = 1
00096 max_link = 0
00097 else
00098 min_link = min(link_add1(1,add1),link_add2(1,add2))
00099 max_link = max(link_add1(2,add1),link_add2(2,add2))
00100 if (min_link == 0) then
00101 min_link = 1
00102 max_link = 0
00103 endif
00104 endif
00105
00106
00107
00108
00109
00110
00111
00112 do nlink=min_link,max_link
00113 if (add1 == Drv_Epios(id_epio_id)%grid1_add_map1(nlink)) then
00114 if (add2 == Drv_Epios(id_epio_id)%grid2_add_map1(nlink)) then
00115
00116 Drv_Epios(id_epio_id)%wts_map1(:,nlink) = &
00117 Drv_Epios(id_epio_id)%wts_map1(:,nlink) + lcl_weights(1:num_wts)
00118 return
00119
00120 endif
00121 endif
00122 end do
00123
00124
00125
00126
00127
00128
00129
00130
00131 il_num_links_map1 = Drv_Epios(id_epio_id)%num_links_map1 + 1
00132 if (il_num_links_map1 > Drv_Epios(id_epio_id)%max_links_map1) &
00133 call prismtrs_resize_remap_vars(Drv_Epios(id_epio_id)%resize_increment, id_epio_id, num_wts)
00134 Drv_Epios(id_epio_id)%grid1_add_map1(il_num_links_map1) = add1
00135 Drv_Epios(id_epio_id)%grid2_add_map1(il_num_links_map1) = add2
00136 Drv_Epios(id_epio_id)%wts_map1 (:,il_num_links_map1) = lcl_weights(1:num_wts)
00137
00138 if (link_add1(1,add1) == 0) link_add1(1,add1) = il_num_links_map1
00139 if (link_add2(1,add2) == 0) link_add2(1,add2) = il_num_links_map1
00140 link_add1(2,add1) = il_num_links_map1
00141 link_add2(2,add2) = il_num_links_map1
00142
00143 Drv_Epios(id_epio_id)%num_links_map1 = il_num_links_map1
00144
00145 #ifdef DEBUGX
00146 PRINT *, '| | | | | | | Quit PRISMTrs_store_link_cnsrv'
00147 call psmile_flushstd
00148 #endif
00149
00150
00151 end subroutine PRISMTrs_store_link_cnsrv
00152