00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011 subroutine prismtrs_apply_weights(il_src_size, &
00012 dda_trans_out, &
00013 il_tgt_size, &
00014 dda_trans_in, &
00015 ida_mask, &
00016 id_nb_neighbors, &
00017 ida_neighbors, &
00018 dda_weights, &
00019 nbr_fields, &
00020 id_err)
00021
00022
00023
00024 USE PRISMDrv, dummy_interface => PRISMTrs_Apply_weights
00025
00026 IMPLICIT NONE
00027
00028
00029
00030
00031
00032 INTEGER, INTENT (IN) :: il_src_size
00033
00034
00035 INTEGER, INTENT (IN) :: nbr_fields
00036
00037
00038 DOUBLE PRECISION, DIMENSION(il_src_size*nbr_fields), INTENT(IN) :: dda_trans_out
00039
00040
00041 INTEGER, INTENT (IN) :: il_tgt_size
00042
00043
00044 INTEGER, DIMENSION(il_tgt_size), INTENT(IN) :: ida_mask
00045
00046
00047
00048 INTEGER, INTENT (IN) :: id_nb_neighbors
00049
00050
00051 INTEGER, DIMENSION(il_tgt_size,id_nb_neighbors), INTENT (IN) ::
00052 ida_neighbors
00053 DOUBLE PRECISION, DIMENSION(il_tgt_size,id_nb_neighbors), INTENT(IN) :: dda_weights
00054
00055
00056
00057
00058
00059 DOUBLE PRECISION, DIMENSION(il_tgt_size*nbr_fields), INTENT (Out) :: dda_trans_in
00060
00061 INTEGER, INTENT (Out) :: id_err
00062
00063
00064
00065
00066
00067
00068
00069
00070
00071
00072
00073
00074
00075
00076
00077
00078
00079
00080
00081 CHARACTER(LEN=len_cvs_string), SAVE :: mycvs =
00082 '$Id: prismtrs_apply_weights.F90 2963 2011-02-17 14:52:53Z coquart $'
00083
00084 INTEGER :: ib, ib_bis, i
00085 INTEGER :: il_count_neighbors
00086 #ifdef SX
00087 INTEGER,DIMENSION (il_tgt_size) :: il_count_neighborsh
00088 #endif
00089
00090
00091
00092 #ifdef VERBOSE
00093 PRINT *, '| | | | | | | Enter PRISMTrs_Apply_weights'
00094 call psmile_flushstd
00095 #endif
00096
00097 id_err = 0
00098
00099
00100
00101
00102 dda_trans_in(:) = 0.
00103 #ifndef SX
00104 DO i = 1, nbr_fields
00105
00106 DO ib = 1, il_tgt_size
00107
00108 il_count_neighbors = 0
00109
00110 DO ib_bis = 1, id_nb_neighbors
00111 IF (ida_neighbors(ib,ib_bis) /= 0) THEN
00112 dda_trans_in(ib+(i-1)*il_tgt_size) = dda_trans_in(ib+(i-1)*il_tgt_size) + &
00113 dda_weights(ib,ib_bis) * &
00114 dda_trans_out(ida_neighbors(ib,ib_bis)+(i-1)*il_src_size)
00115 #ifdef DEBUGX
00116 IF ( (ib == 30) .OR. (ib == 1667) ) THEN
00117 PRINT*, '+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++'
00118 PRINT*, 'weights(ib,ib_bis) for EPIOT number ib, ib_bis :',ib,ib_bis
00119 PRINT*, '+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++'
00120 PRINT*, dda_weights(ib,ib_bis)
00121 ENDIF
00122 #endif
00123 ELSE
00124 il_count_neighbors = il_count_neighbors + 1
00125 END IF
00126
00127 END DO
00128
00129 IF (il_count_neighbors == id_nb_neighbors) &
00130 dda_trans_in(ib+(i-1)*il_tgt_size) = PSMILe_dundef
00131 #ifdef DEBUGX
00132 IF ( (ib == 105) .OR. (ib == 315) ) THEN
00133 PRINT*, '+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++'
00134 PRINT*, 'dda_trans_in for EPIOT number ib :',ib
00135 PRINT*, 'ida_neighbors(ib,:)'
00136 PRINT*, '+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++'
00137 PRINT*, dda_trans_in(ib+(i-1)*il_tgt_size)
00138 PRINT*, ida_neighbors(ib,:)
00139 ENDIF
00140 #endif
00141
00142 END DO
00143
00144 END DO
00145
00146 #else
00147
00148 il_count_neighborsh=0
00149 DO ib = 1, il_tgt_size
00150 DO ib_bis = 1, id_nb_neighbors
00151 IF (ida_mask(ib) == 1 .AND. ida_neighbors(ib,ib_bis) == 0) THEN
00152 il_count_neighborsh(ib) = il_count_neighborsh(ib) + 1
00153 ENDIF
00154 ENDDO
00155 ENDDO
00156
00157 DO i = 1, nbr_fields
00158 DO ib_bis = 1, id_nb_neighbors
00159 DO ib = 1, il_tgt_size
00160
00161 IF (ida_mask(ib) == 1) THEN
00162
00163 IF (ida_neighbors(ib,ib_bis) /= 0) THEN
00164 dda_trans_in(ib+(i-1)*il_tgt_size) = dda_trans_in(ib+(i-1)*il_tgt_size) + &
00165 dda_weights(ib,ib_bis) * &
00166 dda_trans_out(ida_neighbors(ib,ib_bis)+(i-1)*il_src_size)
00167 END IF
00168
00169 ELSE
00170
00171
00172
00173
00174
00175 dda_trans_in(ib+(i-1)*il_tgt_size) = PSMILe_dundef
00176
00177 END IF
00178 END DO
00179 END DO
00180
00181 DO ib = 1, il_tgt_size
00182 IF (il_count_neighborsh(ib) == id_nb_neighbors) &
00183 dda_trans_in(ib+(i-1)*il_tgt_size) = PSMILe_dundef
00184 END DO
00185
00186 ENDDO
00187
00188 #endif
00189
00190 #ifdef VERBOSE
00191 PRINT *, '| | | | | | | Quit PRISMTrs_Apply_weights'
00192 call psmile_flushstd
00193 #endif
00194
00195 END SUBROUTINE PRISMTrs_Apply_weights