[6133] | 1 | ! |
---|
[6928] | 2 | ! $Id$ |
---|
[6133] | 3 | ! |
---|
| 4 | ! AGRIF (Adaptive Grid Refinement In Fortran) |
---|
| 5 | ! |
---|
| 6 | ! Copyright (C) 2003 Laurent Debreu (Laurent.Debreu@imag.fr) |
---|
| 7 | ! Christophe Vouland (Christophe.Vouland@imag.fr) |
---|
| 8 | ! |
---|
| 9 | ! This program is free software; you can redistribute it and/or modify |
---|
| 10 | ! it under the terms of the GNU General Public License as published by |
---|
| 11 | ! the Free Software Foundation; either version 2 of the License, or |
---|
| 12 | ! (at your option) any later version. |
---|
| 13 | ! |
---|
| 14 | ! This program is distributed in the hope that it will be useful, |
---|
| 15 | ! but WITHOUT ANY WARRANTY; without even the implied warranty of |
---|
| 16 | ! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the |
---|
| 17 | ! GNU General Public License for more details. |
---|
| 18 | ! |
---|
| 19 | ! You should have received a copy of the GNU General Public License |
---|
| 20 | ! along with this program; if not, write to the Free Software |
---|
| 21 | ! Foundation, Inc., 59 Temple Place - Suite 330, Boston, MA 02111-1307, USA. |
---|
| 22 | ! |
---|
| 23 | ! |
---|
| 24 | ! |
---|
| 25 | !> Module containing different procedures of update (copy, average, full_weighting) |
---|
| 26 | !! used in the #Agrif_Update module. |
---|
| 27 | !=================================================================================================== |
---|
| 28 | ! |
---|
| 29 | module Agrif_UpdateBasic |
---|
| 30 | ! |
---|
| 31 | use Agrif_Types |
---|
| 32 | |
---|
| 33 | implicit none |
---|
| 34 | |
---|
| 35 | integer, dimension(:,:), allocatable :: indchildcopy |
---|
| 36 | integer, dimension(:,:), allocatable :: indchildaverage |
---|
| 37 | ! |
---|
| 38 | contains |
---|
| 39 | ! |
---|
| 40 | !=================================================================================================== |
---|
| 41 | ! subroutine Agrif_basicupdate_copy1d |
---|
| 42 | ! |
---|
| 43 | !> Carries out a copy on a parent grid (vector x) from its child grid (vector y). |
---|
| 44 | !--------------------------------------------------------------------------------------------------- |
---|
| 45 | subroutine Agrif_basicupdate_copy1d ( x, y, np, nc, s_parent, s_child, ds_parent, ds_child ) |
---|
| 46 | !--------------------------------------------------------------------------------------------------- |
---|
| 47 | real, dimension(np), intent(out) :: x !< Coarse output data to parent |
---|
| 48 | real, dimension(nc), intent(in) :: y !< Fine input data from child |
---|
| 49 | integer, intent(in) :: np !< Length of parent array |
---|
| 50 | integer, intent(in) :: nc !< Length of child array |
---|
| 51 | real, intent(in) :: s_parent !< Parent grid position (s_root = 0) |
---|
| 52 | real, intent(in) :: s_child !< Child grid position (s_root = 0) |
---|
| 53 | real, intent(in) :: ds_parent !< Parent grid dx (ds_root = 1) |
---|
| 54 | real, intent(in) :: ds_child !< Child grid dx (ds_root = 1) |
---|
| 55 | !--------------------------------------------------------------------------------------------------- |
---|
| 56 | integer :: i, locind_child_left, coeffraf |
---|
| 57 | ! |
---|
| 58 | coeffraf = nint(ds_parent/ds_child) |
---|
| 59 | locind_child_left = 1 + nint((s_parent - s_child)/ds_child) |
---|
| 60 | ! |
---|
| 61 | if ( coeffraf == 1 ) then |
---|
| 62 | !CDIR ALTCODE |
---|
| 63 | x(1:np) = y(locind_child_left:locind_child_left+np-1) |
---|
| 64 | return |
---|
| 65 | endif |
---|
| 66 | ! |
---|
| 67 | !CDIR ALTCODE |
---|
| 68 | do i = 1,np |
---|
| 69 | x(i) = y(locind_child_left) |
---|
| 70 | locind_child_left = locind_child_left + coeffraf |
---|
| 71 | enddo |
---|
| 72 | !--------------------------------------------------------------------------------------------------- |
---|
| 73 | end subroutine Agrif_basicupdate_copy1d |
---|
| 74 | !=================================================================================================== |
---|
| 75 | ! |
---|
| 76 | !=================================================================================================== |
---|
| 77 | ! subroutine Agrif_basicupdate_copy1d_before |
---|
| 78 | ! |
---|
| 79 | !> Precomputes index for a copy on a parent grid (vector x) from its child grid (vector y). |
---|
| 80 | !--------------------------------------------------------------------------------------------------- |
---|
| 81 | subroutine Agrif_basicupdate_copy1d_before ( nc2, np, nc, s_parent, s_child, ds_parent, ds_child, dir ) |
---|
| 82 | !--------------------------------------------------------------------------------------------------- |
---|
| 83 | integer, intent(in) :: nc2 !< Length of parent array |
---|
| 84 | integer, intent(in) :: np !< Length of parent array |
---|
| 85 | integer, intent(in) :: nc !< Length of child array |
---|
| 86 | real, intent(in) :: s_parent !< Parent grid position (s_root = 0) |
---|
| 87 | real, intent(in) :: s_child !< Child grid position (s_root = 0) |
---|
| 88 | real, intent(in) :: ds_parent !< Parent grid dx (ds_root = 1) |
---|
| 89 | real, intent(in) :: ds_child !< Child grid dx (ds_root = 1) |
---|
| 90 | integer, intent(in) :: dir !< Direction |
---|
| 91 | !--------------------------------------------------------------------------------------------------- |
---|
| 92 | integer, dimension(:,:), allocatable :: indchildcopy_tmp |
---|
| 93 | integer :: i, n_old, locind_child_left, coeffraf |
---|
| 94 | ! |
---|
| 95 | coeffraf = nint(ds_parent/ds_child) |
---|
| 96 | ! |
---|
| 97 | locind_child_left = 1 + nint((s_parent - s_child)/ds_child) |
---|
| 98 | |
---|
| 99 | if ( .not.allocated(indchildcopy) ) then |
---|
| 100 | allocate(indchildcopy(np*nc2, 3)) |
---|
| 101 | else |
---|
| 102 | n_old = size(indchildcopy,1) |
---|
| 103 | if ( n_old < np*nc2 ) then |
---|
| 104 | allocate( indchildcopy_tmp(1:n_old, 3)) |
---|
| 105 | indchildcopy_tmp = indchildcopy |
---|
| 106 | deallocate(indchildcopy) |
---|
| 107 | allocate(indchildcopy(np*nc2, 3)) |
---|
| 108 | indchildcopy(1:n_old,:) = indchildcopy_tmp |
---|
| 109 | deallocate(indchildcopy_tmp) |
---|
| 110 | endif |
---|
| 111 | endif |
---|
| 112 | ! |
---|
| 113 | do i = 1,np |
---|
| 114 | indchildcopy(i,dir) = locind_child_left |
---|
| 115 | locind_child_left = locind_child_left + coeffraf |
---|
| 116 | enddo |
---|
| 117 | ! |
---|
| 118 | do i = 2,nc2 |
---|
| 119 | indchildcopy(1+(i-1)*np:i*np,dir) = indchildcopy(1+(i-2)*np:(i-1)*np,dir) + nc |
---|
| 120 | enddo |
---|
| 121 | !--------------------------------------------------------------------------------------------------- |
---|
| 122 | end subroutine Agrif_basicupdate_copy1d_before |
---|
| 123 | !=================================================================================================== |
---|
| 124 | ! |
---|
| 125 | !=================================================================================================== |
---|
| 126 | ! subroutine Agrif_basicupdate_copy1d_after |
---|
| 127 | ! |
---|
| 128 | !> Carries out a copy on a parent grid (vector x) from its child grid (vector y) |
---|
| 129 | !! using precomputed index. |
---|
| 130 | !--------------------------------------------------------------------------------------------------- |
---|
| 131 | subroutine Agrif_basicupdate_copy1d_after ( x, y, np, nc, dir ) |
---|
| 132 | !--------------------------------------------------------------------------------------------------- |
---|
| 133 | real, dimension(np), intent(out) :: x !< Coarse output data to parent |
---|
| 134 | real, dimension(nc), intent(in) :: y !< Fine input data from child |
---|
| 135 | integer, intent(in) :: np !< Length of parent array |
---|
| 136 | integer, intent(in) :: nc !< Length of child array |
---|
| 137 | integer, intent(in) :: dir !< Direction |
---|
| 138 | !--------------------------------------------------------------------------------------------------- |
---|
| 139 | integer :: i |
---|
| 140 | ! |
---|
| 141 | !CDIR ALTCODE |
---|
| 142 | do i = 1,np |
---|
| 143 | x(i) = y(indchildcopy(i,dir)) |
---|
| 144 | enddo |
---|
| 145 | !--------------------------------------------------------------------------------------------------- |
---|
| 146 | end subroutine Agrif_basicupdate_copy1d_after |
---|
| 147 | !=================================================================================================== |
---|
| 148 | ! |
---|
| 149 | !=================================================================================================== |
---|
| 150 | ! subroutine Agrif_basicupdate_average1d |
---|
| 151 | ! |
---|
| 152 | !> Carries out an update by average on a parent grid (vector x)from its child grid (vector y). |
---|
| 153 | !--------------------------------------------------------------------------------------------------- |
---|
| 154 | subroutine Agrif_basicupdate_average1d ( x, y, np, nc, s_parent, s_child, ds_parent, ds_child ) |
---|
| 155 | !--------------------------------------------------------------------------------------------------- |
---|
| 156 | REAL, DIMENSION(np), intent(out) :: x |
---|
| 157 | REAL, DIMENSION(nc), intent(in) :: y |
---|
| 158 | INTEGER, intent(in) :: np,nc |
---|
| 159 | REAL, intent(in) :: s_parent, s_child |
---|
| 160 | REAL, intent(in) :: ds_parent, ds_child |
---|
| 161 | ! |
---|
| 162 | INTEGER :: i, ii, locind_child_left, coeffraf |
---|
| 163 | REAL :: xpos, invcoeffraf |
---|
| 164 | INTEGER :: nbnonnuls |
---|
| 165 | INTEGER :: diffmod |
---|
| 166 | ! |
---|
| 167 | coeffraf = nint(ds_parent/ds_child) |
---|
| 168 | invcoeffraf = 1./coeffraf |
---|
| 169 | ! |
---|
| 170 | if (coeffraf == 1) then |
---|
| 171 | locind_child_left = 1 + nint((s_parent - s_child)/ds_child) |
---|
| 172 | x(1:np) = y(locind_child_left:locind_child_left+np-1) |
---|
| 173 | return |
---|
| 174 | endif |
---|
| 175 | ! |
---|
| 176 | xpos = s_parent |
---|
| 177 | x = 0. |
---|
| 178 | ! |
---|
| 179 | diffmod = 0 |
---|
| 180 | ! |
---|
| 181 | IF ( mod(coeffraf,2) == 0 ) diffmod = 1 |
---|
| 182 | ! |
---|
| 183 | locind_child_left = 1 + agrif_int((xpos - s_child)/ds_child) |
---|
| 184 | ! |
---|
| 185 | IF (Agrif_UseSpecialValueInUpdate) THEN |
---|
| 186 | do i = 1,np |
---|
| 187 | nbnonnuls = 0 |
---|
| 188 | !CDIR NOVECTOR |
---|
| 189 | do ii = -coeffraf/2+locind_child_left+diffmod, & |
---|
| 190 | coeffraf/2+locind_child_left |
---|
| 191 | IF (y(ii) /= Agrif_SpecialValueFineGrid) THEN |
---|
| 192 | ! nbnonnuls = nbnonnuls + 1 |
---|
| 193 | x(i) = x(i) + y(ii) |
---|
| 194 | ENDIF |
---|
| 195 | enddo |
---|
| 196 | ! IF (nbnonnuls /= 0) THEN |
---|
| 197 | ! x(i) = x(i)/nbnonnuls |
---|
| 198 | ! ELSE |
---|
| 199 | ! x(i) = Agrif_SpecialValueFineGrid |
---|
| 200 | ! ENDIF |
---|
| 201 | locind_child_left = locind_child_left + coeffraf |
---|
| 202 | enddo |
---|
| 203 | ELSE |
---|
| 204 | ! |
---|
| 205 | !CDIR ALTCODE |
---|
| 206 | do i = 1,np |
---|
| 207 | !CDIR NOVECTOR |
---|
| 208 | do ii = -coeffraf/2+locind_child_left+diffmod, & |
---|
| 209 | coeffraf/2+locind_child_left |
---|
| 210 | x(i) = x(i) + y(ii) |
---|
| 211 | enddo |
---|
| 212 | ! x(i) = x(i)*invcoeffraf |
---|
| 213 | locind_child_left = locind_child_left + coeffraf |
---|
| 214 | enddo |
---|
| 215 | IF (.not.Agrif_Update_Weights) THEN |
---|
| 216 | x = x * invcoeffraf |
---|
| 217 | ENDIF |
---|
| 218 | ENDIF |
---|
| 219 | !--------------------------------------------------------------------------------------------------- |
---|
| 220 | end subroutine Agrif_basicupdate_average1d |
---|
| 221 | !=================================================================================================== |
---|
| 222 | ! |
---|
| 223 | !=================================================================================================== |
---|
| 224 | ! subroutine Average1dPrecompute |
---|
| 225 | ! |
---|
| 226 | !> Carries out an update by average on a parent grid (vector x)from its child grid (vector y). |
---|
| 227 | !--------------------------------------------------------------------------------------------------- |
---|
| 228 | subroutine Average1dPrecompute ( nc2, np, nc, s_parent, s_child, ds_parent, ds_child, dir ) |
---|
| 229 | !--------------------------------------------------------------------------------------------------- |
---|
| 230 | INTEGER, intent(in) :: nc2, np, nc |
---|
| 231 | REAL, intent(in) :: s_parent, s_child |
---|
| 232 | REAL, intent(in) :: ds_parent, ds_child |
---|
| 233 | INTEGER, intent(in) :: dir |
---|
| 234 | ! |
---|
| 235 | INTEGER, DIMENSION(:,:), ALLOCATABLE :: indchildaverage_tmp |
---|
| 236 | INTEGER :: i, locind_child_left, coeffraf |
---|
| 237 | REAL :: xpos |
---|
| 238 | INTEGER :: diffmod |
---|
| 239 | ! |
---|
| 240 | coeffraf = nint(ds_parent/ds_child) |
---|
| 241 | xpos = s_parent |
---|
| 242 | diffmod = 0 |
---|
| 243 | ! |
---|
| 244 | IF ( mod(coeffraf,2) == 0 ) diffmod = 1 |
---|
| 245 | ! |
---|
| 246 | locind_child_left = 1 + agrif_int((xpos - s_child)/ds_child) |
---|
| 247 | ! |
---|
| 248 | if (.not.allocated(indchildaverage)) then |
---|
| 249 | allocate(indchildaverage(np*nc2,3)) |
---|
| 250 | else |
---|
| 251 | if (size(indchildaverage,1)<np*nc2) then |
---|
| 252 | allocate( indchildaverage_tmp(size(indchildaverage,1),size(indchildaverage,2))) |
---|
| 253 | indchildaverage_tmp = indchildaverage |
---|
| 254 | deallocate(indchildaverage) |
---|
| 255 | allocate(indchildaverage(np*nc2,3)) |
---|
| 256 | indchildaverage(1:size(indchildaverage_tmp,1),1:size(indchildaverage_tmp,2)) = indchildaverage_tmp |
---|
| 257 | deallocate(indchildaverage_tmp) |
---|
| 258 | endif |
---|
| 259 | endif |
---|
| 260 | ! |
---|
| 261 | do i = 1,np |
---|
| 262 | indchildaverage(i,dir)= -coeffraf/2+locind_child_left+diffmod |
---|
| 263 | locind_child_left = locind_child_left + coeffraf |
---|
| 264 | enddo |
---|
| 265 | ! |
---|
| 266 | do i = 2,nc2 |
---|
| 267 | indchildaverage(1+(i-1)*np:i*np,dir) = indchildaverage(1+(i-2)*np:(i-1)*np,dir) + nc |
---|
| 268 | enddo |
---|
| 269 | !--------------------------------------------------------------------------------------------------- |
---|
| 270 | end subroutine Average1dPrecompute |
---|
| 271 | !=================================================================================================== |
---|
| 272 | ! |
---|
| 273 | !=================================================================================================== |
---|
| 274 | ! subroutine Average1dAfterCompute |
---|
| 275 | ! |
---|
| 276 | !> Carries out an update by average on a parent grid (vector x) from its child grid (vector y). |
---|
| 277 | !--------------------------------------------------------------------------------------------------- |
---|
| 278 | subroutine Average1dAfterCompute ( x, y, np, nc, s_parent, s_child, ds_parent, ds_child, dir ) |
---|
| 279 | !--------------------------------------------------------------------------------------------------- |
---|
| 280 | REAL, DIMENSION(np), intent(inout) :: x |
---|
| 281 | REAL, DIMENSION(nc), intent(in) :: y |
---|
| 282 | INTEGER, intent(in) :: np, nc |
---|
| 283 | REAL, intent(in) :: s_parent, s_child |
---|
| 284 | REAL, intent(in) :: ds_parent, ds_child |
---|
| 285 | INTEGER, intent(in) :: dir |
---|
| 286 | ! |
---|
| 287 | REAL :: invcoeffraf |
---|
| 288 | INTEGER :: i, j, coeffraf |
---|
| 289 | INTEGER, DIMENSION(np) :: nbnonnuls |
---|
| 290 | REAL, DIMENSION(0:7), parameter :: invcoeff = (/1.,1.,0.5,1./3.,0.25,0.2,1./6.,1./7./) |
---|
| 291 | ! |
---|
| 292 | coeffraf = nint(ds_parent/ds_child) |
---|
| 293 | invcoeffraf = 1./coeffraf |
---|
| 294 | ! |
---|
| 295 | IF (Agrif_UseSpecialValueInUpdate) THEN |
---|
| 296 | ! |
---|
| 297 | ! nbnonnuls = 0 |
---|
| 298 | do j = 1,coeffraf |
---|
| 299 | do i = 1,np |
---|
| 300 | IF (y(indchildaverage(i,dir) + j -1) /= Agrif_SpecialValueFineGrid) THEN |
---|
| 301 | ! nbnonnuls(i) = nbnonnuls(i) + 1 |
---|
| 302 | x(i) = x(i) + y(indchildaverage(i,dir) + j-1 ) |
---|
| 303 | ENDIF |
---|
| 304 | enddo |
---|
| 305 | enddo |
---|
| 306 | do i=1,np |
---|
| 307 | ! x(i) = x(i)*invcoeff(nbnonnuls(i)) |
---|
| 308 | ! if (nbnonnuls(i) == 0) x(i) = Agrif_SpecialValueFineGrid |
---|
| 309 | enddo |
---|
| 310 | ! |
---|
| 311 | ELSE |
---|
| 312 | ! |
---|
| 313 | !CDIR NOLOOPCHG |
---|
| 314 | do j = 1,coeffraf |
---|
| 315 | !CDIR VECTOR |
---|
| 316 | do i= 1,np |
---|
| 317 | x(i) = x(i) + y(indchildaverage(i,dir) + j-1 ) |
---|
| 318 | enddo |
---|
| 319 | enddo |
---|
| 320 | IF (.not.Agrif_Update_Weights) THEN |
---|
| 321 | x = x * invcoeffraf |
---|
| 322 | ENDIF |
---|
| 323 | ! |
---|
| 324 | ENDIF |
---|
| 325 | |
---|
| 326 | !--------------------------------------------------------------------------------------------------- |
---|
| 327 | end subroutine Average1dAfterCompute |
---|
| 328 | !=================================================================================================== |
---|
| 329 | ! |
---|
| 330 | !=================================================================================================== |
---|
| 331 | ! subroutine Agrif_basicupdate_full_weighting1D |
---|
| 332 | ! |
---|
| 333 | !> Carries out an update by full_weighting on a parent grid (vector x) from its child grid (vector y). |
---|
| 334 | !--------------------------------------------------------------------------------------------------- |
---|
| 335 | subroutine Agrif_basicupdate_full_weighting1D ( x, y, np, nc, s_parent, s_child, ds_parent, ds_child ) |
---|
| 336 | !--------------------------------------------------------------------------------------------------- |
---|
| 337 | real, dimension(np), intent(out) :: x |
---|
| 338 | real, dimension(nc), intent(in) :: y |
---|
| 339 | integer, intent(in) :: np, nc |
---|
| 340 | real, intent(in) :: s_parent, s_child |
---|
| 341 | real, intent(in) :: ds_parent, ds_child |
---|
| 342 | !--------------------------------------------------------------------------------------------------- |
---|
| 343 | REAL :: xpos, xposfin |
---|
| 344 | INTEGER :: i, ii, diffmod |
---|
| 345 | INTEGER :: it1, it2 |
---|
| 346 | INTEGER :: i1, i2 |
---|
| 347 | INTEGER :: coeffraf |
---|
| 348 | INTEGER :: locind_child_left |
---|
| 349 | REAL :: sumweight, invsumweight |
---|
| 350 | REAL :: weights(-Agrif_MaxRaff:Agrif_MaxRaff) |
---|
| 351 | |
---|
| 352 | coeffraf = nint(ds_parent/ds_child) |
---|
| 353 | locind_child_left = 1 + agrif_int((s_parent-s_child)/ds_child) |
---|
| 354 | ! |
---|
| 355 | if (coeffraf == 1) then |
---|
| 356 | x(1:np) = y(locind_child_left:locind_child_left+np-1) |
---|
| 357 | return |
---|
| 358 | endif |
---|
| 359 | ! |
---|
| 360 | xpos = s_parent |
---|
| 361 | x = 0. |
---|
| 362 | ! |
---|
| 363 | xposfin = s_child + ds_child * (locind_child_left - 1) |
---|
| 364 | IF (abs(xposfin - xpos) < 0.001) THEN |
---|
| 365 | diffmod = 0 |
---|
| 366 | ELSE |
---|
| 367 | diffmod = 1 |
---|
| 368 | ENDIF |
---|
| 369 | ! |
---|
| 370 | if (diffmod == 1) THEN |
---|
| 371 | invsumweight=1./(2.*coeffraf**2) |
---|
| 372 | do i = -coeffraf,-1 |
---|
| 373 | weights(i) = invsumweight*(2*(coeffraf+i)+1) |
---|
| 374 | enddo |
---|
| 375 | do i = 0,coeffraf-1 |
---|
| 376 | weights(i) = weights(-(i+1)) |
---|
| 377 | enddo |
---|
| 378 | it1 = -coeffraf |
---|
| 379 | i1 = -(coeffraf-1)+locind_child_left |
---|
| 380 | i2 = 2*coeffraf - 1 |
---|
| 381 | |
---|
| 382 | else |
---|
| 383 | invsumweight=1./coeffraf**2 |
---|
| 384 | do i = -(coeffraf-1),0 |
---|
| 385 | weights(i) = invsumweight*(coeffraf + i) |
---|
| 386 | enddo |
---|
| 387 | do i=1,coeffraf-1 |
---|
| 388 | weights(i) = invsumweight*(coeffraf - i) |
---|
| 389 | enddo |
---|
| 390 | it1 = -(coeffraf-1) |
---|
| 391 | i1 = -(coeffraf-1)+locind_child_left |
---|
| 392 | i2 = 2*coeffraf - 2 |
---|
| 393 | |
---|
| 394 | endif |
---|
| 395 | ! |
---|
| 396 | sumweight = 0. |
---|
| 397 | MYLOOP : do i = 1,np |
---|
| 398 | ! |
---|
| 399 | it2 = it1 |
---|
| 400 | |
---|
| 401 | ! sumweight = 0. |
---|
| 402 | |
---|
| 403 | do ii = i1,i1+i2 |
---|
| 404 | ! |
---|
| 405 | IF (Agrif_UseSpecialValueInUpdate) THEN |
---|
| 406 | IF (y(ii) /= Agrif_SpecialValueFineGrid) THEN |
---|
| 407 | x(i) = x(i) + weights(it2)*y(ii) |
---|
| 408 | ! sumweight = sumweight+weights(it2) |
---|
| 409 | ELSE |
---|
| 410 | x(i) = Agrif_SpecialValueFineGrid |
---|
| 411 | i1=i1+coeffraf |
---|
| 412 | CYCLE MYLOOP |
---|
| 413 | ENDIF |
---|
| 414 | ELSE |
---|
| 415 | x(i) = x(i) + weights(it2)*y(ii) |
---|
| 416 | ENDIF |
---|
| 417 | |
---|
| 418 | it2 = it2+1 |
---|
| 419 | ! |
---|
| 420 | enddo |
---|
| 421 | ! |
---|
| 422 | i1 = i1 + coeffraf |
---|
| 423 | ! |
---|
| 424 | enddo MYLOOP |
---|
| 425 | |
---|
| 426 | IF (Agrif_UseSpecialValueInUpdate) THEN |
---|
| 427 | x = x * coeffraf ! x will be divided by coeffraf later in modupdate.F90 |
---|
| 428 | ENDIF |
---|
| 429 | |
---|
| 430 | !--------------------------------------------------------------------------------------------------- |
---|
| 431 | end subroutine Agrif_basicupdate_full_weighting1D |
---|
| 432 | !=================================================================================================== |
---|
| 433 | ! |
---|
| 434 | end module Agrif_UpdateBasic |
---|