remap.f 5.7 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163
  1. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  2. !
  3. ! this routine performs a remapping based on addresses and weights
  4. ! computed in a setup phase
  5. !
  6. !-----------------------------------------------------------------------
  7. !
  8. ! CVS:$Id: remap.f,v 1.5 2000/04/19 21:56:25 pwjones Exp $
  9. !
  10. ! Copyright (c) 1997, 1998 the Regents of the University of
  11. ! California.
  12. !
  13. ! This software and ancillary information (herein called software)
  14. ! called SCRIP is made available under the terms described here.
  15. ! The software has been approved for release with associated
  16. ! LA-CC Number 98-45.
  17. !
  18. ! Unless otherwise indicated, this software has been authored
  19. ! by an employee or employees of the University of California,
  20. ! operator of the Los Alamos National Laboratory under Contract
  21. ! No. W-7405-ENG-36 with the U.S. Department of Energy. The U.S.
  22. ! Government has rights to use, reproduce, and distribute this
  23. ! software. The public may copy and use this software without
  24. ! charge, provided that this Notice and any statement of authorship
  25. ! are reproduced on all copies. Neither the Government nor the
  26. ! University makes any warranty, express or implied, or assumes
  27. ! any liability or responsibility for the use of this software.
  28. !
  29. ! If software is modified to produce derivative works, such modified
  30. ! software should be clearly marked, so as not to confuse it with
  31. ! the version available from Los Alamos National Laboratory.
  32. !
  33. !***********************************************************************
  34. module remap_mod
  35. !-----------------------------------------------------------------------
  36. !
  37. ! this module contains the routines for performing the actual
  38. ! remappings
  39. !
  40. !-----------------------------------------------------------------------
  41. use kinds_mod ! defines common data types
  42. use constants ! defines common constants
  43. implicit none
  44. !***********************************************************************
  45. contains
  46. !***********************************************************************
  47. subroutine remap(dst_array, map_wts, dst_add, src_add,
  48. & src_array, src_grad1, src_grad2, src_grad3)
  49. !-----------------------------------------------------------------------
  50. !
  51. ! performs the remapping based on weights computed elsewhere
  52. !
  53. !-----------------------------------------------------------------------
  54. !-----------------------------------------------------------------------
  55. !
  56. ! input arrays
  57. !
  58. !-----------------------------------------------------------------------
  59. integer (kind=int_kind), dimension(:), intent(in) ::
  60. & dst_add, ! destination address for each link
  61. & src_add ! source address for each link
  62. real (kind=dbl_kind), dimension(:,:), intent(in) ::
  63. & map_wts ! remapping weights for each link
  64. real (kind=dbl_kind), dimension(:), intent(in) ::
  65. & src_array ! array with source field to be remapped
  66. real (kind=dbl_kind), dimension(:), intent(in), optional ::
  67. & src_grad1 ! gradient arrays on source grid necessary for
  68. &, src_grad2 ! higher-order remappings
  69. &, src_grad3
  70. !-----------------------------------------------------------------------
  71. !
  72. ! output variables
  73. !
  74. !-----------------------------------------------------------------------
  75. real (kind=dbl_kind), dimension(:), intent(inout) ::
  76. & dst_array ! array for remapped field on destination grid
  77. !-----------------------------------------------------------------------
  78. !
  79. ! local variables
  80. !
  81. !-----------------------------------------------------------------------
  82. integer (kind=int_kind) :: n, iorder
  83. !-----------------------------------------------------------------------
  84. !
  85. ! check the order of the interpolation
  86. !
  87. !-----------------------------------------------------------------------
  88. if (present(src_grad1)) then
  89. iorder = 2
  90. else
  91. iorder = 1
  92. endif
  93. !-----------------------------------------------------------------------
  94. !
  95. ! first order remapping
  96. !
  97. !-----------------------------------------------------------------------
  98. dst_array = zero
  99. select case (iorder)
  100. case(1)
  101. do n=1,size(dst_add)
  102. dst_array(dst_add(n)) = dst_array(dst_add(n)) +
  103. & src_array(src_add(n))*map_wts(1,n)
  104. end do
  105. !-----------------------------------------------------------------------
  106. !
  107. ! second order remapping
  108. !
  109. !-----------------------------------------------------------------------
  110. case(2)
  111. if (size(map_wts,DIM=1) == 3) then
  112. do n=1,size(dst_add)
  113. dst_array(dst_add(n)) = dst_array(dst_add(n)) +
  114. & src_array(src_add(n))*map_wts(1,n) +
  115. & src_grad1(src_add(n))*map_wts(2,n) +
  116. & src_grad2(src_add(n))*map_wts(3,n)
  117. end do
  118. else if (size(map_wts,DIM=1) == 4) then
  119. do n=1,size(dst_add)
  120. dst_array(dst_add(n)) = dst_array(dst_add(n)) +
  121. & src_array(src_add(n))*map_wts(1,n) +
  122. & src_grad1(src_add(n))*map_wts(2,n) +
  123. & src_grad2(src_add(n))*map_wts(3,n) +
  124. & src_grad3(src_add(n))*map_wts(4,n)
  125. end do
  126. endif
  127. end select
  128. !-----------------------------------------------------------------------
  129. end subroutine remap
  130. !***********************************************************************
  131. end module remap_mod
  132. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!