[2352] | 1 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
| 2 | ! |
---|
| 3 | ! this routine performs a remapping based on addresses and weights |
---|
| 4 | ! computed in a setup phase |
---|
| 5 | ! |
---|
| 6 | !----------------------------------------------------------------------- |
---|
| 7 | ! |
---|
[5420] | 8 | ! CVS:$Id$ |
---|
[2352] | 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 | |
---|
| 35 | module remap_mod |
---|
| 36 | |
---|
| 37 | !----------------------------------------------------------------------- |
---|
| 38 | ! |
---|
| 39 | ! this module contains the routines for performing the actual |
---|
| 40 | ! remappings |
---|
| 41 | ! |
---|
| 42 | !----------------------------------------------------------------------- |
---|
| 43 | |
---|
| 44 | use kinds_mod ! defines common data types |
---|
| 45 | use constants ! defines common constants |
---|
| 46 | |
---|
| 47 | implicit none |
---|
| 48 | |
---|
| 49 | !*********************************************************************** |
---|
| 50 | |
---|
| 51 | contains |
---|
| 52 | |
---|
| 53 | !*********************************************************************** |
---|
| 54 | |
---|
| 55 | subroutine remap(dst_array, map_wts, dst_add, src_add, & |
---|
| 56 | src_array, src_grad1, src_grad2, src_grad3) |
---|
| 57 | |
---|
| 58 | !----------------------------------------------------------------------- |
---|
| 59 | ! |
---|
| 60 | ! performs the remapping based on weights computed elsewhere |
---|
| 61 | ! |
---|
| 62 | !----------------------------------------------------------------------- |
---|
| 63 | |
---|
| 64 | !----------------------------------------------------------------------- |
---|
| 65 | ! |
---|
| 66 | ! input arrays |
---|
| 67 | ! |
---|
| 68 | !----------------------------------------------------------------------- |
---|
| 69 | |
---|
| 70 | integer (kind=int_kind), dimension(:), intent(in) :: & |
---|
| 71 | dst_add, & ! destination address for each link |
---|
| 72 | src_add ! source address for each link |
---|
| 73 | |
---|
| 74 | real (kind=dbl_kind), dimension(:,:), intent(in) :: & |
---|
| 75 | map_wts ! remapping weights for each link |
---|
| 76 | |
---|
| 77 | real (kind=dbl_kind), dimension(:), intent(in) :: & |
---|
| 78 | src_array ! array with source field to be remapped |
---|
| 79 | |
---|
| 80 | real (kind=dbl_kind), dimension(:), intent(in), optional :: & |
---|
| 81 | src_grad1 & ! gradient arrays on source grid necessary for |
---|
| 82 | , src_grad2 & ! higher-order remappings |
---|
| 83 | , src_grad3 |
---|
| 84 | |
---|
| 85 | !----------------------------------------------------------------------- |
---|
| 86 | ! |
---|
| 87 | ! output variables |
---|
| 88 | ! |
---|
| 89 | !----------------------------------------------------------------------- |
---|
| 90 | |
---|
| 91 | real (kind=dbl_kind), dimension(:), intent(inout) :: & |
---|
| 92 | dst_array ! array for remapped field on destination grid |
---|
| 93 | |
---|
| 94 | !----------------------------------------------------------------------- |
---|
| 95 | ! |
---|
| 96 | ! local variables |
---|
| 97 | ! |
---|
| 98 | !----------------------------------------------------------------------- |
---|
| 99 | |
---|
| 100 | integer (kind=int_kind) :: n, iorder |
---|
| 101 | |
---|
| 102 | !----------------------------------------------------------------------- |
---|
| 103 | ! |
---|
| 104 | ! check the order of the interpolation |
---|
| 105 | ! |
---|
| 106 | !----------------------------------------------------------------------- |
---|
| 107 | |
---|
| 108 | if (present(src_grad1)) then |
---|
| 109 | iorder = 2 |
---|
| 110 | else |
---|
| 111 | iorder = 1 |
---|
| 112 | endif |
---|
| 113 | |
---|
| 114 | !----------------------------------------------------------------------- |
---|
| 115 | ! |
---|
| 116 | ! first order remapping |
---|
| 117 | ! |
---|
| 118 | !----------------------------------------------------------------------- |
---|
| 119 | |
---|
| 120 | dst_array = zero |
---|
| 121 | |
---|
| 122 | select case (iorder) |
---|
| 123 | case(1) |
---|
| 124 | |
---|
| 125 | do n=1,size(dst_add) |
---|
| 126 | dst_array(dst_add(n)) = dst_array(dst_add(n)) + & |
---|
| 127 | src_array(src_add(n))*map_wts(1,n) |
---|
| 128 | end do |
---|
| 129 | |
---|
| 130 | !----------------------------------------------------------------------- |
---|
| 131 | ! |
---|
| 132 | ! second order remapping |
---|
| 133 | ! |
---|
| 134 | !----------------------------------------------------------------------- |
---|
| 135 | |
---|
| 136 | case(2) |
---|
| 137 | |
---|
| 138 | if (size(map_wts,DIM=1) == 3) then |
---|
| 139 | do n=1,size(dst_add) |
---|
| 140 | dst_array(dst_add(n)) = dst_array(dst_add(n)) + & |
---|
| 141 | src_array(src_add(n))*map_wts(1,n) + & |
---|
| 142 | src_grad1(src_add(n))*map_wts(2,n) + & |
---|
| 143 | src_grad2(src_add(n))*map_wts(3,n) |
---|
| 144 | end do |
---|
| 145 | else if (size(map_wts,DIM=1) == 4) then |
---|
| 146 | do n=1,size(dst_add) |
---|
| 147 | dst_array(dst_add(n)) = dst_array(dst_add(n)) + & |
---|
| 148 | src_array(src_add(n))*map_wts(1,n) + & |
---|
| 149 | src_grad1(src_add(n))*map_wts(2,n) + & |
---|
| 150 | src_grad2(src_add(n))*map_wts(3,n) + & |
---|
| 151 | src_grad3(src_add(n))*map_wts(4,n) |
---|
| 152 | end do |
---|
| 153 | endif |
---|
| 154 | |
---|
| 155 | end select |
---|
| 156 | |
---|
| 157 | !----------------------------------------------------------------------- |
---|
| 158 | |
---|
| 159 | end subroutine remap |
---|
| 160 | |
---|
| 161 | !*********************************************************************** |
---|
| 162 | |
---|
| 163 | end module remap_mod |
---|
| 164 | |
---|
| 165 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|